ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
mos1d_dg_n.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
16
17
22int main()
23{
24 typedef viennagrid::line_1d_mesh MeshType;
25 typedef viennashe::device<MeshType> DeviceType;
26
27 typedef viennagrid::result_of::const_cell_range<MeshType>::type CellContainer;
28 typedef viennagrid::result_of::iterator<CellContainer>::type CellIterator;
29 typedef viennagrid::result_of::const_facet_range<MeshType>::type FacetContainer;
30
31 std::cout << viennashe::preamble() << std::endl;
32
33 std::cout << "* main(): Creating mesh ..." << std::endl;
34
35 DeviceType device;
36
37 const double len_gate = 1e-9;
38 const double cs_gate = 0.01e-9;
39 const double len_oxide = 1e-9;
40 const double cs_ox = 0.05e-9;
41 const double len_bulk = 100e-9;
42 const double cs_bulk = 1e-9;
43
44 double Vg = 0.2;
45
46 bool ok = false;
47
48 mos1d_mesh_generator mosgen(len_gate, cs_gate, len_oxide, cs_ox, len_bulk, cs_bulk);
49 device.generate_mesh(mosgen);
50
51 // ND NA
52 init_device(device, Vg, 1e8, 3e23 );
53
54 //
55 // DEBUG OUTPUT
56 std::cout << std::endl;
57 CellContainer cells(device.mesh());
58 for (CellIterator cit = cells.begin();
59 cit != cells.end();
60 ++cit)
61 {
62 std::cout << cit->id() << " => " << device.get_material(*cit)
63 << " ## Nd = " << device.get_doping_n(*cit)
64 << " ## Na = " << device.get_doping_p(*cit)
65 << " ## has bnd cond: " << (device.has_contact_potential(*cit)== true ? "TRUE" : "FALSE")
66 << std::endl;
67 } // for cells
68 FacetContainer facets(device.mesh());
69 for (std::size_t i = 0; i < facets.size(); ++i)
70 std::cout << facets[i] << std::endl;
71
72 //
73 // Use a drift-diffusion simulation
74 std::cout << "* main(): Creating DD simulator..." << std::endl;
75
76 viennashe::config dd_cfg;
77 dd_cfg.with_holes(true);
78 dd_cfg.with_electrons(true);
81
82 // We need the dense solver here ...
84 dd_cfg.nonlinear_solver().max_iters(30);
85 dd_cfg.nonlinear_solver().damping(0.6);
86
87 dd_cfg.quantum_correction(true); // solve density gradient equations
88 dd_cfg.with_quantum_correction(true); // apply correction potentials from density gradient
89
90 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
91
92 std::cout << "* main(): Launching simulator..." << std::endl;
93
94 dd_simulator.run();
95
97 gpwriter(device, viennashe::util::any_filter(), dd_simulator.electron_density(), "mos1d_dg_n_dd_p.dat");
98 gpwriter(device, viennashe::util::any_filter(), dd_simulator.hole_density(), "mos1d_dg_n_dd_p.dat");
99
100 std::cout << "* main(): Testing results (DD) ..." << std::endl;
101
102 ok = test_result(dd_simulator.dg_pot_p(), device, "../../tests/data/dg/n/dd_dgpot_holes.dat" );
103 if(!ok) return EXIT_FAILURE;
104 ok = test_result(dd_simulator.dg_pot_n(), device, "../../tests/data/dg/n/dd_dgpot_electrons.dat" );
105 if(!ok) return EXIT_FAILURE;
106 ok = test_result(dd_simulator.potential(), device, "../../tests/data/dg/n/mos1d_dd_pot.dat" );
107 if(!ok) return EXIT_FAILURE;
108 ok = test_result(dd_simulator.electron_density(), device, "../../tests/data/dg/n/mos1d_dd_electrons.dat");
109 if(!ok) return EXIT_FAILURE;
110 ok = test_result(dd_simulator.hole_density(), device, "../../tests/data/dg/n/mos1d_dd_holes.dat");
111 if(!ok) return EXIT_FAILURE;
112
113
114#if 0 // DISABLED DUE TO CONVERGENCY PROBLEMS
115
116 //
117 // Start self-consistent SHE simulations
118 //
119 std::cout << "* main(): Setting up SHE..." << std::endl;
120
121 viennashe::config config;
122 config.with_electrons(true);
124 config.with_holes(true);
126 config.with_traps(false);
127 config.scattering().acoustic_phonon().enabled(true);
128 config.scattering().optical_phonon().enabled(true);
129 config.scattering().ionized_impurity().enabled(true);
130 config.scattering().impact_ionization().enabled(false);
131 config.scattering().electron_electron(false);
132
133 //config.linear_solver().set(viennashe::solvers::linear_solver_ids::dense_linear_solver);
134 config.nonlinear_solver().max_iters(30);
135 config.nonlinear_solver().damping(0.5); //use smaller damping with Newton
136 config.max_expansion_order(1);
137 config.linear_solver().max_iters(2000);
138 config.linear_solver().ilut_drop_tolerance(1e-2);
140 config.energy_spacing(viennashe::physics::constants::q * 0.031); //energy spacing of 31meV
141
142 config.quantum_correction(true); // calculate a correction potential
143 config.with_quantum_correction(true); // use the correction potential for the transport model
144
145 //
146 // Launch solver
147 //
148 std::cout << "* main(): Computing SHE..." << std::endl;
149 viennashe::simulator<DeviceType> she_simulator(device, config);
150 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
151 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
152 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
153 //she_simulator.set_initial_guess(viennashe::quantity::density_gradient_electron_correction(), dd_simulator.dg_pot_n());
154 //she_simulator.set_initial_guess(viennashe::quantity::density_gradient_hole_correction(), dd_simulator.dg_pot_p());
155 she_simulator.run();
156
157 //
158 // TEST
159 //
160 std::cout << "* main(): Testing results (SHE) ..." << std::endl;
161
162 ok = test_result(dd_simulator.dg_pot_p(), device, "../../tests/data/dg/n/she_dgpot_holes.dat" );
163 if(!ok) return EXIT_FAILURE;
164 ok = test_result(dd_simulator.dg_pot_n(), device, "../../tests/data/dg/n/she_dgpot_electrons.dat" );
165 if(!ok) return EXIT_FAILURE;
166 ok = test_result(dd_simulator.potential(), device, "../../tests/data/dg/n/mos1d_she_pot.dat" );
167 if(!ok) return EXIT_FAILURE;
168 ok = test_result(dd_simulator.electron_density(), device, "../../tests/data/dg/n/mos1d_she_electrons.dat");
169 if(!ok) return EXIT_FAILURE;
170 ok = test_result(dd_simulator.hole_density(), device, "../../tests/data/dg/n/mos1d_she_holes.dat");
171 if(!ok) return EXIT_FAILURE;
172#endif
173
174 std::cout << "... \\o/ SUCCESS \\o/ ..." << std::endl;
175
176 return EXIT_SUCCESS;
177}
The main SHE configuration class. To be adjusted by the user for his/her needs.
Definition: config.hpp:124
double min_kinetic_energy_range(viennashe::carrier_type_id ctype) const
Returns the minimum kinetic energy range for the selected carrier.
Definition: config.hpp:388
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 quantum_correction() const
Definition: config.hpp:509
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
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
viennashe::she::scatter_config & scattering()
Returns the configuration object for scattering.
Definition: config.hpp:475
bool with_quantum_correction() const
Definition: config.hpp:512
Defines the physical properties of a device, e.g. doping. This is the implementation for 2d and highe...
Definition: device.hpp:818
ionized_impurity_scattering_parameters const & ionized_impurity() const
Returns the parameters for ionized impurity scattering. Const-version.
Definition: config.hpp:377
impact_ionization_scattering_parameters const & impact_ionization() const
Returns the parameters for impact ionization scattering. Const-version.
Definition: config.hpp:382
acoustic_phonon_scattering_parameters const & acoustic_phonon() const
Returns the parameters for acoustic phonon scattering. Const-version.
Definition: config.hpp:367
optical_phonon_scattering_parameters const & optical_phonon() const
Returns the parameters for optical phonon scattering. Const-version.
Definition: config.hpp:372
bool electron_electron() const
Returns true if electron-electron scattering is activated.
Definition: config.hpp:387
Class for self-consistent SHE simulations.
Definition: simulator.hpp:675
void set(long solver_id)
Sets a new linear solver. Provide IDs defined in.
Definition: config.hpp:84
double ilut_drop_tolerance() const
Returns the drop tolerance for the ILUT preconditioenr.
Definition: config.hpp:177
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
Contains common code for the 1D MOS DG tests.
bool test_result(AccessorType const &quan, DeviceType const &device, std::string filename)
Tests the given quantity against reference data in a file. Uses reference_values.
Definition: mos1d_dg.hpp:182
void init_device(DeviceType &device, double Vg_init, double Nd, double Na)
Initalizes the MOS 1D device for testing.
Definition: mos1d_dg.hpp:94
std::string electron_density()
@ EQUATION_SHE
Definition: forwards.h:118
@ EQUATION_CONTINUITY
Definition: forwards.h:117
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
MOS 1D test-grid generator.
Definition: mos1d_dg.hpp:55
Writes quantities to a file which can be processed by Gnuplot. Works for 1d, 2d and 3d only.
static const double q
Elementary charge.
Definition: constants.hpp:44
A trivial filter, all objects are accepted.
Definition: filter.hpp:41