Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Public Member Functions | Static Public Attributes | Protected Member Functions | Protected Attributes | Friends | List of all members
FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType > Class Template Reference

#include <deal.II/matrix_free/fe_evaluation.h>

Inheritance diagram for FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >:
[legend]

Public Member Functions

 ~FEEvaluationBaseData ()
 
unsigned int get_mapping_data_index_offset () const
 
internal::MatrixFreeFunctions::GeometryType get_cell_type () const
 
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > & get_shape_info () const
 
const internal::MatrixFreeFunctions::DoFInfoget_dof_info () const
 
VectorizedArrayType JxW (const unsigned int q_point) const
 
Tensor< 2, dim, VectorizedArrayType > inverse_jacobian (const unsigned int q_point) const
 
Tensor< 1, dim, VectorizedArrayType > get_normal_vector (const unsigned int q_point) const
 
VectorizedArrayType read_cell_data (const AlignedVector< VectorizedArrayType > &array) const
 
void set_cell_data (AlignedVector< VectorizedArrayType > &array, const VectorizedArrayType &value) const
 
template<typename T >
std::array< T, VectorizedArrayType::size()> read_cell_data (const AlignedVector< std::array< T, VectorizedArrayType::size()> > &array) const
 
template<typename T >
void set_cell_data (AlignedVector< std::array< T, VectorizedArrayType::size()> > &array, const std::array< T, VectorizedArrayType::size()> &value) const
 
std::array< unsigned int, VectorizedArrayType::size()> get_cell_ids () const
 
std::array< unsigned int, VectorizedArrayType::size()> get_cell_or_face_ids () const
 
const std::vector< unsigned int > & get_internal_dof_numbering () const
 
ArrayView< VectorizedArrayType > get_scratch_data () const
 
unsigned int get_quadrature_index () const
 
unsigned int get_current_cell_index () const
 
unsigned int get_active_fe_index () const
 
unsigned int get_active_quadrature_index () const
 
const MatrixFree< dim, Number, VectorizedArrayType > & get_matrix_free () const
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 

Protected Member Functions

 FEEvaluationBaseData (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const unsigned int dof_no, const unsigned int first_selected_component, const unsigned int quad_no, const unsigned int fe_degree, const unsigned int n_q_points, const bool is_interior_face, const unsigned int active_fe_index, const unsigned int active_quad_index, const unsigned int face_type)
 
 FEEvaluationBaseData (const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags, const unsigned int first_selected_component, const FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType > *other)
 
 FEEvaluationBaseData (const FEEvaluationBaseData &other)
 
FEEvaluationBaseDataoperator= (const FEEvaluationBaseData &other)
 

Protected Attributes

AlignedVector< VectorizedArrayType > * scratch_data_array
 
VectorizedArrayType * scratch_data
 
const unsigned int quad_no
 
const MatrixFree< dim, Number, VectorizedArrayType > * matrix_info
 
const internal::MatrixFreeFunctions::DoFInfodof_info
 
const internal::MatrixFreeFunctions::MappingInfoStorage<(is_face ? dim - 1 :dim), dim, Number, VectorizedArrayType > * mapping_data
 
const unsigned int active_fe_index
 
const unsigned int active_quad_index
 
const internal::MatrixFreeFunctions::MappingInfoStorage<(is_face?dim-1:dim), dim, Number, VectorizedArrayType >::QuadratureDescriptor * descriptor
 
const unsigned int n_quadrature_points
 
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > * data
 
const Tensor< 2, dim, VectorizedArrayType > * jacobian
 
const VectorizedArrayType * J_value
 
const Tensor< 1, dim, VectorizedArrayType > * normal_vectors
 
const Tensor< 1, dim, VectorizedArrayType > * normal_x_jacobian
 
const Number * quadrature_weights
 
unsigned int cell
 
bool is_interior_face
 
internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex dof_access_index
 
unsigned int face_no
 
unsigned int face_orientation
 
unsigned int subface_index
 
internal::MatrixFreeFunctions::GeometryType cell_type
 
std::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number, VectorizedArrayType > > mapped_geometry
 

Friends

template<int , int , int , int , typename , typename >
class FEEvaluation
 

Detailed Description

template<int dim, typename Number, bool is_face, typename VectorizedArrayType>
class FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >

This base class of the FEEvaluation and FEFaceEvaluation classes handles mapping-related information independent of the degrees of freedom and finite element in use. This class provides access functionality for user code but is otherwise invisible without any public constructor. The usage is through the class FEEvaluation instead.

This class has four template arguments:

Template Parameters
dimDimension in which this class is to be used
NumberNumber format, usually double or float
is_faceWhether the class is used for a cell integrator (with quadrature dimension the same as the space dimension) or for a face integrator (with quadrature dimension one less)
VectorizedArrayTypeType of array to be woked on in a vectorized fashion, defaults to VectorizedArray<Number>
Note
Currently only VectorizedArray<Number, width> is supported as VectorizedArrayType.

Definition at line 101 of file fe_evaluation.h.

Constructor & Destructor Documentation

◆ ~FEEvaluationBaseData()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::~FEEvaluationBaseData ( )

Destructor.

◆ FEEvaluationBaseData() [1/3]

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::FEEvaluationBaseData ( const MatrixFree< dim, Number, VectorizedArrayType > &  matrix_free,
const unsigned int  dof_no,
const unsigned int  first_selected_component,
const unsigned int  quad_no,
const unsigned int  fe_degree,
const unsigned int  n_q_points,
const bool  is_interior_face,
const unsigned int  active_fe_index,
const unsigned int  active_quad_index,
const unsigned int  face_type 
)
protected

Constructor. Made protected to prevent users from directly using this class. Takes all data stored in MatrixFree. If applied to problems with more than one quadrature formula selected during construction of matrix_free, quad_no allows to select the appropriate formula.

◆ FEEvaluationBaseData() [2/3]

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::FEEvaluationBaseData ( const Mapping< dim > &  mapping,
const FiniteElement< dim > &  fe,
const Quadrature< 1 > &  quadrature,
const UpdateFlags  update_flags,
const unsigned int  first_selected_component,
const FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType > *  other 
)
protected

Constructor that comes with reduced functionality and works similar as FEValues.

◆ FEEvaluationBaseData() [3/3]

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::FEEvaluationBaseData ( const FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType > &  other)
protected

Copy constructor. If FEEvaluationBase was constructed from a mapping, fe, quadrature, and update flags, the underlying geometry evaluation based on FEValues will be deep-copied in order to allow for using in parallel with threads.

Member Function Documentation

◆ get_mapping_data_index_offset()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_mapping_data_index_offset ( ) const

Return the index offset within the geometry fields for the cell the reinit() function has been called for. This index can be used to access an index into a field that has the same compression behavior as the Jacobian of the geometry, e.g., to store an effective coefficient tensors that combines a coefficient with the geometry for lower memory transfer as the available data fields.

◆ get_cell_type()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
internal::MatrixFreeFunctions::GeometryType FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_cell_type ( ) const

Return the type of the cell the reinit() function has been called for. Valid values are cartesian for Cartesian cells (which allows for considerable data compression), affine for cells with affine mappings, and general for general cells without any compressed storage applied.

◆ get_shape_info()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > & FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_shape_info ( ) const

Return a reference to the ShapeInfo object currently in use.

◆ get_dof_info()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const internal::MatrixFreeFunctions::DoFInfo & FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_dof_info ( ) const

Return a reference to the DoFInfo object currently in use.

◆ JxW()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
VectorizedArrayType FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::JxW ( const unsigned int  q_point) const

Return the determinant of the Jacobian from the unit to the real cell times the quadrature weight.

◆ inverse_jacobian()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
Tensor< 2, dim, VectorizedArrayType > FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::inverse_jacobian ( const unsigned int  q_point) const

Return the inverse and transposed version \(J^{-\mathrm T}\) of the Jacobian of the mapping between the unit to the real cell defined as \(J_{ij} = d x_i / d\hat x_j\). The \((i,j)\) entry of the returned tensor contains \(d\hat x_j/dx_i\), i.e., columns refer to reference space coordinates and rows to real cell coordinates. Thus, the returned tensor represents a covariant transformation, which is used in the FEEvaluationBase::get_gradient() function to transform the unit cell gradients to gradients on the real cell by a multiplication \(J^{-\mathrm T} \hat{\nabla} u_h\).

◆ get_normal_vector()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
Tensor< 1, dim, VectorizedArrayType > FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_normal_vector ( const unsigned int  q_point) const

Return the unit normal vector on a face. Note that both sides of a face use the same orientation of the normal vector: For the faces enumerated as interior in FaceToCellTopology and selected with the is_interior_face=true flag of the constructor, this corresponds to the outer normal vector, whereas for faces enumerated as exterior in FaceToCellTopology and selected with the is_interior_face=false flag of the constructor, the normal points into the element as a consequence of the single normal vector.

Note
Only implemented in case is_face == true.

◆ read_cell_data() [1/2]

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
VectorizedArrayType FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::read_cell_data ( const AlignedVector< VectorizedArrayType > &  array) const

Provides a unified interface to access data in a vector of VectorizedArray fields of length MatrixFree::n_cell_batches() + MatrixFree::n_ghost_cell_batches() for both cells (plain read) and faces (indirect addressing).

◆ set_cell_data() [1/2]

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
void FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::set_cell_data ( AlignedVector< VectorizedArrayType > &  array,
const VectorizedArrayType &  value 
) const

Provides a unified interface to set data in a vector of VectorizedArray fields of length MatrixFree::n_cell_batches() + MatrixFree::n_ghost_cell_batches() for both cells (plain read) and faces (indirect addressing).

◆ read_cell_data() [2/2]

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
template<typename T >
std::array< T, VectorizedArrayType::size()> FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::read_cell_data ( const AlignedVector< std::array< T, VectorizedArrayType::size()> > &  array) const

The same as above, just for std::array of length of VectorizedArrayType for arbitrary data type.

◆ set_cell_data() [2/2]

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
template<typename T >
void FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::set_cell_data ( AlignedVector< std::array< T, VectorizedArrayType::size()> > &  array,
const std::array< T, VectorizedArrayType::size()> &  value 
) const

The same as above, just for std::array of length of VectorizedArrayType for arbitrary data type.

◆ get_cell_ids()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
std::array< unsigned int, VectorizedArrayType::size()> FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_cell_ids ( ) const

Return the id of the cells this FEEvaluation or FEFaceEvaluation is associated with.

◆ get_cell_or_face_ids()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
std::array< unsigned int, VectorizedArrayType::size()> FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_cell_or_face_ids ( ) const

Return the id of the cells/faces this FEEvaluation/FEFaceEvaluation is associated with.

◆ get_internal_dof_numbering()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const std::vector< unsigned int > & FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_internal_dof_numbering ( ) const

Return the numbering of local degrees of freedom within the evaluation routines of FEEvaluation in terms of the standard numbering on finite elements.

◆ get_scratch_data()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
ArrayView< VectorizedArrayType > FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_scratch_data ( ) const

Return an ArrayView to internal memory for temporary use. Note that some of this memory is overwritten during evaluate() and integrate() calls so do not assume it to be stable over those calls. The maximum size you can write into is 3*dofs_per_cell+2*n_q_points.

◆ get_quadrature_index()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_quadrature_index ( ) const

Return the number of the quadrature formula of the present cell.

◆ get_current_cell_index()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_current_cell_index ( ) const

Return index of the current cell or face.

◆ get_active_fe_index()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_active_fe_index ( ) const

Return the active FE index for this class for efficient indexing in the hp- case.

◆ get_active_quadrature_index()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_active_quadrature_index ( ) const

Return the active quadrature index for this class for efficient indexing in the hp-case.

◆ get_matrix_free()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const MatrixFree< dim, Number, VectorizedArrayType > & FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::get_matrix_free ( ) const

Return the underlying MatrixFree object.

◆ operator=()

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
FEEvaluationBaseData & FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::operator= ( const FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType > &  other)
protected

Copy assignment operator. If FEEvaluationBase was constructed from a mapping, fe, quadrature, and update flags, the underlying geometry evaluation based on FEValues will be deep-copied in order to allow for using in parallel with threads.

Friends And Related Function Documentation

◆ FEEvaluation

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
template<int , int , int , int , typename , typename >
friend class FEEvaluation
friend

Definition at line 497 of file fe_evaluation.h.

Member Data Documentation

◆ dimension

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
constexpr unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::dimension = dim
staticconstexpr

Definition at line 108 of file fe_evaluation.h.

◆ scratch_data_array

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
AlignedVector<VectorizedArrayType>* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::scratch_data_array
protected

This is the general array for all data fields.

Definition at line 337 of file fe_evaluation.h.

◆ scratch_data

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
VectorizedArrayType* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::scratch_data
protected

This is the user-visible part of scratch_data_array, only showing the last part of scratch_data_array. The first part is consumed by values_dofs, values_quad, etc.

Definition at line 344 of file fe_evaluation.h.

◆ quad_no

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::quad_no
protected

The number of the quadrature formula of the present cell.

Definition at line 349 of file fe_evaluation.h.

◆ matrix_info

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const MatrixFree<dim, Number, VectorizedArrayType>* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::matrix_info
protected

A pointer to the underlying data.

Definition at line 354 of file fe_evaluation.h.

◆ dof_info

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const internal::MatrixFreeFunctions::DoFInfo* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::dof_info
protected

A pointer to the underlying DoF indices and constraint description for the component specified at construction. Also contained in matrix_info, but it simplifies code if we store a reference to it.

Definition at line 361 of file fe_evaluation.h.

◆ mapping_data

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const internal::MatrixFreeFunctions::MappingInfoStorage< (is_face ? dim - 1 : dim), dim, Number, VectorizedArrayType>* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::mapping_data
protected

A pointer to the underlying transformation data from unit to real cells for the given quadrature formula specified at construction. Also contained in matrix_info, but it simplifies code if we store a reference to it.

Definition at line 373 of file fe_evaluation.h.

◆ active_fe_index

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::active_fe_index
protected

The active FE index for this class for efficient indexing in the hp-case.

Definition at line 378 of file fe_evaluation.h.

◆ active_quad_index

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::active_quad_index
protected

The active quadrature index for this class for efficient indexing in the hp-case.

Definition at line 384 of file fe_evaluation.h.

◆ descriptor

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const internal::MatrixFreeFunctions::MappingInfoStorage<(is_face?dim-1:dim),dim,Number,VectorizedArrayType>::QuadratureDescriptor* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::descriptor
protected

A pointer to the underlying quadrature formula specified at construction. Also contained in matrix_info, but it simplifies code if we store a reference to it.

Definition at line 395 of file fe_evaluation.h.

◆ n_quadrature_points

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::n_quadrature_points
protected

The number of quadrature points in the current evaluation context.

Definition at line 400 of file fe_evaluation.h.

◆ data

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const internal::MatrixFreeFunctions::ShapeInfo<VectorizedArrayType>* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::data
protected

A pointer to the unit cell shape data, i.e., values, gradients and Hessians in 1D at the quadrature points that constitute the tensor product. Also contained in matrix_info, but it simplifies code if we store a reference to it.

Definition at line 408 of file fe_evaluation.h.

◆ jacobian

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const Tensor<2, dim, VectorizedArrayType>* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::jacobian
protected

A pointer to the Jacobian information of the present cell. Only set to a useful value if on a non-Cartesian cell.

Definition at line 414 of file fe_evaluation.h.

◆ J_value

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const VectorizedArrayType* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::J_value
protected

A pointer to the Jacobian determinant of the present cell. If on a Cartesian cell or on a cell with constant Jacobian, this is just the Jacobian determinant, otherwise the Jacobian determinant times the quadrature weight.

Definition at line 422 of file fe_evaluation.h.

◆ normal_vectors

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const Tensor<1, dim, VectorizedArrayType>* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::normal_vectors
protected

A pointer to the normal vectors at faces.

Definition at line 427 of file fe_evaluation.h.

◆ normal_x_jacobian

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const Tensor<1, dim, VectorizedArrayType>* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::normal_x_jacobian
protected

A pointer to the normal vectors times the jacobian at faces.

Definition at line 432 of file fe_evaluation.h.

◆ quadrature_weights

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
const Number* FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::quadrature_weights
protected

A pointer to the quadrature weights of the underlying quadrature formula.

Definition at line 437 of file fe_evaluation.h.

◆ cell

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::cell
protected

After a call to reinit(), stores the number of the cell we are currently working with.

Definition at line 443 of file fe_evaluation.h.

◆ is_interior_face

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
bool FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::is_interior_face
protected

Flag holding information whether a face is an interior or exterior face according to the defined direction of the normal. Not used for cells.

Definition at line 449 of file fe_evaluation.h.

◆ dof_access_index

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::dof_access_index
protected

Stores the index an FEFaceEvaluation object is currently pointing into (interior face, exterior face, data associated with cell).

Definition at line 455 of file fe_evaluation.h.

◆ face_no

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::face_no
protected

Stores the current number of a face within the given cell in case is_face==true, using values between 0 and 2*dim.

Definition at line 461 of file fe_evaluation.h.

◆ face_orientation

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::face_orientation
protected

Stores the orientation of the given face with respect to the standard orientation, 0 if in standard orientation.

Definition at line 467 of file fe_evaluation.h.

◆ subface_index

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
unsigned int FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::subface_index
protected

Stores the subface index of the given face. Usually, this variable takes the value numbers::invalid_unsigned_int to indicate integration over the full face, but in case the current physical face has a neighbor that is more refined, it is a subface and must scale the entries in ShapeInfo appropriately.

Definition at line 476 of file fe_evaluation.h.

◆ cell_type

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
internal::MatrixFreeFunctions::GeometryType FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::cell_type
protected

Stores the type of the cell we are currently working with after a call to reinit(). Valid values are cartesian, affine and general, which have different implications on how the Jacobian transformations are stored internally in MappingInfo.

Definition at line 484 of file fe_evaluation.h.

◆ mapped_geometry

template<int dim, typename Number , bool is_face, typename VectorizedArrayType >
std::shared_ptr<internal::MatrixFreeFunctions:: MappingDataOnTheFly<dim, Number, VectorizedArrayType> > FEEvaluationBaseData< dim, Number, is_face, VectorizedArrayType >::mapped_geometry
protected

Geometry data that can be generated FEValues on the fly with the respective constructor.

Definition at line 492 of file fe_evaluation.h.


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