Reference documentation for deal.II version GIT relicensing-218-g1f873f3929 2024-03-28 15:00: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 | Static Public Attributes | Private Attributes | List of all members
PolynomialsRannacherTurek< dim > Class Template Reference

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

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

Public Member Functions

 PolynomialsRannacherTurek ()
 
double compute_value (const unsigned int i, const Point< dim > &p) const override
 
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 override
 
virtual Tensor< 2, dim > compute_2nd_derivative (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 3, dim > compute_3rd_derivative (const unsigned int i, const Point< dim > &p) const override
 
virtual Tensor< 4, dim > compute_4th_derivative (const unsigned int i, const Point< dim > &p) const override
 
Tensor< 1, dim > compute_grad (const unsigned int i, const Point< dim > &p) const override
 
Tensor< 2, dim > compute_grad_grad (const unsigned int i, const Point< dim > &p) const override
 
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 override
 
std::string name () const override
 
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone () const override
 
unsigned int n () const
 
virtual unsigned int degree () const
 
virtual std::size_t memory_consumption () const
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 

Private Attributes

const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Detailed Description

template<int dim>
class PolynomialsRannacherTurek< dim >

Basis for polynomial space on the unit square used for lowest order Rannacher Turek element.

The i-th basis function is the dual basis element corresponding to the dof which evaluates the function's mean value across the i-th face. The numbering can be found in GeometryInfo.

Definition at line 41 of file polynomials_rannacher_turek.h.

Constructor & Destructor Documentation

◆ PolynomialsRannacherTurek()

Constructor, checking that the basis is implemented in this dimension.

Definition at line 25 of file polynomials_rannacher_turek.cc.

Member Function Documentation

◆ compute_value()

template<int dim>
double PolynomialsRannacherTurek< dim >::compute_value ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Value of basis function i at p.

Implements ScalarPolynomialsBase< dim >.

Definition at line 35 of file polynomials_rannacher_turek.cc.

◆ compute_derivative()

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

order-th of basis function i at p.

Consider using evaluate() instead.

Definition at line 239 of file polynomials_rannacher_turek.h.

◆ compute_1st_derivative()

template<int dim>
Tensor< 1, dim > PolynomialsRannacherTurek< dim >::compute_1st_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
inlineoverridevirtual

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

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

Definition at line 250 of file polynomials_rannacher_turek.h.

◆ compute_2nd_derivative()

template<int dim>
Tensor< 2, dim > PolynomialsRannacherTurek< dim >::compute_2nd_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
inlineoverridevirtual

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

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

Definition at line 261 of file polynomials_rannacher_turek.h.

◆ compute_3rd_derivative()

template<int dim>
Tensor< 3, dim > PolynomialsRannacherTurek< dim >::compute_3rd_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
inlineoverridevirtual

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

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

Definition at line 272 of file polynomials_rannacher_turek.h.

◆ compute_4th_derivative()

template<int dim>
Tensor< 4, dim > PolynomialsRannacherTurek< dim >::compute_4th_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
inlineoverridevirtual

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

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

Definition at line 283 of file polynomials_rannacher_turek.h.

◆ compute_grad()

template<int dim>
Tensor< 1, dim > PolynomialsRannacherTurek< dim >::compute_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Gradient of basis function i at p.

Implements ScalarPolynomialsBase< dim >.

Definition at line 68 of file polynomials_rannacher_turek.cc.

◆ compute_grad_grad()

template<int dim>
Tensor< 2, dim > PolynomialsRannacherTurek< dim >::compute_grad_grad ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

Gradient of gradient of basis function i at p.

Implements ScalarPolynomialsBase< dim >.

Definition at line 113 of file polynomials_rannacher_turek.cc.

◆ evaluate()

template<int dim>
void PolynomialsRannacherTurek< 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
overridevirtual

Compute values and derivatives of all basis functions at unit_point.

Size of the vectors must be either equal to the number of polynomials or zero. A size of zero means that we are not computing the vector entries.

Implements ScalarPolynomialsBase< dim >.

Definition at line 162 of file polynomials_rannacher_turek.cc.

◆ name()

template<int dim>
std::string PolynomialsRannacherTurek< dim >::name ( ) const
inlineoverridevirtual

Return the name of the space, which is RannacherTurek.

Implements ScalarPolynomialsBase< dim >.

Definition at line 294 of file polynomials_rannacher_turek.h.

◆ clone()

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

Implements ScalarPolynomialsBase< dim >.

Definition at line 211 of file polynomials_rannacher_turek.cc.

◆ n()

template<int dim>
unsigned int ScalarPolynomialsBase< dim >::n ( ) const
inlineinherited

Return the number of polynomials.

Definition at line 239 of file scalar_polynomials_base.h.

◆ degree()

template<int dim>
unsigned int ScalarPolynomialsBase< dim >::degree ( ) const
inlinevirtualinherited

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

◆ memory_consumption()

template<int dim>
std::size_t ScalarPolynomialsBase< dim >::memory_consumption ( ) const
virtualinherited

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 38 of file scalar_polynomials_base.cc.

Member Data Documentation

◆ dimension

template<int dim>
constexpr unsigned int PolynomialsRannacherTurek< dim >::dimension = dim
staticconstexpr

Dimension we are working in.

Definition at line 47 of file polynomials_rannacher_turek.h.

◆ polynomial_degree

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

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

Definition at line 227 of file scalar_polynomials_base.h.

◆ n_pols

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

The number of polynomials represented by this object.

Definition at line 232 of file scalar_polynomials_base.h.


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