Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/mapping_q.h>
Classes | |
class | InternalData |
Public Member Functions | |
MappingQ (const unsigned int polynomial_degree, const bool use_mapping_q_on_all_cells=false) | |
MappingQ (const MappingQ< dim, spacedim > &mapping) | |
unsigned int | get_degree () const |
virtual bool | preserves_vertex_locations () const override |
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 (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 |
virtual std::unique_ptr< Mapping< dim, spacedim > > | clone () 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 | |
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 |
Protected Member Functions inherited from Mapping< dim, spacedim > | |
virtual CellSimilarity::Similarity | fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const =0 |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const =0 |
virtual void | fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim-1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const =0 |
Protected Attributes | |
const unsigned int | polynomial_degree |
const bool | use_mapping_q_on_all_cells |
std::shared_ptr< const MappingQGeneric< dim, spacedim > > | q1_mapping |
std::shared_ptr< const MappingQGeneric< dim, spacedim > > | qp_mapping |
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) |
A class that implements a polynomial mapping \(Q_p\) of degree \(p\) on cells at the boundary of the domain (or, if requested in the constructor, for all cells) and linear mappings for interior cells.
The class is in fact poorly named since (unless explicitly specified during the construction of the object, see below), it does not actually use mappings of degree \(p\) everywhere, but only on cells at the boundary. This is in contrast to the MappingQGeneric class which indeed does use a polynomial mapping \(Q_p\) of degree \(p\) everywhere. The point of the current class is that 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. This class implements exactly this behavior.
There are a number of special cases worth considering:
For the behavior of the mapping and convergence rates in case of mixing different manifolds, please consult the respective section of MappingQGeneric.
MappingQ< dim, spacedim >::MappingQ | ( | const unsigned int | polynomial_degree, |
const bool | use_mapping_q_on_all_cells = false |
||
) |
Constructor. polynomial_degree
denotes the polynomial degree of the polynomials that are used to map cells boundary.
The second argument determines whether the higher order mapping should also be used on interior cells. If its value is false
(the default), then a lower order mapping is used in the interior. This is sufficient for most cases where higher order mappings are only used to better approximate the boundary. In that case, cells bounded by straight lines are acceptable in the interior. However, there are cases where one would also like to use a higher order mapping in the interior. The MappingQEulerian class is one such case.
The value of use_mapping_q_on_all_cells
is ignored if dim
is not equal to spacedim
, i.e., if we are considering meshes on surfaces embedded into higher dimensional spaces.
Definition at line 58 of file mapping_q.cc.
MappingQ< dim, spacedim >::MappingQ | ( | const MappingQ< dim, spacedim > & | mapping | ) |
Copy constructor.
Definition at line 89 of file mapping_q.cc.
unsigned int MappingQ< dim, spacedim >::get_degree | ( | ) | const |
Return the degree of the mapping, i.e. the value which was passed to the constructor.
Definition at line 118 of file mapping_q.cc.
|
inlineoverridevirtual |
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 >.
Definition at line 128 of file mapping_q.cc.
|
overridevirtual |
Transform the point p
on the unit cell to the point p_real
on the real cell cell
and returns p_real
.
Implements Mapping< dim, spacedim >.
Definition at line 487 of file mapping_q.cc.
|
overridevirtual |
Transform the point p
on the real cell to the point p_unit
on the unit cell cell
and returns p_unit
.
Uses Newton iteration and the transform_unit_to_real_cell
function.
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 return reference coordinates lie inside of outside the reference cell (e.g., using GeometryInfo::is_inside_unit_cell) or whether the exception mentioned above has been thrown. Implements Mapping< dim, spacedim >.
Definition at line 504 of file mapping_q.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 369 of file mapping_q.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 392 of file mapping_q.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 416 of file mapping_q.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 440 of file mapping_q.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 463 of file mapping_q.cc.
|
overridevirtual |
Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.
Implements Mapping< dim, spacedim >.
Reimplemented in MappingQEulerian< dim, VectorType, spacedim >, and MappingC1< dim, spacedim >.
Definition at line 521 of file mapping_q.cc.
|
overrideprotectedvirtual |
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 137 of file mapping_q.cc.
|
overrideprotectedvirtual |
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 148 of file mapping_q.cc.
|
overrideprotectedvirtual |
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 174 of file mapping_q.cc.
|
overrideprotectedvirtual |
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 200 of file mapping_q.cc.
|
protected |
The polynomial degree of the cells to be used on all cells at the boundary of the domain, or everywhere if so specified.
Definition at line 334 of file mapping_q.h.
|
protected |
If this flag is set true
then MappingQ
is used on all cells, not only on boundary cells.
Definition at line 340 of file mapping_q.h.
|
protected |
Pointer to a Q1 mapping. This mapping is used on interior cells unless use_mapping_q_on_all_cells was set in the call to the constructor. The mapping is also used on any cell in the transform_real_to_unit_cell() to compute a cheap initial guess for the position of the point before we employ the more expensive Newton iteration using the full mapping.
Definition at line 359 of file mapping_q.h.
|
protected |
Pointer to a Q_p mapping. This mapping is used on boundary cells unless use_mapping_q_on_all_cells was set in the call to the constructor (in which case it is used for all cells).
Definition at line 375 of file mapping_q.h.