Reference documentation for deal.II version 9.0.0
Classes | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
MappingQEulerian< dim, VectorType, spacedim > Class Template Reference

#include <deal.II/fe/mapping_q_eulerian.h>

Inheritance diagram for MappingQEulerian< dim, VectorType, spacedim >:
[legend]

Classes

class  MappingQEulerianGeneric
 

Public Member Functions

 MappingQEulerian (const unsigned int degree, const DoFHandler< dim, spacedim > &euler_dof_handler, const VectorType &euler_vector, const unsigned int level=numbers::invalid_unsigned_int)
 
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 MappingQ< dim, spacedim >
 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 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
 
- Public Member Functions inherited from Mapping< dim, spacedim >
virtual ~Mapping () override=default
 
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 ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (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)
 

Static Public Member Functions

static ::ExceptionBaseExcInactiveCell ()
 
- Static Public Member Functions inherited from Mapping< dim, spacedim >
static ::ExceptionBaseExcInvalidData ()
 
static ::ExceptionBaseExcTransformationFailed ()
 
static ::ExceptionBaseExcDistortedMappedCell (Point< spacedim > arg1, double arg2, int arg3)
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (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
 
- Protected Member Functions inherited from MappingQ< dim, spacedim >
virtual UpdateFlags requires_update_flags (const UpdateFlags update_flags) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBaseget_data (const UpdateFlags, const Quadrature< dim > &quadrature) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBaseget_face_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) const override
 
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBaseget_subface_data (const UpdateFlags flags, const Quadrature< dim-1 > &quadrature) 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

SmartPointer< const VectorType, MappingQEulerian< dim, VectorType, spacedim > > euler_vector
 
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VectorType, spacedim > > euler_dof_handler
 
- Protected Attributes inherited from MappingQ< dim, spacedim >
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
 

Private Attributes

const unsigned int level
 

Detailed Description

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
class MappingQEulerian< dim, VectorType, spacedim >

This class is an extension of the MappingQ1Eulerian class to higher order \(Q_p\) mappings. It is useful when one wants to calculate shape function information on a domain that is deforming as the computation proceeds.

Usage

The constructor of this class takes three arguments: the polynomial degree of the desired Qp mapping, a reference to the vector that defines the mapping from the initial configuration to the current configuration, and a reference to the DoFHandler. The most common case is to use the solution vector for the problem under consideration as the shift vector. The key requirement is that the number of components of the given vector field must be equal to (or possibly greater than) the number of space dimensions. If there are more components than space dimensions (for example, if one is working with a coupled problem where there are additional solution variables), the first dim components are assumed to represent the displacement field, and the remaining components are ignored. If this assumption does not hold one may need to set up a separate DoFHandler on the triangulation and associate the desired shift vector to it.

Typically, the DoFHandler operates on a finite element that is constructed as a system element (FESystem) from continuous FE_Q objects. An example is shown below:

FESystem<dim> fe(FE_Q<dim>(2), dim, FE_Q<dim>(1), 1);
DoFHandler<dim> dof_handler(triangulation);
dof_handler.distribute_dofs(fe);
Vector<double> displacement_field(dof_handler.n_dofs());
// ... compute displacement field somehow...
MappingQEulerian<dim> q2_mapping(2, dof_handler, displacement_field);

In this example, our element consists of (dim+1) components. Only the first dim components will be used, however, to define the Q2 mapping. The remaining components are ignored.

Note that it is essential to call the distribute_dofs(...) function before constructing a mapping object.

Also 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 MappingQEulerian 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.

Author
Joshua White, 2008

Definition at line 90 of file mapping_q_eulerian.h.

Constructor & Destructor Documentation

◆ MappingQEulerian()

template<int dim, class VectorType , int spacedim>
MappingQEulerian< dim, VectorType, spacedim >::MappingQEulerian ( const unsigned int  degree,
const DoFHandler< dim, spacedim > &  euler_dof_handler,
const VectorType &  euler_vector,
const unsigned int  level = numbers::invalid_unsigned_int 
)

Constructor.

Parameters
[in]degreeThe polynomial degree of the desired \(Q_p\) mapping.
[in]euler_dof_handlerA DoFHandler object that defines a finite element space. This space needs to have at least dim components and the first dim components of the space will be considered displacements relative to the original positions of the cells of the triangulation.
[in]euler_vectorA finite element function in the space defined by the second argument. The first 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.
[in]levelThe multi-grid level at which the mapping will be used. It is mainly used to check if the size of the euler_vector is consistent with the euler_dof_handler .

Definition at line 60 of file mapping_q_eulerian.cc.

Member Function Documentation

◆ get_vertices()

template<int dim, class VectorType , int spacedim>
std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > MappingQEulerian< dim, VectorType, spacedim >::get_vertices ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell) const
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 131 of file mapping_q_eulerian.cc.

◆ clone()

template<int dim, class VectorType , int spacedim>
std::unique_ptr< Mapping< dim, spacedim > > MappingQEulerian< dim, VectorType, spacedim >::clone ( ) const
overridevirtual

Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.

Reimplemented from MappingQ< dim, spacedim >.

Definition at line 84 of file mapping_q_eulerian.cc.

◆ preserves_vertex_locations()

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
virtual bool MappingQEulerian< dim, VectorType, spacedim >::preserves_vertex_locations ( ) const
overridevirtual

Always return false because MappingQEulerian does not in general preserve vertex locations (unless the translation vector happens to provide zero displacements at vertex locations).

Reimplemented from MappingQ< dim, spacedim >.

◆ fill_fe_values()

template<int dim, class VectorType , int spacedim>
CellSimilarity::Similarity MappingQEulerian< dim, VectorType, spacedim >::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
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.

Reimplemented from MappingQ< dim, spacedim >.

Definition at line 246 of file mapping_q_eulerian.cc.

Member Data Documentation

◆ euler_vector

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
SmartPointer<const VectorType, MappingQEulerian<dim,VectorType,spacedim> > MappingQEulerian< dim, VectorType, spacedim >::euler_vector
protected

Reference to the vector of shifts.

Definition at line 165 of file mapping_q_eulerian.h.

◆ euler_dof_handler

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
SmartPointer<const DoFHandler<dim,spacedim>,MappingQEulerian<dim,VectorType,spacedim> > MappingQEulerian< dim, VectorType, spacedim >::euler_dof_handler
protected

Pointer to the DoFHandler to which the mapping vector is associated.

Definition at line 170 of file mapping_q_eulerian.h.

◆ level

template<int dim, typename VectorType = Vector<double>, int spacedim = dim>
const unsigned int MappingQEulerian< dim, VectorType, spacedim >::level
private

Multigrid level at which the mapping is to be used.

Definition at line 178 of file mapping_q_eulerian.h.


The documentation for this class was generated from the following files: