![]() |
Reference documentation for deal.II version 8.1.0
|
#include <mapping_q1_eulerian.h>
Public Member Functions | |
MappingQ1Eulerian (const VECTOR &euler_transform_vectors, const DoFHandler< dim, spacedim > &shiftmap_dof_handler) | |
virtual Mapping< dim, spacedim > * | clone () const |
bool | preserves_vertex_locations () const |
DeclException0 (ExcInactiveCell) | |
![]() | |
MappingQ1 () | |
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 |
virtual void | transform (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, VectorSlice< std::vector< Tensor< 1, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
virtual void | transform (const VectorSlice< const std::vector< DerivativeForm< 1, dim, spacedim > > > input, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
virtual void | transform (const VectorSlice< const std::vector< Tensor< 2, dim > > > input, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) 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, typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, typename std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, typename std::vector< Tensor< 1, spacedim > > &boundary_form, typename std::vector< Point< spacedim > > &normal_vectors) const |
virtual void | fill_fe_subface_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, typename std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, typename std::vector< Tensor< 1, spacedim > > &boundary_form, typename std::vector< Point< spacedim > > &normal_vectors) const |
void | compute_shapes (const std::vector< Point< dim > > &unit_points, InternalData &data) const |
void | compute_data (const UpdateFlags flags, const Quadrature< dim > &quadrature, const unsigned int n_orig_q_points, InternalData &data) const |
void | compute_face_data (const UpdateFlags flags, const Quadrature< dim > &quadrature, const unsigned int n_orig_q_points, InternalData &data) const |
void | compute_fill (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int npts, const DataSetDescriptor data_set, const CellSimilarity::Similarity cell_similarity, InternalData &data, std::vector< Point< spacedim > > &quadrature_points) const |
void | compute_fill_face (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const unsigned int npts, const DataSetDescriptor data_set, const std::vector< double > &weights, InternalData &mapping_data, std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< Tensor< 1, spacedim > > &boundary_form, std::vector< Point< spacedim > > &normal_vectors) const |
virtual void | compute_shapes_virtual (const std::vector< Point< dim > > &unit_points, InternalData &data) const |
Point< spacedim > | transform_unit_to_real_cell_internal (const InternalData &mdata) const |
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, InternalData &mdata) const |
template<> | |
Point< 2 > | transform_real_to_unit_cell_internal (const Triangulation< 2, 3 >::cell_iterator &cell, const Point< 3 > &p, const Point< 2 > &initial_p_unit, InternalData &mdata) const |
template<> | |
Point< 1 > | transform_real_to_unit_cell_internal (const Triangulation< 1, 2 >::cell_iterator &cell, const Point< 2 > &p, const Point< 1 > &initial_p_unit, InternalData &mdata) const |
template<> | |
Point< 1 > | transform_real_to_unit_cell_internal (const Triangulation< 1, 3 >::cell_iterator &cell, const Point< 3 > &p, const Point< 1 > &initial_p_unit, InternalData &mdata) const |
![]() | |
virtual | ~Mapping () |
virtual void | transform (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, VectorSlice< std::vector< Tensor< 1, spacedim > > > output, const InternalDataBase &internal, const MappingType type) const =0 |
virtual void | transform (const VectorSlice< const std::vector< DerivativeForm< 1, dim, spacedim > > > input, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const InternalDataBase &internal, const MappingType type) const =0 |
virtual void | transform (const VectorSlice< const std::vector< Tensor< 2, dim > > > input, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const InternalDataBase &internal, const MappingType type) const =0 |
void | transform_covariant (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 1, spacedim > > > output, const InternalDataBase &internal) const DEAL_II_DEPRECATED |
void | transform_covariant (const VectorSlice< const std::vector< DerivativeForm< 1, dim, spacedim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const InternalDataBase &internal) const DEAL_II_DEPRECATED |
void | transform_contravariant (const VectorSlice< const std::vector< Tensor< 1, dim > > > input, const unsigned int offset, VectorSlice< std::vector< Tensor< 1, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const DEAL_II_DEPRECATED |
void | transform_contravariant (const VectorSlice< const std::vector< DerivativeForm< 1, dim, spacedim > > > input, const unsigned int offset, const VectorSlice< std::vector< Tensor< 2, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const DEAL_II_DEPRECATED |
const Point< spacedim > & | support_point_value (const unsigned int index, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const |
const Tensor< 2, spacedim > & | support_point_gradient (const unsigned int index, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const |
const Tensor< 2, spacedim > & | support_point_inverse_gradient (const unsigned int index, const typename Mapping< dim, spacedim >::InternalDataBase &internal) const |
DeclException0 (ExcInvalidData) | |
DeclException0 (ExcTransformationFailed) | |
DeclException3 (ExcDistortedMappedCell, Point< spacedim >, double, int,<< "The image of the mapping applied to cell with center ["<< arg1<< "] is distorted. The cell geometry or the "<< "mapping are invalid, giving a non-positive volume "<< "fraction of "<< arg2<< " in quadrature point "<< arg3<< ".") | |
![]() | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
DeclException3 (ExcInUse, int, char *, std::string &,<< "Object of class "<< arg2<< " is still used by "<< arg1<< " other objects.\n"<< "(Additional information: "<< arg3<< ")\n"<< "Note the entry in the Frequently Asked Questions of "<< "deal.II (linked to from http://www.dealii.org/) for "<< "more information on what this error means.") | |
DeclException2 (ExcNoSubscriber, char *, char *,<< "No subscriber with identifier \""<< arg2<< "\" did subscribe to this object of class "<< arg1) | |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Protected Member Functions | |
virtual void | fill_fe_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Quadrature< dim > &quadrature, typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, typename std::vector< Point< spacedim > > &quadrature_points, std::vector< double > &JxW_values, std::vector< DerivativeForm< 1, dim, spacedim > > &jacobians, std::vector< DerivativeForm< 2, dim, spacedim > > &jacobian_grads, std::vector< DerivativeForm< 1, spacedim, dim > > &inverse_jacobians, std::vector< Point< spacedim > > &cell_normal_vectors, CellSimilarity::Similarity &cell_similarity) const |
![]() | |
template<int rank> | |
void | transform_fields (const VectorSlice< const std::vector< Tensor< rank, dim > > > input, VectorSlice< std::vector< Tensor< rank, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
template<int rank> | |
void | transform_gradients (const VectorSlice< const std::vector< Tensor< rank, dim > > > input, VectorSlice< std::vector< Tensor< rank, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
template<int rank> | |
void | transform_differential_forms (const VectorSlice< const std::vector< DerivativeForm< rank, dim, spacedim > > > input, VectorSlice< std::vector< DerivativeForm< rank, spacedim, spacedim > > > output, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const MappingType type) const |
template<int dim_> | |
Point< dim_ > | transform_real_to_unit_cell_internal_codim1 (const typename Triangulation< dim_, dim_+1 >::cell_iterator &cell, const Point< dim_+1 > &p, const Point< dim_ > &initial_p_unit, InternalData &mdata) const |
Point< dim > | transform_real_to_unit_cell_initial_guess (const std::vector< Point< spacedim > > &vertex, const Point< spacedim > &p) const |
Protected Attributes | |
SmartPointer< const VECTOR, MappingQ1Eulerian< dim, VECTOR, spacedim > > | euler_transform_vectors |
SmartPointer< const DoFHandler < dim, spacedim > , MappingQ1Eulerian< dim, VECTOR, spacedim > > | shiftmap_dof_handler |
Private Member Functions | |
virtual void | compute_mapping_support_points (const typename Triangulation< dim, spacedim >::cell_iterator &cell, std::vector< Point< spacedim > > &a) const |
Additional Inherited Members | |
![]() | |
typedef QProjector< dim > ::DataSetDescriptor | DataSetDescriptor |
Eulerian mapping of general unit cells by d-linear shape functions. 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 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 wrapper classes, the type of the vector can be specified as template parameter EulerVectorType
Not specifying this template argument in applications using the PETSc vector classes leads to the construction of a copy of the vector which is not acccessible afterwards!
For more information about the spacedim
template parameter check the documentation of FiniteElement or the one of Triangulation.
Definition at line 92 of file mapping_q1_eulerian.h.
MappingQ1Eulerian< dim, VECTOR, spacedim >::MappingQ1Eulerian | ( | const VECTOR & | euler_transform_vectors, |
const DoFHandler< dim, spacedim > & | shiftmap_dof_handler | ||
) |
Constructor. It takes a Vector<double> &
as its first argument to specify the transformation of the whole problem from the reference to the current configuration. The organization of the elements in the Vector
must follow the concept how deal.II stores solutions that are associated to a triangulation. This is automatically the case if the Vector
represents the solution of the previous step of a nonlinear problem. Alternatively, the Vector
can be initialized by DoFAccessor::set_dof_values()
.
|
virtual |
Return a pointer to a copy of the present object. The caller of this copy then assumes ownership of it.
Reimplemented from MappingQ1< dim, spacedim >.
|
virtual |
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 MappingQ1< dim, spacedim >.
MappingQ1Eulerian< dim, VECTOR, spacedim >::DeclException0 | ( | ExcInactiveCell | ) |
Exception.
|
protectedvirtual |
Implementation of the interface in MappingQ1. Overrides the function in the base class, since we cannot use any cell similarity for this class.
Reimplemented from MappingQ1< dim, spacedim >.
|
privatevirtual |
Computes the support points of the mapping. For MappingQ1Eulerian
these are the vertices.
Reimplemented from MappingQ1< dim, spacedim >.
|
protected |
Reference to the vector of shifts.
Definition at line 168 of file mapping_q1_eulerian.h.
|
protected |
Pointer to the DoFHandler to which the mapping vector is associated.
Definition at line 175 of file mapping_q1_eulerian.h.