Reference documentation for deal.II version GIT relicensing-399-g79d89019c5 2024-04-16 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 | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
AnisotropicPolynomials< dim > Class Template Reference

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

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

Public Member Functions

 AnisotropicPolynomials (const std::vector< std::vector< Polynomials::Polynomial< double > > > &base_polynomials)
 
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
 
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
 
virtual std::size_t memory_consumption () const
 

Private 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
 
unsigned int get_n_tensor_pols (const std::vector< std::vector< Polynomials::Polynomial< double > > > &)
 

Static Private Member Functions

static unsigned int get_n_tensor_pols (const std::vector< std::vector< Polynomials::Polynomial< double > > > &pols)
 

Private Attributes

const std::vector< std::vector< Polynomials::Polynomial< double > > > polynomials
 
const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Detailed Description

template<int dim>
class AnisotropicPolynomials< dim >

Anisotropic tensor product of given polynomials.

Given one-dimensional polynomials \(P^x_1(x), P^x_2(x), \ldots\) in \(x\)-direction, \(P^y_1(y), P^y_2(y), \ldots\) in \(y\)-direction, and so on, this class generates polynomials of the form \(Q_{ijk}(x,y,z) = P^x_i(x)P^y_j(y)P^z_k(z)\). (With obvious generalization if dim is in fact only 2. If dim is in fact only 1, then the result is simply the same set of one-dimensional polynomials passed to the constructor.)

If the elements of each set of base polynomials are mutually orthogonal on the interval \([-1,1]\) or \([0,1]\), then the tensor product polynomials are orthogonal on \([-1,1]^d\) or \([0,1]^d\), respectively.

The resulting dim-dimensional tensor product polynomials are ordered as follows: We iterate over the \(x\) coordinates running fastest, then the \(y\) coordinate, etc. For example, for dim==2, the first few polynomials are thus \(P^x_1(x)P^y_1(y)\), \(P^x_2(x)P^y_1(y)\), \(P^x_3(x)P^y_1(y)\), ..., \(P^x_1(x)P^y_2(y)\), \(P^x_2(x)P^y_2(y)\), \(P^x_3(x)P^y_2(y)\), etc.

Definition at line 321 of file tensor_product_polynomials.h.

Constructor & Destructor Documentation

◆ AnisotropicPolynomials()

template<int dim>
AnisotropicPolynomials< dim >::AnisotropicPolynomials ( const std::vector< std::vector< Polynomials::Polynomial< double > > > &  base_polynomials)

Constructor. base_polynomials is a table of one-dimensional polynomials. The number of rows in this table (the first index when indexing into base_polynomials) needs to be equal to the space dimension, with the elements of each row (i.e., the second index) giving the polynomials that shall be used in this particular coordinate direction.

Since we want to build anisotropic polynomials, the dim sets of polynomials passed in as arguments may of course be different, and may also vary in number.

The number of tensor product polynomials is Nx*Ny*Nz, or with terms dropped if the number of space dimensions is less than 3.

Definition at line 642 of file tensor_product_polynomials.cc.

Member Function Documentation

◆ evaluate() [1/2]

template<int dim>
void AnisotropicPolynomials< 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 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 this->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 819 of file tensor_product_polynomials.cc.

◆ compute_value() [1/2]

template<int dim>
double AnisotropicPolynomials< dim >::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 compute function, see above, with values.size()==this->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 701 of file tensor_product_polynomials.cc.

◆ compute_derivative()

template<int dim>
template<int order>
Tensor< order, dim > AnisotropicPolynomials< dim >::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>
virtual Tensor< 1, dim > AnisotropicPolynomials< dim >::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>
virtual Tensor< 2, dim > AnisotropicPolynomials< dim >::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>
virtual Tensor< 3, dim > AnisotropicPolynomials< dim >::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>
virtual Tensor< 4, dim > AnisotropicPolynomials< dim >::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>
Tensor< 1, dim > AnisotropicPolynomials< dim >::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 compute function, see above, with grads.size()==this->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 730 of file tensor_product_polynomials.cc.

◆ compute_grad_grad() [1/2]

template<int dim>
Tensor< 2, dim > AnisotropicPolynomials< dim >::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 compute function, see above, with grad_grads.size()==this->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 771 of file tensor_product_polynomials.cc.

◆ name()

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

Return the name of the space, which is AnisotropicPolynomials.

Implements ScalarPolynomialsBase< dim >.

◆ clone()

template<int dim>
std::unique_ptr< ScalarPolynomialsBase< dim > > AnisotropicPolynomials< 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 948 of file tensor_product_polynomials.cc.

◆ compute_index() [1/2]

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

Each tensor product polynomial \(p_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 661 of file tensor_product_polynomials.cc.

◆ get_n_tensor_pols() [1/2]

template<int dim>
unsigned int AnisotropicPolynomials< dim >::get_n_tensor_pols ( const std::vector< std::vector< Polynomials::Polynomial< double > > > &  pols)
staticprivate

Given the input to the constructor, compute n_pols.

Definition at line 923 of file tensor_product_polynomials.cc.

◆ compute_index() [2/2]

void AnisotropicPolynomials< 0 >::compute_index ( const unsigned int  ,
std::array< unsigned int, 0 > &   
) const
private

Definition at line 691 of file tensor_product_polynomials.cc.

◆ compute_value() [2/2]

double AnisotropicPolynomials< 0 >::compute_value ( const unsigned int  ,
const Point< 0 > &   
) const

Definition at line 718 of file tensor_product_polynomials.cc.

◆ compute_grad() [2/2]

Tensor< 1, 0 > AnisotropicPolynomials< 0 >::compute_grad ( const unsigned int  ,
const Point< 0 > &   
) const

Definition at line 759 of file tensor_product_polynomials.cc.

◆ compute_grad_grad() [2/2]

Tensor< 2, 0 > AnisotropicPolynomials< 0 >::compute_grad_grad ( const unsigned int  ,
const Point< 0 > &   
) const

Definition at line 807 of file tensor_product_polynomials.cc.

◆ evaluate() [2/2]

void AnisotropicPolynomials< 0 >::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 909 of file tensor_product_polynomials.cc.

◆ get_n_tensor_pols() [2/2]

unsigned int AnisotropicPolynomials< 0 >::get_n_tensor_pols ( const std::vector< std::vector< Polynomials::Polynomial< double > > > &  )
private

Definition at line 936 of file tensor_product_polynomials.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

◆ polynomials

template<int dim>
const std::vector<std::vector<Polynomials::Polynomial<double> > > AnisotropicPolynomials< dim >::polynomials
private

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

Definition at line 471 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 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: