Reference documentation for deal.II version 8.5.1
|
#include <deal.II/fe/mapping_fe_field.h>
Classes | |
class | InternalData |
Public Member Functions | |
MappingFEField (const DoFHandlerType &euler_dof_handler, const VectorType &euler_vector, const ComponentMask mask=ComponentMask()) | |
MappingFEField (const MappingFEField< dim, spacedim, VectorType, DoFHandlerType > &mapping) | |
virtual Mapping< dim, spacedim > * | clone () const |
virtual bool | preserves_vertex_locations () const |
unsigned int | get_degree () const |
ComponentMask | get_component_mask () const |
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 |
virtual Point< dim > | transform_real_to_unit_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const |
Functions to transform tensors from reference to real coordinates | |
virtual void | transform (const ArrayView< const Tensor< 1, dim > > &input, const MappingType type, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const |
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 |
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 |
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 |
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 |
Public Member Functions inherited from Mapping< dim, spacedim > | |
virtual | ~Mapping () |
virtual std_cxx11::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 &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInactiveCell () |
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, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Private Member Functions | |
Point< spacedim > | do_transform_unit_to_real_cell (const InternalData &mdata) const |
Point< dim > | do_transform_real_to_unit_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit, InternalData &mdata) const |
void | update_internal_dofs (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const |
virtual void | compute_shapes_virtual (const std::vector< Point< dim > > &unit_points, typename MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData &data) const |
Interface with FEValues | |
virtual UpdateFlags | requires_update_flags (const UpdateFlags update_flags) const |
virtual InternalData * | get_data (const UpdateFlags, const Quadrature< dim > &quadrature) const |
virtual Mapping< dim, spacedim >::InternalDataBase * | get_face_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const |
virtual Mapping< dim, spacedim >::InternalDataBase * | get_subface_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const |
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::FEValues::MappingRelatedData< dim, spacedim > &output_data) 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::FEValues::MappingRelatedData< dim, spacedim > &output_data) const |
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::FEValues::MappingRelatedData< dim, spacedim > &output_data) const |
Private Attributes | |
SmartPointer< const VectorType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > | euler_vector |
SmartPointer< const FiniteElement< dim, spacedim >, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > | fe |
SmartPointer< const DoFHandlerType, MappingFEField< dim, spacedim, VectorType, DoFHandlerType > > | euler_dof_handler |
std::vector< unsigned int > | fe_to_real |
Friends | |
template<int , int , class , class > | |
class | MappingFEField |
Additional Inherited Members | |
Protected Member Functions inherited from Mapping< dim, spacedim > | |
virtual CellSimilarity::Similarity | fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const =0 |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const =0 |
virtual void | fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValues::MappingRelatedData< dim, spacedim > &output_data) const =0 |
The MappingFEField is a generalization of the MappingQEulerian class, for arbitrary vector finite elements. The two main differences are that this class uses a vector of absolute positions instead of a vector of displacements, and it allows for arbitrary FiniteElement types, instead of only FE_Q.
This class effectively decouples the topology from the geometry, by relegating all geometrical information to some components of a FiniteElement vector field. The components that are used for the geometry can be arbitrarily selected at construction time.
The idea is to consider the Triangulation as a parameter configuration space, on which we construct an arbitrary geometrical mapping, using the instruments of the deal.II library: a vector of degrees of freedom, a DoFHandler associated to the geometry of the problem and a ComponentMask that tells us which components of the FiniteElement to use for the mapping.
Typically, the DoFHandler operates on a finite element that is constructed as a system element (FESystem()) from continuous FE_Q() (for iso-parametric discretizations) or FE_Bernstein() (for iso-geometric discretizations) objects. An example is shown below:
Definition at line 77 of file mapping_fe_field.h.
MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::MappingFEField | ( | const DoFHandlerType & | euler_dof_handler, |
const VectorType & | euler_vector, | ||
const ComponentMask | mask = ComponentMask() |
||
) |
Constructor. The first argument is a VectorType that specifies the transformation of the domain from the reference to the current configuration.
In general this class decouples geometry from topology, allowing users to define geometries which are only topologically equivalent to the underlying Triangulation, but which may otherwise be arbitrary. Differently from what happens in MappingQEulerian, the FiniteElement field which is passed to the constructor is interpreted as an absolute geometrical configuration, therefore one has to make sure that the euler_vector actually represents a valid geometry (i.e., one with no inverted cells, or with no zero-volume cells).
If the underlying FiniteElement is a system of FE_Q(), and euler_vector is initialized using VectorTools::get_position_vector(), then this class is in all respects identical to MappingQ().
The optional ComponentMask argument can be used to specify what components of the FiniteElement to use for the geometrical transformation. If no mask is specified at construction time, then a default one is used, which makes this class works in the same way of MappingQEulerian(), i.e., the first spacedim components of the FiniteElement are assumed to represent the geometry of the problem.
Notice that if a mask is specified, it has to match in size the underlying FiniteElement, and it has to have exactly spacedim non-zero elements, indicating the components (in order) of the FiniteElement which will be used for the geometry.
If an incompatible mask is passed, an exception is thrown.
Definition at line 204 of file mapping_fe_field.cc.
MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::MappingFEField | ( | const MappingFEField< dim, spacedim, VectorType, DoFHandlerType > & | mapping | ) |
Copy constructor.
Definition at line 227 of file mapping_fe_field.cc.
|
virtual |
Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.
Implements Mapping< dim, spacedim >.
Definition at line 1969 of file mapping_fe_field.cc.
|
virtual |
Always returns false
.
Implements Mapping< dim, spacedim >.
Definition at line 255 of file mapping_fe_field.cc.
|
virtual |
Maps 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 1765 of file mapping_fe_field.cc.
|
virtual |
Maps 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 1803 of file mapping_fe_field.cc.
|
virtual |
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(\mathbf x)} J(\mathbf x) \hat{\mathbf u}(\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 1648 of file mapping_fe_field.cc.
|
virtual |
Transform a field of differential forms from the reference cell to the physical cell. It is useful to think of \(\mathbf{T} = \nabla \mathbf u\) and \(\hat{\mathbf T} = \hat \nabla \hat{\mathbf u}\), with \(\mathbf u\) a vector field. The mapping types currently implemented by derived classes are:
mapping_covariant:
maps a field of forms on the reference cell to a field of forms on the physical cell. Mathematically, it is the pull back of the differential form
\[ \mathbf T(\mathbf x) = \hat{\mathbf T}(\hat{\mathbf x}) J(\hat{\mathbf x})(J(\hat{\mathbf x})^{T} J(\hat{\mathbf x}))^{-1}. \]
Jacobians of spacedim-vector valued differentiable functions are transformed this way.
In the case when dim=spacedim the previous formula reduces to
\[ \mathbf T(\mathbf x) = \hat{\mathbf u}(\hat{\mathbf x}) J(\hat{\mathbf x})^{-1}. \]
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 1663 of file mapping_fe_field.cc.
|
virtual |
Transform a tensor field from the reference cell to the physical cell. These tensors are usually the Jacobians in the reference cell of vector fields that have been pulled back from the physical cell. The mapping types currently implemented by derived classes are:
mapping_contravariant_gradient:
it assumes \(\mathbf u(\mathbf x) = J \hat{\mathbf u}\) so that \[ \mathbf T(\mathbf x) = J(\hat{\mathbf x}) \hat{\mathbf T}(\hat{\mathbf x}) J(\hat{\mathbf x})^{-1}. \]
mapping_covariant_gradient:
it assumes \(\mathbf u(\mathbf x) = J^{-T} \hat{\mathbf u}\) so that \[ \mathbf T(\mathbf x) = J(\hat{\mathbf x})^{-T} \hat{\mathbf T}(\hat{\mathbf x}) J(\hat{\mathbf x})^{-1}. \]
mapping_piola_gradient:
it assumes \(\mathbf u(\mathbf x) = \frac{1}{\text{det}\;J(\mathbf x)} J(\mathbf x) \hat{\mathbf u}(\mathbf x)\) so that \[ \mathbf T(\mathbf x) = \frac{1}{\text{det}\;J(\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 1678 of file mapping_fe_field.cc.
|
virtual |
Transform a tensor field from the reference cell to the physical cell. This tensors are most of times the hessians in the reference cell of vector fields that have been pulled back from the physical cell.
The mapping types currently implemented by derived classes are:
mapping_covariant_gradient:
maps a field of forms on the reference cell to a field of forms on the physical cell. Mathematically, it is the pull back of the differential form
\[ \mathbf T_{ijk}(\mathbf x) = \hat{\mathbf T}_{iJK}(\hat{\mathbf x}) J_{jJ}^{\dagger} J_{kK}^{\dagger}\]
,
where
\[ J^{\dagger} = J(\hat{\mathbf x})(J(\hat{\mathbf x})^{T} J(\hat{\mathbf x}))^{-1}. \]
Hessians of spacedim-vector valued differentiable functions are transformed this way (After subtraction of the product of the derivative with the Jacobian gradient).
In the case when dim=spacedim the previous formula reduces to
\[J^{\dagger} = J^{-1}\]
[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 1696 of file mapping_fe_field.cc.
|
virtual |
Transform a field of 3-differential forms from the reference cell to the physical cell. It is useful to think of \(\mathbf{T}_{ijk} = D^2_{jk} \mathbf u_i\) and \(\mathbf{\hat T}_{IJK} = \hat D^2_{JK} \mathbf{\hat u}_I\), with \(\mathbf u_i\) a vector field.
The mapping types currently implemented by derived classes are:
mapping_contravariant_hessian:
it assumes \(\mathbf u_i(\mathbf x) = J_{iI} \hat{\mathbf u}_I\) so that \[ \mathbf T_{ijk}(\mathbf x) = J_{iI}(\hat{\mathbf x}) \hat{\mathbf T}_{IJK}(\hat{\mathbf x}) J_{jJ}(\hat{\mathbf x})^{-1} J_{kK}(\hat{\mathbf x})^{-1}. \]
mapping_covariant_hessian:
it assumes \(\mathbf u_i(\mathbf x) = J_{iI}^{-T} \hat{\mathbf u}_I\) so that \[ \mathbf T_{ijk}(\mathbf x) = J_iI(\hat{\mathbf x})^{-1} \hat{\mathbf T}_{IJK}(\hat{\mathbf x}) J_{jJ}(\hat{\mathbf x})^{-1} J_{kK}(\hat{\mathbf x})^{-1}. \]
mapping_piola_hessian:
it assumes \(\mathbf u_i(\mathbf x) = \frac{1}{\text{det}\;J(\mathbf x)} J_{iI}(\mathbf x) \hat{\mathbf u}(\mathbf x)\) so that \[ \mathbf T_{ijk}(\mathbf x) = \frac{1}{\text{det}\;J(\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 1745 of file mapping_fe_field.cc.
unsigned int MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::get_degree | ( | ) | const |
Return the degree of the mapping, i.e. the value which was passed to the constructor.
Definition at line 1953 of file mapping_fe_field.cc.
ComponentMask MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::get_component_mask | ( | ) | const |
Return the ComponentMask of the mapping, i.e. which components to use for the mapping.
Definition at line 1961 of file mapping_fe_field.cc.
|
privatevirtual |
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 297 of file mapping_fe_field.cc.
|
privatevirtual |
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 479 of file mapping_fe_field.cc.
|
privatevirtual |
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 493 of file mapping_fe_field.cc.
|
privatevirtual |
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 508 of file mapping_fe_field.cc.
|
private |
Transforms a point p
on the unit cell to the point p_real
on the real cell cell
and returns p_real
.
This function is called by transform_unit_to_real_cell
and multiple times (through the Newton iteration) by transform_real_to_unit_cell_internal
.
Takes a reference to an InternalData
that must already include the shape values at point p
and the mapping support points of the cell.
This InternalData
argument avoids multiple computations of the shape values at point p
and especially multiple computations of the mapping support points.
Definition at line 1784 of file mapping_fe_field.cc.
|
private |
Transforms the point p
on the real cell to the corresponding point on the unit cell cell
by a Newton iteration.
Takes a reference to an InternalData
that is assumed to be previously created by the get_data
function with UpdateFlags
including update_transformation_values
and update_transformation_gradients
and a one point Quadrature that includes the given initial guess for the transformation initial_p_unit
. Hence this function assumes that mdata
already includes the transformation shape values and gradients computed at initial_p_unit
.
mdata
will be changed by this function.
Definition at line 1847 of file mapping_fe_field.cc.
|
private |
Update internal degrees of freedom.
Definition at line 1978 of file mapping_fe_field.cc.
|
privatevirtual |
See the documentation of the base class for detailed information.
Definition at line 265 of file mapping_fe_field.cc.
|
friend |
Declare other MappingFEField classes friends.
Definition at line 597 of file mapping_fe_field.h.
|
private |
Reference to the vector of shifts.
Definition at line 495 of file mapping_fe_field.h.
|
private |
A FiniteElement 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. We could make this a pointer to prevent construction in 1D and 2D, but since memory and time requirements are not particularly high this seems unnecessary at the moment.
Definition at line 504 of file mapping_fe_field.h.
|
private |
Pointer to the DoFHandler to which the mapping vector is associated.
Definition at line 510 of file mapping_fe_field.h.
|
private |
Mapping between indices in the FE space and the real space. This vector contains one index for each component of the finite element space. If the index is one for which the ComponentMask which is used to construct this element is false, then numbers::invalid_unsigned_int is returned, otherwise the component in real space is returned. For example, if we construct the mapping using ComponentMask(spacedim, true), then this vector contains {0,1,2} in spacedim = 3.
Definition at line 579 of file mapping_fe_field.h.