ViennaSHE 1.3.0
Free open-source semiconductor device simulator using spherical harmonics expansions techniques.
quantity_register.hpp
Go to the documentation of this file.
1#ifndef LIBVIENNASHE_QUANTITY_REGISTER_HPP
2#define LIBVIENNASHE_QUANTITY_REGISTER_HPP
3
4/* ============================================================================
5 Copyright (c) 2011-2022, Institute for Microelectronics,
6 Institute for Analysis and Scientific Computing,
7 TU Wien.
8
9 -----------------
10 ViennaSHE - The Vienna Spherical Harmonics Expansion Boltzmann Solver
11 -----------------
12
13 http://viennashe.sourceforge.net/
14
15 License: MIT (X11), see file LICENSE in the base directory
16=============================================================================== */
17
18
19
25// C++ includes
28
29// C includes
31
32#include <cstring>
33
34namespace libviennashe
35{
36 namespace detail
37 {
38
49 template <typename DeviceType,
50 typename PotentialQuantityType,
51 typename CarrierQuantityType,
52 typename MobilityModel>
54 DeviceType const & device,
55 PotentialQuantityType const & potential,
56 CarrierQuantityType const & carrier,
58 MobilityModel const & mobility_model, std::string name)
59 {
60 typedef typename DeviceType::mesh_type MeshType;
62 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
63
64 CurrentDensityAccessorType Jfield(device, ctype, potential, carrier, mobility_model);
65
67 reg.cell_based.register_quan(quan_cell);
68 }
69
78 template <typename DeviceType, typename SHEQuanT, typename ConfigT>
80 DeviceType const & device,
81 SHEQuanT const & shequan,
82 ConfigT const & conf, std::string name)
83 {
84 typedef typename DeviceType::mesh_type MeshType;
85 typedef typename viennashe::she::current_density_wrapper<DeviceType, SHEQuanT> CurrentDensityAccessorType;
86 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
87
88 CurrentDensityAccessorType Jfield(device, conf, shequan);
89
91 reg.cell_based.register_quan(quan_cell);
92 }
93
102 template <typename DeviceType, typename SHEQuanT, typename ConfigT>
104 DeviceType const & device,
105 SHEQuanT const & shequan,
106 ConfigT const & conf, std::string name)
107 {
108 typedef typename DeviceType::mesh_type MeshType;
110 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
111
112 AccessorType acc(device, conf, shequan);
114 reg.cell_based.register_quan(quan_cell);
115 }
116
125 template <typename DeviceType, typename SHEQuanT, typename ConfigT>
127 DeviceType const & device,
128 SHEQuanT const & shequan,
129 ConfigT const & conf, std::string name)
130 {
131 typedef typename DeviceType::mesh_type MeshType;
132 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
133 typedef typename viennashe::she::carrier_energy_wrapper<SHEQuanT> AccessorType;
134
135 AccessorType acc(conf, shequan);
136
138 reg.cell_based.register_quan(quan_cell);
139 }
140
148 template <typename DeviceType, typename SHEQuanT>
150 DeviceType const & device,
151 SHEQuanT const & shequan,
152 std::string name)
153 {
154 typedef typename DeviceType::mesh_type MeshType;
155 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
157
158 AccessorType acc(shequan, 0.0, 1.0);
159
161 reg.cell_based.register_quan(quan_cell);
162 }
163
164
165
173 template <typename DeviceType, typename PotentialAccessor>
174 void register_electric_field(quan_register_internal & reg, DeviceType const & device,
175 PotentialAccessor const & potential, std::string quantity_name)
176 {
177 typedef typename DeviceType::mesh_type MeshType;
178 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
179 typedef typename viennashe::electric_field_wrapper<DeviceType, PotentialAccessor> ElectricFieldAccessorType;
180 ElectricFieldAccessorType Efield(device, potential);
181
183 reg.cell_based.register_quan(quan_vertex);
184 }
185
193 template <typename DeviceType, typename PotentialAccessor>
194 void register_electric_flux(quan_register_internal & reg, DeviceType const & device,
195 PotentialAccessor const & potential, std::string quantity_name)
196 {
197 typedef typename DeviceType::mesh_type MeshType;
198 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
201
202 AccessorType Dfield(device, potential);
203 WrapperType quan_cell(Dfield, device, quantity_name);
204 reg.cell_based.register_quan(quan_cell);
205 }
206
207 } // namespace detail
208
214 template < typename SimulatorT >
215 void register_quans(SimulatorT const & sim, quan_register_internal & reg)
216 {
217 typedef typename SimulatorT::device_type DeviceType;
218 typedef typename SimulatorT::ResultQuantityType ResultQuantityType;
219
220 typedef typename DeviceType::mesh_type MeshType;
221 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
222
223 DeviceType const & device = sim.device();
224
225 //
226 // Register device quantities
227 //
228
229 // Doping
232 quantity::accessor_based_quantity_wrapper<DeviceType, viennashe::doping_accessor<DeviceType>, CellTag> doping_n_vertex(doping_n, device, "Donor doping concentration");
233 quantity::accessor_based_quantity_wrapper<DeviceType, viennashe::doping_accessor<DeviceType>, CellTag> doping_p_vertex(doping_p, device, "Acceptor doping concentration");
234 reg.cell_based.register_quan(doping_n_vertex);
235 reg.cell_based.register_quan(doping_p_vertex);
236
237 // Built-in Potential
239 quantity::accessor_based_quantity_wrapper<DeviceType, viennashe::built_in_potential_accessor<DeviceType>, CellTag> builtinpot_vt(builtinpot, device, "Built-In potential");
240 reg.cell_based.register_quan(builtinpot_vt);
241
242 //
243 // Register basic quantities
244 //
245
247
248 // Potential
249 ResultQuantityWrapperType pot_vertex(sim.potential(), device, viennashe::quantity::potential());
250 reg.cell_based.register_quan(pot_vertex);
251
252 // Electron Concentration
253 ResultQuantityWrapperType n_vertex(sim.electron_density(), device, viennashe::quantity::electron_density());
254 reg.cell_based.register_quan(n_vertex);
255
256 // Hole Concentration
257 ResultQuantityWrapperType p_vertex(sim.hole_density(), device, viennashe::quantity::hole_density());
258 reg.cell_based.register_quan(p_vertex);
259
260 // DG Electron Correction Potential
261 ResultQuantityWrapperType dgn_vertex(sim.dg_pot_n(), device, viennashe::quantity::density_gradient_electron_correction());
262 reg.cell_based.register_quan(dgn_vertex);
263
264 // DG Hole Correction Potential
265 ResultQuantityWrapperType dgp_vertex(sim.dg_pot_p(), device, viennashe::quantity::density_gradient_hole_correction());
266 reg.cell_based.register_quan(dgp_vertex);
267
268 // Lattice Temperature
269 ResultQuantityWrapperType hde_lat_temp(sim.quantities().lattice_temperature(), device, viennashe::quantity::lattice_temperature());
270 reg.cell_based.register_quan(hde_lat_temp);
271
274 tl_dev(lattice_temp, device, "Transport lattice temperature");
275 reg.cell_based.register_quan(tl_dev);
276
277 //
278 // Register postprocessing quantities
279 //
280
281 // Electric Field (V/m)
282 detail::register_electric_field(reg, device, sim.potential(), "Electric field");
283 // Electric Flux (As/m^2)
284 detail::register_electric_flux(reg, device, sim.potential(), "Electric flux density");
285
286 if (sim.config().with_electrons())
287 {
288 if (sim.config().get_electron_equation() == viennashe::EQUATION_CONTINUITY)
289 {
290 // Current Density Electrons
292 sim.quantities().potential(),
293 sim.quantities().get_unknown_quantity(viennashe::quantity::electron_density()),
295 viennashe::models::create_constant_mobility_model(device, 0.1430), "Electron current density DD");
296 }
297 else if (sim.config().get_electron_equation() == viennashe::EQUATION_SHE)
298 {
299 // Current Density Electrons
300 detail::register_SHE_current_density(reg, device, sim.quantities().electron_distribution_function(),
301 sim.config(), "Electron current density SHE");
302 // Average Carrier Energy
303 detail::register_average_energy(reg, device, sim.quantities().electron_distribution_function(),
304 sim.config(), "Average electron energy SHE");
305 // Average Carrier Velocity
306 detail::register_average_velocity(reg, device, sim.quantities().electron_distribution_function(),
307 sim.config(), "Average electron velocity SHE");
308
309 // Average Expansion Order
310 detail::register_average_expansion_order(reg, device, sim.quantities().electron_distribution_function(), "Average electron expansion order SHE");
311
312 }
313 else
314 {
315 // TODO: throw ?
316 viennashe::log::warn() << "libviennashe: register_quans(): Unknown equation type!" << std::endl;
317 }
318 }
319
320 if (sim.config().with_holes())
321 {
322 if (sim.config().get_hole_equation() == viennashe::EQUATION_CONTINUITY)
323 {
324 // Current Density Holes
326 sim.quantities().potential(),
327 sim.quantities().get_unknown_quantity(viennashe::quantity::hole_density()),
329 viennashe::models::create_constant_mobility_model(device, 0.0460), "Hole current density DD");
330 }
331 else if (sim.config().get_hole_equation() == viennashe::EQUATION_SHE)
332 {
333 // Current Density Holes
334 detail::register_SHE_current_density(reg, device, sim.quantities().hole_distribution_function(),
335 sim.config(), "Hole current density SHE");
336 // Average Carrier Energy
337 detail::register_average_energy(reg, device, sim.quantities().hole_distribution_function(),
338 sim.config(), "Average hole energy SHE");
339 // Average Carrier Velocity
340 detail::register_average_velocity(reg, device, sim.quantities().hole_distribution_function(),
341 sim.config(), "Average hole velocity SHE");
342 // Average Expansion Order
343 detail::register_average_expansion_order(reg, device, sim.quantities().hole_distribution_function(), "Average hole expansion order SHE");
344
345 }
346 else
347 {
348 // TODO: throw ?
349 viennashe::log::warn() << "libviennashe: register_quans(): Unknown equation type!" << std::endl;
350 }
351 }
352
353
354 } // register_quans
355
356
366 template <typename SHEQuanT, typename DFWrapperT, typename CellType >
367 void she_fill_edf_at_cell(SHEQuanT const & quan, DFWrapperT const & edfacc, CellType const & cell, double ** ekin, double ** edf, viennashe_index_type * len)
368 {
369 const std::size_t idx = std::size_t(cell.id().get());
370 len[idx] = static_cast<viennashe_index_type>(quan.get_value_H_size() - 1); // set length
371 for (std::size_t index_H = 1; index_H < quan.get_value_H_size() - 1; ++index_H)
372 {
373 const double energy = quan.get_kinetic_energy(cell, index_H);
374 const double f = edfacc(cell, energy, index_H);
375 ekin[idx][index_H - 1] = energy;
376 edf[idx][index_H - 1] = f;
377 }
378 }
379
380
391 template <typename DOSAccessorT, typename VertexType >
392 void she_fill_dos_at_cell(viennashe_index_type num, double deltaeps, DOSAccessorT const & dosacc, VertexType const & vt, double ** ekin, double ** dos, viennashe_index_type * len)
393 {
394 const std::size_t idx = std::size_t(vt.id().get());
395 len[idx] = num;
396 for (std::size_t index_H = 0; index_H < num; ++index_H)
397 {
398 const double energy = deltaeps * index_H;
399 ekin[idx][index_H - 1] = energy;
400 dos[idx][index_H - 1] = dosacc.density_of_states(energy);
401 }
402 }
403
412 template <typename SimulatorT >
413 void she_fill_edf(SimulatorT const & sim, viennashe::carrier_type_id ctype, double ** ekin, double ** edf, viennashe_index_type * len)
414 {
415 typedef typename SimulatorT::device_type DeviceType;
416 typedef typename DeviceType::mesh_type MeshType;
417
418 typedef typename viennagrid::result_of::const_cell_range<MeshType>::type CellContainer;
419 typedef typename viennagrid::result_of::iterator<CellContainer>::type CellIterator;
420
421 DeviceType const & device = sim.device();
422 CellContainer cells(device.mesh());
423
424 if (ctype == viennashe::ELECTRON_TYPE_ID)
425 {
426 for (CellIterator cit = cells.begin();
427 cit != cells.end();
428 ++cit)
429 {
430 she_fill_edf_at_cell(sim.quantities().electron_distribution_function(), sim.edf(ctype), *cit, ekin, edf, len);
431 }
432 }
433 else
434 {
435 for (CellIterator cit = cells.begin();
436 cit != cells.end();
437 ++cit)
438 {
439 she_fill_edf_at_cell(sim.quantities().hole_distribution_function(), sim.edf(ctype), *cit, ekin, edf, len);
440 }
441 }
442 } // she_fill_edf
443
452 template < typename SimulatorT >
453 void she_fill_dos(SimulatorT const & sim, viennashe::carrier_type_id ctype, double ** ekin, double ** dos, viennashe_index_type * len)
454 {
455 typedef typename SimulatorT::device_type DeviceType;
456 typedef typename DeviceType::mesh_type MeshType;
457
458 typedef typename viennagrid::result_of::const_vertex_range<MeshType>::type VertexContainer;
459 typedef typename viennagrid::result_of::iterator<VertexContainer>::type VertexIterator;
460
461 DeviceType const & device = sim.device();
462 VertexContainer vertices(device.mesh());
463
464 const double ekinmax = viennashe::physics::convert::eV_to_joule(10.0); // TODO: Think about better maximum for the kinetic energy ...
465 const viennashe_index_type num = static_cast<viennashe_index_type>(std::ceil(ekinmax / sim.config().energy_spacing()));
466
467 for (VertexIterator vit = vertices.begin(); vit != vertices.end(); ++vit)
468 {
469 she_fill_dos_at_cell(num, sim.config().energy_spacing(), sim.config().dispersion_relation(ctype), *vit, ekin, dos, len);
470 }
471
472 } // she_fill_dos
473
482 template < typename SimulatorT >
483 void she_fill_group_velocity(SimulatorT const & sim, viennashe::carrier_type_id ctype, double ** ekin, double ** vg, viennashe_index_type * len)
484 {
485 typedef typename SimulatorT::device_type DeviceType;
486 typedef typename DeviceType::mesh_type MeshType;
487
488 typedef typename viennagrid::result_of::const_vertex_range<MeshType>::type VertexContainer;
489 typedef typename viennagrid::result_of::iterator<VertexContainer>::type VertexIterator;
490
491 DeviceType const & device = sim.device();
492 VertexContainer vertices(device.mesh());
493
494 const double ekinmax = viennashe::physics::convert::eV_to_joule(10.0); // TODO: Think about better maximum for the kinetic energy ...
495
496 const viennashe_index_type num = static_cast<viennashe_index_type>(std::ceil(ekinmax / sim.config().energy_spacing()));
497
498 for (VertexIterator vit = vertices.begin(); vit != vertices.end(); ++vit)
499 {
500 const std::size_t idx = std::size_t(vit->id().get());
501 len[idx] = num;
502 for (std::size_t index_H = 0; index_H < num; ++index_H)
503 {
504 const double energy = sim.config().energy_spacing() * index_H;
505 ekin[idx][index_H - 1] = energy;
506 vg[idx][index_H - 1] = sim.config().dispersion_relation(ctype).velocity(energy);
507 }
508 }
509
510 } // she_fill_group_velocity
511
512
513} // namespace libviennashe
514
515
516#endif /* LIBVIENNASHE_QUANTITY_REGISTER_HPP */
517
Implements quantity_wrapper. Wraps vector valued quantities, which are accessible via a device based ...
Implements quantity_wrapper. Wraps scalar quantities, which are accessible via a device based accesso...
Accessor for retrieving the built-in potential in the device.
Definition: accessors.hpp:161
An accessor to the current density on vertices and edges (drift diffusion only!)
Accessor for returning the doping (donator/acceptor doping is defined in the CTOR)
Definition: accessors.hpp:141
Returns the lattice temperature on the device.
Definition: accessors.hpp:107
Accessor class providing the average expansion order inside the device.
An accessor for the average carrier energy at each point inside the device.
Accessor class providing the carrier velocity inside the device.
Accessor class providing the carrier velocity inside the device.
void register_SHE_current_density(quan_register_internal &reg, DeviceType const &device, SHEQuanT const &shequan, ConfigT const &conf, std::string name)
Registers a current density from a SHE simulator.
void register_average_energy(quan_register_internal &reg, DeviceType const &device, SHEQuanT const &shequan, ConfigT const &conf, std::string name)
Registers the average kinetic energy for a SHE simulator.
void register_average_expansion_order(quan_register_internal &reg, DeviceType const &device, SHEQuanT const &shequan, std::string name)
Registers the average expansion order (SHE simulator)
void register_electric_field(quan_register_internal &reg, DeviceType const &device, PotentialAccessor const &potential, std::string quantity_name)
Registers the electric field (vector) for any simulator.
void register_electric_flux(quan_register_internal &reg, DeviceType const &device, PotentialAccessor const &potential, std::string quantity_name)
Registers the electric flux density (vector) for any simulator.
void register_average_velocity(quan_register_internal &reg, DeviceType const &device, SHEQuanT const &shequan, ConfigT const &conf, std::string name)
Registers the average carrier velocity for a SHE simulator.
void register_DD_current_density(quan_register_internal &reg, DeviceType const &device, PotentialQuantityType const &potential, CarrierQuantityType const &carrier, viennashe::carrier_type_id ctype, MobilityModel const &mobility_model, std::string name)
Registers a current density from a DD simulator.
The internal C++ namespace of the library.
void she_fill_dos_at_cell(viennashe_index_type num, double deltaeps, DOSAccessorT const &dosacc, VertexType const &vt, double **ekin, double **dos, viennashe_index_type *len)
Fills the given C-arrays ekin and dos with the DOS at a vertex.
void she_fill_dos(SimulatorT const &sim, viennashe::carrier_type_id ctype, double **ekin, double **dos, viennashe_index_type *len)
Fills the given C-arrays ekin and dos with the DOS for every vertex.
void she_fill_group_velocity(SimulatorT const &sim, viennashe::carrier_type_id ctype, double **ekin, double **vg, viennashe_index_type *len)
Fills the given C-arrays ekin and vg with the group velocity for every vertex.
void register_quans(SimulatorT const &sim, quan_register_internal &reg)
Main quantity regsiter function. Add your simulator quantities here.
void she_fill_edf(SimulatorT const &sim, viennashe::carrier_type_id ctype, double **ekin, double **edf, viennashe_index_type *len)
Fills the given C-arrays with the complete EDF.
void she_fill_edf_at_cell(SHEQuanT const &quan, DFWrapperT const &edfacc, CellType const &cell, double **ekin, double **edf, viennashe_index_type *len)
Fills the given C-arrays with the EDF at a vertex.
reg
This is DD, so we do not need to give any inital guess.
Definition: resistor.py:84
logger< true > warn()
Used to log warnings. The logging level is logWARNING.
Definition: log.hpp:303
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
double eV_to_joule(const double eps)
Definition: physics.hpp:54
std::string lattice_temperature()
std::string density_gradient_hole_correction()
std::string density_gradient_electron_correction()
std::string electron_density()
carrier_type_id
Enumeration type for selecting the carrier type.
Definition: forwards.h:185
@ 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
Routines to wrap ViennaSHE quantities (doping, potential, distribution functions, ....
C++ to C wrapper of the quantity registry.
An accessor to the electric field on vertices and edges. Potential requiered.
An electric flux accessor. Provides the electric flux along edges and on vertices given the electrost...
unsigned long viennashe_index_type
Definition: sys.h:42
Contains all viennashe includes and most of the C++/C-wrappers. Contains macros.