![]() |
Reference documentation for deal.II version 8.1.0
|
#include <mapping_q_eulerian.h>
Classes | |
class | SupportQuadrature |
Public Member Functions | |
MappingQEulerian (const unsigned int degree, const VECTOR &euler_vector, const DoFHandler< dim, spacedim > &euler_dof_handler) | |
virtual Mapping< dim, spacedim > * | clone () const |
bool | preserves_vertex_locations () const |
DeclException0 (ExcInactiveCell) | |
![]() | |
MappingQ (const unsigned int p, const bool use_mapping_q_on_all_cells=false) | |
MappingQ (const MappingQ< dim, spacedim > &mapping) | |
virtual | ~MappingQ () |
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 |
unsigned int | get_degree () const |
![]() | |
MappingQ1 () | |
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 |
![]() | |
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 > > &exterior_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 > > &exterior_form, typename std::vector< Point< spacedim > > &normal_vectors) 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 |
![]() | |
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, MappingQEulerian< dim, VECTOR, spacedim > > | euler_vector |
SmartPointer< const DoFHandler< dim, spacedim >, MappingQEulerian< dim, VECTOR, spacedim > > | euler_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 |
Private Attributes | |
const SupportQuadrature | support_quadrature |
FEValues< dim, spacedim > | fe_values |
Threads::Mutex | fe_values_mutex |
Additional Inherited Members | |
![]() | |
typedef QProjector< dim >::DataSetDescriptor | DataSetDescriptor |
This class is an extension of the MappingQ1Eulerian class to higher order Qp mappings. It is useful when one wants to calculate shape function information on a domain that is deforming as the computation proceeds.
The constructor of this class takes three arguments: the polynomial degree of the desire 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 reqirement is that the number of components of the given vector field 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:
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 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!
Definition at line 94 of file mapping_q_eulerian.h.
MappingQEulerian< dim, VECTOR, spacedim >::MappingQEulerian | ( | const unsigned int | degree, |
const VECTOR & | euler_vector, | ||
const DoFHandler< dim, spacedim > & | euler_dof_handler | ||
) |
Constructor. The first argument is the polynomical degree of the desired Qp mapping. It then takes a Vector<double> &
to specify the transformation of the domain 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 MappingQ< 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 >.
MappingQEulerian< dim, VECTOR, spacedim >::DeclException0 | ( | ExcInactiveCell | ) |
Exception
|
protectedvirtual |
Implementation of the interface in MappingQ. Overrides the function in the base class, since we cannot use any cell similarity for this class.
Reimplemented from MappingQ< dim, spacedim >.
|
privatevirtual |
Compute the positions of the support points in the current configuration
Reimplemented from MappingQ< dim, spacedim >.
|
protected |
Reference to the vector of shifts.
Definition at line 153 of file mapping_q_eulerian.h.
|
protected |
Pointer to the DoFHandler to which the mapping vector is associated.
Definition at line 158 of file mapping_q_eulerian.h.
|
private |
A member variable holding the quadrature points in the right order.
Definition at line 182 of file mapping_q_eulerian.h.
|
mutableprivate |
FEValues object used to query the the given finite element field at the support points in the reference configuration.
The variable is marked as mutable since we have to call FEValues::reinit from compute_mapping_support_points, a function that is 'const'.
Definition at line 191 of file mapping_q_eulerian.h.
|
mutableprivate |
A variable to guard access to the fe_values variable.
Definition at line 196 of file mapping_q_eulerian.h.