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

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

Inheritance diagram for ScalarLagrangePolynomialWedge< dim >:
[legend]

Public Member Functions

 ScalarLagrangePolynomialWedge (const unsigned int degree)
 
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
 
Tensor< 1, dim > compute_1st_derivative (const unsigned int i, const Point< dim > &p) const override
 
Tensor< 2, dim > compute_2nd_derivative (const unsigned int i, const Point< dim > &p) const override
 
Tensor< 3, dim > compute_3rd_derivative (const unsigned int i, const Point< dim > &p) const override
 
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
 
unsigned int n () const
 
virtual unsigned int degree () const
 
virtual std::size_t memory_consumption () const
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 

Private Attributes

const BarycentricPolynomials< 2 > poly_tri
 
const BarycentricPolynomials< 1 > poly_line
 
const unsigned int polynomial_degree
 
const unsigned int n_pols
 

Detailed Description

template<int dim>
class ScalarLagrangePolynomialWedge< dim >

Polynomials defined on wedge entities. This class is basis of FE_WedgeP.

The polynomials are created via a tensor product of a BarycentricPolynomials<2>::get_fe_p_basis(degree) and a BarycentricPolynomials<1>::get_fe_p_basis(degree), however, are re-numerated to better match the definition of FiniteElement.

Definition at line 76 of file polynomials_wedge.h.

Constructor & Destructor Documentation

◆ ScalarLagrangePolynomialWedge()

template<int dim>
ScalarLagrangePolynomialWedge< dim >::ScalarLagrangePolynomialWedge ( const unsigned int  degree)

Definition at line 44 of file polynomials_wedge.cc.

Member Function Documentation

◆ evaluate()

template<int dim>
void ScalarLagrangePolynomialWedge< 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 derivatives of the polynomials 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.

If you need values or derivatives of all 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.

Note
Currently, only the vectors values and grads are filled.

Implements ScalarPolynomialsBase< dim >.

Definition at line 114 of file polynomials_wedge.cc.

◆ compute_value()

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

Compute the value of the ith polynomial at unit point p.

Consider using evaluate() instead.

Implements ScalarPolynomialsBase< dim >.

Definition at line 55 of file polynomials_wedge.cc.

◆ compute_derivative()

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

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

Consider using evaluate() instead.

Template Parameters
orderThe order of the derivative.
Note
Currently, only implemented for first derivative.

Definition at line 187 of file polynomials_wedge.h.

◆ compute_1st_derivative()

template<int dim>
Tensor< 1, dim > ScalarLagrangePolynomialWedge< 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 >.

Definition at line 140 of file polynomials_wedge.cc.

◆ compute_2nd_derivative()

template<int dim>
Tensor< 2, dim > ScalarLagrangePolynomialWedge< 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.

Note
Not implemented yet.

Implements ScalarPolynomialsBase< dim >.

Definition at line 151 of file polynomials_wedge.cc.

◆ compute_3rd_derivative()

template<int dim>
Tensor< 3, dim > ScalarLagrangePolynomialWedge< 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.

Note
Not implemented yet.

Implements ScalarPolynomialsBase< dim >.

Definition at line 167 of file polynomials_wedge.cc.

◆ compute_4th_derivative()

template<int dim>
Tensor< 4, dim > ScalarLagrangePolynomialWedge< 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.

Note
Not implemented yet.

Implements ScalarPolynomialsBase< dim >.

Definition at line 183 of file polynomials_wedge.cc.

◆ compute_grad()

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

Compute the gradient of the ith polynomial at unit point p.

Consider using evaluate() instead.

Note
Not implemented yet.

Implements ScalarPolynomialsBase< dim >.

Definition at line 74 of file polynomials_wedge.cc.

◆ compute_grad_grad()

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

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

Consider using evaluate() instead.

Note
Not implemented yet.

Implements ScalarPolynomialsBase< dim >.

Definition at line 100 of file polynomials_wedge.cc.

◆ name()

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

Return the name of the space.

Implements ScalarPolynomialsBase< dim >.

Definition at line 199 of file polynomials_wedge.cc.

◆ clone()

template<int dim>
std::unique_ptr< ScalarPolynomialsBase< dim > > ScalarLagrangePolynomialWedge< 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 208 of file polynomials_wedge.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.

◆ memory_consumption()

template<int dim>
std::size_t ScalarPolynomialsBase< dim >::memory_consumption
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 39 of file scalar_polynomials_base.cc.

Member Data Documentation

◆ dimension

template<int dim>
constexpr unsigned int ScalarLagrangePolynomialWedge< dim >::dimension = dim
staticconstexpr

Make the dimension available to the outside.

Definition at line 82 of file polynomials_wedge.h.

◆ poly_tri

template<int dim>
const BarycentricPolynomials<2> ScalarLagrangePolynomialWedge< dim >::poly_tri
private

Scalar polynomials defined on a triangle.

Definition at line 174 of file polynomials_wedge.h.

◆ poly_line

template<int dim>
const BarycentricPolynomials<1> ScalarLagrangePolynomialWedge< dim >::poly_line
private

Scalar polynomials defined on a line.

Definition at line 179 of file polynomials_wedge.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: