Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | Public Attributes | List of all members
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > Class Template Reference

#include <deal.II/fe/fe_update_flags.h>

Inheritance diagram for internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >:
[legend]

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
 

Detailed Description

template<int dim, int spacedim = dim>
class internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >

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 403 of file fe_update_flags.h.

Member Function Documentation

◆ initialize()

template<int dim, int spacedim>
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 2704 of file fe_values.cc.

◆ memory_consumption()

template<int dim, int spacedim>
std::size_t internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::memory_consumption

Compute and return an estimate for the memory consumption (in bytes) of this object.

Definition at line 2765 of file fe_values.cc.

Member Data Documentation

◆ JxW_values

template<int dim, int spacedim = dim>
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 433 of file fe_update_flags.h.

◆ jacobians

template<int dim, int spacedim = dim>
std::vector<DerivativeForm<1, dim, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::jacobians

Array of the Jacobian matrices at the quadrature points.

Definition at line 438 of file fe_update_flags.h.

◆ jacobian_grads

template<int dim, int spacedim = dim>
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 444 of file fe_update_flags.h.

◆ inverse_jacobians

template<int dim, int spacedim = dim>
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 449 of file fe_update_flags.h.

◆ jacobian_pushed_forward_grads

template<int dim, int spacedim = dim>
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 455 of file fe_update_flags.h.

◆ jacobian_2nd_derivatives

template<int dim, int spacedim = dim>
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 461 of file fe_update_flags.h.

◆ jacobian_pushed_forward_2nd_derivatives

template<int dim, int spacedim = dim>
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 467 of file fe_update_flags.h.

◆ jacobian_3rd_derivatives

template<int dim, int spacedim = dim>
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 473 of file fe_update_flags.h.

◆ jacobian_pushed_forward_3rd_derivatives

template<int dim, int spacedim = dim>
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 479 of file fe_update_flags.h.

◆ quadrature_points

template<int dim, int spacedim = dim>
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 486 of file fe_update_flags.h.

◆ normal_vectors

template<int dim, int spacedim = dim>
std::vector<Tensor<1, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::normal_vectors

List of outward normal vectors at the quadrature points.

Definition at line 491 of file fe_update_flags.h.

◆ boundary_forms

template<int dim, int spacedim = dim>
std::vector<Tensor<1, spacedim> > internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >::boundary_forms

List of boundary forms at the quadrature points.

Definition at line 496 of file fe_update_flags.h.


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