ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
ushape_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 default configurations:
28#include "viennagrid/config/default_configs.hpp"
29
42template <typename DeviceType>
43void init_device(DeviceType & device)
44{
45 typedef typename DeviceType::segment_type SegmentType;
46
47 SegmentType const & contact_left = device.segment(1);
48 SegmentType const & oxide = device.segment(2);
49 SegmentType const & contact_right = device.segment(3);
50 SegmentType const & body = device.segment(4);
51
52 device.set_material(viennashe::materials::si(), body);
53
54 device.set_material(viennashe::materials::hfo2(), oxide);
55
56 device.set_material(viennashe::materials::metal(), contact_left);
57 device.set_material(viennashe::materials::metal(), contact_right);
58
59 device.set_doping_n(1e24, body);
60 device.set_doping_p(1e8, body);
61
62
63 // Set contact potentials
64 device.set_contact_potential(0.0, contact_left);
65 device.set_contact_potential(0.5, contact_right);
66
67}
68
69int main()
70{
71 typedef viennagrid::triangular_2d_mesh MeshType;
72 typedef viennashe::device<MeshType> DeviceType;
73 typedef DeviceType::segment_type SegmentType;
74
75 std::cout << viennashe::preamble() << std::endl;
76
77 std::cout << "* main(): Creating device..." << std::endl;
78 DeviceType device;
79 device.load_mesh("../../tests/data/ushape2d/ushape125.mesh");
80
81 std::cout << "Vertices: " << viennagrid::vertices(device.mesh()).size() << std::endl;
82 std::cout << "Edges: " << viennagrid::edges(device.mesh()).size() << std::endl;
83
84 // scale device from meter to nanometer
85 device.scale(1e-6);
86
87 init_device(device);
88
89 //
90 // Use a drift-diffusion simulation for obtaining an initial guess of the potential:
91 //
92 std::cout << "* main(): Creating DD simulator..." << std::endl;
93 viennashe::config dd_cfg;
94 dd_cfg.with_electrons(true);
95 dd_cfg.with_holes(true);
96 dd_cfg.nonlinear_solver().max_iters(300);
97 dd_cfg.nonlinear_solver().damping(0.4);
98 viennashe::simulator<DeviceType> dd_simulator(device, dd_cfg);
99
100 std::cout << "* main(): Launching DD simulator..." << std::endl;
101 dd_simulator.run();
102
103 //
104 // Write drift-diffusion quantities to file for inspection (e.g. using ParaView)
105 //
106 viennashe::io::write_quantities_to_VTK_file(dd_simulator, "ushape2d_dd_quan");
107
108 // Check current conservation:
109 std::cout << "#" << std::endl;
110 std::cout << "# Checking electron current: " << std::endl;
111 std::cout << "#" << std::endl;
114 dd_simulator.potential(),
115 dd_simulator.electron_density(),
117
118 std::cout << "#" << std::endl;
119 std::cout << "# Checking hole current: " << std::endl;
120 std::cout << "#" << std::endl;
123 dd_simulator.potential(),
124 dd_simulator.hole_density(),
126
127
128 //
129 // Self-consistent SHE simulations
130 //
131
132 std::cout << "* main(): Setting up SHE..." << std::endl;
133
134 //std::size_t energy_points = 100;
135 viennashe::config config;
137 config.with_electrons(true);
139 config.with_holes(true);
140 config.nonlinear_solver().max_iters(40); //Use higher iteration counts as needed
141 config.nonlinear_solver().damping(0.4);
142 config.max_expansion_order(1);
143 //config.she_boundary_conditions().type(viennashe::BOUNDARY_DIRICHLET);
144 //config.energy_levels(energy_points);
145 config.energy_spacing(31.0 * viennashe::physics::constants::q / 1000.0); //energy spacing of 31 meV
147 //config.linear_solver().set(viennashe::solvers::linear_solver_ids::parallel_linear_solver); //use this if OpenMP is available
148 //config.linear_solver().set(viennashe::solvers::linear_solver_ids::gpu_parallel_linear_solver); //use this for GPU support
149
150 //
151 // Instantiate and launch the SHE simulator
152 //
153 std::cout << "* main(): Computing SHE..." << std::endl;
154 viennashe::simulator<DeviceType> she_simulator(device, config);
155 she_simulator.set_initial_guess(viennashe::quantity::potential(), dd_simulator.potential());
156 she_simulator.set_initial_guess(viennashe::quantity::electron_density(), dd_simulator.electron_density());
157 she_simulator.set_initial_guess(viennashe::quantity::hole_density(), dd_simulator.hole_density());
158 she_simulator.run();
159
160 //
161 // Write output
162 //
163 std::cout << "* main(): Writing SHE result..." << std::endl;
164
165 // Writes EDF to VTK file
167 she_simulator.config(),
168 she_simulator.quantities().electron_distribution_function(),
169 "ushape2d_edf");
170
171 viennashe::io::write_quantities_to_VTK_file(she_simulator, "ushape2d_she_quan");
172
173
174 SegmentType const & contact_left = device.segment(1);
175 SegmentType const & contact_right = device.segment(3);
176 SegmentType const & semiconductor = device.segment(4);
177
178 double current_left = viennashe::get_terminal_current(device,
179 she_simulator.config(),
180 she_simulator.quantities().electron_distribution_function(),
181 semiconductor,
182 contact_left);
183 double current_right = viennashe::get_terminal_current(device,
184 she_simulator.config(),
185 she_simulator.quantities().electron_distribution_function(),
186 semiconductor,
187 contact_right);
188
189 std::cout << "Current left: " << current_left << std::endl;
190 std::cout << "Current right: " << current_right << std::endl;
191
193 she_simulator.config(),
194 she_simulator.quantities().electron_distribution_function());
195
196 if ( std::fabs(current_left + current_right) / std::max(std::fabs(current_left), std::fabs(current_right)) > 1e-5)
197 {
198 std::cerr << "CURRENT MISMATCH DETECTED" << std::endl;
199 std::cerr << "Norm of error: " << std::fabs(current_left - current_right) / std::max(std::fabs(current_left), std::fabs(current_right)) << std::endl;
200 return EXIT_FAILURE;
201 }
202 else
203 std::cout << "SUCCESS: CURRENT MATCHES!" << std::endl;
204
205 std::cout << "* main(): Results can now be viewed with your favorite VTK viewer (e.g. ParaView)." << std::endl;
206 std::cout << "* main(): Don't forget to scale the z-axis by about a factor of 1e11 when examining the distribution function." << std::endl;
207 std::cout << std::endl;
208 std::cout << "****************************************************" << std::endl;
209 std::cout << "* Test finished successfully *" << std::endl;
210 std::cout << "****************************************************" << std::endl;
211
212 return EXIT_SUCCESS;
213}
214
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 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.
viennashe::models::dd::mobility< DeviceType > create_constant_mobility_model(DeviceType const &device, double mu)
Returns a mobility model, which always yields the same mobility.
Definition: mobility.hpp:61
std::string electron_density()
void check_current_conservation(DeviceT const &device, viennashe::config const &conf, SHEQuantity const &quan)
Checks current conservation for SHE. Writes information using log::info().
double get_terminal_current(DeviceT const &device, CurrentAccessorT const &current_accessor, SegmentT const &semiconductor, SegmentT const &terminal)
Returns the terminal current for a number of given vertices. Considers carrier flux and displacement ...
@ HOLE_TYPE_ID
Definition: forwards.h:188
@ ELECTRON_TYPE_ID
Definition: forwards.h:187
@ EQUATION_SHE
Definition: forwards.h:118
@ EQUATION_CONTINUITY
Definition: forwards.h:117
void check_current_conservation(DeviceT const &device, CurrentDensityT const &current_on_facet)
Checks current conservation for SHE. Writes information using log::info().
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
@ serial_linear_solver
Gauss solver (use for systems up to ~1000-by-1000 only)
Definition: config.hpp:44