Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/mapping_q_generic.h>
Classes | |
class | InternalData |
Public Member Functions | |
MappingQGeneric (const unsigned int polynomial_degree) | |
MappingQGeneric (const MappingQGeneric< dim, spacedim > &mapping) | |
virtual std::unique_ptr< Mapping< dim, spacedim > > | clone () const override |
unsigned int | get_degree () const |
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 |
Interface with FEValues | |
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 Quadrature< 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 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 |
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 () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (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) |
Protected Member Functions | |
virtual std::vector< Point< spacedim > > | compute_mapping_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const |
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 |
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 |
Interface with FEValues |
Protected Attributes | |
const unsigned int | polynomial_degree |
const std::unique_ptr< FE_Q< dim > > | fe_q |
std::vector< Table< 2, double > > | support_point_weights_perimeter_to_interior |
Table< 2, double > | support_point_weights_cell |
Friends | |
template<int , int > | |
class | MappingQ |
Additional Inherited Members | |
Static Public Member Functions inherited from Mapping< dim, spacedim > | |
static ::ExceptionBase & | ExcInvalidData () |
static ::ExceptionBase & | ExcTransformationFailed () |
static ::ExceptionBase & | ExcDistortedMappedCell (Point< spacedim > arg1, double arg2, int arg3) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
This class implements the functionality for polynomial mappings \(Q_p\) of polynomial degree \(p\) that will be used on all cells of the mesh. The MappingQ1 and MappingQ classes specialize this behavior slightly.
The class is poorly named. It should really have been called MappingQ because it consistently uses \(Q_p\) mappings on all cells of a triangulation. However, the name MappingQ was already taken when we rewrote the entire class hierarchy for mappings. One might argue that one should always use MappingQGeneric over the existing class MappingQ (which, unless explicitly specified during the construction of the object, only uses mappings of degree \(p\) on cells at the boundary of the domain). On the other hand, there are good reasons to use MappingQ in many situations: in many situations, curved domains are only provided with information about how exactly edges at the boundary are shaped, but we do not know anything about internal edges. Thus, in the absence of other information, we can only assume that internal edges are straight lines, and in that case internal cells may as well be treated is bilinear quadrilaterals or trilinear hexahedra. (An example of how such meshes look is shown in step-1 already, but it is also discussed in the "Results" section of step-6.) Because bi-/trilinear mappings are significantly cheaper to compute than higher order mappings, it is advantageous in such situations to use the higher order mapping only on cells at the boundary of the domain – i.e., the behavior of MappingQ. Of course, MappingQGeneric also uses bilinear mappings for interior cells as long as it has no knowledge about curvature of interior edges, but it implements this the expensive way: as a general \(Q_p\) mapping where the mapping support points just happen to be arranged along linear or bilinear edges or faces.
There are a number of special cases worth considering:
As described above, one often only knows a manifold description of a surface but not the interior of the computational domain. In such a case, a FlatManifold object will be assigned to the interior entities that describes a usual planar coordinate system where the additional points for the higher order mapping are placed exactly according to a bi-/trilinear mapping. When combined with a non-flat manifold on the boundary, for example a circle bulging into the interior of a square cell, the two manifold descriptions are in general incompatible. For example, a FlatManifold defined solely through the cell's vertices would put an interior point located at some small distance epsilon away from the boundary along a straight line and thus in general outside the concave part of a circle. If the polynomial degree of MappingQ is sufficiently high, the transformation from the reference cell to such a cell would in general contain inverted regions close to the boundary.
In order to avoid this situation, this class applies an algorithm for making this transition smooth using a so-called transfinite interpolation that is essentially a linear blend between the descriptions along the surrounding entities. In the algorithm that computes additional points, the compute_mapping_support_points() method, all the entities of the cells are passed through hierarchically, starting from the lines to the quads and finally hexes. Points on objects higher up in the hierarchy are obtained from the manifold associated with that object, taking into account all the points previously computed by the manifolds associated with the lower-dimensional objects, not just the vertices. If only a line is assigned a curved boundary but the adjacent quad is on a flat manifold, the flat manifold on the quad will take the points on the deformed line into account when interpolating the position of the additional points inside the quad and thus always result in a well-defined transformation.
The interpolation scheme used in this class makes sure that curved descriptions can go over to flat descriptions within a single layer of elements, maintaining the overall optimal convergence rates of the finite element interpolation. However, one does often get better solution qualities if the transition between curved boundaries and flat interior domains is spread over a larger range as the mesh is refined. This is provided by the special manifold TransfiniteInterpolationManifold.
Definition at line 128 of file mapping_q_generic.h.
MappingQGeneric< dim, spacedim >::MappingQGeneric | ( | const unsigned int | polynomial_degree | ) |
Constructor. polynomial_degree
denotes the polynomial degree of the polynomials that are used to map cells from the reference to the real cell.
Definition at line 2112 of file mapping_q_generic.cc.
MappingQGeneric< dim, spacedim >::MappingQGeneric | ( | const MappingQGeneric< dim, spacedim > & | mapping | ) |
Copy constructor.
Definition at line 2127 of file mapping_q_generic.cc.
|
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 >.
Reimplemented in MappingQ1Eulerian< dim, VectorType, spacedim >, and MappingQ1< dim, spacedim >.
Definition at line 2141 of file mapping_q_generic.cc.
unsigned int MappingQGeneric< dim, spacedim >::get_degree | ( | ) | const |
Return the degree of the mapping, i.e. the value which was passed to the constructor.
Definition at line 2151 of file mapping_q_generic.cc.
|
overridevirtual |
Always returns true
because the default implementation of functions in this class preserves vertex locations.
Implements Mapping< dim, spacedim >.
Reimplemented in MappingQEulerian< dim, VectorType, spacedim >::MappingQEulerianGeneric, and MappingQ1Eulerian< dim, VectorType, spacedim >.
|
overridevirtual |
Map the point p
on the unit cell to the corresponding point on the real cell cell
.
cell | Iterator to the cell that will be used to define the mapping. |
p | Location of a point on the reference cell. |
Implements Mapping< dim, spacedim >.
Definition at line 2161 of file mapping_q_generic.cc.
|
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
.
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.cell | Iterator to the cell that will be used to define the mapping. |
p | Location of a point on the given cell. |
Implements Mapping< dim, spacedim >.
Definition at line 2371 of file mapping_q_generic.cc.
|
overridevirtual |
Transform a field of vectors or 1-differential forms according to the selected MappingType.
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}). \]
[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | type | The kind of mapping to be applied. |
[in] | internal | A 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] | output | An 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 3467 of file mapping_q_generic.cc.
|
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}. \]
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.[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | type | The kind of mapping to be applied. |
[in] | internal | A 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] | output | An 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 3480 of file mapping_q_generic.cc.
|
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}. \]
[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | type | The kind of mapping to be applied. |
[in] | internal | A 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] | output | An 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 3493 of file mapping_q_generic.cc.
|
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}\]
[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | type | The kind of mapping to be applied. |
[in] | internal | A 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] | output | An 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 3519 of file mapping_q_generic.cc.
|
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}. \]
[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | type | The kind of mapping to be applied. |
[in] | internal | A 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] | output | An array (or part of an array) into which the transformed objects should be placed. |
Implements Mapping< dim, spacedim >.
Definition at line 3568 of file mapping_q_generic.cc.
|
overridevirtual |
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.
Implements Mapping< dim, spacedim >.
Definition at line 2502 of file mapping_q_generic.cc.
|
overridevirtual |
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.
update_flags | A 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. |
quadrature | The quadrature object for which mapping information will have to be computed. This includes the locations and weights of quadrature points. |
Implements Mapping< dim, spacedim >.
Definition at line 2562 of file mapping_q_generic.cc.
|
overridevirtual |
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.
update_flags | A 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. |
quadrature | The quadrature object for which mapping information will have to be computed. This includes the locations and weights of quadrature points. |
Implements Mapping< dim, spacedim >.
Definition at line 2575 of file mapping_q_generic.cc.
|
overridevirtual |
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.
update_flags | A 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. |
quadrature | The quadrature object for which mapping information will have to be computed. This includes the locations and weights of quadrature points. |
Implements Mapping< dim, spacedim >.
Definition at line 2590 of file mapping_q_generic.cc.
|
overridevirtual |
Compute information about the mapping from the reference cell to the real cell indicated by the first argument to this function. Derived classes will have to implement this function based on the kind of mapping they represent. It is called by FEValues::reinit().
Conceptually, this function's represents the application of the mapping \(\mathbf x=\mathbf F_K(\hat {\mathbf x})\) from reference coordinates \(\mathbf\in [0,1]^d\) to real space coordinates \(\mathbf x\) for a given cell \(K\). Its purpose is to compute the following kinds of data:
The information computed by this function is used to fill the various member variables of the output argument of this function. Which of the member variables of that structure should be filled is determined by the update flags stored in the Mapping::InternalDataBase object passed to this function.
An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module.
[in] | cell | The cell of the triangulation for which this function is to compute a mapping from the reference cell to. |
[in] | cell_similarity | Whether or not the cell given as first argument is simply a translation, rotation, etc of the cell for which this function was called the most recent time. This information is computed simply by matching the vertices (as stored by the Triangulation) between the previous and the current cell. The value passed here may be modified by implementations of this function and should then be returned (see the discussion of the return value of this function). |
[in] | quadrature | A 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_data | A 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_data | A 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. |
cell_similarity
argument to this function. The returned value will be used for the corresponding argument when FEValues::reinit() calls FiniteElement::fill_fe_values(). In most cases, derived classes will simply want to return the value passed for cell_similarity
. However, implementations of this function may downgrade the level of cell similarity. This is, for example, the case for classes that take not only into account the locations of the vertices of a cell (as reported by the Triangulation), but also other information specific to the mapping. The purpose is that FEValues::reinit() can compute whether a cell is similar to the previous one only based on the cell's vertices, whereas the mapping may also consider displacement fields (e.g., in the MappingQ1Eulerian and MappingFEField classes). In such cases, the mapping may conclude that the previously computed cell similarity is too optimistic, and invalidate it for subsequent use in FiniteElement::fill_fe_values() by returning a less optimistic cell similarity value.internal_data
and output_data
objects. In other words, if an implementation of this function knows that it has written a piece of data into the output argument in a previous call, then there is no need to copy it there again in a later call if the implementation knows that this is the same value. Implements Mapping< dim, spacedim >.
Definition at line 2606 of file mapping_q_generic.cc.
|
overridevirtual |
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().
[in] | cell | The cell of the triangulation for which this function is to compute a mapping from the reference cell to. |
[in] | face_no | The number of the face of the given cell for which information is requested. |
[in] | quadrature | A 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_data | A 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_data | A 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 3048 of file mapping_q_generic.cc.
|
overridevirtual |
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().
[in] | cell | The cell of the triangulation for which this function is to compute a mapping from the reference cell to. |
[in] | face_no | The number of the face of the given cell for which information is requested. |
[in] | subface_no | The number of the child of a face of the given cell for which information is requested. |
[in] | quadrature | A 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_data | A 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_data | A 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 3093 of file mapping_q_generic.cc.
|
protectedvirtual |
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 in MappingQEulerian< dim, VectorType, spacedim >::MappingQEulerianGeneric, and MappingQ1Eulerian< dim, VectorType, spacedim >.
Definition at line 3820 of file mapping_q_generic.cc.
|
protected |
Transform the point p
on the real cell to the corresponding point on the unit cell cell
by a Newton iteration.
Definition at line 2211 of file mapping_q_generic.cc.
|
protectedvirtual |
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 3590 of file mapping_q_generic.cc.
|
protectedvirtual |
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 3809 of file mapping_q_generic.cc.
Make MappingQ a friend since it needs to call the fill_fe_values() functions on its MappingQGeneric(1) sub-object.
Definition at line 729 of file mapping_q_generic.h.
|
protected |
The degree of the polynomials used as shape functions for the mapping of cells.
Definition at line 591 of file mapping_q_generic.h.
|
protected |
An FE_Q object which is only needed in 3D, since it knows how to reorder shape functions/DoFs on non-standard faces. This is used to reorder support points in the same way.
Definition at line 609 of file mapping_q_generic.h.
|
protected |
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 630 of file mapping_q_generic.h.
|
protected |
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 644 of file mapping_q_generic.h.