ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
mos1d-dg.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
54template <typename DeviceT>
55void generate_mos1d_mesh(DeviceT & mos_device,
56 double len_gate, // gate length
57 double cs_gate, // cell size in the gate
58 double len_oxide, // oxide length
59 double cs_ox, // cell size in the oxide
60 double len_bulk, // length of the silicon bulk
61 double cs_bulk) // cell size in the silicon bulk
62{
64
65 // The gate contact is segment #0:
66 gconf.add_segment(0, // starting location (meter)
67 len_gate, // length (meter)
68 static_cast<unsigned long>(len_gate/cs_gate + 1.0) ); // number of mesh points in the segment
69
70 // The oxide is segment #1
71 gconf.add_segment(len_gate,
72 len_oxide,
73 static_cast<unsigned long>(len_oxide/cs_ox) );
74
75 // The channel is segment #2 and has a fixed length of 20 nm and a cell size of 0.05 nm
76 gconf.add_segment(len_gate + len_oxide,
77 20e-9, // 20nm
78 static_cast<unsigned long>(+20e-9/0.05e-9) );
79
80 // The bulk silicon is segment #3
81 gconf.add_segment(len_gate + len_oxide + 20e-9,
82 len_bulk,
83 static_cast<unsigned long>(len_bulk/cs_bulk) );
84
85 // The last segment #4 is the bulk contact
86 gconf.add_segment(len_gate + len_oxide + 20e-9 + len_bulk,
87 2e-9,
88 static_cast<unsigned long>(2e-9/cs_bulk + 1.0) );
89
90 // Pass the
91 mos_device.generate_mesh(gconf);
92}
93
94
103template <typename DeviceType>
104void init_device(DeviceType & device,
105 double Vg_init, // The initial gate voltage
106 double Nd, // The donor doping (m^-3)
107 double Na) // The acceptor doping (m^-3)
108{
109 typedef typename DeviceType::segment_type SegmentType;
110
111 // Get the segments of the mesh
112 SegmentType const & gate = device.segment(0);
113 SegmentType const & oxide = device.segment(1);
114 SegmentType const & silicon = device.segment(2);
115 SegmentType const & silicon2 = device.segment(3);
116 SegmentType const & bulk = device.segment(4);
117
118 std::cout << "* init_device(): Setting material ..." << std::endl;
119
123 device.set_material(viennashe::materials::sio2(), oxide);
124 device.set_material(viennashe::materials::si(), silicon);
125 device.set_material(viennashe::materials::si(), silicon2);
126 //
127 // Contacts are always conductors/metals (for now)
128 device.set_material(viennashe::materials::metal(), gate);
129 device.set_material(viennashe::materials::metal(), bulk);
130
131
137 std::cout << "* init_device(): Setting doping (per cell) ..." << std::endl;
138 device.set_doping_n(Nd, silicon);
139 device.set_doping_p(Na, silicon);
140 device.set_doping_n(Nd, silicon2);
141 device.set_doping_p(Na, silicon2);
142
143
145 std::cout << "* init_device(): Setting contact potentials (per cell) ..." << std::endl;
146 device.set_contact_potential(Vg_init, gate);
147 device.set_contact_potential(0.0, bulk);
148
149 std::cout << "* init_device(): DONE!" << std::endl;
150
151} // init_device()
152
153
163int main()
164{
166 typedef viennagrid::line_1d_mesh MeshType;
167 typedef viennashe::device<MeshType> DeviceType;
168
169 std::cout << viennashe::preamble() << std::endl;
170
171 // Create the device object
172 DeviceType device;
173
175 std::cout << "* main(): Creating mesh ..." << std::endl;
176
177 // Configuration of the 1d MOS
178 const double len_gate = 1e-9;
179 const double cs_gate = 0.01e-9; // Gate cell size
180 const double len_oxide = 1e-9;
181 const double cs_ox = 0.05e-9; // Oxide cell size
182 const double len_bulk = 100.e-9;
183 const double cs_bulk = 1e-9; // Bulk cell size
184
185 double Vg = -0.2; // Gate voltage
186
187 // Generate the mesh
188 generate_mos1d_mesh(device, len_gate, cs_gate, len_oxide, cs_ox, len_bulk, cs_bulk);
189
191 // ND NA
192 init_device(device, Vg, 3e23, 1e8 ); // configure a p-channel MOS
193
194
201 std::cout << "* main(): Creating DD simulator..." << std::endl;
202
210 viennashe::config dd_cfg;
211 dd_cfg.with_holes(true);
212 dd_cfg.with_electrons(true);
215
216 // Linear solver configuration: We take a dense solver here ...
218
219 // Nonlinear solver: 30 Gummel iterations (Newton is also available)
220 dd_cfg.nonlinear_solver().max_iters(30);
221 dd_cfg.nonlinear_solver().damping(0.6);
222
223 // We will make use of a first-order quantum correction
224 dd_cfg.quantum_correction(true); // solve density gradient equations
225 dd_cfg.with_quantum_correction(true); // apply correction potentials from density gradient
226
227
234 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
235
236 std::cout << "* main(): Launching simulator..." << std::endl;
237 // Run the DD simulation
238 dd_simulator.run();
239
245 viennashe::io::write_cell_quantity_for_gnuplot(dd_simulator.dg_pot_p(), device, "mos1d-dg_dd_dgpot_holes.dat");
246 viennashe::io::write_cell_quantity_for_gnuplot(dd_simulator.dg_pot_n(), device, "mos1d-dg_dd_dgpot_electrons.dat");
247 viennashe::io::write_cell_quantity_for_gnuplot(dd_simulator.potential(), device, "mos1d-dg_dd_pot.dat");
248 viennashe::io::write_cell_quantity_for_gnuplot(dd_simulator.electron_density(), device, "mos1d-dg_dd_electrons.dat");
249 viennashe::io::write_cell_quantity_for_gnuplot(dd_simulator.hole_density(), device, "mos1d-dg_dd_holes.dat");
250
251
258 std::cout << "* main(): Setting up SHE..." << std::endl;
259
267 viennashe::config config;
268 // Configure a first-order SHE for holes only ...
269 config.with_holes(true);
271 // ... and use DD for the electrons
272 config.with_electrons(true);
274 config.max_expansion_order(1);
275
276 // Configure the various scattering mechanisms for SHE
277 config.scattering().acoustic_phonon().enabled(true);
278 config.scattering().optical_phonon().enabled(true);
279 config.scattering().ionized_impurity().enabled(true);
280 config.scattering().impact_ionization().enabled(false);
281 config.scattering().electron_electron(false);
282
283 // Configure the linear solver
284 config.linear_solver().max_iters(2000);
285 config.linear_solver().ilut_drop_tolerance(1e-2);
286
287 // Nonlinear solver configuration: 30 Gummel steps with damping parameter of 0.5
288 config.nonlinear_solver().max_iters(30);
289 config.nonlinear_solver().damping(0.5);
290
291 // We require a minimum kinetic energy range of 0.5 eV at each grid node
293
294 // Set the energy spacing to 31 meV
296
297 // Like for DD use density gradient as a first-order quantum correction model
298 config.quantum_correction(true); // calculate a correction potential
299 config.with_quantum_correction(true); // use the correction potential for the transport model
300
308 std::cout << "* main(): Computing SHE..." << std::endl;
309
310 viennashe::simulator<DeviceType> she_simulator(device, config);
311
312 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
313 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
314 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
315
316 // Run the simulator
317 she_simulator.run();
318
324 std::cout << "* main(): Writing SHE result..." << std::endl;
325
327 {
329 she_simulator.config(),
330 she_simulator.quantities().electron_distribution_function(),
331 "mos1d-dg_edf_n");
332 }
334 {
336 she_simulator.config(),
337 she_simulator.quantities().hole_distribution_function(),
338 "mos1d-dg_edf_p");
339 }
340
342 viennashe::io::write_cell_quantity_for_gnuplot(she_simulator.dg_pot_p(), device, "mos1d-dg_she_dgpot_holes.dat" );
343 viennashe::io::write_cell_quantity_for_gnuplot(she_simulator.dg_pot_n(), device, "mos1d-dg_she_dgpot_electrons.dat" );
344 viennashe::io::write_cell_quantity_for_gnuplot(she_simulator.potential(), device, "mos1d-dg_she_pot.dat" );
345 viennashe::io::write_cell_quantity_for_gnuplot(she_simulator.electron_density(), device, "mos1d-dg_she_electrons.dat");
346 viennashe::io::write_cell_quantity_for_gnuplot(she_simulator.hole_density(), device, "mos1d-dg_she_holes.dat");
347
349 std::cout << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)." << std::endl;
350 std::cout << "* main(): Don't forget to scale the z-axis by about a factor of 1e12 when examining the distribution function." << std::endl;
351 std::cout << std::endl;
352 std::cout << "*********************************************************" << std::endl;
353 std::cout << "* ViennaSHE finished successfully *" << std::endl;
354 std::cout << "*********************************************************" << std::endl;
355
356 return EXIT_SUCCESS;
357}
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
equation_id get_electron_equation() const
Definition: config.hpp:230
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
equation_id get_hole_equation() const
Definition: config.hpp:238
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
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
list Nd
Definition: resistor.py:59
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()
@ 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
A class referring to silicon dioxide.
Definition: all.hpp:111
static const double q
Elementary charge.
Definition: constants.hpp:44