Reference documentation for deal.II version GIT relicensing-249-g48dc7357c7 2024-03-29 12:30: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 Member Functions | Private Attributes | List of all members
PolynomialsRaviartThomas< dim > Class Template Reference

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

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

Public Member Functions

 PolynomialsRaviartThomas (const unsigned int degree_normal, const unsigned int degree_tangential)
 
 PolynomialsRaviartThomas (const unsigned int k)
 
 PolynomialsRaviartThomas (const PolynomialsRaviartThomas &other)=default
 
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
 
std::vector< Point< dim > > get_polynomial_support_points () const
 
unsigned int n () const
 
unsigned int degree () const
 

Static Public Member Functions

static unsigned int n_polynomials (const unsigned int normal_degree, const unsigned int tangential_degree)
 
static unsigned int n_polynomials (const unsigned int degree)
 
static std::vector< unsigned intget_lexicographic_numbering (const unsigned int normal_degree, const unsigned int tangential_degree)
 

Private Attributes

const unsigned int normal_degree
 
const unsigned int tangential_degree
 
const AnisotropicPolynomials< dim > polynomial_space
 
std::vector< unsigned intlexicographic_to_hierarchic
 
std::vector< unsigned inthierarchic_to_lexicographic
 
std::array< std::vector< unsigned int >, dim > renumber_aniso
 
const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Detailed Description

template<int dim>
class PolynomialsRaviartThomas< dim >

This class implements the Hdiv-conforming, vector-valued Raviart-Thomas polynomials as described in the book by Brezzi and Fortin.

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

Definition at line 49 of file polynomials_raviart_thomas.h.

Constructor & Destructor Documentation

◆ PolynomialsRaviartThomas() [1/3]

template<int dim>
PolynomialsRaviartThomas< dim >::PolynomialsRaviartThomas ( const unsigned int  degree_normal,
const unsigned int  degree_tangential 
)

Constructor. Creates all basis functions for Raviart-Thomas polynomials of given degree in normal and tangential directions. The usual FE_RaviartThomas and FE_RaviartThomasNodal classes will use degree + 1 and degree in the two directions, respectively.

Definition at line 63 of file polynomials_raviart_thomas.cc.

◆ PolynomialsRaviartThomas() [2/3]

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

Constructor, using the common Raviart-Thomas space of degree k + 1 in normal direction and k in the tangential directions.

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

Definition at line 125 of file polynomials_raviart_thomas.cc.

◆ PolynomialsRaviartThomas() [3/3]

template<int dim>
PolynomialsRaviartThomas< dim >::PolynomialsRaviartThomas ( const PolynomialsRaviartThomas< dim > &  other)
default

Copy constructor.

Member Function Documentation

◆ evaluate()

template<int dim>
void PolynomialsRaviartThomas< 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 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.

Implements TensorPolynomialsBase< dim >.

Definition at line 133 of file polynomials_raviart_thomas.cc.

◆ name()

template<int dim>
std::string PolynomialsRaviartThomas< dim >::name ( ) const
overridevirtual

Return the name of the space, which is PolynomialsRaviartThomas.

Implements TensorPolynomialsBase< dim >.

Definition at line 226 of file polynomials_raviart_thomas.cc.

◆ n_polynomials() [1/2]

template<int dim>
unsigned int PolynomialsRaviartThomas< dim >::n_polynomials ( const unsigned int  normal_degree,
const unsigned int  tangential_degree 
)
static

Return the number of polynomials in the space without requiring to build an object of PolynomialsRaviartThomas. This is required by the FiniteElement classes.

Definition at line 244 of file polynomials_raviart_thomas.cc.

◆ n_polynomials() [2/2]

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

Variant of the n_polynomials() function taking only a single argument degree, assuming degree + 1 in the normal direction and degree in the tangential directions.

Definition at line 235 of file polynomials_raviart_thomas.cc.

◆ get_lexicographic_numbering()

template<int dim>
std::vector< unsigned int > PolynomialsRaviartThomas< dim >::get_lexicographic_numbering ( const unsigned int  normal_degree,
const unsigned int  tangential_degree 
)
static

Compute the lexicographic to hierarchic numbering underlying this class, computed as a free function.

Definition at line 256 of file polynomials_raviart_thomas.cc.

◆ clone()

template<int dim>
std::unique_ptr< TensorPolynomialsBase< dim > > PolynomialsRaviartThomas< 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_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 321 of file polynomials_raviart_thomas.cc.

◆ get_polynomial_support_points()

template<int dim>
std::vector< Point< dim > > PolynomialsRaviartThomas< dim >::get_polynomial_support_points ( ) const

Compute the generalized support points in the ordering used by the polynomial shape functions. Note that these points are not support points in the classical sense as the Lagrange polynomials of the different components have different points, which need to be combined in terms of Piola transforms.

Definition at line 330 of file polynomials_raviart_thomas.cc.

◆ n()

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

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

Member Data Documentation

◆ normal_degree

template<int dim>
const unsigned int PolynomialsRaviartThomas< dim >::normal_degree
private

The given degree in the normal direction.

Definition at line 142 of file polynomials_raviart_thomas.h.

◆ tangential_degree

template<int dim>
const unsigned int PolynomialsRaviartThomas< dim >::tangential_degree
private

The given degree in the tangential direction.

Definition at line 147 of file polynomials_raviart_thomas.h.

◆ polynomial_space

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

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

Definition at line 153 of file polynomials_raviart_thomas.h.

◆ lexicographic_to_hierarchic

template<int dim>
std::vector<unsigned int> PolynomialsRaviartThomas< dim >::lexicographic_to_hierarchic
private

Renumbering from lexicographic to hierarchic order.

Definition at line 158 of file polynomials_raviart_thomas.h.

◆ hierarchic_to_lexicographic

template<int dim>
std::vector<unsigned int> PolynomialsRaviartThomas< dim >::hierarchic_to_lexicographic
private

Renumbering from hierarchic to lexicographic order. Inverse of lexicographic_to_hierarchic.

Definition at line 164 of file polynomials_raviart_thomas.h.

◆ renumber_aniso

template<int dim>
std::array<std::vector<unsigned int>, dim> PolynomialsRaviartThomas< dim >::renumber_aniso
private

Renumbering from shifted polynomial spaces to lexicographic one.

Definition at line 169 of file polynomials_raviart_thomas.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 139 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 144 of file tensor_polynomials_base.h.


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