ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
np-diode-bipolar.cpp
Go to the documentation of this file.
1/* ============================================================================
2 Copyright (c) 2011-2022, Institute for Microelectronics,
3 Institute for Analysis and Scientific Computing,
4 TU Wien.
5
6 -----------------
7 ViennaSHE - The Vienna Spherical Harmonics Expansion Boltzmann Solver
8 -----------------
9
10 http://viennashe.sourceforge.net/
11
12 License: MIT (X11), see file LICENSE in the base directory
13=============================================================================== */
14
15#if defined(_MSC_VER)
16 // Disable name truncation warning obtained in Visual Studio
17 #pragma warning(disable:4503)
18#endif
19
20#include <iostream>
21#include <cstdlib>
22#include <vector>
23
24// ViennaSHE includes:
25#include "viennashe/core.hpp"
26
27// ViennaGrid mesh configurations and centroid() algorithm:
28#include "viennagrid/config/default_configs.hpp"
29#include "viennagrid/algorithm/centroid.hpp"
30
59template <typename DeviceType>
60void init_device(DeviceType & device)
61{
62 typedef typename DeviceType::segment_type SegmentType;
63
65 SegmentType const & contact_left = device.segment(0);
66 SegmentType const & n_left = device.segment(1);
67 SegmentType const & i_center = device.segment(2);
68 SegmentType const & p_right = device.segment(3);
69 SegmentType const & contact_right = device.segment(4);
70
72 device.set_material(viennashe::materials::metal(), contact_left);
73 device.set_material(viennashe::materials::si(), n_left);
74 device.set_material(viennashe::materials::si(), i_center);
75 device.set_material(viennashe::materials::si(), p_right);
76 device.set_material(viennashe::materials::metal(), contact_right);
77
83 device.set_doping_n(1e22, n_left);
84 device.set_doping_p(1e10, n_left);
85
86 device.set_doping_n(1e16, i_center);
87 device.set_doping_p(1e16, i_center);
88
89 device.set_doping_n(1e10, p_right);
90 device.set_doping_p(1e22, p_right);
91
92
100 trap.collision_cross_section(3.2e-20);
101 trap.density(1e21);
102 trap.energy(0.3 * viennashe::physics::constants::q); // 0.3 eV above center of band gap
103 device.add_trap_level(trap);
104
105 // second trap level with different energy, but same density and collision cross section:
106 trap.energy(-0.3 * viennashe::physics::constants::q); // -0.3 eV below center of band gap
107 device.add_trap_level(trap);
108
109
111 device.set_contact_potential( 0.0, contact_left);
112 device.set_contact_potential(-0.2, contact_right);
113}
114
115
124int main()
125{
130 typedef viennagrid::quadrilateral_2d_mesh MeshType;
131 typedef viennashe::device<MeshType> DeviceType;
132
133 std::cout << viennashe::preamble() << std::endl;
134
136 std::cout << "* main(): Creating device..." << std::endl;
137 DeviceType device;
138
146
158 generator_params.add_segment(0.0, 1.0e-7, 2, //start at x=0, length 100nm, 2 points
159 0.0, 1.0e-8, 3); //start at y=0, length 10nm, 3 points
160 generator_params.add_segment(0.0, 1.0e-7, 2, //start at x=0, length 100nm, 2 points
161 1.0e-8, 1.0e-7, 21); //start at y= 20nm, length 100nm, 21 points
162 generator_params.add_segment(0.0, 1.0e-7, 2, //start at x=0, length 100nm, 2 points
163 1.1e-7, 5.0e-9, 2); //start at y=220nm, length 5nm, 2 points
164 generator_params.add_segment(0.0, 1.0e-7, 2, //start at x=0, length 100nm, 2 points
165 1.15e-7,1.0e-7, 21); //start at y=230nm, length 100nm, 21 points
166 generator_params.add_segment(0.0, 1.0e-7, 2, //start at x=0, length 100nm, 2 points
167 2.15e-7,1.0e-8, 3); //start at y=430, length 10nm, 3 points
168
170 device.generate_mesh(generator_params);
171
175 std::cout << "* main(): Initializing device..." << std::endl;
176 init_device(device);
177
184 std::cout << "* main(): Creating DD simulator..." << std::endl;
185
193 viennashe::config dd_cfg; // Per default the config is set to DD with Gummel iterations
194
195 // Nonlinear solver: Use up to 100 Gummel iterations with a rather strong damping of 0.2
196 dd_cfg.nonlinear_solver().max_iters(100);
197 dd_cfg.nonlinear_solver().damping(0.2);
198
205 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
206
207 std::cout << "* main(): Launching simulator..." << std::endl;
208 dd_simulator.run();
209
215 viennashe::io::write_quantities_to_VTK_file(dd_simulator, "np-diode_dd_quan");
216
217
224 std::cout << "* main(): Setting up SHE..." << std::endl;
225
232 viennashe::config config;
233 // Configure a first-order bipolar SHE
234 config.with_electrons(true);
236
237 config.with_holes(true);
239
240 config.with_traps(true);
241 config.with_trap_selfconsistency(true);
242
243 // Linear solver: Set up to 1000 iterations
244 config.linear_solver().max_iters(1000);
245
246 // Nonlinear solver: Set up to 20 Gummel iterations with a moderate damping of 0.3
247 config.nonlinear_solver().max_iters(20);
248 config.nonlinear_solver().damping(0.3);
249
250 // SHE: Set first-order expansion and use an energy spacing of 12.5 meV
251 config.max_expansion_order(1);
253
254
262 std::cout << "* main(): Computing first-order SHE..." << std::endl;
263
264 viennashe::simulator<DeviceType> she_simulator(device, config);
265
266 // Use the DD solution as initial guess
267 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
268 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
269 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
270
271 she_simulator.run();
272
278 std::cout << "* main(): Writing SHE result..." << std::endl;
280 she_simulator.config(),
281 she_simulator.quantities().electron_distribution_function(),
282 "np-diode-she");
283
286 viennashe::io::write_quantities_to_VTK_file(she_simulator, "np-diode_she_quan");
287
289 std::cout << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)." << std::endl;
290 std::cout << "* main(): Don't forget to scale the z-axis by about a factor of 1e12 when examining the distribution function." << std::endl;
291 std::cout << std::endl;
292 std::cout << "*********************************************************" << std::endl;
293 std::cout << "* ViennaSHE finished successfully *" << std::endl;
294 std::cout << "*********************************************************" << std::endl;
295
296 return EXIT_SUCCESS;
297}
298
The main SHE configuration class. To be adjusted by the user for his/her needs.
Definition: config.hpp:124
long max_expansion_order() const
Returns the current maximum expansion order.
Definition: config.hpp:369
bool with_holes() const
Returns true if holes are considered in the simulation.
Definition: config.hpp:234
bool with_traps() const
Returns true if traps are considered in the simulation.
Definition: config.hpp:243
nonlinear_solver_config_type & nonlinear_solver()
Returns the configuration object for the nonlinear solver.
Definition: config.hpp:495
void set_electron_equation(equation_id equ_id)
Definition: config.hpp:231
linear_solver_config_type & linear_solver()
Returns the configuration object for the linear solver.
Definition: config.hpp:485
double energy_spacing() const
Returns the uniform energy spacing of discrete energies.
Definition: config.hpp:460
bool with_trap_selfconsistency() const
Returns true if traps are considered self-consistently in the simulation.
Definition: config.hpp:247
void set_hole_equation(equation_id equ_id)
Definition: config.hpp:239
bool with_electrons() const
Returns true if electrons are considered in the simulation.
Definition: config.hpp:226
Defines the physical properties of a device, e.g. doping. This is the implementation for 2d and highe...
Definition: device.hpp:818
Class for self-consistent SHE simulations.
Definition: simulator.hpp:675
std::size_t max_iters() const
Returns the maximum number of iterative linear solver iterations.
Definition: config.hpp:155
std::size_t max_iters() const
Returns the maximum number of nonlinear solver iterations.
Definition: config.hpp:317
double damping() const
Returns the damping for the nonlinear solver.
Definition: config.hpp:328
Describes a SRH trap.
Definition: trap_level.hpp:31
void energy(double e)
Returns the trap energy in Joule (zero energy refers to the center of the band gap.
Definition: trap_level.hpp:77
void collision_cross_section(double ccs)
Sets the collision cross section (SI units)
Definition: trap_level.hpp:42
void density(double d)
Returns the trap density for the trap level.
Definition: trap_level.hpp:57
Configuration class for the simple mesh generator.
void add_segment(double start_x, double start_y, double len_x, double len_y, unsigned long points_x, unsigned long points_y)
Convenience header, which includes all core functionality available in ViennaSHE.
void init_device(DeviceType &device, double Vg_init, double Nd, double Na)
Initalizes the MOS 1D device for testing.
Definition: mos1d_dg.hpp:94
void write_quantities_to_VTK_file(viennashe::simulator< DeviceType > const &simulator_obj, std::string filename, bool include_debug_information=false)
Generic interface function for writing simulated quantities to a VTK file.
std::string electron_density()
@ EQUATION_SHE
Definition: forwards.h:118
std::string preamble()
Prints the ViennaSHE preamble (header). Used in all the examples as well as the standalone-applicatio...
Definition: version.hpp:63
int main()
Definition: resistor1d-c.c:108
A class referring to any metal contact.
Definition: all.hpp:44
A class referring to silicon and providing certain material parameters Note that this is the default ...
Definition: all.hpp:50
static const double q
Elementary charge.
Definition: constants.hpp:44