Reference documentation for deal.II version 9.1.1
|
#include <deal.II/fe/mapping_q1_eulerian.h>
Public Member Functions | |
MappingQ1Eulerian (const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector) | |
virtual std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > | get_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override |
virtual std::unique_ptr< Mapping< dim, spacedim > > | clone () const override |
virtual bool | preserves_vertex_locations () const override |
Public Member Functions inherited from MappingQGeneric< dim, spacedim > | |
MappingQGeneric (const unsigned int polynomial_degree) | |
MappingQGeneric (const MappingQGeneric< dim, spacedim > &mapping) | |
unsigned int | get_degree () const |
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 UpdateFlags | requires_update_flags (const UpdateFlags update_flags) const override |
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > | get_data (const UpdateFlags, const Quadrature< dim > &quadrature) const override |
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > | get_face_data (const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override |
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > | get_subface_data (const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override |
virtual CellSimilarity::Similarity | fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override |
virtual void | fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override |
Public Member Functions inherited from Mapping< dim, spacedim > | |
virtual | ~Mapping () override=default |
virtual Point< spacedim > | get_center (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_center_of_reference_cell=true) const |
virtual BoundingBox< spacedim > | get_bounding_box (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 (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 |
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, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Member Functions | |
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 std::vector< Point< spacedim > > | compute_mapping_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const override |
Protected Member Functions inherited from MappingQGeneric< dim, spacedim > | |
Point< dim > | transform_real_to_unit_cell_internal (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p, const Point< dim > &initial_p_unit) const |
virtual void | add_line_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim >> &a) const |
virtual void | add_quad_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim >> &a) const |
Interface with FEValues |
Protected Attributes | |
SmartPointer< const VectorType, MappingQ1Eulerian< dim, VectorType, spacedim > > | euler_transform_vectors |
SmartPointer< const DoFHandler< dim, spacedim >, MappingQ1Eulerian< dim, VectorType, spacedim > > | shiftmap_dof_handler |
Protected Attributes inherited from MappingQGeneric< dim, spacedim > | |
const unsigned int | polynomial_degree |
const std::unique_ptr< FE_Q< dim > > | fe_q |
std::vector< Table< 2, double > > | support_point_weights_perimeter_to_interior |
Table< 2, double > | support_point_weights_cell |
This class provides a mapping that adds to the location of each cell a \(d\)-linear displacement field. (The generalization to higher order polynomials is provided in the MappingQEulerian class.) Each cell is thus shifted in space by values given to the mapping through a finite element field.
The constructor of this class takes two arguments: a reference to the vector that defines the mapping from the reference configuration to the current configuration and a reference to the DoFHandler. The vector should then represent a (flattened out version of a) vector valued field defined at nodes defined by the DoFHandler, where the number of components of the vector field equals the number of space dimensions. Thus, the DoFHandler shall operate on a finite element that has as many components as space dimensions. As an additional requirement, we impose that it have as many degree of freedom per vertex as there are space dimensions; since this object only evaluates the finite element field at the vertices, the values of all other degrees of freedom (not associated to vertices) are ignored. These requirements are met if the finite element which the given DoFHandler operates on is constructed as a system element (FESystem) from dim
continuous FE_Q() objects.
In many cases, the shift vector will also be the solution vector of the problem under investigation. If this is not the case (i.e. the number of components of the solution variable is not equal to the space dimension, e.g. for scalar problems in dim>1
where the Eulerian coordinates only give a background field) or for coupled problems where more variables are computed than just the flow field), then a different DoFHandler has to be set up on the given triangulation, and the shift vector has then to be associated to it.
An example is shown below:
Note that since the vector of shift values and the dof handler are only associated to this object at construction time, you have to make sure that whenever you use this object, the given objects still represent valid data.
To enable the use of the MappingQ1Eulerian class also in the context of parallel codes using the PETSc or Trilinos wrapper classes, the type of the vector can be specified as template parameter VectorType
.
For more information about the spacedim
template parameter check the documentation of FiniteElement or the one of Triangulation.
Definition at line 95 of file mapping_q1_eulerian.h.
MappingQ1Eulerian< dim, VectorType, spacedim >::MappingQ1Eulerian | ( | const DoFHandler< dim, spacedim > & | euler_dof_handler, |
const VectorType & | euler_vector | ||
) |
Constructor.
[in] | euler_dof_handler | A DoFHandler object that defines a finite element space. This space needs to have exactly dim components and these will be considered displacements relative to the original positions of the cells of the triangulation. This DoFHandler must be based on a FESystem(FE_Q(1),dim) finite element. |
[in] | euler_vector | A finite element function in the space defined by the first argument. The dim components of this function will be interpreted as the displacement we use in defining the mapping, relative to the location of cells of the underlying triangulation. |
Definition at line 42 of file mapping_q1_eulerian.cc.
|
overridevirtual |
Return the mapped vertices of the cell. For the current class, this function does not use the support points from the geometry of the current cell but instead evaluates an externally given displacement field in addition to the geometry of the cell.
Reimplemented from Mapping< dim, spacedim >.
Definition at line 54 of file mapping_q1_eulerian.cc.
|
overridevirtual |
Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.
Reimplemented from MappingQGeneric< dim, spacedim >.
Definition at line 122 of file mapping_q1_eulerian.cc.
|
overridevirtual |
Always returns false
because MappingQ1Eulerian does not in general preserve vertex locations (unless the translation vector happens to provide for zero displacements at vertex locations).
Reimplemented from MappingQGeneric< dim, spacedim >.
|
overrideprotectedvirtual |
Compute mapping-related information for a cell. See the documentation of Mapping::fill_fe_values() for a discussion of purpose, arguments, and return value of this function.
This function overrides the function in the base class since we cannot use any cell similarity for this class.
Definition at line 132 of file mapping_q1_eulerian.cc.
|
overrideprotectedvirtual |
Compute the support points of the mapping. For the current class, these are the vertices, as obtained by calling Mapping::get_vertices(). See the documentation of MappingQGeneric::compute_mapping_support_points() for more information.
Reimplemented from MappingQGeneric< dim, spacedim >.
Definition at line 105 of file mapping_q1_eulerian.cc.
|
protected |
Reference to the vector of shifts.
Definition at line 180 of file mapping_q1_eulerian.h.
|
protected |
Pointer to the DoFHandler to which the mapping vector is associated.
Definition at line 187 of file mapping_q1_eulerian.h.