Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Namespaces | Classes | Functions
FETools Namespace Reference

Namespaces

namespace  Compositing
 

Classes

class  FEFactory
 
class  FEFactoryBase
 

Functions

template<int dim, int spacedim>
void compute_component_wise (const FiniteElement< dim, spacedim > &fe, std::vector< unsigned int > &renumbering, std::vector< std::vector< unsigned int > > &start_indices)
 
template<int dim, int spacedim>
void compute_block_renumbering (const FiniteElement< dim, spacedim > &fe, std::vector< types::global_dof_index > &renumbering, std::vector< types::global_dof_index > &block_data, bool return_start_indices=true)
 
template<int dim>
std::vector< unsigned inthierarchic_to_lexicographic_numbering (unsigned int degree)
 
template<int dim>
std::vector< unsigned intlexicographic_to_hierarchic_numbering (unsigned int degree)
 
template<int dim, int spacedim = dim>
std::unique_ptr< FiniteElement< dim, spacedim > > get_fe_by_name (const std::string &name)
 
template<int dim, int spacedim>
void add_fe_name (const std::string &name, const FEFactoryBase< dim, spacedim > *factory)
 
static ::ExceptionBaseExcInvalidFEName (std::string arg1)
 
static ::ExceptionBaseExcInvalidFEDimension (char arg1, int arg2)
 
static ::ExceptionBaseExcInvalidFE ()
 
static ::ExceptionBaseExcFENotPrimitive ()
 
static ::ExceptionBaseExcTriangulationMismatch ()
 
static ::ExceptionBaseExcHangingNodesNotAllowed ()
 
static ::ExceptionBaseExcGridNotRefinedAtLeastOnce ()
 
static ::ExceptionBaseExcMatrixDimensionMismatch (int arg1, int arg2, int arg3, int arg4)
 
static ::ExceptionBaseExcLeastSquaresError (double arg1)
 
static ::ExceptionBaseExcNotGreaterThan (int arg1, int arg2)
 
template std::vector< unsigned inthierarchic_to_lexicographic_numbering< 0 > (unsigned int)
 
template std::vector< unsigned intlexicographic_to_hierarchic_numbering< 0 > (unsigned int)
 
Generation of local matrices
template<int dim, typename number , int spacedim>
void get_interpolation_matrix (const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &interpolation_matrix)
 
template<int dim, typename number , int spacedim>
void get_back_interpolation_matrix (const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &interpolation_matrix)
 
template<int dim, typename number , int spacedim>
void get_interpolation_difference_matrix (const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &difference_matrix)
 
template<int dim, typename number , int spacedim>
void get_projection_matrix (const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &matrix)
 
template<int dim, int spacedim>
FullMatrix< double > compute_node_matrix (const FiniteElement< dim, spacedim > &fe)
 
template<int dim, typename number , int spacedim>
void compute_embedding_matrices (const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false, const double threshold=1.e-12)
 
template<int dim, typename number , int spacedim>
void compute_face_embedding_matrices (const FiniteElement< dim, spacedim > &fe, FullMatrix< number >(&matrices)[GeometryInfo< dim >::max_children_per_face], const unsigned int face_coarse, const unsigned int face_fine, const double threshold=1.e-12)
 
template<int dim, typename number , int spacedim>
void compute_projection_matrices (const FiniteElement< dim, spacedim > &fe, std::vector< std::vector< FullMatrix< number > > > &matrices, const bool isotropic_only=false)
 
template<int dim, int spacedim>
void compute_projection_from_quadrature_points_matrix (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &lhs_quadrature, const Quadrature< dim > &rhs_quadrature, FullMatrix< double > &X)
 
template<int dim, int spacedim>
void compute_interpolation_to_quadrature_points_matrix (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim > &quadrature, FullMatrix< double > &I_q)
 
template<int dim>
void compute_projection_from_quadrature_points (const FullMatrix< double > &projection_matrix, const std::vector< Tensor< 1, dim > > &vector_of_tensors_at_qp, std::vector< Tensor< 1, dim > > &vector_of_tensors_at_nodes)
 
template<int dim>
void compute_projection_from_quadrature_points (const FullMatrix< double > &projection_matrix, const std::vector< SymmetricTensor< 2, dim > > &vector_of_tensors_at_qp, std::vector< SymmetricTensor< 2, dim > > &vector_of_tensors_at_nodes)
 
template<int dim, int spacedim>
void compute_projection_from_face_quadrature_points_matrix (const FiniteElement< dim, spacedim > &fe, const Quadrature< dim - 1 > &lhs_quadrature, const Quadrature< dim - 1 > &rhs_quadrature, const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face, FullMatrix< double > &X)
 
template<int dim, int spacedim, typename number >
void convert_generalized_support_point_values_to_dof_values (const FiniteElement< dim, spacedim > &finite_element, const std::vector< Vector< number > > &support_point_values, std::vector< number > &dof_values)
 
Functions which should be in DoFTools
template<int dim, int spacedim, class InVector , class OutVector >
void interpolate (const DoFHandler< dim, spacedim > &dof1, const InVector &u1, const DoFHandler< dim, spacedim > &dof2, OutVector &u2)
 
template<int dim, int spacedim, class InVector , class OutVector >
void interpolate (const DoFHandler< dim, spacedim > &dof1, const InVector &u1, const DoFHandler< dim, spacedim > &dof2, const AffineConstraints< typename OutVector::value_type > &constraints, OutVector &u2)
 
template<int dim, class InVector , class OutVector , int spacedim>
void back_interpolate (const DoFHandler< dim, spacedim > &dof1, const InVector &u1, const FiniteElement< dim, spacedim > &fe2, OutVector &u1_interpolated)
 
template<int dim, class InVector , class OutVector , int spacedim>
void back_interpolate (const DoFHandler< dim, spacedim > &dof1, const AffineConstraints< typename OutVector::value_type > &constraints1, const InVector &u1, const DoFHandler< dim, spacedim > &dof2, const AffineConstraints< typename OutVector::value_type > &constraints2, OutVector &u1_interpolated)
 
template<int dim, class InVector , class OutVector , int spacedim>
void interpolation_difference (const DoFHandler< dim, spacedim > &dof1, const InVector &z1, const FiniteElement< dim, spacedim > &fe2, OutVector &z1_difference)
 
template<int dim, class InVector , class OutVector , int spacedim>
void interpolation_difference (const DoFHandler< dim, spacedim > &dof1, const AffineConstraints< typename OutVector::value_type > &constraints1, const InVector &z1, const DoFHandler< dim, spacedim > &dof2, const AffineConstraints< typename OutVector::value_type > &constraints2, OutVector &z1_difference)
 
template<int dim, class InVector , class OutVector , int spacedim>
void project_dg (const DoFHandler< dim, spacedim > &dof1, const InVector &u1, const DoFHandler< dim, spacedim > &dof2, OutVector &u2)
 
template<int dim, class InVector , class OutVector , int spacedim>
void extrapolate (const DoFHandler< dim, spacedim > &dof1, const InVector &z1, const DoFHandler< dim, spacedim > &dof2, OutVector &z2)
 
template<int dim, class InVector , class OutVector , int spacedim>
void extrapolate (const DoFHandler< dim, spacedim > &dof1, const InVector &z1, const DoFHandler< dim, spacedim > &dof2, const AffineConstraints< typename OutVector::value_type > &constraints, OutVector &z2)
 

Detailed Description

This namespace offers interpolations and extrapolations of discrete functions of one FiniteElement fe1 to another FiniteElement fe2.

It also provides the local interpolation matrices that interpolate on each cell. Furthermore it provides the difference matrix \(id-I_h\) that is needed for evaluating \((id-I_h)z\) for e.g. the dual solution \(z\).

For more information about the spacedim template parameter check the documentation of FiniteElement or the one of Triangulation.

Function Documentation

◆ compute_component_wise()

template<int dim, int spacedim>
void FETools::compute_component_wise ( const FiniteElement< dim, spacedim > &  fe,
std::vector< unsigned int > &  renumbering,
std::vector< std::vector< unsigned int > > &  start_indices 
)
Warning
In most cases, you will probably want to use compute_base_renumbering().

Compute the vector required to renumber the dofs of a cell by component. Furthermore, compute the vector storing the start indices of each component in the local block vector.

The second vector is organized such that there is a vector for each base element containing the start index for each component served by this base element.

While the first vector is checked to have the correct size, the second one is reinitialized for convenience.

◆ compute_block_renumbering()

template<int dim, int spacedim>
void FETools::compute_block_renumbering ( const FiniteElement< dim, spacedim > &  fe,
std::vector< types::global_dof_index > &  renumbering,
std::vector< types::global_dof_index > &  block_data,
bool  return_start_indices = true 
)

Compute the vector required to renumber the dofs of a cell by block. Furthermore, compute the vector storing either the start indices or the size of each local block vector.

If the bool parameter is true, block_data is filled with the start indices of each local block. If it is false, then the block sizes are returned.

The vector renumbering will be indexed by the standard numbering of local degrees of freedom, namely the first vertex, then the second vertex, after vertices lines, quads, and hexes. For each index, the entry indicates the index which this degree of freedom receives in a numbering scheme, where the first block is numbered completely before the second.

◆ get_interpolation_matrix()

template<int dim, typename number , int spacedim>
void FETools::get_interpolation_matrix ( const FiniteElement< dim, spacedim > &  fe1,
const FiniteElement< dim, spacedim > &  fe2,
FullMatrix< number > &  interpolation_matrix 
)

Compute the interpolation matrix that interpolates a fe1-function to a fe2-function on each cell. The interpolation_matrix needs to be of size (fe2.n_dofs_per_cell(), fe1.n_dofs_per_cell()).

Note, that if the finite element space fe1 is a subset of the finite element space fe2 then the interpolation_matrix is an embedding matrix.

◆ get_back_interpolation_matrix()

template<int dim, typename number , int spacedim>
void FETools::get_back_interpolation_matrix ( const FiniteElement< dim, spacedim > &  fe1,
const FiniteElement< dim, spacedim > &  fe2,
FullMatrix< number > &  interpolation_matrix 
)

Compute the interpolation matrix that interpolates a fe1-function to a fe2-function, and interpolates this to a second fe1-function on each cell. The interpolation_matrix needs to be of size (fe1.n_dofs_per_cell(), fe1.n_dofs_per_cell()).

Note, that this function only makes sense if the finite element space due to fe1 is not a subset of the finite element space due to fe2, as if it were a subset then the interpolation_matrix would be only the unit matrix.

◆ get_interpolation_difference_matrix()

template<int dim, typename number , int spacedim>
void FETools::get_interpolation_difference_matrix ( const FiniteElement< dim, spacedim > &  fe1,
const FiniteElement< dim, spacedim > &  fe2,
FullMatrix< number > &  difference_matrix 
)

Compute the identity matrix minus the back interpolation matrix. The difference_matrix will be of size (fe1.n_dofs_per_cell(), fe1.n_dofs_per_cell()) after this function. Previous content of the argument will be overwritten.

This function computes the matrix that transforms a fe1 function \(z\) to \(z-I_hz\) where \(I_h\) denotes the interpolation operator from the fe1 space to the fe2 space. This matrix hence is useful to evaluate error-representations where \(z\) denotes the dual solution.

◆ get_projection_matrix()

template<int dim, typename number , int spacedim>
void FETools::get_projection_matrix ( const FiniteElement< dim, spacedim > &  fe1,
const FiniteElement< dim, spacedim > &  fe2,
FullMatrix< number > &  matrix 
)

Compute the local \(L^2\)-projection matrix from fe1 to fe2.

◆ compute_node_matrix()

template<int dim, int spacedim>
FullMatrix< double > FETools::compute_node_matrix ( const FiniteElement< dim, spacedim > &  fe)

This is a rather specialized function used during the construction of finite element objects. It is used to build the basis of shape functions for an element, given a set of polynomials and interpolation points. The function is only implemented for finite elements with exactly dim vector components. In particular, this applies to classes derived from the FE_PolyTensor class.

Specifically, the purpose of this function is as follows: FE_PolyTensor receives, from its derived classes, an argument that describes a polynomial space. This space may be parameterized in terms of monomials, or in some other way, but is in general not in the form that we use for finite elements where we typically want to use a basis that is derived from some kind of node functional (e.g., the interpolation at specific points). Concretely, assume that the basis used by the polynomial space is \(\{\tilde\varphi_j(\mathbf x)\}_{j=1}^N\), and that the node functionals of the finite element are \(\{\Psi_i\}_{i=1}^N\). We then want to compute a basis \(\{\varphi_j(\mathbf x)\}_{j=1}^N\) for the finite element space so that \(\Psi_i[\varphi_j] = \delta_{ij}\). To do this, we can set \(\varphi_j(\mathbf x) = \sum_{k=1}^N c_{jk} \tilde\varphi_k(\mathbf x)\) where we need to determine the expansion coefficients \(c_{jk}\). We do this by applying \(\Psi_i\) to both sides of the equation, to obtain

\begin{align*} \Psi_i [\varphi_j] = \sum_{k=1}^N c_{jk} \Psi_i[\tilde\varphi_k], \end{align*}

and we know that the left hand side equals \(\delta_{ij}\). If you think of this as a system of \(N\times N\) equations for the elements of a matrix on the left and on the right, then this can be written as

\begin{align*} I = C X^T \end{align*}

where \(C\) is the matrix of coefficients \(c_{jk}\) and \(X_{ik} = \Psi_i[\tilde\varphi_k]\). Consequently, in order to compute the expansion coefficients \(C=X^{-T}\), we need to apply the node functionals to all functions of the "raw" basis of the polynomial space.

Until the finite element receives this matrix \(X\) back, it describes its shape functions (e.g., in FiniteElement::shape_value()) in the form \(\tilde\varphi_j\). After it calls this function, it has the expansion coefficients and can describe its shape functions as \(\varphi_j\).

This function therefore computes this matrix \(X\), for the following specific circumstances:

  • That the node functionals \(\Psi_i\) are point evaluations at points \(\mathbf x_i\) that the finite element in question describes via its "generalized" support points (through FiniteElement::get_generalized_support_points(), see also this glossary entry). These point evaluations need to necessarily evaluate the value of a shape function at that point (the shape function may be vector-valued, and so the functional may be a linear combination of the individual components of the values); but, in particular, the nodal functions may not be integrals over entire edges or faces, or other non-local functionals. In other words, we assume that \(\Psi_i[\tilde\varphi_j] = f_j(\tilde\varphi_j(\mathbf x_i))\) where \(f_j\) is a function of the (possibly vector-valued) argument that returns a scalar.
  • That the finite element has exactly dim vector components.
  • That the function \(f_j\) is given by whatever the element implements through the FiniteElement::convert_generalized_support_point_values_to_dof_values() function.
Parameters
feThe finite element for which the operations above are to be performed.
Returns
The matrix \(X\) as discussed above.

◆ compute_embedding_matrices()

template<int dim, typename number , int spacedim>
void FETools::compute_embedding_matrices ( const FiniteElement< dim, spacedim > &  fe,
std::vector< std::vector< FullMatrix< number > > > &  matrices,
const bool  isotropic_only = false,
const double  threshold = 1.e-12 
)

For all possible (isotropic and anisotropic) refinement cases compute the embedding matrices from a coarse cell to the child cells. Each column of the resulting matrices contains the representation of a coarse grid basis function by the fine grid basis; the matrices are split such that there is one matrix for every child.

This function computes the coarse grid function in a sufficiently large number of quadrature points and fits the fine grid functions using least squares approximation. Therefore, the use of this function is restricted to the case that the finite element spaces are actually nested.

Note, that matrices[refinement_case-1][child] includes the embedding (or prolongation) matrix of child child for the RefinementCase refinement_case. Here, we use refinement_case-1 instead of refinement_case as for RefinementCase::no_refinement(=0) there are no prolongation matrices available.

Typically this function is called by the various implementations of FiniteElement classes in order to fill the respective FiniteElement::prolongation matrices.

Parameters
feThe finite element class for which we compute the embedding matrices.
matricesA reference to RefinementCase<dim>::isotropic_refinement vectors of FullMatrix objects. Each vector corresponds to one RefinementCase refinement_case and is of the vector size GeometryInfo<dim>::n_children(refinement_case). This is the format used in FiniteElement, where we want to use this function mostly.
isotropic_onlySet to true if you only want to compute matrices for isotropic refinement.
thresholdis the gap allowed in the least squares algorithm computing the embedding.

◆ compute_face_embedding_matrices()

template<int dim, typename number , int spacedim>
void FETools::compute_face_embedding_matrices ( const FiniteElement< dim, spacedim > &  fe,
FullMatrix< number >(&)  matrices[GeometryInfo< dim >::max_children_per_face],
const unsigned int  face_coarse,
const unsigned int  face_fine,
const double  threshold = 1.e-12 
)

Compute the embedding matrices on faces needed for constraint matrices.

Parameters
feThe finite element for which to compute these matrices.
matricesAn array of GeometryInfo<dim>::subfaces_per_face = 2dim-1 FullMatrix objects,holding the embedding matrix for each subface.
face_coarseThe number of the face on the coarse side of the face for which this is computed.
face_fineThe number of the face on the refined side of the face for which this is computed.
thresholdis the gap allowed in the least squares algorithm computing the embedding.
Warning
This function will be used in computing constraint matrices. It is not sufficiently tested yet.

◆ compute_projection_matrices()

template<int dim, typename number , int spacedim>
void FETools::compute_projection_matrices ( const FiniteElement< dim, spacedim > &  fe,
std::vector< std::vector< FullMatrix< number > > > &  matrices,
const bool  isotropic_only = false 
)

For all possible (isotropic and anisotropic) refinement cases compute the L2-projection matrices from the children to a coarse cell.

Note, that matrices[refinement_case-1][child] includes the projection (or restriction) matrix of child child for the RefinementCase refinement_case. Here, we use refinement_case-1 instead of refinement_case as for RefinementCase::no_refinement(=0) there are no projection matrices available.

Typically this function is called by the various implementations of FiniteElement classes in order to fill the respective FiniteElement::restriction matrices.

  • [in] fe The finite element class for which we compute the projection matrices.
  • [in] isotropic_only If set to true, then this function only computes data for the isotropic refinement case. The other elements of the output vector are left untouched (but still exist).

◆ compute_projection_from_quadrature_points_matrix()

template<int dim, int spacedim>
void FETools::compute_projection_from_quadrature_points_matrix ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  lhs_quadrature,
const Quadrature< dim > &  rhs_quadrature,
FullMatrix< double > &  X 
)

Project scalar data defined in quadrature points to a finite element space on a single cell.

What this function does is the following: assume that there is scalar data uq, 0 <= q < Q:=quadrature.size() defined at the quadrature points of a cell, with the points defined by the given rhs_quadrature object. We may then want to ask for that finite element function (on a single cell) vh in the finite- dimensional space defined by the given FE object that is the projection of u in the following sense:

Usually, the projection vh is that function that satisfies (vh,w)=(u,w) for all discrete test functions w. In the present case, we can't evaluate the right hand side, since u is only defined in the quadrature points given by rhs_quadrature, so we replace it by a quadrature approximation. Likewise, the left hand side is approximated using the lhs_quadrature object; if this quadrature object is chosen appropriately, then the integration of the left hand side can be done exactly, without any approximation. The use of different quadrature objects is necessary if the quadrature object for the right hand side has too few quadrature points – for example, if data q is only defined at the cell center, then the corresponding one-point quadrature formula is obviously insufficient to approximate the scalar product on the left hand side by a definite form.

After these quadrature approximations, we end up with a nodal representation Vh of vh that satisfies the following system of linear equations: M Vh = Q U, where Mij=(phi_i,phi_j) is the mass matrix approximated by lhs_quadrature, and Q is the matrix Qiq=phii(xq) wq where wq are quadrature weights; U is the vector of quadrature point data uq.

In order to then get the nodal representation Vh of the projection of U, one computes Vh = X U, X=M-1 Q. The purpose of this function is to compute the matrix X and return it through the last argument of this function.

Note that this function presently only supports scalar data. An extension of the mass matrix is of course trivial, but one has to define the order of data in the vector U if it contains vector valued data in all quadrature points.

A use for this function is described in the introduction to the step-18 example program.

The opposite of this function, interpolation of a finite element function onto quadrature points is essentially what the FEValues::get_function_values functions do; to make things a little simpler, the FETools::compute_interpolation_to_quadrature_points_matrix provides the matrix form of this.

Note that this function works on a single cell, rather than an entire triangulation. In effect, it therefore doesn't matter if you use a continuous or discontinuous version of the finite element.

It is worth noting that there are a few confusing cases of this function. The first one is that it really only makes sense to project onto a finite element that has at most as many degrees of freedom per cell as there are quadrature points; the projection of N quadrature point data into a space with M>N unknowns is well-defined, but often yields funny and non- intuitive results. Secondly, one would think that if the quadrature point data is defined in the support points of the finite element, i.e. the quadrature points of rhs_quadrature equal fe.get_unit_support_points(), then the projection should be the identity, i.e. each degree of freedom of the finite element equals the value of the given data in the support point of the corresponding shape function. However, this is not generally the case: while the matrix Q in that case is the identity matrix, the mass matrix M is not equal to the identity matrix, except for the special case that the quadrature formula lhs_quadrature also has its quadrature points in the support points of the finite element.

Finally, this function only defines a cell wise projection, while one frequently wants to apply it to all cells in a triangulation. However, if it is applied to one cell after the other, the results from later cells may overwrite nodal values computed already from previous cells if degrees of freedom live on the interfaces between cells. The function is therefore most useful for discontinuous elements.

◆ compute_interpolation_to_quadrature_points_matrix()

template<int dim, int spacedim>
void FETools::compute_interpolation_to_quadrature_points_matrix ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim > &  quadrature,
FullMatrix< double > &  I_q 
)

Given a (scalar) local finite element function, compute the matrix that maps the vector of nodal values onto the vector of values of this function at quadrature points as given by the second argument. In a sense, this function does the opposite of the FETools::compute_projection_from_quadrature_points_matrix function.

◆ compute_projection_from_quadrature_points() [1/2]

template<int dim>
void FETools::compute_projection_from_quadrature_points ( const FullMatrix< double > &  projection_matrix,
const std::vector< Tensor< 1, dim > > &  vector_of_tensors_at_qp,
std::vector< Tensor< 1, dim > > &  vector_of_tensors_at_nodes 
)

Compute the projection of tensorial (first-order tensor) data stored at the quadrature points vector_of_tensors_at_qp to data vector_of_tensors_at_nodes at the support points of the cell. The data in vector_of_tensors_at_qp is ordered sequentially following the quadrature point numbering. The size of vector_of_tensors_at_qp must correspond to the number of columns of projection_matrix. The size of vector_of_tensors_at_nodes must correspond to the number of rows of vector_of_tensors_at_nodes . The projection matrix projection_matrix describes the projection of scalar data from the quadrature points and can be obtained from the FETools::compute_projection_from_quadrature_points_matrix function.

◆ compute_projection_from_quadrature_points() [2/2]

template<int dim>
void FETools::compute_projection_from_quadrature_points ( const FullMatrix< double > &  projection_matrix,
const std::vector< SymmetricTensor< 2, dim > > &  vector_of_tensors_at_qp,
std::vector< SymmetricTensor< 2, dim > > &  vector_of_tensors_at_nodes 
)

same as last function but for a SymmetricTensor .

◆ compute_projection_from_face_quadrature_points_matrix()

template<int dim, int spacedim>
void FETools::compute_projection_from_face_quadrature_points_matrix ( const FiniteElement< dim, spacedim > &  fe,
const Quadrature< dim - 1 > &  lhs_quadrature,
const Quadrature< dim - 1 > &  rhs_quadrature,
const typename DoFHandler< dim, spacedim >::active_cell_iterator &  cell,
const unsigned int  face,
FullMatrix< double > &  X 
)

This method implements the FETools::compute_projection_from_quadrature_points_matrix method for faces of a mesh. The matrix that it returns, X, is face specific and its size is fe.n_dofs_per_cell() by rhs_quadrature.size(). The dimension, dim must be larger than 1 for this class, since Quadrature<dim-1> objects are required. See the documentation on the Quadrature class for more information.

◆ convert_generalized_support_point_values_to_dof_values()

template<int dim, int spacedim, typename number >
void FETools::convert_generalized_support_point_values_to_dof_values ( const FiniteElement< dim, spacedim > &  finite_element,
const std::vector< Vector< number > > &  support_point_values,
std::vector< number > &  dof_values 
)

Wrapper around FiniteElement::convert_generalized_support_point_values_to_dof_values() that works with arbitrary number types.

Parameters
[in]finite_elementThe FiniteElement to compute dof values for.
[in]support_point_valuesAn array of size dofs_per_cell (which equals the number of points the get_generalized_support_points() function will return) where each element is a vector with as many entries as the element has vector components. This array should contain the values of a function at the generalized support points of the finite element.
[out]dof_valuesAn array of size dofs_per_cell that contains the node functionals of the element applied to the given function.

◆ interpolate() [1/2]

template<int dim, int spacedim, class InVector , class OutVector >
void FETools::interpolate ( const DoFHandler< dim, spacedim > &  dof1,
const InVector &  u1,
const DoFHandler< dim, spacedim > &  dof2,
OutVector &  u2 
)

Compute the interpolation of a the dof1-function u1 to a dof2-function u2. dof1 and dof2 need to be DoFHandlers based on the same triangulation.

If the elements fe1 and fe2 are either both continuous or both discontinuous then this interpolation is the usual point interpolation. The same is true if fe1 is a continuous and fe2 is a discontinuous finite element. For the case that fe1 is a discontinuous and fe2 is a continuous finite element there is no point interpolation defined at the discontinuities. Therefore the mean value is taken at the DoF values on the discontinuities.

Note that for continuous elements on grids with hanging nodes (i.e. locally refined grids) this function does not give the expected output. Indeed, the resulting output vector does not necessarily respect continuity requirements at hanging nodes: if, for example, you are interpolating a Q2 field to a Q1 field, then at hanging nodes the output field will have the function value of the input field, which however is not usually the mean value of the two adjacent nodes. It is thus not part of the Q1 function space on the whole triangulation, although it is of course Q1 on each cell.

For this case (continuous elements on grids with hanging nodes), please use the interpolate() function with an additional AffineConstraints object as argument, see below, or make the field conforming yourself by calling the distribute function of your hanging node constraints object.

◆ interpolate() [2/2]

template<int dim, int spacedim, class InVector , class OutVector >
void FETools::interpolate ( const DoFHandler< dim, spacedim > &  dof1,
const InVector &  u1,
const DoFHandler< dim, spacedim > &  dof2,
const AffineConstraints< typename OutVector::value_type > &  constraints,
OutVector &  u2 
)

Compute the interpolation of a the dof1-function u1 to a dof2-function u2. dof1 and dof2 need to be DoFHandlers based on the same triangulation. constraints is a hanging node constraints object corresponding to dof2. This object is particular important when interpolating onto continuous elements on grids with hanging nodes (locally refined grids).

If the elements fe1 and fe2 are either both continuous or both discontinuous then this interpolation is the usual point interpolation. The same is true if fe1 is a continuous and fe2 is a discontinuous finite element. For the case that fe1 is a discontinuous and fe2 is a continuous finite element there is no point interpolation defined at the discontinuities. Therefore the mean value is taken at the DoF values at the discontinuities.

◆ back_interpolate() [1/2]

template<int dim, class InVector , class OutVector , int spacedim>
void FETools::back_interpolate ( const DoFHandler< dim, spacedim > &  dof1,
const InVector &  u1,
const FiniteElement< dim, spacedim > &  fe2,
OutVector &  u1_interpolated 
)

Compute the interpolation of the fe1-function u1 to a fe2-function, and interpolates this to a second fe1-function named u1_interpolated.

Note, that this function does not work on continuous elements at hanging nodes. For that case use the back_interpolate function, below, that takes an additional AffineConstraints object.

Furthermore note, that for the specific case when the finite element space corresponding to fe1 is a subset of the finite element space corresponding to fe2, this function is simply an identity mapping.

◆ back_interpolate() [2/2]

template<int dim, class InVector , class OutVector , int spacedim>
void FETools::back_interpolate ( const DoFHandler< dim, spacedim > &  dof1,
const AffineConstraints< typename OutVector::value_type > &  constraints1,
const InVector &  u1,
const DoFHandler< dim, spacedim > &  dof2,
const AffineConstraints< typename OutVector::value_type > &  constraints2,
OutVector &  u1_interpolated 
)

Compute the interpolation of the dof1-function u1 to a dof2-function, and interpolates this to a second dof1-function named u1_interpolated. constraints1 and constraints2 are the hanging node constraints corresponding to dof1 and dof2, respectively. These objects are particular important when continuous elements on grids with hanging nodes (locally refined grids) are involved.

Furthermore note, that for the specific case when the finite element space corresponding to dof1 is a subset of the finite element space corresponding to dof2, this function is simply an identity mapping.

◆ interpolation_difference() [1/2]

template<int dim, class InVector , class OutVector , int spacedim>
void FETools::interpolation_difference ( const DoFHandler< dim, spacedim > &  dof1,
const InVector &  z1,
const FiniteElement< dim, spacedim > &  fe2,
OutVector &  z1_difference 
)

Compute \((Id-I_h)z_1\) for a given dof1-function \(z_1\), where \(I_h\) is the interpolation from fe1 to fe2. The result \((Id-I_h)z_1\) is written into z1_difference.

Note, that this function does not work for continuous elements at hanging nodes. For that case use the interpolation_difference function, below, that takes an additional AffineConstraints object.

◆ interpolation_difference() [2/2]

template<int dim, class InVector , class OutVector , int spacedim>
void FETools::interpolation_difference ( const DoFHandler< dim, spacedim > &  dof1,
const AffineConstraints< typename OutVector::value_type > &  constraints1,
const InVector &  z1,
const DoFHandler< dim, spacedim > &  dof2,
const AffineConstraints< typename OutVector::value_type > &  constraints2,
OutVector &  z1_difference 
)

Compute \((Id-I_h)z_1\) for a given dof1-function \(z_1\), where \(I_h\) is the interpolation from fe1 to fe2. The result \((Id-I_h)z_1\) is written into z1_difference. constraints1 and constraints2 are the hanging node constraints corresponding to dof1 and dof2, respectively. These objects are particular important when continuous elements on grids with hanging nodes (locally refined grids) are involved.

For parallel computations, supply z1 with ghost elements and z1_difference without ghost elements.

◆ project_dg()

template<int dim, class InVector , class OutVector , int spacedim>
void FETools::project_dg ( const DoFHandler< dim, spacedim > &  dof1,
const InVector &  u1,
const DoFHandler< dim, spacedim > &  dof2,
OutVector &  u2 
)

\(L^2\) projection for discontinuous elements. Operates the same direction as interpolate.

The global projection can be computed by local matrices if the finite element spaces are discontinuous. With continuous elements, this is impossible, since a global mass matrix must be inverted.

◆ extrapolate() [1/2]

template<int dim, class InVector , class OutVector , int spacedim>
void FETools::extrapolate ( const DoFHandler< dim, spacedim > &  dof1,
const InVector &  z1,
const DoFHandler< dim, spacedim > &  dof2,
OutVector &  z2 
)

Compute the patchwise extrapolation of a dof1 function z1 to a dof2 function z2. dof1 and dof2 need to be DoFHandler objects based on the same triangulation. This function is used, for example, for extrapolating patchwise a piecewise linear solution to a piecewise quadratic solution.

The function's name is historical and probably not particularly well chosen. The function performs the following operations, one after the other:

  • It interpolates directly from every cell of dof1 to the corresponding cell of dof2 using the interpolation matrix of the finite element spaces used on these cells and provided by the finite element objects involved. This step is done using the FETools::interpolate() function.
  • It then performs a loop over all non-active cells of dof2. If such a non-active cell has at least one active child, then we call the children of this cell a "patch". We then interpolate from the children of this patch to the patch, using the finite element space associated with dof2 and immediately interpolate back to the children. In essence, this information throws away all information in the solution vector that lives on a scale smaller than the patch cell.
  • Since we traverse non-active cells from the coarsest to the finest levels, we may find patches that correspond to child cells of previously treated patches if the mesh had been refined adaptively (this cannot happen if the mesh has been refined globally because there the children of a patch are all active). We also perform the operation described above on these patches, but it is easy to see that on patches that are children of previously treated patches, the operation is now the identity operation (since it interpolates from the children of the current patch a function that had previously been interpolated to these children from an even coarser patch). Consequently, this does not alter the solution vector any more.

The name of the function originates from the fact that it can be used to construct a representation of a function of higher polynomial degree on a once coarser mesh. For example, if you imagine that you start with a \(Q_1\) function on a globally refined mesh, and that dof2 is associated with a \(Q_2\) element, then this function computes the equivalent of the operator \(I_{2h}^{(2)}\) interpolating the original piecewise linear function onto a quadratic function on a once coarser mesh with mesh size \(2h\) (but representing this function on the original mesh with size \(h\)). If the exact solution is sufficiently smooth, then \(u^\ast=I_{2h}^{(2)}u_h\) is typically a better approximation to the exact solution \(u\) of the PDE than \(u_h\) is. In other words, this function provides a postprocessing step that improves the solution in a similar way one often obtains by extrapolating a sequence of solutions, explaining the origin of the function's name.

Note
The resulting field does not satisfy continuity requirements of the given finite elements if the algorithm outlined above is used. When you use continuous elements on grids with hanging nodes, please use the extrapolate function with an additional AffineConstraints argument, see below.
Since this function operates on patches of cells, it requires that the underlying grid is refined at least once for every coarse grid cell. If this is not the case, an exception will be raised.

◆ extrapolate() [2/2]

template<int dim, class InVector , class OutVector , int spacedim>
void FETools::extrapolate ( const DoFHandler< dim, spacedim > &  dof1,
const InVector &  z1,
const DoFHandler< dim, spacedim > &  dof2,
const AffineConstraints< typename OutVector::value_type > &  constraints,
OutVector &  z2 
)

Compute the patchwise extrapolation of a dof1 function z1 to a dof2 function z2. dof1 and dof2 need to be DoFHandler objects based on the same triangulation. constraints is a hanging node constraints object corresponding to dof2. This object is necessary when interpolating onto continuous elements on grids with hanging nodes (locally refined grids).

Otherwise, the function does the same as the other extrapolate function above (for which the documentation provides an extensive description of its operation).

◆ hierarchic_to_lexicographic_numbering()

template<int dim>
std::vector< unsigned int > FETools::hierarchic_to_lexicographic_numbering ( unsigned int  degree)

The numbering of the degrees of freedom in continuous finite elements is hierarchic, i.e. in such a way that we first number the vertex dofs, in the order of the vertices as defined by the triangulation, then the line dofs in the order and respecting the direction of the lines, then the dofs on quads, etc. However, we could have, as well, numbered them in a lexicographic way, i.e. with indices first running in x-direction, then in y-direction and finally in z-direction. Discontinuous elements of class FE_DGQ() are numbered in this way, for example.

This function returns a vector containing information about the lexicographic index each degree of freedom in the hierarchic numbering would have to a given degree of a continuous finite element.

◆ lexicographic_to_hierarchic_numbering()

template<int dim>
std::vector< unsigned int > FETools::lexicographic_to_hierarchic_numbering ( unsigned int  degree)

This is the reverse function to the above one, generating the map from the lexicographic to the hierarchical numbering for a given polynomial degree of a continuous finite element. All the remarks made about the above function are also valid here.

◆ get_fe_by_name()

template<int dim, int spacedim = dim>
std::unique_ptr< FiniteElement< dim, spacedim > > FETools::get_fe_by_name ( const std::string &  name)

Parse the name of a finite element and generate a finite element object accordingly. The parser ignores space characters between words (things matching the regular expression [A-Za-z0-9_]).

The name must be in the form which is returned by the FiniteElement::get_name function, where dimension template parameters <2> etc. can be omitted. Alternatively, the explicit number can be replaced by dim or d. If a number is given, it must match the template parameter of this function.

The names of FESystem elements follow the pattern FESystem[FE_Base1^p1-FE_Base2^p2] The powers p1 etc. may either be numbers or can be replaced by dim or d.

If no finite element can be reconstructed from this string, an exception of type FETools::ExcInvalidFEName is thrown.

The function returns a std::unique_ptr to a newly created finite element meaning the caller obtains ownership over the returned object.

Since the value of the template argument can't be deduced from the (string) argument given to this function, you have to explicitly specify it when you call this function.

This function knows about all the standard elements defined in the library. However, it doesn't by default know about elements that you may have defined in your program. To make your own elements known to this function, use the add_fe_name() function. This function does not work if one wants to get a codimension 1 finite element.

◆ add_fe_name()

template<int dim, int spacedim>
void FETools::add_fe_name ( const std::string &  name,
const FEFactoryBase< dim, spacedim > *  factory 
)

Extend the list of finite elements that can be generated by get_fe_by_name() by the one given as name. If get_fe_by_name() is later called with this name, it will use the object given as second argument to create a finite element object.

The format of the name parameter should include the name of a finite element. However, it is safe to use either the class name alone or to use the result of FiniteElement::get_name (which includes the space dimension as well as the polynomial degree), since everything after the first non- name character will be ignored.

The FEFactory object should be an object newly created with new. FETools will take ownership of this object and delete it once it is not used anymore.

In most cases, if you want objects of type MyFE be created whenever the name my_fe is given to get_fe_by_name, you will want the second argument to this function be of type FEFactory<MyFE>, but you can of course create your custom finite element factory class.

This function takes over ownership of the object given as second argument, i.e. you should never attempt to destroy it later on. The object will be deleted at the end of the program's lifetime.

If the name of the element is already in use, an exception is thrown. Thus, functionality of get_fe_by_name() can only be added, not changed.

Note
This function manipulates a global table (one table for each space dimension). It is thread safe in the sense that every access to this table is secured by a lock. Nevertheless, since each name can be added only once, user code has to make sure that only one thread adds a new element.

Note also that this table exists once for each space dimension. If you have a program that works with finite elements in different space dimensions (for example, step-4 does something like this), then you should call this function for each space dimension for which you want your finite element added to the map.

◆ hierarchic_to_lexicographic_numbering< 0 >()

template std::vector< unsigned int > FETools::hierarchic_to_lexicographic_numbering< 0 > ( unsigned int  )

◆ lexicographic_to_hierarchic_numbering< 0 >()

template std::vector< unsigned int > FETools::lexicographic_to_hierarchic_numbering< 0 > ( unsigned int  )