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 Attributes | Private Attributes | List of all members
TensorProductPolynomialsBubbles< dim > Class Template Reference

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

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

Public Member Functions

template<class Pol >
 TensorProductPolynomialsBubbles (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
 
unsigned int n () const
 
std::string name () const override
 
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone () const override
 
double compute_value (const unsigned int, const Point< 0 > &) const
 
virtual unsigned int degree () const
 
virtual std::size_t memory_consumption () const
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 

Private Attributes

TensorProductPolynomials< dim > tensor_polys
 
std::vector< unsigned intindex_map
 
std::vector< unsigned intindex_map_inverse
 
const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Detailed Description

template<int dim>
class TensorProductPolynomialsBubbles< dim >

A class that represents a space of tensor product polynomials, augmented by \(dim\) (non-normalized) bubble functions of form \(\varphi_j(\mathbf x) = 2^{\text{degree}-1}\left(x_j-frac 12\right)^{\text{degree}-1} \left[\prod_{i=0}^{dim-1}(x_i(1-x_i))\right]\) for \(j=0,\ldots,dim-1\). If degree is one, then the first factor disappears and one receives the usual bubble function centered at the mid-point of the cell.

This class inherits most of its functionality from TensorProductPolynomials. The bubble enrichments are added for the last index.

Definition at line 52 of file tensor_product_polynomials_bubbles.h.

Constructor & Destructor Documentation

◆ TensorProductPolynomialsBubbles()

template<int dim>
template<class Pol >
TensorProductPolynomialsBubbles< dim >::TensorProductPolynomialsBubbles ( const std::vector< Pol > &  pols)

Constructor. pols is a vector of objects that should be derived or otherwise convertible to one-dimensional polynomial objects. It will be copied element by element into a private variable.

Member Function Documentation

◆ output_indices()

template<int dim>
void TensorProductPolynomialsBubbles< dim >::output_indices ( std::ostream &  out) const

Print the list of tensor_polys indices to out.

Definition at line 31 of file tensor_product_polynomials_bubbles.cc.

◆ set_numbering()

template<int dim>
void TensorProductPolynomialsBubbles< dim >::set_numbering ( const std::vector< unsigned int > &  renumber)

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

Definition at line 48 of file tensor_product_polynomials_bubbles.cc.

◆ get_numbering()

template<int dim>
const std::vector< unsigned int > & TensorProductPolynomialsBubbles< dim >::get_numbering ( ) const

Give read access to the renumber vector.

◆ get_numbering_inverse()

template<int dim>
const std::vector< unsigned int > & TensorProductPolynomialsBubbles< dim >::get_numbering_inverse ( ) const

Give read access to the inverse renumber vector.

◆ evaluate()

template<int dim>
void TensorProductPolynomialsBubbles< 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 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 259 of file tensor_product_polynomials_bubbles.cc.

◆ compute_value() [1/2]

template<int dim>
double TensorProductPolynomialsBubbles< 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 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 68 of file tensor_product_polynomials_bubbles.cc.

◆ compute_derivative()

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

Compute the order order 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.

◆ compute_1st_derivative()

template<int dim>
virtual Tensor< 1, dim > TensorProductPolynomialsBubbles< 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 > TensorProductPolynomialsBubbles< 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 > TensorProductPolynomialsBubbles< 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 > TensorProductPolynomialsBubbles< 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()

template<int dim>
Tensor< 1, dim > TensorProductPolynomialsBubbles< 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 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 109 of file tensor_product_polynomials_bubbles.cc.

◆ compute_grad_grad()

template<int dim>
Tensor< 2, dim > TensorProductPolynomialsBubbles< 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 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 155 of file tensor_product_polynomials_bubbles.cc.

◆ n()

template<int dim>
unsigned int TensorProductPolynomialsBubbles< dim >::n ( ) const

Return the number of tensor product polynomials plus the bubble enrichments. For n 1d polynomials this is ndim+1 if the maximum degree of the polynomials is one and ndim+dim otherwise.

◆ name()

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

Return the name of the space, which is TensorProductPolynomialsBubbles.

Implements ScalarPolynomialsBase< dim >.

◆ clone()

template<int dim>
std::unique_ptr< ScalarPolynomialsBase< dim > > TensorProductPolynomialsBubbles< 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 340 of file tensor_product_polynomials_bubbles.cc.

◆ compute_value() [2/2]

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

Definition at line 98 of file tensor_product_polynomials_bubbles.cc.

◆ 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 TensorProductPolynomialsBubbles< dim >::dimension = dim
staticconstexpr

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

Definition at line 59 of file tensor_product_polynomials_bubbles.h.

◆ tensor_polys

template<int dim>
TensorProductPolynomials<dim> TensorProductPolynomialsBubbles< dim >::tensor_polys
private

The TensorProductPolynomials object

Definition at line 230 of file tensor_product_polynomials_bubbles.h.

◆ index_map

template<int dim>
std::vector<unsigned int> TensorProductPolynomialsBubbles< dim >::index_map
private

Index map for reordering the polynomials.

Definition at line 235 of file tensor_product_polynomials_bubbles.h.

◆ index_map_inverse

template<int dim>
std::vector<unsigned int> TensorProductPolynomialsBubbles< dim >::index_map_inverse
private

Index map for reordering the polynomials.

Definition at line 240 of file tensor_product_polynomials_bubbles.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: