ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
simple_impurity_scattering.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#include <stdexcept>
24
25#include "tests/src/common.hpp"
26
27// ViennaSHE includes:
28#include "viennashe/core.hpp"
29
30// ViennaGrid default mesh configuration:
31#include "viennagrid/config/default_configs.hpp"
32
33
45template <typename DeviceType>
46void init_device(DeviceType & device, double voltage)
47{
48 typedef typename DeviceType::mesh_type MeshType;
49
50 // STEP 2: Set doping
51 std::cout << "* init_device(): Setting doping..." << std::endl;
52
53 device.set_doping_n(1e28);
54 device.set_doping_p(1e4);
55 device.set_material(viennashe::materials::si());
56
57 typedef typename viennagrid::result_of::const_cell_range<MeshType>::type CellContainer;
58
59 // STEP 3: Define contacts
60 double gnd = 0.0;
61 double vcc = voltage;// * 2.0;
62
63 CellContainer cells(device.mesh());
64
65 if (cells.size() < 3) throw std::runtime_error("The Mesh is too small. It contains less than 3 cells!");
66
67 device.set_contact_potential(gnd, cells[0]);
68 device.set_material(viennashe::materials::metal(), cells[0]);
69 device.set_contact_potential(vcc, cells[cells.size()-1]);
70 device.set_material(viennashe::materials::metal(), cells[cells.size()-1]);
71
72}
73
74int main()
75{
76 typedef viennagrid::line_1d_mesh MeshType;
77 typedef viennashe::device<MeshType> DeviceType;
78 typedef viennagrid::result_of::const_cell_range<MeshType>::type CellContainer;
79
80 std::cout << viennashe::preamble() << std::endl;
81
82 std::cout << "* main(): Initializing device..." << std::endl;
83 DeviceType device;
84
85 //
86 // Generate device geometry using built-in device generator
87 //
88 {
90 generator_params.add_segment(0.0, 1e-6, 51); //start at x=0, length 1e-6, 51 points
91 device.generate_mesh(generator_params);
92 }
93
94 //
95 // Set up and initialize device
96 //
97
98 init_device(device, 0.1);
99
100 CellContainer cells(device.mesh());
101
102 //
103 // Use a drift-diffusion simulation for obtaining an initial guess of the potential:
104 //
105 std::cout << "* main(): Creating DD simulator..." << std::endl;
106 viennashe::config dd_cfg;
108 dd_cfg.nonlinear_solver().max_iters(50);
109 dd_cfg.nonlinear_solver().damping(0.5);
110 //dd_cfg.nonlinear_solver().set(viennashe::solvers::nonlinear_solver_ids::newton_nonlinear_solver); //uncomment this to use Newton (might require stronger damping)
111
112 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
113
114 std::cout << "* main(): Launching simulator..." << std::endl;
115 dd_simulator.run();
116
117 //
118 // A single SHE postprocessing step (fixed field)
119 //
120 std::cout << "* main(): Setting up SHE..." << std::endl;
121
122 viennashe::config config;
123
125 config.with_electrons(true);
126
127 config.with_holes(false);
129
130 config.nonlinear_solver().max_iters(10);
131 config.nonlinear_solver().damping(0.3);
132
133 config.max_expansion_order(1);
136
137 // configure scattering mechanisms. We use phonon scattering only for now.
138 config.scattering().acoustic_phonon().enabled(true);
139 config.scattering().optical_phonon().enabled(true);
140 config.scattering().ionized_impurity().enabled(false);
141 config.scattering().impact_ionization().enabled(false);
142 config.scattering().electron_electron(false);
143
144 //config.dispersion_relation(viennashe::she::dispersion_relation_ids::ext_vecchi_dispersion); //uncomment this to use the extended Vecchi model
145
146 std::cout << "* main(): Computing SHE (electrons only) without ionized impurity and Modena model ..." << std::endl;
147 viennashe::simulator<DeviceType> she_simulator(device, config);
148
149 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
150 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
151 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
152
153 she_simulator.run();
154
155 std::cout << "* main(): Writing SHE result..." << std::endl;
157
158 double n_mid_without = 0;
159 double Jn_mid_without = 0;
160
161 {
163 she_simulator.config(),
164 she_simulator.quantities().electron_distribution_function(),
165 "simple_impurity_scattering_edf1");
166
167 edfwriter(device, viennashe::util::any_filter(), she_simulator.edf(viennashe::ELECTRON_TYPE_ID), "simple_impurity_scattering_edf1.dat");
168
171 she_simulator.config(),
172 she_simulator.quantities().electron_distribution_function());
173
175 she_simulator.config(),
176 she_simulator.quantities().electron_distribution_function());
177
178 viennashe::io::write_cell_quantity_for_gnuplot(Jn, device, "simple_impurity_scattering_Jn1.dat");
179
180 n_mid_without = n(cells[25]);
181 Jn_mid_without = Jn(cells[25])[0];
182 }
183 // ####################################################################################################################################
184
185 config.scattering().ionized_impurity().enabled(true); // enable ionized impurity scattering
186
187 std::cout << "* main(): Computing SHE (electrons only) WITH ionized impurity and Modena model ..." << std::endl;
188 viennashe::simulator<DeviceType> she_simulator_with(device, config);
189
190 she_simulator_with.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
191 she_simulator_with.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
192 she_simulator_with.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
193
194 she_simulator_with.run();
195
196 std::cout << "* main(): Writing SHE result..." << std::endl;
197
199 she_simulator_with.config(),
200 she_simulator_with.quantities().electron_distribution_function(),
201 "simple_impurity_scattering_edf2");
202
203 //viennashe::io::gnuplot_edf_writer edfwriter;
204 edfwriter(device, viennashe::util::any_filter(), she_simulator_with.edf(viennashe::ELECTRON_TYPE_ID), "simple_impurity_scattering_edf2.dat");
205
208 config,
209 she_simulator_with.quantities().electron_distribution_function());
210
212 config,
213 she_simulator_with.quantities().electron_distribution_function());
214
215 viennashe::io::write_cell_quantity_for_gnuplot(Jn, device, "simple_impurity_scattering_Jn2.dat");
216
217 double n_mid_with = n(cells[25]);
218 double Jn_mid_with = Jn(cells[25])[0];
219
220 if(!viennashe::testing::fuzzy_equal(n_mid_without, n_mid_with))
221 {
222 std::cerr << "ERROR: Ionized impurity scattering should NOT change the electron density in a resistor!" << std::endl;
223 std::cerr << " The densities are not equal." << std::endl;
224 std::cerr << n_mid_without << " != " << n_mid_with << std::endl;
225 return EXIT_FAILURE;
226 }
227
228 if(viennashe::testing::fuzzy_equal(Jn_mid_without, Jn_mid_with, 1e-5))
229 {
230 std::cerr << "ERROR: Ionized impurity scattering should change the electron CURRENT density in a resistor!" << std::endl;
231 std::cerr << " It failed to do so." << std::endl;
232 std::cerr << Jn_mid_without << " == " << Jn_mid_with << std::endl;
233 return EXIT_FAILURE;
234 }
235 else
236 {
237 std::cerr << "BUT: This is ok!" << std::endl;
238 }
239
240 if(Jn_mid_without < Jn_mid_with)
241 {
242 std::cerr << "ERROR: Ionized impurity scattering should DECREASE the electron CURRENT density in a resistor!" << std::endl;
243 std::cerr << " It failed to do so." << std::endl;
244 return EXIT_FAILURE;
245 }
246
247
248 std::cout << std::endl;
249 std::cout << "*********************************************************" << std::endl;
250 std::cout << "* ViennaSHE finished successfully *" << std::endl;
251 std::cout << "*********************************************************" << std::endl;
252
253 return EXIT_SUCCESS;
254}
255
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
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
Defines the physical properties of a device, e.g. doping. This is the implementation for 2d and highe...
Definition: device.hpp:818
An accessor for the carrier density in the device.
Accessor class providing the carrier velocity inside the device.
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
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
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_cell_quantity_for_gnuplot(AccessorType const &quan, DeviceType const &device, std::string filename)
Writes a quantity (on vertices) per point to a text file suitable for gnuplot.
std::string electron_density()
bool fuzzy_equal(double is, double should, double tol=1e-1)
Performs a fuzzy (up to a tolerance) equal. Returns true if is and should are equal within the tolera...
Definition: common.hpp:39
@ ELECTRON_TYPE_ID
Definition: forwards.h:187
@ 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
Writes the energy distribution function to a file which can be processed by Gnuplot....
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
A trivial filter, all objects are accepted.
Definition: filter.hpp:41
Contains common functions, functors and other classes often needed by the tests.