Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
FEValues< dim, spacedim > Class Template Reference

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

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

Public Member Functions

 FEValues (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
 
 FEValues (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim > &quadrature, const UpdateFlags update_flags)
 
 FEValues (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, const UpdateFlags update_flags)
 
 FEValues (const FiniteElement< dim, spacedim > &fe, const hp::QCollection< dim > &quadrature, const UpdateFlags update_flags)
 
template<bool level_dof_access>
void reinit (const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
 
void reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell)
 
const Quadrature< dim > & get_quadrature () const
 
std::size_t memory_consumption () const
 
const FEValues< dim, spacedim > & get_present_fe_values () const
 
void always_allow_check_for_cell_similarity (const bool allow)
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Access to shape function values

These fields are filled by the finite element.

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

Hessian matrices and Laplacians of global finite element fields

template<typename Number >
void get_function_hessians (const ReadVector< Number > &fe_function, std::vector< Tensor< 2, spacedim, Number > > &hessians) const
 
template<typename Number >
void get_function_hessians (const ReadVector< Number > &fe_function, std::vector< std::vector< Tensor< 2, spacedim, Number > > > &hessians, const bool quadrature_points_fastest=false) const
 
template<typename Number >
void get_function_hessians (const ReadVector< Number > &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Tensor< 2, spacedim, Number > > &hessians) const
 
template<typename Number >
void get_function_hessians (const ReadVector< Number > &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< Tensor< 2, spacedim, Number > > > hessians, const bool quadrature_points_fastest=false) const
 
template<typename Number >
void get_function_laplacians (const ReadVector< Number > &fe_function, std::vector< Number > &laplacians) const
 
template<typename Number >
void get_function_laplacians (const ReadVector< Number > &fe_function, std::vector< Vector< Number > > &laplacians) const
 
template<typename Number >
void get_function_laplacians (const ReadVector< Number > &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Number > &laplacians) const
 
template<typename Number >
void get_function_laplacians (const ReadVector< Number > &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Vector< Number > > &laplacians) const
 
template<typename Number >
void get_function_laplacians (const ReadVector< Number > &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< std::vector< Number > > &laplacians, const bool quadrature_points_fastest=false) const
 
Access to third derivatives of global finite element fields
template<typename Number >
void get_function_third_derivatives (const ReadVector< Number > &fe_function, std::vector< Tensor< 3, spacedim, Number > > &third_derivatives) const
 
template<typename Number >
void get_function_third_derivatives (const ReadVector< Number > &fe_function, std::vector< std::vector< Tensor< 3, spacedim, Number > > > &third_derivatives, const bool quadrature_points_fastest=false) const
 
template<typename Number >
void get_function_third_derivatives (const ReadVector< Number > &fe_function, const ArrayView< const types::global_dof_index > &indices, std::vector< Tensor< 3, spacedim, Number > > &third_derivatives) const
 
template<typename Number >
void get_function_third_derivatives (const ReadVector< Number > &fe_function, const ArrayView< const types::global_dof_index > &indices, ArrayView< std::vector< Tensor< 3, spacedim, Number > > > third_derivatives, const 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_point) const
 
const std::vector< Point< spacedim > > & get_quadrature_points () const
 
double JxW (const unsigned int q_point) const
 
const std::vector< double > & get_JxW_values () const
 
const DerivativeForm< 1, dim, spacedim > & jacobian (const unsigned int q_point) const
 
const std::vector< DerivativeForm< 1, dim, spacedim > > & get_jacobians () const
 
const DerivativeForm< 2, dim, spacedim > & jacobian_grad (const unsigned int q_point) const
 
const std::vector< DerivativeForm< 2, dim, spacedim > > & get_jacobian_grads () const
 
const Tensor< 3, spacedim > & jacobian_pushed_forward_grad (const unsigned int q_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 q_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 q_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 q_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 q_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 q_point) const
 
const std::vector< DerivativeForm< 1, spacedim, dim > > & get_inverse_jacobians () const
 
const Tensor< 1, spacedim > & normal_vector (const unsigned int q_point) 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
 
Access to the raw data
const Mapping< dim, spacedim > & get_mapping () const
 
const FiniteElement< dim, spacedim > & get_fe () const
 
UpdateFlags get_update_flags () const
 
Triangulation< dim, spacedim >::cell_iterator get_cell () const
 
CellSimilarity::Similarity get_cell_similarity () const
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcAccessToUninitializedField (std::string arg1)
 
static ::ExceptionBaseExcNotReinited ()
 
static ::ExceptionBaseExcFEDontMatch ()
 
static ::ExceptionBaseExcShapeFunctionNotPrimitive (int arg1)
 
static ::ExceptionBaseExcFENotPrimitive ()
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Public Attributes

const unsigned int n_quadrature_points
 
const unsigned int max_n_quadrature_points
 
const unsigned int dofs_per_cell
 

Static Public Attributes

static constexpr unsigned int integral_dimension = dim
 
static constexpr unsigned int dimension = dim
 
static constexpr unsigned int space_dimension = spacedim
 

Protected Member Functions

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)
 

Protected Attributes

CellIteratorWrapper present_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
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void initialize (const UpdateFlags update_flags)
 
void do_reinit ()
 
void check_no_subscribers () const noexcept
 

Private Attributes

const Quadrature< dim > quadrature
 
::internal::FEValuesViews::Cache< dim, spacedim > fe_values_views_cache
 
bool check_for_cell_similarity_allowed
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Detailed Description

template<int dim, int spacedim = dim>
class FEValues< dim, spacedim >

Finite element evaluated in quadrature points of a cell.

This function implements the initialization routines for FEValuesBase, if values in quadrature points of a cell are needed. For further documentation see this class.

Definition at line 62 of file fe_values.h.

Member Typedef Documentation

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ FEValues() [1/4]

template<int dim, int spacedim = dim>
FEValues< dim, spacedim >::FEValues ( const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags  update_flags 
)

Constructor. Gets cell independent data from mapping and finite element objects, matching the quadrature rule and update flags.

◆ FEValues() [2/4]

template<int dim, int spacedim = dim>
FEValues< dim, spacedim >::FEValues ( const Mapping< dim, spacedim > &  mapping,
const FiniteElement< dim, spacedim > &  fe,
const hp::QCollection< dim > &  quadrature,
const UpdateFlags  update_flags 
)

Like the function above, but taking a collection of quadrature rules.

Note
We require, in contrast to FEFaceValues, that the number of quadrature rules in the collection is one.

◆ FEValues() [3/4]

template<int dim, int spacedim = dim>
FEValues< dim, spacedim >::FEValues ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
const UpdateFlags  update_flags 
)

Constructor. This constructor is equivalent to the other one except that it makes the object use a \(Q_1\) mapping (i.e., an object of type MappingQ(1)) implicitly.

◆ FEValues() [4/4]

template<int dim, int spacedim = dim>
FEValues< dim, spacedim >::FEValues ( const FiniteElement< dim, spacedim > &  fe,
const hp::QCollection< dim > &  quadrature,
const UpdateFlags  update_flags 
)

Like the function above, but taking a collection of quadrature rules.

Note
We require, in contrast to FEFaceValues, that the number of quadrature rules in the collection is one.

Member Function Documentation

◆ reinit() [1/2]

template<int dim, int spacedim = dim>
template<bool level_dof_access>
void FEValues< dim, spacedim >::reinit ( const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &  cell)

Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a DoFHandler object", and the finite element associated with this object. It is assumed that the finite element used by the given cell is also the one used by this FEValues object.

◆ reinit() [2/2]

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::reinit ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell)

Reinitialize the gradients, Jacobi determinants, etc for the given cell of type "iterator into a Triangulation object", and the given finite element. Since iterators into triangulation alone only convey information about the geometry of a cell, but not about degrees of freedom possibly associated with this cell, you will not be able to call some functions of this class if they need information about degrees of freedom. These functions are, above all, the get_function_value/gradients/hessians/laplacians/third_derivatives functions. If you want to call these functions, you have to call the reinit variants that take iterators into DoFHandler or other DoF handler type objects.

◆ get_quadrature()

template<int dim, int spacedim = dim>
const Quadrature< dim > & FEValues< dim, spacedim >::get_quadrature ( ) const

Return a reference to the copy of the quadrature formula stored by this object.

◆ memory_consumption()

template<int dim, int spacedim = dim>
std::size_t FEValues< dim, spacedim >::memory_consumption ( ) const

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

◆ get_present_fe_values()

template<int dim, int spacedim = dim>
const FEValues< dim, spacedim > & FEValues< dim, spacedim >::get_present_fe_values ( ) const

Return a reference to this very object.

Though it seems that it is not very useful, this function is there to provide capability to the hp::FEValues class, in which case it provides the FEValues object for the present cell (remember that for hp-finite elements, the actual FE object used may change from cell to cell, so we also need different FEValues objects for different cells; once you reinitialize the hp::FEValues object for a specific cell, it retrieves the FEValues object for the FE on that cell and returns it through a function of the same name as this one; this function here therefore only provides the same interface so that one can templatize on FEValues and hp::FEValues).

◆ initialize()

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::initialize ( const UpdateFlags  update_flags)
private

Do work common to the two constructors.

◆ do_reinit()

template<int dim, int spacedim = dim>
void FEValues< dim, spacedim >::do_reinit ( )
private

The reinit() functions do only that part of the work that requires knowledge of the type of iterator. After setting present_cell(), they pass on to this function, which does the real work, and which is independent of the actual type of the cell iterator.

◆ always_allow_check_for_cell_similarity()

template<int dim, int spacedim>
void FEValuesBase< dim, spacedim >::always_allow_check_for_cell_similarity ( const bool  allow)
inherited

Explicitly allow to check for cell similarity. The detection of simple geometries with CellSimilarity is sensitive to the first cell detected. When using multiple threads, each thread might get a thread local copy of the FEValues object that is initialized to the first cell the thread sees. As this cell might be different between runs and number of threads used, this slight deviation leads to difference in roundoff errors that propagate through the program. Therefore, deal.II disables the CellSimilarity check by default in programs that use more than one thread. This function can be used to disable this behavior: When called, FEValues objects will always do the similarity check, even in cases where the program uses multiple threads. This substantially accelerates the operations of FEValues because many operations can be avoided if a cell is, for example, just a translation of the previous cell. On the other hand, you might get results that differ by an amount proportional to round-off between the case of using or not using cell similarity information, and because the order of cells assigned to individual threads may differ from run to run, when you call this function you may end up with results that differ by round-off between runs of the same program.

Definition at line 247 of file fe_values_base.cc.

◆ shape_value()

template<int dim, int spacedim>
const double & FEValuesBase< dim, spacedim >::shape_value ( const unsigned int  i,
const unsigned int  q_point 
) const
inherited

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
iNumber of the shape function \(\varphi_i\) to be evaluated. Note that this number runs from zero to dofs_per_cell, even in the case of an FEFaceValues or FESubfaceValues object.
q_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_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  i,
const unsigned int  q_point,
const unsigned int  component 
) const
inherited

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
iNumber of the shape function \(\varphi_i\) to be evaluated.
q_pointNumber 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  i,
const unsigned int  q_point 
) const
inherited

Compute the gradient of the ith 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
iNumber of the shape function \(\varphi_i\) to be evaluated.
q_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  i,
const unsigned int  q_point,
const unsigned int  component 
) const
inherited

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  i,
const unsigned int  q_point 
) const
inherited

Second derivatives of the ith shape function at the q_pointth 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  i,
const unsigned int  q_point,
const unsigned int  component 
) const
inherited

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  i,
const unsigned int  q_point 
) const
inherited

Third derivatives of the ith shape function at the q_pointth 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  i,
const unsigned int  q_point,
const unsigned int  component 
) const
inherited

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<typename Number >
void FEValuesBase< dim, spacedim >::get_function_values ( const ReadVector< Number > &  fe_function,
std::vector< Number > &  values 
) const
inherited

Return the values of a finite element function at the quadrature points of the current cell, face, or subface (selected the last time the reinit() function was called). That is, if the first argument fe_function is a vector of nodal values of a finite element function \(u_h(\mathbf x)\) defined on a DoFHandler object, then the output vector (the second argument, values) is the vector of values \(u_h(\mathbf x_q^K)\) where \(x_q^K\) are the quadrature points on the current cell \(K\). This function is first discussed in the Results section of step-4, and the related get_function_gradients() function is also used in step-15 along with numerous other tutorial programs.

If the current cell is not active (i.e., it has children), then the finite element function is, strictly speaking, defined by shape functions that live on these child cells. Rather than evaluating the shape functions on the child cells, with the quadrature points defined on the current cell, this function first interpolates the finite element function to shape functions defined on the current cell, and then evaluates this interpolated function.

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
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 716 of file fe_values_base.cc.

◆ get_function_values() [2/5]

template<int dim, int spacedim>
template<typename Number >
void FEValuesBase< dim, spacedim >::get_function_values ( const ReadVector< Number > &  fe_function,
std::vector< Vector< Number > > &  values 
) const
inherited

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 763 of file fe_values_base.cc.

◆ get_function_values() [3/5]

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

Generate function values from an arbitrary vector. This function does in essence the same as the first function of this name above, except that it does not make the assumption that the input vector corresponds to a DoFHandler that describes the unknowns of a finite element field (and for which we would then assume that fe_function.size() == dof_handler.n_dofs()). Rather, the nodal values corresponding to the current cell are elements of an otherwise arbitrary vector, and these elements are indexed by the second argument to this function. What the rest of the fe_function input argument corresponds to is of no consequence to this function.

Given this, the function above corresponds to passing fe_function as first argument to the current function, and using the local_dof_indices array that results from the following call as second argument to the current function:

cell->get_dof_indices (local_dof_indices);

(See DoFCellAccessor::get_dof_indices() for more information.)

Likewise, the function above is equivalent to calling

cell->get_dof_values (fe_function, local_dof_values);

and then calling the current function with local_dof_values as first argument, and an array with indices {0,...,fe.dofs_per_cell-1} as second argument.

The point of the current function is that one sometimes wants to evaluate finite element functions at quadrature points with nodal values that are not stored in a global vector – for example, one could modify these local values first, such as by applying a limiter or by ensuring that all nodal values are positive, before evaluating the finite element field that corresponds to these local values on the current cell. Another application is where one wants to postprocess the solution on a cell into a different finite element space on every cell, without actually creating a corresponding DoFHandler – in that case, all one would compute is a local representation of that postprocessed function, characterized by its nodal values; this function then allows the evaluation of that representation at quadrature points.

Parameters
[in]fe_functionA vector of nodal values. This vector can have an arbitrary size, as long as all elements index by indices can actually be accessed.
[in]indicesA vector of indices into fe_function. This vector must have length equal to the number of degrees of freedom on the current cell, and must identify elements in fe_function in the order in which degrees of freedom are indexed on the reference cell.
[out]valuesA vector of values of the given finite element field, at the quadrature points on the current object.
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 740 of file fe_values_base.cc.

◆ get_function_values() [4/5]

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

Generate vector function values from an arbitrary vector.

This function corresponds to the previous one, just for the vector-valued case.

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 789 of file fe_values_base.cc.

◆ get_function_values() [5/5]

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

Generate vector function values from an arbitrary vector. This function is similar to the previous one, but the indices vector may also 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.

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 819 of file fe_values_base.cc.

◆ get_function_gradients() [1/4]

template<int dim, int spacedim>
template<typename Number >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const ReadVector< Number > &  fe_function,
std::vector< Tensor< 1, spacedim, Number > > &  gradients 
) const
inherited

Return the gradients of a finite element function at the quadrature points of the current cell, face, or subface (selected the last time the reinit() function was called). That is, if the first argument fe_function is a vector of nodal values of a finite element function \(u_h(\mathbf x)\) defined on a DoFHandler object, then the output vector (the second argument, values) is the vector of values \(\nabla u_h(\mathbf x_q^K)\) where \(x_q^K\) are the quadrature points on the current cell \(K\). This function is first discussed in the Results section of step-4, and it is also used in step-15 along with numerous other tutorial programs.

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 851 of file fe_values_base.cc.

◆ get_function_gradients() [2/4]

template<int dim, int spacedim>
template<typename Number >
void FEValuesBase< dim, spacedim >::get_function_gradients ( const ReadVector< Number > &  fe_function,
std::vector< std::vector< Tensor< 1, spacedim, Number > > > &  gradients 
) const
inherited

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 898 of file fe_values_base.cc.

◆ get_function_gradients() [3/4]

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

This function relates to the first of the get_function_gradients() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 875 of file fe_values_base.cc.

◆ get_function_gradients() [4/4]

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

This function relates to the first of the get_function_gradients() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 923 of file fe_values_base.cc.

◆ get_function_hessians() [1/4]

template<int dim, int spacedim>
template<typename Number >
void FEValuesBase< dim, spacedim >::get_function_hessians ( const ReadVector< Number > &  fe_function,
std::vector< Tensor< 2, spacedim, Number > > &  hessians 
) const
inherited

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 954 of file fe_values_base.cc.

◆ get_function_hessians() [2/4]

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

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 1001 of file fe_values_base.cc.

◆ get_function_hessians() [3/4]

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

This function relates to the first of the get_function_hessians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 978 of file fe_values_base.cc.

◆ get_function_hessians() [4/4]

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

This function relates to the first of the get_function_hessians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 1028 of file fe_values_base.cc.

◆ get_function_laplacians() [1/5]

template<int dim, int spacedim>
template<typename Number >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const ReadVector< Number > &  fe_function,
std::vector< Number > &  laplacians 
) const
inherited

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 1057 of file fe_values_base.cc.

◆ get_function_laplacians() [2/5]

template<int dim, int spacedim>
template<typename Number >
void FEValuesBase< dim, spacedim >::get_function_laplacians ( const ReadVector< Number > &  fe_function,
std::vector< Vector< Number > > &  laplacians 
) const
inherited

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 1104 of file fe_values_base.cc.

◆ get_function_laplacians() [3/5]

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

This function relates to the first of the get_function_laplacians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 1081 of file fe_values_base.cc.

◆ get_function_laplacians() [4/5]

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

This function relates to the first of the get_function_laplacians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 1129 of file fe_values_base.cc.

◆ get_function_laplacians() [5/5]

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

This function relates to the first of the get_function_laplacians() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 1159 of file fe_values_base.cc.

◆ get_function_third_derivatives() [1/4]

template<int dim, int spacedim>
template<typename Number >
void FEValuesBase< dim, spacedim >::get_function_third_derivatives ( const ReadVector< Number > &  fe_function,
std::vector< Tensor< 3, spacedim, Number > > &  third_derivatives 
) const
inherited

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 1188 of file fe_values_base.cc.

◆ get_function_third_derivatives() [2/4]

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

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 1234 of file fe_values_base.cc.

◆ get_function_third_derivatives() [3/4]

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

This function relates to the first of the get_function_third_derivatives() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 1212 of file fe_values_base.cc.

◆ get_function_third_derivatives() [4/4]

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

This function relates to the first of the get_function_third_derivatives() function above in the same way as the get_function_values() with similar arguments relates to the first of the get_function_values() functions. See there for more information.

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 1261 of file fe_values_base.cc.

◆ dof_indices()

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

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 (...);
FullMatrix<double> cell_matrix (...);
for (auto &cell : dof_handler.active_cell_iterators())
{
cell_matrix = 0;
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
}
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const

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
inherited

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 (...);
FullMatrix<double> cell_matrix (...);
for (auto &cell : dof_handler.active_cell_iterators())
{
cell_matrix = 0;
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
}
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices_starting_at(const unsigned int start_dof_index) const

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
inherited

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 (...);
FullMatrix<double> cell_matrix (...);
for (auto &cell : dof_handler.active_cell_iterators())
{
cell_matrix = 0;
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
}
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices_ending_at(const unsigned int end_dof_index) const

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
inherited

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 ...
}
const Quadrature< dim > quadrature
Definition fe_values.h:172

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 Modern C++ standards

◆ quadrature_point()

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

Return the location of the q_pointth 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
inherited

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  q_point) const
inherited

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 q_pointth 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
inherited

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  q_point) const
inherited

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
inherited

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  q_point) const
inherited

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
inherited

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  q_point) const
inherited

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
inherited

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  q_point) const
inherited

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
inherited

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  q_point) const
inherited

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
inherited

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  q_point) const
inherited

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
inherited

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  q_point) const
inherited

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
inherited

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  q_point) const
inherited

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
inherited

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  q_point) const
inherited

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 q_pointth 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 ( ) const
inherited

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 1298 of file fe_values_base.cc.

◆ operator[]() [1/4]

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

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
inherited

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
inherited

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
inherited

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
inherited

Constant reference to the selected mapping object.

◆ get_fe()

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

Constant reference to the selected finite element object.

◆ get_update_flags()

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

Return the update flags set for this object.

◆ get_cell()

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

Return a triangulation iterator to the current cell.

Definition at line 1289 of file fe_values_base.cc.

◆ get_cell_similarity()

template<int dim, int spacedim>
CellSimilarity::Similarity FEValuesBase< dim, spacedim >::get_cell_similarity ( ) const
inherited

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 1453 of file fe_values_base.cc.

◆ invalidate_present_cell()

template<int dim, int spacedim>
void FEValuesBase< dim, spacedim >::invalidate_present_cell ( )
protectedinherited

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 1351 of file fe_values_base.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)
protectedinherited

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 1369 of file fe_values_base.cc.

◆ compute_update_flags()

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

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 1332 of file fe_values_base.cc.

◆ check_cell_similarity()

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

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 1412 of file fe_values_base.cc.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Member Data Documentation

◆ integral_dimension

template<int dim, int spacedim = dim>
constexpr unsigned int FEValues< dim, spacedim >::integral_dimension = dim
staticconstexpr

Dimension of the object over which we integrate. For the present class, this is equal to dim.

Definition at line 69 of file fe_values.h.

◆ quadrature

template<int dim, int spacedim = dim>
const Quadrature<dim> FEValues< dim, spacedim >::quadrature
private

Store a copy of the quadrature formula here.

Definition at line 172 of file fe_values.h.

◆ dimension

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

Dimension in which this object operates.

Definition at line 160 of file fe_values_base.h.

◆ space_dimension

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

Dimension of the space in which this object operates.

Definition at line 165 of file fe_values_base.h.

◆ n_quadrature_points

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

Number of quadrature points of the current object. Its value is initialized by the value of max_n_quadrature_points and is updated, e.g., if FEFaceValues::reinit() is called for a new cell/face.

Note
The default value equals to the value of max_n_quadrature_points.

Definition at line 174 of file fe_values_base.h.

◆ max_n_quadrature_points

template<int dim, int spacedim>
const unsigned int FEValuesBase< dim, spacedim >::max_n_quadrature_points
inherited

Maximum number of quadrature points. This value might be different from n_quadrature_points, e.g., if a QCollection with different face quadrature rules has been passed to initialize FEFaceValues.

This is mostly useful to initialize arrays to allocate the maximum amount of memory that may be used when re-sizing later on to a the current number of quadrature points given by n_quadrature_points.

Definition at line 185 of file fe_values_base.h.

◆ dofs_per_cell

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

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 192 of file fe_values_base.h.

◆ present_cell

template<int dim, int spacedim>
CellIteratorWrapper FEValuesBase< dim, spacedim >::present_cell
protectedinherited

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 1657 of file fe_values_base.h.

◆ tria_listener_refinement

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

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 1666 of file fe_values_base.h.

◆ tria_listener_mesh_transform

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

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 1675 of file fe_values_base.h.

◆ mapping

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

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

Definition at line 1702 of file fe_values_base.h.

◆ mapping_data

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

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 1710 of file fe_values_base.h.

◆ mapping_output

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

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

Definition at line 1717 of file fe_values_base.h.

◆ fe

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

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

Definition at line 1725 of file fe_values_base.h.

◆ fe_data

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

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 1733 of file fe_values_base.h.

◆ finite_element_output

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

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

Definition at line 1741 of file fe_values_base.h.

◆ update_flags

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

Original update flags handed to the constructor of FEValues.

Definition at line 1747 of file fe_values_base.h.

◆ cell_similarity

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

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 1765 of file fe_values_base.h.

◆ fe_values_views_cache

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

A cache for all possible FEValuesViews objects.

Definition at line 1780 of file fe_values_base.h.

◆ check_for_cell_similarity_allowed

template<int dim, int spacedim>
bool FEValuesBase< dim, spacedim >::check_for_cell_similarity_allowed
privateinherited

Whether checking for cell similarity is allowed.

Definition at line 1785 of file fe_values_base.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


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