FE_DGQ< dim, spacedim > Class Template Reference
[Finite element space descriptions]

Inheritance diagram for FE_DGQ< dim, spacedim >:

Inheritance graph
[legend]

List of all members.

Public Member Functions

 FE_DGQ (const unsigned int p)
virtual std::string get_name () const
virtual void get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual void get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual void get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const
virtual unsigned int memory_consumption () const
Functions to support hp
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_vertex_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_line_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_quad_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const
virtual bool hp_constraints_are_implemented () const
virtual
FiniteElementDomination::Domination 
compare_for_face_domination (const FiniteElement< dim, spacedim > &fe_other) const

Protected Member Functions

 FE_DGQ (const Quadrature< 1 > &points)
virtual FiniteElement< dim,
spacedim > * 
clone () const

Private Member Functions

void rotate_indices (std::vector< unsigned int > &indices, const char direction) const

Static Private Member Functions

static std::vector< unsigned intget_dpo_vector (const unsigned int degree)

Friends

class FE_DGQ
class MappingQ


Detailed Description

template<int dim, int spacedim = dim>
class FE_DGQ< dim, spacedim >

Implementation of scalar, discontinuous tensor product elements based on equidistant support points.

This is a discontinuous finite element based on tensor products of Lagrangian polynomials. The shape functions are Lagrangian interpolants of an equidistant grid of points on the unit cell. The points are numbered in lexicographical order, with x running fastest, then y, then z (if these coordinates are present for a given space dimension at all). For example, these are the node orderings for FE_DGQ(1) in 3d:

 *         6-------7        6-------7
 *        /|       |       /       /|
 *       / |       |      /       / |
 *      /  |       |     /       /  |
 *     4   |       |    4-------5   |
 *     |   2-------3    |       |   3
 *     |  /       /     |       |  /
 *     | /       /      |       | /
 *     |/       /       |       |/
 *     0-------1        0-------1
 *  
and FE_DGQ(2):
 *         24--25--26       24--25--26
 *        /|       |       /       /|
 *      21 |       |     21  22  23 |
 *      /  15  16  17    /       /  17
 *    18   |       |   18--19--20   |
 *     |12 6---7---8    |       |14 8
 *     9  /       /     9  10  11  /
 *     | 3   4   5      |       | 5
 *     |/       /       |       |/
 *     0---1---2        0---1---2
 *  
with node 13 being placed in the interior of the hex.

Note, however, that these are just the Lagrange interpolation points of the shape functions. Even though they may physically be on the surface of the cell, they are logically in the interior since there are no continuity requirements for these shape functions across cell boundaries. This class if partially implemented for the codimension one case (spacedim != dim ), since no passage of information between meshes of different refinement level is possible because the embedding and projection matrices are not computed in the class constructor.

Author:
Ralf Hartmann, Guido Kanschat 2001, 2004

Constructor & Destructor Documentation

template<int dim, int spacedim = dim>
FE_DGQ< dim, spacedim >::FE_DGQ ( const unsigned int  p  ) 

Constructor for tensor product polynomials of degree p. The shape functions created using this constructor correspond to Lagrange interpolation polynomials for equidistantly spaced support points in each coordinate direction.

template<int dim, int spacedim = dim>
FE_DGQ< dim, spacedim >::FE_DGQ ( const Quadrature< 1 > &  points  )  [protected]

Constructor for tensor product polynomials based on Polynomials::Lagrange interpolation of the support points in the quadrature rule points. The degree of these polynomials is points.size()-1.

Note: The FE_DGQ::clone function does not work properly for FE with arbitrary nodes!


Member Function Documentation

template<int dim, int spacedim = dim>
virtual std::string FE_DGQ< dim, spacedim >::get_name (  )  const [virtual]

Return a string that uniquely identifies a finite element. This class returns FE_DGQ<dim>(degree) , with dim and degree replaced by appropriate values.

Implements FiniteElement< dim, spacedim >.

Reimplemented in FE_DGQArbitraryNodes< dim >.

template<int dim, int spacedim = dim>
virtual void FE_DGQ< dim, spacedim >::get_interpolation_matrix ( const FiniteElement< dim, spacedim > &  source,
FullMatrix< double > &  matrix 
) const [virtual]

Return the matrix interpolating from the given finite element to the present one. The size of the matrix is then dofs_per_cell times source.dofs_per_cell.

These matrices are only available if the source element is also a FE_DGQ element. Otherwise, an exception of type FiniteElement<dim>::ExcInterpolationNotImplemented is thrown.

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual void FE_DGQ< dim, spacedim >::get_face_interpolation_matrix ( const FiniteElement< dim, spacedim > &  source,
FullMatrix< double > &  matrix 
) const [virtual]

Return the matrix interpolating from a face of of one element to the face of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.

Derived elements will have to implement this function. They may only provide interpolation matrices for certain source finite elements, for example those from the same family. If they don't implement interpolation from a given element, then they must throw an exception of type FiniteElement<dim>::ExcInterpolationNotImplemented.

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual void FE_DGQ< dim, spacedim >::get_subface_interpolation_matrix ( const FiniteElement< dim, spacedim > &  source,
const unsigned int  subface,
FullMatrix< double > &  matrix 
) const [virtual]

Return the matrix interpolating from a face of of one element to the face of the neighboring element. The size of the matrix is then source.dofs_per_face times this->dofs_per_face.

Derived elements will have to implement this function. They may only provide interpolation matrices for certain source finite elements, for example those from the same family. If they don't implement interpolation from a given element, then they must throw an exception of type FiniteElement<dim>::ExcInterpolationNotImplemented.

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual std::vector<std::pair<unsigned int, unsigned int> > FE_DGQ< dim, spacedim >::hp_vertex_dof_identities ( const FiniteElement< dim, spacedim > &  fe_other  )  const [virtual]

If, on a vertex, several finite elements are active, the hp code first assigns the degrees of freedom of each of these FEs different global indices. It then calls this function to find out which of them should get identical values, and consequently can receive the same global DoF index. This function therefore returns a list of identities between DoFs of the present finite element object with the DoFs of fe_other, which is a reference to a finite element object representing one of the other finite elements active on this particular vertex. The function computes which of the degrees of freedom of the two finite element objects are equivalent, and returns a list of pairs of global dof indices in identities. The first index of each pair denotes one of the vertex dofs of the present element, whereas the second is the corresponding index of the other finite element.

This being a discontinuous element, the set of such constraints is of course empty.

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual std::vector<std::pair<unsigned int, unsigned int> > FE_DGQ< dim, spacedim >::hp_line_dof_identities ( const FiniteElement< dim, spacedim > &  fe_other  )  const [virtual]

Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on lines.

This being a discontinuous element, the set of such constraints is of course empty.

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual std::vector<std::pair<unsigned int, unsigned int> > FE_DGQ< dim, spacedim >::hp_quad_dof_identities ( const FiniteElement< dim, spacedim > &  fe_other  )  const [virtual]

Same as hp_vertex_dof_indices(), except that the function treats degrees of freedom on quads.

This being a discontinuous element, the set of such constraints is of course empty.

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual bool FE_DGQ< dim, spacedim >::hp_constraints_are_implemented (  )  const [virtual]

Return whether this element implements its hanging node constraints in the new way, which has to be used to make elements "hp compatible".

For the FE_DGQ class the result is always true (independent of the degree of the element), as it has no hanging nodes (being a discontinuous element).

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual FiniteElementDomination::Domination FE_DGQ< dim, spacedim >::compare_for_face_domination ( const FiniteElement< dim, spacedim > &  fe_other  )  const [virtual]

Return whether this element dominates the one given as argument when they meet at a common face, whether it is the other way around, whether neither dominates, or if either could dominate.

For a definition of domination, see FiniteElementBase::Domination and in particular the hp paper.

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual bool FE_DGQ< dim, spacedim >::has_support_on_face ( const unsigned int  shape_index,
const unsigned int  face_index 
) const [virtual]

Check for non-zero values on a face.

This function returns true, if the shape function shape_index has non-zero values on the face face_index.

Implementation of the interface in FiniteElement

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual unsigned int FE_DGQ< dim, spacedim >::memory_consumption (  )  const [virtual]

Determine an estimate for the memory consumption (in bytes) of this object.

This function is made virtual, since finite element objects are usually accessed through pointers to their base class, rather than the class itself.

Reimplemented from FiniteElement< dim, spacedim >.

template<int dim, int spacedim = dim>
virtual FiniteElement<dim, spacedim>* FE_DGQ< dim, spacedim >::clone (  )  const [protected, virtual]

clone function instead of a copy constructor.

This function is needed by the constructors of FESystem.

Implements FiniteElement< dim, spacedim >.

Reimplemented in FE_DGQArbitraryNodes< dim >.

template<int dim, int spacedim = dim>
static std::vector<unsigned int> FE_DGQ< dim, spacedim >::get_dpo_vector ( const unsigned int  degree  )  [static, private]

Only for internal use. Its full name is get_dofs_per_object_vector function and it creates the dofs_per_object vector that is needed within the constructor to be passed to the constructor of FiniteElementData.

template<int dim, int spacedim = dim>
void FE_DGQ< dim, spacedim >::rotate_indices ( std::vector< unsigned int > &  indices,
const char  direction 
) const [private]

Compute renumbering for rotation of degrees of freedom.

Rotates a tensor product numbering of degrees of freedom by 90 degrees. It is used to compute the transfer matrices of the children by using only the matrix for the first child.

The direction parameter determines the type of rotation. It is one character of xXyYzZ. The character determines the axis of rotation, case determines the direction. Lower case is counter-clockwise seen in direction of the axis.

Since rotation around the y-axis is not used, it is not implemented either.


Friends And Related Function Documentation

template<int dim, int spacedim = dim>
friend class FE_DGQ [friend]

Allow access from other dimensions.

template<int dim, int spacedim = dim>
friend class MappingQ [friend]

Allows MappingQ class to access to build_renumbering function.


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

deal.II documentation generated on Sat Aug 15 16:51:53 2009 by doxygen 1.5.9