Reference documentation for deal.II version GIT c00406db16 2023-09-28 18:00: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\}}\)
Public Types | Public Member Functions | Static Public Attributes | Private Member Functions | Private Attributes | Static 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 number_type = Number
 
using ScalarNumber = typename internal::VectorizedArrayTrait< Number >::value_type
 
using VectorizedArrayType = typename ::internal::VectorizedArrayTrait< Number >::vectorized_value_type
 
using ETT = typename internal::FEPointEvaluation::EvaluatorTypeTraits< dim, n_components, Number >
 
using value_type = typename ETT::value_type
 
using scalar_value_type = typename ETT::scalar_value_type
 
using vectorized_value_type = typename ETT::vectorized_value_type
 
using gradient_type = typename ETT::gradient_type
 
using interface_vectorized_gradient_type = typename ETT::interface_vectorized_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)
 
 FEPointEvaluation (NonMatching::MappingInfo< dim, spacedim, Number > &mapping_info, const FiniteElement< dim > &fe, const unsigned int first_selected_component=0)
 
 FEPointEvaluation (FEPointEvaluation &other) noexcept
 
 FEPointEvaluation (FEPointEvaluation &&other) noexcept
 
 ~FEPointEvaluation ()
 
void reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim >> &unit_points)
 
void reinit ()
 
void reinit (const unsigned int cell_index)
 
void reinit (const unsigned int cell_index, const unsigned int face_number)
 
template<std::size_t stride_view>
void evaluate (const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
void evaluate (const ArrayView< const ScalarNumber > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
template<std::size_t stride_view>
void integrate (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags)
 
void integrate (const ArrayView< ScalarNumber > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags)
 
template<std::size_t stride_view>
void test_and_sum (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags)
 
void test_and_sum (const ArrayView< ScalarNumber > &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
 
void submit_gradient (const gradient_type &, const unsigned int point_index)
 
DerivativeForm< 1, dim, spacedim, Number > jacobian (const unsigned int point_index) const
 
DerivativeForm< 1, spacedim, dim, Number > inverse_jacobian (const unsigned int point_index) const
 
Number JxW (const unsigned int point_index) const
 
Tensor< 1, spacedim, Number > normal_vector (const unsigned int point_index) const
 
Point< spacedim, Number > real_point (const unsigned int point_index) const
 
Point< dim, Number > unit_point (const unsigned int point_index) const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices () const
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 
static constexpr unsigned int n_components = n_components_
 

Private Member Functions

void setup (const unsigned int first_selected_component)
 
template<bool is_face, bool is_linear>
void do_reinit ()
 
template<bool is_face_path, bool is_linear, std::size_t stride_view>
void prepare_evaluate_fast (const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
template<bool is_face_path, bool is_linear, std::size_t stride_view>
void compute_evaluate_fast (const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags, const unsigned int n_shapes, const unsigned int qb, vectorized_value_type &value, interface_vectorized_gradient_type &gradient)
 
template<bool is_face_path, bool is_linear, std::size_t stride_view>
void evaluate_fast (const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
template<std::size_t stride_view>
void evaluate_slow (const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
template<bool is_face_path, bool is_linear>
void compute_integrate_fast (const EvaluationFlags::EvaluationFlags &integration_flags, const unsigned int n_shapes, const unsigned int qb, const vectorized_value_type value, const interface_vectorized_gradient_type gradient, vectorized_value_type *solution_values_vectorized_linear)
 
template<bool is_face_path, bool is_linear, std::size_t stride_view>
void finish_integrate_fast (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, vectorized_value_type *solution_values_vectorized_linear)
 
template<bool do_JxW, bool is_face_path, bool is_linear, std::size_t stride_view>
void integrate_fast (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags)
 
template<bool do_JxW, std::size_t stride_view>
void integrate_slow (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags)
 
template<bool do_JxW, std::size_t stride_view>
void do_integrate (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags)
 

Private Attributes

const unsigned int n_q_batches
 
const unsigned int n_q_points
 
const unsigned int n_q_points_scalar
 
SmartPointer< const Mapping< dim, spacedim > > mapping
 
SmartPointer< const FiniteElement< dim > > fe
 
std::vector< Polynomials::Polynomial< double > > poly
 
bool use_linear_path
 
std::vector< unsigned int > renumber
 
std::vector< scalar_value_typesolution_renumbered
 
AlignedVector< vectorized_value_typesolution_renumbered_vectorized
 
AlignedVector< ScalarNumberscratch_data_scalar
 
std::vector< value_typevalues
 
std::vector< gradient_typegradients
 
const Point< dim, VectorizedArrayType > * unit_point_ptr
 
const Point< dim - 1, VectorizedArrayType > * unit_point_faces_ptr
 
const Point< spacedim, Number > * real_point_ptr
 
const DerivativeForm< 1, dim, spacedim, Number > * jacobian_ptr
 
const DerivativeForm< 1, spacedim, dim, Number > * inverse_jacobian_ptr
 
const Tensor< 1, spacedim, Number > * normal_ptr
 
const Number * JxW_ptr
 
internal::MatrixFreeFunctions::GeometryType cell_type
 
unsigned int dofs_per_component
 
unsigned int dofs_per_component_face
 
bool use_face_path
 
internal::MatrixFreeFunctions::ShapeInfo< ScalarNumbershape_info
 
unsigned int component_in_base_element
 
std::vector< std::array< bool, n_components > > nonzero_shape_function_component
 
const UpdateFlags update_flags
 
std::shared_ptr< FEValues< dim, spacedim > > fe_values
 
std::unique_ptr< NonMatching::MappingInfo< dim, spacedim, Number > > mapping_info_on_the_fly
 
SmartPointer< NonMatching::MappingInfo< dim, spacedim, Number > > mapping_info
 
unsigned int current_cell_index
 
unsigned int current_face_number
 
bool fast_path
 
boost::signals2::connection connection_is_reinitialized
 
bool is_reinitialized
 
AlignedVector<::ndarray< VectorizedArrayType, 2, dim > > shapes
 
AlignedVector<::ndarray< VectorizedArrayType, 2, dim - 1 > > shapes_faces
 

Static Private Attributes

static constexpr std::size_t n_lanes_user_interface
 
static constexpr std::size_t n_lanes_internal
 
static constexpr std::size_t stride
 

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 MappingQ and MappingCartesian 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 718 of file fe_point_evaluation.h.

Member Typedef Documentation

◆ number_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components_, dim, spacedim, Number >::number_type = Number

Definition at line 724 of file fe_point_evaluation.h.

◆ ScalarNumber

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components_, dim, spacedim, Number >::ScalarNumber = typename internal::VectorizedArrayTrait<Number>::value_type

Definition at line 726 of file fe_point_evaluation.h.

◆ VectorizedArrayType

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components_, dim, spacedim, Number >::VectorizedArrayType = typename ::internal::VectorizedArrayTrait< Number>::vectorized_value_type

Definition at line 728 of file fe_point_evaluation.h.

◆ ETT

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

Definition at line 730 of file fe_point_evaluation.h.

◆ value_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components_, dim, spacedim, Number >::value_type = typename ETT::value_type

Definition at line 732 of file fe_point_evaluation.h.

◆ scalar_value_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components_, dim, spacedim, Number >::scalar_value_type = typename ETT::scalar_value_type

Definition at line 733 of file fe_point_evaluation.h.

◆ vectorized_value_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components_, dim, spacedim, Number >::vectorized_value_type = typename ETT::vectorized_value_type

Definition at line 734 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 ETT::gradient_type

Definition at line 735 of file fe_point_evaluation.h.

◆ interface_vectorized_gradient_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEPointEvaluation< n_components_, dim, spacedim, Number >::interface_vectorized_gradient_type = typename ETT::interface_vectorized_gradient_type

Definition at line 736 of file fe_point_evaluation.h.

Constructor & Destructor Documentation

◆ FEPointEvaluation() [1/4]

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 1430 of file fe_point_evaluation.h.

◆ FEPointEvaluation() [2/4]

template<int n_components_, int dim, int spacedim, typename Number >
FEPointEvaluation< n_components_, dim, spacedim, Number >::FEPointEvaluation ( NonMatching::MappingInfo< dim, spacedim, Number > &  mapping_info,
const FiniteElement< dim > &  fe,
const unsigned int  first_selected_component = 0 
)

Constructor to make the present class able to re-use the geometry data also used by other FEPointEvaluation objects.

Parameters
mapping_infoThe MappingInfo class describes the geometry-related data for evaluating finite-element solutions. This object enables to construct such an object on the outside, possibly re-using it between several objects or between several calls to the same cell and unit points.
feThe FiniteElement object that is used for the evaluation, which is typically the same on all cells to be evaluated.
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 1458 of file fe_point_evaluation.h.

◆ FEPointEvaluation() [3/4]

template<int n_components_, int dim, int spacedim, typename Number >
FEPointEvaluation< n_components_, dim, spacedim, Number >::FEPointEvaluation ( FEPointEvaluation< n_components_, dim, spacedim, Number > &  other)
noexcept

Copy constructor.

Definition at line 1483 of file fe_point_evaluation.h.

◆ FEPointEvaluation() [4/4]

template<int n_components_, int dim, int spacedim, typename Number >
FEPointEvaluation< n_components_, dim, spacedim, Number >::FEPointEvaluation ( FEPointEvaluation< n_components_, dim, spacedim, Number > &&  other)
noexcept

Move constructor.

Definition at line 1525 of file fe_point_evaluation.h.

◆ ~FEPointEvaluation()

template<int n_components_, int dim, int spacedim, typename Number >
FEPointEvaluation< n_components_, dim, spacedim, Number >::~FEPointEvaluation

Destructor.

Definition at line 1562 of file fe_point_evaluation.h.

Member Function Documentation

◆ reinit() [1/4]

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 
)
inline

Set up the mapping information for the given cell, e.g., by computing the Jacobian of the mapping for 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 1657 of file fe_point_evaluation.h.

◆ reinit() [2/4]

template<int n_components_, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components_, dim, spacedim, Number >::reinit
inline

Reinitialize the evaluator to point to the correct precomputed mapping of the single cell in the MappingInfo object.

Definition at line 1687 of file fe_point_evaluation.h.

◆ reinit() [3/4]

template<int n_components_, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components_, dim, spacedim, Number >::reinit ( const unsigned int  cell_index)
inline

Reinitialize the evaluator to point to the correct precomputed mapping of the cell in the MappingInfo object.

Definition at line 1702 of file fe_point_evaluation.h.

◆ reinit() [4/4]

template<int n_components_, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components_, dim, spacedim, Number >::reinit ( const unsigned int  cell_index,
const unsigned int  face_number 
)
inline

Reinitialize the evaluator to point to the correct precomputed mapping of the face in the MappingInfo object.

Definition at line 1737 of file fe_point_evaluation.h.

◆ evaluate() [1/2]

template<int n_components_, int dim, int spacedim, typename Number >
template<std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::evaluate ( const StridedArrayView< const ScalarNumber, stride_view > &  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 read out by FEEvaluation::read_dof_values(global_vector).
[in]evaluation_flagsFlags specifying which quantities should be evaluated at the points.

Definition at line 2233 of file fe_point_evaluation.h.

◆ evaluate() [2/2]

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::evaluate ( const ArrayView< const ScalarNumber > &  solution_values,
const EvaluationFlags::EvaluationFlags evaluation_flags 
)
inline

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 863 of file fe_point_evaluation.h.

◆ integrate() [1/2]

template<int n_components_, int dim, int spacedim, typename Number >
template<std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::integrate ( const StridedArrayView< ScalarNumber, stride_view > &  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 multiplied be the Jacobian determinant times the quadrature weight (JxW).

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

Definition at line 2690 of file fe_point_evaluation.h.

◆ integrate() [2/2]

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::integrate ( const ArrayView< ScalarNumber > &  solution_values,
const EvaluationFlags::EvaluationFlags integration_flags 
)
inline

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 multiplied be the Jacobian determinant times the quadrature weight (JxW).

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 in solution_values not touched by this class will be set to zero.
[in]integration_flagsFlags specifying which quantities should be integrated at the points.

Definition at line 914 of file fe_point_evaluation.h.

◆ test_and_sum() [1/2]

template<int n_components_, int dim, int spacedim, typename Number >
template<std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::test_and_sum ( const StridedArrayView< ScalarNumber, stride_view > &  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 (in contrast to the integrate function of this class). This allows the class to naturally embed point information (e.g. particles) into a finite element formulation.

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

Definition at line 2702 of file fe_point_evaluation.h.

◆ test_and_sum() [2/2]

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::test_and_sum ( const ArrayView< ScalarNumber > &  solution_values,
const EvaluationFlags::EvaluationFlags integration_flags 
)
inline

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 (in contrast to the integrate function of this class). This allows the class to naturally embed point information (e.g. particles) into a finite element formulation.

Parameters
[out]solution_valuesThis array will contain the result of the integral, which can be used 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 in solution_values not touched by this class will be set to zero.
[in]integration_flagsFlags specifying which quantities should be integrated at the points.

Definition at line 973 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::values 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 2714 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 2737 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::gradients 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 2726 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 2749 of file fe_point_evaluation.h.

◆ jacobian()

template<int n_components_, int dim, int spacedim, typename Number >
DerivativeForm< 1, dim, spacedim, Number > 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 2761 of file fe_point_evaluation.h.

◆ inverse_jacobian()

template<int n_components_, int dim, int spacedim, typename Number >
DerivativeForm< 1, spacedim, dim, Number > 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 2779 of file fe_point_evaluation.h.

◆ JxW()

template<int n_components_, int dim, int spacedim, typename Number >
Number FEPointEvaluation< n_components_, dim, spacedim, Number >::JxW ( const unsigned int  point_index) const
inline

Return the Jacobian determinant multiplied by the quadrature weight. This class or the MappingInfo object passed to this function needs to be constructed with UpdateFlags containing update_JxW_values.

Definition at line 2798 of file fe_point_evaluation.h.

◆ normal_vector()

template<int n_components_, int dim, int spacedim, typename Number >
Tensor< 1, spacedim, Number > FEPointEvaluation< n_components_, dim, spacedim, Number >::normal_vector ( const unsigned int  point_index) const
inline

Return the normal vector. This class or the MappingInfo object passed to this function needs to be constructed with UpdateFlags containing update_normal_vectors.

Definition at line 2813 of file fe_point_evaluation.h.

◆ real_point()

template<int n_components_, int dim, int spacedim, typename Number >
Point< spacedim, Number > 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 2828 of file fe_point_evaluation.h.

◆ unit_point()

template<int n_components_, int dim, int spacedim, typename Number >
Point< dim, Number > 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 2843 of file fe_point_evaluation.h.

◆ quadrature_point_indices()

template<int n_components_, int dim, int spacedim, typename Number >
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FEPointEvaluation< n_components_, dim, spacedim, Number >::quadrature_point_indices
inline

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.

Definition at line 2859 of file fe_point_evaluation.h.

◆ setup()

template<int n_components_, int dim, int spacedim, typename Number >
void FEPointEvaluation< n_components_, dim, spacedim, Number >::setup ( const unsigned int  first_selected_component)
private

Common setup function for both constructors. Does the setup for both fast and slow path.

Parameters
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 1571 of file fe_point_evaluation.h.

◆ do_reinit()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool is_face, bool is_linear>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::do_reinit
inlineprivate

Shared functionality of all reinit() functions. Resizes data fields and precomputes the shapes vector, holding the evaluation of 1D basis functions of tensor product polynomials, if necessary.

Definition at line 1755 of file fe_point_evaluation.h.

◆ prepare_evaluate_fast()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool is_face_path, bool is_linear, std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::prepare_evaluate_fast ( const StridedArrayView< const ScalarNumber, stride_view > &  solution_values,
const EvaluationFlags::EvaluationFlags evaluation_flags 
)
inlineprivate

Resizes necessary data fields, reads in and renumbers solution values. Interpolates onto face if face path is selected.

Definition at line 1868 of file fe_point_evaluation.h.

◆ compute_evaluate_fast()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool is_face_path, bool is_linear, std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::compute_evaluate_fast ( const StridedArrayView< const ScalarNumber, stride_view > &  solution_values,
const EvaluationFlags::EvaluationFlags evaluation_flags,
const unsigned int  n_shapes,
const unsigned int  qb,
vectorized_value_type value,
interface_vectorized_gradient_type gradient 
)
inlineprivate

Evaluates the actual interpolation on the cell or face for a quadrature batch.

Definition at line 1938 of file fe_point_evaluation.h.

◆ evaluate_fast()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool is_face_path, bool is_linear, std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::evaluate_fast ( const StridedArrayView< const ScalarNumber, stride_view > &  solution_values,
const EvaluationFlags::EvaluationFlags evaluation_flags 
)
inlineprivate

Fast path of the evaluate function.

Definition at line 2089 of file fe_point_evaluation.h.

◆ evaluate_slow()

template<int n_components_, int dim, int spacedim, typename Number >
template<std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::evaluate_slow ( const StridedArrayView< const ScalarNumber, stride_view > &  solution_values,
const EvaluationFlags::EvaluationFlags evaluation_flags 
)
inlineprivate

Slow path of the evaluate function using FEValues.

Definition at line 2150 of file fe_point_evaluation.h.

◆ compute_integrate_fast()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool is_face_path, bool is_linear>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::compute_integrate_fast ( const EvaluationFlags::EvaluationFlags integration_flags,
const unsigned int  n_shapes,
const unsigned int  qb,
const vectorized_value_type  value,
const interface_vectorized_gradient_type  gradient,
vectorized_value_type solution_values_vectorized_linear 
)
inlineprivate

Integrates the product of the data passed in by submit_value() and submit_gradient() with the values or gradients of test functions on the cell or face for a given quadrature batch.

Definition at line 2276 of file fe_point_evaluation.h.

◆ finish_integrate_fast()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool is_face_path, bool is_linear, std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::finish_integrate_fast ( const StridedArrayView< ScalarNumber, stride_view > &  solution_values,
const EvaluationFlags::EvaluationFlags integration_flags,
vectorized_value_type solution_values_vectorized_linear 
)
inlineprivate

Addition across the lanes of VectorizedArray as accumulated by the compute_integrate_fast_function(), writing the sum into the result vector. Applies face contributions to cell contributions for face path.

Definition at line 2393 of file fe_point_evaluation.h.

◆ integrate_fast()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool do_JxW, bool is_face_path, bool is_linear, std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::integrate_fast ( const StridedArrayView< ScalarNumber, stride_view > &  solution_values,
const EvaluationFlags::EvaluationFlags integration_flags 
)
inlineprivate

Fast path of the integrate function.

Definition at line 2476 of file fe_point_evaluation.h.

◆ integrate_slow()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool do_JxW, std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::integrate_slow ( const StridedArrayView< ScalarNumber, stride_view > &  solution_values,
const EvaluationFlags::EvaluationFlags integration_flags 
)
inlineprivate

Slow path of the integrate function using FEValues.

Definition at line 2558 of file fe_point_evaluation.h.

◆ do_integrate()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool do_JxW, std::size_t stride_view>
void FEPointEvaluation< n_components_, dim, spacedim, Number >::do_integrate ( const StridedArrayView< ScalarNumber, stride_view > &  solution_values,
const EvaluationFlags::EvaluationFlags integration_flags 
)
private

Implementation of the integrate/test_and_sum function.

Definition at line 2635 of file fe_point_evaluation.h.

Member Data Documentation

◆ dimension

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
constexpr unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::dimension = dim
staticconstexpr

Definition at line 721 of file fe_point_evaluation.h.

◆ n_components

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
constexpr unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::n_components = n_components_
staticconstexpr

Definition at line 722 of file fe_point_evaluation.h.

◆ n_lanes_user_interface

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
constexpr std::size_t FEPointEvaluation< n_components_, dim, spacedim, Number >::n_lanes_user_interface
staticconstexprprivate
Initial value:

Definition at line 1081 of file fe_point_evaluation.h.

◆ n_lanes_internal

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
constexpr std::size_t FEPointEvaluation< n_components_, dim, spacedim, Number >::n_lanes_internal
staticconstexprprivate

◆ stride

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
constexpr std::size_t FEPointEvaluation< n_components_, dim, spacedim, Number >::stride
staticconstexprprivate
Initial value:

Definition at line 1085 of file fe_point_evaluation.h.

◆ n_q_batches

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::n_q_batches
private

Number of quadrature batches of the current cell/face.

Definition at line 1210 of file fe_point_evaluation.h.

◆ n_q_points

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::n_q_points
private

Number of quadrature points/batches of the current cell/face.

Definition at line 1215 of file fe_point_evaluation.h.

◆ n_q_points_scalar

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::n_q_points_scalar
private

Number of quadrature points of the current cell/face.

Definition at line 1220 of file fe_point_evaluation.h.

◆ 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 1225 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 1230 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 1236 of file fe_point_evaluation.h.

◆ use_linear_path

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
bool FEPointEvaluation< n_components_, dim, spacedim, Number >::use_linear_path
private

Store whether the linear path should be used.

Definition at line 1241 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 1247 of file fe_point_evaluation.h.

◆ solution_renumbered

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
std::vector<scalar_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, ScalarNumber> type to collect the unknowns for a particular basis function.

Definition at line 1256 of file fe_point_evaluation.h.

◆ solution_renumbered_vectorized

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
AlignedVector<vectorized_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, VectorizedArrayType> format.

Definition at line 1264 of file fe_point_evaluation.h.

◆ scratch_data_scalar

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
AlignedVector<ScalarNumber> FEPointEvaluation< n_components_, dim, spacedim, Number >::scratch_data_scalar
private

Temporary array for the use_face_path path (scalar).

Definition at line 1269 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 1274 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 1279 of file fe_point_evaluation.h.

◆ unit_point_ptr

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const Point<dim, VectorizedArrayType>* FEPointEvaluation< n_components_, dim, spacedim, Number >::unit_point_ptr
private

Pointer to first unit point batch of current cell/face from MappingInfo, set internally during do_reinit().

Definition at line 1285 of file fe_point_evaluation.h.

◆ unit_point_faces_ptr

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const Point<dim - 1, VectorizedArrayType>* FEPointEvaluation< n_components_, dim, spacedim, Number >::unit_point_faces_ptr
private

Pointer to first unit point batch of current face from MappingInfo, set internally during do_reinit(). Needed for use_face_path path.

Definition at line 1291 of file fe_point_evaluation.h.

◆ real_point_ptr

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const Point<spacedim, Number>* FEPointEvaluation< n_components_, dim, spacedim, Number >::real_point_ptr
private

Pointer to real point of first quadrature point of current cell/face from MappingInfo, set internally during do_reinit().

Definition at line 1297 of file fe_point_evaluation.h.

◆ jacobian_ptr

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const DerivativeForm<1, dim, spacedim, Number>* FEPointEvaluation< n_components_, dim, spacedim, Number >::jacobian_ptr
private

Pointer to Jacobian of first quadrature point of current cell/face from MappingInfo, set internally during do_reinit().

Definition at line 1303 of file fe_point_evaluation.h.

◆ inverse_jacobian_ptr

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const DerivativeForm<1, spacedim, dim, Number>* FEPointEvaluation< n_components_, dim, spacedim, Number >::inverse_jacobian_ptr
private

Pointer to inverse Jacobian of first quadrature point of current cell/face from MappingInfo, set internally during do_reinit().

Definition at line 1309 of file fe_point_evaluation.h.

◆ normal_ptr

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const Tensor<1, spacedim, Number>* FEPointEvaluation< n_components_, dim, spacedim, Number >::normal_ptr
private

Pointer to normal vector of first quadrature point of current cell/face from MappingInfo, set internally during do_reinit().

Definition at line 1315 of file fe_point_evaluation.h.

◆ JxW_ptr

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
const Number* FEPointEvaluation< n_components_, dim, spacedim, Number >::JxW_ptr
private

Pointer to Jacobian determinant times quadrature weight of first quadrature point of current cell/face from MappingInfo, set internally during do_reinit().

Definition at line 1322 of file fe_point_evaluation.h.

◆ cell_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
internal::MatrixFreeFunctions::GeometryType FEPointEvaluation< n_components_, dim, spacedim, Number >::cell_type
private

Cell type describing the geometry of the cell and compression of jacobians.

Definition at line 1327 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 1333 of file fe_point_evaluation.h.

◆ dofs_per_component_face

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::dofs_per_component_face
private

Number of unknowns per component, i.e., number of unique basis functions, for a restriction to the face of the chosen FiniteElement (or base element). This means a (dim-1)-dimensional basis.

Definition at line 1340 of file fe_point_evaluation.h.

◆ use_face_path

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
bool FEPointEvaluation< n_components_, dim, spacedim, Number >::use_face_path
private

Bool indicating if use_face_path path should be chosen. Set during do_reinit().

Definition at line 1346 of file fe_point_evaluation.h.

◆ shape_info

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
internal::MatrixFreeFunctions::ShapeInfo<ScalarNumber> FEPointEvaluation< n_components_, dim, spacedim, Number >::shape_info
private

Scalar ShapeInfo object needed for use_face_path path.

Definition at line 1351 of file fe_point_evaluation.h.

◆ component_in_base_element

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::component_in_base_element
private

The first selected component in the active base element.

Definition at line 1356 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 1363 of file fe_point_evaluation.h.

◆ update_flags

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

The desired update flags for the evaluation.

Definition at line 1368 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 1373 of file fe_point_evaluation.h.

◆ mapping_info_on_the_fly

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
std::unique_ptr<NonMatching::MappingInfo<dim, spacedim, Number> > FEPointEvaluation< n_components_, dim, spacedim, Number >::mapping_info_on_the_fly
private

Pointer to mapping info on the fly computed during reinit.

Definition at line 1379 of file fe_point_evaluation.h.

◆ mapping_info

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
SmartPointer<NonMatching::MappingInfo<dim, spacedim, Number> > FEPointEvaluation< n_components_, dim, spacedim, Number >::mapping_info
private

Pointer to currently used mapping info (either on the fly or external precomputed).

Definition at line 1385 of file fe_point_evaluation.h.

◆ current_cell_index

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::current_cell_index
private

The current cell index to access mapping data from mapping info.

Definition at line 1390 of file fe_point_evaluation.h.

◆ current_face_number

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
unsigned int FEPointEvaluation< n_components_, dim, spacedim, Number >::current_face_number
private

The current face number to access mapping data from mapping info.

Definition at line 1395 of file fe_point_evaluation.h.

◆ fast_path

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
bool FEPointEvaluation< n_components_, dim, spacedim, Number >::fast_path
private

Bool indicating if fast path is chosen.

Definition at line 1400 of file fe_point_evaluation.h.

◆ connection_is_reinitialized

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
boost::signals2::connection FEPointEvaluation< n_components_, dim, spacedim, Number >::connection_is_reinitialized
private

Connection to NonMatching::MappingInfo to check whether mapping data has been invalidated.

Definition at line 1406 of file fe_point_evaluation.h.

◆ is_reinitialized

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
bool FEPointEvaluation< n_components_, dim, spacedim, Number >::is_reinitialized
private

Bool indicating if class is reinitialized and data vectors a resized.

Definition at line 1411 of file fe_point_evaluation.h.

◆ shapes

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
AlignedVector<::ndarray<VectorizedArrayType, 2, dim> > FEPointEvaluation< n_components_, dim, spacedim, Number >::shapes
private

Vector containing tensor product shape functions evaluated (during reinit()) at the vectorized unit points.

Definition at line 1417 of file fe_point_evaluation.h.

◆ shapes_faces

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
AlignedVector<::ndarray<VectorizedArrayType, 2, dim - 1> > FEPointEvaluation< n_components_, dim, spacedim, Number >::shapes_faces
private

Vector containing tensor product shape functions evaluated (during reinit()) at the vectorized unit points on faces.

Definition at line 1423 of file fe_point_evaluation.h.


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