hyper.deal
Loading...
Searching...
No Matches
hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType > Class Template Reference

#include <fe_evaluation_face.h>

Inheritance diagram for hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >:
hyperdeal::FEEvaluationBase< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >

Public Types

using PARENT
 
using SpaceType
 
- Public Types inherited from hyperdeal::FEEvaluationBase< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >
using NUMBER_ = Number
 
using VEC_NUMBER_ = VectorizedArrayType
 
using MF = MatrixFree<dim_x, dim_v, Number, VectorizedArrayType>
 
using VectorizedArrayTypeX = typename MF::VectorizedArrayTypeX
 
using VectorizedArrayTypeV = typename MF::VectorizedArrayTypeV
 

Public Member Functions

 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
 
- Public Member Functions inherited from hyperdeal::FEEvaluationBase< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >
 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 Public Attributes

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 Public Attributes inherited from hyperdeal::FEEvaluationBase< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >
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
 

Additional Inherited Members

- Protected Attributes inherited from hyperdeal::FEEvaluationBase< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >
dealii::AlignedVector< VectorizedArrayType > data
 
const MFmatrix_free
 
const dealii::MatrixFree< dim_x, Number, VectorizedArrayTypeX > & matrix_free_x
 
const dealii::MatrixFree< dim_v, Number, VectorizedArrayTypeV > & matrix_free_v
 
unsigned int macro_cell_x
 
unsigned int macro_cell_v
 
unsigned int lane_y
 
unsigned int macro
 
const dealii::AlignedVector< Number > * shape_values
 
const dealii::AlignedVector< Number > * shape_gradients
 

Detailed Description

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.

Member Typedef Documentation

◆ PARENT

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
using hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::PARENT
Initial value:
dim_v,
degree,
n_points,
Number,
VectorizedArrayType>
FEEvaluationBase(const MatrixFree< dim_x, dim_v, Number, VectorizedArrayType > &matrix_free)
Definition fe_evaluation_base.h:144

◆ SpaceType

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
using hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::SpaceType
Initial value:
typename MatrixFree<dim_x, dim_v, Number, VectorizedArrayType>::ID::
SpaceType

Constructor & Destructor Documentation

◆ FEFaceEvaluation()

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::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 )
inline

Constructor.

Parameters
matrix_freeData object that contains all data.
is_minus_faceThis selects which of the two cells of an internal face the current evaluator will be based upon. The interior face is the main face along which the normal vectors are oriented. The exterior face coming from the other side provides the same normal vector as the interior side, so if the outer normal vector to that side is desired, it must be multiplied by -1..
dof_no_xIf x-space matrix_free of matrix_free was set up with multiple DoFHandler objects, this parameter selects to which DoFHandler/AffineConstraints pair the given evaluator should be attached to.
dof_no_vSame as above but for v-space.
quad_no_xIf x-space matrix_free of matrix_free was set up with multiple Quadrature objects, this parameter selects the appropriate number of the quadrature formula.
quad_no_vSame as above but for v-space.

Member Function Documentation

◆ distribute_local_to_global()

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 >::distribute_local_to_global ( dealii::LinearAlgebra::distributed::Vector< Number > & dst) const
inline

Add the face data into a global vector.

◆ distribute_to_buffer()

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 >::distribute_to_buffer ( VectorizedArrayType * dst) const
inline

Add the face data into a cell buffer held e.g. by FEEvaluation.

◆ get_normal_vector_v()

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim_v, VectorizedArrayType > hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::get_normal_vector_v ( const unsigned int qv) const
inline

Get normal for a y-face

◆ get_normal_vector_x()

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
DEAL_II_ALWAYS_INLINE dealii::Tensor< 1, dim_x, VectorizedArrayType > hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::get_normal_vector_x ( const unsigned int qx) const
inline

Get normal for an x-face

◆ get_quadrature_point() [1/2]

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.

◆ get_quadrature_point() [2/2]

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).

◆ n_vectorization_lanes_filled()

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
unsigned int hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::n_vectorization_lanes_filled ( ) const
inline

Return number of filled lanes.

◆ read_dof_values() [1/2]

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 >::read_dof_values ( const dealii::LinearAlgebra::distributed::Vector< Number > & src)
inline

Read data from a global vector into the internal buffer.

◆ read_dof_values() [2/2]

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 >::read_dof_values ( const dealii::LinearAlgebra::distributed::Vector< Number > & src,
VectorizedArrayType * data ) const
inline

Read data from a global vector into a given buffer.

◆ reinit() [1/2]

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 cell_index,
unsigned int face )
inline

Set the view the current face in the context of ECL.

◆ reinit() [2/2]

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

Set the view the current face in the context of FCL.

Reimplemented from hyperdeal::FEEvaluationBase< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >.

◆ submit_value() [1/2]

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).

◆ submit_value() [2/2]

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.

Member Data Documentation

◆ N_Q_POINTS

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
const unsigned int hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::N_Q_POINTS
static
Initial value:
=
dealii::Utilities::pow(n_points, DIM - 1)

◆ N_Q_POINTS_1_CELL

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
const unsigned int hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::N_Q_POINTS_1_CELL
static
Initial value:
=
dealii::Utilities::pow(n_points, dim_x - 0)

◆ N_Q_POINTS_1_FACE

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
const unsigned int hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::N_Q_POINTS_1_FACE
static
Initial value:
=
dealii::Utilities::pow(n_points, dim_x - 1)

◆ N_Q_POINTS_2_CELL

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
const unsigned int hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::N_Q_POINTS_2_CELL
static
Initial value:
=
dealii::Utilities::pow(n_points, dim_v - 0)

◆ N_Q_POINTS_2_FACE

template<int dim_x, int dim_v, int degree, int n_points, typename Number , typename VectorizedArrayType >
const unsigned int hyperdeal::FEFaceEvaluation< dim_x, dim_v, degree, n_points, Number, VectorizedArrayType >::N_Q_POINTS_2_FACE
static
Initial value:
=
dealii::Utilities::pow(n_points, dim_v - 1)

The documentation for this class was generated from the following file: