Namespaces | |
namespace | detail |
Classes | |
class | acoustic_phonon_scattering |
Acoustic phonon scattering process. More... | |
class | acoustic_phonon_scattering_parameters |
Parameters for elastic acoustic phonon scattering. More... | |
class | adaptive_she_not_available_for_this_configuration_exception |
Exception for the case that adaptive SHE cannot deal with the provided configuration (currently holes or ee-scattering) More... | |
class | assembly_exception |
Exception for the case that a scattering term is Inf, NaN or causes an invalid system matrix entry. More... | |
class | average_expansion_order_wrapper |
Accessor class providing the average expansion order inside the device. More... | |
class | carrier_density_wrapper |
An accessor for the carrier density in the device. More... | |
class | carrier_energy_wrapper |
An accessor for the average carrier energy at each point inside the device. More... | |
class | carrier_velocity_wrapper |
Accessor class providing the carrier velocity inside the device. More... | |
class | const_she_quantity |
Common representation of any quantity associated with two objects of a certain type and an index. More... | |
class | coupled_vertices_equal_exception |
Exception for the case that a vertex has a coupling with itself. More... | |
class | current_density_wrapper |
Accessor class providing the carrier velocity inside the device. More... | |
class | df_wrapper |
A convenience wrapper for evaluating the full distribution function. Note: Evaluations are quite costly. If you wish to evaluate integrals over the distribution function, you should definitely consider the use of SHE expansion coefficients as well as the orthogonality relations for spherical harmonics. More... | |
class | division_by_zero |
Exception for the case that any component encounters a divison by 0.0. More... | |
class | edf_wrapper |
A convenience wrapper for accessing the energy distribution function on each vertex of the device. More... | |
class | fixed_charge_scattering |
Trapped charge scattering process. More... | |
class | fixed_charge_scattering_parameters |
class | force_on_facet_accessor |
class | GammaIntegrand_x |
x-component of the Gamma coupling integral (cf. Dissertation Rupp) More... | |
class | GammaIntegrand_y |
y-component of the Gamma coupling integral (cf. Dissertation Rupp) More... | |
class | GammaIntegrand_z |
z-component of the Gamma coupling integral (cf. Dissertation Rupp) More... | |
class | generalized_df_wrapper |
A convenience wrapper for evaluating the generalized distribution function (i.e. f * Z, with density of states Z). Note: Evaluations are quite costly. If you wish to evaluate integrals over the distribution function, you should definitely consider the use of SHE expansion coefficients as well as the orthogonality relations for spherical harmonics. More... | |
class | generalized_edf_wrapper |
A convenience wrapper for accessing the generalized energy distribution function (f * Z, with density of states Z) on each vertex of the device. More... | |
class | impact_ionization_scattering |
Impact ionization scattering process. More... | |
class | impact_ionization_scattering_parameters |
class | interpolated_she_df_wrapper |
A convenience wrapper for accessing the distribution function coefficients over the device on each vertex. Note that odd-order expansion coefficients are averaged, thus one MUST NOT evaluate flux-quantities using this interpolated wrapper! More... | |
class | invalid_expansion_order_exception |
Exception for the case that invalid expansion order is accessed. More... | |
class | invalid_matrixelement_exception |
Exception for the case that invalid expansion order is accessed. More... | |
class | invalid_scattering_term_exception |
Exception for the case that a scattering term is Inf, NaN or causes an invalid system matrix entry. More... | |
class | ionized_impurity_scattering |
Ionized impurity scattering process. More... | |
class | ionized_impurity_scattering_parameters |
Parameters for ionized impurity scattering using an isotropic fit, cf. papers by Jungemann. More... | |
struct | log_acoustic_phonon_scattering |
Configuration class for logging when assembling acoustic phonon scattering. More... | |
struct | log_adaptive_she |
Configuration class for logging when using adaptive SHE. More... | |
struct | log_assemble_all |
Configuration class for logging when assembling the free streaming operator. More... | |
struct | log_assemble_ee_scattering |
Configuration class for logging when assembling electron-electron scattering. More... | |
struct | log_assemble_free_streaming_operator |
Configuration class for logging when assembling the free streaming operator. More... | |
struct | log_assemble_scattering_operator |
Configuration class for logging when assembling any scattering mechanisms. More... | |
struct | log_coupling_matrix_in_direction |
Configuration class for logging when computing the coupling matrix in a particular direction. More... | |
struct | log_fill_coupling_matrices |
Configuration class for logging when filling coupling matrices. More... | |
struct | log_impact_ionization_scattering |
Configuration class for logging when assembling impact ionization scattering. More... | |
struct | log_linear_solver |
Configuration class for logging when running the linear solver. More... | |
struct | log_mapping |
Configuration class for logging when distributing unknown indices over device (mapping) More... | |
struct | log_newton_she |
struct | log_optical_phonon_scattering |
Configuration class for logging when assembling optical phonons scattering. More... | |
struct | log_recover_odd_unknowns |
Configuration class for logging when recovering odd unknowns from the computed even f_{l,m}. More... | |
struct | log_she_solver |
Configuration class for logging during the linear solver phase in a SHE simulation. More... | |
struct | log_smooth_expansion_order |
Configuration class for logging when smoothing expansion orders for mapping. More... | |
struct | log_transfer_to_new_h_space |
Configuration class for logging when dealing with traps. More... | |
struct | log_traps |
Configuration class for logging when dealing with traps. More... | |
class | negative_integration_interval_length_exception |
Exception for the case that neither electrons nor holes are selected for the simulation. More... | |
class | no_carrier_type_id_specified_exception |
Exception for the case that neither electrons nor holes are selected for the simulation. More... | |
class | no_init_guess_found_exception |
Exception for the case that no initial guess was/is specified. More... | |
class | optical_phonon_scattering |
Optical phonon scattering process. More... | |
class | optical_phonon_scattering_parameters |
Parameters for inelastic optical phonon scattering in single valley approximation. More... | |
class | quantity_not_yet_available_exception |
Exception for the case that a macroscopic quantity is accessed, but the simulator has not yet been run. More... | |
class | quantum_correction_with_newton_not_supported_exception |
Exception for the case that a scattering term is Inf, NaN or causes an invalid system matrix entry. More... | |
class | scatter_config |
A configuration class for scattering mechanisms. Enable or disable scattering mechanisms here. More... | |
class | scatter_process_descriptor |
A simple class returning the scattering rate and the energy of a scattered particle. More... | |
class | scattering_base |
class | scattering_parameter_base |
Common base class for all scattering parameter classes. Provides enable/disable interface. More... | |
class | she_df_wrapper |
A convenience wrapper for accessing the distribution function coefficients over the device. Even-order coefficients can be obtained from vertices, odd-order coefficients from edges. More... | |
class | she_simulator_does_not_accept_drift_diffusion_only_exception |
Exception for the case that neither electrons nor holes are selected for the simulation. More... | |
class | she_simulator_requires_bipolar_solution_for_traps |
Exception for the case that traps are enabled without a bipolar SHE simulation. More... | |
class | surface_scattering |
Surface scattering process. More... | |
class | surface_scattering_parameters |
class | timestep_quantities |
The main SHE simulator controller class. Acts as an accessor for all SHE quantities needed for the simulation. More... | |
class | total_energy_too_small_exception |
Exception for the case that the total energy is smaller than the kinetic energy. More... | |
class | trapped_charge_scattering |
Trapped charge scattering process. More... | |
class | trapped_charge_scattering_parameters |
class | unknown_dispersion_relation_exception |
Exception for the case that an invalid dispersion relation is specified. More... | |
class | unknown_she_quantity |
General representation of any solver quantity defined on two different element types (e.g. vertices and edges) in an augmented (x, H) space. More... | |
class | unkown_carrier_type_exception |
Exception thrown in case an unkown or unsupported carrier type is found. More... | |
class | vIntegrand_x |
x-component of the velocity coupling integral \int Y_{l,m} e_r Y_{l',m'} \d \Omega (cf. Dissertation Rupp) More... | |
class | vIntegrand_y |
y-component of the velocity coupling integral \int Y_{l,m} e_r Y_{l',m'} \d \Omega (cf. Dissertation Rupp) More... | |
class | vIntegrand_z |
z-component of the velocity coupling integral \int Y_{l,m} e_r Y_{l',m'} \d \Omega (cf. Dissertation Rupp) More... | |
Functions | |
template<typename DeviceType , typename TimeStepQuantitiesT , typename VertexT , typename EdgeT , typename MatrixType , typename VectorType > | |
void | assemble (DeviceType &device, TimeStepQuantitiesT &old_quantities, TimeStepQuantitiesT &quantities, viennashe::config const &conf, viennashe::she::unknown_she_quantity< VertexT, EdgeT > const &quan, MatrixType &A, VectorType &b, bool use_timedependence, bool quan_valid) |
template<typename DeviceType , typename SHEQuantity , typename MatrixType , typename VectorType , typename CellType , typename CouplingMatrixType > | |
void | assemble_boundary_on_box (DeviceType const &device, viennashe::config const &conf, SHEQuantity const &quan, MatrixType &A, VectorType &b, CellType const &cell, std::size_t index_H, CouplingMatrixType const &coupling_identity) |
Worker function for the assembly of the free streaming operator. Handles the assembly of both even and odd unknowns. More... | |
template<typename DispersionRelation > | |
double | integral_Z_over_hk (double energy_minus, double energy_plus, DispersionRelation const &dispersion) |
Computes \int_a^b Z/(hbar * k) dE, where a and b are the energy boundaries, Z is the density of states and hbar * k is the momentum. More... | |
template<typename MatrixType , typename CellType , typename PolarityTag > | |
MatrixType | coupling_matrix_in_direction_impl (MatrixType const &m1, MatrixType const &m2, MatrixType const &m3, CellType const &c1, CellType const &c2, PolarityTag const &polarity, viennagrid::cartesian_cs< 1 >) |
Returns dot(M, n), where M=(m1, m2, m3) is the vector of coupling matrices, and n = v2 - v1 is the directional vector (i.e. normal vector on box). One spatial dimension. More... | |
template<typename MatrixType , typename CellType , typename PolarityTag > | |
MatrixType | coupling_matrix_in_direction_impl (MatrixType const &m1, MatrixType const &m2, MatrixType const &m3, CellType const &c1, CellType const &c2, PolarityTag const &polarity, viennagrid::cartesian_cs< 2 >) |
Returns dot(M, n), where M=(m1, m2, m3) is the vector of coupling matrices, and n = v2 - v1 is the directional vector (i.e. normal vector on box). Two spatial dimensions. More... | |
template<typename MatrixType , typename CellType , typename PolarityTag > | |
MatrixType | coupling_matrix_in_direction_impl (MatrixType const &m1, MatrixType const &m2, MatrixType const &m3, CellType const &c1, CellType const &c2, PolarityTag const &polarity, viennagrid::cartesian_cs< 3 >) |
Returns dot(M, n), where M=(m1, m2, m3) is the vector of coupling matrices, and n = v2 - v1 is the directional vector (i.e. normal vector on box). Three spatial dimensions. More... | |
template<typename MatrixType , typename VertexType , typename PolarityTag > | |
MatrixType | coupling_matrix_in_direction (MatrixType const &m1, MatrixType const &m2, MatrixType const &m3, VertexType const &v1, VertexType const &v2, PolarityTag const &polarity) |
Returns dot(M, n), where M=(m1, m2, m3) is the vector of coupling matrices, and n = v2 - v1 is the directional vector (i.e. normal vector on box) More... | |
template<typename MatrixType , typename VectorType , typename IteratorType > | |
void | write_boundary (VectorType &rhs, std::size_t row_start, double prefactor, MatrixType const &coupling_matrix, IteratorType row_iter) |
Write Dirichlet boundary conditions to right hand side. More... | |
template<typename DeviceType , typename SHEQuantity , typename FacetType > | |
double | lower_kinetic_energy_at_facet (DeviceType const &device, SHEQuantity const &quan, FacetType const &facet, std::size_t index_H) |
Returns the lower kinetic energy for the discretization box associated with the edge 'edge'. More... | |
template<typename DeviceType , typename SHEQuantity , typename FacetType > | |
double | upper_kinetic_energy_at_facet (DeviceType const &device, SHEQuantity const &quan, FacetType const &facet, std::size_t index_H) |
Returns the upper kinetic energy of the discretization box for the edge 'edge'. More... | |
template<typename SHEQuantity , typename ElementType > | |
double | box_height (SHEQuantity const &quan, ElementType const &elem, std::size_t index_H) |
Returns the height of the control box with respect to energy at the provided element (vertex or edge) and energy. More... | |
template<typename SHEQuantity , typename CellFacetType > | |
double | averaged_density_of_states (SHEQuantity const &quan, viennashe::config::dispersion_relation_type const &dispersion, CellFacetType const &cell_facet, std::size_t index_H) |
Returns the density of states around a vertex or an edge at total energy specified by index_H. Some averaging is applied near the band edge. More... | |
template<typename SHEQuantity , typename CellFacetType > | |
double | integral_v (SHEQuantity const &quan, viennashe::config::dispersion_relation_type const &dispersion, CellFacetType const &cell_facet, std::size_t index_H) |
Computes \int v dH using a Simpson rule, where v is velocity. More... | |
template<typename SHEQuantity , typename CellFacetType > | |
double | integral_vZ (SHEQuantity const &quan, viennashe::config::dispersion_relation_type const &dispersion, CellFacetType const &cell_facet, std::size_t index_H) |
Computes \int v Z dH using a Simpson rule, where v is velocity and Z is the density of states. More... | |
template<typename SHEQuantity , typename CellType > | |
double | averaged_kinetic_energy (SHEQuantity const &quan, CellType const &c, std::size_t index_H) |
Returns the averaged kinetic energy in the box centered at a vertex v with total energy index index_H. More... | |
template<typename SHEQuantity , typename ElementType > | |
long | energy_index_H (SHEQuantity const &quan, ElementType const &elem, long index_H_hint, double kinetic_energy) |
Finds the closest H-index for the provided kinetic energy. More... | |
template<typename ScatterProcessesT , typename DeviceType , typename SHEQuantity , typename MatrixType , typename VectorType , typename ElementType , typename CouplingMatrix > | |
void | assemble_scattering_operator_on_box (ScatterProcessesT const &scatter_processes, DeviceType const &device, viennashe::config const &conf, SHEQuantity const &quan, MatrixType &A, VectorType &b, ElementType const &elem, std::size_t index_H, CouplingMatrix const &coupling_in_scatter, CouplingMatrix const &coupling_out_scatter) |
template<typename DeviceType , typename SHEQuantity , typename MatrixType , typename VectorType , typename CellType , typename FacetType , typename CouplingMatrixType > | |
void | assemble_free_streaming_operator_on_box (DeviceType const &device, viennashe::config const &conf, SHEQuantity const &quan, MatrixType &A, VectorType &, CellType const &cell, FacetType const &facet, std::size_t index_H, CouplingMatrixType const &coupling_matrix_diffusion, CouplingMatrixType const &coupling_matrix_drift, bool odd_assembly) |
Worker function for the assembly of the free streaming operator. Handles the assembly of both even and odd unknowns. More... | |
template<typename DeviceType , typename SHEQuantity , typename MatrixType , typename VectorType , typename ElementType , typename CouplingMatrixType , typename QuantityPotential , typename QuantityPotentialOld > | |
void | assemble_timederivative (DeviceType const &device, viennashe::config const &conf, SHEQuantity const &quan, SHEQuantity const &quan_old, MatrixType &A, VectorType &b, ElementType const &elem, std::size_t index_H, CouplingMatrixType const &identity_matrix, QuantityPotential const &quan_pot, QuantityPotentialOld const &quan_pot_old) |
template<typename DeviceType , typename TimeStepQuantitiesT , typename SHEQuantity , typename MatrixType , typename VectorType > | |
void | assemble_traps_coupling (DeviceType const &device, TimeStepQuantitiesT const &quantities, viennashe::config const &conf, SHEQuantity const &quan, MatrixType &matrix, VectorType &rhs) |
Interface function for the assembly of traps. More... | |
template<typename DeviceType , typename TimeStepQuantitiesT , typename SHEQuantity , typename MatrixType , typename VectorType > | |
void | assemble_traps_solver (DeviceType const &device, TimeStepQuantitiesT &quantities, viennashe::config const &conf, SHEQuantity const &quan, MatrixType &matrix, VectorType &rhs) |
template<typename DeviceType , typename TimeStepQuantitiesT , typename SHEQuantity , typename MatrixType , typename VectorType > | |
void | assemble_traps (DeviceType const &device, TimeStepQuantitiesT &quantities, viennashe::config const &conf, SHEQuantity const &quan, MatrixType &matrix, VectorType &rhs) |
template<typename DeviceType , typename VertexT , typename EdgeT > | |
void | write_boundary_conditions (DeviceType const &device, viennashe::she::unknown_she_quantity< VertexT, EdgeT > &quan, viennashe::config const &conf) |
Writes boundary conditions for SHE to the device. Stores the result using ViennaData. More... | |
template<typename DeviceType > | |
void | write_boundary_conditions (DeviceType const &device, timestep_quantities< DeviceType > &quantities, viennashe::config const &conf) |
template<typename FullMatrixType , typename VectorType > | |
void | diagonalise_odd2odd_coupling_matrix (FullMatrixType &system_matrix, VectorType &rhs, std::size_t num_even) |
Eliminates all odd spherical harmonics expansion coefficients in the off diagonals of S^oo from the system matrix by line operations. More... | |
template<typename FullMatrixType , typename CompressedMatrixType , typename VectorType > | |
void | eliminate_odd_unknowns (FullMatrixType &system_matrix, VectorType const &rhs, CompressedMatrixType &compressed_matrix, VectorType &compressed_rhs) |
Eliminates all odd spherical harmonics expansion coefficients from the system matrix. More... | |
template<typename MatrixT , typename VectorT > | |
VectorT | recover_odd_unknowns (MatrixT const &full_matrix, VectorT const &full_rhs, VectorT const &compressed_result) |
Recovers the odd-order unknowns from the even-order unknowns. More... | |
template<typename DeviceType , typename SHEQuantity > | |
void | lower_expansion_order_near_bandedge (DeviceType const &device, SHEQuantity &quan) |
Assigns the various expansion orders to the device. Uniform expansion order implemented here. More... | |
template<typename DeviceType , typename SHEQuantity > | |
void | smooth_expansion_order (DeviceType const &device, SHEQuantity &quan) |
Smoothes the expansion order over the device such that neighboring vertices differ in their expansion order by less than two. More... | |
template<typename DeviceType , typename SHEQuantity > | |
void | write_expansion_order_to_facets (DeviceType &device, SHEQuantity &quan) |
Writes expansion orders on edges. Edge expansion order is given as max(order(v1), order(v2)), where v1 and v2 denote the cells of the facet. More... | |
template<typename DeviceType , typename VertexT , typename EdgeT > | |
void | distribute_expansion_orders (DeviceType &device, viennashe::she::unknown_she_quantity< VertexT, EdgeT > &quan, viennashe::config const &conf) |
Assigns the various expansion orders to the device. This is the public interface. More... | |
template<typename DeviceType > | |
void | distribute_expansion_orders (DeviceType const &device, timestep_quantities< DeviceType > &quantities, viennashe::config const &conf) |
template<typename IntegrationRule , typename MatrixType > | |
void | fill_coupling_matrices_impl (MatrixType &a_x, MatrixType &a_y, MatrixType &a_z, MatrixType &b_x, MatrixType &b_y, MatrixType &b_z, int L_max) |
Assemble coupling coefficients a_{l,m}^{l',m'} and b_{l,m}^{l',m'} up to order l = L_max. More... | |
template<typename MatrixType > | |
void | fill_coupling_matrices_1 (MatrixType &a_x, MatrixType &a_y, MatrixType &a_z, MatrixType &b_x, MatrixType &b_y, MatrixType &b_z) |
Precomputed coupling matrix up to first order. More... | |
template<typename MatrixType > | |
void | fill_coupling_matrices_3 (MatrixType &a_x, MatrixType &a_y, MatrixType &a_z, MatrixType &b_x, MatrixType &b_y, MatrixType &b_z) |
Precomputed coupling matrix up to third order. More... | |
template<typename MatrixType > | |
void | fill_coupling_matrices_5 (MatrixType &a_x, MatrixType &a_y, MatrixType &a_z, MatrixType &b_x, MatrixType &b_y, MatrixType &b_z) |
Precomputed coupling matrix up to fifth order. More... | |
template<typename MatrixType > | |
void | fill_coupling_matrices_7 (MatrixType &a_x, MatrixType &a_y, MatrixType &a_z, MatrixType &b_x, MatrixType &b_y, MatrixType &b_z) |
Precomputed coupling matrix up to seventh order. More... | |
template<typename MatrixType > | |
void | fill_coupling_matrices (MatrixType &a_x, MatrixType &a_y, MatrixType &a_z, MatrixType &b_x, MatrixType &b_y, MatrixType &b_z, int L_max) |
Public interface for filling coupling matrices. More... | |
template<typename DeviceType , typename SHEQuantity > | |
void | fill_block_indices (DeviceType const &device, SHEQuantity const &quan, std::size_t system_size, std::vector< std::pair< std::size_t, std::size_t > > &indices) |
Writes an array of block indices. More... | |
template<typename DeviceType , typename SHEQuantity , typename MatrixType , typename VectorType > | |
VectorType | solve (DeviceType const &device, SHEQuantity const &quan, viennashe::config const &configuration, MatrixType &full_matrix, VectorType &full_rhs, std::size_t reduced_unknowns) |
Public interface for solving the provided system of discretized SHE equations. More... | |
template<typename DeviceType , typename SHEQuantity , typename VectorType > | |
void | compute_carrier_density_vector (DeviceType const &device, SHEQuantity const &quan, viennashe::config::dispersion_relation_type const &dispersion, VectorType &carrier) |
Computes the carrier density in the device and writes the result to a vector. Used during Gummel iteration. More... | |
template<typename DeviceType , typename SHEQuantity , typename ContainerType > | |
void | write_carrier_density_to_container (DeviceType const &device, SHEQuantity const &quan, viennashe::config::dispersion_relation_type const &dispersion, ContainerType &container, double energy_start=0.0, double energy_end=1.0) |
Convenience function for writing the average expansion order to the container provided. More... | |
template<typename DeviceType , typename SHEQuantity , typename ContainerType > | |
void | write_kinetic_carrier_energy_to_container (DeviceType const &device, viennashe::config const &conf, SHEQuantity const &quan, ContainerType &container) |
Convenience function for writing the average kinetic carrier energy to the container provided. More... | |
template<typename DeviceType , typename SHEQuantity , typename ContainerType > | |
void | write_carrier_velocity_to_container (DeviceType const &device, viennashe::config const &conf, SHEQuantity const &quan, ContainerType &container) |
Convenience function for writing the average expansion order to the container provided. More... | |
template<typename DeviceType , typename SHEQuantity , typename ContainerType > | |
void | write_current_density_to_container (DeviceType const &device, viennashe::config const &conf, SHEQuantity const &quan, ContainerType &container) |
Convenience function for writing the current density to the container provided. More... | |
template<typename DeviceT , typename SHEQuantity > | |
void | check_current_conservation (DeviceT const &device, viennashe::config const &conf, SHEQuantity const &quan) |
Checks current conservation for SHE. Writes information using log::info(). More... | |
template<typename DeviceType , typename SHEQuantity , typename ContainerType > | |
void | write_average_expansion_order_to_container (DeviceType const &device, SHEQuantity const &quan, ContainerType &container, double energy_start=0.0, double energy_end=1.0) |
Convenience function for writing the average expansion order to the container provided. More... | |
template<typename DeviceType , typename SHEQuantity , typename VectorType > | |
void | setup_unknown_scaling (DeviceType const &device, viennashe::config const &conf, SHEQuantity &quan, VectorType &scaling_vector) |
Initializes the unknown scaling. Returns a vector holding the scaling factors. More... | |
template<typename MatrixType , typename VectorType > | |
void | rescale_system (MatrixType &system_matrix, VectorType const &scaling) |
Scales the system matrix (right-preconditioner) with respect to the rescaled unknowns. More... | |
template<typename DispersionRelation > | |
double | ee_scattering_rate (DispersionRelation const &dispersion, double kinetic_energy, double n, double T) |
template<typename DeviceType , typename SHEQuantityT , typename MatrixType , typename VectorType > | |
void | assemble_ee_scattering (DeviceType const &device, viennashe::config const &conf, SHEQuantityT const &quan, SHEQuantityT const &quan_old, MatrixType &matrix, VectorType &rhs) |
Interface function for electron-electron scattering. Differs significantly from ac, op and impurity scattering, thus separate a implementation is used (at least for the moment) More... | |
template<typename DeviceT , typename VertexT , typename EdgeT > | |
void | setup_energies (DeviceT const &device, viennashe::she::unknown_she_quantity< VertexT, EdgeT > &quan, viennashe::config const &conf, viennashe::unknown_quantity< VertexT > const &potential, viennashe::unknown_quantity< VertexT > const &quantum_correction) |
Computes the kinetic energy over the device as obtained from the provided potential. More... | |
template<typename DeviceType > | |
void | setup_energies (DeviceType const &device, timestep_quantities< DeviceType > &quantities, viennashe::config const &conf) |
Computes the kinetic energy over the device as obtained from the provided potential. More... | |
long | even_unknowns_on_node (long L) |
Returns the number of even spherical harmonics up to (including) order L. More... | |
long | odd_unknowns_on_node (long L) |
Returns the number of odd spherical harmonics up to (including) order L. More... | |
template<typename DeviceType > | |
void | transfer_to_new_h_space (DeviceType const &device, viennashe::she::timestep_quantities< DeviceType > const &old_quantities, viennashe::she::timestep_quantities< DeviceType > &new_quantities, viennashe::config const &conf) |
Interface transferring a solution given by 'old_solution' on some other (old) grid to the new grid. More... | |
Definition at line 36 of file common.hpp.
void viennashe::she::assemble | ( | DeviceType & | device, |
TimeStepQuantitiesT & | old_quantities, | ||
TimeStepQuantitiesT & | quantities, | ||
viennashe::config const & | conf, | ||
viennashe::she::unknown_she_quantity< VertexT, EdgeT > const & | quan, | ||
MatrixType & | A, | ||
VectorType & | b, | ||
bool | use_timedependence, | ||
bool | quan_valid | ||
) |
Definition at line 46 of file assemble_all.hpp.
void viennashe::she::assemble_boundary_on_box | ( | DeviceType const & | device, |
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
MatrixType & | A, | ||
VectorType & | b, | ||
CellType const & | cell, | ||
std::size_t | index_H, | ||
CouplingMatrixType const & | coupling_identity | ||
) |
Worker function for the assembly of the free streaming operator. Handles the assembly of both even and odd unknowns.
DeviceType | The device descriptor class |
SHEQuantity | The SHE quantity type giving access to respective data (unknown indices, etc.) |
MatrixType | The system matrix type |
VectorType | The load vector type |
CellType | The topological cell element type for boundary assembly |
CouplingMatrixType | Type of the coupling matrices a_{l,m}^{l',m'} and b_{l,m}^{l',m'} |
Definition at line 69 of file assemble_boundary.hpp.
void viennashe::she::assemble_ee_scattering | ( | DeviceType const & | device, |
viennashe::config const & | conf, | ||
SHEQuantityT const & | quan, | ||
SHEQuantityT const & | quan_old, | ||
MatrixType & | matrix, | ||
VectorType & | rhs | ||
) |
Interface function for electron-electron scattering. Differs significantly from ac, op and impurity scattering, thus separate a implementation is used (at least for the moment)
device | The device description |
conf | The simulator configuration |
matrix | System matrix |
rhs | Load vector |
quan | The current unkown SHE quantity |
quan_old | The last computed SHE quantity (probably from the last non-linear iteration step) |
Definition at line 95 of file assemble_ee_scattering.hpp.
void viennashe::she::assemble_free_streaming_operator_on_box | ( | DeviceType const & | device, |
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
MatrixType & | A, | ||
VectorType & | , | ||
CellType const & | cell, | ||
FacetType const & | facet, | ||
std::size_t | index_H, | ||
CouplingMatrixType const & | coupling_matrix_diffusion, | ||
CouplingMatrixType const & | coupling_matrix_drift, | ||
bool | odd_assembly | ||
) |
Worker function for the assembly of the free streaming operator. Handles the assembly of both even and odd unknowns.
HarmonicsIter1 | Iterator over the indices of the spherical harmonics to be assembled. Even harmonics for even unknowns, odd harmonics for odd unknowns |
HarmonicsIter2 | Iterator over the indices of the coupled spherical harmonics. Opposite parity of HarmonicsIter1 due to even-to-odd and odd-to-even property of the free streaming operator |
DeviceType | The device descriptor class |
SHEController | The SHE controller type giving access to all SHE-specific quantities (unknown indices, etc.) |
MatrixType | The system matrix type |
VectorType | The load vector type |
ElementType1 | The topological element type the assembly is associated with. Vertex for even unknowns, edge for odd unknowns. |
ElementType2 | The topological element type of coupled unknowns. Edge for even unknowns, vertex for odd unknowns. Remember even-to-odd (odd-to-even) property of free streaming operator. |
CouplingMatrixType | Type of the coupling matrices a_{l,m}^{l',m'} and b_{l,m}^{l',m'} |
Definition at line 73 of file assemble_streaming.hpp.
void viennashe::she::assemble_scattering_operator_on_box | ( | ScatterProcessesT const & | scatter_processes, |
DeviceType const & | device, | ||
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
MatrixType & | A, | ||
VectorType & | b, | ||
ElementType const & | elem, | ||
std::size_t | index_H, | ||
CouplingMatrix const & | coupling_in_scatter, | ||
CouplingMatrix const & | coupling_out_scatter | ||
) |
Definition at line 116 of file assemble_scattering.hpp.
void viennashe::she::assemble_timederivative | ( | DeviceType const & | device, |
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
SHEQuantity const & | quan_old, | ||
MatrixType & | A, | ||
VectorType & | b, | ||
ElementType const & | elem, | ||
std::size_t | index_H, | ||
CouplingMatrixType const & | identity_matrix, | ||
QuantityPotential const & | quan_pot, | ||
QuantityPotentialOld const & | quan_pot_old | ||
) |
Definition at line 61 of file assemble_timederivative.hpp.
void viennashe::she::assemble_traps | ( | DeviceType const & | device, |
TimeStepQuantitiesT & | quantities, | ||
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
MatrixType & | matrix, | ||
VectorType & | rhs | ||
) |
Definition at line 401 of file assemble_traps.hpp.
void viennashe::she::assemble_traps_coupling | ( | DeviceType const & | device, |
TimeStepQuantitiesT const & | quantities, | ||
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
MatrixType & | matrix, | ||
VectorType & | rhs | ||
) |
Interface function for the assembly of traps.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quantities | Container of all the SHE quantities |
conf | The simulator configuration |
quan | Quantity describing the trap occupancy |
matrix | System matrix |
rhs | Load vector |
Definition at line 263 of file assemble_traps.hpp.
void viennashe::she::assemble_traps_solver | ( | DeviceType const & | device, |
TimeStepQuantitiesT & | quantities, | ||
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
MatrixType & | matrix, | ||
VectorType & | rhs | ||
) |
Definition at line 330 of file assemble_traps.hpp.
double viennashe::she::averaged_density_of_states | ( | SHEQuantity const & | quan, |
viennashe::config::dispersion_relation_type const & | dispersion, | ||
CellFacetType const & | cell_facet, | ||
std::size_t | index_H | ||
) |
Returns the density of states around a vertex or an edge at total energy specified by index_H. Some averaging is applied near the band edge.
Definition at line 433 of file assemble_common.hpp.
double viennashe::she::averaged_kinetic_energy | ( | SHEQuantity const & | quan, |
CellType const & | c, | ||
std::size_t | index_H | ||
) |
Returns the averaged kinetic energy in the box centered at a vertex v with total energy index index_H.
Definition at line 555 of file assemble_common.hpp.
double viennashe::she::box_height | ( | SHEQuantity const & | quan, |
ElementType const & | elem, | ||
std::size_t | index_H | ||
) |
Returns the height of the control box with respect to energy at the provided element (vertex or edge) and energy.
Definition at line 412 of file assemble_common.hpp.
void viennashe::she::check_current_conservation | ( | DeviceT const & | device, |
viennashe::config const & | conf, | ||
SHEQuantity const & | quan | ||
) |
Checks current conservation for SHE. Writes information using log::info().
device | The device |
conf | The simulator configuration |
quan | The SHE quantities |
Definition at line 253 of file current_density.hpp.
void viennashe::she::compute_carrier_density_vector | ( | DeviceType const & | device, |
SHEQuantity const & | quan, | ||
viennashe::config::dispersion_relation_type const & | dispersion, | ||
VectorType & | carrier | ||
) |
Computes the carrier density in the device and writes the result to a vector. Used during Gummel iteration.
Definition at line 181 of file carrier_density.hpp.
MatrixType viennashe::she::coupling_matrix_in_direction | ( | MatrixType const & | m1, |
MatrixType const & | m2, | ||
MatrixType const & | m3, | ||
VertexType const & | v1, | ||
VertexType const & | v2, | ||
PolarityTag const & | polarity | ||
) |
Returns dot(M, n), where M=(m1, m2, m3) is the vector of coupling matrices, and n = v2 - v1 is the directional vector (i.e. normal vector on box)
m1 | Coupling matrix in x-direction |
m2 | Coupling matrix in y-direction |
m3 | Coupling matrix in z-direction |
v1 | First vertex (of Delaunay triangulation connecting two Voronoi boxes) |
v2 | Second vertex (of Delaunay triangulation connecting two Voronoi boxes) |
polarity | A polarity tag for distinguishing between electrons and holes |
Definition at line 284 of file assemble_common.hpp.
MatrixType viennashe::she::coupling_matrix_in_direction_impl | ( | MatrixType const & | m1, |
MatrixType const & | m2, | ||
MatrixType const & | m3, | ||
CellType const & | c1, | ||
CellType const & | c2, | ||
PolarityTag const & | polarity, | ||
viennagrid::cartesian_cs< 1 > | |||
) |
Returns dot(M, n), where M=(m1, m2, m3) is the vector of coupling matrices, and n = v2 - v1 is the directional vector (i.e. normal vector on box). One spatial dimension.
m1 | Coupling matrix in x-direction |
m2 | Coupling matrix in y-direction |
m3 | Coupling matrix in z-direction |
c1 | First cell (no Delaunay criterion required!) |
c2 | Second cell (no Delaunay criterion required!) |
polarity | A polarity tag for distinguishing between electrons and holes |
Definition at line 177 of file assemble_common.hpp.
MatrixType viennashe::she::coupling_matrix_in_direction_impl | ( | MatrixType const & | m1, |
MatrixType const & | m2, | ||
MatrixType const & | m3, | ||
CellType const & | c1, | ||
CellType const & | c2, | ||
PolarityTag const & | polarity, | ||
viennagrid::cartesian_cs< 2 > | |||
) |
Returns dot(M, n), where M=(m1, m2, m3) is the vector of coupling matrices, and n = v2 - v1 is the directional vector (i.e. normal vector on box). Two spatial dimensions.
m1 | Coupling matrix in x-direction |
m2 | Coupling matrix in y-direction |
m3 | Coupling matrix in z-direction |
c1 | First cell (no Delaunay criterion required!) |
c2 | Second cell (no Delaunay criterion required!) |
polarity | A polarity tag for distinguishing between electrons and holes |
Definition at line 214 of file assemble_common.hpp.
MatrixType viennashe::she::coupling_matrix_in_direction_impl | ( | MatrixType const & | m1, |
MatrixType const & | m2, | ||
MatrixType const & | m3, | ||
CellType const & | c1, | ||
CellType const & | c2, | ||
PolarityTag const & | polarity, | ||
viennagrid::cartesian_cs< 3 > | |||
) |
Returns dot(M, n), where M=(m1, m2, m3) is the vector of coupling matrices, and n = v2 - v1 is the directional vector (i.e. normal vector on box). Three spatial dimensions.
m1 | Coupling matrix in x-direction |
m2 | Coupling matrix in y-direction |
m3 | Coupling matrix in z-direction |
c1 | First cell (no Delaunay criterion required!) |
c2 | Second cell (no Delaunay criterion required!) |
polarity | A polarity tag for distinguishing between electrons and holes |
Definition at line 252 of file assemble_common.hpp.
void viennashe::she::diagonalise_odd2odd_coupling_matrix | ( | FullMatrixType & | system_matrix, |
VectorType & | rhs, | ||
std::size_t | num_even | ||
) |
Eliminates all odd spherical harmonics expansion coefficients in the off diagonals of S^oo from the system matrix by line operations.
FullMatrixType | A the moment this should be the viennashe built-in matrix type (see linalg_util.hpp) |
VectorType | Should be ublas and std::vector compatible |
system_matrix | The full system matrix |
rhs | The full right hand side |
num_even | The number of even unkowns (needed to find S^oo) |
Definition at line 80 of file eliminate.hpp.
void viennashe::she::distribute_expansion_orders | ( | DeviceType & | device, |
viennashe::she::unknown_she_quantity< VertexT, EdgeT > & | quan, | ||
viennashe::config const & | conf | ||
) |
Assigns the various expansion orders to the device. This is the public interface.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The SHE unkown quantities on edges and vertices |
conf | The simulator configuration |
Definition at line 349 of file expansion_order.hpp.
void viennashe::she::distribute_expansion_orders | ( | DeviceType const & | device, |
timestep_quantities< DeviceType > & | quantities, | ||
viennashe::config const & | conf | ||
) |
Definition at line 367 of file expansion_order.hpp.
double viennashe::she::ee_scattering_rate | ( | DispersionRelation const & | dispersion, |
double | kinetic_energy, | ||
double | n, | ||
double | T | ||
) |
Definition at line 55 of file assemble_ee_scattering.hpp.
void viennashe::she::eliminate_odd_unknowns | ( | FullMatrixType & | system_matrix, |
VectorType const & | rhs, | ||
CompressedMatrixType & | compressed_matrix, | ||
VectorType & | compressed_rhs | ||
) |
Eliminates all odd spherical harmonics expansion coefficients from the system matrix.
system_matrix | The full system matrix |
rhs | The full right hand side |
compressed_matrix | The new matrix with only the even harmonics unknowns |
compressed_rhs | The new right hand side vector with only the even harmonics unknowns |
Definition at line 226 of file eliminate.hpp.
long viennashe::she::energy_index_H | ( | SHEQuantity const & | quan, |
ElementType const & | elem, | ||
long | index_H_hint, | ||
double | kinetic_energy | ||
) |
Finds the closest H-index for the provided kinetic energy.
quan | The SHE quantities |
elem | The vertex/edge on which to search |
index_H_hint | Index of the energy from which to scatter |
kinetic_energy | Kinetic energy for which index_H should be found |
Definition at line 65 of file assemble_scattering.hpp.
|
inline |
Returns the number of even spherical harmonics up to (including) order L.
Definition at line 35 of file she_quantity.hpp.
void viennashe::she::fill_block_indices | ( | DeviceType const & | device, |
SHEQuantity const & | quan, | ||
std::size_t | system_size, | ||
std::vector< std::pair< std::size_t, std::size_t > > & | indices | ||
) |
Writes an array of block indices.
Definition at line 44 of file linear_solver.hpp.
void viennashe::she::fill_coupling_matrices | ( | MatrixType & | a_x, |
MatrixType & | a_y, | ||
MatrixType & | a_z, | ||
MatrixType & | b_x, | ||
MatrixType & | b_y, | ||
MatrixType & | b_z, | ||
int | L_max | ||
) |
Public interface for filling coupling matrices.
a_x | x-component of the coupling coefficient a_{l,m}^{l',m'} |
a_y | y-component of the coupling coefficient a_{l,m}^{l',m'} |
a_z | z-component of the coupling coefficient a_{l,m}^{l',m'} |
b_x | x-component of the coupling coefficient b_{l,m}^{l',m'} |
b_y | y-component of the coupling coefficient b_{l,m}^{l',m'} |
b_z | z-component of the coupling coefficient b_{l,m}^{l',m'} |
L_max | Maximum spherical harmonics expansion order |
Definition at line 1050 of file harmonics_coupling.hpp.
void viennashe::she::fill_coupling_matrices_1 | ( | MatrixType & | a_x, |
MatrixType & | a_y, | ||
MatrixType & | a_z, | ||
MatrixType & | b_x, | ||
MatrixType & | b_y, | ||
MatrixType & | b_z | ||
) |
Precomputed coupling matrix up to first order.
Definition at line 282 of file harmonics_coupling.hpp.
void viennashe::she::fill_coupling_matrices_3 | ( | MatrixType & | a_x, |
MatrixType & | a_y, | ||
MatrixType & | a_z, | ||
MatrixType & | b_x, | ||
MatrixType & | b_y, | ||
MatrixType & | b_z | ||
) |
Precomputed coupling matrix up to third order.
Definition at line 309 of file harmonics_coupling.hpp.
void viennashe::she::fill_coupling_matrices_5 | ( | MatrixType & | a_x, |
MatrixType & | a_y, | ||
MatrixType & | a_z, | ||
MatrixType & | b_x, | ||
MatrixType & | b_y, | ||
MatrixType & | b_z | ||
) |
Precomputed coupling matrix up to fifth order.
Definition at line 417 of file harmonics_coupling.hpp.
void viennashe::she::fill_coupling_matrices_7 | ( | MatrixType & | a_x, |
MatrixType & | a_y, | ||
MatrixType & | a_z, | ||
MatrixType & | b_x, | ||
MatrixType & | b_y, | ||
MatrixType & | b_z | ||
) |
Precomputed coupling matrix up to seventh order.
Definition at line 644 of file harmonics_coupling.hpp.
void viennashe::she::fill_coupling_matrices_impl | ( | MatrixType & | a_x, |
MatrixType & | a_y, | ||
MatrixType & | a_z, | ||
MatrixType & | b_x, | ||
MatrixType & | b_y, | ||
MatrixType & | b_z, | ||
int | L_max | ||
) |
Assemble coupling coefficients a_{l,m}^{l',m'} and b_{l,m}^{l',m'} up to order l = L_max.
a_x | x-component of the coupling coefficient a_{l,m}^{l',m'} |
a_y | y-component of the coupling coefficient a_{l,m}^{l',m'} |
a_z | z-component of the coupling coefficient a_{l,m}^{l',m'} |
b_x | x-component of the coupling coefficient b_{l,m}^{l',m'} |
b_y | y-component of the coupling coefficient b_{l,m}^{l',m'} |
b_z | z-component of the coupling coefficient b_{l,m}^{l',m'} |
L_max | Maximum spherical harmonics expansion order |
Definition at line 216 of file harmonics_coupling.hpp.
double viennashe::she::integral_v | ( | SHEQuantity const & | quan, |
viennashe::config::dispersion_relation_type const & | dispersion, | ||
CellFacetType const & | cell_facet, | ||
std::size_t | index_H | ||
) |
Computes \int v dH using a Simpson rule, where v is velocity.
Cf. paper by Hong and Jungemann: A fully coupled scheme..., 2009
dispersion | The dispersion relation used for the SHE simulation |
quan | The SHE quantities |
cell_facet | A topological element which is either a cell or facet |
index_H | The index at which to evaluate the integral |
Definition at line 477 of file assemble_common.hpp.
double viennashe::she::integral_vZ | ( | SHEQuantity const & | quan, |
viennashe::config::dispersion_relation_type const & | dispersion, | ||
CellFacetType const & | cell_facet, | ||
std::size_t | index_H | ||
) |
Computes \int v Z dH using a Simpson rule, where v is velocity and Z is the density of states.
Cf. paper by Hong and Jungemann: A fully coupled scheme..., 2009
dispersion | The dispersion relation used for the SHE simulation |
quan | The SHE quantities |
cell_facet | A topological element which is either a cell or a facet |
index_H | The index at which to evaluate the integral |
Definition at line 521 of file assemble_common.hpp.
double viennashe::she::integral_Z_over_hk | ( | double | energy_minus, |
double | energy_plus, | ||
DispersionRelation const & | dispersion | ||
) |
Computes \int_a^b Z/(hbar * k) dE, where a and b are the energy boundaries, Z is the density of states and hbar * k is the momentum.
Uses the formula Z/(hbar * k) = 0.5 d(Zv)/dE, thus Z and v are evaluated at the energy boundaries a and b only
energy_minus | Lower integration bound (might be smaller than zero -> set to zero) |
energy_plus | Upper integration bound (might be smaller than zero -> set to zero) |
dispersion | The dispersion relation used for the SHE simulation |
Definition at line 153 of file assemble_common.hpp.
void viennashe::she::lower_expansion_order_near_bandedge | ( | DeviceType const & | device, |
SHEQuantity & | quan | ||
) |
Assigns the various expansion orders to the device. Uniform expansion order implemented here.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The SHE quantities |
Definition at line 60 of file expansion_order.hpp.
double viennashe::she::lower_kinetic_energy_at_facet | ( | DeviceType const & | device, |
SHEQuantity const & | quan, | ||
FacetType const & | facet, | ||
std::size_t | index_H | ||
) |
Returns the lower kinetic energy for the discretization box associated with the edge 'edge'.
Definition at line 325 of file assemble_common.hpp.
|
inline |
Returns the number of odd spherical harmonics up to (including) order L.
Definition at line 45 of file she_quantity.hpp.
VectorT viennashe::she::recover_odd_unknowns | ( | MatrixT const & | full_matrix, |
VectorT const & | full_rhs, | ||
VectorT const & | compressed_result | ||
) |
Recovers the odd-order unknowns from the even-order unknowns.
For a system matrix
S = (S^ee, S^eo ) (S^oe, S^oo ),
the odd-order unknowns f^o are recovered from the even unknowns f^e by f^o = (S^oo)^-1 * S^oe * f^e. Note that S^oo is diagonal.
Definition at line 305 of file eliminate.hpp.
void viennashe::she::rescale_system | ( | MatrixType & | system_matrix, |
VectorType const & | scaling | ||
) |
Scales the system matrix (right-preconditioner) with respect to the rescaled unknowns.
system_matrix | The system matrix to be scaled with a diagonal right-preconditioner |
scaling | The values of the diagonal of the right-preconditioner |
Definition at line 91 of file rescaling.hpp.
void viennashe::she::setup_energies | ( | DeviceT const & | device, |
viennashe::she::unknown_she_quantity< VertexT, EdgeT > & | quan, | ||
viennashe::config const & | conf, | ||
viennashe::unknown_quantity< VertexT > const & | potential, | ||
viennashe::unknown_quantity< VertexT > const & | quantum_correction | ||
) |
Computes the kinetic energy over the device as obtained from the provided potential.
device | The device on which simulation is carried out |
quan | The SHE quantities |
conf | The simulator configuration |
potential | The potential as computed with Poisson equation |
quantum_correction | The quantum correction potential |
Definition at line 61 of file setup_energies.hpp.
void viennashe::she::setup_energies | ( | DeviceType const & | device, |
timestep_quantities< DeviceType > & | quantities, | ||
viennashe::config const & | conf | ||
) |
Computes the kinetic energy over the device as obtained from the provided potential.
device | The device on which simulation is carried out |
quantities | The quantities defined on the device |
conf | The simulator configuration |
Definition at line 162 of file setup_energies.hpp.
void viennashe::she::setup_unknown_scaling | ( | DeviceType const & | device, |
viennashe::config const & | conf, | ||
SHEQuantity & | quan, | ||
VectorType & | scaling_vector | ||
) |
Initializes the unknown scaling. Returns a vector holding the scaling factors.
device | The device (includes a ViennaGrid mesh) on which a simulation is carried out |
conf | The simulator config |
quan | The SHE quantities |
scaling_vector | The vector to be filled with the scaling factors |
Definition at line 50 of file rescaling.hpp.
void viennashe::she::smooth_expansion_order | ( | DeviceType const & | device, |
SHEQuantity & | quan | ||
) |
Smoothes the expansion order over the device such that neighboring vertices differ in their expansion order by less than two.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The SHE quantities |
Definition at line 140 of file expansion_order.hpp.
VectorType viennashe::she::solve | ( | DeviceType const & | device, |
SHEQuantity const & | quan, | ||
viennashe::config const & | configuration, | ||
MatrixType & | full_matrix, | ||
VectorType & | full_rhs, | ||
std::size_t | reduced_unknowns | ||
) |
Public interface for solving the provided system of discretized SHE equations.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The SHE quantities |
configuration | The simulator configuration |
full_matrix | The system matrix containing even and odd unknowns |
full_rhs | The load vector containing even and odd unknowns |
reduced_unknowns | The number of reduced unknowns |
Definition at line 100 of file linear_solver.hpp.
void viennashe::she::transfer_to_new_h_space | ( | DeviceType const & | device, |
viennashe::she::timestep_quantities< DeviceType > const & | old_quantities, | ||
viennashe::she::timestep_quantities< DeviceType > & | new_quantities, | ||
viennashe::config const & | conf | ||
) |
Interface transferring a solution given by 'old_solution' on some other (old) grid to the new grid.
device | The device on which simulation is carried out |
conf | The simulator configuration |
new_quantities | The timestep_quantities used for the upcoming simulation, which is filled with the interpolated old values |
old_quantities | The timestep_quantities used for obtaining the old values from another grid (same grid in (x, H)-space, but other potential profile) |
Definition at line 122 of file transfer_to_new_h_space.hpp.
double viennashe::she::upper_kinetic_energy_at_facet | ( | DeviceType const & | device, |
SHEQuantity const & | quan, | ||
FacetType const & | facet, | ||
std::size_t | index_H | ||
) |
Returns the upper kinetic energy of the discretization box for the edge 'edge'.
Definition at line 359 of file assemble_common.hpp.
void viennashe::she::write_average_expansion_order_to_container | ( | DeviceType const & | device, |
SHEQuantity const & | quan, | ||
ContainerType & | container, | ||
double | energy_start = 0.0 , |
||
double | energy_end = 1.0 |
||
) |
Convenience function for writing the average expansion order to the container provided.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The SHE quantities. Passed for compatibility with other write_XYZ_to_domain() routines. |
container | The container to write to |
energy_start | (Optional) Lower bound for the kinetic energy (in eV) range to be considered for the averaging the SHE order |
energy_end | (Optional) Upper bound for the kinetic energy (in eV) range to be considered for the averaging the SHE order |
Definition at line 90 of file expansion_order.hpp.
void viennashe::she::write_boundary | ( | VectorType & | rhs, |
std::size_t | row_start, | ||
double | prefactor, | ||
MatrixType const & | coupling_matrix, | ||
IteratorType | row_iter | ||
) |
Write Dirichlet boundary conditions to right hand side.
rhs | Load vector |
row_start | Load vector index to start |
prefactor | Scalar prefactor for 'coupling_matrix' |
coupling_matrix | SHE coupling matrix |
row_iter | Iterator over row indices in 'coupling_matrix' |
Definition at line 302 of file assemble_common.hpp.
void viennashe::she::write_boundary_conditions | ( | DeviceType const & | device, |
timestep_quantities< DeviceType > & | quantities, | ||
viennashe::config const & | conf | ||
) |
Definition at line 138 of file boundary_conditions.hpp.
void viennashe::she::write_boundary_conditions | ( | DeviceType const & | device, |
viennashe::she::unknown_she_quantity< VertexT, EdgeT > & | quan, | ||
viennashe::config const & | conf | ||
) |
Writes boundary conditions for SHE to the device. Stores the result using ViennaData.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The unknown SHE quantities on edges and vertices |
conf | The simulator configuration |
Definition at line 57 of file boundary_conditions.hpp.
void viennashe::she::write_carrier_density_to_container | ( | DeviceType const & | device, |
SHEQuantity const & | quan, | ||
viennashe::config::dispersion_relation_type const & | dispersion, | ||
ContainerType & | container, | ||
double | energy_start = 0.0 , |
||
double | energy_end = 1.0 |
||
) |
Convenience function for writing the average expansion order to the container provided.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The SHE quantities. Passed for compatibility with other write_XYZ_to_domain() routines. |
dispersion | The dispersion relation |
container | Container for the quantity |
energy_start | (Optional) Lower bound for the kinetic energy (in eV) range to be considered for the averaging the SHE order |
energy_end | (Optional) Upper bound for the kinetic energy (in eV) range to be considered for the averaging the SHE order |
Definition at line 221 of file carrier_density.hpp.
void viennashe::she::write_carrier_velocity_to_container | ( | DeviceType const & | device, |
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
ContainerType & | container | ||
) |
Convenience function for writing the average expansion order to the container provided.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
conf | The simulator configuration |
quan | The SHE quantity (electron or hole distribution function) used for the calculation |
container | The container the values should be written to |
Definition at line 105 of file carrier_velocity.hpp.
void viennashe::she::write_current_density_to_container | ( | DeviceType const & | device, |
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
ContainerType & | container | ||
) |
Convenience function for writing the current density to the container provided.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The SHE quantities. Passed for compatibility with other write_XYZ_to_domain() routines. |
conf | The simulator configuration |
container | Container for the current density vector |
Definition at line 235 of file current_density.hpp.
void viennashe::she::write_expansion_order_to_facets | ( | DeviceType & | device, |
SHEQuantity & | quan | ||
) |
Writes expansion orders on edges. Edge expansion order is given as max(order(v1), order(v2)), where v1 and v2 denote the cells of the facet.
device | The device on which simulation is carried out |
quan | The SHE quantities |
Definition at line 208 of file expansion_order.hpp.
void viennashe::she::write_kinetic_carrier_energy_to_container | ( | DeviceType const & | device, |
viennashe::config const & | conf, | ||
SHEQuantity const & | quan, | ||
ContainerType & | container | ||
) |
Convenience function for writing the average kinetic carrier energy to the container provided.
device | The device (includes a ViennaGrid mesh) on which simulation is carried out |
quan | The SHE quantities. Passed for compatibility with other write_XYZ_to_domain() routines. |
conf | The simulator configuration |
container | The container for storing the average carrier energy values |
Definition at line 164 of file carrier_energy.hpp.