Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/mapping.h>
Public Member Functions | |
InternalDataBase () | |
InternalDataBase (const InternalDataBase &)=delete | |
virtual | ~InternalDataBase ()=default |
virtual std::size_t | memory_consumption () const |
Public Attributes | |
UpdateFlags | update_each |
Base class for internal data of mapping objects. The internal mechanism is that upon construction of a FEValues object, it asks the mapping and finite element classes that are to be used to allocate memory for their own purpose in which they may store data that only needs to be computed once. For example, most finite elements will store the values of the shape functions at the quadrature points in this object, since they do not change from cell to cell and only need to be computed once. The same may be true for Mapping classes that want to only evaluate the shape functions used for mapping once at the quadrature points.
Since different FEValues objects using different quadrature rules might access the same mapping object at the same time, it is necessary to create one such object per FEValues object. FEValues does this by calling Mapping::get_data(), or in reality the implementation of the corresponding function in derived classes. Ownership of the object created by Mapping::get_data() is then transferred to the FEValues object, but a reference to this object is passed to the mapping object every time it is asked to compute information on a concrete cell. This happens when FEValues::reinit() (or the corresponding classes in FEFaceValues and FESubfaceValues) call Mapping::fill_fe_values() (and similarly via Mapping::fill_fe_face_values() and Mapping::fill_fe_subface_values()).
The purpose of this class is for mapping objects to store information that can be computed once at the beginning, on the reference cell, and to access it later when computing information on a concrete cell. As such, the object handed to Mapping::fill_fe_values() is marked as const
, because the assumption is that at the time this information is used, it will not need to modified again. However, classes derived from Mapping can also use such objects for two other purposes:
In both of these cases, the InternalDataBase object being passed around is "morally const", i.e., no external observer can tell whether a scratch array or some intermediate data for Mapping::transform() is being modified by Mapping::fill_fe_values() or not. Consequently, the InternalDataBase objects are always passed around as const
objects. Derived classes that would like to make use of the two additional uses outlined above therefore need to mark the member variables they want to use for these purposes as mutable
to allow for their modification despite the fact that the surrounding object is marked as const
.
Mapping< dim, spacedim >::InternalDataBase::InternalDataBase | ( | ) |
Constructor. Sets update_flags to update_default
and first_cell
to true
.
Definition at line 78 of file mapping.cc.
|
delete |
Copy construction is forbidden.
|
virtualdefault |
Virtual destructor for derived classes
|
virtual |
Return an estimate (in bytes) or the memory consumption of this object.
Reimplemented in MappingQGeneric< dim, spacedim >::InternalData, MappingFEField< dim, spacedim, VectorType, DoFHandlerType >::InternalData, MappingQ< dim, spacedim >::InternalData, MappingManifold< dim, spacedim >::InternalData, and MappingCartesian< dim, spacedim >::InternalData.
Definition at line 86 of file mapping.cc.
UpdateFlags Mapping< dim, spacedim >::InternalDataBase::update_each |
A set of update flags specifying the kind of information that an implementation of the Mapping interface needs to compute on each cell or face, i.e., in Mapping::fill_fe_values() and friends.
This set of flags is stored here by implementations of Mapping::get_data(), Mapping::get_face_data(), or Mapping::get_subface_data(), and is that subset of the update flags passed to those functions that require re-computation on every cell. (The subset of the flags corresponding to information that can be computed once and for all already at the time of the call to Mapping::get_data() – or an implementation of that interface – need not be stored here because it has already been taken care of.)