ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
nin-diode-2d.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
57template <typename DeviceType>
58void init_device(DeviceType & device)
59{
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
81 device.set_contact_potential(0.0, contact_left);
82 device.set_contact_potential(0.5, contact_right);
83}
84
85
94int main()
95{
101
102 std::cout << viennashe::preamble() << std::endl;
103
110 std::cout << "* main(): Creating device..." << std::endl;
111 DeviceType device;
112 device.load_mesh("../examples/data/nin2d.mesh");
113 device.scale(1e-8);
114
118 std::cout << "* main(): Creating device..." << std::endl;
119 init_device(device);
120
121
130 std::cout << "* main(): Creating DD simulator..." << std::endl;
131
132
140 viennashe::config dd_cfg;
141
142 // Nonlinear solver: Use up to 40 Gummel iterations:
143 dd_cfg.nonlinear_solver().max_iters(40);
144
145
152 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
153 std::cout << "* main(): Launching DD simulator..." << std::endl;
154 dd_simulator.run();
155
163 viennashe::io::write_quantities_to_VTK_file(dd_simulator, "nin2d_dd_quan");
164
165
166
180 std::cout << "* main(): Setting up SHE..." << std::endl;
181
182 viennashe::config config;
183
184 // Specify SHE for electrons:
185 config.with_electrons(true);
187
188 // Nonlinear solver: Up to 20 Gummel iterations with moderate damping of 0.5
189 config.nonlinear_solver().max_iters(20);
190 config.nonlinear_solver().damping(0.5);
191
192 // SHE: Maximum expansion order 1 with an energy spacing of 31 meV:
193 config.max_expansion_order(1);
194 config.energy_spacing(31.0 * viennashe::physics::constants::q / 1000.0);
195
196
204 std::cout << "* main(): Computing SHE..." << std::endl;
205 viennashe::simulator<DeviceType> she_simulator(device, config);
206 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
207 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
208 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
209 she_simulator.run();
210
216 std::cout << "* main(): Writing SHE result..." << std::endl;
217
219 she_simulator.config(),
220 she_simulator.quantities().electron_distribution_function(),
221 "nin2d_edf");
222
225 viennashe::io::write_quantities_to_VTK_file(she_simulator, "nin2d_she_quan");
226
228 std::cout << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)." << std::endl;
229 std::cout << "* main(): Don't forget to scale the z-axis by about a factor of 1e11 when examining the distribution function." << std::endl;
230 std::cout << std::endl;
231 std::cout << "*********************************************************" << std::endl;
232 std::cout << "* ViennaSHE finished successfully *" << std::endl;
233 std::cout << "*********************************************************" << std::endl;
234
235 return EXIT_SUCCESS;
236}
237
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
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
double energy_spacing() const
Returns the uniform energy spacing of discrete energies.
Definition: config.hpp:460
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 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
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