E.V.E
v2023.02.15
 
Loading...
Searching...
No Matches

◆ legendre

eve::legendre = functor<legendre_t>
inlineconstexpr
  • The Legendre polynomial of order n is given by \(\displaystyle \mbox{L}_{n}(x) = \frac{e^x}{n!}\frac{d^n}{dx^n}(x^ne^{-x})\).
  • The associated legendre polynomial is given by \(\displaystyle \mbox{L}_{n}^{m} = (-1)^m\frac{d^m}{dx^m}\mbox{L}_{n+m}(x)\).

Defined in header

#include <eve/module/polynomial.hpp>

Callable Signatures

namespace eve
{
// Regular overload
constexpr auto legendre(integral_value auto n, floating_value auto x) noexcept; // 1
// Semantic options
constexpr auto legendre[p_kind](integral_value auto n, floating_value auto x) noexcept; // 1
constexpr auto legendre[q_kind](integral_value auto n, floating_value auto x) noexcept; // 2
constexpr auto legendre[associated](integral_value auto n,
floating_value auto x) noexcept; // 3
constexpr auto legendre[condon_shortley](integral_value auto n,
floating_value auto x) noexcept; // 4
constexpr auto legendre[sph](integral_value auto l, integral_value auto m,
floating_value auto theta) noexcept; // 5
constexpr auto legendre[successor](integral_value auto l,
integral_value auto ln, integral_value lnm1) noexcept; // 6
constexpr auto legendre[successor](integral_value auto l,integral_value auto m,
integral_value auto ln, integral_value lnm1) noexcept; // 6
}
The concept floating_value<T> is satisfied if and only if T satisfies eve::value and the element type...
Definition: value.hpp:95
The concept integral_value<T> is satisfied if and only if T satisfies eve::value and the element type...
Definition: value.hpp:46
constexpr auto legendre
Computes the value of the Legendre and associated Legendre polynomials of order n at x:
Definition: legendre.hpp:121
EVE Main Namespace.
Definition: abi.hpp:18

Parameters

  • n, m, ln, lnm1 : integral positive arguments.
  • x : floating argument.

Return value

  1. The value of the Legendre polynomial of order n at x is returned.
  2. The value of the Legendre polynomial of order n at x of the second kind is returned.
  3. The value of the associated legendre polynomial of orders n, m at x is returned.
  4. multiplies the associated legendre polynomial value by the Condon-Shortley phase \((-1)^m\) to match the definition given by Abramowitz and Stegun (8.6.6). This is currently the version implemented in boost::math.
  5. returns the spherical associated Legendre function of degree l, order m, and polar angle theta in radian (that is the classical spherical harmonic with \(\phi = 0\)), i.e. \(\displaystyle (-1)^m\frac{(2l+1)(l-m)!}{4\pi(l+m)!}\mbox{P}^m_{l}(\cos\theta)\)
  6. The successor option implements the three term recurrence relation for the (associated) Legendre polynomials, \(\displaystyle \mbox{P}^m_{l+1} = \left((2l+1)\mbox{P}^m_{l}(x)-l\mbox{P}^m_{l-1}(x)\right)/(l+m+1)\) ( \(m = 0\) and no \(m\) in call are equivalent here).

External references

Example

#include <eve/module/polynomial.hpp>
#include <eve/wide.hpp>
#include <iostream>
int main()
{
wide_ft xd = {-0.1, -0.2, -0.3, -0.5, 0.0, 0.2, 0.3, 2.0};
wide_it n = {0, 1, 2, 3, 4, 5, 6, 7};
wide_ft x(0.5);
std::cout << "---- simd" << '\n'
<< "<- xd = " << xd << '\n'
<< "<- n = " << n << '\n'
<< "<- x = " << x << '\n'
<< "-> legendre(n, xd) = " << eve::legendre(n, xd) << '\n'
<< "-> legendre(3, xd) = " << eve::legendre(3, xd) << '\n'
<< "-> legendre(n, 0.5) = " << eve::legendre(n, 0.5) << '\n'
<< "-> legendre(n, x) = " << eve::legendre(n, x) << '\n'
;
double xs = 0.1;
std::cout << "---- scalar" << '\n'
<< "<- xs = " << xs << '\n'
<< "-> legendre(4, xs) = " << eve::legendre(4, xs) << '\n';
return 0;
}
Wrapper for SIMD registers.
Definition: wide.hpp:71