ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
mosfet-PETSC.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#include <petscksp.h>
24#include <unistd.h>
25// ViennaSHE includes:
26#include "viennashe/core.hpp"
27
28// ViennaGrid mesh configurations:
29#include "viennagrid/config/default_configs.hpp"
30#define N 16
67template<typename DeviceType>
68void init_device(DeviceType & device)
69{
70 typedef typename DeviceType::segment_type SegmentType;
71
73 SegmentType const & gate_contact = device.segment(1);
74 SegmentType const & source_contact = device.segment(2);
75 SegmentType const & gate_oxide = device.segment(3);
76 SegmentType const & drain_contact = device.segment(4);
77 SegmentType const & source = device.segment(5);
78 SegmentType const & drain = device.segment(6);
79 SegmentType const & body = device.segment(7);
80 SegmentType const & body_contact = device.segment(8);
81
83 std::cout << "* init_device(): Setting materials..." << std::endl;
84 device.set_material(viennashe::materials::metal(), gate_contact);
85 device.set_material(viennashe::materials::metal(), source_contact);
86 device.set_material(viennashe::materials::metal(), drain_contact);
87 device.set_material(viennashe::materials::metal(), body_contact);
88
89 device.set_material(viennashe::materials::hfo2(), gate_oxide);
90
91 device.set_material(viennashe::materials::si(), source);
92 device.set_material(viennashe::materials::si(), drain);
93 device.set_material(viennashe::materials::si(), body);
94
101 std::cout << "* init_device(): Setting doping..." << std::endl;
102 device.set_doping_n(1e24, source);
103 device.set_doping_p(1e8, source);
104
105 device.set_doping_n(1e24, drain);
106 device.set_doping_p(1e8, drain);
107
108 device.set_doping_n(1e17, body);
109 device.set_doping_p(1e15, body);
110
115 device.set_contact_potential(0.8, gate_contact);
116 device.set_contact_potential(0.0, source_contact);
117 device.set_contact_potential(1.0, drain_contact);
118 device.set_contact_potential(0.0, body_contact);
119
120}
121
130int main(int argc, char **argv)
131{
136 int size,opt,s,h;
138 typedef DeviceType::segment_type SegmentType;
139 std::string LSolver("petsc_parallel_linear_solver");
140 std::cout << viennashe::preamble() << std::endl;
141 PetscInitialize(&argc, &argv, (char *) 0, "");
142 MPI_Comm_size(MPI_COMM_WORLD, &size);
143 opt = getopt(argc, argv, "seS");
144 switch(opt){
145 case 's': s = 1;h = 0 ; break; // Space refining
146 case 'e': s = 0;h = 1 ; break; // Energy refining
147 case 'S': s = 0;h = 0 ; break; // Strong Scale
148 default:
149 break;
150 }
151
152 s = 0;h = 0 ;
153// Set MPI here just to fix bug
154// int world_rank;
155// int error = MPI_Init(&argc, &argv);
156// MPI_Comm_rank(MPI_COMM_WORLD, &world_rank);
157
164 std::cout << "* main(): Creating and scaling device..." << std::endl;
165 DeviceType device;
166 device.load_mesh("../examples/data/mosfet840.mesh");
167
168 device.scale(1e-9);
169
173 std::cout << "* main(): Initializing device..." << std::endl;
174 init_device(device);
175 if(s == 1 && size != 1) device.refine(size) ;
176// if (world_rank == 0)
177// {
178
187 std::cout << "* main(): Creating DD simulator..." << std::endl;
188
197 viennashe::config dd_cfg;
198
199 // enable electrons and holes and specify that for each of them a continuity equation should be solved:
200 dd_cfg.with_electrons(true);
201 dd_cfg.with_holes(true);
202 dd_cfg.linear_solver().setArgv(argv);
203 dd_cfg.linear_solver().setArgc(argc);
204 // dd_cfg.linear_solver().set(LSolver);
205 // enable the use of PETSc on the simulator
208
209 // Nonlinear solver parameters: 200 Gummel iterations with rather high damping
210 //dd_cfg.nonlinear_solver().set("newton");
211 dd_cfg.nonlinear_solver().setThreshold(200);
212 dd_cfg.nonlinear_solver().max_iters(200);
213 dd_cfg.nonlinear_solver().damping(0.125);
214
221 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
222 std::cout << "* main(): Launching DD simulator..." << std::endl;
223 dd_simulator.run();
224
232 // viennashe::io::write_quantities_to_VTK_file(dd_simulator,
233 // "mosfet_petsc_dd_quan");
234
239 SegmentType const & drain_contact = device.segment(4);
240 SegmentType const & body = device.segment(6);
241
242// std::cout << "* main(): Drain electron current Id_e = "
243// << viennashe::get_terminal_current(device, viennashe::ELECTRON_TYPE_ID,
244// dd_simulator.potential(), dd_simulator.electron_density(),
245// viennashe::models::create_constant_mobility_model(device, 0.1430),
246// body, drain_contact) * 1e-6 << std::endl;
247// std::cout << "* main(): Drain hole current Id_h = "
248// << viennashe::get_terminal_current(device, viennashe::HOLE_TYPE_ID,
249// dd_simulator.potential(), dd_simulator.hole_density(),
250// viennashe::models::create_constant_mobility_model(device, 0.0460),
251// body, drain_contact) * 1e-6 << std::endl;
252
266 std::cout
267 << "* main(): Setting up first-order SHE (semi-self-consistent using 40 Gummel iterations)..."
268 << std::endl;
269// viennashe::config config (
270// true, true, false, false, false, true, false,
271// viennashe::solvers::linear_solver_ids::petsc_parallel_linear_solver,
272// viennashe::solvers::nonlinear_solver_ids::PETSC_nonlinear_solver);
273// }
274 viennashe::config config;
275 // Set the expansion order to 1 for SHE
276 config.max_expansion_order(1);
277
278 // Use both carrier types
279 config.with_electrons(true);
280 config.with_holes(true);
281
282 // Configure equation
283 config.set_electron_equation(viennashe::EQUATION_SHE); // SHE for electrons
285
286 // Energy spacing of 31 meV
287 h == 0 ?config.energy_spacing(0.031 * viennashe::physics::constants::q/N):
289
290 // Set the scattering mechanisms
291 config.scattering().acoustic_phonon().enabled(true);
292 config.scattering().optical_phonon().enabled(true);
293 config.scattering().ionized_impurity().enabled(false);
294
295 // The linear solver should run for at most 2000 iterations:
296 config.linear_solver().set(LSolver);
297 config.linear_solver().max_iters(2000);
298 config.linear_solver().setArgv(argv);
299 config.linear_solver().setArgc(argc);
300
301 // Configure the nonlinear solver to 40 Gummel iterations with a damping parameter 0.4
302 config.nonlinear_solver().setThreshold(800);
303 config.nonlinear_solver().max_iters(50);
304 // config.nonlinear_solver().set("newton");
305 config.nonlinear_solver().damping(0.4); // 0.2|800 - 0.4|276 - 0.6|182 - 0.8|135 - 1|108 - 1.1|97
306
314 std::cout << "* main(): Computing first-order SHE..." << std::endl;
315 viennashe::simulator<DeviceType> she_simulator(device, config);
316 // Set the previous DD solution as an initial guess
317 she_simulator.set_initial_guess(viennashe::quantity::potential(),
318 dd_simulator.potential());
319 she_simulator.set_initial_guess(viennashe::quantity::electron_density(),
320 dd_simulator.electron_density());
321 she_simulator.set_initial_guess(viennashe::quantity::hole_density(),
322 dd_simulator.hole_density());
323 // }
324// Trigger the actual simulation:
325 she_simulator.run();
326
327// Finalize PETSC
328// viennashe::solvers::PETSC_handler<NumericT,VectorType>::finalize();
329
335// if (world_rank == 0)
336// {
337 PetscFinalize();
338 std::cout
339 << "* main(): Writing energy distribution function from first-order SHE result..."
340 << std::endl;
341 // viennashe::io::she_vtk_writer<DeviceType>()(device, she_simulator.config(),
342 // she_simulator.quantities().electron_distribution_function(),
343 // "mosfet_petsc_she_edf");
344
346// viennashe::io::write_quantity_to_VTK_file(she_simulator.potential(), device,
347 // "mosfet_petsc_she_potential");
348 // viennashe::io::write_quantity_to_VTK_file(she_simulator.electron_density(),
349 // device, "mosfet_petsc_she_electrons");
350
353 // viennashe::io::write_quantities_to_VTK_file<DeviceType>(she_simulator,
354 // "mosfet_petsc_she_quan");
355
360// std::cout << "* main(): Drain hole current Id_e = "
361// << viennashe::get_terminal_current(device, config,
362// she_simulator.quantities().electron_distribution_function(), body,
363// drain_contact) * 1e-6 << std::endl;
364// std::cout << "Id_h = " << viennashe::get_terminal_current(
365// device, config, she_simulator.quantities().hole_distribution_function(), body, drain_contact ) * 1e-6
366// << std::endl;
367
369 std::cout
370 << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)."
371 << std::endl;
372 std::cout
373 << "* main(): Don't forget to scale the z-axis by about a factor of 1e12 when examining the distribution function."
374 << std::endl;
375 std::cout << std::endl;
376 std::cout << "*********************************************************"
377 << std::endl;
378 std::cout << "* ViennaSHE finished successfully *"
379 << std::endl;
380 std::cout << "*********************************************************"
381 << std::endl;
382
383 // }
384 // PetscFinalize();
385 //MPI_Finalize();
386 return EXIT_SUCCESS;
387}
388
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
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
ionized_impurity_scattering_parameters const & ionized_impurity() const
Returns the parameters for ionized impurity scattering. Const-version.
Definition: config.hpp:377
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
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 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
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