|
| FEFaceEvaluation (const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &matrix_free, const bool is_minus_face, const unsigned int dof_no_x, const unsigned int dof_no_v, const unsigned int quad_no_x, const unsigned int quad_no_v) |
|
void | reinit (typename MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::ID face_index) override |
|
void | reinit (typename MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::ID cell_index, unsigned int face) |
|
void | read_dof_values (const dealii::LinearAlgebra::distributed::Vector< Number > &src) |
|
void | read_dof_values (const dealii::LinearAlgebra::distributed::Vector< Number > &src, VectorizedArrayType *data) const |
|
void | read_dof_values_from_buffer (const VectorizedArrayType *src) |
|
void | distribute_to_buffer (VectorizedArrayType *dst) const |
|
void | distribute_local_to_global (dealii::LinearAlgebra::distributed::Vector< Number > &dst) const |
|
unsigned int | n_vectorization_lanes_filled () const |
|
template<SpaceType stype> |
DEAL_II_ALWAYS_INLINE void | submit_value (VectorizedArrayType *__restrict data_ptr, const VectorizedArrayType &value, const unsigned int q, const unsigned int qx, const unsigned int qv) const |
|
template<SpaceType stype> |
DEAL_II_ALWAYS_INLINE void | submit_value (VectorizedArrayType *__restrict data_ptr_1, VectorizedArrayType *__restrict data_ptr_2, const VectorizedArrayType &value, const unsigned int q, const unsigned int q1, const unsigned int q2) const |
|
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim_x, VectorizedArrayType > | get_normal_vector_x (const unsigned int qx) const |
|
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim_v, VectorizedArrayType > | get_normal_vector_v (const unsigned int qv) const |
|
dealii::Point< DIM, VectorizedArrayType > | get_quadrature_point (const unsigned int q) const |
|
template<SpaceType stype> |
dealii::Point< dim, VectorizedArrayType > | get_quadrature_point (const unsigned int qx, const unsigned int qv) const |
|
| FEEvaluationBase (const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &matrix_free) |
|
virtual | ~FEEvaluationBase ()=default |
|
VectorizedArrayType * | get_data_ptr () |
|
const dealii::AlignedVector< Number > * | get_shape_values () const |
|
const dealii::AlignedVector< Number > * | get_shape_gradients () const |
|
|
static const unsigned int | DIM = dim_x + dim_v |
|
static const unsigned int | DIM_ = DIM |
|
static const unsigned int | N_Q_POINTS |
|
static const unsigned int | N_Q_POINTS_1_CELL |
|
static const unsigned int | N_Q_POINTS_2_CELL |
|
static const unsigned int | N_Q_POINTS_1_FACE |
|
static const unsigned int | N_Q_POINTS_2_FACE |
|
static const unsigned int | n_vectors = PARENT::n_vectors |
|
static const unsigned int | n_vectors_v = PARENT::n_vectors_v |
|
static const unsigned int | dim = dim_x + dim_v |
|
static const unsigned int | static_dofs_per_cell = dealii::Utilities::pow(degree + 1, dim) |
|
static const unsigned int | static_dofs_per_cell_x = dealii::Utilities::pow(degree + 1, dim_x) |
|
static const unsigned int | static_dofs_per_cell_v = dealii::Utilities::pow(degree + 1, dim_v) |
|
static const unsigned int | static_dofs_per_face = dealii::Utilities::pow(degree + 1, dim_x + dim_v - 1) |
|
static const int | dim = dim_x + dim_v |
|
static const unsigned int | n_vectors = MF::VectorizedArrayTypeX::size() |
|
static const unsigned int | n_vectors_v = MF::VectorizedArrayTypeV::size() |
|
static const int | static_dofs |
|
template<int dim_x, int dim_v, int degree, int n_points, typename Number, typename VectorizedArrayType>
class hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >
The class that provides all functions necessary to evaluate functions at quadrature points and face integrations in phase space. It delegates the the actual evaluation task to two dealii::FEEvaluation and two dealii::FEFaceEvaluation objects and combines the result on-the-fly.
template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
dealii::Point< DIM, VectorizedArrayType > hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::get_quadrature_point |
( |
const unsigned int | q | ) |
const |
|
inline |
Return position of quadrature point.
TODO: specialize for x-space and v-space.
template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
template<SpaceType stype>
dealii::Point< dim, VectorizedArrayType > hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::get_quadrature_point |
( |
const unsigned int | qx, |
|
|
const unsigned int | qv ) const |
|
inline |
Return coordinate of quadrature point (qx
, qv
).
template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
void hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::reinit |
( |
typename MatrixFree< dim_x, dim_v, Number, VectorizedArrayType >::ID | face_index | ) |
|
|
inlineoverridevirtual |
template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
template<SpaceType stype>
DEAL_II_ALWAYS_INLINE void hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::submit_value |
( |
VectorizedArrayType *__restrict | data_ptr, |
|
|
const VectorizedArrayType & | value, |
|
|
const unsigned int | q, |
|
|
const unsigned int | qx, |
|
|
const unsigned int | qv ) const |
|
inline |
Submit value (i.e. multiply by JxW).
template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
template<SpaceType stype>
DEAL_II_ALWAYS_INLINE void hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::submit_value |
( |
VectorizedArrayType *__restrict | data_ptr_1, |
|
|
VectorizedArrayType *__restrict | data_ptr_2, |
|
|
const VectorizedArrayType & | value, |
|
|
const unsigned int | q, |
|
|
const unsigned int | q1, |
|
|
const unsigned int | q2 ) const |
|
inline |
The same as above but tested from both sides -> FCL.