1#ifndef VIENNASHE_SHE_ASSEMBLE_STREAMING_HPP
2#define VIENNASHE_SHE_ASSEMBLE_STREAMING_HPP
20#include "viennagrid/mesh/mesh.hpp"
21#include "viennagrid/algorithm/norm.hpp"
22#include "viennagrid/algorithm/volume.hpp"
23#include "viennagrid/algorithm/voronoi.hpp"
66 template <
typename DeviceType,
72 typename CouplingMatrixType>
75 SHEQuantity
const & quan,
78 CellType
const & cell, FacetType
const & facet, std::size_t index_H,
79 CouplingMatrixType
const & coupling_matrix_diffusion,
80 CouplingMatrixType
const & coupling_matrix_drift,
83 typedef typename viennagrid::result_of::point<CellType>::type PointType;
90 if (!other_cell_ptr)
return;
92 long row_index = odd_assembly ? quan.get_unknown_index(facet, index_H) : quan.get_unknown_index(cell, index_H);
97 long col_index = odd_assembly ? quan.get_unknown_index(cell, index_H) : quan.get_unknown_index(facet, index_H);
102 PointType cell_connection = viennagrid::centroid(cell) - viennagrid::centroid(*other_cell_ptr);
104 PointType cell_connection_normalized = cell_connection / viennagrid::norm(cell_connection);
106 const double weighted_interface_area = viennagrid::volume(facet) * std::fabs(viennagrid::inner_prod(facet_unit_normal, cell_connection_normalized));
108 long expansion_order_row =
static_cast<long>(odd_assembly ? quan.get_expansion_order(facet, index_H) : quan.get_expansion_order(cell, index_H));
109 long expansion_order_column =
static_cast<long>(odd_assembly ? quan.get_expansion_order(cell, index_H) : quan.get_expansion_order(facet, index_H));
123 double matrix_coeff_diffusion = 0;
127 matrix_coeff_diffusion =
integral_vZ(quan, dispersion, facet, index_H) * weighted_interface_area;
130 matrix_coeff_diffusion =
integral_v(quan, dispersion, facet, index_H) * weighted_interface_area;
132 default:
throw std::runtime_error(
"assemble_free_streaming_operator_on_box(): Unknown SHE discretization type!");
138 log::error() <<
"* assemble_free_streaming_operator(): Invalid column index detected while assembling diffusion operator. " << std::endl;
139 log::error<log_assemble_free_streaming_operator>() <<
" Diagnostics: " << std::endl <<
"index_H: " << index_H << std::endl;
140 log::error<log_assemble_free_streaming_operator>() <<
"cell: " << cell << std::endl;
141 log::error<log_assemble_free_streaming_operator>() <<
"row_index: " << row_index << std::endl;
142 log::error<log_assemble_free_streaming_operator>() <<
"facet: " << facet << std::endl;
143 log::error<log_assemble_free_streaming_operator>() <<
"odd_assembly: " << odd_assembly << std::endl;
144 log::error<log_assemble_free_streaming_operator>() <<
"DOF el2: " << col_index << std::endl;
145 log::error<log_assemble_free_streaming_operator>() <<
"other_cell: " << *other_cell_ptr << std::endl;
149 matrix_coeff_diffusion,
150 coupling_matrix_diffusion,
157 log::debug<log_assemble_free_streaming_operator>() <<
"* assemble_free_streaming_operator_on_box(): diffusion term contribution:" << std::endl;
158 log::debug<log_assemble_free_streaming_operator>() <<
" int_vZ * interface_area" << std::endl;
159 log::debug<log_assemble_free_streaming_operator>() <<
" at (" << row_index <<
", " << col_index <<
")" << std::endl;
160 log::debug<log_assemble_free_streaming_operator>() <<
" matrix_coeff_diffusion = " << matrix_coeff_diffusion << std::endl;
161 log::debug<log_assemble_free_streaming_operator>() <<
"weighted_interface_area = " << weighted_interface_area << std::endl;
162 log::debug<log_assemble_free_streaming_operator>() <<
" coupling_matrix_diffusion: " << coupling_matrix_diffusion << std::endl;
242 upper_kinetic_energy,
247 double matrix_coeff_drift = 0;
251 matrix_coeff_drift = int_Z_over_hk * force * weighted_interface_area * connection_len / 2.0;
254 matrix_coeff_drift = int_Z_over_hk / Z * force * weighted_interface_area * connection_len / 2.0;
256 default:
throw std::runtime_error(
"assemble_free_streaming_operator_on_box(): Unknown SHE discretization type!");
262 log::error() <<
"* assemble_free_streaming_operator(): Invalid column index detected while assembling drift operator. Diagnostics: " << std::endl;
263 log::error<log_assemble_free_streaming_operator>() <<
"index_H: " << index_H << std::endl;
264 log::error<log_assemble_free_streaming_operator>() <<
"cell: " << cell << std::endl;
265 log::error<log_assemble_free_streaming_operator>() <<
"row_index: " << row_index << std::endl;
266 log::error<log_assemble_free_streaming_operator>() <<
"facet: " << facet << std::endl;
267 log::error<log_assemble_free_streaming_operator>() <<
"odd_assembly: " << odd_assembly << std::endl;
268 log::error<log_assemble_free_streaming_operator>() <<
"DOF el2: " << col_index << std::endl;
269 log::error<log_assemble_free_streaming_operator>() <<
"other_cell: " << *other_cell_ptr << std::endl;
274 coupling_matrix_drift,
281 log::debug<log_assemble_free_streaming_operator>() <<
"* assemble_free_streaming_operator_on_box(): drift term contribution:" << std::endl;
282 log::debug<log_assemble_free_streaming_operator>() <<
" int_Z_over_hk / Z * force * interface_area * edge_len / 2.0" << std::endl;
283 log::debug<log_assemble_free_streaming_operator>() <<
" at (" << row_index <<
", " << col_index <<
")" << std::endl;
284 log::debug<log_assemble_free_streaming_operator>() <<
" int_Z_over_hk = " << int_Z_over_hk << std::endl;
285 log::debug<log_assemble_free_streaming_operator>() <<
" force = " << force << std::endl;
286 log::debug<log_assemble_free_streaming_operator>() <<
" weighted_interface_area = " << weighted_interface_area << std::endl;
287 log::debug<log_assemble_free_streaming_operator>() <<
" connection_len = " << connection_len << std::endl;
Common routines for the assembly of SHE equations.
Writes a possibly scaled block matrix to the system matrix. Includes similar functionality for comput...
The main SHE configuration class. To be adjusted by the user for his/her needs.
viennashe::physics::dispersion_proxy dispersion_relation(viennashe::carrier_type_id ctype) const
Returns the dispersion relation for electrons.
she_discretization_type_id she_discretization_type() const
MeshT const & mesh() const
Returns the underlying mesh.
Defines the physical properties of a device, e.g. doping. This is the implementation for 2d and highe...
Iteration over all spherical harmonics indices up to order L, with increasing indices l and m.
A proxy object for a dispersion relation. Does NOT take ownership of the provided pointer!
The SHE configuration class is defined here.
Contains the dispersion relations for different materials and different polarities.
Defines several filter functors for the device. A filter functor returns true if the supplied argumen...
Provides the SHE coupling matrices and computes higher-order coupling matrices if required....
Implementation of numerical integration routines.
A logging facility providing fine-grained control over logging in ViennaSHE.
logger< true > error()
Used to log errors. The logging level is logERROR.
harmonics_iteration_type
An enumeration of spherical harmonics types (all, even, odd) ViennaSHE supports.
@ EVEN_HARMONICS_ITERATION_ID
@ ODD_HARMONICS_ITERATION_ID
VectorType::value_type norm_2(VectorType const &v)
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 state...
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.
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'.
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....
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.
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'.
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 an...
void add_block_matrix(SystemMatrixType &system_matrix, std::size_t row_index, std::size_t col_index, double prefactor, BlockMatrixType const &coupling_matrix, RowIndexIterator row_iter, ColumnIndexIterator const &col_iter_init, bool elastic_scattering_roundoff_error_prevention=false)
Writes a sub-matrix of a block matrix 'coupling_matrix' into the global system matrix,...
viennagrid::result_of::point< viennagrid::element< CellTagT, WrappedConfigT > >::type outer_cell_normal_at_facet(viennagrid::element< CellTagT, WrappedConfigT > const &cell, viennagrid::element< FacetTagT, WrappedConfigT > const &facet)
Returns the unit outer normal of a facet with respect to the provided cell.
CellT const * get_other_cell_of_facet(MeshT const &mesh, FacetT const &facet, CellT const &cell)
Helper function returning a const-pointer to the 'second cell' of a facet, or NULL if there is no sec...
The main ViennaSHE namespace. All functionality resides inside this namespace.
@ SHE_DISCRETIZATION_EVEN_ODD_ORDER_GENERALIZED_DF
@ SHE_DISCRETIZATION_EVEN_ODD_ORDER_DF
Provides a number of fundamental constants. All constants in SI units.
Returns a few helper routines for computing physical quantities. To be replaced in the future.
Provides the exceptions used inside the viennashe::she namespace.
Implementation of spherical harmonics plus helper functions.
Defines the log keys used within the viennashe::she namespace.