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 | Private Attributes | List of all members
ScalarPolynomialsBase< dim > Class Template Referenceabstract

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

Inheritance diagram for ScalarPolynomialsBase< dim >:
[legend]

Public Member Functions

 ScalarPolynomialsBase (const unsigned int deg, const unsigned int n_polynomials)
 
 ScalarPolynomialsBase (ScalarPolynomialsBase< dim > &&)=default
 
 ScalarPolynomialsBase (const ScalarPolynomialsBase< dim > &)=default
 
virtual ~ScalarPolynomialsBase ()=default
 
virtual void evaluate (const Point< dim > &unit_point, std::vector< double > &values, std::vector< Tensor< 1, dim > > &grads, std::vector< Tensor< 2, dim > > &grad_grads, std::vector< Tensor< 3, dim > > &third_derivatives, std::vector< Tensor< 4, dim > > &fourth_derivatives) const =0
 
virtual double compute_value (const unsigned int i, const Point< dim > &p) const =0
 
template<int order>
Tensor< order, dim > compute_derivative (const unsigned int i, const Point< dim > &p) const
 
virtual Tensor< 1, dim > compute_1st_derivative (const unsigned int i, const Point< dim > &p) const =0
 
virtual Tensor< 2, dim > compute_2nd_derivative (const unsigned int i, const Point< dim > &p) const =0
 
virtual Tensor< 3, dim > compute_3rd_derivative (const unsigned int i, const Point< dim > &p) const =0
 
virtual Tensor< 4, dim > compute_4th_derivative (const unsigned int i, const Point< dim > &p) const =0
 
virtual Tensor< 1, dim > compute_grad (const unsigned int, const Point< dim > &) const =0
 
virtual Tensor< 2, dim > compute_grad_grad (const unsigned int, const Point< dim > &) const =0
 
unsigned int n () const
 
virtual unsigned int degree () const
 
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone () const =0
 
virtual std::string name () const =0
 
virtual std::size_t memory_consumption () const
 

Private Attributes

const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Detailed Description

template<int dim>
class ScalarPolynomialsBase< dim >

This class provides a framework for the finite element polynomial classes for use with finite element classes that are derived from FE_Poly. 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_Poly.

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 63 of file scalar_polynomials_base.h.

Constructor & Destructor Documentation

◆ ScalarPolynomialsBase() [1/3]

template<int dim>
ScalarPolynomialsBase< dim >::ScalarPolynomialsBase ( 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 26 of file scalar_polynomials_base.cc.

◆ ScalarPolynomialsBase() [2/3]

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

Move constructor.

◆ ScalarPolynomialsBase() [3/3]

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

Copy constructor.

◆ ~ScalarPolynomialsBase()

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

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

Member Function Documentation

◆ evaluate()

template<int dim>
virtual void ScalarPolynomialsBase< dim >::evaluate ( const Point< dim > &  unit_point,
std::vector< double > &  values,
std::vector< Tensor< 1, dim > > &  grads,
std::vector< Tensor< 2, dim > > &  grad_grads,
std::vector< Tensor< 3, dim > > &  third_derivatives,
std::vector< Tensor< 4, 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.

If you need values or derivatives of all 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.

Implemented in PolynomialSpace< dim >, PolynomialSpace< dim - 1 >, PolynomialsAdini< dim >, BarycentricPolynomials< dim >, BarycentricPolynomials< 1 >, BarycentricPolynomials< 2 >, ScalarLagrangePolynomialPyramid< dim >, PolynomialsRannacherTurek< dim >, ScalarLagrangePolynomialWedge< dim >, TensorProductPolynomials< dim, PolynomialType >, TensorProductPolynomials< dim - 1 >, AnisotropicPolynomials< dim >, TensorProductPolynomialsBubbles< dim >, and TensorProductPolynomialsConst< dim >.

◆ compute_value()

template<int dim>
virtual double ScalarPolynomialsBase< dim >::compute_value ( const unsigned int  i,
const Point< dim > &  p 
) const
pure virtual

◆ compute_derivative()

template<int dim>
template<int order>
Tensor< order, dim > ScalarPolynomialsBase< dim >::compute_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
inline

Compute the orderth derivative of the ith polynomial at unit point p.

Consider using evaluate() instead.

Template Parameters
orderThe order of the derivative.

Definition at line 259 of file scalar_polynomials_base.h.

◆ compute_1st_derivative()

template<int dim>
virtual Tensor< 1, dim > ScalarPolynomialsBase< dim >::compute_1st_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
pure virtual

◆ compute_2nd_derivative()

template<int dim>
virtual Tensor< 2, dim > ScalarPolynomialsBase< dim >::compute_2nd_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
pure virtual

◆ compute_3rd_derivative()

template<int dim>
virtual Tensor< 3, dim > ScalarPolynomialsBase< dim >::compute_3rd_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
pure virtual

◆ compute_4th_derivative()

template<int dim>
virtual Tensor< 4, dim > ScalarPolynomialsBase< dim >::compute_4th_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
pure virtual

◆ compute_grad()

template<int dim>
virtual Tensor< 1, dim > ScalarPolynomialsBase< dim >::compute_grad ( const unsigned int  ,
const Point< dim > &   
) const
pure virtual

◆ compute_grad_grad()

template<int dim>
virtual Tensor< 2, dim > ScalarPolynomialsBase< dim >::compute_grad_grad ( const unsigned int  ,
const Point< dim > &   
) const
pure virtual

◆ n()

template<int dim>
unsigned int ScalarPolynomialsBase< dim >::n
inline

Return the number of polynomials.

Definition at line 240 of file scalar_polynomials_base.h.

◆ degree()

template<int dim>
unsigned int ScalarPolynomialsBase< dim >::degree
inlinevirtual

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.

Reimplemented in PolynomialsP< dim >.

Definition at line 249 of file scalar_polynomials_base.h.

◆ clone()

template<int dim>
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > ScalarPolynomialsBase< 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_Poly, need to make copies of polynomial spaces without knowing their exact type. They do so through this function.

Implemented in PolynomialSpace< dim >, PolynomialSpace< dim - 1 >, PolynomialsAdini< dim >, BarycentricPolynomials< dim >, BarycentricPolynomials< 1 >, BarycentricPolynomials< 2 >, PolynomialsP< dim >, ScalarLagrangePolynomialPyramid< dim >, PolynomialsRannacherTurek< dim >, ScalarLagrangePolynomialWedge< dim >, TensorProductPolynomials< dim, PolynomialType >, TensorProductPolynomials< dim - 1 >, AnisotropicPolynomials< dim >, TensorProductPolynomialsBubbles< dim >, and TensorProductPolynomialsConst< dim >.

◆ name()

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

◆ memory_consumption()

template<int dim>
std::size_t ScalarPolynomialsBase< dim >::memory_consumption
virtual

Return an estimate (in bytes) for the memory consumption of this object.

Reimplemented in BarycentricPolynomials< dim >, BarycentricPolynomials< 1 >, BarycentricPolynomials< 2 >, TensorProductPolynomials< dim, PolynomialType >, and TensorProductPolynomials< dim - 1 >.

Definition at line 39 of file scalar_polynomials_base.cc.

Member Data Documentation

◆ polynomial_degree

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

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

Definition at line 228 of file scalar_polynomials_base.h.

◆ n_pols

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

The number of polynomials represented by this object.

Definition at line 233 of file scalar_polynomials_base.h.


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