ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
half-trigate-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 default configurations:
29#include "viennagrid/config/default_configs.hpp"
30
31#define N 16
32
71template <typename DeviceType>
72void init_device(DeviceType & device)
73{
74 typedef typename DeviceType::segment_type SegmentType;
75
77 SegmentType const & source_segment = device.segmentation()[1];
78 SegmentType const & channel_segment = device.segmentation()[2];
79 SegmentType const & drain_segment = device.segmentation()[3];
80 SegmentType const & oxide_segment = device.segmentation()[4];
81 SegmentType const & gate_segment = device.segmentation()[5];
82 SegmentType const & body_segment = device.segmentation()[6];
83
84 SegmentType const & source_contact_segment = device.segmentation()[7];
85 SegmentType const & drain_contact_segment = device.segmentation()[8];
86 SegmentType const & body_contact_segment = device.segmentation()[9];
87
90 std::cout << "* init_device(): Setting materials..." << std::endl;
91 device.set_material(viennashe::materials::metal(), gate_segment);
92 device.set_material(viennashe::materials::hfo2(), oxide_segment);
93 device.set_material(viennashe::materials::si(), source_segment);
94 device.set_material(viennashe::materials::si(), channel_segment);
95 device.set_material(viennashe::materials::si(), drain_segment);
96 device.set_material(viennashe::materials::si(), body_segment);
97
99 device.set_material(viennashe::materials::metal(), source_contact_segment);
100 device.set_material(viennashe::materials::metal(), drain_contact_segment);
101 device.set_material(viennashe::materials::metal(), body_contact_segment);
102
109 std::cout << "* init_device(): Setting doping..." << std::endl;
110
111 device.set_doping_n(1e18, channel_segment);
112 device.set_doping_p(1e14, channel_segment);
113
114 device.set_doping_n(1e26, source_segment);
115 device.set_doping_p(1e6, source_segment);
116
117 device.set_doping_n(1e26, drain_segment);
118 device.set_doping_p(1e6, drain_segment);
119
120 device.set_doping_n(1e18, body_segment);
121 device.set_doping_p(1e14, body_segment);
122
127 std::cout << "* init_device(): Setting contacts..." << std::endl;
128
129 double gnd = 0.0;
130 double vcc = 0.3;
131 double vgate = 0.8;
132
133 device.set_contact_potential(vgate, gate_segment);
134 device.set_contact_potential(gnd, source_contact_segment);
135 device.set_contact_potential(vcc, drain_contact_segment);
136 device.set_contact_potential(gnd, body_contact_segment);
137}
138
147int main(int argc, char** argv)
148{
154
155 std::cout << viennashe::preamble() << std::endl;
156
163 int size,opt,s,h;
164 std::cout << "* main(): Creating and scaling device..." << std::endl;
165 DeviceType device;
166 device.load_mesh("../examples/data/half-trigate57656.mesh");
167 device.scale(1e-9);
168 PetscInitialize(&argc, &argv, (char *) 0, "");
169 std::string LSolver("petsc_parallel_linear_solver");
170 MPI_Comm_size(MPI_COMM_WORLD, &size);
171 PetscInitialize(&argc, &argv, (char *) 0, "");
175 opt = getopt(argc, argv, "seS");
176 switch(opt){
177 case 's': s = 1;h = 0 ; break; // Space refining
178 case 'e': s = 0;h = 1 ; break; // Energy refining
179 case 'S': s = 0;h = 0 ; break; // Strong Scale
180 default:
181 break;
182 }
183
184 s = 0;h = 0 ;
185 std::cout << "* main(): Initializing device..." << std::endl;
186 init_device(device);
187 if(s == 1 && size != 1) device.refine(size) ;
196 std::cout << "* main(): Creating DD simulator..." << std::endl;
197
206 // Create the configuration object
207 viennashe::config dd_cfg;
208
209 // Specify that the drift-diffusion system for electrons and holes,
210 // both using the continuity equations, should be solved:
211 dd_cfg.with_electrons(true);
212 dd_cfg.with_holes(true);
215
219 // dd_cfg.nonlinear_solver().set(viennashe::solvers::nonlinear_solver_ids::newton_nonlinear_solver);
220
229 dd_cfg.nonlinear_solver().damping(0.4); // The damping factor
230 dd_cfg.nonlinear_solver().max_iters(50); // Maximum of nonlinear iterations
231
238 std::cout << "* main(): Creating and launching DD simulator..." << std::endl;
239 viennashe::simulator<DeviceType> dd_simulator(device, // the device created before
240 dd_cfg); // the configuration
241 dd_simulator.run();
242
250 viennashe::io::write_quantities_to_VTK_file(dd_simulator, // simulator object
251 "half-trigate_dd_quan"); // file name
252
253
275 std::cout << "* main(): Setting up first-order SHE (non-self consistent!)..." << std::endl;
276 viennashe::config config;
277
278 config.with_electrons(true); config.with_holes(true);
279 config.set_electron_equation(viennashe::EQUATION_SHE); // Use SHE for electrons
280 config.set_hole_equation(viennashe::EQUATION_CONTINUITY); // Use DD for holes
281
283 config.max_expansion_order(1);
286 h == 0 ?config.energy_spacing(0.031 * viennashe::physics::constants::q/N):
288 config.linear_solver().set(LSolver);
289 config.linear_solver().max_iters(2000);
290 config.linear_solver().setArgv(argv);
291 config.linear_solver().setArgc(argc);
293 config.nonlinear_solver().max_iters(1);
294
302 std::cout << "* main(): Computing first-order SHE (requires about 5 GB RAM)..." << std::endl;
303 viennashe::simulator<DeviceType> she_simulator(device, config);
304
305 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
306 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
307 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
308
309 she_simulator.run();
310
321 PetscFinalize();
322 viennashe::io::write_quantities_to_VTK_file(she_simulator, // simulator object
323 "half-trigate_she_quan"); // file name
324
327 viennashe::io::write_quantity_to_VTK_file(she_simulator.potential(), device, "trigate_she_potential");
328 viennashe::io::write_quantity_to_VTK_file(she_simulator.electron_density(), device, "trigate_she_electrons");
329
331 std::cout << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)." << std::endl;
332 std::cout << std::endl;
333 std::cout << "*********************************************************" << std::endl;
334 std::cout << "* ViennaSHE finished successfully *" << std::endl;
335 std::cout << "*********************************************************" << std::endl;
336
337 return EXIT_SUCCESS;
338
339}
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
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
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
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