Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Classes | Public Member Functions | Public Attributes | Static Public Attributes | List of all members
FEValuesBase< dim, spacedim > Class Template Reference

#include <deal.II/fe/fe.h>

Inheritance diagram for FEValuesBase< dim, spacedim >:
[legend]

Classes

class  CellIterator
 
class  CellIteratorBase
 
class  TriaCellIterator
 

Public Member Functions

 FEValuesBase (const unsigned int n_q_points, const unsigned int dofs_per_cell, const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe)
 
FEValuesBaseoperator= (const FEValuesBase &)=delete
 
 FEValuesBase (const FEValuesBase &)=delete
 
virtual ~FEValuesBase () override
 
Access to shape function values

These fields are filled by the finite element.

const doubleshape_value (const unsigned int function_no, const unsigned int point_no) const
 
double shape_value_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
 
const Tensor< 1, spacedim > & shape_grad (const unsigned int function_no, const unsigned int quadrature_point) const
 
Tensor< 1, spacedim > shape_grad_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
 
const Tensor< 2, spacedim > & shape_hessian (const unsigned int function_no, const unsigned int point_no) const
 
Tensor< 2, spacedim > shape_hessian_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
 
const Tensor< 3, spacedim > & shape_3rd_derivative (const unsigned int function_no, const unsigned int point_no) const
 
Tensor< 3, spacedim > shape_3rd_derivative_component (const unsigned int function_no, const unsigned int point_no, const unsigned int component) const
 
Access to values of global finite element fields
template<class InputVector >
void get_function_values (const InputVector &fe_function, std::vector< typename InputVector::value_type > &values) const
 
template<class InputVector >
void get_function_values (const InputVector &fe_function, std::vector< Vector< typename InputVector::value_type >> &values) const
 
template<class InputVector >
void get_function_values (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< typename InputVector::value_type > &values) const
 
template<class InputVector >
void get_function_values (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Vector< typename InputVector::value_type >> &values) const
 
template<class InputVector >
void get_function_values (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< typename InputVector::value_type >> values, const bool quadrature_points_fastest) const
 
Access to derivatives of global finite element fields
template<class InputVector >
void get_function_gradients (const InputVector &fe_function, std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &gradients) const
 
template<class InputVector >
void get_function_gradients (const InputVector &fe_function, std::vector< std::vector< Tensor< 1, spacedim, typename InputVector::value_type >>> &gradients) const
 
template<class InputVector >
void get_function_gradients (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &gradients) const
 
template<class InputVector >
void get_function_gradients (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< Tensor< 1, spacedim, typename InputVector::value_type >>> gradients, bool quadrature_points_fastest=false) const
 
Access to second derivatives

Hessian matrices and Laplacians of global finite element fields

template<class InputVector >
void get_function_hessians (const InputVector &fe_function, std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &hessians) const
 
template<class InputVector >
void get_function_hessians (const InputVector &fe_function, std::vector< std::vector< Tensor< 2, spacedim, typename InputVector::value_type >>> &hessians, bool quadrature_points_fastest=false) const
 
template<class InputVector >
void get_function_hessians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &hessians) const
 
template<class InputVector >
void get_function_hessians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< Tensor< 2, spacedim, typename InputVector::value_type >>> hessians, bool quadrature_points_fastest=false) const
 
template<class InputVector >
void get_function_laplacians (const InputVector &fe_function, std::vector< typename InputVector::value_type > &laplacians) const
 
template<class InputVector >
void get_function_laplacians (const InputVector &fe_function, std::vector< Vector< typename InputVector::value_type >> &laplacians) const
 
template<class InputVector >
void get_function_laplacians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< typename InputVector::value_type > &laplacians) const
 
template<class InputVector >
void get_function_laplacians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Vector< typename InputVector::value_type >> &laplacians) const
 
template<class InputVector >
void get_function_laplacians (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< std::vector< typename InputVector::value_type >> &laplacians, bool quadrature_points_fastest=false) const
 
Access to third derivatives of global finite element fields
template<class InputVector >
void get_function_third_derivatives (const InputVector &fe_function, std::vector< Tensor< 3, spacedim, typename InputVector::value_type >> &third_derivatives) const
 
template<class InputVector >
void get_function_third_derivatives (const InputVector &fe_function, std::vector< std::vector< Tensor< 3, spacedim, typename InputVector::value_type >>> &third_derivatives, bool quadrature_points_fastest=false) const
 
template<class InputVector >
void get_function_third_derivatives (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Tensor< 3, spacedim, typename InputVector::value_type >> &third_derivatives) const
 
template<class InputVector >
void get_function_third_derivatives (const InputVector &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< Tensor< 3, spacedim, typename InputVector::value_type >>> third_derivatives, bool quadrature_points_fastest=false) const
 
Cell degrees of freedom
std_cxx20::ranges::iota_view< unsigned int, unsigned intdof_indices () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intdof_indices_starting_at (const unsigned int start_dof_index) const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intdof_indices_ending_at (const unsigned int end_dof_index) const
 
Geometry of the cell
std_cxx20::ranges::iota_view< unsigned int, unsigned intquadrature_point_indices () const
 
const Point< spacedim > & quadrature_point (const unsigned int q) const
 
const std::vector< Point< spacedim > > & get_quadrature_points () const
 
double JxW (const unsigned int quadrature_point) const
 
const std::vector< double > & get_JxW_values () const
 
const DerivativeForm< 1, dim, spacedim > & jacobian (const unsigned int quadrature_point) const
 
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians () const
 
const DerivativeForm< 2, dim, spacedim > & jacobian_grad (const unsigned int quadrature_point) const
 
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads () const
 
const Tensor< 3, spacedim > & jacobian_pushed_forward_grad (const unsigned int quadrature_point) const
 
const std::vector< Tensor< 3, spacedim > > & get_jacobian_pushed_forward_grads () const
 
const DerivativeForm< 3, dim, spacedim > & jacobian_2nd_derivative (const unsigned int quadrature_point) const
 
const std::vector< DerivativeForm< 3, dim, spacedim > > & get_jacobian_2nd_derivatives () const
 
const Tensor< 4, spacedim > & jacobian_pushed_forward_2nd_derivative (const unsigned int quadrature_point) const
 
const std::vector< Tensor< 4, spacedim > > & get_jacobian_pushed_forward_2nd_derivatives () const
 
const DerivativeForm< 4, dim, spacedim > & jacobian_3rd_derivative (const unsigned int quadrature_point) const
 
const std::vector< DerivativeForm< 4, dim, spacedim > > & get_jacobian_3rd_derivatives () const
 
const Tensor< 5, spacedim > & jacobian_pushed_forward_3rd_derivative (const unsigned int quadrature_point) const
 
const std::vector< Tensor< 5, spacedim > > & get_jacobian_pushed_forward_3rd_derivatives () const
 
const DerivativeForm< 1, spacedim, dim > & inverse_jacobian (const unsigned int quadrature_point) const
 
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians () const
 
const Tensor< 1, spacedim > & normal_vector (const unsigned int i) const
 
const std::vector< Tensor< 1, spacedim > > & get_normal_vectors () const
 
Extractors Methods to extract individual components
const FEValuesViews::Scalar< dim, spacedim > & operator[] (const FEValuesExtractors::Scalar &scalar) const
 
const FEValuesViews::Vector< dim, spacedim > & operator[] (const FEValuesExtractors::Vector &vector) const
 
const FEValuesViews::SymmetricTensor< 2, dim, spacedim > & operator[] (const FEValuesExtractors::SymmetricTensor< 2 > &tensor) const
 
const FEValuesViews::Tensor< 2, dim, spacedim > & operator[] (const FEValuesExtractors::Tensor< 2 > &tensor) const
 

Public Attributes

const unsigned int n_quadrature_points
 
const unsigned int dofs_per_cell
 

Static Public Attributes

static const unsigned int dimension = dim
 
static const unsigned int space_dimension = spacedim
 

Access to the raw data

std::unique_ptr< const CellIteratorBasepresent_cell
 
boost::signals2::connection tria_listener_refinement
 
boost::signals2::connection tria_listener_mesh_transform
 
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
 
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > mapping_data
 
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_output
 
const SmartPointer< const FiniteElement< dim, spacedim >, FEValuesBase< dim, spacedim > > fe
 
std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > fe_data
 
::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > finite_element_output
 
UpdateFlags update_flags
 
CellSimilarity::Similarity cell_similarity
 
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
 
template<int , int >
class FEValuesViews::Scalar
 
template<int , int >
class FEValuesViews::Vector
 
template<int , int , int >
class FEValuesViews::SymmetricTensor
 
template<int , int , int >
class FEValuesViews::Tensor
 
const Mapping< dim, spacedim > & get_mapping () const
 
const FiniteElement< dim, spacedim > & get_fe () const
 
UpdateFlags get_update_flags () const
 
const Triangulation< dim, spacedim >::cell_iterator get_cell () const
 
CellSimilarity::Similarity get_cell_similarity () const
 
std::size_t memory_consumption () const
 
static ::ExceptionBaseExcAccessToUninitializedField (std::string arg1)
 
static ::ExceptionBaseExcFEDontMatch ()
 
static ::ExceptionBaseExcShapeFunctionNotPrimitive (int arg1)
 
static ::ExceptionBaseExcFENotPrimitive ()
 
void invalidate_present_cell ()
 
void maybe_invalidate_previous_present_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell)
 
UpdateFlags compute_update_flags (const UpdateFlags update_flags) const
 
void check_cell_similarity (const typename Triangulation< dim, spacedim >::cell_iterator &cell)
 

Detailed Description

template<int dim, int spacedim>
class FEValuesBase< dim, spacedim >

FEValues, FEFaceValues and FESubfaceValues objects are interfaces to finite element and mapping classes on the one hand side, to cells and quadrature rules on the other side. They allow to evaluate values or derivatives of shape functions at the quadrature points of a quadrature formula when projected by a mapping from the unit cell onto a cell in real space. The reason for this abstraction is possible optimization: Depending on the type of finite element and mapping, some values can be computed once on the unit cell. Others must be computed on each cell, but maybe computation of several values at the same time offers ways for optimization. Since this interplay may be complex and depends on the actual finite element, it cannot be left to the applications programmer.

FEValues, FEFaceValues and FESubfaceValues provide only data handling: computations are left to objects of type Mapping and FiniteElement. These provide functions get_*_data and fill_*_values which are called by the constructor and reinit functions of FEValues*, respectively.

General usage

Usually, an object of FEValues* is used in integration loops over all cells of a triangulation (or faces of cells). To take full advantage of the optimization features, it should be constructed before the loop so that information that does not depend on the location and shape of cells can be computed once and for all (this includes, for example, the values of shape functions at quadrature points for the most common elements: we can evaluate them on the unit cell and they will be the same when mapped to the real cell). Then, in the loop over all cells, it must be re-initialized for each grid cell to compute that part of the information that changes depending on the actual cell (for example, the gradient of shape functions equals the gradient on the unit cell – which can be computed once and for all – times the Jacobian matrix of the mapping between unit and real cell, which needs to be recomputed for each cell).

A typical piece of code, adding up local contributions to the Laplace matrix looks like this:

FEValues values (mapping, finite_element, quadrature, flags);
for (const auto &cell : dof_handler.active_cell_iterators())
{
values.reinit(cell);
for (unsigned int q=0; q<quadrature.size(); ++q)
for (unsigned int i=0; i<finite_element.dofs_per_cell; ++i)
for (unsigned int j=0; j<finite_element.dofs_per_cell; ++j)
A(i,j) += fe_values.shape_value(i,q) *
fe_values.shape_value(j,q) *
fe_values.JxW(q);
...
}

The individual functions used here are described below. Note that by design, the order of quadrature points used inside the FEValues object is the same as defined by the quadrature formula passed to the constructor of the FEValues object above.

Member functions

The functions of this class fall into different categories:

Internals about the implementation

The mechanisms by which this class work are discussed on the page on Update flags and about the How Mapping, FiniteElement, and FEValues work together.

Author
Wolfgang Bangerth, 1998, 2003, Guido Kanschat, 2001

Definition at line 36 of file fe.h.

Constructor & Destructor Documentation

◆ FEValuesBase() [1/2]

template<int dim, int spacedim>
FEValuesBase< dim, spacedim >::FEValuesBase ( const unsigned int  n_q_points,
const unsigned int  dofs_per_cell,
const UpdateFlags  update_flags,
const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe 
)

Constructor. Set up the array sizes with n_q_points quadrature points, dofs_per_cell trial functions per cell and with the given pattern to update the fields when the reinit function of the derived classes is called. The fields themselves are not set up, this must happen in the constructor of the derived class.

Definition at line 3076 of file fe_values.cc.

◆ FEValuesBase() [2/2]

template<int dim, int spacedim>
FEValuesBase< dim, spacedim >::FEValuesBase ( const FEValuesBase< dim, spacedim > &  )
delete

The copy constructor is deleted since objects of this class are not copyable.

◆ ~FEValuesBase()

template<int dim, int spacedim>
FEValuesBase< dim, spacedim >::~FEValuesBase
overridevirtual

Destructor.

Definition at line 3099 of file fe_values.cc.

Member Function Documentation

◆ operator=()

template<int dim, int spacedim>
FEValuesBase& FEValuesBase< dim, spacedim >::operator= ( const FEValuesBase< dim, spacedim > &  )
delete

The copy assignment is deleted since objects of this class are not copyable.

◆ shape_value()

template<int dim, int spacedim>
const double& FEValuesBase< dim, spacedim >::shape_value ( const unsigned int  function_no,
const unsigned int  point_no 
) const

Value of a shape function at a quadrature point on the cell, face or subface selected the last time the reinit function of the derived class was called.

If the shape function is vector-valued, then this returns the only non- zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_value_component() function.

Parameters
function_noNumber of the shape function to be evaluated. Note that this number runs from zero to dofs_per_cell, even in the case of an FEFaceValues or FESubfaceValues object.
point_noNumber of the quadrature point at which function is to be evaluated
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ shape_value_component()

template<int dim, int spacedim>
double FEValuesBase< dim, spacedim >::shape_value_component ( const unsigned int  function_no,
const unsigned int  point_no,
const unsigned int  component 
) const

Compute one vector component of the value of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_value() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_value() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.

Parameters
function_noNumber of the shape function to be evaluated.
point_noNumber of the quadrature point at which function is to be evaluated.
componentvector component to be evaluated.
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ shape_grad()

template<int dim, int spacedim>
const Tensor<1, spacedim>& FEValuesBase< dim, spacedim >::shape_grad ( const unsigned int  function_no,
const unsigned int  quadrature_point 
) const

Compute the gradient of the function_noth shape function at the quadrature_pointth quadrature point with respect to real cell coordinates. If you want to get the derivative in one of the coordinate directions, use the appropriate function of the Tensor class to extract one component of the Tensor returned by this function. Since only a reference to the gradient's value is returned, there should be no major performance drawback.

If the shape function is vector-valued, then this returns the only non- zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then it will throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_grad_component() function.

The same holds for the arguments of this function as for the shape_value() function.

Parameters
function_noNumber of the shape function to be evaluated.
quadrature_pointNumber of the quadrature point at which function is to be evaluated.
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ shape_grad_component()

template<int dim, int spacedim>
Tensor<1, spacedim> FEValuesBase< dim, spacedim >::shape_grad_component ( const unsigned int  function_no,
const unsigned int  point_no,
const unsigned int  component 
) const

Return one vector component of the gradient of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_grad() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_grad() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.

The same holds for the arguments of this function as for the shape_value_component() function.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ shape_hessian()

template<int dim, int spacedim>
const Tensor<2, spacedim>& FEValuesBase< dim, spacedim >::shape_hessian ( const unsigned int  function_no,
const unsigned int  point_no 
) const

Second derivatives of the function_noth shape function at the point_noth quadrature point with respect to real cell coordinates. If you want to get the derivatives in one of the coordinate directions, use the appropriate function of the Tensor class to extract one component. Since only a reference to the hessian values is returned, there should be no major performance drawback.

If the shape function is vector-valued, then this returns the only non- zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_hessian_component() function.

The same holds for the arguments of this function as for the shape_value() function.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ shape_hessian_component()

template<int dim, int spacedim>
Tensor<2, spacedim> FEValuesBase< dim, spacedim >::shape_hessian_component ( const unsigned int  function_no,
const unsigned int  point_no,
const unsigned int  component 
) const

Return one vector component of the hessian of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_hessian() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_hessian() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.

The same holds for the arguments of this function as for the shape_value_component() function.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ shape_3rd_derivative()

template<int dim, int spacedim>
const Tensor<3, spacedim>& FEValuesBase< dim, spacedim >::shape_3rd_derivative ( const unsigned int  function_no,
const unsigned int  point_no 
) const

Third derivatives of the function_noth shape function at the point_noth quadrature point with respect to real cell coordinates. If you want to get the 3rd derivatives in one of the coordinate directions, use the appropriate function of the Tensor class to extract one component. Since only a reference to the 3rd derivative values is returned, there should be no major performance drawback.

If the shape function is vector-valued, then this returns the only non- zero component. If the shape function has more than one non-zero component (i.e. it is not primitive), then throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_3rdderivative_component() function.

The same holds for the arguments of this function as for the shape_value() function.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ shape_3rd_derivative_component()

template<int dim, int spacedim>
Tensor<3, spacedim> FEValuesBase< dim, spacedim >::shape_3rd_derivative_component ( const unsigned int  function_no,
const unsigned int  point_no,
const unsigned int  component 
) const

Return one vector component of the third derivative of a shape function at a quadrature point. If the finite element is scalar, then only component zero is allowed and the return value equals that of the shape_3rdderivative() function. If the finite element is vector valued but all shape functions are primitive (i.e. they are non-zero in only one component), then the value returned by shape_3rdderivative() equals that of this function for exactly one component. This function is therefore only of greater interest if the shape function is not primitive, but then it is necessary since the other function cannot be used.

The same holds for the arguments of this function as for the shape_value_component() function.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_function_values() [1/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
std::vector< typename InputVector::value_type > &  values 
) const

Return the values of a finite element function restricted to the current cell, face or subface selected the last time the reinit function of the derived class was called, at the quadrature points.

If the present cell is not active then values are interpolated to the current cell and point values are computed from that.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. To get values of multi- component elements, there is another get_function_values() below, returning a vector of vectors of results.

Parameters
[in]fe_functionA vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell.
[out]valuesThe values of the function specified by fe_function at the quadrature points of the current cell. The object is assume to already have the correct size. The data type stored by this output vector must be what you get when you multiply the values of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument). This happens to be equal to the type of the elements of the solution vector.
Postcondition
values[q] will contain the value of the field described by fe_function at the \(q\)th quadrature point.
Note
The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DoFHandler object with which this FEValues object was last initialized.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3572 of file fe_values.cc.

◆ get_function_values() [2/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
std::vector< Vector< typename InputVector::value_type >> &  values 
) const

This function does the same as the other get_function_values(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.

Postcondition
values[q] is a vector of values of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by values[q] equals the number of components of the finite element, i.e. values[q](c) returns the value of the \(c\)th vector component at the \(q\)th quadrature point.
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3621 of file fe_values.cc.

◆ get_function_values() [3/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
std::vector< typename InputVector::value_type > &  values 
) const

Generate function values from an arbitrary vector.

This function offers the possibility to extract function values in quadrature points from vectors not corresponding to a whole discretization.

The vector indices corresponds to the degrees of freedom on a single cell. Its length may even be a multiple of the number of dofs per cell. Then, the vectors in value should allow for the same multiple of the components of the finite element.

You may want to use this function, if you want to access just a single block from a BlockVector, if you have a multi-level vector or if you already have a local representation of your finite element data.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3597 of file fe_values.cc.

◆ get_function_values() [4/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
std::vector< Vector< typename InputVector::value_type >> &  values 
) const

Generate vector function values from an arbitrary vector.

This function offers the possibility to extract function values in quadrature points from vectors not corresponding to a whole discretization.

The vector indices corresponds to the degrees of freedom on a single cell. Its length may even be a multiple of the number of dofs per cell. Then, the vectors in value should allow for the same multiple of the components of the finite element.

You may want to use this function, if you want to access just a single block from a BlockVector, if you have a multi-level vector or if you already have a local representation of your finite element data.

Since this function allows for fairly general combinations of argument sizes, be aware that the checks on the arguments may not detect errors.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3649 of file fe_values.cc.

◆ get_function_values() [5/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_values ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
ArrayView< std::vector< typename InputVector::value_type >>  values,
const bool  quadrature_points_fastest 
) const

Generate vector function values from an arbitrary vector.

This function offers the possibility to extract function values in quadrature points from vectors not corresponding to a whole discretization.

The vector indices corresponds to the degrees of freedom on a single cell. Its length may even be a multiple of the number of dofs per cell. Then, the vectors in value should allow for the same multiple of the components of the finite element.

Depending on the value of the last argument, the outer vector of values has either the length of the quadrature rule (quadrature_points_fastest == false) or the length of components to be filled quadrature_points_fastest == true. If p is the current quadrature point number and i is the vector component of the solution desired, the access to values is values[p][i] if quadrature_points_fastest == false, and values[i][p] otherwise.

You may want to use this function, if you want to access just a single block from a BlockVector, if you have a multi-level vector or if you already have a local representation of your finite element data.

Since this function allows for fairly general combinations of argument sizes, be aware that the checks on the arguments may not detect errors.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3680 of file fe_values.cc.

◆ get_function_gradients() [1/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &  gradients 
) const

Compute the gradients of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's gradient instead of its value.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.

Parameters
[in]fe_functionA vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell.
[out]gradientsThe gradients of the function specified by fe_function at the quadrature points of the current cell. The gradients are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. The data type stored by this output vector must be what you get when you multiply the gradients of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).
Postcondition
gradients[q] will contain the gradient of the field described by fe_function at the \(q\)th quadrature point. gradients[q][d] represents the derivative in coordinate direction \(d\) at quadrature point \(q\).
Note
The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DoFHandler object with which this FEValues object was last initialized.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3713 of file fe_values.cc.

◆ get_function_gradients() [2/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
std::vector< std::vector< Tensor< 1, spacedim, typename InputVector::value_type >>> &  gradients 
) const

This function does the same as the other get_function_gradients(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.

Postcondition
gradients[q] is a vector of gradients of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by gradients[q] equals the number of components of the finite element, i.e. gradients[q][c] returns the gradient of the \(c\)th vector component at the \(q\)th quadrature point. Consequently, gradients[q][c][d] is the derivative in coordinate direction \(d\) of the \(c\)th vector component of the vector field at quadrature point \(q\) of the current cell.
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3764 of file fe_values.cc.

◆ get_function_gradients() [3/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
std::vector< Tensor< 1, spacedim, typename InputVector::value_type >> &  gradients 
) const

Function gradient access with more flexibility. See get_function_values() with corresponding arguments.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3739 of file fe_values.cc.

◆ get_function_gradients() [4/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
ArrayView< std::vector< Tensor< 1, spacedim, typename InputVector::value_type >>>  gradients,
bool  quadrature_points_fastest = false 
) const

Function gradient access with more flexibility. See get_function_values() with corresponding arguments.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_gradients flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3793 of file fe_values.cc.

◆ get_function_hessians() [1/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &  hessians 
) const

Compute the tensor of second derivatives of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's second derivatives instead of its value.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.

Parameters
[in]fe_functionA vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell.
[out]hessiansThe Hessians of the function specified by fe_function at the quadrature points of the current cell. The Hessians are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. The data type stored by this output vector must be what you get when you multiply the Hessians of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).
Postcondition
hessians[q] will contain the Hessian of the field described by fe_function at the \(q\)th quadrature point. hessians[q][i][j] represents the \((i,j)\)th component of the matrix of second derivatives at quadrature point \(q\).
Note
The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DoFHandler object with which this FEValues object was last initialized.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3826 of file fe_values.cc.

◆ get_function_hessians() [2/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
std::vector< std::vector< Tensor< 2, spacedim, typename InputVector::value_type >>> &  hessians,
bool  quadrature_points_fastest = false 
) const

This function does the same as the other get_function_hessians(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.

Postcondition
hessians[q] is a vector of Hessians of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by hessians[q] equals the number of components of the finite element, i.e. hessians[q][c] returns the Hessian of the \(c\)th vector component at the \(q\)th quadrature point. Consequently, hessians[q][c][i][j] is the \((i,j)\)th component of the matrix of second derivatives of the \(c\)th vector component of the vector field at quadrature point \(q\) of the current cell.
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3877 of file fe_values.cc.

◆ get_function_hessians() [3/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
std::vector< Tensor< 2, spacedim, typename InputVector::value_type >> &  hessians 
) const

Access to the second derivatives of a function with more flexibility. See get_function_values() with corresponding arguments.

Definition at line 3852 of file fe_values.cc.

◆ get_function_hessians() [4/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_hessians ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
ArrayView< std::vector< Tensor< 2, spacedim, typename InputVector::value_type >>>  hessians,
bool  quadrature_points_fastest = false 
) const

Access to the second derivatives of a function with more flexibility. See get_function_values() with corresponding arguments.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3908 of file fe_values.cc.

◆ get_function_laplacians() [1/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
std::vector< typename InputVector::value_type > &  laplacians 
) const

Compute the (scalar) Laplacian (i.e. the trace of the tensor of second derivatives) of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's second derivatives instead of its value.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.

Parameters
[in]fe_functionA vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell.
[out]laplaciansThe Laplacians of the function specified by fe_function at the quadrature points of the current cell. The Laplacians are computed in real space (as opposed to on the unit cell). The object is assume to already have the correct size. The data type stored by this output vector must be what you get when you multiply the Laplacians of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument). This happens to be equal to the type of the elements of the input vector.
Postcondition
laplacians[q] will contain the Laplacian of the field described by fe_function at the \(q\)th quadrature point.
For each component of the output vector, there holds laplacians[q]=trace(hessians[q]), where hessians would be the output of the get_function_hessians() function.
Note
The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DoFHandler object with which this FEValues object was last initialized.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3939 of file fe_values.cc.

◆ get_function_laplacians() [2/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
std::vector< Vector< typename InputVector::value_type >> &  laplacians 
) const

This function does the same as the other get_function_laplacians(), but applied to multi-component (vector-valued) elements. The meaning of the arguments is as explained there.

Postcondition
laplacians[q] is a vector of Laplacians of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by laplacians[q] equals the number of components of the finite element, i.e. laplacians[q][c] returns the Laplacian of the \(c\)th vector component at the \(q\)th quadrature point.
For each component of the output vector, there holds laplacians[q][c]=trace(hessians[q][c]), where hessians would be the output of the get_function_hessians() function.
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3988 of file fe_values.cc.

◆ get_function_laplacians() [3/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
std::vector< typename InputVector::value_type > &  laplacians 
) const

Access to the second derivatives of a function with more flexibility. See get_function_values() with corresponding arguments.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 3964 of file fe_values.cc.

◆ get_function_laplacians() [4/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
std::vector< Vector< typename InputVector::value_type >> &  laplacians 
) const

Access to the second derivatives of a function with more flexibility. See get_function_values() with corresponding arguments.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 4015 of file fe_values.cc.

◆ get_function_laplacians() [5/5]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
std::vector< std::vector< typename InputVector::value_type >> &  laplacians,
bool  quadrature_points_fastest = false 
) const

Access to the second derivatives of a function with more flexibility. See get_function_values() with corresponding arguments.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_hessians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 4046 of file fe_values.cc.

◆ get_function_third_derivatives() [1/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_third_derivatives ( const InputVector &  fe_function,
std::vector< Tensor< 3, spacedim, typename InputVector::value_type >> &  third_derivatives 
) const

Compute the tensor of third derivatives of a finite element at the quadrature points of a cell. This function is the equivalent of the corresponding get_function_values() function (see there for more information) but evaluates the finite element field's third derivatives instead of its value.

This function may only be used if the finite element in use is a scalar one, i.e. has only one vector component. There is a corresponding function of the same name for vector-valued finite elements.

Parameters
[in]fe_functionA vector of values that describes (globally) the finite element function that this function should evaluate at the quadrature points of the current cell.
[out]third_derivativesThe third derivatives of the function specified by fe_function at the quadrature points of the current cell. The third derivatives are computed in real space (as opposed to on the unit cell). The object is assumed to already have the correct size. The data type stored by this output vector must be what you get when you multiply the third derivatives of shape function times the type used to store the values of the unknowns \(U_j\) of your finite element vector \(U\) (represented by the fe_function argument).
Postcondition
third_derivatives[q] will contain the third derivatives of the field described by fe_function at the \(q\)th quadrature point. third_derivatives[q][i][j][k] represents the \((i,j,k)\)th component of the 3rd order tensor of third derivatives at quadrature point \(q\).
Note
The actual data type of the input vector may be either a Vector<T>, BlockVector<T>, or one of the PETSc or Trilinos vector wrapper classes. It represents a global vector of DoF values associated with the DoFHandler object with which this FEValues object was last initialized.
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 4076 of file fe_values.cc.

◆ get_function_third_derivatives() [2/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_third_derivatives ( const InputVector &  fe_function,
std::vector< std::vector< Tensor< 3, spacedim, typename InputVector::value_type >>> &  third_derivatives,
bool  quadrature_points_fastest = false 
) const

This function does the same as the other get_function_third_derivatives(), but applied to multi-component (vector- valued) elements. The meaning of the arguments is as explained there.

Postcondition
third_derivatives[q] is a vector of third derivatives of the field described by fe_function at the \(q\)th quadrature point. The size of the vector accessed by third_derivatives[q] equals the number of components of the finite element, i.e. third_derivatives[q][c] returns the third derivative of the \(c\)th vector component at the \(q\)th quadrature point. Consequently, third_derivatives[q][c][i][j][k] is the \((i,j,k)\)th component of the tensor of third derivatives of the \(c\)th vector component of the vector field at quadrature point \(q\) of the current cell.
Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 4129 of file fe_values.cc.

◆ get_function_third_derivatives() [3/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_third_derivatives ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
std::vector< Tensor< 3, spacedim, typename InputVector::value_type >> &  third_derivatives 
) const

Access to the third derivatives of a function with more flexibility. See get_function_values() with corresponding arguments.

Definition at line 4103 of file fe_values.cc.

◆ get_function_third_derivatives() [4/4]

template<int dim, int spacedim>
template<class InputVector >
void FEValuesBase< dim, spacedim >::get_function_third_derivatives ( const InputVector &  fe_function,
const ArrayView< const types::global_dof_index > &  indices,
ArrayView< std::vector< Tensor< 3, spacedim, typename InputVector::value_type >>>  third_derivatives,
bool  quadrature_points_fastest = false 
) const

Access to the third derivatives of a function with more flexibility. See get_function_values() with corresponding arguments.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 4160 of file fe_values.cc.

◆ dof_indices()

template<int dim, int spacedim>
std_cxx20::ranges::iota_view<unsigned int, unsigned int> FEValuesBase< dim, spacedim >::dof_indices ( ) const

Return an object that can be thought of as an array containing all indices from zero (inclusive) to dofs_per_cell (exclusive). This allows one to write code using range-based for loops of the following kind:

FEValues<dim> fe_values (...);
for (auto &cell : dof_handler.active_cell_iterators())
{
fe_values.reinit(cell);
for (const auto q : fe_values.quadrature_point_indices())
for (const auto i : fe_values.dof_indices())
for (const auto j : fe_values.dof_indices())
cell_matrix(i,j) += ...; // Do something for DoF indices (i,j)
// at quadrature point q
}

Here, we are looping over all degrees of freedom on all cells, with i and j taking on all valid indices for cell degrees of freedom, as defined by the finite element passed to fe_values.

◆ dof_indices_starting_at()

template<int dim, int spacedim>
std_cxx20::ranges::iota_view<unsigned int, unsigned int> FEValuesBase< dim, spacedim >::dof_indices_starting_at ( const unsigned int  start_dof_index) const

Return an object that can be thought of as an array containing all indices from start_dof_index (inclusive) to dofs_per_cell (exclusive). This allows one to write code using range-based for loops of the following kind:

FEValues<dim> fe_values (...);
for (auto &cell : dof_handler.active_cell_iterators())
{
fe_values.reinit(cell);
for (const auto q : fe_values.quadrature_point_indices())
for (const auto i : fe_values.dof_indices())
for (const auto j : fe_values.dof_indices_starting_at(i))
cell_matrix(i,j) += ...; // Do something for DoF indices (i,j)
// at quadrature point q
}

Here, we are looping over all local degrees of freedom on all cells, with i taking on all valid indices for cell degrees of freedom, as defined by the finite element passed to fe_values, and j taking on a specified subset of i's range, starting at i itself and ending at the number of cell degrees of freedom. In this way, we can construct the upper half and the diagonal of a stiffness matrix contribution (assuming it is symmetric, and that only one half of it needs to be computed), for example.

Note
If the start_dof_index is equal to the number of DoFs in the cell, then the returned index range is empty.

◆ dof_indices_ending_at()

template<int dim, int spacedim>
std_cxx20::ranges::iota_view<unsigned int, unsigned int> FEValuesBase< dim, spacedim >::dof_indices_ending_at ( const unsigned int  end_dof_index) const

Return an object that can be thought of as an array containing all indices from zero (inclusive) to end_dof_index (inclusive). This allows one to write code using range-based for loops of the following kind:

FEValues<dim> fe_values (...);
for (auto &cell : dof_handler.active_cell_iterators())
{
fe_values.reinit(cell);
for (const auto q : fe_values.quadrature_point_indices())
for (const auto i : fe_values.dof_indices())
for (const auto j : fe_values.dof_indices_ending_at(i))
cell_matrix(i,j) += ...; // Do something for DoF indices (i,j)
// at quadrature point q
}

Here, we are looping over all local degrees of freedom on all cells, with i taking on all valid indices for cell degrees of freedom, as defined by the finite element passed to fe_values, and j taking on a specified subset of i's range, starting at zero and ending at i itself. In this way, we can construct the lower half and the diagonal of a stiffness matrix contribution (assuming it is symmetric, and that only one half of it needs to be computed), for example.

Note
If the end_dof_index is equal to zero, then the returned index range is empty.

◆ quadrature_point_indices()

template<int dim, int spacedim>
std_cxx20::ranges::iota_view<unsigned int, unsigned int> FEValuesBase< dim, spacedim >::quadrature_point_indices ( ) const

Return an object that can be thought of as an array containing all indices from zero to n_quadrature_points. This allows to write code using range-based for loops of the following kind:

FEValues<dim> fe_values (...);
for (auto &cell : dof_handler.active_cell_iterators())
{
fe_values.reinit(cell);
for (const auto q_point : fe_values.quadrature_point_indices())
... do something at the quadrature point ...
}

Here, we are looping over all quadrature points on all cells, with q_point taking on all valid indices for quadrature points, as defined by the quadrature rule passed to fe_values.

See also
deal.II and the C++11 standard

◆ quadrature_point()

template<int dim, int spacedim>
const Point<spacedim>& FEValuesBase< dim, spacedim >::quadrature_point ( const unsigned int  q) const

Position of the qth quadrature point in real space.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_quadrature_points flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_quadrature_points()

template<int dim, int spacedim>
const std::vector<Point<spacedim> >& FEValuesBase< dim, spacedim >::get_quadrature_points ( ) const

Return a reference to the vector of quadrature points in real space.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_quadrature_points flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ JxW()

template<int dim, int spacedim>
double FEValuesBase< dim, spacedim >::JxW ( const unsigned int  quadrature_point) const

Mapped quadrature weight. If this object refers to a volume evaluation (i.e. the derived class is of type FEValues), then this is the Jacobi determinant times the weight of the *ith unit quadrature point.

For surface evaluations (i.e. classes FEFaceValues or FESubfaceValues), it is the mapped surface element times the weight of the quadrature point.

You can think of the quantity returned by this function as the volume or surface element \(dx, ds\) in the integral that we implement here by quadrature.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_JxW_values flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_JxW_values()

template<int dim, int spacedim>
const std::vector<double>& FEValuesBase< dim, spacedim >::get_JxW_values ( ) const

Return a reference to the array holding the values returned by JxW().

◆ jacobian()

template<int dim, int spacedim>
const DerivativeForm<1, dim, spacedim>& FEValuesBase< dim, spacedim >::jacobian ( const unsigned int  quadrature_point) const

Return the Jacobian of the transformation at the specified quadrature point, i.e. \(J_{ij}=dx_i/d\hat x_j\)

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jacobians()

template<int dim, int spacedim>
const std::vector<DerivativeForm<1, dim, spacedim> >& FEValuesBase< dim, spacedim >::get_jacobians ( ) const

Return a reference to the array holding the values returned by jacobian().

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ jacobian_grad()

template<int dim, int spacedim>
const DerivativeForm<2, dim, spacedim>& FEValuesBase< dim, spacedim >::jacobian_grad ( const unsigned int  quadrature_point) const

Return the second derivative of the transformation from unit to real cell, i.e. the first derivative of the Jacobian, at the specified quadrature point, i.e. \(G_{ijk}=dJ_{jk}/d\hat x_i\).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_grads flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jacobian_grads()

template<int dim, int spacedim>
const std::vector<DerivativeForm<2, dim, spacedim> >& FEValuesBase< dim, spacedim >::get_jacobian_grads ( ) const

Return a reference to the array holding the values returned by jacobian_grads().

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_grads flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ jacobian_pushed_forward_grad()

template<int dim, int spacedim>
const Tensor<3, spacedim>& FEValuesBase< dim, spacedim >::jacobian_pushed_forward_grad ( const unsigned int  quadrature_point) const

Return the second derivative of the transformation from unit to real cell, i.e. the first derivative of the Jacobian, at the specified quadrature point, pushed forward to the real cell coordinates, i.e. \(G_{ijk}=dJ_{iJ}/d\hat x_K (J_{jJ})^{-1} (J_{kK})^{-1}\).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_pushed_forward_grads flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jacobian_pushed_forward_grads()

template<int dim, int spacedim>
const std::vector<Tensor<3, spacedim> >& FEValuesBase< dim, spacedim >::get_jacobian_pushed_forward_grads ( ) const

Return a reference to the array holding the values returned by jacobian_pushed_forward_grads().

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_pushed_forward_grads flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ jacobian_2nd_derivative()

template<int dim, int spacedim>
const DerivativeForm<3, dim, spacedim>& FEValuesBase< dim, spacedim >::jacobian_2nd_derivative ( const unsigned int  quadrature_point) const

Return the third derivative of the transformation from unit to real cell, i.e. the second derivative of the Jacobian, at the specified quadrature point, i.e. \(G_{ijkl}=\frac{d^2J_{ij}}{d\hat x_k d\hat x_l}\).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_2nd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jacobian_2nd_derivatives()

template<int dim, int spacedim>
const std::vector<DerivativeForm<3, dim, spacedim> >& FEValuesBase< dim, spacedim >::get_jacobian_2nd_derivatives ( ) const

Return a reference to the array holding the values returned by jacobian_2nd_derivatives().

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_2nd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ jacobian_pushed_forward_2nd_derivative()

template<int dim, int spacedim>
const Tensor<4, spacedim>& FEValuesBase< dim, spacedim >::jacobian_pushed_forward_2nd_derivative ( const unsigned int  quadrature_point) const

Return the third derivative of the transformation from unit to real cell, i.e. the second derivative of the Jacobian, at the specified quadrature point, pushed forward to the real cell coordinates, i.e. \(G_{ijkl}=\frac{d^2J_{iJ}}{d\hat x_K d\hat x_L} (J_{jJ})^{-1} (J_{kK})^{-1}(J_{lL})^{-1}\).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_pushed_forward_2nd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jacobian_pushed_forward_2nd_derivatives()

template<int dim, int spacedim>
const std::vector<Tensor<4, spacedim> >& FEValuesBase< dim, spacedim >::get_jacobian_pushed_forward_2nd_derivatives ( ) const

Return a reference to the array holding the values returned by jacobian_pushed_forward_2nd_derivatives().

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_pushed_forward_2nd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ jacobian_3rd_derivative()

template<int dim, int spacedim>
const DerivativeForm<4, dim, spacedim>& FEValuesBase< dim, spacedim >::jacobian_3rd_derivative ( const unsigned int  quadrature_point) const

Return the fourth derivative of the transformation from unit to real cell, i.e. the third derivative of the Jacobian, at the specified quadrature point, i.e. \(G_{ijklm}=\frac{d^2J_{ij}}{d\hat x_k d\hat x_l d\hat x_m}\).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jacobian_3rd_derivatives()

template<int dim, int spacedim>
const std::vector<DerivativeForm<4, dim, spacedim> >& FEValuesBase< dim, spacedim >::get_jacobian_3rd_derivatives ( ) const

Return a reference to the array holding the values returned by jacobian_3rd_derivatives().

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ jacobian_pushed_forward_3rd_derivative()

template<int dim, int spacedim>
const Tensor<5, spacedim>& FEValuesBase< dim, spacedim >::jacobian_pushed_forward_3rd_derivative ( const unsigned int  quadrature_point) const

Return the fourth derivative of the transformation from unit to real cell, i.e. the third derivative of the Jacobian, at the specified quadrature point, pushed forward to the real cell coordinates, i.e. \(G_{ijklm}=\frac{d^3J_{iJ}}{d\hat x_K d\hat x_L d\hat x_M} (J_{jJ})^{-1} (J_{kK})^{-1} (J_{lL})^{-1} (J_{mM})^{-1}\).

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_pushed_forward_3rd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_jacobian_pushed_forward_3rd_derivatives()

template<int dim, int spacedim>
const std::vector<Tensor<5, spacedim> >& FEValuesBase< dim, spacedim >::get_jacobian_pushed_forward_3rd_derivatives ( ) const

Return a reference to the array holding the values returned by jacobian_pushed_forward_3rd_derivatives().

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_jacobian_pushed_forward_2nd_derivatives flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ inverse_jacobian()

template<int dim, int spacedim>
const DerivativeForm<1, spacedim, dim>& FEValuesBase< dim, spacedim >::inverse_jacobian ( const unsigned int  quadrature_point) const

Return the inverse Jacobian of the transformation at the specified quadrature point, i.e. \(J_{ij}=d\hat x_i/dx_j\)

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_inverse_jacobians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_inverse_jacobians()

template<int dim, int spacedim>
const std::vector<DerivativeForm<1, spacedim, dim> >& FEValuesBase< dim, spacedim >::get_inverse_jacobians ( ) const

Return a reference to the array holding the values returned by inverse_jacobian().

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_inverse_jacobians flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ normal_vector()

template<int dim, int spacedim>
const Tensor<1, spacedim>& FEValuesBase< dim, spacedim >::normal_vector ( const unsigned int  i) const

Return the normal vector at a quadrature point. If you call this function for a face (i.e., when using a FEFaceValues or FESubfaceValues object), then this function returns the outward normal vector to the cell at the ith quadrature point of the face.

In contrast, if you call this function for a cell of codimension one (i.e., when using a FEValues<dim,spacedim> object with spacedim>dim), then this function returns the normal vector to the cell – in other words, an approximation to the normal vector to the manifold in which the triangulation is embedded. There are of course two normal directions to a manifold in that case, and this function returns the "up" direction as induced by the numbering of the vertices.

The length of the vector is normalized to one.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_normal_vectors flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

◆ get_normal_vectors()

template<int dim, int spacedim>
const std::vector< Tensor< 1, spacedim > > & FEValuesBase< dim, spacedim >::get_normal_vectors

Return the normal vectors at all quadrature points represented by this object. See the normal_vector() function for what the normal vectors represent.

Note
For this function to work properly, the underlying FEValues, FEFaceValues, or FESubfaceValues object on which you call it must have computed the information you are requesting. To do so, the update_normal_vectors flag must be an element of the list of UpdateFlags that you passed to the constructor of this object. See The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues for more information.

Definition at line 4199 of file fe_values.cc.

◆ operator[]() [1/4]

template<int dim, int spacedim>
const FEValuesViews::Scalar<dim, spacedim>& FEValuesBase< dim, spacedim >::operator[] ( const FEValuesExtractors::Scalar scalar) const

Create a view of the current FEValues object that represents a particular scalar component of the possibly vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.

◆ operator[]() [2/4]

template<int dim, int spacedim>
const FEValuesViews::Vector<dim, spacedim>& FEValuesBase< dim, spacedim >::operator[] ( const FEValuesExtractors::Vector vector) const

Create a view of the current FEValues object that represents a set of dim scalar components (i.e. a vector) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.

◆ operator[]() [3/4]

template<int dim, int spacedim>
const FEValuesViews::SymmetricTensor<2, dim, spacedim>& FEValuesBase< dim, spacedim >::operator[] ( const FEValuesExtractors::SymmetricTensor< 2 > &  tensor) const

Create a view of the current FEValues object that represents a set of (dim*dim + dim)/2 scalar components (i.e. a symmetric 2nd order tensor) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.

◆ operator[]() [4/4]

template<int dim, int spacedim>
const FEValuesViews::Tensor<2, dim, spacedim>& FEValuesBase< dim, spacedim >::operator[] ( const FEValuesExtractors::Tensor< 2 > &  tensor) const

Create a view of the current FEValues object that represents a set of (dim*dim) scalar components (i.e. a 2nd order tensor) of the vector-valued finite element. The concept of views is explained in the documentation of the namespace FEValuesViews and in particular in the Handling vector valued problems module.

◆ get_mapping()

template<int dim, int spacedim>
const Mapping<dim, spacedim>& FEValuesBase< dim, spacedim >::get_mapping ( ) const

Constant reference to the selected mapping object.

◆ get_fe()

template<int dim, int spacedim>
const FiniteElement<dim, spacedim>& FEValuesBase< dim, spacedim >::get_fe ( ) const

Constant reference to the selected finite element object.

◆ get_update_flags()

template<int dim, int spacedim>
UpdateFlags FEValuesBase< dim, spacedim >::get_update_flags ( ) const

Return the update flags set for this object.

◆ get_cell()

template<int dim, int spacedim>
const Triangulation< dim, spacedim >::cell_iterator FEValuesBase< dim, spacedim >::get_cell

Return a triangulation iterator to the current cell.

Definition at line 4190 of file fe_values.cc.

◆ get_cell_similarity()

template<int dim, int spacedim>
CellSimilarity::Similarity FEValuesBase< dim, spacedim >::get_cell_similarity

Return the relation of the current cell to the previous cell. This allows re-use of some cell data (like local matrices for equations with constant coefficients) if the result is CellSimilarity::translation.

Definition at line 4367 of file fe_values.cc.

◆ memory_consumption()

template<int dim, int spacedim>
std::size_t FEValuesBase< dim, spacedim >::memory_consumption

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 4212 of file fe_values.cc.

◆ invalidate_present_cell()

template<int dim, int spacedim>
void FEValuesBase< dim, spacedim >::invalidate_present_cell
protected

A function that is connected to the triangulation in order to reset the stored 'present_cell' iterator to an invalid one whenever the triangulation is changed and the iterator consequently becomes invalid.

Definition at line 4251 of file fe_values.cc.

◆ maybe_invalidate_previous_present_cell()

template<int dim, int spacedim>
void FEValuesBase< dim, spacedim >::maybe_invalidate_previous_present_cell ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell)
protected

This function is called by the various reinit() functions in derived classes. Given the cell indicated by the argument, test whether we have to throw away the previously stored present_cell argument because it would require us to compare cells from different triangulations. In checking all this, also make sure that we have tria_listener connected to the triangulation to which we will set present_cell right after calling this function.

Definition at line 4269 of file fe_values.cc.

◆ compute_update_flags()

template<int dim, int spacedim>
UpdateFlags FEValuesBase< dim, spacedim >::compute_update_flags ( const UpdateFlags  update_flags) const
protected

Initialize some update flags. Called from the initialize functions of derived classes, which are in turn called from their constructors.

Basically, this function finds out using the finite element and mapping object already stored which flags need to be set to compute everything the user wants, as expressed through the flags passed as argument.

Definition at line 4232 of file fe_values.cc.

◆ check_cell_similarity()

template<int dim, int spacedim>
void FEValuesBase< dim, spacedim >::check_cell_similarity ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell)
inlineprotected

A function that checks whether the new cell is similar to the one previously used. Then, a significant amount of the data can be reused, e.g. the derivatives of the basis functions in real space, shape_grad.

Definition at line 4312 of file fe_values.cc.

Friends And Related Function Documentation

◆ FEValuesViews::Scalar

template<int dim, int spacedim>
template<int , int >
friend class FEValuesViews::Scalar
friend

Definition at line 3594 of file fe_values.h.

◆ FEValuesViews::Vector

template<int dim, int spacedim>
template<int , int >
friend class FEValuesViews::Vector
friend

Definition at line 3596 of file fe_values.h.

◆ FEValuesViews::SymmetricTensor

template<int dim, int spacedim>
template<int , int , int >
friend class FEValuesViews::SymmetricTensor
friend

Definition at line 3598 of file fe_values.h.

◆ FEValuesViews::Tensor

template<int dim, int spacedim>
template<int , int , int >
friend class FEValuesViews::Tensor
friend

Definition at line 3600 of file fe_values.h.

Member Data Documentation

◆ dimension

template<int dim, int spacedim>
const unsigned int FEValuesBase< dim, spacedim >::dimension = dim
static

Dimension in which this object operates.

Definition at line 2091 of file fe_values.h.

◆ space_dimension

template<int dim, int spacedim>
const unsigned int FEValuesBase< dim, spacedim >::space_dimension = spacedim
static

Dimension of the space in which this object operates.

Definition at line 2096 of file fe_values.h.

◆ n_quadrature_points

template<int dim, int spacedim>
const unsigned int FEValuesBase< dim, spacedim >::n_quadrature_points

Number of quadrature points.

Definition at line 2101 of file fe_values.h.

◆ dofs_per_cell

template<int dim, int spacedim>
const unsigned int FEValuesBase< dim, spacedim >::dofs_per_cell

Number of shape functions per cell. If we use this base class to evaluate a finite element on faces of cells, this is still the number of degrees of freedom per cell, not per face.

Definition at line 2108 of file fe_values.h.

◆ present_cell

template<int dim, int spacedim>
std::unique_ptr<const CellIteratorBase> FEValuesBase< dim, spacedim >::present_cell
protected

Store the cell selected last time the reinit() function was called. This is necessary for the get_function_* functions as well as the functions of same name in the extractor classes.

Definition at line 3458 of file fe_values.h.

◆ tria_listener_refinement

template<int dim, int spacedim>
boost::signals2::connection FEValuesBase< dim, spacedim >::tria_listener_refinement
protected

A signal connection we use to ensure we get informed whenever the triangulation changes by refinement. We need to know about that because it invalidates all cell iterators and, as part of that, the 'present_cell' iterator we keep around between subsequent calls to reinit() in order to compute the cell similarity.

Definition at line 3474 of file fe_values.h.

◆ tria_listener_mesh_transform

template<int dim, int spacedim>
boost::signals2::connection FEValuesBase< dim, spacedim >::tria_listener_mesh_transform
protected

A signal connection we use to ensure we get informed whenever the triangulation changes by mesh transformations. We need to know about that because it invalidates all cell iterators and, as part of that, the 'present_cell' iterator we keep around between subsequent calls to reinit() in order to compute the cell similarity.

Definition at line 3483 of file fe_values.h.

◆ mapping

template<int dim, int spacedim>
const SmartPointer<const Mapping<dim, spacedim>, FEValuesBase<dim, spacedim> > FEValuesBase< dim, spacedim >::mapping
protected

A pointer to the mapping object associated with this FEValues object.

Definition at line 3510 of file fe_values.h.

◆ mapping_data

template<int dim, int spacedim>
std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> FEValuesBase< dim, spacedim >::mapping_data
protected

A pointer to the internal data object of mapping, obtained from Mapping::get_data(), Mapping::get_face_data(), or Mapping::get_subface_data().

Definition at line 3518 of file fe_values.h.

◆ mapping_output

template<int dim, int spacedim>
::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> FEValuesBase< dim, spacedim >::mapping_output
protected

An object into which the Mapping::fill_fe_values() and similar functions place their output.

Definition at line 3525 of file fe_values.h.

◆ fe

template<int dim, int spacedim>
const SmartPointer<const FiniteElement<dim, spacedim>, FEValuesBase<dim, spacedim> > FEValuesBase< dim, spacedim >::fe
protected

A pointer to the finite element object associated with this FEValues object.

Definition at line 3534 of file fe_values.h.

◆ fe_data

template<int dim, int spacedim>
std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase> FEValuesBase< dim, spacedim >::fe_data
protected

A pointer to the internal data object of finite element, obtained from FiniteElement::get_data(), Mapping::get_face_data(), or FiniteElement::get_subface_data().

Definition at line 3542 of file fe_values.h.

◆ finite_element_output

template<int dim, int spacedim>
::internal::FEValuesImplementation::FiniteElementRelatedData<dim, spacedim> FEValuesBase< dim, spacedim >::finite_element_output
protected

An object into which the FiniteElement::fill_fe_values() and similar functions place their output.

Definition at line 3550 of file fe_values.h.

◆ update_flags

template<int dim, int spacedim>
UpdateFlags FEValuesBase< dim, spacedim >::update_flags
protected

Original update flags handed to the constructor of FEValues.

Definition at line 3556 of file fe_values.h.

◆ cell_similarity

template<int dim, int spacedim>
CellSimilarity::Similarity FEValuesBase< dim, spacedim >::cell_similarity
protected

An enum variable that can store different states of the current cell in comparison to the previously visited cell. If wanted, additional states can be checked here and used in one of the methods used during reinit.

Definition at line 3574 of file fe_values.h.

◆ fe_values_views_cache

template<int dim, int spacedim>
::internal::FEValuesViews::Cache<dim, spacedim> FEValuesBase< dim, spacedim >::fe_values_views_cache
private

A cache for all possible FEValuesViews objects.

Definition at line 3589 of file fe_values.h.


The documentation for this class was generated from the following files:
LocalIntegrators::Advection::cell_matrix
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double >> &velocity, const double factor=1.)
Definition: advection.h:80
OpenCASCADE::point
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
Definition: utilities.cc:188
FEValues
Definition: fe.h:38
FEValuesBase::mapping
const SmartPointer< const Mapping< dim, spacedim >, FEValuesBase< dim, spacedim > > mapping
Definition: fe_values.h:3510
LAPACKSupport::A
static const char A
Definition: lapack_support.h:155
FullMatrix< double >