ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
half-trigate.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 configurations:
28#include "viennagrid/config/default_configs.hpp"
29
30
69template <typename DeviceType>
70void init_device(DeviceType & device)
71{
72 typedef typename DeviceType::segment_type SegmentType;
73
75 SegmentType const & source_segment = device.segmentation()[1];
76 SegmentType const & channel_segment = device.segmentation()[2];
77 SegmentType const & drain_segment = device.segmentation()[3];
78 SegmentType const & oxide_segment = device.segmentation()[4];
79 SegmentType const & gate_segment = device.segmentation()[5];
80 SegmentType const & body_segment = device.segmentation()[6];
81
82 SegmentType const & source_contact_segment = device.segmentation()[7];
83 SegmentType const & drain_contact_segment = device.segmentation()[8];
84 SegmentType const & body_contact_segment = device.segmentation()[9];
85
88 std::cout << "* init_device(): Setting materials..." << std::endl;
89 device.set_material(viennashe::materials::metal(), gate_segment);
90 device.set_material(viennashe::materials::hfo2(), oxide_segment);
91 device.set_material(viennashe::materials::si(), source_segment);
92 device.set_material(viennashe::materials::si(), channel_segment);
93 device.set_material(viennashe::materials::si(), drain_segment);
94 device.set_material(viennashe::materials::si(), body_segment);
95
97 device.set_material(viennashe::materials::metal(), source_contact_segment);
98 device.set_material(viennashe::materials::metal(), drain_contact_segment);
99 device.set_material(viennashe::materials::metal(), body_contact_segment);
100
107 std::cout << "* init_device(): Setting doping..." << std::endl;
108
109 device.set_doping_n(1e18, channel_segment);
110 device.set_doping_p(1e14, channel_segment);
111
112 device.set_doping_n(1e26, source_segment);
113 device.set_doping_p(1e6, source_segment);
114
115 device.set_doping_n(1e26, drain_segment);
116 device.set_doping_p(1e6, drain_segment);
117
118 device.set_doping_n(1e18, body_segment);
119 device.set_doping_p(1e14, body_segment);
120
125 std::cout << "* init_device(): Setting contacts..." << std::endl;
126
127 double gnd = 0.0;
128 double vcc = 0.3;
129 double vgate = 0.8;
130
131 device.set_contact_potential(vgate, gate_segment);
132 device.set_contact_potential(gnd, source_contact_segment);
133 device.set_contact_potential(vcc, drain_contact_segment);
134 device.set_contact_potential(gnd, body_contact_segment);
135}
136
145int main()
146{
152
153 std::cout << viennashe::preamble() << std::endl;
154
161 std::cout << "* main(): Creating and scaling device..." << std::endl;
162 DeviceType device;
163 try
164 {
165 device.load_mesh("../examples/data/half-trigate57656.mesh");
166 }
167 catch (std::runtime_error const & e)
168 {
169 std::cerr << "-----------------------------------------------------------" << std::endl;
170 std::cerr << "--- NOTE: Please download the mesh file from:" << std::endl;
171 std::cerr << "--- http://viennashe.sourceforge.net/half-trigate57656.mesh" << std::endl;
172 std::cerr << "--- and place it in folder ../examples/data/" << std::endl;
173 std::cerr << "--- (this way the repository remains light-weight)" << std::endl;
174 std::cerr << "-----------------------------------------------------------" << std::endl;
175 throw e;
176 }
177
178 device.scale(1e-9);
179
180
184 std::cout << "* main(): Initializing device..." << std::endl;
185 init_device(device);
186
195 std::cout << "* main(): Creating DD simulator..." << std::endl;
196
205 // Create the configuration object
206 viennashe::config dd_cfg;
207
208 // Specify that the drift-diffusion system for electrons and holes,
209 // both using the continuity equations, should be solved:
210 dd_cfg.with_electrons(true);
211 dd_cfg.with_holes(true);
214
218 // dd_cfg.nonlinear_solver().set(viennashe::solvers::nonlinear_solver_ids::newton_nonlinear_solver);
219
228 dd_cfg.nonlinear_solver().damping(0.4); // The damping factor
229 dd_cfg.nonlinear_solver().max_iters(50); // Maximum of nonlinear iterations
230
237 std::cout << "* main(): Creating and launching DD simulator..." << std::endl;
238 viennashe::simulator<DeviceType> dd_simulator(device, // the device created before
239 dd_cfg); // the configuration
240 dd_simulator.run();
241
249 viennashe::io::write_quantities_to_VTK_file(dd_simulator, // simulator object
250 "half-trigate_dd_quan"); // file name
251
252
274 std::cout << "* main(): Setting up first-order SHE (non-self consistent!)..." << std::endl;
275 viennashe::config config;
276
277 config.with_electrons(true); config.with_holes(true);
278 config.set_electron_equation(viennashe::EQUATION_SHE); // Use SHE for electrons
279 config.set_hole_equation(viennashe::EQUATION_CONTINUITY); // Use DD for holes
280
282 config.max_expansion_order(1);
286
288 config.nonlinear_solver().max_iters(1);
289
297 std::cout << "* main(): Computing first-order SHE (requires about 5 GB RAM)..." << std::endl;
298 viennashe::simulator<DeviceType> she_simulator(device, config);
299
300 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
301 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
302 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
303
304 she_simulator.run();
305
316 viennashe::io::write_quantities_to_VTK_file(she_simulator, // simulator object
317 "half-trigate_she_quan"); // file name
318
321 viennashe::io::write_quantity_to_VTK_file(she_simulator.potential(), device, "trigate_she_potential");
322 viennashe::io::write_quantity_to_VTK_file(she_simulator.electron_density(), device, "trigate_she_electrons");
323
325 std::cout << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)." << std::endl;
326 std::cout << std::endl;
327
328 std::cout << "*********************************************************" << std::endl;
329 std::cout << "* ViennaSHE finished successfully *" << std::endl;
330 std::cout << "*********************************************************" << std::endl;
331
332 return EXIT_SUCCESS;
333
334}
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
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
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 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.
void write_quantity_to_VTK_file(QuantityType const &quantity, DeviceType const &device, std::string filename, std::string name_in_file="viennashe_quantity")
Generic interface function for writing a quantity to a VTK file. Automatically dispatches between ver...
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 hafnium dioxide.
Definition: all.hpp:119
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