Reference documentation for deal.II version 9.6.0
|
#include <deal.II/matrix_free/fe_point_evaluation.h>
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 | |
FEPointEvaluation (const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const UpdateFlags update_flags, const unsigned int first_selected_component=0) | |
FEPointEvaluation (NonMatching::MappingInfo< dim, spacedim, Number > &mapping_info, const FiniteElement< dim, spacedim > &fe, const unsigned int first_selected_component=0) | |
void | reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points) |
void | reinit () |
void | reinit (const unsigned int cell_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) |
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_type & | get_value (const unsigned int point_index) const |
void | submit_value (const value_type &value, const unsigned int point_index) |
const gradient_type & | get_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 int > | quadrature_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 () |
Private Member Functions | |
template<bool is_linear, std::size_t stride_view> | |
void | prepare_evaluate_fast (const StridedArrayView< const ScalarNumber, stride_view > &solution_values) |
template<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_unit_gradient_type &gradient) |
template<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_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_unit_gradient_type gradient, vectorized_value_type *solution_values_vectorized_linear) |
template<bool is_linear, std::size_t stride_view> | |
void | finish_integrate_fast (const StridedArrayView< ScalarNumber, stride_view > &solution_values, vectorized_value_type *solution_values_vectorized_linear, const bool sum_into_values) |
template<bool do_JxW, bool is_linear, std::size_t stride_view> | |
void | integrate_fast (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values) |
template<bool do_JxW, std::size_t stride_view> | |
void | integrate_slow (const StridedArrayView< ScalarNumber, stride_view > &solution_values, const EvaluationFlags::EvaluationFlags &integration_flags, const bool sum_into_values) |
template<bool do_JxW, 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) |
void | internal_reinit_single_cell_state_mapping_info () |
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 |
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 topic. 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 1124 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::number_type = Number |
Definition at line 1131 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::ScalarNumber |
Definition at line 1133 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::VectorizedArrayType |
Definition at line 1135 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::ETT |
Definition at line 1137 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::value_type = typename ETT::value_type |
Definition at line 1139 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::scalar_value_type = typename ETT::scalar_value_type |
Definition at line 1140 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::vectorized_value_type = typename ETT::vectorized_value_type |
Definition at line 1141 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::unit_gradient_type = typename ETT::unit_gradient_type |
Definition at line 1142 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::gradient_type = typename ETT::real_gradient_type |
Definition at line 1143 of file fe_point_evaluation.h.
using FEPointEvaluation< n_components_, dim, spacedim, Number >::interface_vectorized_unit_gradient_type |
Definition at line 1144 of file fe_point_evaluation.h.
FEPointEvaluation< n_components_, dim, spacedim, Number >::FEPointEvaluation | ( | const Mapping< dim, spacedim > & | mapping, |
const FiniteElement< dim, spacedim > & | fe, | ||
const UpdateFlags | update_flags, | ||
const unsigned int | first_selected_component = 0 ) |
Constructor.
mapping | The Mapping class describing the actual geometry of a cell passed to the evaluate() function. |
fe | The FiniteElement object that is used for the evaluation, which is typically the same on all cells to be evaluated. |
update_flags | Specify 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_component | For multi-component FiniteElement objects, this parameter allows to select a range of n_components components starting from this parameter. |
Definition at line 2354 of file fe_point_evaluation.h.
FEPointEvaluation< n_components_, dim, spacedim, Number >::FEPointEvaluation | ( | NonMatching::MappingInfo< dim, spacedim, Number > & | mapping_info, |
const FiniteElement< dim, spacedim > & | 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.
mapping_info | The 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. |
fe | The FiniteElement object that is used for the evaluation, which is typically the same on all cells to be evaluated. |
first_selected_component | For multi-component FiniteElement objects, this parameter allows to select a range of n_components components starting from this parameter. |
Definition at line 2341 of file fe_point_evaluation.h.
|
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.
[in] | cell | An iterator to the current cell |
[in] | unit_points | List 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 2396 of file fe_point_evaluation.h.
|
inline |
Reinitialize the evaluator to point to the correct precomputed mapping of the single cell in the MappingInfo object.
Definition at line 2386 of file fe_point_evaluation.h.
|
inline |
Reinitialize the evaluator to point to the correct precomputed mapping of the cell in the MappingInfo object.
Definition at line 2428 of file fe_point_evaluation.h.
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().
[in] | solution_values | This array is supposed to contain the unknown values on the element read out by FEEvaluation::read_dof_values(global_vector) . |
[in] | evaluation_flags | Flags specifying which quantities should be evaluated at the points. |
Definition at line 2466 of file fe_point_evaluation.h.
void FEPointEvaluation< 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().
[in] | solution_values | This array is supposed to contain the unknown values on the element as returned by cell->get_dof_values(global_vector,
solution_values) . |
[in] | evaluation_flags | Flags specifying which quantities should be evaluated at the points. |
Definition at line 2498 of file fe_point_evaluation.h.
void FEPointEvaluation< 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).
[out] | solution_values | This 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_flags | Flags specifying which quantities should be integrated at the points. |
[in] | sum_into_values | Flag 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 2512 of file fe_point_evaluation.h.
void FEPointEvaluation< 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).
[out] | solution_values | This 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_flags | Flags specifying which quantities should be integrated at the points. |
[in] | sum_into_values | Flag 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 2524 of file fe_point_evaluation.h.
void FEPointEvaluation< 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. 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.
[out] | solution_values | This 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_flags | Flags specifying which quantities should be integrated at the points. |
[in] | sum_into_values | Flag 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 2579 of file fe_point_evaluation.h.
void FEPointEvaluation< 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. 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.
[out] | solution_values | This 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_flags | Flags specifying which quantities should be integrated at the points. |
[in] | sum_into_values | Flag 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 2591 of file fe_point_evaluation.h.
|
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 3164 of file fe_point_evaluation.h.
|
inline |
Return the normal derivative in real coordinates at the point with index point_index
after a call to FEPointEvaluation::evaluate() with EvaluationFlags::gradients set.
Definition at line 3182 of file fe_point_evaluation.h.
|
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 3193 of file fe_point_evaluation.h.
|
inlineprivate |
Resizes necessary data fields, reads in and renumbers solution values. Interpolates onto face if face path is selected.
Definition at line 2607 of file fe_point_evaluation.h.
|
inlineprivate |
Evaluates the actual interpolation on the cell or face for a quadrature batch.
Definition at line 2641 of file fe_point_evaluation.h.
|
inlineprivate |
Fast path of the evaluate function.
Definition at line 2714 of file fe_point_evaluation.h.
|
inlineprivate |
Slow path of the evaluate function using FEValues.
Definition at line 2775 of file fe_point_evaluation.h.
|
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 2862 of file fe_point_evaluation.h.
|
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 2903 of file fe_point_evaluation.h.
|
inlineprivate |
Fast path of the integrate function.
Definition at line 2955 of file fe_point_evaluation.h.
|
inlineprivate |
Slow path of the integrate function using FEValues.
Definition at line 3037 of file fe_point_evaluation.h.
|
private |
Implementation of the integrate/test_and_sum function.
Definition at line 3118 of file fe_point_evaluation.h.
|
inlineprivate |
Internal function to initialize the pointers of this class when an external MappingInfo has already queried the mapping for the relevant information.
Definition at line 2370 of file fe_point_evaluation.h.
|
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.
|
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.
|
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.
|
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.
|
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.
|
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.
Definition at line 774 of file fe_point_evaluation.h.
|
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.
|
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.
|
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.
|
inlineinherited |
Return the position in real coordinates of the given point index among the points passed to reinit().
Definition at line 808 of file fe_point_evaluation.h.
|
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.
|
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.
|
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.
|
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.
|
inherited |
Returns how many lanes of a quadrature batch are active.
Definition at line 846 of file fe_point_evaluation.h.
|
protectedinherited |
Common setup function for both constructors. Does the setup for both fast and slow path.
first_selected_component | For 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.
|
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.
|
staticconstexpr |
Definition at line 1128 of file fe_point_evaluation.h.
|
staticconstexpr |
Definition at line 1129 of file fe_point_evaluation.h.
|
staticconstexprprivate |
Definition at line 1397 of file fe_point_evaluation.h.
|
staticconstexprprivate |
Definition at line 1399 of file fe_point_evaluation.h.
|
staticconstexprprivate |
Definition at line 1401 of file fe_point_evaluation.h.
|
protectedinherited |
Number of quadrature batches of the current cell/face.
Definition at line 879 of file fe_point_evaluation.h.
|
protectedinherited |
Number of quadrature points/batches of the current cell/face.
Definition at line 884 of file fe_point_evaluation.h.
|
protectedinherited |
Number of quadrature points of the current cell/face.
Definition at line 889 of file fe_point_evaluation.h.
|
protectedinherited |
Pointer to the Mapping object passed to the constructor.
Definition at line 894 of file fe_point_evaluation.h.
|
protectedinherited |
Pointer to the FiniteElement object passed to the constructor.
Definition at line 899 of file fe_point_evaluation.h.
|
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.
|
protectedinherited |
Store whether the linear path should be used.
Definition at line 910 of file fe_point_evaluation.h.
|
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.
|
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.
|
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.
|
protectedinherited |
Temporary array for the face path (scalar).
Definition at line 938 of file fe_point_evaluation.h.
|
protectedinherited |
Temporary array to store the values at the points.
Definition at line 943 of file fe_point_evaluation.h.
|
protectedinherited |
Temporary array to store the gradients in real coordinates at the points.
Definition at line 948 of file fe_point_evaluation.h.
|
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.
|
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.
|
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.
|
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.
|
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.
|
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.
|
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.
|
protectedinherited |
Cell type describing the geometry of the cell and compression of jacobians.
Definition at line 996 of file fe_point_evaluation.h.
|
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.
|
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.
|
protectedinherited |
Scalar ShapeInfo object needed for face path.
Definition at line 1014 of file fe_point_evaluation.h.
|
protectedinherited |
The first selected component in the active base element.
Definition at line 1019 of file fe_point_evaluation.h.
|
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.
|
protectedinherited |
The desired update flags for the evaluation.
Definition at line 1031 of file fe_point_evaluation.h.
|
protectedinherited |
The FEValues object underlying the slow evaluation path.
Definition at line 1036 of file fe_point_evaluation.h.
|
protectedinherited |
Pointer to mapping info on the fly computed during reinit.
Definition at line 1042 of file fe_point_evaluation.h.
|
protectedinherited |
Pointer to currently used mapping info (either on the fly or external precomputed).
Definition at line 1048 of file fe_point_evaluation.h.
|
protectedinherited |
The current cell index to access mapping data from mapping info.
Definition at line 1053 of file fe_point_evaluation.h.
|
protectedinherited |
The current face number to access mapping data from mapping info.
Definition at line 1058 of file fe_point_evaluation.h.
|
protectedinherited |
Bool indicating if fast path is chosen.
Definition at line 1063 of file fe_point_evaluation.h.
|
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.
|
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.
|
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.
|
protectedinherited |
Definition at line 1085 of file fe_point_evaluation.h.