Class Laguerre

java.lang.Object
com.polytechnik.utils.Laguerre
All Implemented Interfaces:
BasisFunctionsCalculatable, BasisFunctionsMultipliable, BasisPolynomials, ConfederateMatrixCalculatable, SimpleBasisPolynomials

public class Laguerre extends Object implements BasisPolynomials
Laguerre polynomials.
  • Field Details

    • BASIS

      private static final Laguerre BASIS
  • Constructor Details

    • Laguerre

      public Laguerre()
  • Method Details

    • getNext

      public double[] getNext(int n, double[] Lnm1, double[] Lnm2, boolean flag_type2, BasisFunctionsCalculatable PBASIS)
      Calculates n-th Laguerre polynomial using \( L_n=(x*(-1/n)-(1-2*n)/n)*L_{n-1}-((n-1)/n)*L_{n-2} \).
      Specified by:
      getNext in interface BasisPolynomials
      Parameters:
      n - L index.
      Lnm1 - Input \( L_{n-1} \).
      Lnm2 - Input \( L_{n-2} \).
      flag_type2 - Whether to use type 2 polynomials.
      PBASIS - The basis in which to calculte.
      Returns:
      L_{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 (not PBASIS!) basis. This method may not be stable for a large dimensions of pol.length. Equivalent to numQl applications of SimpleBasisPolynomials.getNext(int, double[], double[], com.polytechnik.utils.BasisFunctionsCalculatable) to pol 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
      
      new 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 }
      }
      
      means
      \( 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 interface SimpleBasisPolynomials
      Parameters:
      numQl - The number of basis functions to calculate.
      pol - The polynomial to use in \( Q_l(pol(x))\), it is given in the PBASIS basis.
      PBASIS - the basis in which pol 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 Laguerre polynomial at given x.
      Specified by:
      calculate in interface BasisFunctionsCalculatable
      Parameters:
      n - L index.
      x - Argument.
      Returns:
      The value of the Laguerre polynomial Ln 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 interface BasisFunctionsCalculatable
      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[] Lcoefs, double x)
      Calculates \( \sum_k L_k(x) Lcoefs[k] \) using Clenshaw recurrence formula.
      Specified by:
      getD0 in interface BasisFunctionsCalculatable
      Parameters:
      Lcoefs - Coefficients.
      x - Point to estimate the sum.
      Returns:
      The sum.
    • getD1

      public double getD1(double[] Lcoefs, double x)
      Calculates \( d/dx \sum_k L_k(x) Lcoefs[k] \) using Clenshaw recurrence formula.
      Specified by:
      getD1 in interface BasisFunctionsCalculatable
      Parameters:
      Lcoefs - Coefficients.
      x - Pont to estimate the sum.
      Returns:
      The derivative.
    • getD2

      public double getD2(double[] Lcoefs, double x)
      Calculates \( d^2/dx^2 \sum_k L_k(x) Lcoefs[k] \) using Clenshaw recurrence formula.
      Specified by:
      getD2 in interface BasisFunctionsCalculatable
      Parameters:
      Lcoefs - Coefficients.
      x - Pont to estimate the sum.
      Returns:
      The derivative.
    • getDifferentiated

      public double[] getDifferentiated(double[] Lcoefs)
      Differentiate a function in Laguerre basis, formulae 8.971.2, page 1051, and 8.974.3, page 1052, Gradshtein & Ryzhik 1963.
      Specified by:
      getDifferentiated in interface BasisFunctionsCalculatable
      Parameters:
      Lcoefs - input, \( F(x)=\sum_k L_k(x)*Lcoefs[k] \)
      Returns:
      \( \sum_k L_k(x)'*Lcoefs[k] \).
    • getIntegrated

      public double[] getIntegrated(double[] Lcoefs)
      Integrate a function in Laguerre basis, formulae Prudnikov, Brychkov, Marichev, volume 2, 2003, formulae 1.14.3.1 page 47.
      Specified by:
      getIntegrated in interface BasisFunctionsCalculatable
      Parameters:
      Lcoefs - Input: \( F(x)=\sum_k L_k(x)*Lcoefs[k] \).
      Returns:
      \( \int_0^x F(x)dx=\sum_k \int_0^x L_k(x)dx*Lcoefs[k] \).
    • setQlQmExpansionCoefs

      public void setQlQmExpansionCoefs(int l, int m, double[] c)
      Expand \( L_l(z)L_m(z)=\sum_{k=0}^{l+m} c_k L_k(z) \), G.N. Watson London Math. Soc. Jn., v 13, 1938, p.29. See also: Joseph Gillis and George Weiss, Products of Laguerre Polynomials.
      Specified by:
      setQlQmExpansionCoefs in interface BasisFunctionsMultipliable
      Parameters:
      l - l.
      m - m.
      c - coefs.
    • getaxbP

      public double[] getaxbP(double[] p, double ma, double mb)
      Multiply p by ma*x+mb. Calculations in Laguerre. Use Laguerre recurrence, see formulae 8.971.6, MO109, page 1051, Gradshtein & Ryzhik 1963. Do not mistake with a and b recurrence coefficiemts of getAB().
      Specified by:
      getaxbP in interface BasisFunctionsCalculatable
      Parameters:
      p - Input polynomial.
      Returns:
      p*(ma*x+mb), an array of the length 1+p.length.
    • sdiv1

      public SDivRes sdiv1(double[] p, double x0)
      Divide p by (x-x0), returns quotient and remainder. Calculations in Laguerre basis. Use Laguerre recurrence, see formulae 8.971.6, MO109, page 1051, Gradshtein & Ryzhik 1963.
      Specified by:
      sdiv1 in interface BasisFunctionsCalculatable
      Parameters:
      p - Input polynomial.
      x0 - To divide on (x-x0)
      Returns:
      Quotient and remainer.
    • setNewtonBinomialLikeCoefs

      public void setNewtonBinomialLikeCoefs(int n, double[] b, double[] q, double a, double x0)
      Calculates \( L_n(a(x+x0)) \) using 3 term recurrence. Calculations in Laguerre basis. The \( b_k^{n} \) are defined as: \( L_n(a(x+x0))=\sum_{k=0}^{k=n}b_k^{n} L_k(x) \) The method should be applied consequently to get \( b_k^{n} \) for n=0,1,2,3,... Uses Laguerre recurrence, see formulae 8.971.6, MO109, page 1051, Gradshtein & Ryzhik 1963.
      Specified by:
      setNewtonBinomialLikeCoefs in interface BasisPolynomials
      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|xLk> moments for k=0..n-1 from <F|Lk>, k=0..n
      Specified by:
      getXQkMomentsFromQkMoments in interface BasisFunctionsCalculatable
      Parameters:
      moments - <F|Lk>.
      Returns:
      <F|xLk>.
    • getConfederateMatrix

      public double[] getConfederateMatrix(double[] coefs)
      Calculates confederate matrix to solve polynomial equations in Laguerre basis directly. We use the three term recurrence: \[ L_k=(x*(-1/k)-(1-2*k)/k)*L_{k-1}-((k-1)/k)*L_{k-2} \] Use ABASISLAGUERRE to calculate roots.
      Specified by:
      getConfederateMatrix in interface ConfederateMatrixCalculatable
    • getxP

      public static double[] getxP(double[] p)
      Multiply p by x. Calculations in Laguerre. Use Laguerre recurrence, see formulae 8.971.6, MO109, page 1051, Gradshtein & Ryzhik 1963.
      Parameters:
      p - Input polynomial.
      Returns:
      x*p
    • get2LaguerreAveraged

      public static double get2LaguerreAveraged(double[] p1, double[] p2)
    • getL1FromL0

      public static double[] getL1FromL0(double[] l0)
      Calculate <F|L^1_n> moments.
      Parameters:
      l0 - Input, l0[k]=<F|L^0_k>.
      Returns:
      <F|L^1_k>
    • getLaguerreExpFIntegrated

      public static double[] getLaguerreExpFIntegrated(double[] Lcoefs)
      \( \int_x^{\infty} f(x) \exp(-x) dx = i(x) \exp(-x) \). Prudnikov, Brychkov, Marichev, volume 2, 2003, formulae 1.14.3.6 page 47.
      Returns:
      \( i(x) \).
    • getRoots

      static double[] getRoots(int n, double x0, boolean flag_radau)
      A different way to calculate roots. Use Jacobi matrix.
      Parameters:
      n - The Laguerre polynomial order.
      x0 - Radau edge, ignored if flag_radau=false;
      flag_radau - Whether to use Radau polynomial.
      Returns:
      The roots.
    • getRoots

      public static double[] getRoots(int n)
    • getRadauRoots

      public static double[] getRadauRoots(int n)
    • getWeights

      public static double[] getWeights(double[] xx)
      Uses Abramowitz and Stegun 1972, p. 890.
    • getNodesPolynomial

      public static double[] getNodesPolynomial(int n)
    • getRadauNodesPolynomial

      public static double[] getRadauNodesPolynomial(int n)
    • getRadauWeights

      public static double[] getRadauWeights(double[] xr)
      Radau exact weights from Walter Gautschi, "Gauss-Radau formulae for Jacobi and Laguerre weight functions". Mathematics and Computers in Simulation 54 (2000) 403-412. https://www.cs.purdue.edu/homes/wxg/selected_works/section_08/164.pdf
    • getChristoffelKernel

      public static double[] getChristoffelKernel(int n)
      Christoffel function.
    • setQlQmExpansionCoefs_test

      public static void setQlQmExpansionCoefs_test(int l, int m, double[] s)
      The same as setQlQmExpansionCoefs, but uses 3 term recursion instead of hypergeometric function expansion. It is rather slow at large l,m. Used for unit test.
    • sumQlQmExpansionCoefs

      private static void sumQlQmExpansionCoefs(int l, int m, double[] s, double c)
    • main

      public static void main(String[] args)
      A unit test.