Defines Legendre polynomials in the 1D reference element (the interval \([-1.1]\))
More...
#include <ptems/spaces/legendrepolynomials.hpp>
|
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...
|
|
Defines Legendre polynomials in the 1D reference element (the interval \([-1.1]\))
◆ EstimateAnalyticity()
template<typename T = double>
Estimates the analyticity of a function with the current basis and specified coefficients.
- Template Parameters
-
T | The type of the coefficients (default=double) |
- Parameters
-
refElement | The reference element |
polydeg | The polynomial degree of the function |
coefficients | Array 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()
Computes the value of the first \(p+1\) Legendre polynomials, and their first order derivatives, on the reference interval \(4[-1.1]\).
- Template Parameters
-
- Parameters
-
pt | The point to compute the values of the Legendre polynomial at |
result | Array of at least polydeg+1 entries to store the result into |
polydeg | The polynomial degree of the Legendre basis functions |
- Returns
- Pointer to the entry AFTER the last basis filled (
=result+polydeg+1
)
◆ Hessian()
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
-
- Parameters
-
pt | The point to compute the values of the Legendre polynomial at |
result | Array of at least polydeg+1 entries to store the result into |
polydeg | The 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>
Returns the value of integrating over the reference element the specified Legendre polynomial.
- Note
- This is usually an analytical value
- Template Parameters
-
T | The type of the result (default=double) |
- Parameters
-
polydeg | The maximum polynomial degree of the Legendre polynomials |
i | The 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>
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
-
T | The type of the result (default=double) |
- Parameters
-
polydeg | The maximum polynomial degree of the basis |
i | The index of the first Legendre polynomial \(L_i\) |
j | The index of the second Legendre polynomial \(L_j\) |
- Returns
- The integration over the reference interval \(\int_{-1}^1 L_i L_j \dx\)
◆ NumberDoFs()
|
inlinestaticconstexprinherited |
Gets the number of degrees of freedom required for a polynomial of the specified degree on the specified reference element.
- Parameters
-
refElement | The reference interval (parameter is ignored) |
polydeg | The degree of the polynomial |
- Returns
- The number of degrees of freedom
◆ Polynomials()
Computes the value of the first \(p+1\) Legendre polynomials on the reference interval \([-1,1]\).
- Template Parameters
-
- Parameters
-
pt | The point to compute the values of the Legendre polynomial at |
result | Array of at least polydeg+1 entries to store the result into |
polydeg | The 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:
- ptems/spaces/legendrepolynomials.hpp