1#ifndef LIBVIENNASHE_QUANTITY_REGISTER_HPP
2#define LIBVIENNASHE_QUANTITY_REGISTER_HPP
49 template <
typename DeviceType,
50 typename PotentialQuantityType,
51 typename CarrierQuantityType,
52 typename MobilityModel>
54 DeviceType
const & device,
56 CarrierQuantityType
const & carrier,
58 MobilityModel
const & mobility_model, std::string name)
60 typedef typename DeviceType::mesh_type MeshType;
62 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
64 CurrentDensityAccessorType Jfield(device, ctype,
potential, carrier, mobility_model);
67 reg.cell_based.register_quan(quan_cell);
78 template <
typename DeviceType,
typename SHEQuanT,
typename ConfigT>
80 DeviceType
const & device,
81 SHEQuanT
const & shequan,
82 ConfigT
const & conf, std::string name)
84 typedef typename DeviceType::mesh_type MeshType;
86 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
88 CurrentDensityAccessorType Jfield(device, conf, shequan);
91 reg.cell_based.register_quan(quan_cell);
102 template <
typename DeviceType,
typename SHEQuanT,
typename ConfigT>
104 DeviceType
const & device,
105 SHEQuanT
const & shequan,
106 ConfigT
const & conf, std::string name)
108 typedef typename DeviceType::mesh_type MeshType;
110 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
112 AccessorType acc(device, conf, shequan);
114 reg.cell_based.register_quan(quan_cell);
125 template <
typename DeviceType,
typename SHEQuanT,
typename ConfigT>
127 DeviceType
const & device,
128 SHEQuanT
const & shequan,
129 ConfigT
const & conf, std::string name)
131 typedef typename DeviceType::mesh_type MeshType;
132 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
135 AccessorType acc(conf, shequan);
138 reg.cell_based.register_quan(quan_cell);
148 template <
typename DeviceType,
typename SHEQuanT>
150 DeviceType
const & device,
151 SHEQuanT
const & shequan,
154 typedef typename DeviceType::mesh_type MeshType;
155 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
158 AccessorType acc(shequan, 0.0, 1.0);
161 reg.cell_based.register_quan(quan_cell);
173 template <
typename DeviceType,
typename PotentialAccessor>
175 PotentialAccessor
const &
potential, std::string quantity_name)
177 typedef typename DeviceType::mesh_type MeshType;
178 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
180 ElectricFieldAccessorType Efield(device,
potential);
183 reg.cell_based.register_quan(quan_vertex);
193 template <
typename DeviceType,
typename PotentialAccessor>
195 PotentialAccessor
const &
potential, std::string quantity_name)
197 typedef typename DeviceType::mesh_type MeshType;
198 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
203 WrapperType quan_cell(Dfield, device, quantity_name);
204 reg.cell_based.register_quan(quan_cell);
214 template <
typename SimulatorT >
217 typedef typename SimulatorT::device_type DeviceType;
218 typedef typename SimulatorT::ResultQuantityType ResultQuantityType;
220 typedef typename DeviceType::mesh_type MeshType;
221 typedef typename viennagrid::result_of::cell_tag<MeshType>::type CellTag;
223 DeviceType
const & device = sim.device();
234 reg.cell_based.register_quan(doping_n_vertex);
235 reg.cell_based.register_quan(doping_p_vertex);
240 reg.cell_based.register_quan(builtinpot_vt);
250 reg.cell_based.register_quan(pot_vertex);
254 reg.cell_based.register_quan(n_vertex);
258 reg.cell_based.register_quan(p_vertex);
262 reg.cell_based.register_quan(dgn_vertex);
266 reg.cell_based.register_quan(dgp_vertex);
270 reg.cell_based.register_quan(hde_lat_temp);
274 tl_dev(lattice_temp, device,
"Transport lattice temperature");
275 reg.cell_based.register_quan(tl_dev);
286 if (sim.config().with_electrons())
292 sim.quantities().potential(),
301 sim.config(),
"Electron current density SHE");
304 sim.config(),
"Average electron energy SHE");
307 sim.config(),
"Average electron velocity SHE");
316 viennashe::log::warn() <<
"libviennashe: register_quans(): Unknown equation type!" << std::endl;
320 if (sim.config().with_holes())
326 sim.quantities().potential(),
335 sim.config(),
"Hole current density SHE");
338 sim.config(),
"Average hole energy SHE");
341 sim.config(),
"Average hole velocity SHE");
349 viennashe::log::warn() <<
"libviennashe: register_quans(): Unknown equation type!" << std::endl;
366 template <
typename SHEQuanT,
typename DFWrapperT,
typename CellType >
369 const std::size_t idx = std::size_t(cell.id().get());
371 for (std::size_t index_H = 1; index_H < quan.get_value_H_size() - 1; ++index_H)
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;
391 template <
typename DOSAccessorT,
typename VertexType >
394 const std::size_t idx = std::size_t(vt.id().get());
396 for (std::size_t index_H = 0; index_H < num; ++index_H)
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);
412 template <
typename SimulatorT >
415 typedef typename SimulatorT::device_type DeviceType;
416 typedef typename DeviceType::mesh_type MeshType;
418 typedef typename viennagrid::result_of::const_cell_range<MeshType>::type CellContainer;
419 typedef typename viennagrid::result_of::iterator<CellContainer>::type CellIterator;
421 DeviceType
const & device = sim.device();
422 CellContainer cells(device.mesh());
426 for (CellIterator cit = cells.begin();
430 she_fill_edf_at_cell(sim.quantities().electron_distribution_function(), sim.edf(ctype), *cit, ekin, edf, len);
435 for (CellIterator cit = cells.begin();
439 she_fill_edf_at_cell(sim.quantities().hole_distribution_function(), sim.edf(ctype), *cit, ekin, edf, len);
452 template <
typename SimulatorT >
455 typedef typename SimulatorT::device_type DeviceType;
456 typedef typename DeviceType::mesh_type MeshType;
458 typedef typename viennagrid::result_of::const_vertex_range<MeshType>::type VertexContainer;
459 typedef typename viennagrid::result_of::iterator<VertexContainer>::type VertexIterator;
461 DeviceType
const & device = sim.device();
462 VertexContainer vertices(device.mesh());
467 for (VertexIterator vit = vertices.begin(); vit != vertices.end(); ++vit)
469 she_fill_dos_at_cell(num, sim.config().energy_spacing(), sim.config().dispersion_relation(ctype), *vit, ekin, dos, len);
482 template <
typename SimulatorT >
485 typedef typename SimulatorT::device_type DeviceType;
486 typedef typename DeviceType::mesh_type MeshType;
488 typedef typename viennagrid::result_of::const_vertex_range<MeshType>::type VertexContainer;
489 typedef typename viennagrid::result_of::iterator<VertexContainer>::type VertexIterator;
491 DeviceType
const & device = sim.device();
492 VertexContainer vertices(device.mesh());
498 for (VertexIterator vit = vertices.begin(); vit != vertices.end(); ++vit)
500 const std::size_t idx = std::size_t(vit->id().get());
502 for (std::size_t index_H = 0; index_H < num; ++index_H)
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);
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.
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)
Returns the lattice temperature on the device.
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 ®, 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 ®, 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 ®, DeviceType const &device, SHEQuanT const &shequan, std::string name)
Registers the average expansion order (SHE simulator)
void register_electric_field(quan_register_internal ®, 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 ®, 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 ®, 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 ®, 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 ®)
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.
logger< true > warn()
Used to log warnings. The logging level is logWARNING.
viennashe::models::dd::mobility< DeviceType > create_constant_mobility_model(DeviceType const &device, double mu)
Returns a mobility model, which always yields the same mobility.
double eV_to_joule(const double eps)
std::string lattice_temperature()
std::string hole_density()
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.
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
Contains all viennashe includes and most of the C++/C-wrappers. Contains macros.