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
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
PolynomialsABF< dim > Class Template Reference

#include <deal.II/base/polynomials_abf.h>

Inheritance diagram for PolynomialsABF< dim >:
[legend]

Public Member Functions

 PolynomialsABF (const unsigned int k)
 
void evaluate (const Point< dim > &unit_point, std::vector< Tensor< 1, dim > > &values, std::vector< Tensor< 2, dim > > &grads, std::vector< Tensor< 3, dim > > &grad_grads, std::vector< Tensor< 4, dim > > &third_derivatives, std::vector< Tensor< 5, dim > > &fourth_derivatives) const override
 
std::string name () const override
 
virtual std::unique_ptr< TensorPolynomialsBase< dim > > clone () const override
 
unsigned int n () const
 
unsigned int degree () const
 

Static Public Member Functions

static unsigned int n_polynomials (const unsigned int degree)
 

Private Attributes

const AnisotropicPolynomials< dim > polynomial_space
 
Threads::Mutex mutex
 
std::vector< double > p_values
 
std::vector< Tensor< 1, dim > > p_grads
 
std::vector< Tensor< 2, dim > > p_grad_grads
 
std::vector< Tensor< 3, dim > > p_third_derivatives
 
std::vector< Tensor< 4, dim > > p_fourth_derivatives
 
const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Detailed Description

template<int dim>
class PolynomialsABF< dim >

This class implements the Hdiv-conforming, vector-valued Arnold-Boffi-Falk polynomials as described in the article by Arnold-Boffi- Falk: Quadrilateral H(div) finite elements, SIAM J. Numer. Anal. Vol.42, No.6, pp.2429-2451

The ABF polynomials are constructed such that the divergence is in the tensor product polynomial space Qk. Therefore, the polynomial order of each component must be two orders higher in the corresponding direction, yielding the polynomial spaces (Qk+2,k, Qk,k+2) and (Qk+2,k,k, Qk,k+2,k, Qk,k,k+2) in 2d and 3d, resp.

Definition at line 53 of file polynomials_abf.h.

Constructor & Destructor Documentation

◆ PolynomialsABF()

template<int dim>
PolynomialsABF< dim >::PolynomialsABF ( const unsigned int  k)

Constructor. Creates all basis functions for Raviart-Thomas polynomials of given degree.

  • k: the degree of the Raviart-Thomas-space, which is the degree of the largest tensor product polynomial space Qk contained.

Definition at line 50 of file polynomials_abf.cc.

Member Function Documentation

◆ evaluate()

template<int dim>
void PolynomialsABF< dim >::evaluate ( const Point< dim > &  unit_point,
std::vector< Tensor< 1, dim > > &  values,
std::vector< Tensor< 2, dim > > &  grads,
std::vector< Tensor< 3, dim > > &  grad_grads,
std::vector< Tensor< 4, dim > > &  third_derivatives,
std::vector< Tensor< 5, dim > > &  fourth_derivatives 
) const
overridevirtual

Compute the value and the first and second derivatives of each Raviart- Thomas polynomial at unit_point.

The size of the vectors must either be zero or equal n(). In the first case, the function will not compute these values.

If you need values or derivatives of all tensor product polynomials then use this function, rather than using any of the compute_value, compute_grad or compute_grad_grad functions, see below, in a loop over all tensor product polynomials.

Implements TensorPolynomialsBase< dim >.

Definition at line 64 of file polynomials_abf.cc.

◆ name()

template<int dim>
std::string PolynomialsABF< dim >::name
inlineoverridevirtual

Return the name of the space, which is ABF.

Implements TensorPolynomialsBase< dim >.

Definition at line 148 of file polynomials_abf.h.

◆ n_polynomials()

template<int dim>
unsigned int PolynomialsABF< dim >::n_polynomials ( const unsigned int  degree)
static

Return the number of polynomials in the space RT(degree) without requiring to build an object of PolynomialsABF. This is required by the FiniteElement classes.

Definition at line 158 of file polynomials_abf.cc.

◆ clone()

template<int dim>
std::unique_ptr< TensorPolynomialsBase< dim > > PolynomialsABF< dim >::clone
overridevirtual

A sort of virtual copy constructor, this function returns a copy of the polynomial space object. Derived classes need to override the function here in this base class and return an object of the same type as the derived class.

Some places in the library, for example the constructors of FE_PolyTensor, need to make copies of polynomial spaces without knowing their exact type. They do so through this function.

Implements TensorPolynomialsBase< dim >.

Definition at line 186 of file polynomials_abf.cc.

◆ n()

template<int dim>
unsigned int TensorPolynomialsBase< dim >::n
inlineinherited

Return the number of polynomials.

Definition at line 157 of file tensor_polynomials_base.h.

◆ degree()

template<int dim>
unsigned int TensorPolynomialsBase< dim >::degree
inlineinherited

Return the highest polynomial degree of polynomials represented by this class. A derived class may override this if its value is different from my_degree.

Definition at line 166 of file tensor_polynomials_base.h.

Member Data Documentation

◆ polynomial_space

template<int dim>
const AnisotropicPolynomials<dim> PolynomialsABF< dim >::polynomial_space
private

An object representing the polynomial space for a single component. We can re-use it for the other vector components by rotating the coordinates of the evaluation point.

Definition at line 112 of file polynomials_abf.h.

◆ mutex

template<int dim>
Threads::Mutex PolynomialsABF< dim >::mutex
mutableprivate

A mutex that guards the following scratch arrays.

Definition at line 117 of file polynomials_abf.h.

◆ p_values

template<int dim>
std::vector<double> PolynomialsABF< dim >::p_values
mutableprivate

Auxiliary memory.

Definition at line 122 of file polynomials_abf.h.

◆ p_grads

template<int dim>
std::vector<Tensor<1, dim> > PolynomialsABF< dim >::p_grads
mutableprivate

Auxiliary memory.

Definition at line 127 of file polynomials_abf.h.

◆ p_grad_grads

template<int dim>
std::vector<Tensor<2, dim> > PolynomialsABF< dim >::p_grad_grads
mutableprivate

Auxiliary memory.

Definition at line 132 of file polynomials_abf.h.

◆ p_third_derivatives

template<int dim>
std::vector<Tensor<3, dim> > PolynomialsABF< dim >::p_third_derivatives
mutableprivate

Auxiliary memory.

Definition at line 137 of file polynomials_abf.h.

◆ p_fourth_derivatives

template<int dim>
std::vector<Tensor<4, dim> > PolynomialsABF< dim >::p_fourth_derivatives
mutableprivate

Auxiliary memory.

Definition at line 142 of file polynomials_abf.h.

◆ polynomial_degree

template<int dim>
const unsigned int TensorPolynomialsBase< dim >::polynomial_degree
privateinherited

The highest polynomial degree of this functions represented by this object.

Definition at line 145 of file tensor_polynomials_base.h.

◆ n_pols

template<int dim>
const unsigned int TensorPolynomialsBase< dim >::n_pols
privateinherited

The number of polynomials represented by this object.

Definition at line 150 of file tensor_polynomials_base.h.


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