Reference documentation for deal.II version 9.1.0-pre
Protected Member Functions | Protected Attributes | Private Member Functions | Friends | List of all members
FEEvaluationBase< dim, n_components_, Number, is_face > Class Template Reference

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

Inheritance diagram for FEEvaluationBase< dim, n_components_, Number, is_face >:
[legend]

Public Member Functions

1: General operations
 ~FEEvaluationBase ()
 
unsigned int get_cell_data_number () const
 
unsigned int get_mapping_data_index_offset () const
 
internal::MatrixFreeFunctions::GeometryType get_cell_type () const
 
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArray< Number > > & get_shape_info () const
 
2: Reading from and writing to vectors
template<typename VectorType >
void read_dof_values (const VectorType &src, const unsigned int first_index=0)
 
template<typename VectorType >
void read_dof_values_plain (const VectorType &src, const unsigned int first_index=0)
 
template<typename VectorType >
void distribute_local_to_global (VectorType &dst, const unsigned int first_index=0) const
 
template<typename VectorType >
void set_dof_values (VectorType &dst, const unsigned int first_index=0) const
 
3: Data access
value_type get_dof_value (const unsigned int dof) const
 
void submit_dof_value (const value_type val_in, const unsigned int dof)
 
value_type get_value (const unsigned int q_point) const
 
void submit_value (const value_type val_in, const unsigned int q_point)
 
gradient_type get_gradient (const unsigned int q_point) const
 
value_type get_normal_derivative (const unsigned int q_point) const
 
void submit_gradient (const gradient_type grad_in, const unsigned int q_point)
 
void submit_normal_derivative (const value_type grad_in, const unsigned int q_point)
 
Tensor< 1, n_components_, Tensor< 2, dim, VectorizedArray< Number > > > get_hessian (const unsigned int q_point) const
 
gradient_type get_hessian_diagonal (const unsigned int q_point) const
 
value_type get_laplacian (const unsigned int q_point) const
 
VectorizedArray< Number > get_divergence (const unsigned int q_point) const
 
SymmetricTensor< 2, dim, VectorizedArray< Number > > get_symmetric_gradient (const unsigned int q_point) const
 
Tensor< 1,(dim==2?1:dim), VectorizedArray< Number > > get_curl (const unsigned int q_point) const
 
void submit_divergence (const VectorizedArray< Number > div_in, const unsigned int q_point)
 
void submit_symmetric_gradient (const SymmetricTensor< 2, dim, VectorizedArray< Number >> grad_in, const unsigned int q_point)
 
void submit_curl (const Tensor< 1, dim==2?1:dim, VectorizedArray< Number >> curl_in, const unsigned int q_point)
 
value_type integrate_value () const
 
VectorizedArray< Number > JxW (const unsigned int q_index) const
 
void fill_JxW_values (AlignedVector< VectorizedArray< Number >> &JxW_values) const
 
Tensor< 2, dim, VectorizedArray< Number > > inverse_jacobian (const unsigned int q_index) const
 
Tensor< 1, dim, VectorizedArray< Number > > get_normal_vector (const unsigned int q_point) const
 
VectorizedArray< Number > read_cell_data (const AlignedVector< VectorizedArray< Number >> &array) const
 
4: Access to internal data
const VectorizedArray< Number > * begin_dof_values () const
 
VectorizedArray< Number > * begin_dof_values ()
 
const VectorizedArray< Number > * begin_values () const
 
VectorizedArray< Number > * begin_values ()
 
const VectorizedArray< Number > * begin_gradients () const
 
VectorizedArray< Number > * begin_gradients ()
 
const VectorizedArray< Number > * begin_hessians () const
 
VectorizedArray< Number > * begin_hessians ()
 
const std::vector< unsigned int > & get_internal_dof_numbering () const
 
ArrayView< VectorizedArray< Number > > get_scratch_data () const
 

Protected Member Functions

 FEEvaluationBase (const MatrixFree< dim, Number > &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)
 
template<int n_components_other>
 FEEvaluationBase (const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags, const unsigned int first_selected_component, const FEEvaluationBase< dim, n_components_other, Number > *other)
 
 FEEvaluationBase (const FEEvaluationBase &other)
 
FEEvaluationBaseoperator= (const FEEvaluationBase &other)
 
template<typename VectorType , typename VectorOperation >
void read_write_operation (const VectorOperation &operation, VectorType *vectors[], const bool apply_constraints=true) const
 
template<typename VectorType , typename VectorOperation >
void read_write_operation_contiguous (const VectorOperation &operation, VectorType *vectors[]) const
 
template<typename VectorType , typename VectorOperation >
void read_write_operation_global (const VectorOperation &operation, VectorType *vectors[]) const
 

Protected Attributes

AlignedVector< VectorizedArray< Number > > * scratch_data_array
 
VectorizedArray< Number > * scratch_data
 
VectorizedArray< Number > * values_dofs [n_components]
 
VectorizedArray< Number > * values_quad [n_components]
 
VectorizedArray< Number > * gradients_quad [n_components][dim]
 
VectorizedArray< Number > * hessians_quad [n_components][(dim *(dim+1))/2]
 
const unsigned int quad_no
 
const unsigned int n_fe_components
 
const unsigned int active_fe_index
 
const unsigned int active_quad_index
 
const unsigned int n_quadrature_points
 
const MatrixFree< dim, Number > * matrix_info
 
const internal::MatrixFreeFunctions::DoFInfodof_info
 
const internal::MatrixFreeFunctions::MappingInfoStorage<(is_face?dim-1:dim), dim, Number > * mapping_data
 
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArray< Number > > * data
 
const Tensor< 2, dim, VectorizedArray< Number > > * jacobian
 
const VectorizedArray< Number > * J_value
 
const Tensor< 1, dim, VectorizedArray< Number > > * normal_vectors
 
const Tensor< 1, dim, VectorizedArray< Number > > * 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
 
bool dof_values_initialized
 
bool values_quad_initialized
 
bool gradients_quad_initialized
 
bool hessians_quad_initialized
 
bool values_quad_submitted
 
bool gradients_quad_submitted
 
std::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly< dim, Number > > mapped_geometry
 
const unsigned int first_selected_component
 
std::vector< types::global_dof_indexlocal_dof_indices
 

Private Member Functions

void set_data_pointers ()
 

Friends

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

Detailed Description

template<int dim, int n_components_, typename Number, bool is_face = false>
class FEEvaluationBase< dim, n_components_, Number, is_face >

This is the base class for the FEEvaluation classes. This class is a base class and needs usually not be called in user code. It does not have any public constructor. The usage is through the class FEEvaluation instead. It implements a reinit method that is used to set pointers so that operations on quadrature points can be performed quickly, access functions to vectors for the read_dof_values, set_dof_values, and distribute_local_to_global functions, as well as methods to access values and gradients of finite element functions.

This class has three template arguments:

Template Parameters
dimDimension in which this class is to be used
n_componentsNumber of vector components when solving a system of PDEs. If the same operation is applied to several components of a PDE (e.g. a vector Laplace equation), they can be applied simultaneously with one call (and often more efficiently)
NumberNumber format, usually double or float
Author
Katharina Kormann and Martin Kronbichler, 2010-2018

Definition at line 91 of file fe_evaluation.h.

Constructor & Destructor Documentation

template<int dim, int n_components_, typename Number, bool is_face = false>
FEEvaluationBase< dim, n_components_, Number, is_face >::~FEEvaluationBase ( )

Destructor.

template<int dim, int n_components_, typename Number, bool is_face = false>
FEEvaluationBase< dim, n_components_, Number, is_face >::FEEvaluationBase ( const MatrixFree< dim, Number > &  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 
)
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 finite element or more than one quadrature formula selected during construction of matrix_free, dof_no, first_selected_component and quad_no allow to select the appropriate components.

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

Constructor that comes with reduced functionality and works similar as FEValues. The arguments are similar to the ones passed to the constructor of FEValues, with the notable difference that FEEvaluation expects a one- dimensional quadrature formula, Quadrature<1>, instead of a dim dimensional one. The finite element can be both scalar or vector valued, but this method always only selects a scalar base element at a time (with n_components copies as specified by the class template argument). For vector-valued elements, the optional argument first_selected_component allows to specify the index of the base element to be used for evaluation. Note that the internal data structures always assume that the base element is primitive, non-primitive are not supported currently.

As known from FEValues, a call to the reinit method with a Triangulation::cell_iterator is necessary to make the geometry and degrees of freedom of the current class known. If the iterator includes DoFHandler information (i.e., it is a DoFHandler::cell_iterator or similar), the initialization allows to also read from or write to vectors in the standard way for DoFHandler::active_cell_iterator types for one cell at a time. However, this approach is much slower than the path with MatrixFree with MPI since index translation has to be done. As only one cell at a time is used, this method does not vectorize over several elements (which is most efficient for vector operations), but only possibly within the element if the evaluate/integrate routines are combined inside user code (e.g. for computing cell matrices).

The optional FEEvaluationBase object allows several FEEvaluation objects to share the geometry evaluation, i.e., the underlying mapping and quadrature points do only need to be evaluated once. This only works if the quadrature formulas are the same. Otherwise, a new evaluation object is created. Make sure to not pass an optional object around when you intend to use the FEEvaluation object in parallel with another one because otherwise the intended sharing may create race conditions.

template<int dim, int n_components_, typename Number, bool is_face = false>
FEEvaluationBase< dim, n_components_, Number, is_face >::FEEvaluationBase ( const FEEvaluationBase< dim, n_components_, Number, is_face > &  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

template<int dim, int n_components_, typename Number, bool is_face = false>
unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::get_cell_data_number ( ) const
template<int dim, int n_components_, typename Number, bool is_face = false>
unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::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.

template<int dim, int n_components_, typename Number, bool is_face = false>
internal::MatrixFreeFunctions::GeometryType FEEvaluationBase< dim, n_components_, Number, is_face >::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.

template<int dim, int n_components_, typename Number, bool is_face = false>
const internal::MatrixFreeFunctions::ShapeInfo<VectorizedArray<Number> >& FEEvaluationBase< dim, n_components_, Number, is_face >::get_shape_info ( ) const

Return a reference to the ShapeInfo object currently in use.

template<int dim, int n_components_, typename Number, bool is_face = false>
template<typename VectorType >
void FEEvaluationBase< dim, n_components_, Number, is_face >::read_dof_values ( const VectorType &  src,
const unsigned int  first_index = 0 
)

For the vector src, read out the values on the degrees of freedom of the current cell, and store them internally. Similar functionality as the function DoFAccessor::get_interpolated_dof_values when no constraints are present, but it also includes constraints from hanging nodes, so one can see it as a similar function to AffineConstraints::read_dof_values as well. Note that if vectorization is enabled, the DoF values for several cells are set.

If some constraints on the vector are inhomogeneous, use the function read_dof_values_plain instead and provide the vector with useful data also in constrained positions by calling AffineConstraints::distribute. When accessing vector entries during the solution of linear systems, the temporary solution should always have homogeneous constraints and this method is the correct one.

If this class was constructed without a MatrixFree object and the information is acquired on the fly through a DoFHandler<dim>::cell_iterator, only one single cell is used by this class and this function extracts the values of the underlying components on the given cell. This call is slower than the ones done through a MatrixFree object and lead to a structure that does not effectively use vectorization in the evaluate routines based on these values (instead, VectorizedArray::n_array_elements same copies are worked on).

If the given vector template class is a block vector (determined through the template function 'IsBlockVector<VectorType>::value', which checks for vectors derived from BlockVectorBase) or an std::vector<VectorType> or std::vector<VectorType *>, this function reads n_components blocks from the block vector starting at the index first_index. For non-block vectors, first_index is ignored.

template<int dim, int n_components_, typename Number, bool is_face = false>
template<typename VectorType >
void FEEvaluationBase< dim, n_components_, Number, is_face >::read_dof_values_plain ( const VectorType &  src,
const unsigned int  first_index = 0 
)

For the vector src, read out the values on the degrees of freedom of the current cell, and store them internally. Similar functionality as the function DoFAccessor::get_interpolated_dof_values. As opposed to the read_dof_values function, this function reads out the plain entries from vectors, without taking stored constraints into account. This way of access is appropriate when the constraints have been distributed on the vector by a call to AffineConstraints::distribute previously. This function is also necessary when inhomogeneous constraints are to be used, as MatrixFree can only handle homogeneous constraints. Note that if vectorization is enabled, the DoF values for several cells are set.

If this class was constructed without a MatrixFree object and the information is acquired on the fly through a DoFHandler<dim>::cell_iterator, only one single cell is used by this class and this function extracts the values of the underlying components on the given cell. This call is slower than the ones done through a MatrixFree object and lead to a structure that does not effectively use vectorization in the evaluate routines based on these values (instead, VectorizedArray::n_array_elements same copies are worked on).

If the given vector template class is a block vector (determined through the template function 'IsBlockVector<VectorType>::value', which checks for vectors derived from BlockVectorBase) or an std::vector<VectorType> or std::vector<VectorType *>, this function reads n_components blocks from the block vector starting at the index first_index. For non-block vectors, first_index is ignored.

template<int dim, int n_components_, typename Number, bool is_face = false>
template<typename VectorType >
void FEEvaluationBase< dim, n_components_, Number, is_face >::distribute_local_to_global ( VectorType &  dst,
const unsigned int  first_index = 0 
) const

Takes the values stored internally on dof values of the current cell and sums them into the vector dst. The function also applies constraints during the write operation. The functionality is hence similar to the function AffineConstraints::distribute_local_to_global. If vectorization is enabled, the DoF values for several cells are used.

If this class was constructed without a MatrixFree object and the information is acquired on the fly through a DoFHandler<dim>::cell_iterator, only one single cell is used by this class and this function extracts the values of the underlying components on the given cell. This call is slower than the ones done through a MatrixFree object and lead to a structure that does not effectively use vectorization in the evaluate routines based on these values (instead, VectorizedArray::n_array_elements same copies are worked on).

If the given vector template class is a block vector (determined through the template function 'IsBlockVector<VectorType>::value', which checks for vectors derived from BlockVectorBase) or an std::vector<VectorType> or std::vector<VectorType *>, this function writes to n_components blocks of the block vector starting at the index first_index. For non-block vectors, first_index is ignored.

template<int dim, int n_components_, typename Number, bool is_face = false>
template<typename VectorType >
void FEEvaluationBase< dim, n_components_, Number, is_face >::set_dof_values ( VectorType &  dst,
const unsigned int  first_index = 0 
) const

Takes the values stored internally on dof values of the current cell and writes them into the vector dst. The function skips the degrees of freedom which are constrained. As opposed to the distribute_local_to_global method, the old values at the position given by the current cell are overwritten. Thus, if a degree of freedom is associated to more than one cell (as usual in continuous finite elements), the values will be overwritten and only the value written last is retained.

If this class was constructed without a MatrixFree object and the information is acquired on the fly through a DoFHandler<dim>::cell_iterator, only one single cell is used by this class and this function extracts the values of the underlying components on the given cell. This call is slower than the ones done through a MatrixFree object and lead to a structure that does not effectively use vectorization in the evaluate routines based on these values (instead, VectorizedArray::n_array_elements same copies are worked on).

If the given vector template class is a block vector (determined through the template function 'IsBlockVector<VectorType>::value', which checks for vectors derived from BlockVectorBase) or an std::vector<VectorType> or std::vector<VectorType *>, this function writes to n_components blocks of the block vector starting at the index first_index. For non-block vectors, first_index is ignored.

template<int dim, int n_components_, typename Number, bool is_face = false>
value_type FEEvaluationBase< dim, n_components_, Number, is_face >::get_dof_value ( const unsigned int  dof) const

Return the value stored for the local degree of freedom with index dof. If the object is vector-valued, a vector-valued return argument is given. Thus, the argument dof can at most run until dofs_per_component rather than dofs_per_cell since the different components of a vector-valued FE are return together. Note that when vectorization is enabled, values from several cells are grouped together. If set_dof_values was called last, the value corresponds to the one set there. If integrate was called last, it instead corresponds to the value of the integrated function with the test function of the given index.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::submit_dof_value ( const value_type  val_in,
const unsigned int  dof 
)

Write a value to the field containing the degrees of freedom with component dof. Writes to the same field as is accessed through get_dof_value. Therefore, the original data that was read from a vector is overwritten as soon as a value is submitted.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
value_type FEEvaluationBase< dim, n_components_, Number, is_face >::get_value ( const unsigned int  q_point) const

Return the value of a finite element function at quadrature point number q_point after a call to evaluate(true,...), or the value that has been stored there with a call to submit_value. If the object is vector-valued, a vector-valued return argument is given. Note that when vectorization is enabled, values from several cells are grouped together.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::submit_value ( const value_type  val_in,
const unsigned int  q_point 
)

Write a value to the field containing the values on quadrature points with component q_point. Access to the same field as through get_value. If applied before the function integrate(true,...) is called, this specifies the value which is tested by all basis function on the current cell and integrated over.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
gradient_type FEEvaluationBase< dim, n_components_, Number, is_face >::get_gradient ( const unsigned int  q_point) const

Return the gradient of a finite element function at quadrature point number q_point after a call to evaluate(...,true,...), or the value that has been stored there with a call to submit_gradient.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
value_type FEEvaluationBase< dim, n_components_, Number, is_face >::get_normal_derivative ( const unsigned int  q_point) const

Return the derivative of a finite element function at quadrature point number q_point after a call to evaluate(...,true,...) in the direction normal to the face: \(\boldsymbol \nabla u(\mathbf x_q) \cdot \mathbf n(\mathbf x_q)\)

This call is equivalent to calling get_gradient() * get_normal_vector() but will use a more efficient internal representation of data.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::submit_gradient ( const gradient_type  grad_in,
const unsigned int  q_point 
)

Write a contribution that is tested by the gradient to the field containing the values on quadrature points with component q_point. Access to the same field as through get_gradient(). If applied before the function integrate(...,true) is called, this specifies what is tested by all basis function gradients on the current cell and integrated over.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::submit_normal_derivative ( const value_type  grad_in,
const unsigned int  q_point 
)

Write a contribution that is tested by the gradient to the field containing the values on quadrature points with component q_point. Access to the same field as through get_gradient() or get_normal_derivative(). If applied before the function integrate(...,true) is called, this specifies what is tested by all basis function gradients on the current cell and integrated over.

Note
This operation writes the data to the same field as submit_gradient(). As a consequence, only one of these two can be used. Usually, the contribution of a potential call to this function must be added into the contribution for submit_gradient().
The derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).
template<int dim, int n_components_, typename Number, bool is_face = false>
Tensor<1, n_components_, Tensor<2, dim, VectorizedArray<Number> > > FEEvaluationBase< dim, n_components_, Number, is_face >::get_hessian ( const unsigned int  q_point) const

Return the Hessian of a finite element function at quadrature point number q_point after a call to evaluate(...,true). If only the diagonal or even the trace of the Hessian, the Laplacian, is needed, use the other functions below.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
gradient_type FEEvaluationBase< dim, n_components_, Number, is_face >::get_hessian_diagonal ( const unsigned int  q_point) const

Return the diagonal of the Hessian of a finite element function at quadrature point number q_point after a call to evaluate(...,true).

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
value_type FEEvaluationBase< dim, n_components_, Number, is_face >::get_laplacian ( const unsigned int  q_point) const

Return the Laplacian (i.e., the trace of the Hessian) of a finite element function at quadrature point number q_point after a call to evaluate(...,true). Compared to the case when computing the full Hessian, some operations can be saved when only the Laplacian is requested.

Note that the derived class FEEvaluationAccess overloads this operation with specializations for the scalar case (n_components == 1) and for the vector-valued case (n_components == dim).

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number> FEEvaluationBase< dim, n_components_, Number, is_face >::get_divergence ( const unsigned int  q_point) const

Return the divergence of a vector-valued finite element at quadrature point number q_point after a call to evaluate(...,true,...).

Note
Only available for n_components_==dim.
template<int dim, int n_components_, typename Number, bool is_face = false>
SymmetricTensor<2, dim, VectorizedArray<Number> > FEEvaluationBase< dim, n_components_, Number, is_face >::get_symmetric_gradient ( const unsigned int  q_point) const

Return the symmetric gradient of a vector-valued finite element at quadrature point number q_point after a call to evaluate(...,true,...). It corresponds to 0.5 (grad+gradT).

Note
Only available for n_components_==dim.
template<int dim, int n_components_, typename Number, bool is_face = false>
Tensor<1, (dim == 2 ? 1 : dim), VectorizedArray<Number> > FEEvaluationBase< dim, n_components_, Number, is_face >::get_curl ( const unsigned int  q_point) const

Return the curl of the vector field, \(\nabla \times v\) after a call to evaluate(...,true,...).

Note
Only available for n_components_==dim.
template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::submit_divergence ( const VectorizedArray< Number >  div_in,
const unsigned int  q_point 
)

Write a contribution that is tested by the divergence to the field containing the values on quadrature points with component q_point. Access to the same field as through get_gradient. If applied before the function integrate(...,true) is called, this specifies what is tested by all basis function gradients on the current cell and integrated over.

Note
Only available for n_components_==dim.
This operation writes the data to the same field as submit_gradient(). As a consequence, only one of these two can be used. Usually, the contribution of a potential call to this function must be added into the diagonal of the contribution for submit_gradient().
template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::submit_symmetric_gradient ( const SymmetricTensor< 2, dim, VectorizedArray< Number >>  grad_in,
const unsigned int  q_point 
)

Write a contribution that is tested by the symmetric gradient to the field containing the values on quadrature points with component q_point. Access to the same field as through get_symmetric_gradient. If applied before the function integrate(...,true) is called, this specifies the symmetric gradient which is tested by all basis function symmetric gradients on the current cell and integrated over.

Note
Only available for n_components_==dim.
This operation writes the data to the same field as submit_gradient(). As a consequence, only one of these two can be used. Usually, the contribution of a potential call to this function must be added to the respective entries of the rank-2 tensor for submit_gradient().
template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::submit_curl ( const Tensor< 1, dim==2?1:dim, VectorizedArray< Number >>  curl_in,
const unsigned int  q_point 
)

Write the components of a curl containing the values on quadrature point q_point. Access to the same data field as through get_gradient.

Note
Only available for n_components_==dim.
This operation writes the data to the same field as submit_gradient(). As a consequence, only one of these two can be used. Usually, the contribution of a potential call to this function must be added to the respective entries of the rank-2 tensor for submit_gradient().
template<int dim, int n_components_, typename Number, bool is_face = false>
value_type FEEvaluationBase< dim, n_components_, Number, is_face >::integrate_value ( ) const

Takes values on quadrature points, multiplies by the Jacobian determinant and quadrature weights (JxW) and sums the values for all quadrature points on the cell. The result is a scalar, representing the integral over the function over the cell. If a vector-element is used, the resulting components are still separated. Moreover, if vectorization is enabled, the integral values of several cells are contained in the slots of the returned VectorizedArray field.

Note
In case the FEEvaluation object is initialized with a batch of cells where not all lanes in the SIMD vector VectorizedArray are representing actual data, this method performs computations on dummy data (that is copied from the last valid lane) that will not make sense. Thus, the user needs to make sure that it is not used in any computation explicitly, like when summing the results of several cells.
template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number> FEEvaluationBase< dim, n_components_, Number, is_face >::JxW ( const unsigned int  q_index) const

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

template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::fill_JxW_values ( AlignedVector< VectorizedArray< Number >> &  JxW_values) const

Fills the JxW values currently used into the given array.

template<int dim, int n_components_, typename Number, bool is_face = false>
Tensor<2, dim, VectorizedArray<Number> > FEEvaluationBase< dim, n_components_, Number, is_face >::inverse_jacobian ( const unsigned int  q_index) const

Return the inverse and transposed version of Jacobian of the mapping between the unit to the real cell (representing the covariant transformation). This is exactly the matrix used internally to transform the unit cell gradients to gradients on the real cell.

template<int dim, int n_components_, typename Number, bool is_face = false>
Tensor<1, dim, VectorizedArray<Number> > FEEvaluationBase< dim, n_components_, Number, is_face >::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.
template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number> FEEvaluationBase< dim, n_components_, Number, is_face >::read_cell_data ( const AlignedVector< VectorizedArray< Number >> &  array) const

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

template<int dim, int n_components_, typename Number, bool is_face = false>
const VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::begin_dof_values ( ) const

Return a read-only pointer to the first field of the dof values. This is the data field the read_dof_values() functions write into. First come the dof values for the first component, then all values for the second component, and so on. This is related to the internal data structures used in this class. In general, it is safer to use the get_dof_value() function instead.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::begin_dof_values ( )

Return a read and write pointer to the first field of the dof values. This is the data field the read_dof_values() functions write into. First come the dof values for the first component, then all values for the second component, and so on. This is related to the internal data structures used in this class. In general, it is safer to use the get_dof_value() function instead.

template<int dim, int n_components_, typename Number, bool is_face = false>
const VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::begin_values ( ) const

Return a read-only pointer to the first field of function values on quadrature points. First come the function values on all quadrature points for the first component, then all values for the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_value() function instead, which does all the transformation internally.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::begin_values ( )

Return a read and write pointer to the first field of function values on quadrature points. First come the function values on all quadrature points for the first component, then all values for the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_value() function instead, which does all the transformation internally.

template<int dim, int n_components_, typename Number, bool is_face = false>
const VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::begin_gradients ( ) const

Return a read-only pointer to the first field of function gradients on quadrature points. First comes the x-component of the gradient for the first component on all quadrature points, then the y-component, and so on. Next comes the x-component of the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_gradient() function instead, which does all the transformation internally.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::begin_gradients ( )

Return a read and write pointer to the first field of function gradients on quadrature points. First comes the x-component of the gradient for the first component on all quadrature points, then the y-component, and so on. Next comes the x-component of the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_gradient() function instead, which does all the transformation internally.

template<int dim, int n_components_, typename Number, bool is_face = false>
const VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::begin_hessians ( ) const

Return a read-only pointer to the first field of function hessians on quadrature points. First comes the xx-component of the hessian for the first component on all quadrature points, then the yy-component, zz- component in (3D), then the xy-component, and so on. Next comes the xx- component of the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_laplacian() or get_hessian() functions instead, which does all the transformation internally.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::begin_hessians ( )

Return a read and write pointer to the first field of function hessians on quadrature points. First comes the xx-component of the hessian for the first component on all quadrature points, then the yy-component, zz- component in (3D), then the xy-component, and so on. Next comes the xx- component of the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_laplacian() or get_hessian() functions instead, which does all the transformation internally.

template<int dim, int n_components_, typename Number, bool is_face = false>
const std::vector<unsigned int>& FEEvaluationBase< dim, n_components_, Number, is_face >::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.

template<int dim, int n_components_, typename Number, bool is_face = false>
ArrayView<VectorizedArray<Number> > FEEvaluationBase< dim, n_components_, Number, is_face >::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.

template<int dim, int n_components_, typename Number, bool is_face = false>
FEEvaluationBase& FEEvaluationBase< dim, n_components_, Number, is_face >::operator= ( const FEEvaluationBase< dim, n_components_, Number, is_face > &  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.

template<int dim, int n_components_, typename Number, bool is_face = false>
template<typename VectorType , typename VectorOperation >
void FEEvaluationBase< dim, n_components_, Number, is_face >::read_write_operation ( const VectorOperation operation,
VectorType *  vectors[],
const bool  apply_constraints = true 
) const
protected

A unified function to read from and write into vectors based on the given template operation. It can perform the operation for read_dof_values, distribute_local_to_global, and set_dof_values. It performs the operation for several vectors at a time.

template<int dim, int n_components_, typename Number, bool is_face = false>
template<typename VectorType , typename VectorOperation >
void FEEvaluationBase< dim, n_components_, Number, is_face >::read_write_operation_contiguous ( const VectorOperation operation,
VectorType *  vectors[] 
) const
protected

A unified function to read from and write into vectors based on the given template operation for DG-type schemes where all degrees of freedom on cells are contiguous. It can perform the operation for read_dof_values(), distribute_local_to_global(), and set_dof_values() for several vectors at a time, depending on n_components.

template<int dim, int n_components_, typename Number, bool is_face = false>
template<typename VectorType , typename VectorOperation >
void FEEvaluationBase< dim, n_components_, Number, is_face >::read_write_operation_global ( const VectorOperation operation,
VectorType *  vectors[] 
) const
protected

A unified function to read from and write into vectors based on the given template operation for the case when we do not have an underlying MatrixFree object. It can perform the operation for read_dof_values, distribute_local_to_global, and set_dof_values. It performs the operation for several vectors at a time, depending on n_components.

template<int dim, int n_components_, typename Number, bool is_face = false>
void FEEvaluationBase< dim, n_components_, Number, is_face >::set_data_pointers ( )
private

Sets the pointers for values, gradients, hessians to the central scratch_data_array.

Friends And Related Function Documentation

template<int dim, int n_components_, typename Number, bool is_face = false>
template<int , int , typename , bool >
friend class FEEvaluationBase
friend

Make other FEEvaluationBase as well as FEEvaluation objects friends.

Definition at line 1091 of file fe_evaluation.h.

Member Data Documentation

template<int dim, int n_components_, typename Number, bool is_face = false>
AlignedVector<VectorizedArray<Number> >* FEEvaluationBase< dim, n_components_, Number, is_face >::scratch_data_array
protected

This is the general array for all data fields.

Definition at line 820 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::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 827 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::values_dofs[n_components]
protected

This field stores the values for local degrees of freedom (e.g. after reading out from a vector but before applying unit cell transformations or before distributing them into a result vector). The methods get_dof_value() and submit_dof_value() read from or write to this field.

The values of this array are stored in the start section of scratch_data_array. Due to its access as a thread local memory, the memory can get reused between different calls. As opposed to requesting memory on the stack, this approach allows for very large polynomial degrees.

Definition at line 841 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::values_quad[n_components]
protected

This field stores the values of the finite element function on quadrature points after applying unit cell transformations or before integrating. The methods get_value() and submit_value() access this field.

The values of this array are stored in the start section of scratch_data_array. Due to its access as a thread local memory, the memory can get reused between different calls. As opposed to requesting memory on the stack, this approach allows for very large polynomial degrees.

Definition at line 854 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::gradients_quad[n_components][dim]
protected

This field stores the gradients of the finite element function on quadrature points after applying unit cell transformations or before integrating. The methods get_gradient() and submit_gradient() (as well as some specializations like get_symmetric_gradient() or get_divergence()) access this field.

The values of this array are stored in the start section of scratch_data_array. Due to its access as a thread local memory, the memory can get reused between different calls. As opposed to requesting memory on the stack, this approach allows for very large polynomial degrees.

Definition at line 869 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::hessians_quad[n_components][(dim *(dim+1))/2]
protected

This field stores the Hessians of the finite element function on quadrature points after applying unit cell transformations. The methods get_hessian(), get_laplacian(), get_hessian_diagonal() access this field.

The values of this array are stored in the start section of scratch_data_array. Due to its access as a thread local memory, the memory can get reused between different calls. As opposed to requesting memory on the stack, this approach allows for very large polynomial degrees.

Definition at line 882 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::quad_no
protected

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

Definition at line 887 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::n_fe_components
protected

Stores the number of components in the finite element as detected in the MatrixFree storage class for comparison with the template argument.

Definition at line 893 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::active_fe_index
protected

Stores the active fe index for this class for efficient indexing in the hp case.

Definition at line 899 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::active_quad_index
protected

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

Definition at line 905 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::n_quadrature_points
protected

Stores the number of quadrature points in the current evaluation context.

Definition at line 910 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const MatrixFree<dim, Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::matrix_info
protected

Stores a pointer to the underlying data.

Definition at line 915 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const internal::MatrixFreeFunctions::DoFInfo* FEEvaluationBase< dim, n_components_, Number, is_face >::dof_info
protected

Stores 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 922 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const internal::MatrixFreeFunctions:: MappingInfoStorage<(is_face ? dim - 1 : dim), dim, Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::mapping_data
protected

Stores 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 931 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const internal::MatrixFreeFunctions::ShapeInfo<VectorizedArray<Number> >* FEEvaluationBase< dim, n_components_, Number, is_face >::data
protected

Stores 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 939 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const Tensor<2, dim, VectorizedArray<Number> >* FEEvaluationBase< dim, n_components_, Number, is_face >::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 945 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const VectorizedArray<Number>* FEEvaluationBase< dim, n_components_, Number, is_face >::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 953 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const Tensor<1, dim, VectorizedArray<Number> >* FEEvaluationBase< dim, n_components_, Number, is_face >::normal_vectors
protected

A pointer to the normal vectors at faces.

Definition at line 958 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const Tensor<1, dim, VectorizedArray<Number> >* FEEvaluationBase< dim, n_components_, Number, is_face >::normal_x_jacobian
protected

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

Definition at line 963 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const Number* FEEvaluationBase< dim, n_components_, Number, is_face >::quadrature_weights
protected

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

Definition at line 968 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::cell
protected

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

Definition at line 974 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
bool FEEvaluationBase< dim, n_components_, Number, is_face >::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 980 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex FEEvaluationBase< dim, n_components_, Number, is_face >::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 986 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::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 992 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::face_orientation
protected

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

Definition at line 998 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::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 1007 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
internal::MatrixFreeFunctions::GeometryType FEEvaluationBase< dim, n_components_, Number, is_face >::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 1015 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
bool FEEvaluationBase< dim, n_components_, Number, is_face >::dof_values_initialized
protected

Debug information to track whether dof values have been initialized before accessed. Used to control exceptions when uninitialized data is used.

Definition at line 1022 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
bool FEEvaluationBase< dim, n_components_, Number, is_face >::values_quad_initialized
protected

Debug information to track whether values on quadrature points have been initialized before accessed. Used to control exceptions when uninitialized data is used.

Definition at line 1029 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
bool FEEvaluationBase< dim, n_components_, Number, is_face >::gradients_quad_initialized
protected

Debug information to track whether gradients on quadrature points have been initialized before accessed. Used to control exceptions when uninitialized data is used.

Definition at line 1036 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
bool FEEvaluationBase< dim, n_components_, Number, is_face >::hessians_quad_initialized
protected

Debug information to track whether Hessians on quadrature points have been initialized before accessed. Used to control exceptions when uninitialized data is used.

Definition at line 1043 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
bool FEEvaluationBase< dim, n_components_, Number, is_face >::values_quad_submitted
protected

Debug information to track whether values on quadrature points have been submitted for integration before the integration is actually stared. Used to control exceptions when uninitialized data is used.

Definition at line 1050 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
bool FEEvaluationBase< dim, n_components_, Number, is_face >::gradients_quad_submitted
protected

Debug information to track whether gradients on quadrature points have been submitted for integration before the integration is actually stared. Used to control exceptions when uninitialized data is used.

Definition at line 1057 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
std::shared_ptr< internal::MatrixFreeFunctions::MappingDataOnTheFly<dim, Number> > FEEvaluationBase< dim, n_components_, Number, is_face >::mapped_geometry
protected

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

Definition at line 1065 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
const unsigned int FEEvaluationBase< dim, n_components_, Number, is_face >::first_selected_component
protected

For a FiniteElement with more than one base element, select at which component this data structure should start.

Definition at line 1071 of file fe_evaluation.h.

template<int dim, int n_components_, typename Number, bool is_face = false>
std::vector<types::global_dof_index> FEEvaluationBase< dim, n_components_, Number, is_face >::local_dof_indices
mutableprotected

A temporary data structure necessary to read degrees of freedom when no MatrixFree object was given at initialization.

Definition at line 1077 of file fe_evaluation.h.


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