Go to the documentation of this file. 1 #pragma once // Ensure that file is included only once in a single compilation.
24 template <
unsigned int poly_deg>
30 static constexpr
unsigned int n_fun() {
return degree() + 1; }
34 static constexpr
unsigned int dim() {
return 1U; }
38 static constexpr
unsigned int degree() {
return poly_deg; }
43 template <
typename return_t,
typename input_t>
44 static constexpr return_t
clenshaw_alpha(
const unsigned int k,
const input_t& x)
46 return (((return_t)(2 * k + 1)) / ((return_t)(k + 1))) * ((return_t)x);
51 template <
typename return_t,
typename input_t>
52 static constexpr return_t
clenshaw_beta(
const unsigned int k,
const input_t&)
54 return -((return_t)k) / ((return_t)(k + 1));
59 template <
typename return_t,
typename input_t>
60 static constexpr return_t
clenshaw_b(
const unsigned int k,
const unsigned int n,
const input_t& x)
67 return clenshaw_alpha<return_t>(k, x) * clenshaw_b<return_t>(k + 1, n, x) +
68 clenshaw_beta<return_t>(k + 1, x) * clenshaw_b<return_t>(k + 2, n, x);
84 template <
typename return_t,
typename input_t>
87 const bool normalized =
true)
90 const return_t x = 2. * (return_t)x_val - 1.;
93 return_t legendre_val;
104 x * clenshaw_b<return_t>(1,
index, x) - 0.5 * clenshaw_b<return_t>(2,
index, x);
126 template <
typename return_t,
typename input_t>
128 const input_t& x_val,
129 const bool normalized =
true)
135 const return_t x = 2. * (return_t)x_val - 1.;
140 return_t legendre_val = 0., x_pot = 1.;
141 for (
unsigned int p =
index; p > 0; --p)
143 legendre_val += ((return_t)p) * x_pot * fct_val<return_t>(p - 1, x_val,
false);
151 return 2. * legendre_val;
168 template <
unsigned int poly_deg>
174 static constexpr
unsigned int n_fun() {
return degree() + 1; }
178 static constexpr
unsigned int dim() {
return 1U; }
182 static constexpr
unsigned int degree() {
return poly_deg; }
197 template <
typename return_t,
typename input_t>
198 static constexpr return_t fct_val(
const unsigned int index,
199 const input_t& x_val,
200 const bool normalized =
true)
208 const return_t x = 2. * (return_t)x_val - 1.;
211 return_t lobatto_val =
212 (x * Legendre<poly_deg>::template fct_val<return_t>(
index - 1, x_val,
false) -
213 Legendre<poly_deg>::template fct_val<return_t>(
index - 2, x_val,
false)) /
220 return 0.5 * lobatto_val;
235 template <
typename return_t,
typename input_t>
236 static constexpr return_t der_val(
const unsigned int index,
237 const input_t& x_val,
238 const bool normalized =
true)
243 return Legendre<poly_deg>::template fct_val<return_t>(
index - 1, x_val, normalized);
static constexpr return_t clenshaw_beta(const unsigned int k, const input_t &)
Evaluate variable beta of the Clenshaw algorithm.
Definition: one_dimensional.hxx:52
index
Definition: check_push_test.py:10
Struct that handles the evaluation of one-dimensional Legendre polynomials.
Definition: one_dimensional.hxx:25
static constexpr return_t clenshaw_alpha(const unsigned int k, const input_t &x)
Evaluate variable alpha of the Clenshaw algorithm.
Definition: one_dimensional.hxx:44
static constexpr float_t heron_root(const float_t square)
Calculate the non-negative square-root of a non-negative number at compile time.
Definition: compile_time_tricks.hxx:17
This file provides the function tpp_assert.
Contains the functionalities to evaluate and integrate tensor product type polynomials.
Definition: compile_time_tricks.hxx:7
static constexpr unsigned int degree()
Maximum degree of all polynomials.
Definition: one_dimensional.hxx:38
static constexpr return_t fct_val(const unsigned int index, const input_t &x_val, const bool normalized=true)
Evaluate value of one shape function.
Definition: one_dimensional.hxx:85
static constexpr return_t der_val(const unsigned int index, const input_t &x_val, const bool normalized=true)
Evaluate value of the derivative of orthonormal shape function.
Definition: one_dimensional.hxx:127
static constexpr unsigned int n_fun()
Number of shape functions that span the local polynomial space.
Definition: one_dimensional.hxx:30
static constexpr return_t clenshaw_b(const unsigned int k, const unsigned int n, const input_t &x)
Evaluate variable b of the Clenshaw algorithm.
Definition: one_dimensional.hxx:60
static constexpr unsigned int dim()
Dimension of abscissas of the polynomials.
Definition: one_dimensional.hxx:34