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 Attributes | Protected Member Functions | Protected Attributes | Private Attributes | Friends | List of all members
TensorProductPolynomials< dim, PolynomialType > Class Template Reference

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

Inheritance diagram for TensorProductPolynomials< dim, PolynomialType >:
[legend]

Public Member Functions

template<class Pol >
 TensorProductPolynomials (const std::vector< Pol > &pols)
 
void output_indices (std::ostream &out) const
 
void set_numbering (const std::vector< unsigned int > &renumber)
 
const std::vector< unsigned int > & get_numbering () const
 
const std::vector< unsigned int > & get_numbering_inverse () const
 
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
 
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
 
std::string name () const override
 
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone () const override
 
virtual std::size_t memory_consumption () const override
 
std::vector< PolynomialType > get_underlying_polynomials () const
 
void output_indices (std::ostream &) const
 
void set_numbering (const std::vector< unsigned int > &)
 
double compute_value (const unsigned int, const Point< 0 > &) const
 
Tensor< 1, 0 > compute_grad (const unsigned int, const Point< 0 > &) const
 
Tensor< 2, 0 > compute_grad_grad (const unsigned int, const Point< 0 > &) const
 
void evaluate (const Point< 0 > &, std::vector< double > &, std::vector< Tensor< 1, 0 > > &, std::vector< Tensor< 2, 0 > > &, std::vector< Tensor< 3, 0 > > &, std::vector< Tensor< 4, 0 > > &) const
 
unsigned int n () const
 
virtual unsigned int degree () const
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 

Protected Member Functions

void compute_index (const unsigned int i, std::array< unsigned int, dim > &indices) const
 
void compute_index (const unsigned int, std::array< unsigned int, 0 > &) const
 

Protected Attributes

std::vector< PolynomialType > polynomials
 
std::vector< unsigned intindex_map
 
std::vector< unsigned intindex_map_inverse
 

Private Attributes

const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Friends

class TensorProductPolynomialsBubbles< dim >
 
class TensorProductPolynomialsConst< dim >
 

Detailed Description

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
class TensorProductPolynomials< dim, PolynomialType >

Tensor product of given polynomials.

Given a vector of n one-dimensional polynomials P1 to Pn, this class generates ndim polynomials of the form Qijk(x,y,z) = Pi(x)Pj(y)Pk(z). If the base polynomials are mutually orthogonal on the interval [-1,1] or [0,1], then the tensor product polynomials are orthogonal on [-1,1]dim or [0,1]dim, respectively.

Indexing is as follows: the order of dim-dimensional polynomials is x-coordinates running fastest, then y-coordinate, etc. The first few polynomials are thus P1(x)P1(y), P2(x)P1(y), P3(x)P1(y), ..., P1(x)P2(y), P2(x)P2(y), P3(x)P2(y), ... and likewise in 3d.

The output_indices() function prints the ordering of the dim-dimensional polynomials, i.e. for each polynomial in the polynomial space it gives the indices i,j,k of the one-dimensional polynomials in x,y and z direction. The ordering of the dim-dimensional polynomials can be changed by using the set_numbering() function.

Template Parameters
PolynomialTypeA class that satisfies the required interface for computing tensor products. Typical choices for this template argument are Polynomials::Polynomial and Polynomials::PiecewisePolynomial.

Definition at line 76 of file tensor_product_polynomials.h.

Constructor & Destructor Documentation

◆ TensorProductPolynomials()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
template<class Pol >
TensorProductPolynomials< dim, PolynomialType >::TensorProductPolynomials ( const std::vector< Pol > &  pols)

Constructor. pols is a vector of objects that should be derived or otherwise convertible to one-dimensional polynomial objects of type PolynomialType (template argument of class). It will be copied element by element into a protected member variable.

Member Function Documentation

◆ output_indices() [1/2]

template<int dim, typename PolynomialType >
void TensorProductPolynomials< dim, PolynomialType >::output_indices ( std::ostream &  out) const

Print the list of the indices to out.

Definition at line 112 of file tensor_product_polynomials.cc.

◆ set_numbering() [1/2]

template<int dim, typename PolynomialType >
void TensorProductPolynomials< dim, PolynomialType >::set_numbering ( const std::vector< unsigned int > &  renumber)

Set the ordering of the polynomials. Requires renumber.size()==n(). Stores a copy of renumber.

Definition at line 141 of file tensor_product_polynomials.cc.

◆ get_numbering()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
const std::vector< unsigned int > & TensorProductPolynomials< dim, PolynomialType >::get_numbering ( ) const

Give read access to the renumber vector.

◆ get_numbering_inverse()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
const std::vector< unsigned int > & TensorProductPolynomials< dim, PolynomialType >::get_numbering_inverse ( ) const

Give read access to the inverse renumber vector.

◆ evaluate() [1/2]

template<int dim, typename PolynomialType >
void TensorProductPolynomials< dim, PolynomialType >::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 the value and the first and second derivatives of each tensor product polynomial at unit_point.

The size of the vectors must either be equal 0 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 ScalarPolynomialsBase< dim >.

Definition at line 308 of file tensor_product_polynomials.cc.

◆ compute_value() [1/2]

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

Compute the value of the ith tensor product polynomial at unit_point. Here i is given in tensor product numbering.

Note, that using this function within a loop over all tensor product polynomials is not efficient, because then each point value of the underlying (one-dimensional) polynomials is (unnecessarily) computed several times. Instead use the evaluate() function with values.size()==n() to get the point values of all tensor polynomials all at once and in a much more efficient way.

Implements ScalarPolynomialsBase< dim >.

Definition at line 167 of file tensor_product_polynomials.cc.

◆ compute_derivative()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
template<int order>
Tensor< order, dim > TensorProductPolynomials< dim, PolynomialType >::compute_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const

Compute the orderth derivative of the ith tensor product polynomial at unit_point. Here i is given in tensor product numbering.

Note, that using this function within a loop over all tensor product polynomials is not efficient, because then each derivative value of the underlying (one-dimensional) polynomials is (unnecessarily) computed several times. Instead use the evaluate() function, see above, with the size of the appropriate parameter set to n() to get the point value of all tensor polynomials all at once and in a much more efficient way.

Template Parameters
orderThe derivative order.

◆ compute_1st_derivative()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
virtual Tensor< 1, dim > TensorProductPolynomials< dim, PolynomialType >::compute_1st_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

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

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

◆ compute_2nd_derivative()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
virtual Tensor< 2, dim > TensorProductPolynomials< dim, PolynomialType >::compute_2nd_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

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

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

◆ compute_3rd_derivative()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
virtual Tensor< 3, dim > TensorProductPolynomials< dim, PolynomialType >::compute_3rd_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

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

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

◆ compute_4th_derivative()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
virtual Tensor< 4, dim > TensorProductPolynomials< dim, PolynomialType >::compute_4th_derivative ( const unsigned int  i,
const Point< dim > &  p 
) const
overridevirtual

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

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

◆ compute_grad() [1/2]

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

Compute the grad of the ith tensor product polynomial at unit_point. Here i is given in tensor product numbering.

Note, that using this function within a loop over all tensor product polynomials is not efficient, because then each derivative value of the underlying (one-dimensional) polynomials is (unnecessarily) computed several times. Instead use the evaluate() function, see above, with grads.size()==n() to get the point value of all tensor polynomials all at once and in a much more efficient way.

Implements ScalarPolynomialsBase< dim >.

Definition at line 201 of file tensor_product_polynomials.cc.

◆ compute_grad_grad() [1/2]

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

Compute the second derivative (grad_grad) of the ith tensor product polynomial at unit_point. Here i is given in tensor product numbering.

Note, that using this function within a loop over all tensor product polynomials is not efficient, because then each derivative value of the underlying (one-dimensional) polynomials is (unnecessarily) computed several times. Instead use the evaluate() function, see above, with grad_grads.size()==n() to get the point value of all tensor polynomials all at once and in a much more efficient way.

Implements ScalarPolynomialsBase< dim >.

Definition at line 252 of file tensor_product_polynomials.cc.

◆ name()

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
std::string TensorProductPolynomials< dim, PolynomialType >::name ( ) const
overridevirtual

Return the name of the space, which is TensorProductPolynomials.

Implements ScalarPolynomialsBase< dim >.

◆ clone()

template<int dim, typename PolynomialType >
std::unique_ptr< ScalarPolynomialsBase< dim > > TensorProductPolynomials< dim, PolynomialType >::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_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 487 of file tensor_product_polynomials.cc.

◆ memory_consumption()

template<int dim, typename PolynomialType >
std::size_t TensorProductPolynomials< dim, PolynomialType >::memory_consumption
overridevirtual

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

Reimplemented from ScalarPolynomialsBase< dim >.

Definition at line 496 of file tensor_product_polynomials.cc.

◆ get_underlying_polynomials()

template<int dim, typename PolynomialType >
std::vector< PolynomialType > TensorProductPolynomials< dim, PolynomialType >::get_underlying_polynomials

Return a copy of the underlying one-dimensional polynomials given to the constructor of this class.

Definition at line 507 of file tensor_product_polynomials.cc.

◆ compute_index() [1/2]

template<int dim, typename PolynomialType >
void TensorProductPolynomials< dim, PolynomialType >::compute_index ( const unsigned int  i,
std::array< unsigned int, dim > &  indices 
) const
inlineprotected

Each tensor product polynomial i is a product of one-dimensional polynomials in each space direction. Compute the indices of these one- dimensional polynomials for each space direction, given the index i.

Definition at line 84 of file tensor_product_polynomials.cc.

◆ compute_index() [2/2]

void TensorProductPolynomials< 0, Polynomials::Polynomial< double > >::compute_index ( const unsigned int  ,
std::array< unsigned int, 0 > &   
) const
inlineprotected

Definition at line 100 of file tensor_product_polynomials.cc.

◆ output_indices() [2/2]

void TensorProductPolynomials< 0, Polynomials::Polynomial< double > >::output_indices ( std::ostream &  ) const

Definition at line 130 of file tensor_product_polynomials.cc.

◆ set_numbering() [2/2]

void TensorProductPolynomials< 0, Polynomials::Polynomial< double > >::set_numbering ( const std::vector< unsigned int > &  )

Definition at line 156 of file tensor_product_polynomials.cc.

◆ compute_value() [2/2]

double TensorProductPolynomials< 0, Polynomials::Polynomial< double > >::compute_value ( const unsigned int  ,
const Point< 0 > &   
) const

Definition at line 187 of file tensor_product_polynomials.cc.

◆ compute_grad() [2/2]

Tensor< 1, 0 > TensorProductPolynomials< 0, Polynomials::Polynomial< double > >::compute_grad ( const unsigned int  ,
const Point< 0 > &   
) const

Definition at line 238 of file tensor_product_polynomials.cc.

◆ compute_grad_grad() [2/2]

Tensor< 2, 0 > TensorProductPolynomials< 0, Polynomials::Polynomial< double > >::compute_grad_grad ( const unsigned int  ,
const Point< 0 > &   
) const

Definition at line 297 of file tensor_product_polynomials.cc.

◆ evaluate() [2/2]

void TensorProductPolynomials< 0, Polynomials::Polynomial< double > >::evaluate ( const Point< 0 > &  ,
std::vector< double > &  ,
std::vector< Tensor< 1, 0 > > &  ,
std::vector< Tensor< 2, 0 > > &  ,
std::vector< Tensor< 3, 0 > > &  ,
std::vector< Tensor< 4, 0 > > &   
) const

Definition at line 471 of file tensor_product_polynomials.cc.

◆ n()

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

Return the number of polynomials.

Definition at line 240 of file scalar_polynomials_base.h.

◆ degree()

template<int dim>
unsigned int ScalarPolynomialsBase< dim >::degree
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 249 of file scalar_polynomials_base.h.

Friends And Related Symbol Documentation

◆ TensorProductPolynomialsBubbles< dim >

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
friend class TensorProductPolynomialsBubbles< dim >
friend

TensorProductPolynomialsBubbles has a TensorProductPolynomials class so we declare it as a friend class.

Definition at line 278 of file tensor_product_polynomials.h.

◆ TensorProductPolynomialsConst< dim >

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
friend class TensorProductPolynomialsConst< dim >
friend

TensorProductPolynomialsConst has a TensorProductPolynomials class so we declare it as a friend class.

Definition at line 278 of file tensor_product_polynomials.h.

Member Data Documentation

◆ dimension

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
constexpr unsigned int TensorProductPolynomials< dim, PolynomialType >::dimension = dim
staticconstexpr

Access to the dimension of this object, for checking and automatic setting of dimension in other classes.

Definition at line 83 of file tensor_product_polynomials.h.

◆ polynomials

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
std::vector<PolynomialType> TensorProductPolynomials< dim, PolynomialType >::polynomials
protected

Copy of the vector pols of polynomials given to the constructor.

Definition at line 259 of file tensor_product_polynomials.h.

◆ index_map

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
std::vector<unsigned int> TensorProductPolynomials< dim, PolynomialType >::index_map
protected

Index map for reordering the polynomials.

Definition at line 264 of file tensor_product_polynomials.h.

◆ index_map_inverse

template<int dim, typename PolynomialType = Polynomials::Polynomial<double>>
std::vector<unsigned int> TensorProductPolynomials< dim, PolynomialType >::index_map_inverse
protected

Index map for reordering the polynomials.

Definition at line 269 of file tensor_product_polynomials.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 228 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 233 of file scalar_polynomials_base.h.


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