Reference documentation for deal.II version 9.0.0
Classes | Public Member Functions | Friends | List of all members
Mapping< dim, spacedim > Class Template Referenceabstract

Abstract base class for mapping classes. More...

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

Inheritance diagram for Mapping< dim, spacedim >:
[legend]

Classes

class  InternalDataBase
 

Public Member Functions

virtual ~Mapping () override=default
 
virtual std::unique_ptr< Mapping< dim, spacedim > > clone () const =0
 
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
 
virtual bool preserves_vertex_locations () const =0
 
Mapping points between reference and real cells
virtual Point< spacedim > transform_unit_to_real_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const =0
 
virtual Point< dim > transform_real_to_unit_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
 
Point< dim-1 > project_real_point_to_unit_point_on_face (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int &face_no, const Point< spacedim > &p) const
 
Functions to transform tensors from reference to real coordinates
virtual void transform (const ArrayView< const Tensor< 1, dim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const =0
 
virtual void transform (const ArrayView< const DerivativeForm< 1, dim, spacedim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 2, spacedim > > &output) const =0
 
virtual void transform (const ArrayView< const Tensor< 2, dim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 2, spacedim > > &output) const =0
 
virtual void transform (const ArrayView< const DerivativeForm< 2, dim, spacedim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 3, spacedim > > &output) const =0
 
virtual void transform (const ArrayView< const Tensor< 3, dim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 3, spacedim > > &output) const =0
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (const char *identifier=nullptr) const
 
void unsubscribe (const char *identifier=nullptr) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Static Public Member Functions

Exceptions
static ::ExceptionBaseExcInvalidData ()
 
static ::ExceptionBaseExcTransformationFailed ()
 
static ::ExceptionBaseExcDistortedMappedCell (Point< spacedim > arg1, double arg2, int arg3)
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Protected Member Functions

Interface with FEValues
virtual UpdateFlags requires_update_flags (const UpdateFlags update_flags) const =0
 
virtual std::unique_ptr< InternalDataBaseget_data (const UpdateFlags update_flags, const Quadrature< dim > &quadrature) const =0
 
virtual std::unique_ptr< InternalDataBaseget_face_data (const UpdateFlags update_flags, const Quadrature< dim-1 > &quadrature) const =0
 
virtual std::unique_ptr< InternalDataBaseget_subface_data (const UpdateFlags update_flags, const Quadrature< dim-1 > &quadrature) const =0
 
virtual CellSimilarity::Similarity fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const =0
 
virtual void fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const =0
 
virtual void fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const =0
 

Friends

class FEValuesBase< dim, spacedim >
 
class FEValues< dim, spacedim >
 
class FEFaceValues< dim, spacedim >
 
class FESubfaceValues< dim, spacedim >
 

Detailed Description

template<int dim, int spacedim = dim>
class Mapping< dim, spacedim >

Abstract base class for mapping classes.

This class declares the interface for the functionality to describe mappings from the reference (unit) cell to a cell in real space, as well as for filling the information necessary to use the FEValues, FEFaceValues, and FESubfaceValues classes. Concrete implementations of these interfaces are provided in derived classes.

Mathematics of the mapping

The mapping is a transformation \(\mathbf x = \mathbf F_K(\hat{\mathbf x})\) which maps points \(\hat{\mathbf x}\) in the reference cell \([0,1]^\text{dim}\) to points \(\mathbf x\) in the actual grid cell \(K\subset{\mathbb R}^\text{spacedim}\). Many of the applications of such mappings require the Jacobian of this mapping, \(J(\hat{\mathbf x}) = \hat\nabla {\mathbf F}_K(\hat{\mathbf x})\). For instance, if dim=spacedim=2, we have

\[ J(\hat{\mathbf x}) = \left(\begin{matrix} \frac{\partial x}{\partial \hat x} & \frac{\partial x}{\partial \hat y} \\ \frac{\partial y}{\partial \hat x} & \frac{\partial y}{\partial \hat y} \end{matrix}\right) \]

Mapping of scalar functions

The shape functions of scalar finite elements are typically defined on a reference cell and are then simply mapped according to the rule

\[ \varphi(\mathbf x) = \varphi\bigl(\mathbf F_K(\hat{\mathbf x})\bigr) = \hat \varphi(\hat{\mathbf x}). \]

Mapping of integrals

Using simply a change of variables, integrals of scalar functions over a cell \(K\) can be expressed as an integral over the reference cell \(\hat K\). Specifically, The volume form \(d\hat x\) is transformed so that

\[ \int_K u(\mathbf x)\,dx = \int_{\hat K} \hat u(\hat{\mathbf x}) \left|\text{det}J(\hat{\mathbf x})\right| \,d\hat x. \]

In expressions where such integrals are approximated by quadrature, this then leads to terms of the form

\[ \int_K u(\mathbf x)\,dx \approx \sum_{q} \hat u(\hat{\mathbf x}_q) \underbrace{\left|\text{det}J(\hat{\mathbf x}_q)\right| w_q}_{=: \text{JxW}_q}. \]

Here, the weights \(\text{JxW}_q\) of each quadrature point (where JxW mnemonically stands for Jacobian times Quadrature Weights) take the role of the \(dx\) in the original integral. Consequently, they appear in all code that computes integrals approximated by quadrature, and are accessed by FEValues::JxW().

Todo:
Document what happens in the codimension-1 case.

Mapping of vector fields, differential forms and gradients of vector fields

The transformation of vector fields or differential forms (gradients of scalar functions) \(\mathbf v\), and gradients of vector fields \(\mathbf T\) follows the general form

\[ \mathbf v(\mathbf x) = \mathbf A(\hat{\mathbf x}) \hat{\mathbf v}(\hat{\mathbf x}), \qquad \mathbf T(\mathbf x) = \mathbf A(\hat{\mathbf x}) \hat{\mathbf T}(\hat{\mathbf x}) \mathbf B(\hat{\mathbf x}). \]

The differential forms A and B are determined by the kind of object being transformed. These transformations are performed through the transform() functions, and the type of object being transformed is specified by their MappingType argument. See the documentation there for possible choices.

Derivatives of the mapping

Some applications require the derivatives of the mapping, of which the first order derivative is the mapping Jacobian, \(J_{iJ}(\hat{\mathbf x})=\frac{\partial x_i}{\partial \hat x_J}\), described above. Higher order derivatives of the mapping are similarly defined, for example the Jacobian derivative, \(\hat H_{iJK}(\hat{\mathbf x}) = \frac{\partial^2 x_i}{\partial \hat x_J \partial \hat x_K}\), and the Jacobian second derivative, \(\hat K_{iJKL}(\hat{\mathbf x}) = \frac{\partial^3 x_i}{\partial \hat x_J \partial \hat x_K \partial \hat x_L}\). It is also useful to define the "pushed-forward" versions of the higher order derivatives: the Jacobian pushed-forward derivative, \(H_{ijk}(\hat{\mathbf x}) = \frac{\partial^2 x_i}{\partial \hat x_J \partial \hat x_K}(J_{jJ})^{-1}(J_{kK})^{-1}\), and the Jacobian pushed-forward second derivative, \(K_{ijkl}(\hat{\mathbf x}) = \frac{\partial^3 x_i}{\partial \hat x_J \partial \hat x_K \partial \hat x_L}(J_{jJ})^{-1}(J_{kK})^{-1}(J_{lL})^{-1}\). These pushed-forward versions can be used to compute the higher order derivatives of functions defined on the reference cell with respect to the real cell coordinates. For instance, the Jacobian derivative with respect to the real cell coordinates is given by:

\[ \frac{\partial}{\partial x_j}\left[J_{iJ}(\hat{\mathbf x})\right] = H_{ikn}(\hat{\mathbf x})J_{nJ}(\hat{\mathbf x}), \]

and the derivative of the Jacobian inverse with respect to the real cell coordinates is similarly given by:

\[ \frac{\partial}{\partial x_j}\left[\left(J_{iJ}(\hat{\mathbf x})\right)^{-1}\right] = -H_{nik}(\hat{\mathbf x})\left(J_{nJ}(\hat{\mathbf x})\right)^{-1}. \]

In a similar fashion, higher order derivatives, with respect to the real cell coordinates, of functions defined on the reference cell can be defined using the Jacobian pushed-forward higher-order derivatives. For example, the derivative, with respect to the real cell coordinates, of the Jacobian pushed-forward derivative is given by:

\[ \frac{\partial}{\partial x_l}\left[H_{ijk}(\hat{\mathbf x})\right] = K_{ijkl}(\hat{\mathbf x}) -H_{mjl}(\hat{\mathbf x})H_{imk}(\hat{\mathbf x})-H_{mkl}(\hat{\mathbf x})H_{imj}(\hat{\mathbf x}). \]

References

A general publication on differential geometry and finite elements is the survey

The description of the Piola transform has been taken from the lecture notes by Ronald H. W. Hoppe, University of Houston, Chapter 7.

Author
Guido Kanschat, Ralf Hartmann 2000, 2001

Definition at line 46 of file dof_tools.h.

Constructor & Destructor Documentation

◆ ~Mapping()

template<int dim, int spacedim = dim>
virtual Mapping< dim, spacedim >::~Mapping ( )
overridevirtualdefault

Virtual destructor.

Member Function Documentation

◆ clone()

template<int dim, int spacedim = dim>
virtual std::unique_ptr<Mapping<dim,spacedim> > Mapping< dim, spacedim >::clone ( ) const
pure virtual

Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.

The function is declared abstract virtual in this base class, and derived classes will have to implement it.

This function is mainly used by the hp::MappingCollection class.

Implemented in MappingQ< dim, spacedim >, MappingQGeneric< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ1Eulerian< dim, VectorType, spacedim >, MappingManifold< dim, spacedim >, MappingCartesian< dim, spacedim >, MappingQ1< dim, spacedim >, and MappingC1< dim, spacedim >.

◆ get_vertices()

template<int dim, int spacedim>
std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > Mapping< dim, spacedim >::get_vertices ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell) const
virtual

Return the mapped vertices of a cell.

Most of the time, these values will simply be the coordinates of the vertices of a cell as returned by cell->vertex(v) for vertex v, i.e., information stored by the triangulation. However, there are also mappings that add displacements or choose completely different locations, e.g., MappingQEulerian, MappingQ1Eulerian, or MappingFEField.

The default implementation of this function simply returns the information stored by the triangulation, i.e., cell->vertex(v).

Reimplemented in MappingQEulerian< dim, VectorType, spacedim >::MappingQEulerianGeneric, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQEulerian< dim, VectorType, spacedim >, and MappingQ1Eulerian< dim, VectorType, spacedim >.

Definition at line 25 of file mapping.cc.

◆ preserves_vertex_locations()

template<int dim, int spacedim = dim>
virtual bool Mapping< dim, spacedim >::preserves_vertex_locations ( ) const
pure virtual

Return whether the mapping preserves vertex locations. In other words, this function returns whether the mapped location of the reference cell vertices (given by GeometryInfo::unit_cell_vertex()) equals the result of cell->vertex() (i.e., information stored by the triangulation).

For example, implementations in derived classes return true for MappingQ, MappingQGeneric, MappingCartesian, but false for MappingQEulerian, MappingQ1Eulerian, and MappingFEField.

Implemented in MappingQEulerian< dim, VectorType, spacedim >::MappingQEulerianGeneric, MappingQGeneric< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ1Eulerian< dim, VectorType, spacedim >, MappingQ< dim, spacedim >, MappingManifold< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ transform_unit_to_real_cell()

template<int dim, int spacedim = dim>
virtual Point<spacedim> Mapping< dim, spacedim >::transform_unit_to_real_cell ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const Point< dim > &  p 
) const
pure virtual

Map the point p on the unit cell to the corresponding point on the real cell cell.

Parameters
cellIterator to the cell that will be used to define the mapping.
pLocation of a point on the reference cell.
Returns
The location of the reference point mapped to real space using the mapping defined by the class derived from the current one that implements the mapping, and the coordinates of the cell identified by the first argument.

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ< dim, spacedim >, MappingManifold< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ transform_real_to_unit_cell()

template<int dim, int spacedim = dim>
virtual Point<dim> Mapping< dim, spacedim >::transform_real_to_unit_cell ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const Point< spacedim > &  p 
) const
pure virtual

Map the point p on the real cell to the corresponding point on the unit cell, and return its coordinates. This function provides the inverse of the mapping provided by transform_unit_to_real_cell().

In the codimension one case, this function returns the normal projection of the real point p on the curve or surface identified by the cell.

Note
Polynomial mappings from the reference (unit) cell coordinates to the coordinate system of a real cell are not always invertible if the point for which the inverse mapping is to be computed lies outside the cell's boundaries. In such cases, the current function may fail to compute a point on the reference cell whose image under the mapping equals the given point p. If this is the case then this function throws an exception of type Mapping::ExcTransformationFailed . Whether the given point p lies outside the cell can therefore be determined by checking whether the returned reference coordinates lie inside or outside the reference cell (e.g., using GeometryInfo::is_inside_unit_cell()) or whether the exception mentioned above has been thrown.
Parameters
cellIterator to the cell that will be used to define the mapping.
pLocation of a point on the given cell.
Returns
The reference cell location of the point that when mapped to real space equals the coordinates given by the second argument. This mapping uses the mapping defined by the class derived from the current one that implements the mapping, and the coordinates of the cell identified by the first argument.

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ< dim, spacedim >, MappingManifold< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ project_real_point_to_unit_point_on_face()

template<int dim, int spacedim>
Point< dim-1 > Mapping< dim, spacedim >::project_real_point_to_unit_point_on_face ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int &  face_no,
const Point< spacedim > &  p 
) const

Transform the point p on the real cell to the corresponding point on the unit cell, and then projects it to a dim-1 point on the face with the given face number face_no. Ideally the point p is near the face face_no, but any point in the cell can technically be projected.

This function does not make physical sense when dim=1, so it throws an exception in this case.

Definition at line 40 of file mapping.cc.

◆ requires_update_flags()

template<int dim, int spacedim = dim>
virtual UpdateFlags Mapping< dim, spacedim >::requires_update_flags ( const UpdateFlags  update_flags) const
protectedpure virtual

Given a set of update flags, compute which other quantities also need to be computed in order to satisfy the request by the given flags. Then return the combination of the original set of flags and those just computed.

As an example, if update_flags contains update_JxW_values (i.e., the product of the determinant of the Jacobian and the weights provided by the quadrature formula), a mapping may require the computation of the full Jacobian matrix in order to compute its determinant. They would then return not just update_JxW_values, but also update_jacobians. (This is not how it is actually done internally in the derived classes that compute the JxW values – they set update_contravariant_transformation instead, from which the determinant can also be computed – but this does not take away from the instructiveness of the example.)

An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module.

See also
UpdateFlags

Implemented in MappingQGeneric< dim, spacedim >, MappingManifold< dim, spacedim >, MappingQ< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, and MappingCartesian< dim, spacedim >.

◆ get_data()

template<int dim, int spacedim = dim>
virtual std::unique_ptr<InternalDataBase> Mapping< dim, spacedim >::get_data ( const UpdateFlags  update_flags,
const Quadrature< dim > &  quadrature 
) const
protectedpure virtual

Create and return a pointer to an object into which mappings can store data that only needs to be computed once but that can then be used whenever the mapping is applied to a concrete cell (e.g., in the various transform() functions, as well as in the fill_fe_values(), fill_fe_face_values() and fill_fe_subface_values() that form the interface of mappings with the FEValues class).

Derived classes will return pointers to objects of a type derived from Mapping::InternalDataBase (see there for more information) and may pre- compute some information already (in accordance with what will be asked of the mapping in the future, as specified by the update flags) and for the given quadrature object. Subsequent calls to transform() or fill_fe_values() and friends will then receive back the object created here (with the same set of update flags and for the same quadrature object). Derived classes can therefore pre-compute some information in their get_data() function and store it in the internal data object.

The mapping classes do not keep track of the objects created by this function. Ownership will therefore rest with the caller.

An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module.

Parameters
update_flagsA set of flags that define what is expected of the mapping class in future calls to transform() or the fill_fe_values() group of functions. This set of flags may contain flags that mappings do not know how to deal with (e.g., for information that is in fact computed by the finite element classes, such as UpdateFlags::update_values). Derived classes will need to store these flags, or at least that subset of flags that will require the mapping to perform any actions in fill_fe_values(), in InternalDataBase::update_each.
quadratureThe quadrature object for which mapping information will have to be computed. This includes the locations and weights of quadrature points.
Returns
A pointer to a newly created object of type InternalDataBase (or a derived class). Ownership of this object passes to the calling function.
Note
C++ allows that virtual functions in derived classes may return pointers to objects not of type InternalDataBase but in fact pointers to objects of classes derived from InternalDataBase. (This feature is called "covariant return types".) This is useful in some contexts where the calling is within the derived class and will immediately make use of the returned object, knowing its real (derived) type.

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingManifold< dim, spacedim >, MappingQ< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ get_face_data()

template<int dim, int spacedim = dim>
virtual std::unique_ptr<InternalDataBase> Mapping< dim, spacedim >::get_face_data ( const UpdateFlags  update_flags,
const Quadrature< dim-1 > &  quadrature 
) const
protectedpure virtual

Like get_data(), but in preparation for later calls to transform() or fill_fe_face_values() that will need information about mappings from the reference face to a face of a concrete cell.

Parameters
update_flagsA set of flags that define what is expected of the mapping class in future calls to transform() or the fill_fe_values() group of functions. This set of flags may contain flags that mappings do not know how to deal with (e.g., for information that is in fact computed by the finite element classes, such as UpdateFlags::update_values). Derived classes will need to store these flags, or at least that subset of flags that will require the mapping to perform any actions in fill_fe_values(), in InternalDataBase::update_each.
quadratureThe quadrature object for which mapping information will have to be computed. This includes the locations and weights of quadrature points.
Returns
A pointer to a newly created object of type InternalDataBase (or a derived class). Ownership of this object passes to the calling function.
Note
C++ allows that virtual functions in derived classes may return pointers to objects not of type InternalDataBase but in fact pointers to objects of classes derived from InternalDataBase. (This feature is called "covariant return types".) This is useful in some contexts where the calling is within the derived class and will immediately make use of the returned object, knowing its real (derived) type.

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingManifold< dim, spacedim >, MappingQ< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ get_subface_data()

template<int dim, int spacedim = dim>
virtual std::unique_ptr<InternalDataBase> Mapping< dim, spacedim >::get_subface_data ( const UpdateFlags  update_flags,
const Quadrature< dim-1 > &  quadrature 
) const
protectedpure virtual

Like get_data() and get_face_data(), but in preparation for later calls to transform() or fill_fe_subface_values() that will need information about mappings from the reference face to a child of a face (i.e., subface) of a concrete cell.

Parameters
update_flagsA set of flags that define what is expected of the mapping class in future calls to transform() or the fill_fe_values() group of functions. This set of flags may contain flags that mappings do not know how to deal with (e.g., for information that is in fact computed by the finite element classes, such as UpdateFlags::update_values). Derived classes will need to store these flags, or at least that subset of flags that will require the mapping to perform any actions in fill_fe_values(), in InternalDataBase::update_each.
quadratureThe quadrature object for which mapping information will have to be computed. This includes the locations and weights of quadrature points.
Returns
A pointer to a newly created object of type InternalDataBase (or a derived class). Ownership of this object passes to the calling function.
Note
C++ allows that virtual functions in derived classes may return pointers to objects not of type InternalDataBase but in fact pointers to objects of classes derived from InternalDataBase. (This feature is called "covariant return types".) This is useful in some contexts where the calling is within the derived class and will immediately make use of the returned object, knowing its real (derived) type.

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingManifold< dim, spacedim >, MappingQ< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ fill_fe_values()

template<int dim, int spacedim = dim>
virtual CellSimilarity::Similarity Mapping< dim, spacedim >::fill_fe_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const CellSimilarity::Similarity  cell_similarity,
const Quadrature< dim > &  quadrature,
const typename Mapping< dim, spacedim >::InternalDataBase internal_data,
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  output_data 
) const
protectedpure virtual

Compute information about the mapping from the reference cell to the real cell indicated by the first argument to this function. Derived classes will have to implement this function based on the kind of mapping they represent. It is called by FEValues::reinit().

Conceptually, this function's represents the application of the mapping \(\mathbf x=\mathbf F_K(\hat {\mathbf x})\) from reference coordinates \(\mathbf\in [0,1]^d\) to real space coordinates \(\mathbf x\) for a given cell \(K\). Its purpose is to compute the following kinds of data:

  • Data that results from the application of the mapping itself, e.g., computing the location \(\mathbf x_q = \mathbf F_K(\hat{\mathbf x}_q)\) of quadrature points on the real cell, and that is directly useful to users of FEValues, for example during assembly.
  • Data that is necessary for finite element implementations to compute their shape functions on the real cell. To this end, the FEValues::reinit() function calls FiniteElement::fill_fe_values() after the current function, and the output of this function serves as input to FiniteElement::fill_fe_values(). Examples of information that needs to be computed here for use by the finite element classes is the Jacobian of the mapping, \(\hat\nabla \mathbf F_K(\hat{\mathbf x})\) or its inverse, for example to transform the gradients of shape functions on the reference cell to the gradients of shape functions on the real cell.

The information computed by this function is used to fill the various member variables of the output argument of this function. Which of the member variables of that structure should be filled is determined by the update flags stored in the Mapping::InternalDataBase object passed to this function.

An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module.

Parameters
[in]cellThe cell of the triangulation for which this function is to compute a mapping from the reference cell to.
[in]cell_similarityWhether or not the cell given as first argument is simply a translation, rotation, etc of the cell for which this function was called the most recent time. This information is computed simply by matching the vertices (as stored by the Triangulation) between the previous and the current cell. The value passed here may be modified by implementations of this function and should then be returned (see the discussion of the return value of this function).
[in]quadratureA reference to the quadrature formula in use for the current evaluation. This quadrature object is the same as the one used when creating the internal_data object. The object is used both to map the location of quadrature points, as well as to compute the JxW values for each quadrature point (which involves the quadrature weights).
[in]internal_dataA reference to an object previously created by get_data() and that may be used to store information the mapping can compute once on the reference cell. See the documentation of the Mapping::InternalDataBase class for an extensive description of the purpose of these objects.
[out]output_dataA reference to an object whose member variables should be computed. Not all of the members of this argument need to be filled; which ones need to be filled is determined by the update flags stored inside the internal_data object.
Returns
An updated value of the cell_similarity argument to this function. The returned value will be used for the corresponding argument when FEValues::reinit() calls FiniteElement::fill_fe_values(). In most cases, derived classes will simply want to return the value passed for cell_similarity. However, implementations of this function may downgrade the level of cell similarity. This is, for example, the case for classes that take not only into account the locations of the vertices of a cell (as reported by the Triangulation), but also other information specific to the mapping. The purpose is that FEValues::reinit() can compute whether a cell is similar to the previous one only based on the cell's vertices, whereas the mapping may also consider displacement fields (e.g., in the MappingQ1Eulerian and MappingFEField classes). In such cases, the mapping may conclude that the previously computed cell similarity is too optimistic, and invalidate it for subsequent use in FiniteElement::fill_fe_values() by returning a less optimistic cell similarity value.
Note
FEValues ensures that this function is always called with the same pair of internal_data and output_data objects. In other words, if an implementation of this function knows that it has written a piece of data into the output argument in a previous call, then there is no need to copy it there again in a later call if the implementation knows that this is the same value.

Implemented in MappingQGeneric< dim, spacedim >, and MappingManifold< dim, spacedim >.

◆ fill_fe_face_values()

template<int dim, int spacedim = dim>
virtual void Mapping< dim, spacedim >::fill_fe_face_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const Quadrature< dim-1 > &  quadrature,
const typename Mapping< dim, spacedim >::InternalDataBase internal_data,
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  output_data 
) const
protectedpure virtual

This function is the equivalent to Mapping::fill_fe_values(), but for faces of cells. See there for an extensive discussion of its purpose. It is called by FEFaceValues::reinit().

Parameters
[in]cellThe cell of the triangulation for which this function is to compute a mapping from the reference cell to.
[in]face_noThe number of the face of the given cell for which information is requested.
[in]quadratureA reference to the quadrature formula in use for the current evaluation. This quadrature object is the same as the one used when creating the internal_data object. The object is used both to map the location of quadrature points, as well as to compute the JxW values for each quadrature point (which involves the quadrature weights).
[in]internal_dataA reference to an object previously created by get_data() and that may be used to store information the mapping can compute once on the reference cell. See the documentation of the Mapping::InternalDataBase class for an extensive description of the purpose of these objects.
[out]output_dataA reference to an object whose member variables should be computed. Not all of the members of this argument need to be filled; which ones need to be filled is determined by the update flags stored inside the internal_data object.

Implemented in MappingQGeneric< dim, spacedim >, and MappingManifold< dim, spacedim >.

◆ fill_fe_subface_values()

template<int dim, int spacedim = dim>
virtual void Mapping< dim, spacedim >::fill_fe_subface_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  subface_no,
const Quadrature< dim-1 > &  quadrature,
const typename Mapping< dim, spacedim >::InternalDataBase internal_data,
::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  output_data 
) const
protectedpure virtual

This function is the equivalent to Mapping::fill_fe_values(), but for subfaces (i.e., children of faces) of cells. See there for an extensive discussion of its purpose. It is called by FESubfaceValues::reinit().

Parameters
[in]cellThe cell of the triangulation for which this function is to compute a mapping from the reference cell to.
[in]face_noThe number of the face of the given cell for which information is requested.
[in]subface_noThe number of the child of a face of the given cell for which information is requested.
[in]quadratureA reference to the quadrature formula in use for the current evaluation. This quadrature object is the same as the one used when creating the internal_data object. The object is used both to map the location of quadrature points, as well as to compute the JxW values for each quadrature point (which involves the quadrature weights).
[in]internal_dataA reference to an object previously created by get_data() and that may be used to store information the mapping can compute once on the reference cell. See the documentation of the Mapping::InternalDataBase class for an extensive description of the purpose of these objects.
[out]output_dataA reference to an object whose member variables should be computed. Not all of the members of this argument need to be filled; which ones need to be filled is determined by the update flags stored inside the internal_data object.

Implemented in MappingQGeneric< dim, spacedim >, and MappingManifold< dim, spacedim >.

◆ transform() [1/5]

template<int dim, int spacedim = dim>
virtual void Mapping< dim, spacedim >::transform ( const ArrayView< const Tensor< 1, dim > > &  input,
const MappingType  type,
const typename Mapping< dim, spacedim >::InternalDataBase internal,
const ArrayView< Tensor< 1, spacedim > > &  output 
) const
pure virtual

Transform a field of vectors or 1-differential forms according to the selected MappingType.

Note
Normally, this function is called by a finite element, filling FEValues objects. For this finite element, there should be an alias MappingType like mapping_bdm, mapping_nedelec, etc. This alias should be preferred to using the types below.

The mapping types currently implemented by derived classes are:

  • mapping_contravariant: maps a vector field on the reference cell to the physical cell through the Jacobian:

    \[ \mathbf u(\mathbf x) = J(\hat{\mathbf x})\hat{\mathbf u}(\hat{\mathbf x}). \]

    In physics, this is usually referred to as the contravariant transformation. Mathematically, it is the push forward of a vector field.

  • mapping_covariant: maps a field of one-forms on the reference cell to a field of one-forms on the physical cell. (Theoretically this would refer to a DerivativeForm<1,dim,1> but we canonically identify this type with a Tensor<1,dim>). Mathematically, it is the pull back of the differential form

    \[ \mathbf u(\mathbf x) = J(\hat{\mathbf x})(J(\hat{\mathbf x})^{T} J(\hat{\mathbf x}))^{-1}\hat{\mathbf u}(\hat{\mathbf x}). \]

    Gradients of scalar differentiable functions are transformed this way.

    In the case when dim=spacedim the previous formula reduces to

    \[ \mathbf u(\mathbf x) = J(\hat{\mathbf x})^{-T}\hat{\mathbf u}(\hat{\mathbf x}) \]

    because we assume that the mapping \(\mathbf F_K\) is always invertible, and consequently its Jacobian \(J\) is an invertible matrix.

  • mapping_piola: A field of dim-1-forms on the reference cell is also represented by a vector field, but again transforms differently, namely by the Piola transform

    \[ \mathbf u(\mathbf x) = \frac{1}{\text{det}\;J(\hat{\mathbf x})} J(\hat{\mathbf x}) \hat{\mathbf u}(\hat{\mathbf x}). \]

Parameters
[in]inputAn array (or part of an array) of input objects that should be mapped.
[in]typeThe kind of mapping to be applied.
[in]internalA pointer to an object of type Mapping::InternalDataBase that contains information previously stored by the mapping. The object pointed to was created by the get_data(), get_face_data(), or get_subface_data() function, and will have been updated as part of a call to fill_fe_values(), fill_fe_face_values(), or fill_fe_subface_values() for the current cell, before calling the current function. In other words, this object also represents with respect to which cell the transformation should be applied to.
[out]outputAn array (or part of an array) into which the transformed objects should be placed. (Note that the array view is const, but the tensors it points to are not.)

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ< dim, spacedim >, MappingManifold< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ transform() [2/5]

template<int dim, int spacedim = dim>
virtual void Mapping< dim, spacedim >::transform ( const ArrayView< const DerivativeForm< 1, dim, spacedim > > &  input,
const MappingType  type,
const typename Mapping< dim, spacedim >::InternalDataBase internal,
const ArrayView< Tensor< 2, spacedim > > &  output 
) const
pure virtual

Transform a field of differential forms from the reference cell to the physical cell. It is useful to think of \(\mathbf{T} = \nabla \mathbf u\) and \(\hat{\mathbf T} = \hat \nabla \hat{\mathbf u}\), with \(\mathbf u\) a vector field. The mapping types currently implemented by derived classes are:

  • mapping_covariant: maps a field of forms on the reference cell to a field of forms on the physical cell. Mathematically, it is the pull back of the differential form

    \[ \mathbf T(\mathbf x) = \hat{\mathbf T}(\hat{\mathbf x}) J(\hat{\mathbf x})(J(\hat{\mathbf x})^{T} J(\hat{\mathbf x}))^{-1}. \]

    Jacobians of spacedim-vector valued differentiable functions are transformed this way.

    In the case when dim=spacedim the previous formula reduces to

    \[ \mathbf T(\mathbf x) = \hat{\mathbf u}(\hat{\mathbf x}) J(\hat{\mathbf x})^{-1}. \]

Note
It would have been more reasonable to make this transform a template function with the rank in DerivativeForm<1, dim, rank>. Unfortunately C++ does not allow templatized virtual functions. This is why we identify DerivativeForm<1, dim, 1> with a Tensor<1,dim> when using mapping_covariant() in the function transform() above this one.
Parameters
[in]inputAn array (or part of an array) of input objects that should be mapped.
[in]typeThe kind of mapping to be applied.
[in]internalA pointer to an object of type Mapping::InternalDataBase that contains information previously stored by the mapping. The object pointed to was created by the get_data(), get_face_data(), or get_subface_data() function, and will have been updated as part of a call to fill_fe_values(), fill_fe_face_values(), or fill_fe_subface_values() for the current cell, before calling the current function. In other words, this object also represents with respect to which cell the transformation should be applied to.
[out]outputAn array (or part of an array) into which the transformed objects should be placed. (Note that the array view is const, but the tensors it points to are not.)

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ< dim, spacedim >, MappingManifold< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ transform() [3/5]

template<int dim, int spacedim = dim>
virtual void Mapping< dim, spacedim >::transform ( const ArrayView< const Tensor< 2, dim > > &  input,
const MappingType  type,
const typename Mapping< dim, spacedim >::InternalDataBase internal,
const ArrayView< Tensor< 2, spacedim > > &  output 
) const
pure virtual

Transform a tensor field from the reference cell to the physical cell. These tensors are usually the Jacobians in the reference cell of vector fields that have been pulled back from the physical cell. The mapping types currently implemented by derived classes are:

  • mapping_contravariant_gradient: it assumes \(\mathbf u(\mathbf x) = J \hat{\mathbf u}\) so that

    \[ \mathbf T(\mathbf x) = J(\hat{\mathbf x}) \hat{\mathbf T}(\hat{\mathbf x}) J(\hat{\mathbf x})^{-1}. \]

  • mapping_covariant_gradient: it assumes \(\mathbf u(\mathbf x) = J^{-T} \hat{\mathbf u}\) so that

    \[ \mathbf T(\mathbf x) = J(\hat{\mathbf x})^{-T} \hat{\mathbf T}(\hat{\mathbf x}) J(\hat{\mathbf x})^{-1}. \]

  • mapping_piola_gradient: it assumes \(\mathbf u(\mathbf x) = \frac{1}{\text{det}\;J(\hat{\mathbf x})} J(\hat{\mathbf x}) \hat{\mathbf u}(\hat{\mathbf x})\) so that

    \[ \mathbf T(\mathbf x) = \frac{1}{\text{det}\;J(\hat{\mathbf x})} J(\hat{\mathbf x}) \hat{\mathbf T}(\hat{\mathbf x}) J(\hat{\mathbf x})^{-1}. \]

Todo:
The formulas for mapping_covariant_gradient, mapping_contravariant_gradient and mapping_piola_gradient are only true as stated for linear mappings. If, for example, the mapping is bilinear (or has a higher order polynomial degree) then there is a missing term associated with the derivative of \(J\).
Parameters
[in]inputAn array (or part of an array) of input objects that should be mapped.
[in]typeThe kind of mapping to be applied.
[in]internalA pointer to an object of type Mapping::InternalDataBase that contains information previously stored by the mapping. The object pointed to was created by the get_data(), get_face_data(), or get_subface_data() function, and will have been updated as part of a call to fill_fe_values(), fill_fe_face_values(), or fill_fe_subface_values() for the current cell, before calling the current function. In other words, this object also represents with respect to which cell the transformation should be applied to.
[out]outputAn array (or part of an array) into which the transformed objects should be placed. (Note that the array view is const, but the tensors it points to are not.)

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ< dim, spacedim >, MappingManifold< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ transform() [4/5]

template<int dim, int spacedim = dim>
virtual void Mapping< dim, spacedim >::transform ( const ArrayView< const DerivativeForm< 2, dim, spacedim > > &  input,
const MappingType  type,
const typename Mapping< dim, spacedim >::InternalDataBase internal,
const ArrayView< Tensor< 3, spacedim > > &  output 
) const
pure virtual

Transform a tensor field from the reference cell to the physical cell. This tensors are most of times the hessians in the reference cell of vector fields that have been pulled back from the physical cell.

The mapping types currently implemented by derived classes are:

  • mapping_covariant_gradient: maps a field of forms on the reference cell to a field of forms on the physical cell. Mathematically, it is the pull back of the differential form

    \[ \mathbf T_{ijk}(\mathbf x) = \hat{\mathbf T}_{iJK}(\hat{\mathbf x}) J_{jJ}^{\dagger} J_{kK}^{\dagger}\]

    ,

    where

    \[ J^{\dagger} = J(\hat{\mathbf x})(J(\hat{\mathbf x})^{T} J(\hat{\mathbf x}))^{-1}. \]

Hessians of spacedim-vector valued differentiable functions are transformed this way (After subtraction of the product of the derivative with the Jacobian gradient).

In the case when dim=spacedim the previous formula reduces to

\[J^{\dagger} = J^{-1}\]

Parameters
[in]inputAn array (or part of an array) of input objects that should be mapped.
[in]typeThe kind of mapping to be applied.
[in]internalA pointer to an object of type Mapping::InternalDataBase that contains information previously stored by the mapping. The object pointed to was created by the get_data(), get_face_data(), or get_subface_data() function, and will have been updated as part of a call to fill_fe_values(), fill_fe_face_values(), or fill_fe_subface_values() for the current cell, before calling the current function. In other words, this object also represents with respect to which cell the transformation should be applied to.
[out]outputAn array (or part of an array) into which the transformed objects should be placed. (Note that the array view is const, but the tensors it points to are not.)

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ< dim, spacedim >, MappingManifold< dim, spacedim >, and MappingCartesian< dim, spacedim >.

◆ transform() [5/5]

template<int dim, int spacedim = dim>
virtual void Mapping< dim, spacedim >::transform ( const ArrayView< const Tensor< 3, dim > > &  input,
const MappingType  type,
const typename Mapping< dim, spacedim >::InternalDataBase internal,
const ArrayView< Tensor< 3, spacedim > > &  output 
) const
pure virtual

Transform a field of 3-differential forms from the reference cell to the physical cell. It is useful to think of \(\mathbf{T}_{ijk} = D^2_{jk} \mathbf u_i\) and \(\mathbf{\hat T}_{IJK} = \hat D^2_{JK} \mathbf{\hat u}_I\), with \(\mathbf u_i\) a vector field.

The mapping types currently implemented by derived classes are:

  • mapping_contravariant_hessian: it assumes \(\mathbf u_i(\mathbf x) = J_{iI} \hat{\mathbf u}_I\) so that

    \[ \mathbf T_{ijk}(\mathbf x) = J_{iI}(\hat{\mathbf x}) \hat{\mathbf T}_{IJK}(\hat{\mathbf x}) J_{jJ}(\hat{\mathbf x})^{-1} J_{kK}(\hat{\mathbf x})^{-1}. \]

  • mapping_covariant_hessian: it assumes \(\mathbf u_i(\mathbf x) = J_{iI}^{-T} \hat{\mathbf u}_I\) so that

    \[ \mathbf T_{ijk}(\mathbf x) = J_iI(\hat{\mathbf x})^{-1} \hat{\mathbf T}_{IJK}(\hat{\mathbf x}) J_{jJ}(\hat{\mathbf x})^{-1} J_{kK}(\hat{\mathbf x})^{-1}. \]

  • mapping_piola_hessian: it assumes \(\mathbf u_i(\mathbf x) = \frac{1}{\text{det}\;J(\hat{\mathbf x})} J_{iI}(\hat{\mathbf x}) \hat{\mathbf u}(\hat{\mathbf x})\) so that

    \[ \mathbf T_{ijk}(\mathbf x) = \frac{1}{\text{det}\;J(\hat{\mathbf x})} J_{iI}(\hat{\mathbf x}) \hat{\mathbf T}_{IJK}(\hat{\mathbf x}) J_{jJ}(\hat{\mathbf x})^{-1} J_{kK}(\hat{\mathbf x})^{-1}. \]

Parameters
[in]inputAn array (or part of an array) of input objects that should be mapped.
[in]typeThe kind of mapping to be applied.
[in]internalA pointer to an object of type Mapping::InternalDataBase that contains information previously stored by the mapping. The object pointed to was created by the get_data(), get_face_data(), or get_subface_data() function, and will have been updated as part of a call to fill_fe_values(), fill_fe_face_values(), or fill_fe_subface_values() for the current cell, before calling the current function. In other words, this object also represents with respect to which cell the transformation should be applied to.
[out]outputAn array (or part of an array) into which the transformed objects should be placed.

Implemented in MappingQGeneric< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >, MappingQ< dim, spacedim >, MappingManifold< dim, spacedim >, and MappingCartesian< dim, spacedim >.

Friends And Related Function Documentation

◆ FEValuesBase< dim, spacedim >

template<int dim, int spacedim = dim>
friend class FEValuesBase< dim, spacedim >
friend

Give class FEValues access to the private get_...data and fill_fe_...values functions.

Definition at line 1176 of file mapping.h.


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