Polynomials.
Constructor | Description |
Full Usage:
Polynomials()
|
|
Static member | Description |
Full Usage:
Polynomials.ChebyshevT(n, x, y, dy, d2y)
Parameters:
int
-
Degree of the polynomial >= 0.
x : float
-
Point in which the computation is performed, -1 <= x <= 1.
y : byref<float>
-
Output: value of the polynomial in x.
dy : byref<float>
-
Output: value of the first derivative in x.
d2y : byref<float>
-
Output: value of the second derivative in x.
|
Computes the value of the Chebyshev polynomial of degree n and its first and second derivatives at a given point.
Note: This C++ implementation is based on the Fortran function VACHPO from "Fortran routines for spectral methods" by Daniele Funaro Department of Mathematics University of Pavia Via Abbiategrasso 209, 27100 Pavia, Italy e-mails: fun18@ipvian.ian.pv.cnr.it funaro@dragon.ian.pv.cnr.it
|
Full Usage:
Polynomials.HermiteH(n, x, y, dy, d2y)
Parameters:
int
-
Degree of the polynomial >= 0.
x : float
-
Point in which the computation is performed.
y : byref<float>
-
Output: value of the polynomial in x.
dy : byref<float>
-
Output: value of the first derivative in x.
d2y : byref<float>
-
Output: value of the second derivative in x.
|
Computes the value of the Hermite polynomial of degree n and its first and second derivatives at a given point.
Note: This C++ implementation is based on the Fortran function VAHEPO from "Fortran routines for spectral methods" by Daniele Funaro Department of Mathematics University of Pavia Via Abbiategrasso 209, 27100 Pavia, Italy e-mails: fun18@ipvian.ian.pv.cnr.it funaro@dragon.ian.pv.cnr.it
|
Full Usage:
Polynomials.JacobiP(n, a, b, x, y, dy, d2y)
Parameters:
int
-
Degree of the polynomial >= 0.
a : float
-
Parameter > -1.
b : float
-
Parameter > -1.
x : float
-
Point in which the computation is performed, -1 <= x <= 1.
y : byref<float>
-
Output: value of the polynomial in x.
dy : byref<float>
-
Output: value of the first derivative in x.
d2y : byref<float>
-
Output: value of the second derivative in x.
|
Computes the value of the Jacobi polynomial of degree n and its first and second derivatives at a given point.
Note: This C++ implementation is based on the Fortran function VAJAPO from "Fortran routines for spectral methods" by Daniele Funaro Department of Mathematics University of Pavia Via Abbiategrasso 209, 27100 Pavia, Italy e-mails: fun18@ipvian.ian.pv.cnr.it funaro@dragon.ian.pv.cnr.it
|
Full Usage:
Polynomials.LaguerreL(n, a, x, y, dy, d2y)
Parameters:
int
-
Degree of the polynomial
a : float
-
Parameter > -1
x : float
-
Point in which the computation is performed, x >= 0
y : byref<float>
-
Output: value of the polynomial in x
dy : byref<float>
-
Output: value of the first derivative in x
d2y : byref<float>
-
Output: value of the second derivative in x
|
Computes the value of the Laguerre polynomial of degree n and its first and second derivatives at a given point.
Note: This C++ implementation is based on the Fortran function VALAPO from "Fortran routines for spectral methods" by Daniele Funaro Department of Mathematics University of Pavia Via Abbiategrasso 209, 27100 Pavia, Italy e-mails: fun18@ipvian.ian.pv.cnr.it funaro@dragon.ian.pv.cnr.it
|
Full Usage:
Polynomials.LegendreP(n, x, y, dy, d2y)
Parameters:
int
-
Degree of the polynomial >= 0.
x : float
-
Point in which the computation is performed, -1 <= x <= 1.
y : byref<float>
-
Output: value of the polynomial in x.
dy : byref<float>
-
Output: value of the first derivative in x.
d2y : byref<float>
-
Output: value of the second derivative in x.
|
Computes the value of the Legendre polynomial of degree n and its first and second derivatives at a given point.
Note: This C++ implementation is based on the Fortran function VALEPO from "Fortran routines for spectral methods" by Daniele Funaro Department of Mathematics University of Pavia Via Abbiategrasso 209, 27100 Pavia, Italy e-mails: fun18@ipvian.ian.pv.cnr.it funaro@dragon.ian.pv.cnr.it
|
Full Usage:
Polynomials.LegendreP(l, m, x)
Parameters:
int
m : int
x : float
Returns: float
|
|
Full Usage:
Polynomials.SphericalHarmonicY(l, m, theta, phi)
Parameters:
int
-
First integer.
m : int
-
Second integer, must be in the range -l <= m <= l.
theta : float
-
First angle.
phi : float
-
Second angle.
Returns: Complex
The spherical harmonics Y_lm(theta,phi).
|
Computes the spherical harmonics Y_lm(theta,phi) with l and m integers satisfying -l <= m <= l and arbitrary angles theta and phi.
|