Package com.polytechnik.utils
Class Legendre
java.lang.Object
com.polytechnik.utils.Legendre
- All Implemented Interfaces:
BasisFunctionsCalculatable
,BasisFunctionsMultipliable
,BasisPolynomials
,ConfederateMatrixCalculatable
,SimpleBasisPolynomials
Legendre polynomials in [-1:1] interval.
-
Field Summary
Fields -
Constructor Summary
Constructors -
Method Summary
Modifier and TypeMethodDescriptiondouble
calculate
(int n, double x) Calculation of the value of n-th Legendre polynomial at givenx
.double[]
getaxbP
(double[] p, double ma, double mb) Multiply p by ma*x+mb.double[][]
getBasisFunctionsOnPolynomialArgument
(int numQl, double[] pol, BasisFunctionsCalculatable PBASIS) This method calculates \( Q_l(pol(x))\), where pol(x) is a given polynomial and \( Q_l \) is basis finction of this (notPBASIS
!) basis.double[]
getConfederateMatrix
(double[] coefs) Calculates confederate matrix to solve polynomial equations in Legendre basis directly.double
getD0
(double[] Pcoefs, double x) Calculates \( \sum_k P_k(x) Pcoefs[k] \) using Clenshaw recurrence formula.double
getD1
(double[] Pcoefs, double x) Calculates \( d/dx \sum_k P_k(x) Pcoefs[k] \) using Clenshaw recurrence formula.double
getD2
(double[] Pcoefs, double x) Calculates \( d^2/dx^2 \sum_k P_k(x) Pcoefs[k] \) using Clenshaw recurrence formula.double[]
getDifferentiated
(double[] Pcoefs) Differentiate in Legendre basis.double[]
getIntegrated
(double[] Pcoefs) Integrate Legendre polynomial.double[]
getMomentsFromSample
(int nmoments, double[] x, double[] w, int nelsinsample) Given observations x and weights w calculate the 0..nmoms-1 moments in the basis.double[]
getNext
(int n, double[] Pnm1, double[] Pnm2, boolean flag_type2, BasisFunctionsCalculatable PBASIS) Calculates n-th Legendre polynomial.double[]
getXQkMomentsFromQkMoments
(double[] moments) Calculate <F|xPk> moments for k=0..n-1 from <F|Pk>, k=0..nstatic void
A unit test.sdiv1
(double[] p, double x0) Divide p by (x-x0), returns quotient and remainder.void
setNewtonBinomialLikeCoefs
(int n, double[] b, double[] q, double a, double x0) Calculates P_n(a(x+x0)) using 3 term recurrence Calculations in Legendre basis.void
setQlQmExpansionCoefs
(int l, int m, double[] s) Expand \( P_l(z)P_m(z)=\sum_{k=0}^{l+m} s_k P_k(z) \), formulae 8.915.5, A(9036), page 1040, Gradshtein & Ryzhik 1963.Methods inherited from class java.lang.Object
clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
Methods inherited from interface com.polytechnik.utils.BasisFunctionsCalculatable
convertFunctionTo_Basis_from_Monomials, convertFunctionTo_Basis_from_PBASIS, convertFunctionTo_Monomials_from_Basis, getMomentsFromSample, getMomentsFromSample, getONE, getX0Moments
Methods inherited from interface com.polytechnik.utils.BasisFunctionsMultipliable
getBasisFunctionsMultipliedByPolynomial, getKK, getMomentsOfMeasureProducingPolynomialInKK, getMomentsOfMeasureProducingPolynomialInKK_MQQM, getPQkMomentsFromQkMoments, getQQMatr, getTwoQuadraticFormsProductAsQuadraticForm, mult2Pol, sdiv
Methods inherited from interface com.polytechnik.utils.BasisPolynomials
getNext
Methods inherited from interface com.polytechnik.utils.SimpleBasisPolynomials
convertBasisToPBASIS, convertFunctionTo_PBASIS_from_Basis, getBasisFunctionsMultipliedByPolynomial, getBasisFunctionsOnPolynomialArgument
-
Field Details
-
BASIS
-
-
Constructor Details
-
Legendre
public Legendre()
-
-
Method Details
-
getNext
public double[] getNext(int n, double[] Pnm1, double[] Pnm2, boolean flag_type2, BasisFunctionsCalculatable PBASIS) Calculates n-th Legendre polynomial. Uses: \( nP_n=x*(2n-1)*P_{n-1}-(n-1)*P_{n-2} \).- Specified by:
getNext
in interfaceBasisPolynomials
- Parameters:
n
- P index.Pnm1
- Input \( P_{n-1} \).Pnm2
- Input \( P_{n-2} \).flag_type2
- Whether to use type 2 polynomials.PBASIS
- The basis in which to calculte.- Returns:
- \( P_{n} \) polynomial, array dimension n+1.
-
getBasisFunctionsOnPolynomialArgument
public double[][] getBasisFunctionsOnPolynomialArgument(int numQl, double[] pol, BasisFunctionsCalculatable PBASIS) Description copied from interface:SimpleBasisPolynomials
This method calculates \( Q_l(pol(x))\), where pol(x) is a given polynomial and \( Q_l \) is basis finction of this (notPBASIS
!) basis. This method may not be stable for a large dimensions ofpol.length
. Equivalent to numQl applications ofSimpleBasisPolynomials.getNext(int, double[], double[], com.polytechnik.utils.BasisFunctionsCalculatable)
topol
array. The method can be used to calculate e.g. \( T_n(P_k(x)) \) or \( T_n(P_k(L_m(x^j))) \). For example
meansnew Legendre().getBasisFunctionsOnPolynomialArgument(4,new double []{0,0,0,0,1},new Chebyshev()) ={ { 1.0 }, { 0.0, 0.0, 0.0, 0.0, 1.0 }, { 0.25, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.75 }, { 0.0, 0.0, 0.0, 0.0, 0.375, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.0, 0.625 } }
\( P_0(T_4(x))=T_0 \)
\( P_1(T_4(x))=T_4(x) \)
\( P_2(T_4(x))=0.25*T_0+0.75*T_8(x) \)
\( P_3(T_4(x))=0.375*T_4(x)+0.625*T_{12}(x) \).- Specified by:
getBasisFunctionsOnPolynomialArgument
in interfaceSimpleBasisPolynomials
- Parameters:
numQl
- The number of basis functions to calculate.pol
- The polynomial to use in \( Q_l(pol(x))\), it is given in thePBASIS
basis.PBASIS
- the basis in whichpol
is given and the result is returned.- Returns:
- An array of numQl dimension, elements: \( Q_0(pol(x)) \),
\( Q_1(pol(x)) \), ... \( Q_{numQl-1}(pol(x)) \) in the basis
PBASIS
.
-
calculate
public double calculate(int n, double x) Calculation of the value of n-th Legendre polynomial at givenx
. Uses: \( nP_n=x *(2n-1)P_{n-1}-(n-1)P_{n-2} \).- Specified by:
calculate
in interfaceBasisFunctionsCalculatable
- Parameters:
n
- P index.x
- Argument.- Returns:
- The value of the Legendre polynomial \( P_n \) at
x
.
-
getMomentsFromSample
public double[] getMomentsFromSample(int nmoments, double[] x, double[] w, int nelsinsample) Description copied from interface:BasisFunctionsCalculatable
Given observations x and weights w calculate the 0..nmoms-1 moments in the basis.- Specified by:
getMomentsFromSample
in interfaceBasisFunctionsCalculatable
- Parameters:
nmoments
- The number of moments to calculate.x
- observation samples.w
- weights.nelsinsample
- How many elements to sample.- Returns:
- nelsinsample generalzied moments, \( k=0\dots nmoms-1 \): \( \sum_{j=0}^{nelsinsample-1} Q_k(x[j])*w[j] \).
-
getD0
public double getD0(double[] Pcoefs, double x) Calculates \( \sum_k P_k(x) Pcoefs[k] \) using Clenshaw recurrence formula.- Specified by:
getD0
in interfaceBasisFunctionsCalculatable
- Parameters:
Pcoefs
- Coefficients.x
- Pont to estimate the sum.- Returns:
- The sum.
-
getD1
public double getD1(double[] Pcoefs, double x) Calculates \( d/dx \sum_k P_k(x) Pcoefs[k] \) using Clenshaw recurrence formula.- Specified by:
getD1
in interfaceBasisFunctionsCalculatable
- Parameters:
Pcoefs
- Coefficients.x
- Pont to estimate the sum.- Returns:
- The derivative.
-
getD2
public double getD2(double[] Pcoefs, double x) Calculates \( d^2/dx^2 \sum_k P_k(x) Pcoefs[k] \) using Clenshaw recurrence formula.- Specified by:
getD2
in interfaceBasisFunctionsCalculatable
- Parameters:
Pcoefs
- Coefficients.x
- Pont to estimate the sum.- Returns:
- The second derivative.
-
getDifferentiated
public double[] getDifferentiated(double[] Pcoefs) Differentiate in Legendre basis. formulae 8.915.2, A(9036), page 1040, Gradshtein & Ryzhik 1963.- Specified by:
getDifferentiated
in interfaceBasisFunctionsCalculatable
- Parameters:
Pcoefs
- : \( F(x)=\sum_k Q_k(x) Qcoefs[k] \).- Returns:
- \( dF/dx=\sum_k Q'_k(x) Qcoefs[k] \).
-
getIntegrated
public double[] getIntegrated(double[] Pcoefs) Integrate Legendre polynomial. Prudnikov, Brychkov, Marichev, volume 2, 2003, formulae 1.14.1.1 page 46.- Specified by:
getIntegrated
in interfaceBasisFunctionsCalculatable
- Parameters:
Pcoefs
- : \( F(x)=\sum_k Q_k(x) Qcoefs[k] \).- Returns:
- \( \int F(x)dx=\sum_k Qcoefs[k] \int Q_k(x) dx \). The integral is typically normalized as I(0)=0.
-
setQlQmExpansionCoefs
public void setQlQmExpansionCoefs(int l, int m, double[] s) Expand \( P_l(z)P_m(z)=\sum_{k=0}^{l+m} s_k P_k(z) \), formulae 8.915.5, A(9036), page 1040, Gradshtein & Ryzhik 1963. The formulae of Neumann-Adams, see W. N. Bailey, Generalized hypergeometric series, p.101, Cambridge, 1962. Second option is to use Gaunt's formula: see John C. Slater Quantum Theory of Atomic Structure, McGraw-Hill (New York, 1960), Volume I, page 309, which cites the original work of J. A. Gaunt, Philosophical Transactions of the Royal Society of London, A228:151 (1929).- Specified by:
setQlQmExpansionCoefs
in interfaceBasisFunctionsMultipliable
- Parameters:
l
- l.m
- m.s
- coefs.
-
getaxbP
public double[] getaxbP(double[] p, double ma, double mb) Multiply p by ma*x+mb. Calculations in Legendre basis. Use Legendre recurrence, see formulae 8.914.1, UVII99, page 1040, Gradshtein & Ryzhik 1963. Do not mistake with a and b recurrence coefficiemts ofgetAB()
.- Specified by:
getaxbP
in interfaceBasisFunctionsCalculatable
- Parameters:
p
- Input polynomial.- Returns:
- p*(ma*x+mb), an array of the length 1+p.length.
-
sdiv1
Divide p by (x-x0), returns quotient and remainder. Calculations in Legendre basis. Use Legendre recurrence, see formulae 8.914.1, UVII99, page 1040, Gradshtein & Ryzhik 1963.- Specified by:
sdiv1
in interfaceBasisFunctionsCalculatable
- Parameters:
p
- Input polynomial.x0
- To divide on (x-x0)- Returns:
- Quotient and remainder.
-
setNewtonBinomialLikeCoefs
public void setNewtonBinomialLikeCoefs(int n, double[] b, double[] q, double a, double x0) Calculates P_n(a(x+x0)) using 3 term recurrence Calculations in Legendre basis. The b_k are defined as: \( P_n(a(x+x0))=\sum_{k=0}^{k=n}b_k^{n} P_k(x) \) The method should be applied consequently to get \( b_k^{n} \) for n=0,1,2,3,... Uses Legendre recurrence, see formulae 8.914.1, UVII99, page 1040, Gradshtein & Ryzhik 1963.- Specified by:
setNewtonBinomialLikeCoefs
in interfaceBasisPolynomials
- Parameters:
n
- The max dimension of b.b
- On input: \( b_k^{n-1} \), on output \( b_k^{n} \).q
- On input: \( b_k^{n-2} \), on output: \( b_k^{n-1} \).
-
getXQkMomentsFromQkMoments
public double[] getXQkMomentsFromQkMoments(double[] moments) Calculate <F|xPk> moments for k=0..n-1 from <F|Pk>, k=0..n- Specified by:
getXQkMomentsFromQkMoments
in interfaceBasisFunctionsCalculatable
- Parameters:
moments
- <F|Pk>.- Returns:
- <F|xPk>
-
getConfederateMatrix
public double[] getConfederateMatrix(double[] coefs) Calculates confederate matrix to solve polynomial equations in Legendre basis directly. We use the three term recurrence: \[ P_k=(x*(2*k-1)/k)*P_{k-1}-((k-1)/k)*P_{k-2} \] UseABASISLEGENDRE
to calculate roots.- Specified by:
getConfederateMatrix
in interfaceConfederateMatrixCalculatable
-
main
A unit test.
-