deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/matrix_free/fe_evaluation.h>
Public Types | |
using | BaseClass = FEEvaluationBase< dim, n_components_, Number, true, VectorizedArrayType > |
using | number_type = Number |
using | value_type = typename BaseClass::value_type |
using | gradient_type = typename BaseClass::gradient_type |
using | hessian_type = std::conditional_t< n_components_==1, Tensor< 2, dim, VectorizedArrayType >, std::conditional_t< n_components_==dim, Tensor< 3, dim, VectorizedArrayType >, Tensor< 1, n_components_, Tensor< 2, dim, VectorizedArrayType > > > > |
using | NumberType = VectorizedArrayType |
using | ScalarNumber = typename internal::VectorizedArrayTrait< VectorizedArrayType >::value_type |
using | shape_info_number_type = ScalarNumber |
Public Member Functions | |
FEFaceEvaluation (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const bool is_interior_face=true, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0, const unsigned int active_fe_index=numbers::invalid_unsigned_int, const unsigned int active_quad_index=numbers::invalid_unsigned_int, const unsigned int face_type=numbers::invalid_unsigned_int) | |
FEFaceEvaluation (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const std::pair< unsigned int, unsigned int > &range, const bool is_interior_face=true, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0) | |
void | reinit (const unsigned int face_batch_number) |
void | reinit (const unsigned int cell_batch_number, const unsigned int face_number) |
void | evaluate (const EvaluationFlags::EvaluationFlags evaluation_flag) |
void | evaluate (const VectorizedArrayType *values_array, const EvaluationFlags::EvaluationFlags evaluation_flag) |
void | project_to_face (const EvaluationFlags::EvaluationFlags evaluation_flag) |
void | project_to_face (const VectorizedArrayType *values_array, const EvaluationFlags::EvaluationFlags evaluation_flag) |
void | evaluate_in_face (const EvaluationFlags::EvaluationFlags evaluation_flag) |
template<typename VectorType > | |
void | gather_evaluate (const VectorType &input_vector, const EvaluationFlags::EvaluationFlags evaluation_flag) |
void | integrate (const EvaluationFlags::EvaluationFlags integration_flag, const bool sum_into_values=false) |
void | integrate (const EvaluationFlags::EvaluationFlags integration_flag, VectorizedArrayType *values_array, const bool sum_into_values=false) |
void | integrate_in_face (const EvaluationFlags::EvaluationFlags integration_flag) |
void | collect_from_face (const EvaluationFlags::EvaluationFlags integration_flag, const bool sum_into_values=false) |
void | collect_from_face (const EvaluationFlags::EvaluationFlags integration_flag, VectorizedArrayType *values_array, const bool sum_into_values=false) |
template<typename VectorType > | |
void | integrate_scatter (const EvaluationFlags::EvaluationFlags integration_flag, VectorType &output_vector) |
template<typename VectorType > | |
void | integrate_scatter (const bool integrate_values, const bool integrate_gradients, VectorType &output_vector) |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | dof_indices () const |
bool | at_boundary () const |
types::boundary_id | boundary_id () const |
unsigned int | get_dofs_per_component_projected_to_face () |
unsigned int | get_dofs_projected_to_face () |
const MatrixFree< dim, Number, VectorizedArrayType > & | get_matrix_free () const |
void | set_data_pointers (AlignedVector< VectorizedArrayType > *scratch_data, const unsigned int n_components) |
void | reinit_face (const internal::MatrixFreeFunctions::FaceToCellTopology< n_lanes > &face) |
Reading from and writing to vectors | |
template<typename VectorType > | |
void | read_dof_values (const VectorType &src, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip()) |
template<typename VectorType > | |
void | read_dof_values_plain (const VectorType &src, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip()) |
template<typename VectorType > | |
void | distribute_local_to_global (VectorType &dst, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip()) const |
template<typename VectorType > | |
void | set_dof_values (VectorType &dst, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip()) const |
template<typename VectorType > | |
void | set_dof_values_plain (VectorType &dst, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip()) const |
Access to data at quadrature points or the data at cell DoFs | |
value_type | get_dof_value (const unsigned int dof) const |
void | submit_dof_value (const value_type val_in, const unsigned int dof) |
value_type | get_value (const unsigned int q_point) const |
void | submit_value (const value_type val_in, const unsigned int q_point) |
template<int n_components_local = n_components, typename = std::enable_if_t<n_components == n_components_local>> | |
void | submit_value (const Tensor< 1, 1, VectorizedArrayType > val_in, const unsigned int q_point) |
gradient_type | get_gradient (const unsigned int q_point) const |
value_type | get_normal_derivative (const unsigned int q_point) const |
void | submit_gradient (const gradient_type grad_in, const unsigned int q_point) |
template<int dim_ = dim, typename = std::enable_if_t<dim_ == 1 && n_components == dim_>> | |
void | submit_gradient (const Tensor< 2, 1, VectorizedArrayType > val_in, const unsigned int q_point) |
void | submit_normal_derivative (const value_type grad_in, const unsigned int q_point) |
hessian_type | get_hessian (const unsigned int q_point) const |
gradient_type | get_hessian_diagonal (const unsigned int q_point) const |
value_type | get_laplacian (const unsigned int q_point) const |
value_type | get_normal_hessian (const unsigned int q_point) const |
void | submit_hessian (const hessian_type hessian_in, const unsigned int q_point) |
void | submit_normal_hessian (const value_type normal_hessian_in, const unsigned int q_point) |
template<int dim_ = dim, typename = std::enable_if_t<n_components_ == dim_>> | |
VectorizedArrayType | get_divergence (const unsigned int q_point) const |
template<int dim_ = dim, typename = std::enable_if_t<n_components_ == dim_>> | |
void | submit_divergence (const VectorizedArrayType div_in, const unsigned int q_point) |
template<int dim_ = dim, typename = std::enable_if_t<n_components_ == dim_>> | |
SymmetricTensor< 2, dim, VectorizedArrayType > | get_symmetric_gradient (const unsigned int q_point) const |
template<int dim_ = dim, typename = std::enable_if_t<n_components_ == dim_>> | |
void | submit_symmetric_gradient (const SymmetricTensor< 2, dim, VectorizedArrayType > grad_in, const unsigned int q_point) |
template<int dim_ = dim, typename = std::enable_if_t<n_components_ == dim_ && dim_ != 1>> | |
Tensor< 1,(dim==2 ? 1 :dim), VectorizedArrayType > | get_curl (const unsigned int q_point) const |
template<int dim_ = dim, typename = std::enable_if_t<n_components_ == dim_ && dim != 1>> | |
void | submit_curl (const Tensor< 1, dim==2 ? 1 :dim, VectorizedArrayType > curl_in, const unsigned int q_point) |
value_type | integrate_value () const |
Access to geometry data at quadrature points | |
ArrayView< VectorizedArrayType > | get_scratch_data () const |
VectorizedArrayType | JxW (const unsigned int q_point) const |
Point< dim, VectorizedArrayType > | quadrature_point (const unsigned int q) const |
Tensor< 2, dim, VectorizedArrayType > | inverse_jacobian (const unsigned int q_point) const |
Tensor< 1, dim, VectorizedArrayType > | normal_vector (const unsigned int q_point) const |
Tensor< 1, dim, VectorizedArrayType > | get_normal_vector (const unsigned int q_point) const |
Access to internal data arrays | |
const VectorizedArrayType * | begin_dof_values () const |
VectorizedArrayType * | begin_dof_values () |
const VectorizedArrayType * | begin_values () const |
VectorizedArrayType * | begin_values () |
const VectorizedArrayType * | begin_gradients () const |
VectorizedArrayType * | begin_gradients () |
const VectorizedArrayType * | begin_hessians () const |
VectorizedArrayType * | begin_hessians () |
Information about the current cell this class operates on | |
unsigned int | get_mapping_data_index_offset () const |
internal::MatrixFreeFunctions::GeometryType | get_cell_type () const |
const ShapeInfoType & | get_shape_info () const |
const internal::MatrixFreeFunctions::DoFInfo & | get_dof_info () const |
const std::vector< unsigned int > & | get_internal_dof_numbering () const |
unsigned int | get_quadrature_index () const |
unsigned int | get_current_cell_index () const |
unsigned int | get_active_fe_index () const |
unsigned int | get_active_quadrature_index () const |
unsigned int | get_first_selected_component () const |
std::uint8_t | get_face_no (const unsigned int v=0) const |
unsigned int | get_subface_index () const |
std::uint8_t | get_face_orientation (const unsigned int v=0) const |
internal::MatrixFreeFunctions::DoFInfo::DoFAccessIndex | get_dof_access_index () const |
bool | is_interior_face () const |
const std::array< unsigned int, n_lanes > & | get_cell_ids () const |
const std::array< unsigned int, n_lanes > & | get_face_ids () const |
unsigned int | get_cell_or_face_batch_id () const |
const std::array< unsigned int, n_lanes > & | get_cell_or_face_ids () const |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | quadrature_point_indices () const |
Functions to access cell- and face-data vectors. | |
T | read_cell_data (const AlignedVector< T > &array) const |
void | set_cell_data (AlignedVector< T > &array, const T &value) const |
T | read_face_data (const AlignedVector< T > &array) const |
void | set_face_data (AlignedVector< T > &array, const T &value) const |
Static Public Member Functions | |
static bool | fast_evaluation_supported (const unsigned int given_degree, const unsigned int given_n_q_points_1d) |
Public Attributes | |
const unsigned int | dofs_per_component |
const unsigned int | dofs_per_cell |
const unsigned int | n_q_points |
Static Public Attributes | |
static constexpr unsigned int | dimension = dim |
static constexpr unsigned int | n_components = n_components_ |
static constexpr unsigned int | n_lanes = VectorizedArrayType::size() |
static constexpr unsigned int | static_n_q_points |
static constexpr unsigned int | static_n_q_points_cell |
static constexpr unsigned int | static_dofs_per_component |
static constexpr unsigned int | tensor_dofs_per_cell |
static constexpr unsigned int | static_dofs_per_cell |
Protected Member Functions | |
template<typename VectorType , typename VectorOperation > | |
void | read_write_operation (const VectorOperation &operation, const std::array< VectorType *, n_components_ > &vectors, const std::array< const std::vector< ArrayView< const typename VectorType::value_type > > *, n_components_ > &vectors_sm, const std::bitset< n_lanes > &mask, const bool apply_constraints=true) const |
template<typename VectorType , typename VectorOperation > | |
void | read_write_operation_contiguous (const VectorOperation &operation, const std::array< VectorType *, n_components_ > &vectors, const std::array< const std::vector< ArrayView< const typename VectorType::value_type > > *, n_components_ > &vectors_sm, const std::bitset< n_lanes > &mask) const |
template<typename VectorType , typename VectorOperation > | |
void | read_write_operation_global (const VectorOperation &operation, const std::array< VectorType *, n_components_ > &vectors) const |
void | apply_hanging_node_constraints (const bool transpose) const |
Private Types | |
using | ShapeInfoType = internal::MatrixFreeFunctions::ShapeInfo< typename internal::VectorizedArrayTrait< VectorizedArrayType >::value_type > |
using | MappingInfoStorageType = internal::MatrixFreeFunctions::MappingInfoStorage<(is_face ? dim - 1 :dim), dim, VectorizedArrayType > |
using | DoFInfo = internal::MatrixFreeFunctions::DoFInfo |
The class that provides all functions necessary to evaluate functions at quadrature points and face integrations. The design of the class is similar to FEEvaluation and most of the interfaces are shared with that class, in particular most access functions that come from the common base classes FEEvaluationData and FEEvaluationBase. Furthermore, the relation of this class to FEEvaluation is similar to the relation between FEValues and FEFaceValues.
dim | Dimension in which this class is to be used |
fe_degree | Degree of the tensor product finite element with fe_degree+1 degrees of freedom per coordinate direction. If set to -1, the degree of the underlying element will be used, which acts as a run time constant rather than a compile time constant that slows down the execution. |
n_q_points_1d | Number of points in the quadrature formula in 1d, usually chosen as fe_degree+1 |
n_components | Number of vector components when solving a system of PDEs. If the same operation is applied to several components of a PDE (e.g. a vector Laplace equation), they can be applied simultaneously with one call (and often more efficiently) |
Number | Number format, usually double or float |
VectorizedArrayType | Type of array to be woked on in a vectorized fashion, defaults to VectorizedArray<Number> |
Definition at line 1844 of file fe_evaluation.h.
using FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::BaseClass = FEEvaluationBase<dim, n_components_, Number, true, VectorizedArrayType> |
An alias to the base class.
Definition at line 1858 of file fe_evaluation.h.
using FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::number_type = Number |
A underlying number type specified as template argument.
Definition at line 1864 of file fe_evaluation.h.
using FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::value_type = typename BaseClass::value_type |
The type of function values, e.g. VectorizedArrayType
for n_components=1
or Tensor<1,dim,VectorizedArrayType >
for n_components=dim
.
Definition at line 1871 of file fe_evaluation.h.
using FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::gradient_type = typename BaseClass::gradient_type |
The type of gradients, e.g. Tensor<1,dim,VectorizedArrayType>
for n_components=1
or Tensor<2,dim,VectorizedArrayType >
for n_components=dim
.
Definition at line 1878 of file fe_evaluation.h.
|
inherited |
Definition at line 109 of file fe_evaluation.h.
|
privateinherited |
Definition at line 115 of file fe_evaluation_data.h.
|
privateinherited |
Definition at line 117 of file fe_evaluation_data.h.
|
privateinherited |
Definition at line 119 of file fe_evaluation_data.h.
|
inherited |
Definition at line 124 of file fe_evaluation_data.h.
|
inherited |
Definition at line 125 of file fe_evaluation_data.h.
|
inherited |
Definition at line 127 of file fe_evaluation_data.h.
FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::FEFaceEvaluation | ( | const MatrixFree< dim, Number, VectorizedArrayType > & | matrix_free, |
const bool | is_interior_face = true , |
||
const unsigned int | dof_no = 0 , |
||
const unsigned int | quad_no = 0 , |
||
const unsigned int | first_selected_component = 0 , |
||
const unsigned int | active_fe_index = numbers::invalid_unsigned_int , |
||
const unsigned int | active_quad_index = numbers::invalid_unsigned_int , |
||
const unsigned int | face_type = numbers::invalid_unsigned_int |
||
) |
Constructor. Takes all data stored in MatrixFree. If applied to problems with more than one finite element or more than one quadrature formula selected during construction of matrix_free
, the appropriate component can be selected by the optional arguments.
matrix_free | Data object that contains all data |
is_interior_face | This selects which of the two cells of an internal face the current evaluator will be based upon. The interior face is the main face along which the normal vectors are oriented. The exterior face coming from the other side provides the same normal vector as the interior side, so if the outer normal vector to that side is desired, it must be multiplied by -1. |
dof_no | If matrix_free was set up with multiple DoFHandler objects, this parameter selects to which DoFHandler/AffineConstraints pair the given evaluator should be attached to. |
quad_no | If matrix_free was set up with multiple Quadrature objects, this parameter selects the appropriate number of the quadrature formula. |
first_selected_component | If the dof_handler selected by dof_no uses an FESystem consisting of more than one base element, this parameter selects the number of the base element in FESystem. Note that this does not directly relate to the component of the respective element due to the possibility for a multiplicity in the element. |
active_fe_index | If matrix_free was set up with DoFHandler objects with hp::FECollections, this parameter selects to which DoFHandler/AffineConstraints pair the given evaluator should be attached to. |
face_type | In the case of a face, indicate its reference-cell type (0 for line or quadrilateral 1 for triangle). |
active_quad_index | If matrix_free was set up with hp::Collection objects, this parameter selects the appropriate number of the quadrature formula. |
FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::FEFaceEvaluation | ( | const MatrixFree< dim, Number, VectorizedArrayType > & | matrix_free, |
const std::pair< unsigned int, unsigned int > & | range, | ||
const bool | is_interior_face = true , |
||
const unsigned int | dof_no = 0 , |
||
const unsigned int | quad_no = 0 , |
||
const unsigned int | first_selected_component = 0 |
||
) |
Constructor. Takes all data stored in MatrixFree for a given face range, which allows to automatically identify the active_fe_index and active_quad_index in case of a p-adaptive strategy.
The rest of the arguments are the same as in the constructor above.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::reinit | ( | const unsigned int | face_batch_number | ) |
Initializes the operation pointer to the current face. This method is the default choice for face integration as the data stored in MappingInfo is stored according to this numbering. Unlike the reinit functions taking a cell iterator as argument below and the FEValues::reinit() methods, where the information related to a particular cell is generated in the reinit call, this function is very cheap since all data is pre-computed in matrix_free
, and only a few indices and pointers have to be set appropriately.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::reinit | ( | const unsigned int | cell_batch_number, |
const unsigned int | face_number | ||
) |
As opposed to the reinit() method from the base class, this reinit() method initializes for a given number of cells and a face number. This method is less efficient than the other reinit() method taking a numbering of the faces because it needs to copy the data associated with the faces to the cells in this call.
|
static |
Check if fast evaluation/integration is supported.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::evaluate | ( | const EvaluationFlags::EvaluationFlags | evaluation_flag | ) |
Evaluates the function values, the gradients, and the Laplacians of the FE function given at the DoF values stored in the internal data field dof_values
(that is usually filled by the read_dof_values() method) at the quadrature points on the unit cell. The function arguments specify which parts shall actually be computed. Needs to be called before the functions get_value(), get_gradient() or get_normal_derivative() give useful information (unless these values have been set manually by accessing the internal data pointers).
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::evaluate | ( | const VectorizedArrayType * | values_array, |
const EvaluationFlags::EvaluationFlags | evaluation_flag | ||
) |
Evaluates the function values, the gradients, and the Laplacians of the FE function given at the DoF values in the input array values_array
at the quadrature points on the unit cell. If multiple components are involved in the current FEEvaluation object, the sorting in values_array is such that all degrees of freedom for the first component come first, then all degrees of freedom for the second, and so on. The function arguments specify which parts shall actually be computed. Needs to be called before the functions get_value(), get_gradient(), or get_normal_derivative() give useful information (unless these values have been set manually).
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::project_to_face | ( | const EvaluationFlags::EvaluationFlags | evaluation_flag | ) |
Projects the values, the gradients, and the Hessians into the face DoFs of the current face using the internally stored cell DoFs.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::project_to_face | ( | const VectorizedArrayType * | values_array, |
const EvaluationFlags::EvaluationFlags | evaluation_flag | ||
) |
Projects the values, the gradients, and the Hessians into the face DoFs of the current face using the cell DoFs provided via values_array
.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::evaluate_in_face | ( | const EvaluationFlags::EvaluationFlags | evaluation_flag | ) |
Evaluates the values, the gradients, and the Hessians in-face, interpolating into the face quadrature points.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::gather_evaluate | ( | const VectorType & | input_vector, |
const EvaluationFlags::EvaluationFlags | evaluation_flag | ||
) |
Reads from the input vector and evaluates the function values, the gradients, and the Laplacians of the FE function at the quadrature points on the unit cell. The function arguments specify which parts shall actually be computed. Needs to be called before the functions get_value(), get_gradient(), or get_normal_derivative() give useful information.
This call is equivalent to calling read_dof_values() followed by evaluate(), but might internally use some additional optimizations.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::integrate | ( | const EvaluationFlags::EvaluationFlags | integration_flag, |
const bool | sum_into_values = false |
||
) |
This function takes the values and/or gradients that are stored on quadrature points, tests them by all the basis functions/gradients on the cell and performs the cell integration. The function argument integration_flag
is used to control which of the submitted contributions are used during the integration. The result is written into the internal data field dof_values
(that is usually written into the result vector by the distribute_local_to_global() or set_dof_values() methods).
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::integrate | ( | const EvaluationFlags::EvaluationFlags | integration_flag, |
VectorizedArrayType * | values_array, | ||
const bool | sum_into_values = false |
||
) |
This function takes the values and/or gradients that are stored on quadrature points, tests them by all the basis functions/gradients on the cell and performs the cell integration. The function argument integration_flag
is used to control which of the submitted contributions are used during the integration. As opposed to the other integrate() method, this call stores the result of the testing in the given array values_array
.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::integrate_in_face | ( | const EvaluationFlags::EvaluationFlags | integration_flag | ) |
This function tests the values, gradients and Hessians submitted on the face quadrature points by multiplying with the in-face basis function values, gradients and Hessians and accumulating to the respective face DoFs.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::collect_from_face | ( | const EvaluationFlags::EvaluationFlags | integration_flag, |
const bool | sum_into_values = false |
||
) |
Collects the contributions from the face DoFs of values, normal gradients and normal Hessians to the internal cell DoFs.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::collect_from_face | ( | const EvaluationFlags::EvaluationFlags | integration_flag, |
VectorizedArrayType * | values_array, | ||
const bool | sum_into_values = false |
||
) |
Collects the contributions from the face DoFs of values, normal gradients and normal Hessians to the cell DoFs specified via values_array
.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::integrate_scatter | ( | const EvaluationFlags::EvaluationFlags | integration_flag, |
VectorType & | output_vector | ||
) |
This function takes the values and/or gradients that are stored on quadrature points, tests them by all the basis functions/gradients on the cell and performs the cell integration. The function argument integration_flag
is used to control which of the submitted contributions are used during the integration.
This call is equivalent to calling integrate() followed by distribute_local_to_global(), but might internally use some additional optimizations.
void FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::integrate_scatter | ( | const bool | integrate_values, |
const bool | integrate_gradients, | ||
VectorType & | output_vector | ||
) |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::dof_indices | ( | ) | const |
Return an object that can be thought of as an array containing all indices from zero to dofs_per_cell
. This allows to write code using range-based for loops.
bool FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::at_boundary | ( | ) | const |
Return whether the face associated to this FEFaceEvaluation object is at the boundary.
types::boundary_id FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::boundary_id | ( | ) | const |
Return the boundary indicator of the face associated to this FEFaceEvaluation object.
If the return value is the special value numbers::internal_face_boundary_id, then the face is in the interior of the domain.
unsigned int FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::get_dofs_per_component_projected_to_face | ( | ) |
Get the number of degrees of freedom of a single component which are projected onto a face.
unsigned int FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::get_dofs_projected_to_face | ( | ) |
Get the number of degrees of freedom accumulated over all components which are projected onto a face.
|
inherited |
For the vector src
, read out the values on the degrees of freedom of the current cell, and store them internally. Similar functionality as the function DoFAccessor::get_interpolated_dof_values when no constraints are present, but it also includes constraints from hanging nodes, so one can see it as a similar function to AffineConstraints::read_dof_values as well. Note that if vectorization is enabled, the DoF values for several cells are set.
If some constraints on the vector are inhomogeneous, use the function read_dof_values_plain instead and provide the vector with useful data also in constrained positions by calling AffineConstraints::distribute. When accessing vector entries during the solution of linear systems, the temporary solution should always have homogeneous constraints and this method is the correct one.
If the given vector template class is a block vector (determined through the template function 'IsBlockVector<VectorType>::value', which checks for vectors derived from BlockVectorBase) or an std::vector<VectorType> or std::vector<VectorType *>, this function reads n_components
blocks from the block vector starting at the index first_index
. For non-block vectors, first_index
is ignored.
|
inherited |
For the vector src
, read out the values on the degrees of freedom of the current cell, and store them internally. Similar functionality as the function DoFAccessor::get_interpolated_dof_values. As opposed to the read_dof_values function, this function reads out the plain entries from vectors, without taking stored constraints into account. This way of access is appropriate when the constraints have been distributed on the vector by a call to AffineConstraints::distribute previously. This function is also necessary when inhomogeneous constraints are to be used, as MatrixFree can only handle homogeneous constraints. Note that if vectorization is enabled, the DoF values for several cells are set.
If the given vector template class is a block vector (determined through the template function 'IsBlockVector<VectorType>::value', which checks for vectors derived from BlockVectorBase) or an std::vector<VectorType> or std::vector<VectorType *>, this function reads n_components
blocks from the block vector starting at the index first_index
. For non-block vectors, first_index
is ignored.
|
inherited |
Takes the values stored internally on dof values of the current cell and sums them into the vector dst
. The function also applies constraints during the write operation. The functionality is hence similar to the function AffineConstraints::distribute_local_to_global. If vectorization is enabled, the DoF values for several cells are used.
If the given vector template class is a block vector (determined through the template function 'IsBlockVector<VectorType>::value', which checks for vectors derived from BlockVectorBase) or an std::vector<VectorType> or std::vector<VectorType *>, this function writes to n_components
blocks of the block vector starting at the index first_index
. For non-block vectors, first_index
is ignored.
The mask
can be used to suppress the write access for some of the cells contained in the current cell vectorization batch, e.g. in case of local time stepping, where some cells are excluded from a call. A value of true
in the bitset means that the respective lane index will be processed, whereas a value of false
skips this index. The default setting is a bitset that contains all ones, which will write the accumulated integrals to all cells in the batch.
|
inherited |
Takes the values stored internally on dof values of the current cell and writes them into the vector dst
. The function skips the degrees of freedom which are constrained. As opposed to the distribute_local_to_global method, the old values at the position given by the current cell are overwritten. Thus, if a degree of freedom is associated to more than one cell (as usual in continuous finite elements), the values will be overwritten and only the value written last is retained. Please note that in a parallel context this function might also touch degrees of freedom owned by other MPI processes, so that a subsequent update or accumulation of ghost values as done by MatrixFree::loop() might invalidate the degrees of freedom set by this function.
If the given vector template class is a block vector (determined through the template function 'IsBlockVector<VectorType>::value', which checks for vectors derived from BlockVectorBase) or an std::vector<VectorType> or std::vector<VectorType *>, this function writes to n_components
blocks of the block vector starting at the index first_index
. For non-block vectors, first_index
is ignored.
The mask
can be used to suppress the write access for some of the cells contained in the current cell vectorization batch, e.g. in case of local time stepping, where some cells are excluded from a call. A value of true
in the bitset means that the respective lane index will be processed, whereas a value of false
skips this index. The default setting is a bitset that contains all ones, which will write the accumulated integrals to all cells in the batch.
|
inherited |
Same as set_dof_values(), but without resolving constraints.
|
inherited |
Return the value stored in the array of coefficients for the elemental finite element solution expansion for the local degree of freedom with index dof
. If the object is vector-valued, a vector-valued return argument is given. Thus, the argument dof
can at most run until dofs_per_component
rather than dofs_per_cell
since the different components of a vector-valued FE are returned together. Note that when vectorization is enabled, values from several cells are grouped together in the inner VectorizedArray argument of value_type
. If submit_dof_value
was called last, the value corresponds to the data set there for the respective index. If integrate
was called last, it instead corresponds to the value of the integrated function with the test function of the given index that gets written into the global vector by distribute_local_to_global().
|
inherited |
Write a value to the field containing coefficients associated with the cell-local degrees of freedom with index dof
. This function writes to the same field as is accessed through get_dof_value
. Therefore, the original data that is stored in this location, e.g. from reading a global vector via read_dof_values(), is overwritten as soon as a value is submitted for a particular DoF index.
|
inherited |
Return the value of a finite element function interpolated to the quadrature point with index q_point
after a call to evaluate()
with EvaluationFlags::values
set, or the value that has been stored there with a call to FEEvaluationBase::submit_value(). If the object is vector-valued, a vector-valued return argument is given. In case vectorization is enabled, values from several cells are grouped together as a VectorizedArray for each component in value_type
.
|
inherited |
Write a contribution that gets multiplied by the value of the test function to the field containing the values at quadrature points with index q_point
. As part of this function, the data gets multiplied by the quadrature weight and possible Jacobian determinants. When this function has been called for all quadrature points and a subsequent call to the function integrate(EvaluationFlags::values)
has been made, the result is an array of entries, each representing the result of the integral of product of the test function multiplied by the data passed to this function.
|
inherited |
For scalar elements, the value_type and gradient_type can be unintentionally mixed up because FEEvaluationBase does not distinguish between scalar accessors and vector-valued accessors and the respective types, but solely in terms of the number of components and dimension. Thus, enable the use of submit_value() also for tensors with a single component.
|
inherited |
Return the gradient of the finite element function evaluated at quadrature point with index q_point
after a call to evaluate()
with EvaluationFlags::gradients
set, collecting all components in a vector-valued problem as the outer tensor index and all partial derivatives as the inner (second) tensor index. For scalar problems, gradient_type
is overloaded as a Tensor<1, dim>. For further information, see the get_value() function.
|
inherited |
Return the derivative of a finite element function interpolated to the quadrature point with index q_point
after a call to evaluate(EvaluationFlags::gradients)
the direction normal to the face: \(\boldsymbol \nabla u(\mathbf x_q) \cdot \mathbf n(\mathbf
x_q)\).
This call is equivalent to calling get_gradient() * normal_vector() but will use a more efficient internal representation of data.
|
inherited |
Write a contribution that gets multiplied by the gradient of the test function to the field containing the gradients at quadrature points with index q_point
. When this function has queued information for all quadrature points and followed by a call to the function integrate(EvaluationFlags::gradients)
, the result is an array of entries, each representing the result of the integral of product of the test function gradient multiplied by the values passed to this function.
|
inherited |
In 1D, the value_type and gradient_type can be unintentionally mixed up because FEEvaluationBase does not distinguish between scalar accessors and vector-valued accessors and the respective types, but solely in terms of the number of components and dimension. Thus, enable the use of submit_gradient() also for rank-2 tensors with a single component.
|
inherited |
Write a contribution that gets multiplied by the gradient of the test function times the normal vector to the field containing the gradients at quadrature points with index q_point
, see submit_gradient() for further information.
|
inherited |
Return the Hessian of the finite element function interpolated to quadrature point with index q_point
after a call to evaluate(EvaluationFlags::hessians)
. If only the diagonal or even the trace of the Hessian, the Laplacian, is needed, use the other functions below.
|
inherited |
Return the diagonal of the Hessian of the finite element function interpolated to the quadrature point with index q_point
after a call to evaluate(EvaluationFlags::hessians)
.
|
inherited |
Return the Laplacian (i.e., the trace of the Hessian) of the finite element function interpolated to the quadrature point with index q_point
after a call to evaluate(EvaluationFlags::hessians)
. Compared to the case when computing the full Hessian, some operations can be saved when only the Laplacian is requested.
|
inherited |
Return the second derivative along the normal direction \(\partial_{n}^2
u_h\) (i.e., the Hessian of the function \(u_h\) contracted twice with the direction of the normal vector) of the finite element function interpolated to the quadrature point with index q_point
after a call to evaluate(EvaluationFlags::hessians)
. Compared to the case when computing the full Hessian, some operations can be saved when only the normal Hessian is requested.
|
inherited |
Write a contribution that gets multiplied by the Hessian of the test function to the field containing the Hessians at quadrature points with index q_point
. When this function has queued information for all quadrature points and followed by a call to the function integrate(EvaluationFlags::hessians)
, the result is an array of entries, each representing the result of the integral of product of the test function Hessian multiplied by the values passed to this function.
|
inherited |
Write a contribution that gets multiplied by the Hessian of the test function times the normal projector to the field containing the Hessians at quadrature points with index q_point
. When this function has queued information for all quadrature points and followed by a call to the function integrate(EvaluationFlags::hessians)
, the result is an array of entries, each representing the result of the integral of product of the test function Hessian multiplied by the values times the normal projector passed to this function.
|
inherited |
Return the divergence of a vector-valued finite element at quadrature point number q_point
after a call to evaluate(EvaluationFlags::gradients)
.
n_components == dim
).
|
inherited |
Write a contribution that is multiplied by the divergence of the test function to the field containing the gradients at quadrature points with index q_point
. See submit_gradient() for further information.
|
inherited |
Return the symmetric gradient of the vector-valued finite element function interpolated at the quadrature point with index q_point
after a call to evaluate(EvaluationFlags::gradients)
. It corresponds to 0.5 (grad+gradT)
.
|
inherited |
Write a contribution that is multiplied by the symmetric gradient of the test function to the field containing the gradients at quadrature points with index q_point
. See submit_gradient() for further information.
|
inherited |
Return the curl of the vector field, \(\nabla \times v\) interpolated to the quadrature point index after calling evaluate(EvaluationFlags::gradients)
.
|
inherited |
Write the components of a curl containing the values on quadrature point q_point
. Access to the same data field as through get_gradient().
|
inherited |
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 over the function over the cell. If a vector-element is used, the resulting components are still separated. Moreover, if vectorization is enabled, the integral values of several cells are contained in the slots of the returned VectorizedArray field.
|
inherited |
Return the underlying MatrixFree object.
|
protectedinherited |
A unified function to read from and write into vectors based on the given template operation. It can perform the operation for read_dof_values
, distribute_local_to_global
, and set_dof_values
. It performs the operation for several vectors at a time.
|
protectedinherited |
A unified function to read from and write into vectors based on the given template operation for DG-type schemes where all degrees of freedom on cells are contiguous. It can perform the operation for read_dof_values(), distribute_local_to_global(), and set_dof_values() for several vectors at a time, depending on n_components.
|
protectedinherited |
A unified function to read from and write into vectors based on the given template operation for the case when we do not have an underlying MatrixFree object. It can perform the operation for read_dof_values
, distribute_local_to_global
, and set_dof_values
. It performs the operation for several vectors at a time, depending on n_components.
|
protectedinherited |
Apply hanging-node constraints.
|
inherited |
Sets the pointers for values, gradients, hessians to the central scratch_data_array inside the given scratch array, for a given number of components as provided by one of the derived classes.
|
inherited |
Initialize the indices related to the given face description. Used during the reinit() call of the FEFaceEvaluation class.
|
inherited |
Return an ArrayView to internal memory for temporary use. Note that some of this memory is overwritten during evaluate() and integrate() calls so do not assume it to be stable over those calls. The maximum size you can write into is 3*dofs_per_cell+2*n_q_points.
|
inherited |
Return the determinant of the Jacobian from the unit to the real cell times the quadrature weight.
|
inherited |
Returns the q-th quadrature point on the face in real coordinates stored in MappingInfo.
|
inherited |
Return the inverse and transposed version \(J^{-\mathrm T}\) of the Jacobian of the mapping between the unit to the real cell defined as \(J_{ij} = d x_i / d\hat x_j\). The \((i,j)\) entry of the returned tensor contains \(d\hat x_j/dx_i\), i.e., columns refer to reference space coordinates and rows to real cell coordinates. Thus, the returned tensor represents a covariant transformation, which is used in the FEEvaluationBase::get_gradient() function to transform the unit cell gradients to gradients on the real cell by a multiplication \(J^{-\mathrm T} \hat{\nabla} u_h\).
|
inherited |
Return the unit normal vector on a face. Note that both sides of a face use the same orientation of the normal vector: For the faces enumerated as interior
in FaceToCellTopology and selected with the is_interior_face=true
flag of the constructor, this corresponds to the outer normal vector, whereas for faces enumerated as exterior
in FaceToCellTopology and selected with the is_interior_face=false
flag of the constructor, the normal points into the element as a consequence of the single normal vector.
is_face == true
.
|
inherited |
Same as normal_vector(const unsigned int q_point)
.
|
inherited |
Return a read-only pointer to the first field of the dof values. This is the data field the read_dof_values() functions write into. First come the dof values for the first component, then all values for the second component, and so on. This is related to the internal data structures used in this class. In general, it is safer to use the get_dof_value() function instead.
|
inherited |
Return a read and write pointer to the first field of the dof values. This is the data field the read_dof_values() functions write into. First come the dof values for the first component, then all values for the second component, and so on. This is related to the internal data structures used in this class. In general, it is safer to use the get_dof_value() function instead.
|
inherited |
Return a read-only pointer to the first field of function values on quadrature points. First come the function values on all quadrature points for the first component, then all values for the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate
only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_value() function instead, which does all the transformation internally.
|
inherited |
Return a read and write pointer to the first field of function values on quadrature points. First come the function values on all quadrature points for the first component, then all values for the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate
only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_value() function instead, which does all the transformation internally.
|
inherited |
Return a read-only pointer to the first field of function gradients on quadrature points. First comes the x-component of the gradient for the first component on all quadrature points, then the y-component, and so on. Next comes the x-component of the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate
only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_gradient() function instead, which does all the transformation internally.
|
inherited |
Return a read and write pointer to the first field of function gradients on quadrature points. First comes the x-component of the gradient for the first component on all quadrature points, then the y-component, and so on. Next comes the x-component of the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate
only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_gradient() function instead, which does all the transformation internally.
|
inherited |
Return a read-only pointer to the first field of function hessians on quadrature points. First comes the xx-component of the hessian for the first component on all quadrature points, then the yy-component, zz-component in (3d), then the xy-component, and so on. Next comes the xx-component of the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate
only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_laplacian() or get_hessian() functions instead, which does all the transformation internally.
|
inherited |
Return a read and write pointer to the first field of function hessians on quadrature points. First comes the xx-component of the hessian for the first component on all quadrature points, then the yy-component, zz-component in (3d), then the xy-component, and so on. Next comes the xx-component of the second component, and so on. This is related to the internal data structures used in this class. The raw data after a call to evaluate
only contains unit cell operations, so possible transformations, quadrature weights etc. must be applied manually. In general, it is safer to use the get_laplacian() or get_hessian() functions instead, which does all the transformation internally.
|
inherited |
Return the index offset within the geometry fields for the cell the reinit()
function has been called for. This index can be used to access an index into a field that has the same compression behavior as the Jacobian of the geometry, e.g., to store an effective coefficient tensors that combines a coefficient with the geometry for lower memory transfer as the available data fields.
|
inherited |
Return the type of the cell the reinit()
function has been called for. Valid values are cartesian
for Cartesian cells (which allows for considerable data compression), affine
for cells with affine mappings, and general
for general cells without any compressed storage applied.
|
inherited |
Return a reference to the ShapeInfo object currently in use.
|
inherited |
Return a reference to the DoFInfo object currently in use.
|
inherited |
Return the numbering of local degrees of freedom within the evaluation routines of FEEvaluation in terms of the standard numbering on finite elements.
|
inherited |
Return the number of the quadrature formula of the present cell.
|
inherited |
Return index of the current cell or face.
|
inherited |
Return the active FE index for this class for efficient indexing in the hp-case.
|
inherited |
Return the active quadrature index for this class for efficient indexing in the hp-case.
|
inherited |
Return the first selected component in a multi-component system.
|
inherited |
If is_face is true, this function returns the face number of the given lane within a cell for the face this object was initialized to. On cells where the face makes no sense, an exception is thrown.
dof_access_index == dof_access_cell
and is_interior_face == false
.
|
inherited |
If is_face is true, this function returns the index of a subface along a face in case this object was initialized to a face with hanging nodes. On faces On cells where the face makes no sense, an exception is thrown.
|
inherited |
If is_face is true, this function returns for a given lane the orientation index within an array of orientations as stored in ShapeInfo for unknowns and quadrature points.
dof_access_index == dof_access_cell
and is_interior_face == false
.
|
inherited |
Return the current index in the access to compressed indices.
|
inherited |
Return whether this object was set up to work on what is called "interior" faces in the evaluation context.
|
inlineinherited |
Return the id of the cells this FEEvaluation or FEFaceEvaluation is associated with.
Definition at line 496 of file fe_evaluation_data.h.
|
inlineinherited |
Return the id of the faces this FEFaceEvaluation is associated with.
Definition at line 510 of file fe_evaluation_data.h.
|
inlineinherited |
Return the id of the cell/face batch this FEEvaluation/FEFaceEvaluation is associated with.
Definition at line 524 of file fe_evaluation_data.h.
|
inlineinherited |
Return the id of the cells/faces this FEEvaluation/FEFaceEvaluation is associated with.
Definition at line 539 of file fe_evaluation_data.h.
|
inherited |
Return an object that can be thought of as an array containing all indices from zero to n_quadrature_points
. This allows to write code using range-based for loops.
|
inherited |
Provides a unified interface to access data in a vector of VectorizedArray fields of length MatrixFree::n_cell_batches() + MatrixFree::n_ghost_cell_batches() for cell data. It is implemented both for cells and faces (access data to the associated cell).
The underlying type can be both the Number type as given by the class template (i.e., VectorizedArrayType) as well as std::array with as many entries as n_lanes.
|
inherited |
Provides a unified interface to set data in a vector of VectorizedArray fields of length MatrixFree::n_cell_batches() + MatrixFree::n_ghost_cell_batches() for cell data. It is implemented both for cells and faces (access data to the associated cell).
|
inherited |
Provides a unified interface to access data in a vector of VectorizedArray fields of length MatrixFree::n_inner_face_batches() + MatrixFree::n_boundary_face_batches() + MatrixFree::n_ghost_inner_face_batches() for face data.
|
inherited |
Provides a unified interface to set data in a vector of VectorizedArray fields of length MatrixFree::n_inner_face_batches() + MatrixFree::n_boundary_face_batches() + MatrixFree::n_ghost_inner_face_batches() for face data.
|
staticconstexpr |
The dimension given as template argument.
Definition at line 1883 of file fe_evaluation.h.
|
staticconstexpr |
The number of solution components of the evaluator given as template argument.
Definition at line 1889 of file fe_evaluation.h.
|
staticconstexpr |
The number of lanes of the template argument VectorizedArrayType.
Definition at line 1894 of file fe_evaluation.h.
|
staticconstexpr |
The static number of quadrature points determined from the given template argument n_q_points_1d
taken to the power of dim-1.
n_q_points
, can be different if fe_degree=-1
is given and run-time loop lengths are used rather than compile time ones. Definition at line 1905 of file fe_evaluation.h.
|
staticconstexpr |
The static number of quadrature points on a cell with the same quadrature formula.
static_n_q_points
variable. Definition at line 1916 of file fe_evaluation.h.
|
staticconstexpr |
The static number of degrees of freedom of a scalar component determined from the given template argument fe_degree
.
dofs_per_component
can be different if fe_degree=-1
is given. Definition at line 1927 of file fe_evaluation.h.
|
staticconstexpr |
The static number of degrees of freedom of all components determined from the given template argument fe_degree
.
dofs_per_cell
can be different if fe_degree=-1
is given. Definition at line 1938 of file fe_evaluation.h.
|
staticconstexpr |
The static number of degrees of freedom of all components determined from the given template argument fe_degree
.
dofs_per_cell
can be different if fe_degree=-1
is given. Definition at line 1949 of file fe_evaluation.h.
const unsigned int FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::dofs_per_component |
The number of degrees of freedom of a single component on the cell for the underlying evaluation object. Usually close to static_dofs_per_component, but the number depends on the actual element selected and is thus not static.
Definition at line 2243 of file fe_evaluation.h.
const unsigned int FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::dofs_per_cell |
The number of degrees of freedom on the cell accumulated over all components in the current evaluation object. Usually close to static_dofs_per_cell = static_dofs_per_component*n_components, but the number depends on the actual element selected and is thus not static.
Definition at line 2251 of file fe_evaluation.h.
const unsigned int FEFaceEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::n_q_points |
The number of quadrature points in use. If the number of quadrature points in 1d is given as a template, this number is simply the dim-1
-th power of that value. If the element degree is set to -1 (dynamic selection of element degree), the static value of quadrature points is inaccurate and this value must be used instead.
Definition at line 2260 of file fe_evaluation.h.
|
protectedinherited |
This is the general array for all data fields.
Definition at line 787 of file fe_evaluation.h.
|
protectedinherited |
A pointer to the underlying data.
Definition at line 792 of file fe_evaluation.h.
|
mutableprotectedinherited |
A temporary data structure necessary to read degrees of freedom when no MatrixFree object was given at initialization.
Definition at line 798 of file fe_evaluation.h.
|
protectedinherited |
A pointer to the unit cell shape data, i.e., values, gradients and Hessians in 1d at the quadrature points that constitute the tensor product. Also contained in matrix_info, but it simplifies code if we store a reference to it.
Definition at line 661 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the underlying DoF indices and constraint description for the component specified at construction. Also contained in matrix_info, but it simplifies code if we store a reference to it.
Definition at line 668 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the underlying transformation data from unit to real cells for the given quadrature formula specified at construction. Also contained in matrix_info, but it simplifies code if we store a reference to it.
Definition at line 676 of file fe_evaluation_data.h.
|
protectedinherited |
The number of the quadrature formula of the present cell among all quadrature formulas available in the MatrixFree objects pass to derived classes.
Definition at line 683 of file fe_evaluation_data.h.
|
protectedinherited |
Stores the number of components in the finite element as detected in the MatrixFree storage class for comparison with the template argument.
Definition at line 689 of file fe_evaluation_data.h.
|
protectedinherited |
For a FiniteElement with more than one base element, select at which component this data structure should start.
Definition at line 695 of file fe_evaluation_data.h.
|
protectedinherited |
The active FE index for this class for efficient indexing in the hp-case.
Definition at line 700 of file fe_evaluation_data.h.
|
protectedinherited |
The active quadrature index for this class for efficient indexing in the hp-case.
Definition at line 706 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the underlying quadrature formula specified at construction. Also contained in mapping_data, but it simplifies code if we store a reference to it.
Definition at line 713 of file fe_evaluation_data.h.
|
protectedinherited |
The number of quadrature points in the current evaluation context.
Definition at line 718 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the quadrature-point information of the present cell. Only set to a useful value if on a non-Cartesian cell.
Definition at line 724 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the inverse transpose Jacobian information of the present cell. Only the first inverse transpose Jacobian (q_point = 0) is set for Cartesian/affine cells, and the actual Jacobian is stored at index 1 instead. For faces on hypercube elements, the derivatives are reorder s.t. the derivative orthogonal to the face is stored last, i.e for dim = 3 and face_no = 0 or 1, the derivatives are ordered as [dy, dz, dx], face_no = 2 or 3: [dz, dx, dy], and face_no = 5 or 6: [dx, dy, dz]. If the Jacobian also is stored, the components are instead reordered in the same way. Filled from MappingInfoStorage.jacobians in include/deal.II/matrix_free/mapping_info.templates.h
Definition at line 738 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the gradients of the inverse Jacobian transformation of the present cell.
Definition at line 745 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the gradients of the Jacobian transformation of the present cell.
Definition at line 752 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the Jacobian determinant of the present cell. If on a Cartesian cell or on a cell with constant Jacobian, this is just the Jacobian determinant, otherwise the Jacobian determinant times the quadrature weight.
Definition at line 760 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the normal vectors at faces.
Definition at line 765 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the normal vectors times the Jacobian at faces.
Definition at line 770 of file fe_evaluation_data.h.
|
protectedinherited |
A pointer to the quadrature weights of the underlying quadrature formula.
Definition at line 775 of file fe_evaluation_data.h.
|
mutableprotectedinherited |
This is the user-visible part of FEEvaluationBase::scratch_data_array, only showing the part that can be consumed by various users. It is set during the allocation of the internal data structures in FEEvaluationBase.
Definition at line 783 of file fe_evaluation_data.h.
|
protectedinherited |
This field stores the values for local degrees of freedom (e.g. after reading out from a vector but before applying unit cell transformations or before distributing them into a result vector). The methods get_dof_value() and submit_dof_value() read from or write to this field.
The values of this array are stored in the start section of scratch_data_array
. Due to its access as a thread local memory, the memory can get reused between different calls.
Definition at line 795 of file fe_evaluation_data.h.
|
protectedinherited |
This field stores the values of the finite element function on quadrature points after applying unit cell transformations or before integrating. The methods get_value() and submit_value() access this field.
The values of this array are stored in the start section of scratch_data_array
. Due to its access as a thread local memory, the memory can get reused between different calls.
Definition at line 806 of file fe_evaluation_data.h.
|
protectedinherited |
This field stores the values of the finite element function on quadrature points after applying unit cell transformations or before integrating. This field is accessed when performing the contravariant Piola transform for gradients on general cells. This is done by the functions get_gradient() and submit_gradient() when using a H(div)- conforming finite element such as FE_RaviartThomasNodal.
The values of this array are stored in the start section of scratch_data_array
. Due to its access as a thread local memory, the memory can get reused between different calls.
Definition at line 820 of file fe_evaluation_data.h.
|
protectedinherited |
This field stores the gradients of the finite element function on quadrature points after applying unit cell transformations or before integrating. The methods get_gradient() and submit_gradient() (as well as some specializations like get_symmetric_gradient() or get_divergence()) access this field.
The values of this array are stored in the start section of scratch_data_array
. Due to its access as a thread local memory, the memory can get reused between different calls.
Definition at line 833 of file fe_evaluation_data.h.
|
protectedinherited |
This field stores the gradients of the finite element function on quadrature points after applying unit cell transformations or before integrating. The methods get_hessian() and submit_hessian() (as well as some specializations like get_hessian_diagonal() or get_laplacian()) access this field for general cell/face types.
The values of this array are stored in the start section of scratch_data_array
. Due to its access as a thread local memory, the memory can get reused between different calls.
Definition at line 846 of file fe_evaluation_data.h.
|
protectedinherited |
This field stores the Hessians of the finite element function on quadrature points after applying unit cell transformations. The methods get_hessian(), get_laplacian(), get_hessian_diagonal() access this field.
The values of this array are stored in the start section of scratch_data_array
. Due to its access as a thread local memory, the memory can get reused between different calls.
Definition at line 857 of file fe_evaluation_data.h.
|
protectedinherited |
Debug information to track whether the reinit() function of derived classes was called.
Definition at line 863 of file fe_evaluation_data.h.
|
protectedinherited |
Debug information to track whether dof values have been initialized before accessed. Used to control exceptions when uninitialized data is used.
Definition at line 870 of file fe_evaluation_data.h.
|
protectedinherited |
Debug information to track whether values on quadrature points have been initialized before accessed. Used to control exceptions when uninitialized data is used.
Definition at line 877 of file fe_evaluation_data.h.
|
protectedinherited |
Debug information to track whether gradients on quadrature points have been initialized before accessed. Used to control exceptions when uninitialized data is used.
Definition at line 884 of file fe_evaluation_data.h.
|
protectedinherited |
Debug information to track whether Hessians on quadrature points have been initialized before accessed. Used to control exceptions when uninitialized data is used.
Definition at line 891 of file fe_evaluation_data.h.
|
protectedinherited |
Debug information to track whether values on quadrature points have been submitted for integration before the integration is actually stared. Used to control exceptions when uninitialized data is used.
Definition at line 898 of file fe_evaluation_data.h.
|
protectedinherited |
Debug information to track whether gradients on quadrature points have been submitted for integration before the integration is actually stared. Used to control exceptions when uninitialized data is used.
Definition at line 905 of file fe_evaluation_data.h.
|
protectedinherited |
Debug information to track whether hessians on quadrature points have been submitted for integration before the integration is actually stared. Used to control exceptions when uninitialized data is used.
Definition at line 912 of file fe_evaluation_data.h.
|
protectedinherited |
After a call to reinit(), stores the number of the cell we are currently working with.
Definition at line 918 of file fe_evaluation_data.h.
|
protectedinherited |
Flag holding information whether a face is an interior or exterior face according to the defined direction of the normal. For cells it defines if the dof values should be read from the actual cell corresponding to the interior face or the neighboring cell corresponding to the exterior face.
Definition at line 926 of file fe_evaluation_data.h.
|
protectedinherited |
Stores the index an FEFaceEvaluation object is currently pointing into (interior face, exterior face, data associated with cell).
Definition at line 932 of file fe_evaluation_data.h.
|
protectedinherited |
Stores the (non-vectorized) number of faces within cells in case of ECL for the exterior cells where the single number face_no
is not enough to cover all cases.
dof_access_index == dof_access_cell
and is_interior_face == false
. Definition at line 942 of file fe_evaluation_data.h.
|
protectedinherited |
Store the orientation of the neighbor's faces with respect to the current cell for the case of exterior faces on ECL with possibly different orientations behind different cells.
dof_access_index == dof_access_cell
and is_interior_face == false
. Definition at line 952 of file fe_evaluation_data.h.
|
protectedinherited |
Stores the subface index of the given face. Usually, this variable takes the value numbers::invalid_unsigned_int to indicate integration over the full face, but in case the current physical face has a neighbor that is more refined, it is a subface and must scale the entries in ShapeInfo appropriately.
Definition at line 961 of file fe_evaluation_data.h.
|
protectedinherited |
Stores the type of the cell we are currently working with after a call to reinit(). Valid values are cartesian
, affine
and general
, which have different implications on how the Jacobian transformations are stored internally in MappingInfo.
Definition at line 969 of file fe_evaluation_data.h.
|
protectedinherited |
Stores the (non-vectorized) id of the cells or faces this object is initialized to. Relevant for ECL.
Definition at line 975 of file fe_evaluation_data.h.
|
protectedinherited |
Stores the (non-vectorized) id of the cells or faces this object is initialized to. Relevant for ECL.
Definition at line 981 of file fe_evaluation_data.h.
|
protectedinherited |
Geometry data that can be generated FEValues on the fly with the respective constructor, as an alternative to the entry point with MatrixFree.
Definition at line 990 of file fe_evaluation_data.h.
|
protectedinherited |
Bool indicating if the divergence is requested. Used internally in the case of the Piola transform.
Definition at line 996 of file fe_evaluation_data.h.