FE_ABF< dim > Class Template Reference
[Finite element space descriptions]

Inheritance diagram for FE_ABF< dim >:

Inheritance graph
[legend]

List of all members.

Classes

class  InternalData

Public Member Functions

 FE_ABF (const unsigned int p)
virtual std::string get_name () const
virtual unsigned int n_base_elements () const
virtual const FiniteElement
< dim > & 
base_element (const unsigned int index) const
virtual unsigned int element_multiplicity (const unsigned int index) const
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const
virtual void interpolate (std::vector< double > &local_dofs, const std::vector< double > &values) const
virtual void interpolate (std::vector< double > &local_dofs, const std::vector< Vector< double > > &values, unsigned int offset=0) const
virtual void interpolate (std::vector< double > &local_dofs, const VectorSlice< const std::vector< std::vector< double > > > &values) const
virtual unsigned int memory_consumption () const
virtual FiniteElement< dim > * clone () const

Private Member Functions

void initialize_support_points (const unsigned int rt_degree)
void initialize_restriction ()
virtual UpdateFlags update_once (const UpdateFlags flags) const
virtual UpdateFlags update_each (const UpdateFlags flags) const

Static Private Member Functions

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

Private Attributes

const unsigned int rt_order
Table< 2, doubleboundary_weights
Table< 3, doubleinterior_weights
Table< 2, doubleboundary_weights_abf
Table< 3, doubleinterior_weights_abf

Friends

class FE_ABF


Detailed Description

template<int dim>
class FE_ABF< dim >

Implementation of Arnold-Boffi-Falk (ABF) elements, conforming with the space Hdiv. These elements generate vector fields with normal components continuous between mesh cells.

These elements are based on an article from Arnold, Boffi and Falk: Quadrilateral H(div) finite elements, SIAM J. Numer. Anal. Vol.42, No.6, pp.2429-2451

In this article, the authors demonstrate that the usual RT elements and also BDM and other proposed finite dimensional subspaces of H(div) do not work properly on arbitrary FE grids. I.e. the convergence rates deteriorate on these meshes. As a solution the authors propose the ABF elements, which are implemented in this module.

This class is not implemented for the codimension one case (spacedim != dim).

Todo:
Even if this element is implemented for two and three space dimensions, the definition of the node values relies on consistently oriented faces in 3D. Therefore, care should be taken on complicated meshes.

Interpolation

The interpolation operators associated with the RT element are constructed such that interpolation and computing the divergence are commuting operations. We require this from interpolating arbitrary functions as well as the restriction matrices. It can be achieved by two interpolation schemes, the simplified one in FE_RaviartThomasNodal and the original one here:

Node values on edges/faces

On edges or faces, the node values are the moments of the normal component of the interpolated function with respect to the traces of the RT polynomials. Since the normal trace of the RT space of degree k on an edge/face is the space Qk, the moments are taken with respect to this space.

Interior node values

Higher order RT spaces have interior nodes. These are moments taken with respect to the gradient of functions in Qk on the cell (this space is the matching space for RTk in a mixed formulation).

Generalized support points

The node values above rely on integrals, which will be computed by quadrature rules themselves. The generalized support points are a set of points such that this quadrature can be performed with sufficient accuracy. The points needed are those of QGaussk+1 on each face as well as QGaussk in the interior of the cell (or none for RT0).

Author:
Oliver Kayser-Herold, 2006, based on previous work by Guido Kanschat and Wolfgang Bangerth

Constructor & Destructor Documentation

template<int dim>
FE_ABF< dim >::FE_ABF ( const unsigned int  p  ) 

Constructor for the ABF element of degree p.


Member Function Documentation

template<int dim>
virtual std::string FE_ABF< dim >::get_name (  )  const [virtual]

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

Implements FiniteElement< dim, dim >.

template<int dim>
virtual unsigned int FE_ABF< dim >::n_base_elements (  )  const [virtual]

Number of base elements in a mixed discretization. Here, this is of course equal to one.

Reimplemented from FE_PolyTensor< PolynomialsABF< dim >, dim >.

template<int dim>
virtual const FiniteElement<dim>& FE_ABF< dim >::base_element ( const unsigned int  index  )  const [virtual]

Access to base element objects. Since this element is atomic, base_element(0) is this, and all other indices throw an error.

Reimplemented from FE_PolyTensor< PolynomialsABF< dim >, dim >.

template<int dim>
virtual unsigned int FE_ABF< dim >::element_multiplicity ( const unsigned int  index  )  const [virtual]

Multiplicity of base element index. Since this is an atomic element, element_multiplicity(0) returns one, and all other indices will throw an error.

Reimplemented from FE_PolyTensor< PolynomialsABF< dim >, dim >.

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

Check whether a shape function may be non-zero on a face.

Right now, this is only implemented for RT0 in 1D. Otherwise, returns always true.

Reimplemented from FiniteElement< dim, dim >.

template<int dim>
virtual void FE_ABF< dim >::interpolate ( std::vector< double > &  local_dofs,
const std::vector< double > &  values 
) const [virtual]

Interpolate a set of scalar values, computed in the generalized support points.

Note:
This function is implemented in FiniteElement for the case that the element has support points. In this case, the resulting coefficients are just the values in the suport points. All other elements must reimplement it.

Reimplemented from FiniteElement< dim, dim >.

template<int dim>
virtual void FE_ABF< dim >::interpolate ( std::vector< double > &  local_dofs,
const std::vector< Vector< double > > &  values,
unsigned int  offset = 0 
) const [virtual]

Interpolate a set of vector values, computed in the generalized support points.

Since a finite element often only interpolates part of a vector, offset is used to determine the first component of the vector to be interpolated. Maybe consider changing your data structures to use the next function.

Reimplemented from FiniteElement< dim, dim >.

template<int dim>
virtual void FE_ABF< dim >::interpolate ( std::vector< double > &  local_dofs,
const VectorSlice< const std::vector< std::vector< double > > > &  values 
) const [virtual]

Interpolate a set of vector values, computed in the generalized support points.

Reimplemented from FiniteElement< dim, dim >.

template<int dim>
virtual unsigned int FE_ABF< dim >::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, dim >.

template<int dim>
virtual FiniteElement<dim>* FE_ABF< dim >::clone (  )  const [virtual]

A sort of virtual copy constructor. Some places in the library, for example the constructors of FESystem as well as the hp::FECollection class, need to make copied of finite elements without knowing their exact type. They do so through this function.

Implements FiniteElement< dim, dim >.

template<int dim>
static std::vector<unsigned int> FE_ABF< dim >::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>
void FE_ABF< dim >::initialize_support_points ( const unsigned int  rt_degree  )  [private]

Initialize the generalized_support_points field of the FiniteElement class and fill the tables with interpolation weights (boundary_weights and interior_weights). Called from the constructor.

template<int dim>
void FE_ABF< dim >::initialize_restriction (  )  [private]

Initialize the interpolation from functions on refined mesh cells onto the father cell. According to the philosophy of the Raviart-Thomas element, this restriction operator preserves the divergence of a function weakly.

template<int dim>
virtual UpdateFlags FE_ABF< dim >::update_once ( const UpdateFlags  flags  )  const [private, virtual]

Given a set of flags indicating what quantities are requested from a FEValues object, return which of these can be precomputed once and for all. Often, the values of shape function at quadrature points can be precomputed, for example, in which case the return value of this function would be the logical and of the input flags and update_values.

For the present kind of finite element, this is exactly the case.

Reimplemented from FE_PolyTensor< PolynomialsABF< dim >, dim >.

template<int dim>
virtual UpdateFlags FE_ABF< dim >::update_each ( const UpdateFlags  flags  )  const [private, virtual]

This is the opposite to the above function: given a set of flags indicating what we want to know, return which of these need to be computed each time we visit a new cell.

If for the computation of one quantity something else is also required (for example, we often need the covariant transformation when gradients need to be computed), include this in the result as well.

Reimplemented from FE_PolyTensor< PolynomialsABF< dim >, dim >.


Friends And Related Function Documentation

template<int dim>
friend class FE_ABF [friend]

Allow access from other dimensions.


Member Data Documentation

template<int dim>
const unsigned int FE_ABF< dim >::rt_order [private]

The order of the ABF element. The lowest order elements are usually referred to as RT0, even though their shape functions are piecewise quadratics.

template<int dim>
Table<2, double> FE_ABF< dim >::boundary_weights [private]

These are the factors multiplied to a function in the generalized_face_support_points when computing the integration. They are organized such that there is one row for each generalized face support point and one column for each degree of freedom on the face.

template<int dim>
Table<3, double> FE_ABF< dim >::interior_weights [private]

Precomputed factors for interpolation of interior degrees of freedom. The rationale for this Table is the same as for boundary_weights. Only, this table has a third coordinate for the space direction of the component evaluated.

template<int dim>
Table<2, double> FE_ABF< dim >::boundary_weights_abf [private]

These are the factors multiplied to a function in the generalized_face_support_points when computing the integration. They are organized such that there is one row for each generalized face support point and one column for each degree of freedom on the face.

template<int dim>
Table<3, double> FE_ABF< dim >::interior_weights_abf [private]

Precomputed factors for interpolation of interior degrees of freedom. The rationale for this Table is the same as for boundary_weights. Only, this table has a third coordinate for the space direction of the component evaluated.


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

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