MappingQEulerian< dim, EulerVectorType, spacedim > Class Template Reference
[Mappings between reference and real cell]

Inheritance diagram for MappingQEulerian< dim, EulerVectorType, spacedim >:

Inheritance graph
[legend]

List of all members.

Classes

class  ExcInactiveCell
class  ExcWrongNoOfComponents
class  ExcWrongVectorSize
class  SupportQuadrature

Public Member Functions

 MappingQEulerian (const unsigned int degree, const EulerVectorType &euler_vector, const DoFHandler< dim > &euler_dof_handler)
virtual Mapping< dim, spacedim > * clone () const

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< Tensor< 2, spacedim > > &jacobians, std::vector< Tensor< 3, spacedim > > &jacobian_grads, std::vector< Tensor< 2, spacedim > > &inverse_jacobians, std::vector< Point< spacedim > > &cell_normal_vectors, enum CellSimilarity::Similarity &cell_similarity) const

Protected Attributes

const EulerVectorType & euler_vector
const SmartPointer< const
DoFHandler< dim > > 
euler_dof_handler

Private Member Functions

virtual void compute_mapping_support_points (const typename Triangulation< dim >::cell_iterator &cell, std::vector< Point< dim > > &a) const

Private Attributes

const SupportQuadrature support_quadrature
FEValues< dim > fe_values
Threads::ThreadMutex fe_values_mutex


Detailed Description

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

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.

Usage

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:

 *    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> soln_vector(dof_handler.n_dofs());
 *    MappingQEulerian<dim> q2_mapping(2,soln_vector,dof_handler);
 * 

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!

Author:
Joshua White, 2008

Constructor & Destructor Documentation

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
MappingQEulerian< dim, EulerVectorType, spacedim >::MappingQEulerian ( const unsigned int  degree,
const EulerVectorType &  euler_vector,
const DoFHandler< dim > &  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().


Member Function Documentation

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
virtual Mapping<dim,spacedim>* MappingQEulerian< dim, EulerVectorType, spacedim >::clone (  )  const [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 >.

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
virtual void MappingQEulerian< dim, EulerVectorType, spacedim >::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< Tensor< 2, spacedim > > &  jacobians,
std::vector< Tensor< 3, spacedim > > &  jacobian_grads,
std::vector< Tensor< 2, spacedim > > &  inverse_jacobians,
std::vector< Point< spacedim > > &  cell_normal_vectors,
enum CellSimilarity::Similarity cell_similarity 
) const [protected, virtual]

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

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
virtual void MappingQEulerian< dim, EulerVectorType, spacedim >::compute_mapping_support_points ( const typename Triangulation< dim >::cell_iterator &  cell,
std::vector< Point< dim > > &  a 
) const [private, virtual]

Compute the positions of the support points in the current configuration


Member Data Documentation

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
const EulerVectorType& MappingQEulerian< dim, EulerVectorType, spacedim >::euler_vector [protected]

Reference to the vector of shifts.

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
const SmartPointer<const DoFHandler<dim> > MappingQEulerian< dim, EulerVectorType, spacedim >::euler_dof_handler [protected]

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

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
const SupportQuadrature MappingQEulerian< dim, EulerVectorType, spacedim >::support_quadrature [private]

A member variable holding the quadrature points in the right order.

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
FEValues<dim> MappingQEulerian< dim, EulerVectorType, spacedim >::fe_values [mutable, private]

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'.

template<int dim, class EulerVectorType = Vector<double>, int spacedim = dim>
Threads::ThreadMutex MappingQEulerian< dim, EulerVectorType, spacedim >::fe_values_mutex [mutable, private]

A variable to guard access to the fe_values variable.


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

deal.II documentation generated on Sat Aug 15 16:52:06 2009 by doxygen 1.5.9