deal.II version GIT relicensing-2289-g1e5549a87a 2024-12-21 21:30:00+00:00
|
#include <deal.II/fe/mapping_related_data.h>
Public Member Functions | |
void | initialize (const unsigned int n_quadrature_points, const UpdateFlags flags) |
std::size_t | memory_consumption () const |
Public Attributes | |
std::vector< double > | JxW_values |
std::vector< DerivativeForm< 1, dim, spacedim > > | jacobians |
std::vector< DerivativeForm< 2, dim, spacedim > > | jacobian_grads |
std::vector< DerivativeForm< 1, spacedim, dim > > | inverse_jacobians |
std::vector< Tensor< 3, spacedim > > | jacobian_pushed_forward_grads |
std::vector< DerivativeForm< 3, dim, spacedim > > | jacobian_2nd_derivatives |
std::vector< Tensor< 4, spacedim > > | jacobian_pushed_forward_2nd_derivatives |
std::vector< DerivativeForm< 4, dim, spacedim > > | jacobian_3rd_derivatives |
std::vector< Tensor< 5, spacedim > > | jacobian_pushed_forward_3rd_derivatives |
std::vector< Point< spacedim > > | quadrature_points |
std::vector< Tensor< 1, spacedim > > | normal_vectors |
std::vector< Tensor< 1, spacedim > > | boundary_forms |
A class that stores all of the mapping related data used in FEValues, FEFaceValues, and FESubfaceValues objects. Objects of this kind will be given as output argument when FEValues::reinit() calls Mapping::fill_fe_values() for a given cell, face, or subface.
The data herein will then be provided as input argument in the following call to FiniteElement::fill_fe_values().
Definition at line 51 of file mapping_related_data.h.
void internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::initialize | ( | const unsigned int | n_quadrature_points, |
const UpdateFlags | flags | ||
) |
Initialize all vectors to correct size.
Definition at line 29 of file mapping_related_data.cc.
std::size_t internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::memory_consumption | ( | ) | const |
Compute and return an estimate for the memory consumption (in bytes) of this object.
Definition at line 90 of file mapping_related_data.cc.
std::vector<double> internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::JxW_values |
Store an array of weights times the Jacobi determinant at the quadrature points. This function is reset each time reinit() is called. The Jacobi determinant is actually the reciprocal value of the Jacobi matrices stored in this class, see the general documentation of this class for more information.
However, if this object refers to an FEFaceValues or FESubfaceValues object, then the JxW_values correspond to the Jacobian of the transformation of the face, not the cell, i.e. the dimensionality is that of a surface measure, not of a volume measure. In this case, it is computed from the boundary forms, rather than the Jacobian matrix.
Definition at line 81 of file mapping_related_data.h.
std::vector<DerivativeForm<1, dim, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::jacobians |
Array of the Jacobian matrices at the quadrature points.
Definition at line 86 of file mapping_related_data.h.
std::vector<DerivativeForm<2, dim, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::jacobian_grads |
Array of the derivatives of the Jacobian matrices at the quadrature points.
Definition at line 92 of file mapping_related_data.h.
std::vector<DerivativeForm<1, spacedim, dim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::inverse_jacobians |
Array of the inverse Jacobian matrices at the quadrature points.
Definition at line 97 of file mapping_related_data.h.
std::vector<Tensor<3, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::jacobian_pushed_forward_grads |
Array of the derivatives of the Jacobian matrices at the quadrature points, pushed forward to the real cell coordinates.
Definition at line 103 of file mapping_related_data.h.
std::vector<DerivativeForm<3, dim, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::jacobian_2nd_derivatives |
Array of the second derivatives of the Jacobian matrices at the quadrature points.
Definition at line 109 of file mapping_related_data.h.
std::vector<Tensor<4, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::jacobian_pushed_forward_2nd_derivatives |
Array of the second derivatives of the Jacobian matrices at the quadrature points, pushed forward to the real cell coordinates.
Definition at line 115 of file mapping_related_data.h.
std::vector<DerivativeForm<4, dim, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::jacobian_3rd_derivatives |
Array of the third derivatives of the Jacobian matrices at the quadrature points.
Definition at line 121 of file mapping_related_data.h.
std::vector<Tensor<5, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::jacobian_pushed_forward_3rd_derivatives |
Array of the third derivatives of the Jacobian matrices at the quadrature points, pushed forward to the real cell coordinates.
Definition at line 127 of file mapping_related_data.h.
std::vector<Point<spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::quadrature_points |
Array of quadrature points. This array is set up upon calling reinit() and contains the quadrature points on the real element, rather than on the reference element.
Definition at line 134 of file mapping_related_data.h.
std::vector<Tensor<1, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::normal_vectors |
List of outward normal vectors at the quadrature points.
Definition at line 139 of file mapping_related_data.h.
std::vector<Tensor<1, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::boundary_forms |
List of boundary forms at the quadrature points.
Definition at line 144 of file mapping_related_data.h.