A class that performs interpolation from a face E-vector to quadrature point values and/or derivatives (Q-vectors) on the faces.
More...
#include <quadinterpolator_face.hpp>
|
| | FaceQuadratureInterpolator (const FiniteElementSpace &fes, const IntegrationRule &ir, FaceType type) |
| |
| void | DisableTensorProducts (bool disable=true) const |
| | Disable the use of tensor product evaluations, for tensor-product elements, e.g. quads and hexes.
|
| |
| QVectorLayout | GetOutputLayout () const |
| | Query the current output Q-vector layout. The default value is QVectorLayout::byNODES.
|
| |
| void | SetOutputLayout (QVectorLayout layout) const |
| | Set the desired output Q-vector layout. The default value is QVectorLayout::byNODES.
|
| |
| void | Mult (const Vector &e_vec, unsigned eval_flags, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor) const |
| | Interpolate the E-vector e_vec to quadrature points.
|
| |
| void | Values (const Vector &e_vec, Vector &q_val) const |
| | Interpolate the values of the E-vector e_vec at quadrature points.
|
| |
|
| template<const int T_VDIM = 0, const int T_ND = 0, const int T_NQ = 0> |
| static void | Eval2D (const int NF, const int vdim, const QVectorLayout q_layout, const DofToQuad &maps, const Array< bool > &signs, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor, const int eval_flags) |
| | Template compute kernel for 2D.
|
| |
| template<const int T_VDIM = 0, const int T_ND = 0, const int T_NQ = 0> |
| static void | Eval3D (const int NF, const int vdim, const QVectorLayout q_layout, const DofToQuad &maps, const Array< bool > &signs, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor, const int eval_flags) |
| | Template compute kernel for 3D.
|
| |
| template<const int T_VDIM = 0, const int T_ND = 0, const int T_NQ = 0> |
| static void | SmemEval3D (const int NF, const int vdim, const QVectorLayout q_layout, const DofToQuad &maps, const Array< bool > &signs, const Vector &e_vec, Vector &q_val, Vector &q_der, Vector &q_det, Vector &q_nor, const int eval_flags) |
| |
A class that performs interpolation from a face E-vector to quadrature point values and/or derivatives (Q-vectors) on the faces.
A face E-vector represents the face-wise discontinuous version of the trace FE space and can be obtained, for example, from a GridFunction using the Operator returned by FiniteElementSpace::GetFaceRestriction().
Definition at line 25 of file quadinterpolator_face.hpp.
◆ FaceEvalFlags
| Enumerator |
|---|
| VALUES | Evaluate the values at quadrature points.
|
| DERIVATIVES | Evaluate the derivatives at quadrature points.
|
| DETERMINANTS | Assuming the derivative at quadrature points form a matrix, this flag can be used to compute and store their determinants. This flag can only be used in Mult().
|
| NORMALS | |
Definition at line 53 of file quadinterpolator_face.hpp.
◆ FaceQuadratureInterpolator()
◆ DisableTensorProducts()
| void mfem::FaceQuadratureInterpolator::DisableTensorProducts |
( |
bool | disable = true | ) |
const |
|
inline |
Disable the use of tensor product evaluations, for tensor-product elements, e.g. quads and hexes.
Currently, tensor product evaluations are not implemented and this method has no effect.
Definition at line 71 of file quadinterpolator_face.hpp.
◆ Eval2D()
template<const int T_VDIM, const int T_ND1D, const int T_NQ1D>
| void mfem::FaceQuadratureInterpolator::Eval2D |
( |
const int | NF, |
|
|
const int | vdim, |
|
|
const QVectorLayout | q_layout, |
|
|
const DofToQuad & | maps, |
|
|
const Array< bool > & | signs, |
|
|
const Vector & | e_vec, |
|
|
Vector & | q_val, |
|
|
Vector & | q_der, |
|
|
Vector & | q_det, |
|
|
Vector & | q_nor, |
|
|
const int | eval_flags ) |
|
static |
◆ Eval3D()
template<const int T_VDIM, const int T_ND1D, const int T_NQ1D>
| void mfem::FaceQuadratureInterpolator::Eval3D |
( |
const int | NF, |
|
|
const int | vdim, |
|
|
const QVectorLayout | q_layout, |
|
|
const DofToQuad & | maps, |
|
|
const Array< bool > & | signs, |
|
|
const Vector & | e_vec, |
|
|
Vector & | q_val, |
|
|
Vector & | q_der, |
|
|
Vector & | q_det, |
|
|
Vector & | q_nor, |
|
|
const int | eval_flags ) |
|
static |
◆ GetOutputLayout()
| QVectorLayout mfem::FaceQuadratureInterpolator::GetOutputLayout |
( |
| ) |
const |
|
inline |
◆ Mult()
| void mfem::FaceQuadratureInterpolator::Mult |
( |
const Vector & | e_vec, |
|
|
unsigned | eval_flags, |
|
|
Vector & | q_val, |
|
|
Vector & | q_der, |
|
|
Vector & | q_det, |
|
|
Vector & | q_nor ) const |
Interpolate the E-vector e_vec to quadrature points.
The eval_flags are a bitwise mask of constants from the FaceEvalFlags enumeration. When the VALUES flag is set, the values at quadrature points are computed and stored in the Vector q_val. Similarly, when the flag DERIVATIVES is set, the derivatives are computed and stored in q_der. When the DETERMINANTS flags is set, it is assumed that the derivatives form a matrix at each quadrature point (i.e. the associated FiniteElementSpace is a vector space) and their determinants are computed and stored in q_det.
Definition at line 642 of file quadinterpolator_face.cpp.
◆ SetOutputLayout()
| void mfem::FaceQuadratureInterpolator::SetOutputLayout |
( |
QVectorLayout | layout | ) |
const |
|
inline |
◆ SmemEval3D()
template<const int T_VDIM, const int T_ND1D, const int T_NQ1D>
| void mfem::FaceQuadratureInterpolator::SmemEval3D |
( |
const int | NF, |
|
|
const int | vdim, |
|
|
const QVectorLayout | q_layout, |
|
|
const DofToQuad & | maps, |
|
|
const Array< bool > & | signs, |
|
|
const Vector & | e_vec, |
|
|
Vector & | q_val, |
|
|
Vector & | q_der, |
|
|
Vector & | q_det, |
|
|
Vector & | q_nor, |
|
|
const int | eval_flags ) |
|
static |
◆ Values()
| void mfem::FaceQuadratureInterpolator::Values |
( |
const Vector & | e_vec, |
|
|
Vector & | q_val ) const |
◆ FiniteElementSpace
◆ fespace
◆ IntRule
◆ MAX_ND1D
| const int mfem::FaceQuadratureInterpolator::MAX_ND1D = 10 |
|
staticprotected |
◆ MAX_ND2D
| const int mfem::FaceQuadratureInterpolator::MAX_ND2D = 100 |
|
staticprotected |
◆ MAX_ND3D
| const int mfem::FaceQuadratureInterpolator::MAX_ND3D = 1000 |
|
staticprotected |
◆ MAX_NQ1D
| const int mfem::FaceQuadratureInterpolator::MAX_NQ1D = 10 |
|
staticprotected |
◆ MAX_NQ2D
| const int mfem::FaceQuadratureInterpolator::MAX_NQ2D = 100 |
|
staticprotected |
◆ MAX_NQ3D
| const int mfem::FaceQuadratureInterpolator::MAX_NQ3D = 1000 |
|
staticprotected |
◆ MAX_VDIM1D
| const int mfem::FaceQuadratureInterpolator::MAX_VDIM1D = 1 |
|
staticprotected |
◆ MAX_VDIM2D
| const int mfem::FaceQuadratureInterpolator::MAX_VDIM2D = 2 |
|
staticprotected |
◆ MAX_VDIM3D
| const int mfem::FaceQuadratureInterpolator::MAX_VDIM3D = 3 |
|
staticprotected |
◆ q_layout
◆ use_tensor_products
| bool mfem::FaceQuadratureInterpolator::use_tensor_products |
|
mutableprotected |
The documentation for this class was generated from the following files: