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

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

Public Types

using value_type = typename internal::FEPointEvaluation::EvaluatorTypeTraits< dim, n_components, Number >::value_type
 
using gradient_type = typename internal::FEPointEvaluation::EvaluatorTypeTraits< dim, n_components, Number >::gradient_type
 

Public Member Functions

 FEPointEvaluation (const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const UpdateFlags update_flags, const unsigned int first_selected_component=0)
 
void reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points)
 
void evaluate (const ArrayView< const Number > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
void integrate (const ArrayView< Number > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags)
 
const value_typeget_value (const unsigned int point_index) const
 
void submit_value (const value_type &value, const unsigned int point_index)
 
const gradient_typeget_gradient (const unsigned int point_index) const
 
const gradient_typeget_unit_gradient (const unsigned int point_index) const
 
void submit_gradient (const gradient_type &, const unsigned int point_index)
 
DerivativeForm< 1, dim, spacedim > jacobian (const unsigned int point_index) const
 
DerivativeForm< 1, spacedim, dim > inverse_jacobian (const unsigned int point_index) const
 
Point< spacedim > real_point (const unsigned int point_index) const
 
Point< dim > unit_point (const unsigned int point_index) const
 

Private Attributes

SmartPointer< const Mapping< dim, spacedim > > mapping
 
const MappingQGeneric< dim, spacedim > * mapping_q_generic
 
SmartPointer< const FiniteElement< dim > > fe
 
std::vector< Polynomials::Polynomial< double > > poly
 
bool polynomials_are_hat_functions
 
std::vector< unsigned intrenumber
 
std::vector< value_typesolution_renumbered
 
AlignedVector< typename internal::FEPointEvaluation::EvaluatorTypeTraits< dim, n_components, VectorizedArray< Number > >::value_typesolution_renumbered_vectorized
 
std::vector< value_typevalues
 
std::vector< gradient_typeunit_gradients
 
std::vector< gradient_typegradients
 
unsigned int dofs_per_component
 
std::vector< std::array< bool, n_components > > nonzero_shape_function_component
 
UpdateFlags update_flags
 
UpdateFlags update_flags_mapping
 
std::shared_ptr< FEValues< dim, spacedim > > fe_values
 
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > mapping_data
 
std::vector< Point< dim > > unit_points
 

Detailed Description

template<int n_components, int dim, int spacedim = dim, typename Number = double>
class FEPointEvaluation< n_components, dim, spacedim, Number >

This class provides an interface to the evaluation of interpolated solution values and gradients on cells on arbitrary reference point positions. These points can change from cell to cell, both with respect to their quantity as well to the location. The two typical use cases are evaluations on non-matching grids and particle simulations.

The use of this class is similar to FEValues or FEEvaluation: The class is first initialized to a cell by calling FEPointEvaluation::reinit(cell, unit_points), with the main difference to the other concepts that the underlying points in reference coordinates need to be passed along. Then, upon call to evaluate() or integrate(), the user can compute information at the give points. Eventually, the access functions get_value() or get_gradient() allow to query this information at a specific point index.

The functionality is similar to creating an FEValues object with a Quadrature object on the unit_points on every cell separately and then calling FEValues::get_function_values or FEValues::get_function_gradients, and for some elements and mappings this is what actually happens internally. For specific combinations of Mapping and FiniteElement realizations, however, there is a much more efficient implementation that avoids the memory allocation and other expensive start-up cost of FEValues. Currently, the functionality is specialized for mappings derived from MappingQGeneric and for finite elements with tensor product structure that work with the Matrix-free infrastructure module. In those cases, the cost implied by this class is similar (or sometimes even somewhat lower) than using FEValues::reinit(cell) followed by FEValues::get_function_gradients.

Definition at line 393 of file fe_point_evaluation.h.

Member Typedef Documentation

◆ value_type

template<int n_components, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components, dim, spacedim, Number >::value_type = typename internal::FEPointEvaluation:: EvaluatorTypeTraits<dim, n_components, Number>::value_type

Definition at line 396 of file fe_point_evaluation.h.

◆ gradient_type

template<int n_components, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components, dim, spacedim, Number >::gradient_type = typename internal::FEPointEvaluation:: EvaluatorTypeTraits<dim, n_components, Number>::gradient_type

Definition at line 398 of file fe_point_evaluation.h.

Constructor & Destructor Documentation

◆ FEPointEvaluation()

template<int n_components, int dim, int spacedim, typename Number >
FEPointEvaluation< n_components, dim, spacedim, Number >::FEPointEvaluation ( const Mapping< dim > &  mapping,
const FiniteElement< dim > &  fe,
const UpdateFlags  update_flags,
const unsigned int  first_selected_component = 0 
)

Constructor.

Parameters
mappingThe Mapping class describing the actual geometry of a cell passed to the evaluate() function.
feThe FiniteElement object that is used for the evaluation, which is typically the same on all cells to be evaluated.
update_flagsSpecify the quantities to be computed by the mapping during the call of reinit(). During evaluate() or integrate(), this data is queried to produce the desired result (e.g., the gradient of a finite element solution).
first_selected_componentFor multi-component FiniteElement objects, this parameter allows to select a range of n_components components starting from this parameter.

Definition at line 682 of file fe_point_evaluation.h.

Member Function Documentation

◆ reinit()

template<int n_components, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components, dim, spacedim, Number >::reinit ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const ArrayView< const Point< dim > > &  unit_points 
)

Set up the mapping information for the given cell, e.g., by computing the Jacobian of the mapping the given points if gradients of the functions are requested.

Parameters
[in]cellAn iterator to the current cell
[in]unit_pointsList of points in the reference locations of the current cell where the FiniteElement object should be evaluated/integrated in the evaluate() and integrate() functions.

Definition at line 759 of file fe_point_evaluation.h.

◆ evaluate()

template<int n_components, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components, dim, spacedim, Number >::evaluate ( const ArrayView< const Number > &  solution_values,
const EvaluationFlags::EvaluationFlags evaluation_flags 
)

This function interpolates the finite element solution, represented by solution_values, on the cell and unit_points passed to reinit().

Parameters
[in]solution_valuesThis array is supposed to contain the unknown values on the element as returned by cell->get_dof_values(global_vector, solution_values).
[in]evaluation_flagsFlags specifying which quantities should be evaluated at the points.

Definition at line 801 of file fe_point_evaluation.h.

◆ integrate()

template<int n_components, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components, dim, spacedim, Number >::integrate ( const ArrayView< Number > &  solution_values,
const EvaluationFlags::EvaluationFlags integration_flags 
)

This function multiplies the quantities passed in by previous submit_value() or submit_gradient() calls by the value or gradient of the test functions, and performs summation over all given points. This is similar to the integration of a bilinear form in terms of the test function, with the difference that this formula does not include a JxW factor. This allows the class to naturally embed point information (e.g. particles) into a finite element formulation. Of course, by multiplication of a JxW information of the data given to submit_value(), the integration can also be represented by this class.

Parameters
[out]solution_valuesThis array will contain the result of the integral, which can be used to during cell->set_dof_values(solution_values, global_vector) or cell->distribute_local_to_global(solution_values, global_vector). Note that for multi-component systems where only some of the components are selected by the present class, the entries not touched by this class will be zeroed out.
[in]integration_flagsFlags specifying which quantities should be integrated at the points.

Definition at line 935 of file fe_point_evaluation.h.

◆ get_value()

template<int n_components, int dim, int spacedim, typename Number >
const FEPointEvaluation< n_components, dim, spacedim, Number >::value_type & FEPointEvaluation< n_components, dim, spacedim, Number >::get_value ( const unsigned int  point_index) const
inline

Return the value at quadrature point number point_index after a call to FEPointEvaluation::evaluate() with EvaluationFlags::value set, or the value that has been stored there with a call to FEPointEvaluation::submit_value(). If the object is vector-valued, a vector-valued return argument is given.

Definition at line 1094 of file fe_point_evaluation.h.

◆ submit_value()

template<int n_components, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components, dim, spacedim, Number >::submit_value ( const value_type value,
const unsigned int  point_index 
)
inline

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

Definition at line 1133 of file fe_point_evaluation.h.

◆ get_gradient()

template<int n_components, int dim, int spacedim, typename Number >
const FEPointEvaluation< n_components, dim, spacedim, Number >::gradient_type & FEPointEvaluation< n_components, dim, spacedim, Number >::get_gradient ( const unsigned int  point_index) const
inline

Return the gradient in real coordinates at the point with index point_index after a call to FEPointEvaluation::evaluate() with EvaluationFlags::gradient set, or the gradient that has been stored there with a call to FEPointEvaluation::submit_gradient(). The gradient in real coordinates is obtained by taking the unit gradient (also accessible via get_unit_gradient()) and applying the inverse Jacobian of the mapping. If the object is vector-valued, a vector-valued return argument is given.

Definition at line 1106 of file fe_point_evaluation.h.

◆ get_unit_gradient()

template<int n_components, int dim, int spacedim, typename Number >
const FEPointEvaluation< n_components, dim, spacedim, Number >::gradient_type & FEPointEvaluation< n_components, dim, spacedim, Number >::get_unit_gradient ( const unsigned int  point_index) const
inline

Return the gradient in unit coordinates at the point with index point_index after a call to FEPointEvaluation::evaluate() with EvaluationFlags::gradient set, or the gradient that has been stored there with a call to FEPointEvaluation::submit_gradient(). If the object is vector-valued, a vector-valued return argument is given. Note that when vectorization is enabled, values from several points are grouped together.

Definition at line 1118 of file fe_point_evaluation.h.

◆ submit_gradient()

template<int n_components, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components, dim, spacedim, Number >::submit_gradient ( const gradient_type gradient,
const unsigned int  point_index 
)
inline

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

Definition at line 1145 of file fe_point_evaluation.h.

◆ jacobian()

template<int n_components, int dim, int spacedim, typename Number >
DerivativeForm< 1, dim, spacedim > FEPointEvaluation< n_components, dim, spacedim, Number >::jacobian ( const unsigned int  point_index) const
inline

Return the Jacobian of the transformation on the current cell with the given point index. Prerequisite: This class needs to be constructed with UpdateFlags containing update_jacobian.

Definition at line 1157 of file fe_point_evaluation.h.

◆ inverse_jacobian()

template<int n_components, int dim, int spacedim, typename Number >
DerivativeForm< 1, spacedim, dim > FEPointEvaluation< n_components, dim, spacedim, Number >::inverse_jacobian ( const unsigned int  point_index) const
inline

Return the inverse of the Jacobian of the transformation on the current cell with the given point index. Prerequisite: This class needs to be constructed with UpdateFlags containing update_inverse_jacobian or update_gradients.

Definition at line 1169 of file fe_point_evaluation.h.

◆ real_point()

template<int n_components, int dim, int spacedim, typename Number >
Point< spacedim > FEPointEvaluation< n_components, dim, spacedim, Number >::real_point ( const unsigned int  point_index) const
inline

Return the position in real coordinates of the given point index among the points passed to reinit().

Definition at line 1183 of file fe_point_evaluation.h.

◆ unit_point()

template<int n_components, int dim, int spacedim, typename Number >
Point< dim > FEPointEvaluation< n_components, dim, spacedim, Number >::unit_point ( const unsigned int  point_index) const
inline

Return the position in unit/reference coordinates of the given point index, i.e., the respective point passed to the reinit() function.

Definition at line 1195 of file fe_point_evaluation.h.

Member Data Documentation

◆ mapping

template<int n_components, int dim, int spacedim = dim, typename Number = double>
SmartPointer<const Mapping<dim, spacedim> > FEPointEvaluation< n_components, dim, spacedim, Number >::mapping
private

Pointer to the Mapping object passed to the constructor.

Definition at line 572 of file fe_point_evaluation.h.

◆ mapping_q_generic

template<int n_components, int dim, int spacedim = dim, typename Number = double>
const MappingQGeneric<dim, spacedim>* FEPointEvaluation< n_components, dim, spacedim, Number >::mapping_q_generic
private

Pointer to MappingQGeneric class that enables the fast path of this class.

Definition at line 578 of file fe_point_evaluation.h.

◆ fe

template<int n_components, int dim, int spacedim = dim, typename Number = double>
SmartPointer<const FiniteElement<dim> > FEPointEvaluation< n_components, dim, spacedim, Number >::fe
private

Pointer to the FiniteElement object passed to the constructor.

Definition at line 583 of file fe_point_evaluation.h.

◆ poly

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::vector<Polynomials::Polynomial<double> > FEPointEvaluation< n_components, dim, spacedim, Number >::poly
private

Description of the 1D polynomial basis for tensor product elements used for the fast path of this class using tensor product evaluators.

Definition at line 589 of file fe_point_evaluation.h.

◆ polynomials_are_hat_functions

template<int n_components, int dim, int spacedim = dim, typename Number = double>
bool FEPointEvaluation< n_components, dim, spacedim, Number >::polynomials_are_hat_functions
private

Store whether the polynomials are linear with nodes at 0 and 1.

Definition at line 594 of file fe_point_evaluation.h.

◆ renumber

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::vector<unsigned int> FEPointEvaluation< n_components, dim, spacedim, Number >::renumber
private

Renumbering between the unknowns of unknowns implied by the FiniteElement class and a lexicographic numbering used for the tensorized code path.

Definition at line 600 of file fe_point_evaluation.h.

◆ solution_renumbered

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::vector<value_type> FEPointEvaluation< n_components, dim, spacedim, Number >::solution_renumbered
private

Temporary array to store the solution_values passed to the evaluate() function in a format compatible with the tensor product evaluators. For vector-valued setups, this array uses a Tensor<1, n_components> type to collect the unknowns for a particular basis function.

Definition at line 608 of file fe_point_evaluation.h.

◆ solution_renumbered_vectorized

template<int n_components, int dim, int spacedim = dim, typename Number = double>
AlignedVector<typename internal::FEPointEvaluation::EvaluatorTypeTraits< dim, n_components, VectorizedArray<Number> >::value_type> FEPointEvaluation< n_components, dim, spacedim, Number >::solution_renumbered_vectorized
private

Temporary array to store a vectorized version of the solution_values computed during integrate() in a format compatible with the tensor product evaluators. For vector-valued setups, this array uses a Tensor<1, n_components, VectorizedArray<Number>> format.

Definition at line 620 of file fe_point_evaluation.h.

◆ values

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::vector<value_type> FEPointEvaluation< n_components, dim, spacedim, Number >::values
private

Temporary array to store the values at the points.

Definition at line 625 of file fe_point_evaluation.h.

◆ unit_gradients

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::vector<gradient_type> FEPointEvaluation< n_components, dim, spacedim, Number >::unit_gradients
private

Temporary array to store the gradients in unit coordinates at the points.

Definition at line 630 of file fe_point_evaluation.h.

◆ gradients

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::vector<gradient_type> FEPointEvaluation< n_components, dim, spacedim, Number >::gradients
private

Temporary array to store the gradients in real coordinates at the points.

Definition at line 635 of file fe_point_evaluation.h.

◆ dofs_per_component

template<int n_components, int dim, int spacedim = dim, typename Number = double>
unsigned int FEPointEvaluation< n_components, dim, spacedim, Number >::dofs_per_component
private

Number of unknowns per component, i.e., number of unique basis functions, for the chosen FiniteElement (or base element).

Definition at line 641 of file fe_point_evaluation.h.

◆ nonzero_shape_function_component

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::vector<std::array<bool, n_components> > FEPointEvaluation< n_components, dim, spacedim, Number >::nonzero_shape_function_component
private

For complicated FiniteElement objects this variable informs us about which unknowns actually carry degrees of freedom in the selected components.

Definition at line 648 of file fe_point_evaluation.h.

◆ update_flags

template<int n_components, int dim, int spacedim = dim, typename Number = double>
UpdateFlags FEPointEvaluation< n_components, dim, spacedim, Number >::update_flags
private

The desired update flags for the evaluation.

Definition at line 653 of file fe_point_evaluation.h.

◆ update_flags_mapping

template<int n_components, int dim, int spacedim = dim, typename Number = double>
UpdateFlags FEPointEvaluation< n_components, dim, spacedim, Number >::update_flags_mapping
private

The update flags specific for the mapping in the fast evaluation path.

Definition at line 658 of file fe_point_evaluation.h.

◆ fe_values

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::shared_ptr<FEValues<dim, spacedim> > FEPointEvaluation< n_components, dim, spacedim, Number >::fe_values
private

The FEValues object underlying the slow evaluation path.

Definition at line 663 of file fe_point_evaluation.h.

◆ mapping_data

template<int n_components, int dim, int spacedim = dim, typename Number = double>
internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> FEPointEvaluation< n_components, dim, spacedim, Number >::mapping_data
private

Array to store temporary data computed by the mapping for the fast evaluation path.

Definition at line 670 of file fe_point_evaluation.h.

◆ unit_points

template<int n_components, int dim, int spacedim = dim, typename Number = double>
std::vector<Point<dim> > FEPointEvaluation< n_components, dim, spacedim, Number >::unit_points
private

The reference points specified at reinit().

Definition at line 675 of file fe_point_evaluation.h.


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