deal.II version GIT relicensing-2287-g6548a49e0a 2024-12-20 18:30:00+00:00
|
#include <deal.II/matrix_free/fe_evaluation.h>
Public Types | |
using | BaseClass = FEEvaluationBase< dim, n_components_, Number, false, 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 | |
FEEvaluation (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, 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) | |
FEEvaluation (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free, const std::pair< unsigned int, unsigned int > &range, const unsigned int dof_no=0, const unsigned int quad_no=0, const unsigned int first_selected_component=0) | |
FEEvaluation (const Mapping< dim > &mapping, const FiniteElement< dim > &fe, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags, const unsigned int first_selected_component=0) | |
FEEvaluation (const FiniteElement< dim > &fe, const Quadrature< 1 > &quadrature, const UpdateFlags update_flags, const unsigned int first_selected_component=0) | |
FEEvaluation (const FiniteElement< dim > &fe, const FEEvaluationData< dim, VectorizedArrayType, false > &other, const unsigned int first_selected_component=0) | |
FEEvaluation (const FEEvaluation &other) | |
FEEvaluation & | operator= (const FEEvaluation &other) |
void | reinit (const unsigned int cell_batch_index) |
void | reinit (const std::array< unsigned int, n_lanes > &cell_ids) |
template<bool level_dof_access> | |
void | reinit (const TriaIterator< DoFCellAccessor< dim, dim, level_dof_access > > &cell) |
void | reinit (const typename Triangulation< dim >::cell_iterator &cell) |
void | evaluate (const EvaluationFlags::EvaluationFlags evaluation_flag) |
void | evaluate (const VectorizedArrayType *values_array, 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) |
void | integrate (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) |
std_cxx20::ranges::iota_view< unsigned int, unsigned int > | dof_indices () const |
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 | |
void | read_dof_values (const VectorType &src, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip()) |
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()) |
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 |
void | set_dof_values (VectorType &dst, const unsigned int first_index=0, const std::bitset< n_lanes > &mask=std::bitset< n_lanes >().flip()) const |
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) |
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) |
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) |
VectorizedArrayType | get_divergence (const unsigned int q_point) const |
void | submit_divergence (const VectorizedArrayType div_in, const unsigned int q_point) |
SymmetricTensor< 2, dim, VectorizedArrayType > | get_symmetric_gradient (const unsigned int q_point) const |
void | submit_symmetric_gradient (const SymmetricTensor< 2, dim, VectorizedArrayType > grad_in, const unsigned int q_point) |
Tensor< 1,(dim==2 ? 1 :dim), VectorizedArrayType > | get_curl (const unsigned int q_point) const |
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_dofs_per_component |
static constexpr unsigned int | tensor_dofs_per_cell |
static constexpr unsigned int | static_dofs_per_cell |
Protected Member Functions | |
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 |
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 |
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 |
Private Member Functions | |
void | check_template_arguments (const unsigned int fe_no, const unsigned int first_selected_component) |
The class that provides all functions necessary to evaluate functions at quadrature points and cell integrations. In functionality, this class is similar to FEValues, however, it includes a lot of specialized functions that make it much faster (between 5 and 500, depending on the polynomial degree). For evaluation of face terms in DG, see the class FEFaceEvaluation.
The first and foremost way of usage is to initialize this class from a MatrixFree object that caches everything related to the degrees of freedom and the mapping information. This way, it is possible to use vectorization for applying a differential operator for several cells at once.
The capabilities of FEEvaluation span a large spectrum of integration tasks for weak forms. In general, there are two classes of tasks that get done. One is the evaluate
path that interpolates from a solution vector to quadrature points:
Likewise, a gradient of the finite element solution represented by vector
can be interpolated to the quadrature points by fe_eval.get_gradient(q)
. The combination of read_dof_values(), evaluate() and get_value() is similar to what FEValues::get_function_values or FEValues::get_function_gradients does, but it is in general much faster because it makes use of the tensor product, see the description of the evaluation routines below, and can do this operation for several cells at once through vectorization.
The second class of tasks done by FEEvaluation are integration tasks for right hand sides. In finite element computations, these typically consist of multiplying a quantity on quadrature points (a function value, or a field interpolated by the finite element space itself) by a set of test functions and integrating over the cell through summation of the values in each quadrature point, multiplied by the quadrature weight and the Jacobian determinant of the transformation. If a generic Function object is given and we want to compute \(v_i = \int_\Omega \varphi_i f dx\), this is done by the following cell-wise integration:
In this code, the call to fe_eval.submit_value()
prepares for the multiplication by the test function prior to the actual integration (inside the submit call, the value to be tested is also multiplied by the determinant of the Jacobian and the quadrature weight). In the integrate()
call, an integral contribution tested by each basis function underlying the FEEvaluation object (e.g. the four linear shape functions of FE_Q<2>(1) in 2d) is computed, which gives the vector entries to be summed into the dst
vector. Note that the above code needs to explicitly loop over the components in the vectorized array for evaluating the function, which is necessary for interfacing with a generic Function object with double arguments. Simple functions can also be implemented in VectorizedArray form directly as VectorizedArray provides the basic math operations.
For evaluating a bilinear form, the evaluation on a source vector is combined with the integration involving test functions that get written into a result vector. This setting is the context of matrix-free operator evaluation and explained in the step-37 and step-48 tutorial programs.
Note that the two vector accesses through FEEvaluation::read_dof_values and FEEvaluation::distribute_local_to_global resolve constraints on the fly, based on the AffineConstraints object specified at the MatrixFree::reinit() call. In case the values in the degrees of freedom are of interest (usually only the values in quadrature points are necessary), these can be accessed through FEEvaluation::get_dof_value(i), where i is the index of the basis function. Note that the numbering of the degrees of freedom for continuous elements in FEEvaluation is different from the ordering in FE_Q (or FEValues) because FEEvaluation needs to access them in lexicographic order, which is the ordering used in FE_DGQ, for instance. Re-indexing would be too expensive because the access inside evaluate() and integrate() is on the critical path in the tensorial evaluation parts. An alternative to filling the DoF values by read_dof_values() before an evaluate() call is to manually assign a value by a set_dof_value() call. Likewise, if the local result of integration should be further processed rather than scattered into a vector by distribute_local_to_global(), one can access it by get_dof_value() after an integrate() call. An example for using the values of an integral in a different context is fast assembly of matrices as shown in the next subsection.
For most operator evaluation tasks that repeatedly go through the mesh, the realization by MatrixFree that combines pre-computed data for the mapping (Jacobian transformations for the geometry description) with on-the-fly evaluation of basis functions is the most efficient way of doing things. In other words, the framework selects a trade-off between memory usage and initialization of objects that is suitable for replacement of matrix-vector products or explicit time integration in a matrix-free way.
The second form of usage is to initialize FEEvaluation from geometry information generated by FEValues. This allows to apply the integration loops on the fly without prior initialization of MatrixFree objects. This can be useful when the memory and initialization cost of MatrixFree is not acceptable, e.g. when a different number of quadrature points should be used for one single evaluation in error computation. Also, when using the routines of this class to assemble matrices the trade-off implied by the MatrixFree class may not be desired. In such a case, the cost to initialize the necessary geometry data on the fly is comparably low and thus avoiding a global object MatrixFree can be useful. When used in this way, reinit methods reminiscent from FEValues with a cell iterator are used. However, note that this model results in working on a single cell at a time, with geometry data duplicated in all components of the vectorized array. Thus, vectorization is only useful when it can apply the same operation on different data, e.g. when performing matrix assembly.
As an example, consider the following code to assemble the contributions to the Laplace matrix:
This code generates the columns of the cell matrix with the loop over i
above. The way this is done is the following: FEEvaluation's routines focus on the evaluation of finite element operators, so for computing a cell matrix out of an operator evaluation it is applied to all the unit vectors on the cell. Applying the operator on a unit vector might seem inefficient but the evaluation routines used here are so quick that they still work much faster than what is possible with FEValues. In particular, the complexity is (fe_degree+1)2*dim+1
rather than (fe_degree+1)3*dim
.
Due to vectorization, we can generate matrix columns for several unit vectors at a time (e.g. 4). The variable n_items
make sure that we do the last iteration where the number of cell dofs is not divisible by the vectorization length correctly. Also note that we need to get the internal dof numbering applied by fe_eval because FEEvaluation internally uses a lexicographic numbering of degrees of freedom as explained above.
The temporary data for holding the solution values on the local degrees of freedom as well as the interpolated values, gradients, and Hessians on quadrature points is a scratch array provided by MatrixFree::acquire_scratch_data() that is re-used between different calls to FEEvaluation. Therefore, constructing an FEEvaluation object is typically cheap and does not involve any expensive operation. Only a few dozen pointers to the actual data fields are set during construction. Therefore, no negative performance impact arises when creating an FEEvaluation several times per loop, such as at the top of a local_cell_operation
operation that is split in small chunks for a parallel for loop, obviating a separate scratch data field for parallel loops as necessary in the loop of WorkStream
.
When using the FEEvaluation class in multithreaded mode, the thread local storage of the scratch data in MatrixFree automatically makes sure that each thread gets it private data array. Note, however, that deal.II must be compiled with thread support also when all the thread parallelization is provided externally and not done via deal.II's routines, such as OpenMP. This is because deal.II needs to know the notation of thread local storage. The FEEvaluation kernels have been verified to work within OpenMP loops.
This class is designed to perform all arithmetics on single-instruction multiple-data (SIMD) instructions present on modern CPUs by explicit vectorization, which are made available in deal.II through the class VectorizedArray, using the widest vector width available at configure/compile time. In order to keep programs flexible, FEEvaluation always applies vectorization over several elements. This is often the best compromise because computations on different elements are usually independent in the finite element method (except of course the process of adding an integral contribution to a global residual vector), also in more complicated scenarios: Stabilization parameter can e.g. be defined as the maximum of some quantities on all quadrature points of a cell divided by the cell's volume, but without locally mixing the results with neighbors. Using the terminology from computer architecture, the design of FEEvaluation relies on not doing any cross-lane data exchange when operating on the cell in typical integration scenarios.
When the number of cells in the problem is not a multiple of the number of array elements in the SIMD vector, the implementation of FEEvaluation fills in some dummy entries in the unused SIMD lanes and carries them around nonetheless, a choice made necessary since the length of VectorizedArray is fixed at compile time. Yet, this approach most often results in superior code as compared to an auto-vectorization setup where an alternative unvectorized code path would be necessary next to the vectorized version to be used on fully populated lanes, together with a dispatch mechanism. In read_dof_values
, the empty lanes resulting from a reinit() call to an incomplete batch of cells are set to zero, whereas distribute_local_to_global
or set_dof_values
simply ignores the content in the empty lanes. The number of actually filled SIMD lanes can by queried by MatrixFree::n_active_entries_per_cell_batch() or MatrixFree::n_active_entries_per_face_batch().
Obviously, the computations performed on the artificial lanes (without real data) should never be mixed with valid results. The contract in using this class is that the user makes sure that lanes are not crossed in user code, in particular since it is not clear a priori which cells are going to be put together in vectorization. For example, results on an element should not be added to results on other elements except through the global vector access methods or by access that is masked by MatrixFree::n_active_entries_per_cell_batch(). No guarantee can be made that results on artificial lanes will always be zero that can safely be added to other results: The data on JxW or Jacobians is copied from the last valid lane in order to avoid division by zero that could trigger floating point exceptions or trouble in other situations.
This class contains specialized evaluation routines for elements based on tensor-product quadrature formulas and tensor-product-like shape functions, including standard FE_Q or FE_DGQ elements and quadrature points symmetric around 0.5 (like Gauss quadrature), FE_DGP elements based on truncated tensor products as well as the faster case of Gauss-Lobatto elements with Gauss-Lobatto quadrature which give diagonal mass matrices and quicker evaluation internally. The main benefit of this class is the evaluation of all shape functions in all quadrature or integration over all shape functions in dim (fe_degree+1)dim+1
operations instead of the slower (fe_degree+1)2*dim
complexity in the evaluation routines of FEValues. This is done by an algorithm called sum factorization which factors out constant factors during the evaluation along a coordinate direction. This algorithm is the basis of many spectral element algorithms.
Note that many of the operations available through this class are inherited from the base class FEEvaluationBase, in particular reading from and writing to vectors. Furthermore, functionality to access to values, gradients and Hessians of the finite element function at quadrature points is inherited.
This class assumes that the shape functions of the FiniteElement under consideration do not depend on the geometry of the cells in real space. Currently, other finite elements cannot be treated with the matrix-free concept.
The class FEEvaluation as two usage models. The first usage model is to specify the polynomial degree as a template parameter. This guarantees maximum efficiency: The evaluation with sum factorization performs a number of nested short 1d loops of length equal to the polynomial degree plus one. If the loop bounds are known at compile time, the compiler can unroll loops as deemed most efficient by its heuristics. At least the innermost loop is almost always completely unrolled, avoiding the loop overhead.
However, carrying the polynomial degree (and the number of quadrature points) as a template parameter makes things more complicated in codes where different polynomial degrees should be considered, e.g. in application codes where the polynomial degree is given through an input file. The second usage model is to rely on pre-compiled code for polynomial degrees. While a user code can use different functions for the cells (that get e.g. invoked by some dynamic dispatch mechanism for the various degree templates), deal.II also supports usage of this class based on the information in the element passed to the initialization. For this usage model, set the template parameter for the polynomial degree to -1 and choose an arbitrary number for the number of quadrature points. That code part contains pre-compiled templated code for polynomial degrees between 1 and 6 and common quadrature formulas, which runs almost as fast as the templated version. In case the chosen degree is not pre-compiled, an evaluator object with template specialization for -1 is invoked that runs according to run-time bounds.
An overview of the performance of FEEvaluation is given in the following figure. It considers the time spent per degree of freedom for evaluating the Laplacian with continuous finite elements using a code similar to the step-37 tutorial program for single-precision arithmetics. The time is based on an experiment on a single core of an Intel Xeon E5-2687W v4, running at 3.4 GHz and measured at problem sizes of around 10 million. The plot lists the computational time (around 0.1 seconds) divided by the number of degrees freedom.
The figure shows that the templated computational kernels are between 2.5 and 3 times faster than the non-templated ones. The fastest turnaround on this setup is for polynomial degree 5 at 7.4e-9 seconds per degree of freedom or 134 million degrees of freedom per second - on a single core. The non-templated version is also fastest at polynomial degree 5 with 2.1e-9 seconds per degree of freedom or 48 million degrees of freedom per second. Note that using FEEvaluation with template degree=-1
selects the fast path for degrees between one and six, and the slow path for other degrees.
It is also possible to pre-compile the code in FEEvaluation for a different maximal polynomial degree. This is controlled by the class internal::FEEvaluationFactory and the implementation in include/deal.II/matrix_free/evaluation_template_factory.templates.h
. By setting the macro FE_EVAL_FACTORY_DEGREE_MAX
to the desired integer and instantiating the classes FEEvaluationFactory and FEFaceEvaluationFactory (the latter for FEFaceEvaluation) creates paths to templated functions for a possibly larger set of degrees. This can both be set when configuring deal.II by passing the flag -D FE_EVAL_FACTORY_DEGREE_MAX=8
(in case you want to compile all degrees up to eight; recommended setting) or by compiling evaluation_template_factory.templates.h
and evaluation_template_face_factory.templates.h
with the FE_EVAL_FACTORY_DEGREE_MAX
overridden to the desired value. In the second option, symbols will be available twice, and it depends on your linker and dynamic library loader whether the user-specified setting takes precedence; use LD_PRELOAD
to select the desired library. You can check if fast evaluation/integration for a given degree/n_quadrature_points pair by calling FEEvaluation::fast_evaluation_supported() or FEFaceEvaluation::fast_evaluation_supported().
FEEvaluation also allows for treating vector-valued problems through a template parameter on the number of components:
If used this way, the components can be gathered from several components of an std::vector<VectorType>
through the call
where the 0 means that the vectors starting from the zeroth vector in the std::vector
should be used, src[0], src[1], ..., src[n_components-1]
.
An alternative way for reading multi-component systems is possible if the DoFHandler underlying the MatrixFree data is based on an FESystem of n_components
entries. In that case, a single vector is provided for the read_dof_values() and distribute_local_to_global() calls.
An important property of FEEvaluation in multi-component systems is the layout of multiple components in the get_value(), get_gradient(), or get_dof_value() calls. In this case, instead of a scalar return field VectorizedArray<double> a tensor is returned,
In a similar vein, the submit_value() and submit_gradient() calls take tensors of values. Note that there exist specializations of these types for n_components=1
and n_components=dim
. In the scalar case, these provide the scalar return types described above. In the vector-valued case, the gradient is converted from Tensor<1,dim,Tensor<1,dim,VectorizedArray<double> > >
to Tensor<2,dim,VectorizedArray<double> >
. Furthermore, additional operations such as the divergence or curl are available.
In case different shape functions are combined, for example mixed finite element formulations in Stokes flow, two FEEvaluation objects are created, one for the velocity and one for the pressure. Those are then combined on quadrature points:
This code assumes that a BlockVector of two components describes the velocity and pressure components, respectively. For identifying the different DoFHandler objects for velocity and pressure, the second argument to the FEEvaluation objects specify the respective component 0 for velocity and 1 for pressure. For further examples of vector-valued problems, the deal.II test suite includes a few additional examples as well, e.g. the Stokes operator described above is found at https://github.com/dealii/dealii/blob/master/tests/matrix_free/matrix_vector_stokes_noflux.cc
The design of FEEvaluation and MatrixFree separates the geometry from the basis functions. Therefore, several DoFHandler objects (or the same DoFHandler equipped with different constraint objects) can share the same geometry information like in the Stokes example above. All geometry is cached once in MatrixFree, so FEEvaluation does not need to do expensive initialization calls and rather sets a few pointers. This realization is based on the idea that the geometry information is needed only once also when several fields are evaluated, in a departure from FEValues which sets up the internal mapping data for each field. If for example a multi-component PDE involves the shape values on one component and the shape gradient on the other, no efficiency is lost if both are based on the same MatrixFree object where the update flags specify that both update_values
, update_gradients
, and update_jxw_values
are given. The selection of desired quantities of shape values is through the flags in the evaluate() or integrate calls and the access at quadrature points:
In the loop over quadrature points, one can ask any of the two FEEvaluation objects — it does not really matter which one because they only keep pointers to the quadrature point data — to provide the quadrature point location.
This observation also translates to the case when different differential operators are implemented in a program, for example the action of a mass matrix for one phase of the algorithm and the action of a stiffness matrix in another one. Only a single MatrixFree object is necessary, maintaining full efficiency by using different local functions with the respective implementation in separate FEEvaluation objects. In other words, a user does not need to bother about being conservative when providing update_flags to the initialization of MatrixFree for efficiency reasons - no overhead incurs inside FEEvaluation, except for at most one or two more if
statements inside the FEEvaluation::reinit() call. Rather, the largest set of flags necessary among all calls is perfectly fine from an efficiency point of view.
For the combination of different fields, including different solution vectors that come from different time steps, it is mandatory that all FEEvaluation objects share the same MatrixFree object. This is because the way cells are looped by MatrixFree::cell_loop() can be different for different DoFHandler or AffineConstraints arguments. More precisely, even though the layout is going to be the same in serial, there is no guarantee about the ordering for different DoFHandler/AffineConstraints in the MPI case. The reason is that the algorithm detects cells that need data exchange with MPI and those can change for different elements — FE_Q with hanging node constraints connects to more neighbors than a FE_DGQ element, for instance, and cells which need data exchange are put in different positions inside the cell loop. Of course, if the exact same DoFHandler, AffineConstraints, and options (such as the setting for thread parallelism) are set, then the order is going to be the same because the algorithm is deterministic.
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. Can be set to -1 if the degree is not known at compile time, but performance will usually be worse by a factor of 2-3. |
n_q_points_1d | Number of points in the quadrature formula in 1d, defaults to 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). Defaults to 1. |
Number | Number format, usually double or float . Defaults to double |
Definition at line 1381 of file fe_evaluation.h.
using FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::BaseClass = FEEvaluationBase<dim, n_components_, Number, false, VectorizedArrayType> |
An alias to the base class.
Definition at line 1395 of file fe_evaluation.h.
using FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::number_type = Number |
An underlying number type specified as template argument.
Definition at line 1401 of file fe_evaluation.h.
using FEEvaluation< 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 1408 of file fe_evaluation.h.
using FEEvaluation< 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 1415 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.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::FEEvaluation | ( | const MatrixFree< dim, Number, VectorizedArrayType > & | matrix_free, |
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 |
||
) |
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 |
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 component, this parameter allows for selecting the component where the current evaluation routine should start. Note that one evaluator does not support combining different shape functions in different components. In other words, the same base element of a FESystem needs to be set for the components between first_selected_component and first_selected_component+n_components_ . |
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. |
active_quad_index | If matrix_free was set up with hp::Collection objects, this parameter selects the appropriate number of the quadrature formula. |
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::FEEvaluation | ( | const MatrixFree< dim, Number, VectorizedArrayType > & | matrix_free, |
const std::pair< unsigned int, unsigned int > & | range, | ||
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 cell 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.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::FEEvaluation | ( | const Mapping< dim > & | mapping, |
const FiniteElement< dim > & | fe, | ||
const Quadrature< 1 > & | quadrature, | ||
const UpdateFlags | update_flags, | ||
const unsigned int | first_selected_component = 0 |
||
) |
Constructor that comes with reduced functionality and works similar as FEValues. The arguments are similar to the ones passed to the constructor of FEValues, with the notable difference that FEEvaluation expects a one-dimensional quadrature formula, Quadrature<1>, instead of a dim
dimensional one. The finite element can be both scalar or vector valued, but this method always only selects a scalar base element at a time (with n_components
copies as specified by the class template). For vector-valued elements, the optional argument first_selected_component
allows to specify the index of the base element to be used for evaluation. Note that the internal data structures always assume that the base element is primitive, non-primitive are not supported currently.
As known from FEValues, a call to the reinit method with a Triangulation<dim>::cell_iterator is necessary to make the geometry and degrees of freedom of the current class known. If the iterator includes DoFHandler information (i.e., it is a DoFHandler<dim>::cell_iterator or similar), the initialization allows to also read from or write to vectors in the standard way for DoFHandler<dim>::active_cell_iterator types for one cell at a time. However, this approach is much slower than the path with MatrixFree with MPI since index translation has to be done. As only one cell at a time is used, this method does not vectorize over several elements (which is most efficient for vector operations), but only possibly within the element if the evaluate/integrate routines are combined inside user code (e.g. for computing cell matrices).
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::FEEvaluation | ( | const FiniteElement< dim > & | fe, |
const Quadrature< 1 > & | quadrature, | ||
const UpdateFlags | update_flags, | ||
const unsigned int | first_selected_component = 0 |
||
) |
Constructor for the reduced functionality. This constructor is equivalent to the other one except that it makes the object use a \(Q_1\) mapping (i.e., an object of type MappingQ(1)) implicitly.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::FEEvaluation | ( | const FiniteElement< dim > & | fe, |
const FEEvaluationData< dim, VectorizedArrayType, false > & | other, | ||
const unsigned int | first_selected_component = 0 |
||
) |
Constructor for the reduced functionality. Similar to the other constructor with FiniteElement argument but using another FEEvaluationBase object to provide information about the geometry. This allows several FEEvaluation objects to share the geometry evaluation, i.e., the underlying mapping and quadrature points do only need to be evaluated once. Make sure to not pass an optional object around when you intend to use the FEEvaluation object in parallel to the given one because otherwise the intended sharing may create race conditions.
FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::FEEvaluation | ( | const FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType > & | other | ) |
Copy constructor. If FEEvaluationBase was constructed from a mapping, fe, quadrature, and update flags, the underlying geometry evaluation based on FEValues will be deep-copied in order to allow for using in parallel with threads.
FEEvaluation & FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::operator= | ( | const FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType > & | other | ) |
Copy assignment operator. If FEEvaluationBase was constructed from a mapping, fe, quadrature, and update flags, the underlying geometry evaluation based on FEValues will be deep-copied in order to allow for using in parallel with threads.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::reinit | ( | const unsigned int | cell_batch_index | ) |
Initialize the operation pointer to the current cell batch index. 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 have to be set appropriately.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::reinit | ( | const std::array< unsigned int, n_lanes > & | cell_ids | ) |
Similar as the above function but allowing to define customized cell batches on the fly. A cell batch is defined by the (matrix-free) index of its cells: see also the documentation of get_cell_ids () or get_cell_or_face_ids ().
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::reinit | ( | const TriaIterator< DoFCellAccessor< dim, dim, level_dof_access > > & | cell | ) |
Initialize the data to the current cell using a TriaIterator object as usual in FEValues. The argument is either of type DoFHandler::active_cell_iterator or DoFHandler::level_cell_iterator. This option is only available if the FEEvaluation object was created with a finite element, quadrature formula and correct update flags and without a MatrixFree object. This initialization method loses the ability to use vectorization, see also the description of the FEEvaluation class. When this reinit method is used, FEEvaluation can also read from vectors (but less efficient than with data coming from MatrixFree).
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::reinit | ( | const typename Triangulation< dim >::cell_iterator & | cell | ) |
Initialize the data to the current cell using a TriaIterator object as usual in FEValues. This option is only available if the FEEvaluation object was created with a finite element, quadrature formula and correct update flags and without a MatrixFree object. This initialization method loses the ability to use vectorization, see also the description of the FEEvaluation class. When this reinit method is used, FEEvaluation can not read from vectors because no DoFHandler information is available.
|
static |
Check if fast evaluation/integration is supported.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::evaluate | ( | const EvaluationFlags::EvaluationFlags | evaluation_flag | ) |
Evaluate the function values, the gradients, and the Hessians of the polynomial interpolation from the DoF values in the input vector to the quadrature points on the unit cell. The function arguments specify which parts shall actually be computed. This function has to be called first so that the access functions get_value(), get_gradient() or get_laplacian() give useful information (unless these values have been set manually).
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::evaluate | ( | const VectorizedArrayType * | values_array, |
const EvaluationFlags::EvaluationFlags | evaluation_flag | ||
) |
Evaluate the function values, the gradients, and the Hessians of the polynomial interpolation from the DoF values in the input array values_array
to 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. This function has to be called first so that the access functions get_value(), get_gradient() or get_laplacian() give useful information (unless these values have been set manually).
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::gather_evaluate | ( | const VectorType & | input_vector, |
const EvaluationFlags::EvaluationFlags | evaluation_flag | ||
) |
Read from the input vector and evaluates the function values, the gradients, and the Hessians of the polynomial interpolation of the vector entries from input_vector
associated with the current cell to the quadrature points on the unit cell. The function arguments specify which parts shall actually be computed. This function has to be called first so that the access functions get_value(), get_gradient() or get_laplacian() give useful information (unless these values have been set manually).
This call is equivalent to calling read_dof_values() followed by evaluate(), but might internally use some additional optimizations.
void FEEvaluation< dim, fe_degree, n_q_points_1d, n_components_, Number, VectorizedArrayType >::integrate | ( | const EvaluationFlags::EvaluationFlags | integration_flag | ) |
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 FEEvaluation< 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
, whose previous results is overwritten, rather than writing it on the internal data structures behind begin_dof_values().
void FEEvaluation< 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, performs the cell integration, and adds the result into the global vector output_vector
on the degrees of freedom associated with the present cell index. 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.
std_cxx20::ranges::iota_view< unsigned int, unsigned int > FEEvaluation< 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.
|
private |
Checks if the template arguments regarding degree of the element corresponds to the actual element used at initialization.
|
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 1420 of file fe_evaluation.h.
|
staticconstexpr |
The number of solution components of the evaluator given as template argument.
Definition at line 1426 of file fe_evaluation.h.
|
staticconstexpr |
The number of lanes of the template argument VectorizedArrayType.
Definition at line 1431 of file fe_evaluation.h.
|
staticconstexpr |
The static number of quadrature points determined from the given template argument n_q_points_1d
.
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 1441 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 or if the underlying is of more complicated type than the usual FE_Q or FE_DGQ ones, such as FE_DGP. Definition at line 1453 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 or if the underlying is of more complicated type than the usual FE_Q or FE_DGQ ones, such as FE_DGP. Definition at line 1465 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 or if the underlying is of more complicated type than the usual FE_Q or FE_DGQ ones, such as FE_DGP. Definition at line 1477 of file fe_evaluation.h.
const unsigned int FEEvaluation< 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 1772 of file fe_evaluation.h.
const unsigned int FEEvaluation< 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 1780 of file fe_evaluation.h.
const unsigned int FEEvaluation< 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
-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 1789 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.