Classes | |
struct | Matrices |
Public Member Functions | |
FE_DGP (const unsigned int p) | |
virtual std::string | get_name () 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 | |
virtual FiniteElement< dim, spacedim > * | clone () const |
Static Private Member Functions | |
static std::vector< unsigned int > | get_dpo_vector (const unsigned int degree) |
This finite element implements complete polynomial spaces, that is, dim-dimensional polynomials of degree p. For example, in 2d the element FE_DGP(1) would represent the span of the functions , which is in contrast to the element FE_DGQ(1) that is formed by the span of
. Since the DGP space has only three unknowns for each quadrilateral, it is immediately clear that this element can not be continuous.
The basis functions for this element are chosen to form a Legendre basis on the unit square. Thus, the mass matrix is diagonal, if the grid cells are parallelograms. Note that this is in contrast to the FE_DGPMonomial class that actually uses the monomial basis listed above as basis functions.
The shape functions are defined in the class PolynomialSpace. The polynomials used inside PolynomialSpace are Polynomials::Legendre up to degree p
given in FE_DGP. For the ordering of the basis functions, refer to PolynomialSpace, remebering that the Legendre polynomials are ordered by ascending degree.
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.
It is worth noting that under a (bi-, tri-)linear mapping, the space described by this element does not contain , even if we use a basis of polynomials of degree
. Consequently, for example, on meshes with non-affine cells, a linear function can not be exactly represented by elements of type FE_DGP(1) or FE_DGPMonomial(1).
This can be understood by the following 2-d example: consider the cell with vertices at :
For this cell, a bilinear transformation produces the relations
and
that correlate reference coordinates
and coordinates in real space
. Under this mapping, the constant function is clearly mapped onto itself, but the two other shape functions of the
space, namely
and
are mapped onto
where
.
For the simple case that , i.e. if the real cell is the unit square, the expressions can be simplified to
and
. However, for all other cases, the functions
are not linear any more, and neither is any linear combincation of them. Consequently, the linear functions are not within the range of the mapped
polynomials.
Constructor for tensor product polynomials of degree p
.
virtual std::string FE_DGP< dim, spacedim >::get_name | ( | ) | const [virtual] |
Return a string that uniquely identifies a finite element. This class returns FE_DGP<dim>(degree)
, with dim
and degree
replaced by appropriate values.
Implements FiniteElement< dim, spacedim >.
virtual std::vector<std::pair<unsigned int, unsigned int> > FE_DGP< 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 >.
virtual std::vector<std::pair<unsigned int, unsigned int> > FE_DGP< 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 >.
virtual std::vector<std::pair<unsigned int, unsigned int> > FE_DGP< 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 >.
virtual bool FE_DGP< 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_DGP 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 >.
virtual FiniteElementDomination::Domination FE_DGP< 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 >.
virtual void FE_DGP< 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,spacedim>::ExcInterpolationNotImplemented.
Reimplemented from FiniteElement< dim, spacedim >.
virtual void FE_DGP< 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,spacedim>::ExcInterpolationNotImplemented.
Reimplemented from FiniteElement< dim, spacedim >.
virtual bool FE_DGP< 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 >.
virtual unsigned int FE_DGP< 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 >.
virtual FiniteElement<dim,spacedim>* FE_DGP< 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 >.
static std::vector<unsigned int> FE_DGP< 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
.