Reference documentation for deal.II version 9.5.0
\(\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
Classes | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
MappingQEulerian< dim, VectorType, spacedim > Class Template Reference

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

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

Classes

class  SupportQuadrature
 

Public Member Functions

 MappingQEulerian (const unsigned int degree, const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector, const unsigned int level=numbers::invalid_unsigned_int)
 
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > get_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
 
virtual std::unique_ptr< Mapping< dim, spacedim > > clone () const override
 
virtual bool preserves_vertex_locations () const override
 
virtual std::vector< Point< spacedim > > compute_mapping_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
 
unsigned int get_degree () const
 
virtual BoundingBox< spacedim > get_bounding_box (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override
 
virtual bool is_compatible_with (const ReferenceCell &reference_cell) const override
 
void fill_mapping_data_for_generic_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points, const UpdateFlags update_flags, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
 
void fill_mapping_data_for_face_quadrature (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_number, const Quadrature< dim - 1 > &face_quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
 
boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_face > get_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no) const
 
virtual Point< spacedim > get_center (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_barycenter_of_reference_cell=true) const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
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
 
virtual void transform_points_real_to_unit_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const override
 
Functions to transform tensors from reference to real coordinates
virtual void transform (const ArrayView< const Tensor< 1, dim > > &input, const MappingKind kind, 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 MappingKind kind, 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 MappingKind kind, 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 MappingKind kind, 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 MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 3, spacedim > > &output) const override
 
Mapping points between reference and real cells
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
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcInactiveCell ()
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
Exceptions
static ::ExceptionBaseExcInvalidData ()
 
static ::ExceptionBaseExcTransformationFailed ()
 
static ::ExceptionBaseExcDistortedMappedCell (Point< spacedim > arg1, double arg2, int arg3)
 

Protected Member Functions

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
 
Point< dim > transform_real_to_unit_cell_internal (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit) const
 
Point< 1 > transform_real_to_unit_cell_internal (const Triangulation< 1, 1 >::cell_iterator &cell, const Point< 1 > &p, const Point< 1 > &initial_p_unit) const
 
Point< 2 > transform_real_to_unit_cell_internal (const Triangulation< 2, 2 >::cell_iterator &cell, const Point< 2 > &p, const Point< 2 > &initial_p_unit) const
 
Point< 3 > transform_real_to_unit_cell_internal (const Triangulation< 3, 3 >::cell_iterator &cell, const Point< 3 > &p, const Point< 3 > &initial_p_unit) const
 
Point< 1 > transform_real_to_unit_cell_internal (const Triangulation< 1, 2 >::cell_iterator &cell, const Point< 2 > &p, const Point< 1 > &initial_p_unit) const
 
Point< 2 > transform_real_to_unit_cell_internal (const Triangulation< 2, 3 >::cell_iterator &cell, const Point< 3 > &p, const Point< 2 > &initial_p_unit) const
 
Point< 1 > transform_real_to_unit_cell_internal (const Triangulation< 1, 3 >::cell_iterator &, const Point< 3 > &, const Point< 1 > &) const
 
virtual void add_line_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
 
virtual void add_quad_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const
 
void add_quad_support_points (const Triangulation< 3, 3 >::cell_iterator &cell, std::vector< Point< 3 > > &a) const
 
void add_quad_support_points (const Triangulation< 2, 3 >::cell_iterator &cell, std::vector< Point< 3 > > &a) const
 
Interface with FEValues and friends
virtual UpdateFlags requires_update_flags (const UpdateFlags update_flags) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data (const UpdateFlags, const Quadrature< dim > &quadrature) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data (const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data (const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
 
virtual void fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< 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
 
virtual void fill_fe_immersed_surface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const NonMatching::ImmersedSurfaceQuadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
 
Interface with FEValues
virtual std::unique_ptr< InternalDataBaseget_face_data (const UpdateFlags update_flags, const Quadrature< dim - 1 > &quadrature) const
 
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
 

Protected Attributes

SmartPointer< const VectorType, MappingQEulerian< dim, VectorType, spacedim > > euler_vector
 
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim > > euler_dof_handler
 
const unsigned int polynomial_degree
 
const std::vector< Point< 1 > > line_support_points
 
const std::vector< Polynomials::Polynomial< double > > polynomials_1d
 
const std::vector< unsigned intrenumber_lexicographic_to_hierarchic
 
const std::vector< Point< dim > > unit_cell_support_points
 
const std::vector< Table< 2, double > > support_point_weights_perimeter_to_interior
 
const Table< 2, double > support_point_weights_cell
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

void check_no_subscribers () const noexcept
 

Private Attributes

const unsigned int level
 
const SupportQuadrature support_quadrature
 
FEValues< dim, spacedim > fe_values
 
Threads::Mutex fe_values_mutex
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Detailed Description

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
class MappingQEulerian< dim, VectorType, spacedim >

This class is an extension of the MappingQ1Eulerian class to higher order \(Q_p\) mappings. It is useful when one wants to calculate shape function information on a domain that is deforming as the computation proceeds.

Usage

The constructor of this class takes three arguments: the polynomial degree of the desired Qp mapping, a reference to the vector that defines the mapping from the initial configuration to the current configuration, and a reference to the DoFHandler. The most common case is to use the solution vector for the problem under consideration as the shift vector. The key requirement is that the number of components of the given vector field must be equal to (or possibly greater than) the number of space dimensions. If there are more components than space dimensions (for example, if one is working with a coupled problem where there are additional solution variables), the first dim components are assumed to represent the displacement field, and the remaining components are ignored. If this assumption does not hold one may need to set up a separate DoFHandler on the triangulation and associate the desired shift vector to it.

Typically, the DoFHandler operates on a finite element that is constructed as a system element (FESystem) from continuous FE_Q objects. An example is shown below:

FESystem<dim> fe(FE_Q<dim>(2), dim, FE_Q<dim>(1), 1);
dof_handler.distribute_dofs(fe);
Vector<double> displacement_field(dof_handler.n_dofs());
// ... compute displacement field somehow...
MappingQEulerian<dim> q2_mapping(2, dof_handler, displacement_field);
Definition fe_q.h:551
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation

In this example, our element consists of (dim+1) components. Only the first dim components will be used, however, to define the Q2 mapping. The remaining components are ignored.

Note that it is essential to call the distribute_dofs(...) function before constructing a mapping object.

Also note that since the vector of shift values and the dof handler are only associated to this object at construction time, you have to make sure that whenever you use this object, the given objects still represent valid data.

To enable the use of the MappingQEulerian class also in the context of parallel codes using the PETSc or Trilinos wrapper classes, the type of the vector can be specified as template parameter VectorType.

Definition at line 95 of file mapping_q_eulerian.h.

Member Typedef Documentation

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 230 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 235 of file subscriptor.h.

Constructor & Destructor Documentation

◆ MappingQEulerian()

template<int dim, class VectorType , int spacedim>
MappingQEulerian< dim, VectorType, spacedim >::MappingQEulerian ( const unsigned int  degree,
const DoFHandler< dim, spacedim > &  euler_dof_handler,
const VectorType &  euler_vector,
const unsigned int  level = numbers::invalid_unsigned_int 
)

Constructor.

Parameters
[in]degreeThe polynomial degree of the desired \(Q_p\) mapping.
[in]euler_dof_handlerA DoFHandler object that defines a finite element space. This space needs to have at least dim components and the first dim components of the space will be considered displacements relative to the original positions of the cells of the triangulation.
[in]euler_vectorA finite element function in the space defined by the second argument. The first dim components of this function will be interpreted as the displacement we use in defining the mapping, relative to the location of cells of the underlying triangulation.
[in]levelThe multi-grid level at which the mapping will be used. It is mainly used to check if the size of the euler_vector is consistent with the euler_dof_handler .

Definition at line 48 of file mapping_q_eulerian.cc.

Member Function Documentation

◆ get_vertices() [1/2]

template<int dim, class VectorType , int spacedim>
boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > MappingQEulerian< dim, VectorType, spacedim >::get_vertices ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell) const
overridevirtual

Return the mapped vertices of the cell. For the current class, this function does not use the support points from the geometry of the current cell but instead evaluates an externally given displacement field in addition to the geometry of the cell.

Reimplemented from Mapping< dim, spacedim >.

Definition at line 105 of file mapping_q_eulerian.cc.

◆ clone()

template<int dim, class VectorType , int spacedim>
std::unique_ptr< Mapping< dim, spacedim > > MappingQEulerian< dim, VectorType, spacedim >::clone
overridevirtual

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

Reimplemented from MappingQ< dim, spacedim >.

Definition at line 67 of file mapping_q_eulerian.cc.

◆ preserves_vertex_locations()

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
virtual bool MappingQEulerian< dim, VectorType, spacedim >::preserves_vertex_locations ( ) const
overridevirtual

Always return false because MappingQEulerian does not in general preserve vertex locations (unless the translation vector happens to provide zero displacements at vertex locations).

Reimplemented from MappingQ< dim, spacedim >.

◆ compute_mapping_support_points()

template<int dim, class VectorType , int spacedim>
std::vector< Point< spacedim > > MappingQEulerian< dim, VectorType, spacedim >::compute_mapping_support_points ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell) const
overridevirtual

Return the locations of support points for the mapping. For example, for \(Q_1\) mappings these are the vertices, and for higher order polynomial mappings they are the vertices plus interior points on edges, faces, and the cell interior that are placed in consultation with the Manifold description of the domain and its boundary. However, other classes may override this function differently. In particular, the MappingQ1Eulerian class does exactly this by not computing the support points from the geometry of the current cell but instead evaluating an externally given displacement field in addition to the geometry of the cell.

The default implementation of this function is appropriate for most cases. It takes the locations of support points on the boundary of the cell from the underlying manifold. Interior support points (ie. support points in quads for 2d, in hexes for 3d) are then computed using an interpolation from the lower-dimensional entities (lines, quads) in order to make the transformation as smooth as possible without introducing additional boundary layers within the cells due to the placement of support points.

The function works its way from the vertices (which it takes from the given cell) via the support points on the line (for which it calls the add_line_support_points() function) and the support points on the quad faces (in 3d, for which it calls the add_quad_support_points() function). It then adds interior support points that are either computed by interpolation from the surrounding points using weights for transfinite interpolation, or if dim<spacedim, it asks the underlying manifold for the locations of interior points.

Reimplemented from MappingQ< dim, spacedim >.

Definition at line 123 of file mapping_q_eulerian.cc.

◆ fill_fe_values()

template<int dim, class VectorType , int spacedim>
CellSimilarity::Similarity MappingQEulerian< dim, VectorType, 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
overrideprotectedvirtual

Compute mapping-related information for a cell. See the documentation of Mapping::fill_fe_values() for a discussion of purpose, arguments, and return value of this function.

This function overrides the function in the base class since we cannot use any cell similarity for this class.

Reimplemented from MappingQ< dim, spacedim >.

Definition at line 196 of file mapping_q_eulerian.cc.

◆ get_degree()

template<int dim, int spacedim>
unsigned int MappingQ< dim, spacedim >::get_degree
inherited

Return the degree of the mapping, i.e. the value which was passed to the constructor.

Definition at line 281 of file mapping_q.cc.

◆ get_bounding_box()

template<int dim, int spacedim>
BoundingBox< spacedim > MappingQ< dim, spacedim >::get_bounding_box ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell) const
overridevirtualinherited

Return the bounding box of a mapped cell.

If you are using a (bi-,tri-)linear mapping that preserves vertex locations, this function simply returns the value also produced by cell->bounding_box(). However, there are also mappings that add displacements or choose completely different locations, e.g., MappingQEulerian, MappingQ1Eulerian, or MappingFEField.

For linear mappings, this function returns the bounding box containing all the vertices of the cell, as returned by the get_vertices() method. For higher order mappings defined through support points, the bounding box is only guaranteed to contain all the support points, and it is, in general, only an approximation of the true bounding box, which may be larger.

Parameters
[in]cellThe cell for which you want to compute the bounding box

Reimplemented from Mapping< dim, spacedim >.

Definition at line 1731 of file mapping_q.cc.

◆ is_compatible_with()

template<int dim, int spacedim>
bool MappingQ< dim, spacedim >::is_compatible_with ( const ReferenceCell reference_cell) const
overridevirtualinherited

Returns if this instance of Mapping is compatible with the type of cell in reference_cell.

Implements Mapping< dim, spacedim >.

Definition at line 1741 of file mapping_q.cc.

◆ transform_unit_to_real_cell()

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

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 290 of file mapping_q.cc.

◆ transform_real_to_unit_cell()

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

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 473 of file mapping_q.cc.

◆ transform_points_real_to_unit_cell()

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::transform_points_real_to_unit_cell ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const ArrayView< const Point< spacedim > > &  real_points,
const ArrayView< Point< dim > > &  unit_points 
) const
overridevirtualinherited

Map multiple points from the real point locations to points in reference locations. The functionality is essentially the same as looping over all points and calling the Mapping::transform_real_to_unit_cell() function for each point individually, but it can be much faster for certain mappings that implement a more specialized version such as MappingQ. The only difference in behavior is that this function will never throw an ExcTransformationFailed() exception. If the transformation fails for real_points[i], the returned unit_points[i] contains std::numeric_limits<double>::infinity() as the first entry.

Reimplemented from Mapping< dim, spacedim >.

Definition at line 593 of file mapping_q.cc.

◆ transform() [1/5]

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::transform ( const ArrayView< const Tensor< 1, dim > > &  input,
const MappingKind  kind,
const typename Mapping< dim, spacedim >::InternalDataBase &  internal,
const ArrayView< Tensor< 1, spacedim > > &  output 
) const
overridevirtualinherited

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

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

The mapping kinds 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]kindThe 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 1314 of file mapping_q.cc.

◆ transform() [2/5]

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::transform ( const ArrayView< const DerivativeForm< 1, dim, spacedim > > &  input,
const MappingKind  kind,
const typename Mapping< dim, spacedim >::InternalDataBase &  internal,
const ArrayView< Tensor< 2, spacedim > > &  output 
) const
overridevirtualinherited

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 kinds 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]kindThe 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 1330 of file mapping_q.cc.

◆ transform() [3/5]

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::transform ( const ArrayView< const Tensor< 2, dim > > &  input,
const MappingKind  kind,
const typename Mapping< dim, spacedim >::InternalDataBase &  internal,
const ArrayView< Tensor< 2, spacedim > > &  output 
) const
overridevirtualinherited

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 kinds 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]kindThe 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 1346 of file mapping_q.cc.

◆ transform() [4/5]

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::transform ( const ArrayView< const DerivativeForm< 2, dim, spacedim > > &  input,
const MappingKind  kind,
const typename Mapping< dim, spacedim >::InternalDataBase &  internal,
const ArrayView< Tensor< 3, spacedim > > &  output 
) const
overridevirtualinherited

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 kinds 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]kindThe 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 1378 of file mapping_q.cc.

◆ transform() [5/5]

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::transform ( const ArrayView< const Tensor< 3, dim > > &  input,
const MappingKind  kind,
const typename Mapping< dim, spacedim >::InternalDataBase &  internal,
const ArrayView< Tensor< 3, spacedim > > &  output 
) const
overridevirtualinherited

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 kinds 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]kindThe 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 1430 of file mapping_q.cc.

◆ fill_mapping_data_for_generic_points()

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::fill_mapping_data_for_generic_points ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const ArrayView< const Point< dim > > &  unit_points,
const UpdateFlags  update_flags,
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  output_data 
) const
inherited

As opposed to the other fill_fe_values() and fill_fe_face_values() functions that rely on pre-computed information of InternalDataBase, this function chooses the flexible evaluation path on the cell and points passed in to the current function.

Parameters
[in]cellThe cell where to evaluate the mapping
[in]unit_pointsThe points in reference coordinates where the transformation (Jacobians, positions) should be computed.
[in]update_flagsThe kind of information that should be computed.
[out]output_dataA struct containing the evaluated quantities such as the Jacobian resulting from application of the mapping on the given cell with its underlying manifolds.

Definition at line 1238 of file mapping_q.cc.

◆ fill_mapping_data_for_face_quadrature()

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::fill_mapping_data_for_face_quadrature ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_number,
const Quadrature< dim - 1 > &  face_quadrature,
const typename Mapping< dim, spacedim >::InternalDataBase &  internal_data,
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  output_data 
) const
inherited

As opposed to the fill_fe_face_values() function that relies on pre-computed information of InternalDataBase, this function chooses the flexible evaluation path on the cell and points passed in to the current function.

Parameters
[in]cellThe cell where to evaluate the mapping.
[in]face_numberThe face number where to evaluate the mapping.
[in]face_quadratureThe quadrature points where the transformation (Jacobians, positions) should be computed.
[in]internal_dataA reference to an object previously created 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 struct containing the evaluated quantities such as the Jacobian resulting from application of the mapping on the given cell with its underlying manifolds.

Definition at line 1278 of file mapping_q.cc.

◆ requires_update_flags()

template<int dim, int spacedim>
UpdateFlags MappingQ< dim, spacedim >::requires_update_flags ( const UpdateFlags  update_flags) const
overrideprotectedvirtualinherited

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
The interplay of UpdateFlags, Mapping, and FiniteElement in FEValues

Implements Mapping< dim, spacedim >.

Definition at line 675 of file mapping_q.cc.

◆ get_data()

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

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 731 of file mapping_q.cc.

◆ get_face_data() [1/2]

template<int dim, int spacedim>
std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > MappingQ< dim, spacedim >::get_face_data ( const UpdateFlags  update_flags,
const hp::QCollection< dim - 1 > &  quadrature 
) const
overrideprotectedvirtualinherited

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.

Reimplemented from Mapping< dim, spacedim >.

Definition at line 746 of file mapping_q.cc.

◆ get_face_data() [2/2]

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
protectedvirtualinherited
Deprecated:
Use the version taking a hp::QCollection argument.

◆ get_subface_data()

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

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 767 of file mapping_q.cc.

◆ fill_fe_face_values() [1/2]

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::fill_fe_face_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no,
const hp::QCollection< dim - 1 > &  quadrature,
const typename Mapping< dim, spacedim >::InternalDataBase &  internal_data,
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtualinherited

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.

Reimplemented from Mapping< dim, spacedim >.

Definition at line 1002 of file mapping_q.cc.

◆ fill_fe_face_values() [2/2]

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
protectedvirtualinherited
Deprecated:
Use the version taking a hp::QCollection argument.

◆ fill_fe_subface_values()

template<int dim, int spacedim>
void MappingQ< 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
overrideprotectedvirtualinherited

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.

Implements Mapping< dim, spacedim >.

Definition at line 1054 of file mapping_q.cc.

◆ fill_fe_immersed_surface_values()

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::fill_fe_immersed_surface_values ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const NonMatching::ImmersedSurfaceQuadrature< dim > &  quadrature,
const typename Mapping< dim, spacedim >::InternalDataBase &  internal_data,
internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &  output_data 
) const
overrideprotectedvirtualinherited

The equivalent of Mapping::fill_fe_values(), but for the case that the quadrature is an ImmersedSurfaceQuadrature. See there for a comprehensive description of the input parameters. This function is called by FEImmersedSurfaceValues::reinit().

Reimplemented from Mapping< dim, spacedim >.

Definition at line 1107 of file mapping_q.cc.

◆ transform_real_to_unit_cell_internal() [1/7]

template<int dim, int spacedim>
Point< dim > MappingQ< dim, spacedim >::transform_real_to_unit_cell_internal ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const Point< spacedim > &  p,
const Point< dim > &  initial_p_unit 
) const
protectedinherited

Transform the point p on the real cell to the corresponding point on the unit cell cell by a Newton iteration.

Definition at line 324 of file mapping_q.cc.

◆ transform_real_to_unit_cell_internal() [2/7]

Point< 1 > MappingQ< 1, 1 >::transform_real_to_unit_cell_internal ( const Triangulation< 1, 1 >::cell_iterator &  cell,
const Point< 1 > &  p,
const Point< 1 > &  initial_p_unit 
) const
protectedinherited

Definition at line 338 of file mapping_q.cc.

◆ transform_real_to_unit_cell_internal() [3/7]

Point< 2 > MappingQ< 2, 2 >::transform_real_to_unit_cell_internal ( const Triangulation< 2, 2 >::cell_iterator &  cell,
const Point< 2 > &  p,
const Point< 2 > &  initial_p_unit 
) const
protectedinherited

Definition at line 358 of file mapping_q.cc.

◆ transform_real_to_unit_cell_internal() [4/7]

Point< 3 > MappingQ< 3, 3 >::transform_real_to_unit_cell_internal ( const Triangulation< 3, 3 >::cell_iterator &  cell,
const Point< 3 > &  p,
const Point< 3 > &  initial_p_unit 
) const
protectedinherited

Definition at line 376 of file mapping_q.cc.

◆ transform_real_to_unit_cell_internal() [5/7]

Point< 1 > MappingQ< 1, 2 >::transform_real_to_unit_cell_internal ( const Triangulation< 1, 2 >::cell_iterator &  cell,
const Point< 2 > &  p,
const Point< 1 > &  initial_p_unit 
) const
protectedinherited

Definition at line 394 of file mapping_q.cc.

◆ transform_real_to_unit_cell_internal() [6/7]

Point< 2 > MappingQ< 2, 3 >::transform_real_to_unit_cell_internal ( const Triangulation< 2, 3 >::cell_iterator &  cell,
const Point< 3 > &  p,
const Point< 2 > &  initial_p_unit 
) const
protectedinherited

Definition at line 427 of file mapping_q.cc.

◆ transform_real_to_unit_cell_internal() [7/7]

Point< 1 > MappingQ< 1, 3 >::transform_real_to_unit_cell_internal ( const Triangulation< 1, 3 >::cell_iterator &  ,
const Point< 3 > &  ,
const Point< 1 > &   
) const
protectedinherited

Definition at line 460 of file mapping_q.cc.

◆ add_line_support_points()

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::add_line_support_points ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
std::vector< Point< spacedim > > &  a 
) const
protectedvirtualinherited

Append the support points of all shape functions located on bounding lines of the given cell to the vector a. Points located on the vertices of a line are not included.

This function uses the underlying manifold object of the line (or, if none is set, of the cell) for the location of the requested points. This function is usually called by compute_mapping_support_points() function.

This function is made virtual in order to allow derived classes to choose shape function support points differently than the present class, which chooses the points as interpolation points on the boundary.

Definition at line 1455 of file mapping_q.cc.

◆ add_quad_support_points() [1/3]

template<int dim, int spacedim>
void MappingQ< dim, spacedim >::add_quad_support_points ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
std::vector< Point< spacedim > > &  a 
) const
protectedvirtualinherited

Append the support points of all shape functions located on bounding faces (quads in 3d) of the given cell to the vector a. This function is only defined for dim=3. Points located on the vertices or lines of a quad are not included.

This function uses the underlying manifold object of the quad (or, if none is set, of the cell) for the location of the requested points. This function is usually called by compute_mapping_support_points().

This function is made virtual in order to allow derived classes to choose shape function support points differently than the present class, which chooses the points as interpolation points on the boundary.

Definition at line 1621 of file mapping_q.cc.

◆ add_quad_support_points() [2/3]

void MappingQ< 3, 3 >::add_quad_support_points ( const Triangulation< 3, 3 >::cell_iterator &  cell,
std::vector< Point< 3 > > &  a 
) const
protectedinherited

Definition at line 1522 of file mapping_q.cc.

◆ add_quad_support_points() [3/3]

void MappingQ< 2, 3 >::add_quad_support_points ( const Triangulation< 2, 3 >::cell_iterator &  cell,
std::vector< Point< 3 > > &  a 
) const
protectedinherited

Definition at line 1591 of file mapping_q.cc.

◆ get_vertices() [2/2]

template<int dim, int spacedim = dim>
boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_face > Mapping< dim, spacedim >::get_vertices ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const unsigned int  face_no 
) const
inherited

Return the mapped vertices of a face.

Same as above but working on a given face of a cell.

Parameters
[in]cellThe cell containing the face.
[in]face_noThe number of the face within the cell.

◆ get_center()

template<int dim, int spacedim = dim>
virtual Point< spacedim > Mapping< dim, spacedim >::get_center ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const bool  map_barycenter_of_reference_cell = true 
) const
virtualinherited

Return one of two possible mapped centers of a cell.

If you are using a (bi-,tri-)linear mapping that preserves vertex locations, this function simply returns the value also produced by cell->center(). However, there are also mappings that add displacements or choose completely different locations, e.g., MappingQEulerian, MappingQ1Eulerian, or MappingFEField, and mappings based on high order polynomials, for which the center may not coincide with the average of the vertex locations.

By default, this function returns the push forward of the barycenter of the reference cell. If the parameter map_barycenter_of_reference_cell is set to false, then the returned value will be the average of the vertex locations, as returned by the get_vertices() method.

Parameters
[in]cellThe cell for which you want to compute the center
[in]map_barycenter_of_reference_cellA flag that switches the algorithm for the computation of the cell center from transform_unit_to_real_cell() applied to the center of the reference cell to computing the vertex averages.

◆ project_real_point_to_unit_point_on_face()

template<int dim, int spacedim = dim>
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
inherited

Transform the point p on the real cell to the corresponding point on the reference cell, and then project this point to a (dim-1)-dimensional point in the coordinate system of 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.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 136 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 156 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 204 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 53 of file subscriptor.cc.

Member Data Documentation

◆ euler_vector

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
SmartPointer<const VectorType, MappingQEulerian<dim, VectorType, spacedim> > MappingQEulerian< dim, VectorType, spacedim >::euler_vector
protected

Reference to the vector of shifts.

Definition at line 179 of file mapping_q_eulerian.h.

◆ euler_dof_handler

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
SmartPointer<const DoFHandler<dim, spacedim>, MappingQEulerian<dim, VectorType, spacedim> > MappingQEulerian< dim, VectorType, spacedim >::euler_dof_handler
protected

Pointer to the DoFHandler to which the mapping vector is associated.

Definition at line 186 of file mapping_q_eulerian.h.

◆ level

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
const unsigned int MappingQEulerian< dim, VectorType, spacedim >::level
private

Multigrid level at which the mapping is to be used.

Definition at line 192 of file mapping_q_eulerian.h.

◆ support_quadrature

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
const SupportQuadrature MappingQEulerian< dim, VectorType, spacedim >::support_quadrature
private

A member variable holding the quadrature points in the right order.

Definition at line 210 of file mapping_q_eulerian.h.

◆ fe_values

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
FEValues<dim, spacedim> MappingQEulerian< dim, VectorType, spacedim >::fe_values
mutableprivate

FEValues object used to query the given finite element field at the support points in the reference configuration.

The variable is marked as mutable since we have to call FEValues::reinit from compute_mapping_support_points, a function that is 'const'.

Definition at line 220 of file mapping_q_eulerian.h.

◆ fe_values_mutex

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
Threads::Mutex MappingQEulerian< dim, VectorType, spacedim >::fe_values_mutex
mutableprivate

A variable to guard access to the fe_values variable.

Definition at line 225 of file mapping_q_eulerian.h.

◆ polynomial_degree

template<int dim, int spacedim = dim>
const unsigned int MappingQ< dim, spacedim >::polynomial_degree
protectedinherited

The degree of the polynomials used as shape functions for the mapping of cells.

Definition at line 534 of file mapping_q.h.

◆ line_support_points

template<int dim, int spacedim = dim>
const std::vector<Point<1> > MappingQ< dim, spacedim >::line_support_points
protectedinherited

Definition at line 544 of file mapping_q.h.

◆ polynomials_1d

template<int dim, int spacedim = dim>
const std::vector<Polynomials::Polynomial<double> > MappingQ< dim, spacedim >::polynomials_1d
protectedinherited

Definition at line 551 of file mapping_q.h.

◆ renumber_lexicographic_to_hierarchic

template<int dim, int spacedim = dim>
const std::vector<unsigned int> MappingQ< dim, spacedim >::renumber_lexicographic_to_hierarchic
protectedinherited

Definition at line 558 of file mapping_q.h.

◆ unit_cell_support_points

template<int dim, int spacedim = dim>
const std::vector<Point<dim> > MappingQ< dim, spacedim >::unit_cell_support_points
protectedinherited

Definition at line 570 of file mapping_q.h.

◆ support_point_weights_perimeter_to_interior

template<int dim, int spacedim = dim>
const std::vector<Table<2, double> > MappingQ< dim, spacedim >::support_point_weights_perimeter_to_interior
protectedinherited

A vector of tables of weights by which we multiply the locations of the support points on the perimeter of an object (line, quad, hex) to get the location of interior support points.

Access into this table is by [structdim-1], i.e., use 0 to access the support point weights on a line (i.e., the interior points of the GaussLobatto quadrature), use 1 to access the support point weights from to perimeter to the interior of a quad, and use 2 to access the support point weights from the perimeter to the interior of a hex.

The table itself contains as many columns as there are surrounding points to a particular object (2 for a line, 4 + 4*(degree-1) for a quad, 8 + 12*(degree-1) + 6*(degree-1)*(degree-1) for a hex) and as many rows as there are strictly interior points.

For the definition of this table see equation (8) of the ‘mapping’ report.

Definition at line 592 of file mapping_q.h.

◆ support_point_weights_cell

template<int dim, int spacedim = dim>
const Table<2, double> MappingQ< dim, spacedim >::support_point_weights_cell
protectedinherited

A table of weights by which we multiply the locations of the vertex points of the cell to get the location of all additional support points, both on lines, quads, and hexes (as appropriate). This data structure is used when we fill all support points at once, which is the case if the same manifold is attached to all sub-entities of a cell. This way, we can avoid some of the overhead in transforming data for mappings.

The table has as many rows as there are vertices to the cell (2 in 1d, 4 in 2d, 8 in 3d), and as many rows as there are additional support points in the mapping, i.e., (degree+1)^dim - 2^dim.

Definition at line 606 of file mapping_q.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 219 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 225 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 241 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 249 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when printing out the list of subscribers.

Definition at line 271 of file subscriptor.h.


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