Reference documentation for deal.II version 9.0.0
Classes | Public Member Functions | Private Member Functions | Static Private Attributes | List of all members
MappingCartesian< dim, spacedim > Class Template Reference

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

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

Classes

class  InternalData
 

Public Member Functions

virtual std::unique_ptr< Mapping< dim, spacedim > > clone () const override
 
virtual bool preserves_vertex_locations () const override
 
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 override
 
virtual Point< dim > transform_real_to_unit_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
 
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 override
 
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 override
 
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 override
 
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 override
 
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 override
 
- Public Member Functions inherited from Mapping< dim, spacedim >
virtual ~Mapping () override=default
 
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
 
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
 
- 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)
 

Private Member Functions

void compute_fill (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const CellSimilarity::Similarity cell_similarity, const InternalData &data, std::vector< Point< dim > > &quadrature_points, std::vector< Tensor< 1, dim > > &normal_vectors) const
 
Interface with FEValues
virtual UpdateFlags requires_update_flags (const UpdateFlags update_flags) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBaseget_data (const UpdateFlags, const Quadrature< dim > &quadrature) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBaseget_face_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBaseget_subface_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
 
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 override
 
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 override
 
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 override
 

Static Private Attributes

static const unsigned int invalid_face_number = numbers::invalid_unsigned_int
 

Additional Inherited Members

- Static Public Member Functions inherited from Mapping< dim, spacedim >
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 inherited from Mapping< dim, spacedim >
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
 

Detailed Description

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

A class providing a mapping from the reference cell to cells that are axiparallel.

This class maps the unit cell to a grid cell with surfaces parallel to the coordinate lines/planes. It is specifically developed for Cartesian meshes. In other words, the mapping is meant for cells for which the mapping from the reference to the real cell is a scaling along the coordinate directions: The transformation from reference coordinates \(\hat {\mathbf x}\) to real coordinates \(\mathbf x\) on each cell is of the form

\begin{align*} {\mathbf x}(\hat {\mathbf x}) = \begin{pmatrix} h_x & 0 \\ 0 & h_y \end{pmatrix} \hat{\mathbf x} + {\mathbf v}_0 \end{align*}

in 2d, and

\begin{align*} {\mathbf x}(\hat {\mathbf x}) = \begin{pmatrix} h_x & 0 & 0 \\ 0 & h_y & 0 \\ 0 & 0 & h_z \end{pmatrix} \hat{\mathbf x} + {\mathbf v}_0 \end{align*}

in 3d, where \({\mathbf v}_0\) is the bottom left vertex and \(h_x,h_y,h_z\) are the extents of the cell along the axes.

The class is intended for efficiency, and it does not do a whole lot of error checking. If you apply this mapping to a cell that does not conform to the requirements above, you will get strange results.

Author
Guido Kanschat, 2001; Ralf Hartmann, 2005

Definition at line 61 of file mapping_cartesian.h.

Member Function Documentation

◆ clone()

template<int dim, int spacedim>
std::unique_ptr< Mapping< dim, spacedim > > MappingCartesian< dim, spacedim >::clone ( ) const
overridevirtual

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.

Implements Mapping< dim, spacedim >.

Definition at line 1066 of file mapping_cartesian.cc.

◆ preserves_vertex_locations()

template<int dim, int spacedim>
bool MappingCartesian< dim, spacedim >::preserves_vertex_locations ( ) const
overridevirtual

Return true because MappingCartesian preserves vertex locations.

Implements Mapping< dim, spacedim >.

Definition at line 66 of file mapping_cartesian.cc.

◆ transform_unit_to_real_cell()

template<int dim, int spacedim>
Point< spacedim > MappingCartesian< dim, spacedim >::transform_unit_to_real_cell ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const Point< dim > &  p 
) const
overridevirtual

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.

Implements Mapping< dim, spacedim >.

Definition at line 997 of file mapping_cartesian.cc.

◆ transform_real_to_unit_cell()

template<int dim, int spacedim>
Point< dim > MappingCartesian< dim, spacedim >::transform_real_to_unit_cell ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const Point< spacedim > &  p 
) const
overridevirtual

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.

Implements Mapping< dim, spacedim >.

Definition at line 1032 of file mapping_cartesian.cc.

◆ transform() [1/5]

template<int dim, int spacedim>
void MappingCartesian< 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
overridevirtual

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.)

Implements Mapping< dim, spacedim >.

Definition at line 618 of file mapping_cartesian.cc.

◆ transform() [2/5]

template<int dim, int spacedim>
void MappingCartesian< 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
overridevirtual

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.)

Implements Mapping< dim, spacedim >.

Definition at line 673 of file mapping_cartesian.cc.

◆ transform() [3/5]

template<int dim, int spacedim>
void MappingCartesian< 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
overridevirtual

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.)

Implements Mapping< dim, spacedim >.

Definition at line 774 of file mapping_cartesian.cc.

◆ transform() [4/5]

template<int dim, int spacedim>
void MappingCartesian< 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
overridevirtual

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.)

Implements Mapping< dim, spacedim >.

Definition at line 875 of file mapping_cartesian.cc.

◆ transform() [5/5]

template<int dim, int spacedim>
void MappingCartesian< 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
overridevirtual

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.

Implements Mapping< dim, spacedim >.

Definition at line 913 of file mapping_cartesian.cc.

◆ requires_update_flags()

template<int dim, int spacedim>
UpdateFlags MappingCartesian< dim, spacedim >::requires_update_flags ( const UpdateFlags  update_flags) const
overrideprivatevirtual

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

Implements Mapping< dim, spacedim >.

Definition at line 75 of file mapping_cartesian.cc.

◆ get_data()

template<int dim, int spacedim>
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > MappingCartesian< dim, spacedim >::get_data ( const UpdateFlags  update_flags,
const Quadrature< dim > &  quadrature 
) const
overrideprivatevirtual

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.

Implements Mapping< dim, spacedim >.

Definition at line 93 of file mapping_cartesian.cc.

◆ get_face_data()

template<int dim, int spacedim>
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > MappingCartesian< dim, spacedim >::get_face_data ( const UpdateFlags  update_flags,
const Quadrature< dim-1 > &  quadrature 
) const
overrideprivatevirtual

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.

Implements Mapping< dim, spacedim >.

Definition at line 110 of file mapping_cartesian.cc.

◆ get_subface_data()

template<int dim, int spacedim>
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > MappingCartesian< dim, spacedim >::get_subface_data ( const UpdateFlags  update_flags,
const Quadrature< dim-1 > &  quadrature 
) const
overrideprivatevirtual

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.

Implements Mapping< dim, spacedim >.

Definition at line 131 of file mapping_cartesian.cc.

◆ compute_fill()

template<int dim, int spacedim>
void MappingCartesian< dim, spacedim >::compute_fill ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const unsigned int  sub_no,
const CellSimilarity::Similarity  cell_similarity,
const InternalData data,
std::vector< Point< dim > > &  quadrature_points,
std::vector< Tensor< 1, dim > > &  normal_vectors 
) const
private

Do the computation for the fill_* functions.

Definition at line 153 of file mapping_cartesian.cc.

Member Data Documentation

◆ invalid_face_number

template<int dim, int spacedim = dim>
const unsigned int MappingCartesian< dim, spacedim >::invalid_face_number = numbers::invalid_unsigned_int
staticprivate

Value to indicate that a given face or subface number is invalid.

Definition at line 267 of file mapping_cartesian.h.


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