ARTS  2.4.0(git:4fb77825)
legendre.cc File Reference

Contains the code to calculate Legendre polynomials. More...

#include <cmath>
#include <limits>
#include <sstream>
#include "exceptions.h"
#include "legendre.h"
#include "math_funcs.h"

Go to the source code of this file.

Functions

Numeric legendre_poly (Index l, Index m, Numeric x)
 legendre_poly More...
 
Numeric legendre_poly_norm_schmidt (Index l, Index m, Numeric x)
 legendre_poly_norm_schmidt More...
 
Numeric legendre_poly_deriv (Index l, Index m, Numeric x)
 legendre_poly_deriv More...
 
Numeric legendre_poly_norm_schmidt_deriv (Index l, Index m, Numeric x)
 legendre_poly_norm_schmidt_deriv More...
 
Numeric g_legendre_poly (Index l, Index m, Numeric x)
 g_legendre_poly More...
 
Numeric g_legendre_poly_norm_schmidt (Index l, Index m, Numeric x)
 g_legendre_poly_norm_schmidt More...
 
Numeric g_legendre_poly_deriv (Index l, Index m, Numeric x)
 g_legendre_poly_deriv More...
 
Numeric g_legendre_poly_norm_schmidt_deriv (Index l, Index m, Numeric x)
 g_legendre_poly_norm_schmidt_deriv More...
 
Numeric g_legendre_poly_norm_schmidt_deriv1 (Index l, Index m, Numeric x)
 g_legendre_poly_norm_schmidt_deriv1 More...
 
Numeric g_legendre_poly_norm_schmidt_deriv2 (Index l, Index m, Numeric x)
 g_legendre_poly_norm_schmidt_deriv2 More...
 
Numeric g_legendre_poly_norm_schmidt_deriv3 (Index l, Index m, Numeric x)
 g_legendre_poly_norm_schmidt_deriv3 More...
 
Numeric g_legendre_poly_norm_schmidt_deriv4 (Index l, Index m, Numeric x)
 g_legendre_poly_norm_schmidt_deriv4 More...
 
bool gsl_integration_glfixed_table_alloc (Vector &x, Vector &w, Index n)
 gsl_integration_glfixed_table_alloc More...
 

Detailed Description

Contains the code to calculate Legendre polynomials.

Author
Oliver Lemke
Date
2003-08-14

Definition in file legendre.cc.

Function Documentation

◆ g_legendre_poly()

Numeric g_legendre_poly ( Index  l,
Index  m,
Numeric  x 
)

g_legendre_poly

Returns the associated Legendre polynomial Plm(x) without the factor (-1)^m.

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| <= 1

The code is based on the Numerical recipes. Results were compared to the Legendre calculations from the GNU Scientific library and found to be identical.

Returns
Plm
Parameters
lIndex
mIndex
xValue
Author
Nikolay Koulev
Date
2003-09-02

Definition at line 290 of file legendre.cc.

References abs, ll, sqrt(), and ARTS::Var::x().

Referenced by g_legendre_poly_deriv(), g_legendre_poly_norm_schmidt(), g_legendre_poly_norm_schmidt_deriv(), g_legendre_poly_norm_schmidt_deriv1(), g_legendre_poly_norm_schmidt_deriv2(), g_legendre_poly_norm_schmidt_deriv3(), and g_legendre_poly_norm_schmidt_deriv4().

◆ g_legendre_poly_deriv()

Numeric g_legendre_poly_deriv ( Index  l,
Index  m,
Numeric  x 
)

g_legendre_poly_deriv

Returns the derivative of the associated Legendre polynomial Plm(x)) without the factor (-1)^m..

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| < 1

Returns
dPlm
Parameters
lIndex
mIndex
xValue
Author
Nikolay Koulev
Date
2003-09-02

Definition at line 385 of file legendre.cc.

References g_legendre_poly(), sqrt(), and ARTS::Var::x().

◆ g_legendre_poly_norm_schmidt()

Numeric g_legendre_poly_norm_schmidt ( Index  l,
Index  m,
Numeric  x 
)

g_legendre_poly_norm_schmidt

Returns the Schmidt quasi-normalized associated Legendre polynomial Plm(x)) without the factor (-1)^m..

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| <= 1

The code is based on the Numerical recipes. Results were compared to the Legendre calculations from the GNU Scientific library and found to be identical.

Returns
Plm
Parameters
lIndex
mIndex
xValue
Author
Nikolay koulev
Date
2003-09-02

Definition at line 355 of file legendre.cc.

References fac(), g_legendre_poly(), sqrt(), and ARTS::Var::x().

◆ g_legendre_poly_norm_schmidt_deriv()

Numeric g_legendre_poly_norm_schmidt_deriv ( Index  l,
Index  m,
Numeric  x 
)

g_legendre_poly_norm_schmidt_deriv

Returns the derivative of the Schmidt quasi-normalized associated Legendre polynomial Plm(x)) without the factor (-1)^m.

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| < 1

Returns
dPlm
Parameters
lIndex
mIndex
xValue
Author
Nikolay Koulev
Date
2003-09-02

Definition at line 451 of file legendre.cc.

References fac(), g_legendre_poly(), sqrt(), and ARTS::Var::x().

◆ g_legendre_poly_norm_schmidt_deriv1()

Numeric g_legendre_poly_norm_schmidt_deriv1 ( Index  l,
Index  m,
Numeric  x 
)

g_legendre_poly_norm_schmidt_deriv1

Returns the derivative of the Schmidt quasi-normalized associated Legendre polynomial Plm(x)) without the factor (-1)^m. Utilizes the simplest recurrence scheme.

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| < 1

Returns
dPlm
Parameters
lIndex
mIndex
xValue
Author
Nikolay Koulev
Date
2003-09-02

Definition at line 522 of file legendre.cc.

References fac(), g_legendre_poly(), sqrt(), and ARTS::Var::x().

◆ g_legendre_poly_norm_schmidt_deriv2()

Numeric g_legendre_poly_norm_schmidt_deriv2 ( Index  l,
Index  m,
Numeric  x 
)

g_legendre_poly_norm_schmidt_deriv2

Returns the derivative of the Schmidt quasi-normalized associated Legendre polynomial Plm(x)) without the factor (-1)^m.

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| < 1

Returns
dPlm
Parameters
lIndex
mIndex
xValue
Author
Nikolay Koulev
Date
2003-09-02

Definition at line 594 of file legendre.cc.

References fac(), g_legendre_poly(), sqrt(), and ARTS::Var::x().

◆ g_legendre_poly_norm_schmidt_deriv3()

Numeric g_legendre_poly_norm_schmidt_deriv3 ( Index  l,
Index  m,
Numeric  x 
)

g_legendre_poly_norm_schmidt_deriv3

Returns the derivative of the Schmidt quasi-normalized associated Legendre polynomial Plm(x)) without the factor (-1)^m.

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| < 1

Returns
dPlm
Parameters
lIndex
mIndex
xValue
Author
Nikolay Koulev
Date
2003-09-02

Definition at line 663 of file legendre.cc.

References fac(), g_legendre_poly(), sqrt(), and ARTS::Var::x().

◆ g_legendre_poly_norm_schmidt_deriv4()

Numeric g_legendre_poly_norm_schmidt_deriv4 ( Index  l,
Index  m,
Numeric  x 
)

g_legendre_poly_norm_schmidt_deriv4

Returns the derivative of the Schmidt quasi-normalized associated Legendre polynomial Plm(x)) without the factor (-1)^m.

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| < 1

Returns
dPlm
Parameters
lIndex
mIndex
xValue
Author
Nikolay Koulev
Date
2003-09-02

Definition at line 732 of file legendre.cc.

References fac(), g_legendre_poly(), sqrt(), and ARTS::Var::x().

◆ gsl_integration_glfixed_table_alloc()

bool gsl_integration_glfixed_table_alloc ( Vector x,
Vector w,
Index  n 
)

gsl_integration_glfixed_table_alloc

This function determines the Gauss-Legendre abscissae and weights necessary for an n-point fixed order integration scheme. If possible, high precision precomputed coefficients are used. If precomputed weights are not available, lower precision coefficients are computed on the fly.

Code is taken from the GNU Scientific Library under GPLv3 license. Adapted to return ARTS datatypes.

Parameters
[out]xVector with Gauss-Legendre abscissae.
[out]wVector with Gauss-Legendre weights.
[in]nIntegration scheme order.
Returns
true if precomputed values have been used.
Author
Oliver Lemke
Date
2018-05-29

Definition at line 2685 of file legendre.cc.

References w(), and ARTS::Var::x().

Referenced by AngularGridsSetFluxCalc(), and test_gsl_int().

◆ legendre_poly()

Numeric legendre_poly ( Index  l,
Index  m,
Numeric  x 
)

legendre_poly

Returns the associated Legendre polynomial Plm(x).

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| <= 1

The code is based on the Numerical recipes. Results were compared to the Legendre calculations from the GNU Scientific library and found to be identical.

Returns
Plm
Parameters
lIndex
mIndex
xValue
Author
Oliver Lemke
Date
2003-08-14

Definition at line 58 of file legendre.cc.

References abs, ll, sqrt(), and ARTS::Var::x().

Referenced by legendre_poly_deriv(), legendre_poly_norm_schmidt(), and main().

◆ legendre_poly_deriv()

Numeric legendre_poly_deriv ( Index  l,
Index  m,
Numeric  x 
)

legendre_poly_deriv

Returns the derivative of the associated Legendre polynomial Plm(x).

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| < 1

Returns
dPlm
Parameters
lIndex
mIndex
xValue
Author
Oliver Lemke
Date
2003-08-18

Definition at line 151 of file legendre.cc.

References legendre_poly(), sqrt(), and ARTS::Var::x().

Referenced by main().

◆ legendre_poly_norm_schmidt()

Numeric legendre_poly_norm_schmidt ( Index  l,
Index  m,
Numeric  x 
)

legendre_poly_norm_schmidt

Returns the Schmidt quasi-normalized associated Legendre polynomial Plm(x).

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| <= 1

The code is based on the Numerical recipes. Results were compared to the Legendre calculations from the GNU Scientific library and found to be identical.

Returns
Plm
Parameters
lIndex
mIndex
xValue
Author
Oliver Lemke
Date
2003-08-15

Definition at line 122 of file legendre.cc.

References fac(), legendre_poly(), sqrt(), and ARTS::Var::x().

Referenced by legendre_poly_norm_schmidt_deriv(), and main().

◆ legendre_poly_norm_schmidt_deriv()

Numeric legendre_poly_norm_schmidt_deriv ( Index  l,
Index  m,
Numeric  x 
)

legendre_poly_norm_schmidt_deriv

Returns the derivative of the Schmidt quasi-normalized associated Legendre polynomial Plm(x).

The input parameters must fulfill the following conditions: 0 <= m <= l and |x| < 1

Returns
dPlm
Parameters
lIndex
mIndex
xValue
Author
Nikolay Koulev
Date
2003-08-18

Definition at line 216 of file legendre.cc.

References fac(), legendre_poly_norm_schmidt(), sqrt(), and ARTS::Var::x().

Referenced by main().