ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
nin-diode-2d-hde.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 default mesh configurations:
28#include "viennagrid/config/default_configs.hpp"
29
30
56template <typename DeviceType>
57void init_device(DeviceType & device)
58{
59 typedef typename DeviceType::mesh_type MeshType;
60 typedef typename DeviceType::segment_type SegmentType;
61
63 SegmentType const & contact_left = device.segment(1);
64 SegmentType const & i_center = device.segment(3);
65 SegmentType const & contact_right = device.segment(5);
66
68 device.set_material(viennashe::materials::si());
69 device.set_doping_n(1e24);
70 device.set_doping_p(1e8);
71
73 device.set_doping_n(1e21, i_center);
74 device.set_doping_p(1e11, i_center);
75
77 device.set_material(viennashe::materials::metal(), contact_left);
78 device.set_material(viennashe::materials::metal(), contact_right);
79
80
82 device.set_contact_potential(0.0, contact_left);
83 device.set_contact_potential(0.5, contact_right);
84
89 typedef typename viennagrid::result_of::const_cell_range<MeshType>::type CellContainer;
90 typedef typename viennagrid::result_of::iterator<CellContainer>::type CellIterator;
91
92 CellContainer cells(device.mesh());
93 for (CellIterator cit = cells.begin();
94 cit != cells.end();
95 ++cit)
96 {
97 device.set_lattice_temperature(300.0, *cit);
98 }
99}
100
101
110int main()
111{
116 typedef viennagrid::triangular_2d_mesh MeshType;
117 typedef viennashe::device<MeshType> DeviceType;
118
119 std::cout << viennashe::preamble() << std::endl;
120
127 std::cout << "* main(): Creating device..." << std::endl;
128 DeviceType device;
129 device.load_mesh("../examples/data/nin2d.mesh");
130 device.scale(1e-8);
131
135 std::cout << "* main(): Creating device..." << std::endl;
136 init_device(device);
137
146 std::cout << "* main(): Creating DD simulator..." << std::endl;
147
155 viennashe::config dd_cfg;
156
157 // Nonlinear solver: Up to 80 Gummel iterations:
158 dd_cfg.nonlinear_solver().max_iters(80);
159
160 // Here we enable the heat equation to solve for temperature explicitly:
161 dd_cfg.with_hde(true);
162
163
170 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
171
172 std::cout << "* main(): Launching DD simulator..." << std::endl;
173 dd_simulator.run();
174
182 viennashe::io::write_quantities_to_VTK_file(dd_simulator, "nin2d_dd_quan");
183
184
198 std::cout << "* main(): Setting up SHE..." << std::endl;
199 viennashe::config config;
200
201 // SHE for electrons
202 config.with_electrons(true);
204
205 // SHE for holes
206 config.with_holes(true);
208
209 // Linear solver: Up to 1000 iterations:
210 config.linear_solver().max_iters(1000);
211
212 // Nonlinear solver: Up to 10 Gummel iterations with a damping parameter of 0.5:
213 config.nonlinear_solver().max_iters(10); //Use higher iteration counts as needed
214 config.nonlinear_solver().damping(0.5);
215
216 // Explicitly enable the Heat Diffusion Equation:
217 config.with_hde(true);
218
219 // SHE: Specify a first-order simulation with an energy spacing of 31 meV
220 config.max_expansion_order(1);
221 config.energy_spacing(31.0 * viennashe::physics::constants::q / 1000.0);
222
223
231 std::cout << "* main(): Computing SHE..." << std::endl;
232 viennashe::simulator<DeviceType> she_simulator(device, config);
233 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
234 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
235 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
236 she_simulator.run();
237
243 std::cout << "* main(): Writing SHE result..." << std::endl;
244
246 she_simulator.config(),
247 she_simulator.quantities().electron_distribution_function(),
248 "nin2d-hde_edf");
249
252 viennashe::io::write_quantities_to_VTK_file(she_simulator, "nin2d-hde_she_quan");
253
255 std::cout << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)." << std::endl;
256 std::cout << "* main(): Don't forget to scale the z-axis by about a factor of 1e11 when examining the distribution function." << std::endl;
257 std::cout << std::endl;
258 std::cout << "*********************************************************" << std::endl;
259 std::cout << "* ViennaSHE finished successfully *" << std::endl;
260 std::cout << "*********************************************************" << std::endl;
261
262 return EXIT_SUCCESS;
263}
264
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_hde() const
Definition: config.hpp:502
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
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
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
@ 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
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