FiniteElement< dim > Class Template Reference
[Base classesFinite element space descriptions]

Inheritance diagram for FiniteElement< dim >:

Inheritance graph
[legend]

List of all members.

Public Member Functions

 FiniteElement (const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< std::vector< bool > > &nonzero_components)
virtual ~FiniteElement ()
virtual std::string get_name () const =0
const FiniteElement< dim > & operator[] (const unsigned int fe_index) const
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_vertex_dof_identities (const FiniteElement< dim > &fe_other) const
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_line_dof_identities (const FiniteElement< dim > &fe_other) const
virtual std::vector< std::pair
< unsigned int, unsigned int > > 
hp_quad_dof_identities (const FiniteElement< dim > &fe_other) const
virtual
FiniteElementDomination::Domination 
compare_for_face_domination (const FiniteElement< dim > &fe_other) const
bool operator== (const FiniteElement< dim > &) const
virtual unsigned int memory_consumption () const
Shape function access
virtual double shape_value (const unsigned int i, const Point< dim > &p) const
virtual double shape_value_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual Tensor< 1, dim > shape_grad (const unsigned int i, const Point< dim > &p) const
virtual Tensor< 1, dim > shape_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual Tensor< 2, dim > shape_grad_grad (const unsigned int i, const Point< dim > &p) const
virtual Tensor< 2, dim > shape_grad_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual bool has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const
Transfer and constraint matrices
const FullMatrix< double > & get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
const FullMatrix< double > & get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
bool prolongation_is_implemented () const
bool isotropic_prolongation_is_implemented () const
bool restriction_is_implemented () const
bool isotropic_restriction_is_implemented () const
bool restriction_is_additive (const unsigned int index) const
const FullMatrix< double > & constraints (const internal::SubfaceCase< dim > &subface_case=internal::SubfaceCase< dim >::case_isotropic) const
bool constraints_are_implemented (const internal::SubfaceCase< dim > &subface_case=internal::SubfaceCase< dim >::case_isotropic) const
virtual bool hp_constraints_are_implemented () const
virtual void get_interpolation_matrix (const FiniteElement< dim > &source, FullMatrix< double > &matrix) const
Functions to support hp
virtual void get_face_interpolation_matrix (const FiniteElement< dim > &source, FullMatrix< double > &matrix) const
virtual void get_subface_interpolation_matrix (const FiniteElement< dim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Index computations
std::pair< unsigned int,
unsigned int
system_to_component_index (const unsigned int index) const
unsigned int component_to_system_index (const unsigned int component, const unsigned int index) const
std::pair< unsigned int,
unsigned int
face_system_to_component_index (const unsigned int index) const
unsigned int adjust_quad_dof_index_for_face_orientation (const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const
unsigned int adjust_line_dof_index_for_line_orientation (const unsigned int index, const bool line_orientation) const
const std::vector< bool > & get_nonzero_components (const unsigned int i) const
unsigned int n_nonzero_components (const unsigned int i) const
bool is_primitive (const unsigned int i) const
bool is_primitive () const
virtual unsigned int n_base_elements () const =0
virtual const FiniteElement
< dim > & 
base_element (const unsigned int index) const =0
virtual unsigned int element_multiplicity (const unsigned int index) const =0
std::pair< std::pair< unsigned
int, unsigned int >, unsigned
int
system_to_base_index (const unsigned int index) const
std::pair< std::pair< unsigned
int, unsigned int >, unsigned
int
face_system_to_base_index (const unsigned int index) const
unsigned int first_block_of_base (const unsigned int b) const
std::pair< unsigned int,
unsigned int
component_to_base_index (const unsigned int component) const
std::pair< unsigned int,
unsigned int
block_to_base_index (const unsigned int block) const
std::pair< unsigned int,
unsigned int
system_to_block_index (const unsigned int component) const
unsigned int component_to_block_index (const unsigned int component) const
Support points and interpolation
const std::vector< Point< dim > > & get_unit_support_points () const
bool has_support_points () const
virtual Point< dim > unit_support_point (const unsigned int index) const
const std::vector< Point< dim-1 > > & get_unit_face_support_points () const
bool has_face_support_points () const
virtual Point< dim-1 > unit_face_support_point (const unsigned int index) const
const std::vector< Point< dim > > & get_generalized_support_points () const
bool has_generalized_support_points () const
const std::vector< Point< dim-1 > > & get_generalized_face_support_points () const
bool has_generalized_face_support_points () 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

Protected Member Functions

void reinit_restriction_and_prolongation_matrices (const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
TableIndices< 2 > interface_constraints_size () const
void compute_2nd (const Mapping< dim > &mapping, const typename Triangulation< dim >::cell_iterator &cell, const unsigned int offset, typename Mapping< dim >::InternalDataBase &mapping_internal, InternalDataBase &fe_internal, FEValuesData< dim > &data) const
virtual UpdateFlags update_once (const UpdateFlags flags) const =0
virtual UpdateFlags update_each (const UpdateFlags flags) const =0
virtual FiniteElement< dim > * clone () const =0

Static Protected Member Functions

static std::vector< unsigned intcompute_n_nonzero_components (const std::vector< std::vector< bool > > &nonzero_components)

Protected Attributes

const bool cached_primitivity
std::vector< std::vector
< FullMatrix< double > > > 
restriction
std::vector< std::vector
< FullMatrix< double > > > 
prolongation
FullMatrix< doubleinterface_constraints
std::vector< Point< dim > > unit_support_points
std::vector< Point< dim-1 > > unit_face_support_points
std::vector< Point< dim > > generalized_support_points
std::vector< Point< dim-1 > > generalized_face_support_points
Table< 2, intadjust_quad_dof_index_for_face_orientation_table
std::vector< intadjust_line_dof_index_for_line_orientation_table

Private Member Functions

virtual Mapping< dim >
::InternalDataBase
get_data (const UpdateFlags flags, const Mapping< dim > &mapping, const Quadrature< dim > &quadrature) const =0
virtual Mapping< dim >
::InternalDataBase
get_face_data (const UpdateFlags flags, const Mapping< dim > &mapping, const Quadrature< dim-1 > &quadrature) const
virtual Mapping< dim >
::InternalDataBase
get_subface_data (const UpdateFlags flags, const Mapping< dim > &mapping, const Quadrature< dim-1 > &quadrature) const
virtual void fill_fe_values (const Mapping< dim > &mapping, const typename Triangulation< dim >::cell_iterator &cell, const Quadrature< dim > &quadrature, typename Mapping< dim >::InternalDataBase &mapping_internal, typename Mapping< dim >::InternalDataBase &fe_internal, FEValuesData< dim > &data) const =0
virtual void fill_fe_face_values (const Mapping< dim > &mapping, const typename Triangulation< dim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, typename Mapping< dim >::InternalDataBase &mapping_internal, typename Mapping< dim >::InternalDataBase &fe_internal, FEValuesData< dim > &data) const =0
virtual void fill_fe_subface_values (const Mapping< dim > &mapping, const typename Triangulation< dim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, typename Mapping< dim >::InternalDataBase &mapping_internal, typename Mapping< dim >::InternalDataBase &fe_internal, FEValuesData< dim > &data) const =0

Private Attributes

std::vector< std::pair
< unsigned int, unsigned int > > 
system_to_component_table
std::vector< std::pair
< unsigned int, unsigned int > > 
face_system_to_component_table
std::vector< std::pair
< std::pair< unsigned int,
unsigned int >, unsigned int > > 
system_to_base_table
std::vector< std::pair
< std::pair< unsigned int,
unsigned int >, unsigned int > > 
face_system_to_base_table
std::vector< unsigned intfirst_block_of_base_table
std::vector< std::pair
< std::pair< unsigned int,
unsigned int >, unsigned int > > 
component_to_base_table
const std::vector< boolrestriction_is_additive_flags
const std::vector< std::vector
< bool > > 
nonzero_components
const std::vector< unsigned intn_nonzero_components_table

Static Private Attributes

static const double fd_step_length

Friends

class InternalDataBase
class FEValuesBase< dim >
class FEValues< dim >
class FEFaceValues< dim >
class FESubfaceValues< dim >
class FESystem
class FE_PolyTensor
class hp::FECollection< dim >

Classes

class  ExcBoundaryFaceUsed
class  ExcComponentIndexInvalid
class  ExcConstraintsVoid
class  ExcEmbeddingVoid
class  ExcFEHasNoSupportPoints
class  ExcFENotPrimitive
class  ExcInterpolationNotImplemented
class  ExcJacobiDeterminantHasWrongSign
class  ExcProjectionVoid
class  ExcShapeFunctionNotPrimitive
class  ExcUnitShapeValuesDoNotExist
class  ExcWrongInterfaceMatrixSize
class  InternalDataBase


Detailed Description

template<int dim>
class FiniteElement< dim >

Base class for finite elements in arbitrary dimensions. This class provides several fields which describe a specific finite element and which are filled by derived classes. It more or less only offers the fields and access functions which makes it possible to copy finite elements without knowledge of the actual type (linear, quadratic, etc). In particular, the functions to fill the data fields of FEValues and its derived classes are declared.

The interface of this class is very restrictive. The reason is that finite element values should be accessed only by use of FEValues objects. These, together with FiniteElement are responsible to provide an optimized implementation.

This class declares the shape functions and their derivatives on the unit cell $[0,1]^d$. The means to transform them onto a given cell in physical space is provided by the FEValues class with a Mapping object.

The different matrices are initialized with the correct size, such that in the derived (concrete) finite element classes, their entries only have to be filled in; no resizing is needed. If the matrices are not defined by a concrete finite element, they should be resized to zero. This way functions using them can find out, that they are missing. On the other hand, it is possible to use finite element classes without implementation of the full functionality, if only part of it is needed. The functionality under consideration here is hanging nodes constraints and grid transfer, respectively.

Components and blocks

For vector valued elements shape functions may have nonzero entries in one or several components of the vector valued function. If the element is primitive, there is indeed a single component with a nonzero entry for each shape function. This component can be determined by system_to_component_index(), the number of components is FiniteElementData::n_components().

Furthermore, you may want to split your linear system into blocks for the use in BlockVector, BlockSparseMatrix, BlockMatrixArray and so on. If you use non-primitive elements, you cannot determine the block number by system_to_component_index(). Instead, you can use system_to_block_index(), which will automatically take care of the additional components occupied by vector valued elements. The number of generated blocks can be determined by FiniteElementData::n_blocks().

If you decide to operate by base element and multiplicity, the function first_block_of_base() will be helpful.

Support points

Since a FiniteElement does not have information on the actual grid cell, it can only provide support points on the unit cell. Support points on the actual grid cell must be computed by mapping these points. The class used for this kind of operation is FEValues. In most cases, code of the following type will serve to provide the mapped support points.

 Quadrature<dim> dummy_quadrature (fe.get_unit_support_points());
 FEValues<dim>   fe_values (mapping, fe, dummy_quadrature,
                            update_quadrature_points);
 fe_values.reinit (cell);
 Point<dim>& mapped_point = fe_values.quadrature_point (i);

Alternatively, the points can be transformed one-by-one:

 const vector<Point<dim> >& unit_points =
    fe.get_unit_support_points();

 Point<dim> mapped_point =
    mapping.transform_unit_to_real_cell (cell, unit_points[i]);
This is a shortcut, and as all shortcuts should be used cautiously. If the mapping of all support points is needed, the first variant should be preferred for efficiency.

Notes on the implementation of derived classes

The following sections list the information to be provided by derived classes, depending on the space dimension. They are followed by a list of functions helping to generate these values.

Finite elements in one dimension

Finite elements in one dimension need only set the restriction and prolongation matrices. The constructor of this class in one dimension presets the interface_constraints matrix to have dimension zero. Changing this behaviour in derived classes is generally not a reasonable idea and you risk getting into trouble.

Finite elements in two dimensions

In addition to the fields already present in 1D, a constraint matrix is needed, if the finite element has node values located on edges or vertices. These constraints are represented by an $m\times n$-matrix interface_constraints, where m is the number of degrees of freedom on the refined side without the corner vertices (those dofs on the middle vertex plus those on the two lines), and n is that of the unrefined side (those dofs on the two vertices plus those on the line). The matrix is thus a rectangular one. The $m\times n$ size of the interface_constraints matrix can also be accessed through the interface_constraints_size() function.

The mapping of the dofs onto the indices of the matrix on the unrefined side is as follows: let $d_v$ be the number of dofs on a vertex, $d_l$ that on a line, then $n=0...d_v-1$ refers to the dofs on vertex zero of the unrefined line, $n=d_v...2d_v-1$ to those on vertex one, $n=2d_v...2d_v+d_l-1$ to those on the line.

Similarly, $m=0...d_v-1$ refers to the dofs on the middle vertex of the refined side (vertex one of child line zero, vertex zero of child line one), $m=d_v...d_v+d_l-1$ refers to the dofs on child line zero, $m=d_v+d_l...d_v+2d_l-1$ refers to the dofs on child line one. Please note that we do not need to reserve space for the dofs on the end vertices of the refined lines, since these must be mapped one-to-one to the appropriate dofs of the vertices of the unrefined line.

It should be noted that it is not possible to distribute a constrained degree of freedom to other degrees of freedom which are themselves constrained. Only one level of indirection is allowed. It is not known at the time of this writing whether this is a constraint itself.

Finite elements in three dimensions

For the interface constraints, almost the same holds as for the 2D case. The numbering for the indices $n$ on the mother face is obvious and keeps to the usual numbering of degrees of freedom on quadrilaterals.

The numbering of the degrees of freedom on the interior of the refined faces for the index $m$ is as follows: let $d_v$ and $d_l$ be as above, and $d_q$ be the number of degrees of freedom per quadrilateral (and therefore per face), then $m=0...d_v-1$ denote the dofs on the vertex at the center, $m=d_v...5d_v-1$ for the dofs on the vertices at the center of the bounding lines of the quadrilateral, $m=5d_v..5d_v+4*d_l-1$ are for the degrees of freedom on the four lines connecting the center vertex to the outer boundary of the mother face, $m=5d_v+4*d_l...5d_v+4*d_l+8*d_l-1$ for the degrees of freedom on the small lines surrounding the quad, and $m=5d_v+12*d_l...5d_v+12*d_l+4*d_q-1$ for the dofs on the four child faces. Note the direction of the lines at the boundary of the quads, as shown below.

The order of the twelve lines and the four child faces can be extracted from the following sketch, where the overall order of the different dof groups is depicted:

 *    *--15--4--16--*
 *    |      |      |
 *    10 19  6  20  12
 *    |      |      |
 *    1--7---0--8---2
 *    |      |      |
 *    9  17  5  18  11
 *    |      |      |
 *    *--13--3--14--*
 * 
The numbering of vertices and lines, as well as the numbering of children within a line is consistent with the one described in Triangulation. Therefore, this numbering is seen from the outside and inside, respectively, depending on the face.

The three-dimensional case has a few pitfalls available for derived classes that want to implement constraint matrices. Consider the following case:

 *          *-------*
 *         /       /|
 *        /       / |
 *       /       /  |
 *      *-------*   |
 *      |       |   *-------*
 *      |       |  /       /|
 *      |   1   | /       / |
 *      |       |/       /  |
 *      *-------*-------*   |
 *      |       |       |   *
 *      |       |       |  /
 *      |   2   |   3   | /
 *      |       |       |/
 *      *-------*-------*
 * 
Now assume that we want to refine cell 2. We will end up with two faces with hanging nodes, namely the faces between cells 1 and 2, as well as between cells 2 and 3. Constraints have to be applied to the degrees of freedom on both these faces. The problem is that there is now an edge (the top right one of cell 2) which is part of both faces. The hanging node(s) on this edge are therefore constrained twice, once from both faces. To be meaningful, these constraints of course have to be consistent: both faces have to constrain the hanging nodes on the edge to the same nodes on the coarse edge (and only on the edge, as there can then be no constraints to nodes on the rest of the face), and they have to do so with the same weights. This is sometimes tricky since the nodes on the edge may have different local numbers.

For the constraint matrix this means the following: if a degree of freedom on one edge of a face is constrained by some other nodes on the same edge with some weights, then the weights have to be exactly the same as those for constrained nodes on the three other edges with respect to the corresponding nodes on these edges. If this isn't the case, you will get into trouble with the ConstraintMatrix class that is the primary consumer of the constraint information: while that class is able to handle constraints that are entered more than once (as is necessary for the case above), it insists that the weights are exactly the same.

Helper functions

Construction of a finite element and computation of the matrices described above may be a tedious task, in particular if it has to be performed for several space dimensions. Therefore, some functions in FETools have been provided to help with these tasks.

Computing the correct basis from "raw" basis functions

First, aready the basis of the shape function space may be difficult to implement for arbitrary order and dimension. On the other hand, if the node values are given, then the duality relation between node functionals and basis functions defines the basis. As a result, the shape function space may be defined with arbitrary "raw" basis functions, such that the actual finite element basis is computed from linear combinations of them. The coefficients of these combinations are determined by the duality of node values.

Using this matrix allows the construction of the basis of shape functions in two steps.

  1. Define the space of shape functions using an arbitrary basis wj and compute the matrix M of node functionals Ni applied to these basis functions, such that its entries are mij = Ni(wj).

  2. Compute the basis vj of the finite element shape function space by applying M-1 to the basis wj.

The function computing the matrix M for you is FETools::compute_node_matrix(). It relies on the existence of generalized_support_points and implementation of interpolate() with VectorSlice argument.

The piece of code in the constructor of a finite element responsible for this looks like

  FullMatrix<double> M(this->dofs_per_cell, this->dofs_per_cell);
  FETools::compute_node_matrix(M, *this);
  this->inverse_node_matrix.reinit(this->dofs_per_cell, this->dofs_per_cell);
  this->inverse_node_matrix.invert(M);
Don't forget to make sure that unit_support_points or generalized_support_points are initialized before this!

Computing the prolongation matrices for multigrid

Once the shape functions are set up, the grid transfer matrices for Multigrid accessed by get_prolongation_matrix() can be computed automatically, using FETools::compute_embedding_matrices().

This can be achieved by

  for (unsigned int i=0; i<GeometryInfo<dim>::children_per_cell; ++i)
    this->prolongation[i].reinit (this->dofs_per_cell,
                                  this->dofs_per_cell);
  FETools::compute_embedding_matrices (*this, this->prolongation);

Computing the restriction matrices for error estimators

missing...

Computing interface_constraints

Constraint matrices can be computed semi-automatically using FETools::compute_face_embedding_matrices(). This function computes the representation of the coarse mesh functions by fine mesh functions for each child of a face separately. These matrices must be convoluted into a single rectangular constraint matrix, eliminating degrees of freedom on common vertices and edges as well as on the coarse grid vertices. See the discussion above for details.

Author:
Wolfgang Bangerth, Guido Kanschat, Ralf Hartmann, 1998, 2000, 2001, 2005

Constructor & Destructor Documentation

template<int dim>
FiniteElement< dim >::FiniteElement ( const FiniteElementData< dim > &  fe_data,
const std::vector< bool > &  restriction_is_additive_flags,
const std::vector< std::vector< bool > > &  nonzero_components 
)

Constructor

template<int dim>
virtual FiniteElement< dim >::~FiniteElement (  )  [virtual]

Virtual destructor. Makes sure that pointers to this class are deleted properly.


Member Function Documentation

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

Return a string that uniquely identifies a finite element. The general convention is that this is the class name, followed by the space dimension in angle brackets, and the polynomial degree and whatever else is necessary in parentheses. For example, FE_Q<2>(3) is the value returned for a cubic element in 2d.

Systems of elements have their own naming convention, see the FESystem class.

Implemented in FE_ABF< dim >, FE_DGP< dim >, FE_DGPMonomial< dim >, FE_DGPNonparametric< dim >, FE_DGQ< dim >, FE_DGQArbitraryNodes< dim >, FE_Nedelec< dim >, FE_Q< dim >, FE_Q_Hierarchical< dim >, FE_RaviartThomas< dim >, FE_RaviartThomasNodal< dim >, and FESystem< dim >.

template<int dim>
const FiniteElement< dim > & FiniteElement< dim >::operator[] ( const unsigned int  fe_index  )  const [inline]

This operator returns a reference to the present object if the argument given equals to zero. While this does not seem particularly useful, it is helpful in writing code that works with both DoFHandler and the hp version hp::DoFHandler, since one can then write code like this:

				      *   dofs_per_cell
				      *     = dof_handler->get_fe()[cell->active_fe_index()].dofs_per_cell;
				      * 

This code doesn't work in both situations without the present operator because DoFHandler::get_fe() returns a finite element, whereas hp::DoFHandler::get_fe() returns a collection of finite elements that doesn't offer a dofs_per_cell member variable: one first has to select which finite element to work on, which is done using the operator[]. Fortunately, cell->active_fe_index() also works for non-hp classes and simply returns zero in that case. The present operator[] accepts this zero argument, by returning the finite element with index zero within its collection (that, of course, consists only of the present finite element anyway).

References Assert.

template<int dim>
virtual double FiniteElement< dim >::shape_value ( const unsigned int  i,
const Point< dim > &  p 
) const [virtual]

Return the value of the ith shape function at the point p. p is a point on the reference element. If the finite element is vector-valued, then return the value of the only non-zero component of the vector value of this shape function. If the shape function has more than one non-zero component (which we refer to with the term non-primitive), then derived classes implementing this function should throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_value_component() function.

An ExcUnitShapeValuesDoNotExist is thrown if the shape values of the FiniteElement under consideration depends on the shape of the cell in real space.

Reimplemented in FE_DGPNonparametric< dim >, FE_Poly< POLY, dim >, FE_PolyTensor< POLY, dim >, FESystem< dim >, FE_Poly< PolynomialSpace< dim >, dim >, FE_Poly< PolynomialsP< dim >, dim >, FE_Poly< TensorProductPolynomials< dim >, dim >, FE_PolyTensor< PolynomialsRaviartThomas< dim >, dim >, and FE_PolyTensor< PolynomialsABF< dim >, dim >.

template<int dim>
virtual double FiniteElement< dim >::shape_value_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const [virtual]

Just like for shape_value(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the value of the component-th vector component of the ith shape function at point p.

Reimplemented in FE_DGPNonparametric< dim >, FE_Nedelec< dim >, FE_Poly< POLY, dim >, FE_PolyTensor< POLY, dim >, FESystem< dim >, FE_Poly< PolynomialSpace< dim >, dim >, FE_Poly< PolynomialsP< dim >, dim >, FE_Poly< TensorProductPolynomials< dim >, dim >, FE_PolyTensor< PolynomialsRaviartThomas< dim >, dim >, and FE_PolyTensor< PolynomialsABF< dim >, dim >.

template<int dim>
virtual Tensor<1,dim> FiniteElement< dim >::shape_grad ( const unsigned int  i,
const Point< dim > &  p 
) const [virtual]

Return the gradient of the ith shape function at the point p. p is a point on the reference element, and likewise the gradient is the gradient on the unit cell with respect to unit cell coordinates. If the finite element is vector-valued, then return the value of the only non-zero component of the vector value of this shape function. If the shape function has more than one non-zero component (which we refer to with the term non-primitive), then derived classes implementing this function should throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_grad_component() function.

An ExcUnitShapeValuesDoNotExist is thrown if the shape values of the FiniteElement under consideration depends on the shape of the cell in real space.

Reimplemented in FE_DGPNonparametric< dim >, FE_Poly< POLY, dim >, FE_PolyTensor< POLY, dim >, FESystem< dim >, FE_Poly< PolynomialSpace< dim >, dim >, FE_Poly< PolynomialsP< dim >, dim >, FE_Poly< TensorProductPolynomials< dim >, dim >, FE_PolyTensor< PolynomialsRaviartThomas< dim >, dim >, and FE_PolyTensor< PolynomialsABF< dim >, dim >.

template<int dim>
virtual Tensor<1,dim> FiniteElement< dim >::shape_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const [virtual]

Just like for shape_grad(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the gradient of the component-th vector component of the ith shape function at point p.

Reimplemented in FE_DGPNonparametric< dim >, FE_Nedelec< dim >, FE_Poly< POLY, dim >, FE_PolyTensor< POLY, dim >, FESystem< dim >, FE_Poly< PolynomialSpace< dim >, dim >, FE_Poly< PolynomialsP< dim >, dim >, FE_Poly< TensorProductPolynomials< dim >, dim >, FE_PolyTensor< PolynomialsRaviartThomas< dim >, dim >, and FE_PolyTensor< PolynomialsABF< dim >, dim >.

template<int dim>
virtual Tensor<2,dim> FiniteElement< dim >::shape_grad_grad ( const unsigned int  i,
const Point< dim > &  p 
) const [virtual]

Return the tensor of second derivatives of the ith shape function at point p on the unit cell. The derivatives are derivatives on the unit cell with respect to unit cell coordinates. If the finite element is vector-valued, then return the value of the only non-zero component of the vector value of this shape function. If the shape function has more than one non-zero component (which we refer to with the term non-primitive), then derived classes implementing this function should throw an exception of type ExcShapeFunctionNotPrimitive. In that case, use the shape_grad_grad_component() function.

An ExcUnitShapeValuesDoNotExist is thrown if the shape values of the FiniteElement under consideration depends on the shape of the cell in real space.

Reimplemented in FE_DGPNonparametric< dim >, FE_Poly< POLY, dim >, FE_PolyTensor< POLY, dim >, FESystem< dim >, FE_Poly< PolynomialSpace< dim >, dim >, FE_Poly< PolynomialsP< dim >, dim >, FE_Poly< TensorProductPolynomials< dim >, dim >, FE_PolyTensor< PolynomialsRaviartThomas< dim >, dim >, and FE_PolyTensor< PolynomialsABF< dim >, dim >.

template<int dim>
virtual Tensor<2,dim> FiniteElement< dim >::shape_grad_grad_component ( const unsigned int  i,
const Point< dim > &  p,
const unsigned int  component 
) const [virtual]

Just like for shape_grad_grad(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the gradient of the component-th vector component of the ith shape function at point p.

Reimplemented in FE_DGPNonparametric< dim >, FE_Nedelec< dim >, FE_Poly< POLY, dim >, FE_PolyTensor< POLY, dim >, FESystem< dim >, FE_Poly< PolynomialSpace< dim >, dim >, FE_Poly< PolynomialsP< dim >, dim >, FE_Poly< TensorProductPolynomials< dim >, dim >, FE_PolyTensor< PolynomialsRaviartThomas< dim >, dim >, and FE_PolyTensor< PolynomialsABF< dim >, dim >.

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

Check for non-zero values on a face in order to optimize out matrix elements.

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

A default implementation is provided in this basis class which always returns true. This is the safe way to go.

Reimplemented in FE_ABF< dim >, FE_DGP< dim >, FE_DGPMonomial< dim >, FE_DGPNonparametric< dim >, FE_DGQ< dim >, FE_Nedelec< dim >, FE_Q< dim >, FE_Q_Hierarchical< dim >, FE_RaviartThomas< dim >, FESystem< dim >, and FE_Q_Hierarchical< dim >.

template<int dim>
const FullMatrix<double>& FiniteElement< dim >::get_restriction_matrix ( const unsigned int  child,
const RefinementCase< dim > &  refinement_case = RefinementCase< dim >::isotropic_refinement 
) const

Projection from a fine grid space onto a coarse grid space. If this projection operator is associated with a matrix P, then the restriction of this matrix P_i to a single child cell is returned here.

The matrix P is the concatenation or the sum of the cell matrices P_i, depending on the restriction_is_additive_flags. This distinguishes interpolation (concatenation) and projection with respect to scalar products (summation).

Row and column indices are related to coarse grid and fine grid spaces, respectively, consistent with the definition of the associated operator.

If projection matrices are not implemented in the derived finite element class, this function aborts with ExcProjectionVoid. You can check whether this is the case by calling the restriction_is_implemented() or the isotropic_restriction_is_implemented() function.

template<int dim>
const FullMatrix<double>& FiniteElement< dim >::get_prolongation_matrix ( const unsigned int  child,
const RefinementCase< dim > &  refinement_case = RefinementCase< dim >::isotropic_refinement 
) const

Embedding matrix between grids.

The identity operator from a coarse grid space into a fine grid space is associated with a matrix P. The restriction of this matrix P_i to a single child cell is returned here.

The matrix P is the concatenation, not the sum of the cell matrices P_i. That is, if the same non-zero entry j,k exists in in two different child matrices P_i, the value should be the same in both matrices and it is copied into the matrix P only once.

Row and column indices are related to fine grid and coarse grid spaces, respectively, consistent with the definition of the associated operator.

These matrices are used by routines assembling the prolongation matrix for multi-level methods. Upon assembling the transfer matrix between cells using this matrix array, zero elements in the prolongation matrix are discarded and will not fill up the transfer matrix.

If projection matrices are not implemented in the derived finite element class, this function aborts with ExcEmbeddingVoid. You can check whether this is the case by calling the prolongation_is_implemented() or the isotropic_prolongation_is_implemented() function.

template<int dim>
bool FiniteElement< dim >::prolongation_is_implemented (  )  const

Return whether this element implements its prolongation matrices. The return value also indicates whether a call to the get_prolongation_matrix() function will generate an error or not.

Note, that this function returns true only if the prolongation matrices of the isotropic and all anisotropic refinement cases are implemented. If you are interested in the prolongation matrices for isotropic refinement only, use the isotropic_prolongation_is_implemented function instead.

This function is mostly here in order to allow us to write more efficient test programs which we run on all kinds of weird elements, and for which we simply need to exclude certain tests in case something is not implemented. It will in general probably not be a great help in applications, since there is not much one can do if one needs these features and they are not implemented. This function could be used to check whether a call to get_prolongation_matrix() will succeed; however, one then still needs to cope with the lack of information this just expresses.

template<int dim>
bool FiniteElement< dim >::isotropic_prolongation_is_implemented (  )  const

Return whether this element implements its prolongation matrices for isotropic children. The return value also indicates whether a call to the get_prolongation_matrix function will generate an error or not.

This function is mostly here in order to allow us to write more efficient test programs which we run on all kinds of weird elements, and for which we simply need to exclude certain tests in case something is not implemented. It will in general probably not be a great help in applications, since there is not much one can do if one needs these features and they are not implemented. This function could be used to check whether a call to get_prolongation_matrix() will succeed; however, one then still needs to cope with the lack of information this just expresses.

template<int dim>
bool FiniteElement< dim >::restriction_is_implemented (  )  const

Return whether this element implements its restriction matrices. The return value also indicates whether a call to the get_restriction_matrix() function will generate an error or not.

Note, that this function returns true only if the restriction matrices of the isotropic and all anisotropic refinement cases are implemented. If you are interested in the restriction matrices for isotropic refinement only, use the isotropic_restriction_is_implemented function instead.

This function is mostly here in order to allow us to write more efficient test programs which we run on all kinds of weird elements, and for which we simply need to exclude certain tests in case something is not implemented. It will in general probably not be a great help in applications, since there is not much one can do if one needs these features and they are not implemented. This function could be used to check whether a call to get_restriction_matrix() will succeed; however, one then still needs to cope with the lack of information this just expresses.

template<int dim>
bool FiniteElement< dim >::isotropic_restriction_is_implemented (  )  const

Return whether this element implements its restriction matrices for isotropic children. The return value also indicates whether a call to the get_restriction_matrix function will generate an error or not.

This function is mostly here in order to allow us to write more efficient test programs which we run on all kinds of weird elements, and for which we simply need to exclude certain tests in case something is not implemented. It will in general probably not be a great help in applications, since there is not much one can do if one needs these features and they are not implemented. This function could be used to check whether a call to get_restriction_matrix() will succeed; however, one then still needs to cope with the lack of information this just expresses.

template<int dim>
bool FiniteElement< dim >::restriction_is_additive ( const unsigned int  index  )  const [inline]

Access the restriction_is_additive_flags field. See there for more information on its contents.

The index must be between zero and the number of shape functions of this element.

References Assert, FiniteElementData< dim >::dofs_per_cell, and FiniteElement< dim >::restriction_is_additive_flags.

template<int dim>
const FullMatrix<double>& FiniteElement< dim >::constraints ( const internal::SubfaceCase< dim > &  subface_case = internal::SubfaceCase< dim >::case_isotropic  )  const

Return a readonly reference to the matrix which describes the constraints at the interface between a refined and an unrefined cell.

The matrix is obviously empty in only one space dimension, since there are no constraints then.

Note that some finite elements do not (yet) implement hanging node constraints. If this is the case, then this function will generate an exception, since no useful return value can be generated. If you should have a way to live with this, then you might want to use the constraints_are_implemented() function to check up front whethehr this function will succeed or generate the exception.

template<int dim>
bool FiniteElement< dim >::constraints_are_implemented ( const internal::SubfaceCase< dim > &  subface_case = internal::SubfaceCase< dim >::case_isotropic  )  const

Return whether this element implements its hanging node constraints. The return value also indicates whether a call to the constraints() function will generate an error or not.

This function is mostly here in order to allow us to write more efficient test programs which we run on all kinds of weird elements, and for which we simply need to exclude certain tests in case hanging node constraints are not implemented. It will in general probably not be a great help in applications, since there is not much one can do if one needs hanging node constraints and they are not implemented. This function could be used to check whether a call to constraints() will succeed; however, one then still needs to cope with the lack of information this just expresses.

template<int dim>
virtual bool FiniteElement< dim >::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". That means, the element properly implements the get_face_interpolation_matrix and get_subface_interpolation_matrix methods. Therefore the return value also indicates whether a call to the get_face_interpolation_matrix() method and the get_subface_interpolation_matrix() method will generate an error or not.

Currently the main purpose of this function is to allow the make_hanging_node_constraints method to decide whether the new procedures, which are supposed to work in the hp framework can be used, or if the old well verified but not hp capable functions should be used. Once the transition to the new scheme for computing the interface constraints is complete, this function will be superfluous and will probably go away.

Derived classes should implement this function accordingly. The default assumption is that a finite element does not provide hp capable face interpolation, and the default implementation therefore returns false.

Reimplemented in FE_DGP< dim >, FE_DGPMonomial< dim >, FE_DGPNonparametric< dim >, FE_DGQ< dim >, FE_Q< dim >, FE_RaviartThomasNodal< dim >, and FESystem< dim >.

template<int dim>
virtual void FiniteElement< dim >::get_interpolation_matrix ( const FiniteElement< dim > &  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.

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

Reimplemented in FE_DGPMonomial< dim >, FE_DGQ< dim >, FE_Q< dim >, and FESystem< dim >.

template<int dim>
virtual void FiniteElement< dim >::get_face_interpolation_matrix ( const FiniteElement< dim > &  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 ExcInterpolationNotImplemented.

Reimplemented in FE_DGP< dim >, FE_DGPMonomial< dim >, FE_DGPNonparametric< dim >, FE_DGQ< dim >, FE_Q< dim >, FE_RaviartThomasNodal< dim >, and