ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
resistor.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:
28#include "viennagrid/config/default_configs.hpp"
29
30
55template <typename DeviceType>
56void init_device(DeviceType & device, double vcc)
57{
58 // ViennaGrid types for mesh handling:
59 typedef typename DeviceType::mesh_type MeshType;
60
61 typedef typename viennagrid::result_of::const_cell_range<MeshType>::type CellContainer;
62
63 // Container of all cells in the mesh:
64 CellContainer cells(device.mesh());
65
67 device.set_material(viennashe::materials::si()); // Silicon everywhere
68
69 // overwrite the material of the two cells on the left and right:
70 device.set_material(viennashe::materials::metal(), cells[0]);
71 device.set_material(viennashe::materials::metal(), cells[cells.size()-1]);
72
78 std::cout << "* init_device(): Setting doping..." << std::endl;
79 device.set_doping_n(1e28);
80 device.set_doping_p(1e4);
81
83 device.set_contact_potential(0.0, cells[0]);
84 device.set_contact_potential(vcc, cells[cells.size()-1]);
85}
86
87
96int main()
97{
102 typedef viennagrid::line_1d_mesh MeshType;
103 typedef viennashe::device<MeshType> DeviceType;
104
105 std::cout << viennashe::preamble() << std::endl;
106
108 DeviceType device;
109
110
117
124 generator_params.add_segment(0.0, 1e-6, 11); //start at x=0, length 1um, 11 points
125
127 device.generate_mesh(generator_params);
128
129
133 std::cout << "* main(): Initializing device..." << std::endl;
134 init_device(device, 1.0);
135
136
143 std::cout << "* main(): Creating DD simulator..." << std::endl;
144
152 viennashe::config dd_cfg;
153
154 // Linear solver: Use a dense linear solver, because the mesh is small (no preconditioning troubles):
156
157 // Nonlinear solver: 20 Gummel iterations with a damping parameter of 0.25:
158 dd_cfg.nonlinear_solver().max_iters(20);
159 dd_cfg.nonlinear_solver().damping(0.25);
160
167 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
168
169 std::cout << "* main(): Launching simulator..." << std::endl;
170 dd_simulator.run();
171
172
178 viennashe::io::write_quantities_to_VTK_file(dd_simulator, "resistor_dd_quan");
179
187 std::cout << "* main(): Setting up SHE..." << std::endl;
188
195 viennashe::config config;
196
197 config.with_electrons(true);
199
200 config.with_holes(true);
202
203 // Nonlinear solver: Only one Gummel iteration with very high damping (to avoid changes to the potential)
204 config.nonlinear_solver().max_iters(1);
205 config.nonlinear_solver().damping(0.001);
206
207 // SHE: First-order simulation using an energy spacing of 15.5meV and a kinetic energy range of at least 1eV:
208 config.max_expansion_order(1);
211
212 // configure scattering mechanisms. We use phonon scattering only for now.
213 config.scattering().acoustic_phonon().enabled(true);
214 config.scattering().optical_phonon().enabled(true);
215 config.scattering().ionized_impurity().enabled(false);
216 config.scattering().impact_ionization().enabled(false);
217 config.scattering().electron_electron(false);
218
219
227 std::cout << "* main(): Computing SHE (electrons only) ..." << std::endl;
228 viennashe::simulator<DeviceType> she_simulator(device, config);
229
230 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
231 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
232 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
233
234 she_simulator.run();
235
236 std::cout << "* main(): Writing SHE result..." << std::endl;
237
244 edfwriter(device, viennashe::util::any_filter(), she_simulator.edf(viennashe::ELECTRON_TYPE_ID), "resistor_edf.dat");
245
249 she_simulator.quantities().electron_distribution_function());
250
251 viennashe::io::write_cell_quantity_for_gnuplot(Jn, device, "resistor_Jn.dat");
252 viennashe::io::write_cell_quantity_for_gnuplot(she_simulator.electron_density(), device, "resistor_n.dat");
253
255 viennashe::io::write_quantities_to_VTK_file(she_simulator, "resistor_she_quan");
256
257
271 std::cout << "* main(): Computing SHE (electrons only) with EE scattering ..." << std::endl;
272 config.nonlinear_solver().max_iters(10);
273 config.scattering().electron_electron(true);
274
279 viennashe::simulator<DeviceType> she_simulator_ee(device, config);
280
281 she_simulator_ee.set_initial_guess(viennashe::quantity::potential(), she_simulator.potential());
282 she_simulator_ee.set_initial_guess(viennashe::quantity::electron_density(), she_simulator.electron_density());
283 she_simulator_ee.set_initial_guess(viennashe::quantity::hole_density(), she_simulator.hole_density());
284
285 she_simulator_ee.run();
286
292 std::cout << "* main(): Writing SHE result..." << std::endl;
293
294 edfwriter(device, viennashe::util::any_filter(), she_simulator_ee.edf(viennashe::ELECTRON_TYPE_ID), "resistor_edf_ee.dat");
295
298 viennashe::io::write_cell_quantity_for_gnuplot(she_simulator_ee.electron_density(), device, "resistor_n_ee.dat");
299
300 viennashe::she::carrier_energy_wrapper<she_quan_type> avge (config, she_simulator_ee.quantities().electron_distribution_function());
301 viennashe::io::write_cell_quantity_for_gnuplot(avge, device, "resistor_avg_ee.dat");
302
304 std::cout << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)." << std::endl;
305 std::cout << "* main(): Don't forget to scale the z-axis by about a factor of 1e12 when examining the distribution function." << std::endl;
306 std::cout << std::endl;
307 std::cout << "*********************************************************" << std::endl;
308 std::cout << "* ViennaSHE finished successfully *" << std::endl;
309 std::cout << "*********************************************************" << std::endl;
310
311 return EXIT_SUCCESS;
312}
313
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 average carrier energy at each point inside 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.
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()
@ 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