PTEMS  0.1.0-dev+git.81fd0e4
PolyTopic Element Method Solver
ptems::LegendrePolynomials< 1 > Struct Reference

Defines Legendre polynomials in the 1D reference element (the interval \([-1.1]\)) More...

#include <ptems/spaces/legendrepolynomials.hpp>

Inheritance diagram for ptems::LegendrePolynomials< 1 >:
Collaboration diagram for ptems::LegendrePolynomials< 1 >:

Static Public Member Functions

template<typename T >
static T * Polynomials (const ReferenceElement< 1 > &, const Vector< 1 > &pt, T *result, std::size_t polydeg)
 Computes the value of the first \(p+1\) Legendre polynomials on the reference interval \([-1,1]\). More...
 
template<typename T >
static FuncAndGradData< 1, T > * Gradient (const ReferenceElement< 1 > &, const Vector< 1 > &pt, FuncAndGradData< 1, T > *result, std::size_t polydeg)
 Computes the value of the first \(p+1\) Legendre polynomials, and their first order derivatives, on the reference interval \(4[-1.1]\). More...
 
template<typename T >
static FuncGradAndHessianData< 1, T > * Hessian (const ReferenceElement< 1 > &, const Vector< 1 > &pt, FuncGradAndHessianData< 1, T > *result, std::size_t polydeg)
 Computes the value of the first \(p+1\) Legendre polynomials, and their first and second order derivatives, on the reference interval \(4[-1.1]\). More...
 
template<typename T = double>
constexpr static T Integrate (const ReferenceElement< 1 > &, [[maybe_unused]] std::size_t polydeg, std::size_t i)
 Returns the value of integrating over the reference element the specified Legendre polynomial. More...
 
template<typename T = double>
constexpr static T Integrate (const ReferenceElement< 1 > &, [[maybe_unused]] std::size_t polydeg, std::size_t i, std::size_t j)
 Returns the value of integrating over the reference element the multiplication of the specified Legendre polynomials. More...
 
template<typename T = double>
static Vector< 1 > EstimateAnalyticity ([[maybe_unused]] const ReferenceElement< 1 > &refElement, std::size_t polydeg, const T *coefficients)
 Estimates the analyticity of a function with the current basis and specified coefficients. More...
 
static constexpr std::size_t NumberDoFs ([[maybe_unused]] const ReferenceElement< 1 > &refElement, std::size_t polydeg)
 Gets the number of degrees of freedom required for a polynomial of the specified degree on the specified reference element. More...
 

Detailed Description

Defines Legendre polynomials in the 1D reference element (the interval \([-1.1]\))

Member Function Documentation

◆ EstimateAnalyticity()

template<typename T = double>
static Vector<1> ptems::LegendrePolynomials< 1 >::EstimateAnalyticity ( [[maybe_unused] ] const ReferenceElement< 1 > &  refElement,
std::size_t  polydeg,
const T *  coefficients 
)
inlinestatic

Estimates the analyticity of a function with the current basis and specified coefficients.

Template Parameters
TThe type of the coefficients (default=double)
Parameters
refElementThe reference element
polydegThe polynomial degree of the function
coefficientsArray of coefficients for each basis. Must be length of NumberDoFs
Returns
Vector of analyticity for each dimension in range [0,1], or -1 if analyticity computation invalid

◆ Gradient()

template<typename T >
static FuncAndGradData<1,T>* ptems::LegendrePolynomials< 1 >::Gradient ( const ReferenceElement< 1 > &  ,
const Vector< 1 > &  pt,
FuncAndGradData< 1, T > *  result,
std::size_t  polydeg 
)
inlinestatic

Computes the value of the first \(p+1\) Legendre polynomials, and their first order derivatives, on the reference interval \(4[-1.1]\).

Template Parameters
TThe type of the result
Parameters
ptThe point to compute the values of the Legendre polynomial at
resultArray of at least polydeg+1 entries to store the result into
polydegThe polynomial degree of the Legendre basis functions
Returns
Pointer to the entry AFTER the last basis filled (=result+polydeg+1)

◆ Hessian()

template<typename T >
static FuncGradAndHessianData<1,T>* ptems::LegendrePolynomials< 1 >::Hessian ( const ReferenceElement< 1 > &  ,
const Vector< 1 > &  pt,
FuncGradAndHessianData< 1, T > *  result,
std::size_t  polydeg 
)
inlinestatic

Computes the value of the first \(p+1\) Legendre polynomials, and their first and second order derivatives, on the reference interval \(4[-1.1]\).

Template Parameters
TThe type of the result
Parameters
ptThe point to compute the values of the Legendre polynomial at
resultArray of at least polydeg+1 entries to store the result into
polydegThe polynomial degree of the Legendre basis functions
Returns
Pointer to the entry AFTER the last basis filled (=result+polydeg+1)

◆ Integrate() [1/2]

template<typename T = double>
constexpr static T ptems::LegendrePolynomials< 1 >::Integrate ( const ReferenceElement< 1 > &  ,
[[maybe_unused] ] std::size_t  polydeg,
std::size_t  i 
)
inlinestaticconstexpr

Returns the value of integrating over the reference element the specified Legendre polynomial.

Note
This is usually an analytical value
Template Parameters
TThe type of the result (default=double)
Parameters
polydegThe maximum polynomial degree of the Legendre polynomials
iThe index of the Legendre polynomial \(L_i\)
Returns
The integration over the reference interval \(\int_{-1}^1 L_i\dx\)

◆ Integrate() [2/2]

template<typename T = double>
constexpr static T ptems::LegendrePolynomials< 1 >::Integrate ( const ReferenceElement< 1 > &  ,
[[maybe_unused] ] std::size_t  polydeg,
std::size_t  i,
std::size_t  j 
)
inlinestaticconstexpr

Returns the value of integrating over the reference element the multiplication of the specified Legendre polynomials.

Note
This is usually an analytical value
Template Parameters
TThe type of the result (default=double)
Parameters
polydegThe maximum polynomial degree of the basis
iThe index of the first Legendre polynomial \(L_i\)
jThe index of the second Legendre polynomial \(L_j\)
Returns
The integration over the reference interval \(\int_{-1}^1 L_i L_j \dx\)

◆ NumberDoFs()

static constexpr std::size_t ptems::Polynomials< 1 >::NumberDoFs ( [[maybe_unused] ] const ReferenceElement< 1 > &  refElement,
std::size_t  polydeg 
)
inlinestaticconstexprinherited

Gets the number of degrees of freedom required for a polynomial of the specified degree on the specified reference element.

Parameters
refElementThe reference interval (parameter is ignored)
polydegThe degree of the polynomial
Returns
The number of degrees of freedom

◆ Polynomials()

template<typename T >
static T* ptems::LegendrePolynomials< 1 >::Polynomials ( const ReferenceElement< 1 > &  ,
const Vector< 1 > &  pt,
T *  result,
std::size_t  polydeg 
)
inlinestatic

Computes the value of the first \(p+1\) Legendre polynomials on the reference interval \([-1,1]\).

Template Parameters
TThe type of the result
Parameters
ptThe point to compute the values of the Legendre polynomial at
resultArray of at least polydeg+1 entries to store the result into
polydegThe polynomial degree of the Legendre basis functions
Returns
Pointer to the entry AFTER the last basis filled (=result+polydeg+1)

The documentation for this struct was generated from the following file: