deal.II version GIT relicensing-2055-gd703fda1ad 2024-10-30 20:00:00+00:00
|
#include <deal.II/fe/mapping_cartesian.h>
Classes | |
class | InternalData |
Public Member Functions | |
virtual std::unique_ptr< Mapping< dim, spacedim > > | clone () const override |
virtual bool | preserves_vertex_locations () 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 |
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > | get_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell) 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 |
virtual BoundingBox< spacedim > | get_bounding_box (const typename Triangulation< dim, spacedim >::cell_iterator &cell) 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 |
EnableObserverPointer functionality | |
Classes derived from EnableObserverPointer provide a facility to subscribe to this object. This is mostly used by the ObserverPointer 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 ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Exceptions | |
static ::ExceptionBase & | ExcInvalidData () |
static ::ExceptionBase & | ExcTransformationFailed () |
static ::ExceptionBase & | ExcDistortedMappedCell (Point< spacedim > arg1, double arg2, int arg3) |
Protected Member Functions | |
Interface with FEValues | |
virtual std::unique_ptr< InternalDataBase > | get_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 |
Private Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
void | update_cell_extents (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const InternalData &data) const |
void | maybe_update_cell_quadrature_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const ArrayView< const Point< dim > > &unit_quadrature_points, std::vector< Point< dim > > &quadrature_points) const |
void | maybe_update_face_quadrature_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const InternalData &data, std::vector< Point< dim > > &quadrature_points) const |
void | maybe_update_subface_quadrature_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const InternalData &data, std::vector< Point< dim > > &quadrature_points) const |
void | maybe_update_normal_vectors (const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, dim > > &normal_vectors) const |
void | maybe_update_jacobian_derivatives (const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const |
void | maybe_update_volume_elements (const InternalData &data) const |
void | maybe_update_jacobians (const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const |
void | maybe_update_inverse_jacobians (const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const |
void | check_no_subscribers () const noexcept |
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 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 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 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 |
Private Attributes | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
Static Private Attributes | |
static std::mutex | mutex |
A class providing a mapping from the reference cell to cells that are axiparallel, i.e., that have the shape of rectangles (in 2d) or boxes (in 3d) with edges parallel to the coordinate directions. The class therefore provides functionality that is equivalent to what, for example, MappingQ would provide for such cells. However, knowledge of the shape of cells allows this class to be substantially more efficient.
Specifically, the mapping is meant for cells for which the mapping from the reference to the real cell is a scaling along the coordinate directions: The transformation from reference coordinates \(\hat {\mathbf x}\) to real coordinates \(\mathbf x\) on each cell is of the form
\begin{align*} {\mathbf x}(\hat {\mathbf x}) = \begin{pmatrix} h_x & 0 \\ 0 & h_y \end{pmatrix} \hat{\mathbf x} + {\mathbf v}_0 \end{align*}
in 2d, and
\begin{align*} {\mathbf x}(\hat {\mathbf x}) = \begin{pmatrix} h_x & 0 & 0 \\ 0 & h_y & 0 \\ 0 & 0 & h_z \end{pmatrix} \hat{\mathbf x} + {\mathbf v}_0 \end{align*}
in 3d, where \({\mathbf v}_0\) is the bottom left vertex and \(h_x,h_y,h_z\) are the extents of the cell along the axes.
The class is intended for efficiency, and it does not do a whole lot of error checking. If you apply this mapping to a cell that does not conform to the requirements above, you will get strange results.
Definition at line 78 of file mapping_cartesian.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 238 of file enable_observer_pointer.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 243 of file enable_observer_pointer.h.
|
overridevirtual |
Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.
The function is declared abstract virtual in this base class, and derived classes will have to implement it.
This function is mainly used by the hp::MappingCollection class.
Implements Mapping< dim, spacedim >.
Definition at line 1256 of file mapping_cartesian.cc.
|
overridevirtual |
Return true
because MappingCartesian preserves vertex locations.
Implements Mapping< dim, spacedim >.
Definition at line 129 of file mapping_cartesian.cc.
|
overridevirtual |
Returns if this instance of Mapping is compatible with the type of cell in reference_cell
.
Implements Mapping< dim, spacedim >.
Definition at line 138 of file mapping_cartesian.cc.
|
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 1188 of file mapping_cartesian.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 1208 of file mapping_cartesian.cc.
|
overridevirtual |
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 1229 of file mapping_cartesian.cc.
|
overridevirtual |
Transform a field of vectors or 1-differential forms according to the selected MappingKind.
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}). \]
[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | kind | 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 772 of file mapping_cartesian.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 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}. \]
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] | kind | 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 832 of file mapping_cartesian.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 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}. \]
[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | kind | 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 945 of file mapping_cartesian.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 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}\]
[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | kind | 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 1058 of file mapping_cartesian.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 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}. \]
[in] | input | An array (or part of an array) of input objects that should be mapped. |
[in] | kind | 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 1097 of file mapping_cartesian.cc.
void MappingCartesian< 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 |
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.
[in] | cell | The cell where to evaluate the mapping |
[in] | unit_points | The points in reference coordinates where the transformation (Jacobians, positions) should be computed. |
[in] | update_flags | The kind of information that should be computed. |
[out] | output_data | A 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 557 of file mapping_cartesian.cc.
|
overrideprivatevirtual |
Given a set of update flags, compute which other quantities also need to be computed in order to satisfy the request by the given flags. Then return the combination of the original set of flags and those just computed.
As an example, if update_flags
contains update_JxW_values (i.e., the product of the determinant of the Jacobian and the weights provided by the quadrature formula), a mapping may require the computation of the full Jacobian matrix in order to compute its determinant. They would then return not just update_JxW_values, but also update_jacobians. (This is not how it is actually done internally in the derived classes that compute the JxW values – they set update_contravariant_transformation instead, from which the determinant can also be computed – but this does not take away from the instructiveness of the example.)
An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation topic.
Implements Mapping< dim, spacedim >.
Definition at line 155 of file mapping_cartesian.cc.
|
overrideprivatevirtual |
Create and return a pointer to an object into which mappings can store data that only needs to be computed once but that can then be used whenever the mapping is applied to a concrete cell (e.g., in the various transform() functions, as well as in the fill_fe_values(), fill_fe_face_values() and fill_fe_subface_values() that form the interface of mappings with the FEValues class).
Derived classes will return pointers to objects of a type derived from Mapping::InternalDataBase (see there for more information) and may pre-compute some information already (in accordance with what will be asked of the mapping in the future, as specified by the update flags) and for the given quadrature object. Subsequent calls to transform() or fill_fe_values() and friends will then receive back the object created here (with the same set of update flags and for the same quadrature object). Derived classes can therefore pre-compute some information in their get_data() function and store it in the internal data object.
The mapping classes do not keep track of the objects created by this function. Ownership will therefore rest with the caller.
An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation topic.
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 174 of file mapping_cartesian.cc.
|
overrideprivatevirtual |
Like get_data(), but in preparation for later calls to transform() or fill_fe_face_values() that will need information about mappings from the reference face to a face of a concrete cell.
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. |
Reimplemented from Mapping< dim, spacedim >.
Definition at line 188 of file mapping_cartesian.cc.
|
overrideprivatevirtual |
Like get_data() and get_face_data(), but in preparation for later calls to transform() or fill_fe_subface_values() that will need information about mappings from the reference face to a child of a face (i.e., subface) of a concrete cell.
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 215 of file mapping_cartesian.cc.
|
overrideprivatevirtual |
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 topic.
[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 507 of file mapping_cartesian.cc.
|
overrideprivatevirtual |
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. |
Reimplemented from Mapping< dim, spacedim >.
Definition at line 594 of file mapping_cartesian.cc.
|
overrideprivatevirtual |
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 647 of file mapping_cartesian.cc.
|
overrideprivatevirtual |
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 706 of file mapping_cartesian.cc.
|
private |
Update the cell_extents field of the incoming InternalData object with the size of the incoming cell.
Definition at line 240 of file mapping_cartesian.cc.
|
private |
Compute the quadrature points if the UpdateFlags of the incoming InternalData object say that they should be updated.
Called from fill_fe_values.
Definition at line 288 of file mapping_cartesian.cc.
|
private |
Compute the quadrature points if the UpdateFlags of the incoming InternalData object say that they should be updated.
Called from fill_fe_face_values.
Definition at line 310 of file mapping_cartesian.cc.
|
private |
Compute the quadrature points if the UpdateFlags of the incoming InternalData object say that they should be updated.
Called from fill_fe_subface_values.
Definition at line 339 of file mapping_cartesian.cc.
|
private |
Compute the normal vectors if the UpdateFlags of the incoming InternalData object say that they should be updated.
Definition at line 375 of file mapping_cartesian.cc.
|
private |
Since the Jacobian is constant for this mapping all derivatives of the Jacobian are identically zero. Fill these quantities with zeros if the corresponding update flags say that they should be updated.
Definition at line 394 of file mapping_cartesian.cc.
|
private |
Compute the volume elements if the UpdateFlags of the incoming InternalData object say that they should be updated.
Definition at line 446 of file mapping_cartesian.cc.
|
private |
Compute the Jacobians if the UpdateFlags of the incoming InternalData object say that they should be updated.
Definition at line 462 of file mapping_cartesian.cc.
|
private |
Compute the inverse Jacobians if the UpdateFlags of the incoming InternalData object say that they should be updated.
Definition at line 484 of file mapping_cartesian.cc.
|
virtualinherited |
Return the mapped vertices of a cell.
Most of the time, these values will simply be the coordinates of the vertices of a cell as returned by cell->vertex(v)
for vertex v
, i.e., information stored by the triangulation. However, there are also mappings that add displacements or choose completely different locations, e.g., MappingQEulerian, MappingQ1Eulerian, or MappingFEField.
The default implementation of this function simply returns the information stored by the triangulation, i.e., cell->vertex(v)
.
Reimplemented in MappingFEField< dim, spacedim, VectorType >, MappingQ1Eulerian< dim, VectorType, spacedim >, MappingQCache< dim, spacedim >, and MappingQEulerian< dim, VectorType, spacedim >.
|
inherited |
Return the mapped vertices of a face.
Same as above but working on a given face of a cell.
[in] | cell | The cell containing the face. |
[in] | face_no | The number of the face within the cell. |
|
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.
[in] | cell | The cell for which you want to compute the center |
[in] | map_barycenter_of_reference_cell | A 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. |
|
virtualinherited |
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.
[in] | cell | The cell for which you want to compute the bounding box |
Reimplemented in MappingFE< dim, spacedim >, MappingQ< dim, spacedim >, and MappingQ< dim, dim >.
|
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.
|
protectedvirtualinherited |
|
protectedvirtualinherited |
|
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 131 of file enable_observer_pointer.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 151 of file enable_observer_pointer.cc.
|
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 322 of file enable_observer_pointer.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 339 of file enable_observer_pointer.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 199 of file enable_observer_pointer.cc.
|
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 331 of file enable_observer_pointer.h.
|
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.
Definition at line 53 of file enable_observer_pointer.cc.
|
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 227 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 233 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the ObserverPointer objects that subscribe to this class.
Definition at line 249 of file enable_observer_pointer.h.
|
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 257 of file enable_observer_pointer.h.
|
staticprivateinherited |
A mutex used to ensure data consistency when accessing the mutable
members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers()
.
Definition at line 280 of file enable_observer_pointer.h.