Reference documentation for deal.II version GIT relicensing-245-g36f19064f7 2024-03-29 07:20:02+00:00
\(\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 | Private Attributes | List of all members
TensorPolynomialsBase< dim > Class Template Referenceabstract

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

Inheritance diagram for TensorPolynomialsBase< dim >:
Inheritance graph
[legend]

Public Member Functions

 TensorPolynomialsBase (const unsigned int deg, const unsigned int n_polynomials)
 
 TensorPolynomialsBase (TensorPolynomialsBase< dim > &&)=default
 
 TensorPolynomialsBase (const TensorPolynomialsBase< dim > &)=default
 
virtual ~TensorPolynomialsBase ()=default
 
virtual 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 =0
 
unsigned int n () const
 
unsigned int degree () const
 
virtual std::unique_ptr< TensorPolynomialsBase< dim > > clone () const =0
 
virtual std::string name () const =0
 

Private Attributes

const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Detailed Description

template<int dim>
class TensorPolynomialsBase< dim >

This class provides a framework for the finite element polynomial classes for use with finite element classes that are derived from FE_PolyTensor. An object of this type (or rather of a type derived from this class) is stored as a member variable in each object of type FE_PolyTensor.

Deriving classes

Any derived class must provide the most basic properties for shape functions evaluated on the reference cell. This includes, but is not limited to, implementing the evaluate(), name(), and clone() member functions. These functions are necessary to store the most basic information of how the polynomials in the derived class evaluate at a given point on the reference cell. More information on each function can be found in the corresponding function's documentation.

Some classes that derive from this class include

Definition at line 61 of file tensor_polynomials_base.h.

Constructor & Destructor Documentation

◆ TensorPolynomialsBase() [1/3]

template<int dim>
TensorPolynomialsBase< dim >::TensorPolynomialsBase ( const unsigned int  deg,
const unsigned int  n_polynomials 
)

Constructor. This takes the degree of the space, deg from the finite element class, and n, the number of polynomials for the space.

Definition at line 25 of file tensor_polynomials_base.cc.

◆ TensorPolynomialsBase() [2/3]

template<int dim>
TensorPolynomialsBase< dim >::TensorPolynomialsBase ( TensorPolynomialsBase< dim > &&  )
default

Move constructor.

◆ TensorPolynomialsBase() [3/3]

template<int dim>
TensorPolynomialsBase< dim >::TensorPolynomialsBase ( const TensorPolynomialsBase< dim > &  )
default

Copy constructor.

◆ ~TensorPolynomialsBase()

template<int dim>
virtual TensorPolynomialsBase< dim >::~TensorPolynomialsBase ( )
virtualdefault

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

Member Function Documentation

◆ evaluate()

template<int dim>
virtual void TensorPolynomialsBase< 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
pure virtual

Compute the value and the derivatives of the polynomials 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.

Implemented in PolynomialsABF< dim >, PolynomialsBDM< dim >, PolynomialsBernardiRaugel< dim >, PolynomialsNedelec< dim >, PolynomialsRaviartThomas< dim >, and PolynomialsRT_Bubbles< dim >.

◆ n()

template<int dim>
unsigned int TensorPolynomialsBase< dim >::n ( ) const
inline

Return the number of polynomials.

Definition at line 151 of file tensor_polynomials_base.h.

◆ degree()

template<int dim>
unsigned int TensorPolynomialsBase< dim >::degree ( ) const
inline

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 160 of file tensor_polynomials_base.h.

◆ clone()

template<int dim>
virtual std::unique_ptr< TensorPolynomialsBase< dim > > TensorPolynomialsBase< dim >::clone ( ) const
pure virtual

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.

Implemented in PolynomialsABF< dim >, PolynomialsBDM< dim >, PolynomialsBernardiRaugel< dim >, PolynomialsNedelec< dim >, PolynomialsRaviartThomas< dim >, and PolynomialsRT_Bubbles< dim >.

◆ name()

template<int dim>
virtual std::string TensorPolynomialsBase< dim >::name ( ) const
pure virtual

Member Data Documentation

◆ polynomial_degree

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

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

Definition at line 139 of file tensor_polynomials_base.h.

◆ n_pols

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

The number of polynomials represented by this object.

Definition at line 144 of file tensor_polynomials_base.h.


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