Reference documentation for deal.II version 9.6.0
|
#include <deal.II/base/tensor_product_polynomials.h>
Public Member Functions | |
AnisotropicPolynomials (const std::vector< std::vector< Polynomials::Polynomial< double > > > &base_polynomials) | |
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 |
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 |
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 |
std::vector< unsigned int > | index_map |
std::vector< unsigned int > | index_map_inverse |
const unsigned int | polynomial_degree |
const unsigned int | n_pols |
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.
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 684 of file tensor_product_polynomials.cc.
void AnisotropicPolynomials< dim >::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 156 of file tensor_product_polynomials.cc.
|
inline |
Give read access to the renumber vector.
Definition at line 138 of file tensor_product_polynomials.cc.
|
inline |
Give read access to the inverse renumber vector.
Definition at line 147 of file tensor_product_polynomials.cc.
|
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 871 of file tensor_product_polynomials.cc.
|
overridevirtual |
Compute the value of the i
th 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 753 of file tensor_product_polynomials.cc.
Tensor< order, dim > AnisotropicPolynomials< dim >::compute_derivative | ( | const unsigned int | i, |
const Point< dim > & | p ) const |
Compute the order
th derivative of the i
th 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.
order | The derivative order. |
|
overridevirtual |
Compute the first derivative of the i
th polynomial at unit point p
.
Consider using evaluate() instead.
Implements ScalarPolynomialsBase< dim >.
|
overridevirtual |
Compute the second derivative of the i
th polynomial at unit point p
.
Consider using evaluate() instead.
Implements ScalarPolynomialsBase< dim >.
|
overridevirtual |
Compute the third derivative of the i
th polynomial at unit point p
.
Consider using evaluate() instead.
Implements ScalarPolynomialsBase< dim >.
|
overridevirtual |
Compute the fourth derivative of the i
th polynomial at unit point p
.
Consider using evaluate() instead.
Implements ScalarPolynomialsBase< dim >.
|
overridevirtual |
Compute the grad of the i
th 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 782 of file tensor_product_polynomials.cc.
|
overridevirtual |
Compute the second derivative (grad_grad) of the i
th 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 823 of file tensor_product_polynomials.cc.
|
overridevirtual |
Return the name of the space, which is AnisotropicPolynomials
.
Implements ScalarPolynomialsBase< dim >.
|
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 1000 of file tensor_product_polynomials.cc.
|
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 713 of file tensor_product_polynomials.cc.
|
staticprivate |
Given the input to the constructor, compute n_pols
.
Definition at line 975 of file tensor_product_polynomials.cc.
void AnisotropicPolynomials< 0 >::set_numbering | ( | const std::vector< unsigned int > & | ) |
Definition at line 171 of file tensor_product_polynomials.cc.
|
private |
Definition at line 743 of file tensor_product_polynomials.cc.
double AnisotropicPolynomials< 0 >::compute_value | ( | const unsigned int | , |
const Point< 0 > & | ) const |
Definition at line 770 of file tensor_product_polynomials.cc.
Tensor< 1, 0 > AnisotropicPolynomials< 0 >::compute_grad | ( | const unsigned int | , |
const Point< 0 > & | ) const |
Definition at line 811 of file tensor_product_polynomials.cc.
Tensor< 2, 0 > AnisotropicPolynomials< 0 >::compute_grad_grad | ( | const unsigned int | , |
const Point< 0 > & | ) const |
Definition at line 859 of file tensor_product_polynomials.cc.
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 961 of file tensor_product_polynomials.cc.
|
private |
Definition at line 988 of file tensor_product_polynomials.cc.
|
inlineinherited |
Return the number of polynomials.
Definition at line 239 of file scalar_polynomials_base.h.
|
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.
|
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.
|
private |
Copy of the vector pols
of polynomials given to the constructor.
Definition at line 490 of file tensor_product_polynomials.h.
|
private |
Index map for reordering the polynomials.
Definition at line 495 of file tensor_product_polynomials.h.
|
private |
Index map for reordering the polynomials.
Definition at line 500 of file tensor_product_polynomials.h.
|
privateinherited |
The highest polynomial degree of this functions represented by this object.
Definition at line 227 of file scalar_polynomials_base.h.
|
privateinherited |
The number of polynomials represented by this object.
Definition at line 232 of file scalar_polynomials_base.h.