Reference documentation for deal.II version GIT relicensing-399-g79d89019c5 2024-04-16 15:00:02+00:00
\(\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
MappingManifold< dim, spacedim >::InternalData Class Reference

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

Inheritance diagram for MappingManifold< dim, spacedim >::InternalData:
Inheritance graph
[legend]

Public Member Functions

 InternalData ()=default
 
virtual void reinit (const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
 
void initialize_face (const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
 
void compute_manifold_quadrature_weights (const Quadrature< dim > &quadrature)
 
void store_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
 
virtual std::size_t memory_consumption () const override
 

Public Attributes

std::vector< Point< spacedim > > vertices
 
Triangulation< dim, spacedim >::cell_iterator cell
 
Quadrature< dim > quad
 
std::vector< std::vector< double > > cell_manifold_quadrature_weights
 
std::vector< double > vertex_weights
 
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
 
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
 
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
 
std::vector< std::vector< Tensor< 1, spacedim > > > aux
 
std::vector< double > volume_elements
 
SmartPointer< const Manifold< dim, spacedim > > manifold
 
UpdateFlags update_each
 

Detailed Description

template<int dim, int spacedim = dim>
class MappingManifold< dim, spacedim >::InternalData

Storage for internal data of polynomial mappings. See Mapping::InternalDataBase for an extensive description.

For the current class, the InternalData class stores data that is computed once when the object is created (in get_data()) as well as data the class wants to store from between the call to fill_fe_values(), fill_fe_face_values(), or fill_fe_subface_values() until possible later calls from the finite element to functions such as transform(). The latter class of member variables are marked as 'mutable'.

Definition at line 162 of file mapping_manifold.h.

Constructor & Destructor Documentation

◆ InternalData()

template<int dim, int spacedim = dim>
MappingManifold< dim, spacedim >::InternalData::InternalData ( )
default

Constructor.

Member Function Documentation

◆ reinit()

template<int dim, int spacedim>
void MappingManifold< dim, spacedim >::InternalData::reinit ( const UpdateFlags  update_flags,
const Quadrature< dim > &  quadrature 
)
overridevirtual

Definition at line 69 of file mapping_manifold.cc.

◆ initialize_face()

template<int dim, int spacedim>
void MappingManifold< dim, spacedim >::InternalData::initialize_face ( const UpdateFlags  update_flags,
const Quadrature< dim > &  quadrature,
const unsigned int  n_original_q_points 
)

Initialize the object's member variables related to cell and face data based on the given arguments. In order to initialize cell data, this function calls initialize().

Definition at line 105 of file mapping_manifold.cc.

◆ compute_manifold_quadrature_weights()

template<int dim, int spacedim>
void MappingManifold< dim, spacedim >::InternalData::compute_manifold_quadrature_weights ( const Quadrature< dim > &  quadrature)

Compute the weights associated to the Manifold object, that need to be passed when computing the location of the quadrature points.

Definition at line 170 of file mapping_manifold.cc.

◆ store_vertices()

template<int dim, int spacedim>
void MappingManifold< dim, spacedim >::InternalData::store_vertices ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell) const

Store vertices internally.

Definition at line 157 of file mapping_manifold.cc.

◆ memory_consumption()

template<int dim, int spacedim>
std::size_t MappingManifold< dim, spacedim >::InternalData::memory_consumption ( ) const
overridevirtual

Return an estimate (in bytes) for the memory consumption of this object.

Definition at line 48 of file mapping_manifold.cc.

Member Data Documentation

◆ vertices

template<int dim, int spacedim = dim>
std::vector<Point<spacedim> > MappingManifold< dim, spacedim >::InternalData::vertices
mutable

The current cell vertices.

Computed each.

Definition at line 211 of file mapping_manifold.h.

◆ cell

template<int dim, int spacedim = dim>
Triangulation<dim,spacedim>::cell_iterator MappingManifold< dim, spacedim >::InternalData::cell
mutable

The current cell.

Computed each.

Definition at line 218 of file mapping_manifold.h.

◆ quad

template<int dim, int spacedim = dim>
Quadrature<dim> MappingManifold< dim, spacedim >::InternalData::quad

The actual quadrature on the reference cell.

Computed once.

Definition at line 225 of file mapping_manifold.h.

◆ cell_manifold_quadrature_weights

template<int dim, int spacedim = dim>
std::vector<std::vector<double> > MappingManifold< dim, spacedim >::InternalData::cell_manifold_quadrature_weights

Values of quadrature weights for manifold quadrature formulas.

The Manifold class has a function (Manifold::get_new_point()) that returns new points according to a weighted average of some surrounding points on the Manifold. For each quadrature point, we call this function with a Quadrature formula constructed using the vertices of the current cell, and the values of the basis functions of an FE_Q(1) finite element evaluated at the quadrature point itself. While the vertices of the cell change for every cell, the weights can be computed once for each quadrature point. We store this information in the following variable, where the first index runs through the quadrature points, and the second index runs through the vertex indices.

Computed once.

Definition at line 245 of file mapping_manifold.h.

◆ vertex_weights

template<int dim, int spacedim = dim>
std::vector<double> MappingManifold< dim, spacedim >::InternalData::vertex_weights
mutable

A vector of weights for use in Manifold::get_new_point(). For each point (interior to a cell), we compute the weight each vertex has for this point. If the point lies at a vertex, then this vertex has weight one and all others have weight zero. If the point lies interior to a cell, then the weight every vertex has is just the \(d\)-linear shape functions associated with each vertex evaluated at that point.

This array has size GeometryInfo<dim>::vertices_per_cell, but it can't be converted into a fixed size array because it is used as input for Manifold::get_new_point() which wants to see a std::vector<double> for the weights.

Definition at line 261 of file mapping_manifold.h.

◆ unit_tangentials

template<int dim, int spacedim = dim>
std::array<std::vector<Tensor<1, dim> >, GeometryInfo<dim>::faces_per_cell *(dim - 1)> MappingManifold< dim, spacedim >::InternalData::unit_tangentials

Unit tangential vectors. Used for the computation of boundary forms and normal vectors.

This array has (dim-1) * GeometryInfo::faces_per_cell entries. The first GeometryInfo::faces_per_cell contain the vectors in the first tangential direction for each face; the second set of GeometryInfo<dim>::faces_per_cell entries contain the vectors in the second tangential direction (only in 3d, since there we have 2 tangential directions per face), etc.

Filled once.

Definition at line 278 of file mapping_manifold.h.

◆ covariant

template<int dim, int spacedim = dim>
std::vector<DerivativeForm<1, dim, spacedim> > MappingManifold< dim, spacedim >::InternalData::covariant
mutable

Tensors of covariant transformation at each of the quadrature points. The matrix stored is the Jacobian * G^{-1}, where G = Jacobian^{t} * Jacobian, is the first fundamental form of the map; if dim=spacedim then it reduces to the transpose of the inverse of the Jacobian matrix, which itself is stored in the contravariant field of this structure.

Computed on each cell.

Definition at line 289 of file mapping_manifold.h.

◆ contravariant

template<int dim, int spacedim = dim>
std::vector<DerivativeForm<1, dim, spacedim> > MappingManifold< dim, spacedim >::InternalData::contravariant
mutable

Tensors of contravariant transformation at each of the quadrature points. The contravariant matrix is the Jacobian of the transformation, i.e. \(J_{ij}=dx_i/d\hat x_j\).

Computed on each cell.

Definition at line 298 of file mapping_manifold.h.

◆ aux

template<int dim, int spacedim = dim>
std::vector<std::vector<Tensor<1, spacedim> > > MappingManifold< dim, spacedim >::InternalData::aux
mutable

Auxiliary vectors for internal use.

Definition at line 303 of file mapping_manifold.h.

◆ volume_elements

template<int dim, int spacedim = dim>
std::vector<double> MappingManifold< dim, spacedim >::InternalData::volume_elements
mutable

The determinant of the Jacobian in each quadrature point. Filled if update_volume_elements.

Definition at line 309 of file mapping_manifold.h.

◆ manifold

template<int dim, int spacedim = dim>
SmartPointer<const Manifold<dim, spacedim> > MappingManifold< dim, spacedim >::InternalData::manifold
mutable

A pointer to the Manifold in use.

Updated each.

Definition at line 316 of file mapping_manifold.h.

◆ update_each

UpdateFlags Mapping< dim, spacedim >::InternalDataBase::update_each
inherited

A set of update flags specifying the kind of information that an implementation of the Mapping interface needs to compute on each cell or face, i.e., in Mapping::fill_fe_values() and friends.

This set of flags is stored here by implementations of Mapping::get_data(), Mapping::get_face_data(), or Mapping::get_subface_data(), and is that subset of the update flags passed to those functions that require re-computation on every cell. (The subset of the flags corresponding to information that can be computed once and for all already at the time of the call to Mapping::get_data() – or an implementation of that interface – need not be stored here because it has already been taken care of.)

Definition at line 691 of file mapping.h.


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