Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
FEFacePointEvaluation< n_components_, dim, spacedim, Number > Class Template Reference

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

Inheritance diagram for FEFacePointEvaluation< n_components_, dim, spacedim, Number >:

Public Types

using number_type = Number
 
using ScalarNumber
 
using VectorizedArrayType
 
using ETT
 
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 unit_gradient_type = typename ETT::unit_gradient_type
 
using gradient_type = typename ETT::real_gradient_type
 
using interface_vectorized_unit_gradient_type
 

Public Member Functions

 FEFacePointEvaluation (NonMatching::MappingInfo< dim, spacedim, Number > &mapping_info, const FiniteElement< dim, spacedim > &fe, const bool is_interior=true, const unsigned int first_selected_component=0)
 
void reinit (const unsigned int cell_index, const unsigned int face_number)
 
void reinit (const unsigned int face_index)
 
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, const bool sum_into_values=false)
 
void integrate (const ArrayView< ScalarNumber > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
 
template<std::size_t stride_view>
void test_and_sum (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
 
void test_and_sum (const ArrayView< ScalarNumber > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
 
template<int stride_face_dof = VectorizedArrayType::size()>
void evaluate_in_face (const ScalarNumber *face_dof_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
template<int stride_face_dof = VectorizedArrayType::size()>
void integrate_in_face (ScalarNumber *face_dof_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values=false)
 
Tensor< 1, spacedim, Number > normal_vector (const unsigned int point_index) const
 
const value_type get_normal_derivative (const unsigned int point_index) const
 
void submit_normal_derivative (const value_type &, const unsigned int point_index)
 
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)
 
double get_divergence (const unsigned int point_index) const
 
void submit_divergence (const double &value, const unsigned int point_index)
 
DerivativeForm< 1, dim, spacedim, double > jacobian (const unsigned int point_index) const
 
DerivativeForm< 1, spacedim, dim, double > inverse_jacobian (const unsigned int point_index) const
 
double JxW (const unsigned int point_index) const
 
Point< spacedim, double > real_point (const unsigned int point_index) const
 
Point< spacedim, double > quadrature_point (const unsigned int point_index) const
 
Point< dim, double > unit_point (const unsigned int point_index) const
 
scalar_value_type integrate_value () const
 
std_cxx20::ranges::iota_view< unsigned int, unsigned intquadrature_point_indices () const
 
unsigned int n_active_entries_per_quadrature_batch (unsigned int q)
 

Static Public Attributes

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

Protected Member Functions

void setup (const unsigned int first_selected_component)
 
void do_reinit ()
 

Protected 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, spacedim > > fe
 
std::vector< Polynomials::Polynomial< double > > poly
 
bool use_linear_path
 
std::vector< unsigned intrenumber
 
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, double > * real_point_ptr
 
const DerivativeForm< 1, dim, spacedim, double > * jacobian_ptr
 
const DerivativeForm< 1, spacedim, dim, double > * inverse_jacobian_ptr
 
const Tensor< 1, spacedim, double > * normal_ptr
 
const double * JxW_ptr
 
internal::MatrixFreeFunctions::GeometryType cell_type
 
unsigned int dofs_per_component
 
unsigned int dofs_per_component_face
 
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, double > > mapping_info_on_the_fly
 
SmartPointer< NonMatching::MappingInfo< dim, spacedim, double > > mapping_info
 
unsigned int current_cell_index
 
unsigned int current_face_number
 
bool fast_path
 
bool must_reinitialize_pointers
 
AlignedVector<::ndarray< VectorizedArrayType, 2, dim > > shapes
 
AlignedVector<::ndarray< VectorizedArrayType, 2, dim - 1 > > shapes_faces
 
const bool is_interior
 

Private Member Functions

template<bool is_linear, std::size_t stride_view>
void do_evaluate (const StridedArrayView< const ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
template<bool do_JxW, bool is_linear, std::size_t stride_view>
void do_integrate (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values)
 
template<bool is_linear, int stride_face_dof>
void do_evaluate_in_face (const ScalarNumber *face_dof_values, const EvaluationFlags::EvaluationFlags &evaluation_flags)
 
template<bool do_JxW, bool is_linear, int stride_face_dof>
void do_integrate_in_face (ScalarNumber *face_dof_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values)
 

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 FEFacePointEvaluation< n_components_, dim, spacedim, Number >

This class provides an interface to the evaluation of interpolated solution values and gradients on faces on arbitrary reference point positions. These points can change from face to face, both with respect to their quantity as well to the location. A typical use case is evaluations on non-matching grids.

The use of this class is similar to FEEvaluation: In the constructor, a reference to a NonMatching::MappingInfo object is passed, where the quadrature points in reference position is stored together with the mapping information. The class is then reinitialized to a cell by calling FEFacePointEvaluation::reinit(face_index) or FEFacePointEvaluation::reinit(cell_index, face_number). Then, upon call to evaluate() or integrate(), the user can compute information at the given points. Eventually, the access functions get_value() or get_gradient() allow to query this information at a specific point index.

Definition at line 1534 of file fe_point_evaluation.h.

Member Typedef Documentation

◆ number_type

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

Definition at line 1541 of file fe_point_evaluation.h.

◆ ScalarNumber

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEFacePointEvaluation< n_components_, dim, spacedim, Number >::ScalarNumber

◆ VectorizedArrayType

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

Definition at line 1545 of file fe_point_evaluation.h.

◆ ETT

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

Definition at line 1547 of file fe_point_evaluation.h.

◆ value_type

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

Definition at line 1549 of file fe_point_evaluation.h.

◆ scalar_value_type

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

Definition at line 1550 of file fe_point_evaluation.h.

◆ vectorized_value_type

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

Definition at line 1551 of file fe_point_evaluation.h.

◆ unit_gradient_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEFacePointEvaluation< n_components_, dim, spacedim, Number >::unit_gradient_type = typename ETT::unit_gradient_type

Definition at line 1552 of file fe_point_evaluation.h.

◆ gradient_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEFacePointEvaluation< n_components_, dim, spacedim, Number >::gradient_type = typename ETT::real_gradient_type

Definition at line 1553 of file fe_point_evaluation.h.

◆ interface_vectorized_unit_gradient_type

template<int n_components_, int dim, int spacedim = dim, typename Number = double>
using FEFacePointEvaluation< n_components_, dim, spacedim, Number >::interface_vectorized_unit_gradient_type
Initial value:
typename ETT::interface_vectorized_unit_gradient_type

Definition at line 1554 of file fe_point_evaluation.h.

Constructor & Destructor Documentation

◆ FEFacePointEvaluation()

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

Constructor. Allows to select if interior or exterior face is selected.

Definition at line 3208 of file fe_point_evaluation.h.

Member Function Documentation

◆ reinit() [1/2]

template<int n_components_, int dim, int spacedim, typename Number >
void FEFacePointEvaluation< 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. Used in element-centric loops (ECL).

Definition at line 3225 of file fe_point_evaluation.h.

◆ reinit() [2/2]

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

Reinitialize the evaluator to point to the correct precomputed mapping of the face in the MappingInfo object. Used in face-centric loops (FCL).

Definition at line 3243 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 FEFacePointEvaluation< 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 3262 of file fe_point_evaluation.h.

◆ evaluate() [2/2]

template<int n_components_, int dim, int spacedim, typename Number >
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::evaluate ( const ArrayView< const ScalarNumber > & 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 3290 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 FEFacePointEvaluation< n_components_, dim, spacedim, Number >::integrate ( const StridedArrayView< ScalarNumber, stride_view > & solution_values,
const EvaluationFlags::EvaluationFlags & integration_flags,
const bool sum_into_values = false )

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.
[in]sum_into_valuesFlag specifying if the integrated values should be summed into the solution values. For the default value sum_into_values=false every value of solution_values is zeroed out.

Definition at line 3358 of file fe_point_evaluation.h.

◆ integrate() [2/2]

template<int n_components_, int dim, int spacedim, typename Number >
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::integrate ( const ArrayView< ScalarNumber > & solution_values,
const EvaluationFlags::EvaluationFlags & integration_flags,
const bool sum_into_values = false )

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.
[in]sum_into_valuesFlag specifying if the integrated values should be summed into the solution values. For the default value sum_into_values=false every value of solution_values is zeroed out.

Definition at line 3395 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 FEFacePointEvaluation< n_components_, dim, spacedim, Number >::test_and_sum ( const StridedArrayView< ScalarNumber, stride_view > & solution_values,
const EvaluationFlags::EvaluationFlags & integration_flags,
const bool sum_into_values = false )

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.
[in]sum_into_valuesFlag specifying if the integrated values should be summed into the solution values. For the default value sum_into_values=false every value of solution_values is zeroed out.

Definition at line 3411 of file fe_point_evaluation.h.

◆ test_and_sum() [2/2]

template<int n_components_, int dim, int spacedim, typename Number >
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::test_and_sum ( const ArrayView< ScalarNumber > & solution_values,
const EvaluationFlags::EvaluationFlags & integration_flags,
const bool sum_into_values = false )

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.
[in]sum_into_valuesFlag specifying if the integrated values should be summed into the solution values. For the default value sum_into_values=false every value of solution_values is zeroed out.

Definition at line 3448 of file fe_point_evaluation.h.

◆ evaluate_in_face()

template<int n_components_, int dim, int spacedim, typename Number >
template<int stride_face_dof>
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::evaluate_in_face ( const ScalarNumber * face_dof_values,
const EvaluationFlags::EvaluationFlags & evaluation_flags )

Evaluate values and gradients in face for the selected face (lane) of the batch. Default stride into the face dofs is width of VectorizedArray<selected_floating_point_type> which is the default vectorization over faces for FEFaceEvaluation.

Definition at line 3556 of file fe_point_evaluation.h.

◆ integrate_in_face()

template<int n_components_, int dim, int spacedim, typename Number >
template<int stride_face_dof>
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::integrate_in_face ( ScalarNumber * face_dof_values,
const EvaluationFlags::EvaluationFlags & integration_flags,
const bool sum_into_values = false )

Integrate values and gradients in face for the selected face (lane) of the batch. Default stride into the face dofs is width of VectorizedArray<selected_floating_point_type> which is the default vectorization over faces for FEFaceEvaluation.

Definition at line 3729 of file fe_point_evaluation.h.

◆ normal_vector()

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

◆ get_normal_derivative()

template<int n_components_, int dim, int spacedim, typename Number >
const FEFacePointEvaluation< n_components_, dim, spacedim, Number >::value_type FEFacePointEvaluation< n_components_, dim, spacedim, Number >::get_normal_derivative ( const unsigned int point_index) const
inline

Return the normal derivative in real coordinates at the point with index point_index after a call to FEFacePointEvaluation::evaluate() with EvaluationFlags::gradients set.

Definition at line 3940 of file fe_point_evaluation.h.

◆ submit_normal_derivative()

template<int n_components_, int dim, int spacedim, typename Number >
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::submit_normal_derivative ( const value_type & value,
const unsigned int point_index )
inline

Write a contribution that is tested by the normal derivative to the field containing the values on points with the given point_index. Access to the same field as through set_gradient()/get_gradient.

Definition at line 3951 of file fe_point_evaluation.h.

◆ do_evaluate()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool is_linear, std::size_t stride_view>
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::do_evaluate ( const StridedArrayView< const ScalarNumber, stride_view > & solution_values,
const EvaluationFlags::EvaluationFlags & evaluation_flags )
private

Definition at line 3304 of file fe_point_evaluation.h.

◆ do_integrate()

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

Definition at line 3464 of file fe_point_evaluation.h.

◆ do_evaluate_in_face()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool is_linear, int stride_face_dof>
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::do_evaluate_in_face ( const ScalarNumber * face_dof_values,
const EvaluationFlags::EvaluationFlags & evaluation_flags )
inlineprivate

Actually does the evaluation templated on the chosen code path (linear or higher order).

Definition at line 3573 of file fe_point_evaluation.h.

◆ do_integrate_in_face()

template<int n_components_, int dim, int spacedim, typename Number >
template<bool do_JxW, bool is_linear, int stride_face_dof>
void FEFacePointEvaluation< n_components_, dim, spacedim, Number >::do_integrate_in_face ( ScalarNumber * face_dof_values,
const EvaluationFlags::EvaluationFlags & integration_flags,
const bool sum_into_values )
inlineprivate

Actually does the integration templated on the chosen code path (linear or higher order).

Definition at line 3749 of file fe_point_evaluation.h.

◆ get_value()

const FEPointEvaluationBase< n_components_, dim, spacedim, double >::value_type & FEPointEvaluationBase< n_components_, dim, spacedim, double >::get_value ( const unsigned int point_index) const
inlineinherited

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 FEPointEvaluationBase::submit_value(). If the object is vector-valued, a vector-valued return argument is given.

Definition at line 712 of file fe_point_evaluation.h.

◆ submit_value()

void FEPointEvaluationBase< n_components_, dim, spacedim, double >::submit_value ( const value_type & value,
const unsigned int point_index )
inlineinherited

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

◆ get_gradient()

const FEPointEvaluationBase< n_components_, dim, spacedim, double >::gradient_type & FEPointEvaluationBase< n_components_, dim, spacedim, double >::get_gradient ( const unsigned int point_index) const
inlineinherited

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 FEPointEvaluationBase::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 735 of file fe_point_evaluation.h.

◆ submit_gradient()

void FEPointEvaluationBase< n_components_, dim, spacedim, double >::submit_gradient ( const gradient_type & gradient,
const unsigned int point_index )
inlineinherited

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

◆ get_divergence()

double FEPointEvaluationBase< n_components_, dim, spacedim, double >::get_divergence ( const unsigned int point_index) const
inlineinherited

Return the divergence in real coordinates at the point with index point_index after a call to FEPointEvaluation::evaluate() with EvaluationFlags::gradients set, or the divergence that has been stored there with a call to FEPointEvaluationBase::submit_divergence(). This functions only makes sense for a vector field with dim components.

Definition at line 756 of file fe_point_evaluation.h.

◆ submit_divergence()

void FEPointEvaluationBase< n_components_, dim, spacedim, double >::submit_divergence ( const double & value,
const unsigned int point_index )
inlineinherited

Write a contribution that is tested by the divergence to the field containing the values on points with the given point_index. Access to the same field as through get_divergence(). 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.

Note
This operation writes the data to the same field as submit_gradient(). As a consequence, only one of these functions can be used. In case several terms of this kind appear in a weak form, the contribution of a potential call to this function must be added into the diagonal of the rank-2 tensor contribution passed to submit_gradient(), in order not to overwrite information.

Definition at line 774 of file fe_point_evaluation.h.

◆ jacobian()

DerivativeForm< 1, dim, spacedim, double > FEPointEvaluationBase< n_components_, dim, spacedim, double >::jacobian ( const unsigned int point_index) const
inlineinherited

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

◆ inverse_jacobian()

DerivativeForm< 1, spacedim, dim, double > FEPointEvaluationBase< n_components_, dim, spacedim, double >::inverse_jacobian ( const unsigned int point_index) const
inlineinherited

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

◆ JxW()

double FEPointEvaluationBase< n_components_, dim, spacedim, double >::JxW ( const unsigned int point_index) const
inlineinherited

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

◆ real_point()

Point< spacedim, double > FEPointEvaluationBase< n_components_, dim, spacedim, double >::real_point ( const unsigned int point_index) const
inlineinherited

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

Deprecated
Use the function quadrature_point() instead.

Definition at line 808 of file fe_point_evaluation.h.

◆ quadrature_point()

Point< spacedim, double > FEPointEvaluationBase< n_components_, dim, spacedim, double >::quadrature_point ( const unsigned int point_index) const
inlineinherited

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

Definition at line 815 of file fe_point_evaluation.h.

◆ unit_point()

Point< dim, double > FEPointEvaluationBase< n_components_, dim, spacedim, double >::unit_point ( const unsigned int point_index) const
inlineinherited

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

◆ integrate_value()

FEPointEvaluationBase< n_components_, dim, spacedim, double >::scalar_value_type FEPointEvaluationBase< n_components_, dim, spacedim, double >::integrate_value ( ) const
inlineinherited

 * Take values collected at quadrature points via the submit_value() function, multiply by the Jacobian determinant and quadrature weights (JxW) and sum the values for all quadrature  * points on the cell. The result is a scalar, representing the integral  * of the function over the cell. 

Definition at line 832 of file fe_point_evaluation.h.

◆ quadrature_point_indices()

std_cxx20::ranges::iota_view< unsigned int, unsigned int > FEPointEvaluationBase< n_components_, dim, spacedim, double >::quadrature_point_indices ( ) const
inlineinherited

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

◆ n_active_entries_per_quadrature_batch()

unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::n_active_entries_per_quadrature_batch ( unsigned int q)
inherited

Returns how many lanes of a quadrature batch are active.

Definition at line 846 of file fe_point_evaluation.h.

◆ setup()

void FEPointEvaluationBase< n_components_, dim, spacedim, double >::setup ( const unsigned int first_selected_component)
protectedinherited

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

◆ do_reinit()

void FEPointEvaluationBase< n_components_, dim, spacedim, double >::do_reinit ( )
inlineprotectedinherited

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

Member Data Documentation

◆ dimension

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

Definition at line 1538 of file fe_point_evaluation.h.

◆ n_components

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

Definition at line 1539 of file fe_point_evaluation.h.

◆ n_lanes_user_interface

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

Definition at line 1771 of file fe_point_evaluation.h.

◆ n_lanes_internal

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

◆ stride

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

Definition at line 1775 of file fe_point_evaluation.h.

◆ n_q_batches

const unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::n_q_batches
protectedinherited

Number of quadrature batches of the current cell/face.

Definition at line 879 of file fe_point_evaluation.h.

◆ n_q_points

const unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::n_q_points
protectedinherited

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

Definition at line 884 of file fe_point_evaluation.h.

◆ n_q_points_scalar

const unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::n_q_points_scalar
protectedinherited

Number of quadrature points of the current cell/face.

Definition at line 889 of file fe_point_evaluation.h.

◆ mapping

SmartPointer<const Mapping<dim, spacedim> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::mapping
protectedinherited

Pointer to the Mapping object passed to the constructor.

Definition at line 894 of file fe_point_evaluation.h.

◆ fe

SmartPointer<const FiniteElement<dim, spacedim> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::fe
protectedinherited

Pointer to the FiniteElement object passed to the constructor.

Definition at line 899 of file fe_point_evaluation.h.

◆ poly

std::vector<Polynomials::Polynomial<double> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::poly
protectedinherited

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

◆ use_linear_path

bool FEPointEvaluationBase< n_components_, dim, spacedim, double >::use_linear_path
protectedinherited

Store whether the linear path should be used.

Definition at line 910 of file fe_point_evaluation.h.

◆ renumber

std::vector<unsigned int> FEPointEvaluationBase< n_components_, dim, spacedim, double >::renumber
protectedinherited

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

Definition at line 916 of file fe_point_evaluation.h.

◆ solution_renumbered

std::vector<scalar_value_type> FEPointEvaluationBase< n_components_, dim, spacedim, double >::solution_renumbered
protectedinherited

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

◆ solution_renumbered_vectorized

AlignedVector<vectorized_value_type> FEPointEvaluationBase< n_components_, dim, spacedim, double >::solution_renumbered_vectorized
protectedinherited

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

◆ scratch_data_scalar

AlignedVector<ScalarNumber> FEPointEvaluationBase< n_components_, dim, spacedim, double >::scratch_data_scalar
protectedinherited

Temporary array for the face path (scalar).

Definition at line 938 of file fe_point_evaluation.h.

◆ values

std::vector<value_type> FEPointEvaluationBase< n_components_, dim, spacedim, double >::values
protectedinherited

Temporary array to store the values at the points.

Definition at line 943 of file fe_point_evaluation.h.

◆ gradients

std::vector<gradient_type> FEPointEvaluationBase< n_components_, dim, spacedim, double >::gradients
protectedinherited

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

Definition at line 948 of file fe_point_evaluation.h.

◆ unit_point_ptr

const Point<dim, VectorizedArrayType>* FEPointEvaluationBase< n_components_, dim, spacedim, double >::unit_point_ptr
protectedinherited

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

Definition at line 954 of file fe_point_evaluation.h.

◆ unit_point_faces_ptr

const Point<dim - 1, VectorizedArrayType>* FEPointEvaluationBase< n_components_, dim, spacedim, double >::unit_point_faces_ptr
protectedinherited

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

Definition at line 960 of file fe_point_evaluation.h.

◆ real_point_ptr

const Point<spacedim, double>* FEPointEvaluationBase< n_components_, dim, spacedim, double >::real_point_ptr
protectedinherited

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

Definition at line 966 of file fe_point_evaluation.h.

◆ jacobian_ptr

const DerivativeForm<1, dim, spacedim, double>* FEPointEvaluationBase< n_components_, dim, spacedim, double >::jacobian_ptr
protectedinherited

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

Definition at line 972 of file fe_point_evaluation.h.

◆ inverse_jacobian_ptr

const DerivativeForm<1, spacedim, dim, double>* FEPointEvaluationBase< n_components_, dim, spacedim, double >::inverse_jacobian_ptr
protectedinherited

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

Definition at line 978 of file fe_point_evaluation.h.

◆ normal_ptr

const Tensor<1, spacedim, double>* FEPointEvaluationBase< n_components_, dim, spacedim, double >::normal_ptr
protectedinherited

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

Definition at line 984 of file fe_point_evaluation.h.

◆ JxW_ptr

const double* FEPointEvaluationBase< n_components_, dim, spacedim, double >::JxW_ptr
protectedinherited

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

◆ cell_type

internal::MatrixFreeFunctions::GeometryType FEPointEvaluationBase< n_components_, dim, spacedim, double >::cell_type
protectedinherited

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

Definition at line 996 of file fe_point_evaluation.h.

◆ dofs_per_component

unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::dofs_per_component
protectedinherited

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

Definition at line 1002 of file fe_point_evaluation.h.

◆ dofs_per_component_face

unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::dofs_per_component_face
protectedinherited

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

◆ shape_info

internal::MatrixFreeFunctions::ShapeInfo<ScalarNumber> FEPointEvaluationBase< n_components_, dim, spacedim, double >::shape_info
protectedinherited

Scalar ShapeInfo object needed for face path.

Definition at line 1014 of file fe_point_evaluation.h.

◆ component_in_base_element

unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::component_in_base_element
protectedinherited

The first selected component in the active base element.

Definition at line 1019 of file fe_point_evaluation.h.

◆ nonzero_shape_function_component

std::vector<std::array<bool, n_components> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::nonzero_shape_function_component
protectedinherited

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

Definition at line 1026 of file fe_point_evaluation.h.

◆ update_flags

const UpdateFlags FEPointEvaluationBase< n_components_, dim, spacedim, double >::update_flags
protectedinherited

The desired update flags for the evaluation.

Definition at line 1031 of file fe_point_evaluation.h.

◆ fe_values

std::shared_ptr<FEValues<dim, spacedim> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::fe_values
protectedinherited

The FEValues object underlying the slow evaluation path.

Definition at line 1036 of file fe_point_evaluation.h.

◆ mapping_info_on_the_fly

std::unique_ptr<NonMatching::MappingInfo<dim, spacedim, double> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::mapping_info_on_the_fly
protectedinherited

Pointer to mapping info on the fly computed during reinit.

Definition at line 1042 of file fe_point_evaluation.h.

◆ mapping_info

SmartPointer<NonMatching::MappingInfo<dim, spacedim, double> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::mapping_info
protectedinherited

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

Definition at line 1048 of file fe_point_evaluation.h.

◆ current_cell_index

unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::current_cell_index
protectedinherited

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

Definition at line 1053 of file fe_point_evaluation.h.

◆ current_face_number

unsigned int FEPointEvaluationBase< n_components_, dim, spacedim, double >::current_face_number
protectedinherited

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

Definition at line 1058 of file fe_point_evaluation.h.

◆ fast_path

bool FEPointEvaluationBase< n_components_, dim, spacedim, double >::fast_path
protectedinherited

Bool indicating if fast path is chosen.

Definition at line 1063 of file fe_point_evaluation.h.

◆ must_reinitialize_pointers

bool FEPointEvaluationBase< n_components_, dim, spacedim, double >::must_reinitialize_pointers
protectedinherited

Bool indicating if class needs to call reinit() inside evaluate()/integrate() functions, which is the case when the present class does not own the MappingInfo object but shares evaluation points with another object.

Definition at line 1071 of file fe_point_evaluation.h.

◆ shapes

AlignedVector<::ndarray<VectorizedArrayType, 2, dim> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::shapes
protectedinherited

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

Definition at line 1077 of file fe_point_evaluation.h.

◆ shapes_faces

AlignedVector<::ndarray<VectorizedArrayType, 2, dim - 1> > FEPointEvaluationBase< n_components_, dim, spacedim, double >::shapes_faces
protectedinherited

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

Definition at line 1083 of file fe_point_evaluation.h.

◆ is_interior

const bool FEPointEvaluationBase< n_components_, dim, spacedim, double >::is_interior
protectedinherited

Definition at line 1085 of file fe_point_evaluation.h.


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