# -*- coding: utf-8 -*-
"""
=========================================
DABI Spectral 1D
=========================================
Created on Thu Mar 29 11:33:32 2012
@author: Daniele Bigoni (dabi@imm.dtu.dk)
Implementation of Spectral Methods in 1 dimension.
Available polynomials:
* Jacobi
* Hermite Physicist
* Hermite Function
* Hermite Probabilistic
* Laguerre Polynomial
* Laguerre Function
* ORTHPOL package (generation of recursion coefficients using [4]_)
Some of the useful references used while developing this toolbox are listed in the following
.. [1] "Implemenenting Spectral Methods for Partial Differential Equations" by David A. Kopriva, Springer, 2009
.. [2] J. Shen and L.L. Wang, “Some recent advances on spectral methods for unbounded domains”. Communications in Computational Physics, vol. 5, no. 2–4, pp. 195–241, 2009
.. [3] W.H. Press, Numerical Recipes: "The Art of Scientific Computing". Cambridge University Press, 2007
.. [4] W. Gautschi, "Algorithm 726: ORTHPOL -- a package of routines for generating orthogonal polynomials and Gauss-type quadrature rules". ACM Trans. Math. Softw., vol. 20, issue 1, pp. 21-62, 1994
"""
__revision__ = filter(str.isdigit, "$Revision: 93 $")
__author__ = "Daniele Bigoni"
__copyright__ = """Copyright 2012, Daniele Bigoni"""
__credits__ = ["Daniele Bigoni"]
__maintainer__ = "Daniele Bigoni"
__email__ = "dabi@imm.dtu.dk"
__status__ = "Production"
import sys
import numpy as np
from numpy import linalg as LA
from numpy import fft as FFT
from scipy.special import gamma as gammaF
from scipy import factorial
import orthpol
JACOBI = 'Jacobi'
HERMITEP = 'HermiteP'
HERMITEF = 'HermiteF'
HERMITEP_PROB = 'HermitePprob'
LAGUERREP = 'LaguerreP'
LAGUERREF = 'LaguerreF'
ORTHPOL = 'ORTHPOL'
AVAIL_POLY = [JACOBI, HERMITEP, HERMITEF, HERMITEP_PROB, LAGUERREP, LAGUERREF, ORTHPOL]
[docs]class Poly1D(object):
poly = []
params = []
def __init__(self,poly,params):
"""
Initialization of the Polynomial instance.
Syntax:
``p = Poly1D(poly,params)``
Input:
* ``poly`` = The orthogonal polynomial type desired
* ``params`` = The parameters needed by the selected polynomial
Description:
This method generates an instance of the Poly1D class, to be used in order to generate
orthogonal basis of the polynomial type selected. Avaliable polynomial types can be
selected using their string name or by predefined attributes
* 'Jacobi' or ``DABISpectral1D.JACOBI``
* 'HermiteP' or ``DABISpectral1D.HERMITEP``
* 'HermiteF' or ``DABISpectral1D.HERMITEF``
* 'HermitePprob' or ``DABISpectral1D.HERMITEP_PROB``
* 'LaguerreP' or ``DABISpectral1D.LAGUERREP``
* 'LaguerreF' or ``DABISpectral1D.LAGUERREF``
* 'ORTHPOL' or ``DABISpectral1D.ORTHPOL``
Additional parameters are required for some polynomials.
+--------------+--------------+
| Polynomial | Parameters |
+==============+==============+
| Jacobi | (alpha,beta) |
+--------------+--------------+
| HermiteP | None |
+--------------+--------------+
| HermiteF | None |
+--------------+--------------+
| HermitePprob | None |
+--------------+--------------+
| LaguerreP | alpha |
+--------------+--------------+
| LaguerreF | alpha |
+--------------+--------------+
| ORTHPOL | see notes |
+--------------+--------------+
.. note:: The ORTHPOL polynomials are built up using the "Multiple-Component Discretization Procedure" described in [4]_. The following parameters describing the measure function are required in order to use the procedure for finding the recursion coefficients (alpha,beta) and have to be provided at construction time:
* ``ncapm``: (int) maximum integer N0 (default = 500)
* ``mc``: (int) number of component intervals in the continuous part of the spectrum
* ``mp``: (int) number of points in the discrete part of the spectrum. If the measure has no discrete part, set mp=0
* ``xp``, ``yp``: (Numpy 1d-array) of dimension mp, containing the abscissas and the jumps of the point spectrum
* ``mu``: (function) measure function that returns the mass (float) with arguments: ``x`` (float) absissa, ``i`` (int) interval number in the continuous part
* ``eps``: (float) desired relative accuracy of the nonzero recurion coefficents
* ``irout``: (int) selects the routine for generating the recursion coefficients from the discrete inner product; ``irout=1`` selects the routine ``sti``, ``irout!=1`` selects the routine ``lancz``
* ``finl``, ``finr``: (bool) specify whether the extreme left/right interval is finite (false for infinite)
* ``endl``, ``endr``: (Numpy 1d-array) of dimension ``mc`` containing the left and right endpoints of the component intervals. If the first of these extends to -infinity, endl[0] is not being used by the routine.
Parameters ``iq``, ``quad``, ``idelta`` in [4]_ are suppressed. Instead the routine ``qgp`` of ORTHPOL [4]_ is used by default (``iq=1`` and ``idelta=2``)
"""
# Check consistency of polynomial types and parameters
if (poly in AVAIL_POLY):
if (poly == JACOBI):
if len(params) != 2:
print "The number of parameters inserted for the polynomial of type '%s' is not correct" % poly
return
if ((poly == LAGUERREP) or (poly == LAGUERREF)):
if len(params) != 1:
print "The number of parameters inserted for the polynomial of type '%s' is not correct" % poly
return
if ((poly == ORTHPOL)):
if len(params) != 12:
print "The number of parameters inserted for the polynomial of type '%s' is not correct" % poly
return
else:
print "The inserted type of polynomial is not available."
return
self.poly = poly
self.params = params
def __AlmostEqual(self,a,b):
"""
b = AlmostEqual(a,b)
Test equality of two floating point numbers. Returns a boolean.
"""
eps = np.finfo(np.float64).eps
if ((a == 0) or (b == 0)):
if (abs(a-b) <= 2*eps):
return True
else:
return False
else:
if ( (abs(a-b) <= eps*abs(a)) and (abs(a-b) <= eps*abs(b)) ):
return True
else:
return False
def __JacobiGQ(self,N):
"""
Purpose: Compute the N'th order Gauss quadrature points, x,
and weights, w, associated with the Jacobi
polynomial, of type (alpha,beta) > -1 ( <> -0.5).
"""
if (self.poly != JACOBI):
print "The method __JacobiGQ cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,-0.5) or self.__AlmostEqual(beta,-0.5) ):
return self.__JacobiCGQ(N)
x = np.zeros((N+1))
w = np.zeros((N+1))
if (N == 0):
x[0] = -(alpha-beta)/(alpha+beta+2)
w[0] = 2
return (x,w)
J = np.zeros((N+1,N+1))
h1 = 2.* np.arange(0.,N+1)+alpha+beta
J = np.diag(-1/2*(alpha**2-beta**2)/(h1+2)/h1) + np.diag( 2/(h1[range(0,N)] + 2) * np.sqrt( np.arange(1.,N+1) * (np.arange(1.,N+1)+alpha+beta) * (np.arange(1.,N+1)+alpha) * (np.arange(1.,N+1)+beta) / (h1[range(0,N)] + 1) / (h1[range(0,N)] + 3) ) , 1)
if (alpha + beta < 10.*np.finfo(np.float64).eps): J[0,0] = 0.0
J = J + np.transpose(J)
# Compute quadrature by eigenvalue solve
vals,vecs = LA.eigh(J)
perm = np.argsort(vals)
x = vals[perm]
vecs = vecs[:,perm]
w = np.power(np.transpose(vecs[0,:]),2) * 2**(alpha+beta+1)/(alpha+beta+1)*gammaF(alpha+1)*gammaF(beta+1)/gammaF(alpha+beta+1)
return (np.asarray(x),np.asarray(w))
def __JacobiGL(self,N):
"""
x = JacobiGL
Purpose: Compute the N'th order Gauss Lobatto quadrature points, x,
and weights, w, associated with the Jacobi
polynomial, of type (alpha,beta) > -1 ( <> -0.5).
"""
if (self.poly != JACOBI):
print "The method __JacobiGL cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,0.0) or self.__AlmostEqual(beta,0.0) ):
return self.__JacobiLGL(N)
elif ( self.__AlmostEqual(alpha,-0.5) or self.__AlmostEqual(beta,-0.5) ):
return self.__JacobiCGL(N)
x = np.mat(np.zeros((N+1)))
if ( N == 1 ):
x[0] = -1.
x[1] = 1.
return x
[xint,w] = self.JacobiGQ(alpha+1, beta+1, N-2)
x = np.concatenate(([-1.], xint, [1.]))
return x
def __JacobiLGL(self,N):
"""
x,w = JacobiLGL(N)
Compute the Legendre Gauss Lobatto points and weights for polynomials up
to degree N
Algorithm (25) taken from [1]_
"""
if (self.poly != JACOBI):
print "The method __JacobiLGL cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( not self.__AlmostEqual(alpha,0.0) or not self.__AlmostEqual(beta,0.0) ):
print "The method can be called only for Legendre Polynomials. Actual values of alpha and beta: %f, %f" % (alpha,beta)
else:
maxNewtonIter = 1000
NewtonTOL = 1e-12
x = np.zeros((N+1))
w = np.zeros((N+1))
if ( N == 1 ):
x[0] = -1.
x[1] = 1.
w[0] = 1.
w[1] = 1.
else:
x[0] = -1.
w[0] = 2./(N * (N + 1.))
x[N] = 1.
w[N] = w[0]
for j in range(1,int(np.floor((N+1.)/2.)-1) + 1):
x[j] = -np.cos( (j + 1./4.)*np.pi/N - 3./(8.*N*np.pi) * 1./(j+1./4.) )
# Newton iteratio for getting the point
k = 0
delta = 10. * NewtonTOL
while ( (k < maxNewtonIter) and (abs(delta) > NewtonTOL * abs(x[j])) ):
k = k + 1
q,dq,LN = self.__qLEvaluation(N,x[j])
delta = -q/dq
x[j] = x[j] + delta
q,dq,LN = self.__qLEvaluation(N,x[j])
x[N-j] = -x[j]
w[j] = 2./(N*(N+1)*LN**2)
w[N-j] = w[j]
if ( np.remainder(N,2) == 0 ):
q,dq,LN = self.__qLEvaluation(N,0.)
x[N/2] = 0
w[N/2] = 2./(N*(N+1)*LN**2)
return (np.asarray(x),np.asarray(w))
def __JacobiCGQ(self,N):
"""
Compute the Chebyshev Gauss points and weights for polynomials up
to degree N
Algorithm (26) taken from [1]_
"""
if (self.poly != JACOBI):
print "The method __JacobiCGQ cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( not self.__AlmostEqual(alpha,-0.5) or not self.__AlmostEqual(beta,-0.5) ):
print "The method can be called only for Chebyshev Polynomials. Actual values of alpha and beta: %f, %f" % (alpha,beta)
else:
x = np.zeros((N+1))
w = np.zeros((N+1))
for j in range(0,N+1):
x[j] = -np.cos( (2.*j + 1.)/(2.*N + 2.) * np.pi )
w[j] = np.pi / (N + 1.)
return (np.asarray(x),np.asarray(w))
def __JacobiCGL(self,N):
"""
x,w = JacobiCL(N)
Compute the Chebyshev Gauss Lobatto points and weights for polynomials up
to degree N
Algorithm (27) taken from [1]_
"""
if (self.poly != JACOBI):
print "The method __JacobiCGL cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( not self.__AlmostEqual(alpha,-0.5) or not self.__AlmostEqual(beta,-0.5) ):
print "The method can be called only for Chebyshev Polynomials. Actual values of alpha and beta: %f, %f" % (alpha,beta)
else:
x = np.zeros((N+1))
w = np.zeros((N+1))
for j in range(0,N+1):
x[j] = -np.cos(float(j)/float(N) * np.pi)
w[j] = np.pi/float(N)
w[0] = w[0]/2.
w[N] = w[N]/2.
return (np.asarray(x),np.asarray(w))
def __HermitePGQ(self,N):
"""
Compute the Hermite-Gauss quadrature points and weights for Hermite Physicists
polynomials up to degree N
For further details see [2]_
"""
if (self.poly != HERMITEP):
print "The method __HermitePGQ cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
j = np.asarray(range(1,N+1))
b = np.sqrt(j / 2.)
D = np.diag(b,1)
D = D + D.T
x,vec = np.linalg.eig(D)
x = np.sort(x)
hp = self.__HermiteP(x,N)
w = np.sqrt(np.pi) * 2.**N * factorial(N) / ((N+1) * hp**2.)
return (x,w)
def __HermiteFGQ(self,N):
"""
Compute the Hermite-Gauss quadrature points and weights for the Hermite
functions up to degree N
For further details see [2]_
"""
if (self.poly != HERMITEF):
print "The method __HermiteFGQ cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
j = np.asarray(range(1,N+1))
b = np.sqrt(j / 2.)
D = np.diag(b,1)
D = D + D.T
x,vec = np.linalg.eig(D)
hf = self.__HermiteF(x,N)
w = np.sqrt(np.pi) / ((N+1) * hf**2.)
return (x,w)
def __HermiteP_Prob_GQ(self,N):
"""
Compute the Hermite-Gauss quadrature points and weights for Hermite
Probabilistic polynomials up to degree N
For further details see Golub-Welsh algorithm in [3]
"""
if (self.poly != HERMITEP_PROB):
print "The method __HermiteP_Prob_GQ cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
j = np.asarray(range(1,N+1))
b = np.sqrt(j)
D = np.diag(b,1)
D = D + D.T
x,vec = np.linalg.eig(D)
x = np.sort(x)
hp = self.__HermiteP_Prob(x,N)
w = factorial(N) / ((N+1) * hp[:,0]**2.)
return (x,w)
def __LaguerrePGQ(self,N,alpha=None):
"""
__LaguerrePGQ(): Compute the Laguerre-Gauss quadrature points and weights for Laguerre polynomials up to degree N
Syntax:
``(x,w) = __LaguerrePGQ(N,[alpha=None])``
Input:
* ``N`` = (int) Order of polynomial accuracy
* ``alpha`` = (optional,float) Laguerre constant
Output:
* ``x`` = set of Laguerre-Gauss quadrature points
* ``w`` = set of Laguerre-Gauss quadrature weights
.. note:: For further details see [2]_
"""
if (self.poly != LAGUERREP):
print "The method __LaguerrePGQ cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
if (alpha is None):
(alpha,) = self.params
# Compute points and weights
j = np.asarray(range(0,N+1))
a = 2. * j + alpha + 1.
b = - np.sqrt( j[1:] * (j[1:] + alpha) )
D = np.diag(b,1)
D = D + D.T
D = D + np.diag(a,0)
x,vec = np.linalg.eig(D)
x = np.sort(x)
lp = self.__LaguerreP(x,N).ravel()
w = gammaF(N+alpha+1.)/ ((N+alpha+1.) * factorial(N+1)) * ( x / lp**2. )
return (x,w)
def __LaguerreFGQ(self,N,alpha=None):
"""
__LaguerreFGQ(): Compute the Laguerre-Gauss quadrature points and weights for Laguerre functions up to degree N
Syntax:
``(x,w) = __LaguerreFGQ(N,[alpha=None])``
Input:
* ``N`` = (int) Order of polynomial accuracy
* ``alpha`` = (optional,float) Laguerre constant
Output:
* ``x`` = set of Laguerre-Gauss quadrature points
* ``w`` = set of Laguerre-Gauss quadrature weights
.. note:: For further details see [2]_
"""
if (self.poly != LAGUERREF):
print "The method __LaguerreFGQ cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
if (alpha is None):
(alpha,) = self.params
# Compute points and weights
j = np.asarray(range(0,N+1))
a = 2. * j + alpha + 1.
b = - np.sqrt( j[1:] * (j[1:] + alpha) )
D = np.diag(b,1)
D = D + D.T
D = D + np.diag(a,0)
x,vec = np.linalg.eig(D)
x = np.sort(x)
lf = self.__LaguerreF(x,N).ravel()
w = gammaF(N+alpha+1.)/ ((N+alpha+1.) * factorial(N+1)) * ( x / lf**2. )
return (x,w)
def __LaguerrePGR(self,N,alpha=None):
"""
__LaguerrePGR(): Compute the Laguerre-Gauss-Radau quadrature points and weights for Laguerre polynomials up to degree N
Syntax:
``(x,w) = __LaguerrePGR(N,[alpha=None])``
Input:
* ``N`` = (int) Order of polynomial accuracy
* ``alpha`` = (optional,float) Laguerre constant
Output:
* ``x`` = set of Laguerre-Gauss-Radau quadrature points
* ``w`` = set of Laguerre-Gauss-Radau quadrature weights
.. note:: For further details see [2]_
"""
if (self.poly != LAGUERREP):
print "The method __LaguerrePGR cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
if (alpha is None):
(alpha,) = self.params
# Compute points and weights x1...xN, w1...wN
j = np.asarray(range(0,N))
a = 2. * j + (alpha+1.) + 1.
b = - np.sqrt( j[1:] * (j[1:] + (alpha+1.) ) )
D = np.diag(b,1)
D = D + D.T
D = D + np.diag(a,0)
x,vec = np.linalg.eig(D)
x = np.sort(x)
lp = self.__LaguerreP(x,N).ravel()
w = gammaF(N+alpha+1.)/ ((N+alpha+1.) * factorial(N)) * ( 1. / lp**2. )
# Add x0 and w0
x = np.hstack((0.0,x))
w0 = (alpha+1.) * gammaF(alpha+1.)**2. * gammaF(N+1) / gammaF(N+alpha+2.)
w = np.hstack((w0,w))
return (x,w)
def __LaguerreFGR(self,N,alpha=None):
"""
__LaguerreFGR(): Compute the Laguerre-Gauss-Radau quadrature points and weights for Laguerre functions up to degree N
Syntax:
``(x,w) = __LaguerreFGR(N,[alpha=None])``
Input:
* ``N`` = (int) Order of polynomial accuracy
* ``alpha`` = (optional,float) Laguerre constant
Output:
* ``x`` = set of Laguerre-Gauss-Radau quadrature points
* ``w`` = set of Laguerre-Gauss-Radau quadrature weights
.. note:: For further details see [2]_
"""
if (self.poly != LAGUERREF):
print "The method __LaguerreFGR cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
if (alpha is None):
(alpha,) = self.params
# Compute points and weights x1...xN, w1...wN
j = np.asarray(range(0,N))
a = 2. * j + (alpha+1.) + 1.
b = - np.sqrt( j[1:] * (j[1:] + (alpha+1.) ) )
D = np.diag(b,1)
D = D + D.T
D = D + np.diag(a,0)
x,vec = LA.eig(D)
x = np.sort(x)
lp = self.__LaguerreF(x,N).ravel()
w = gammaF(N+alpha+1.)/ ((N+alpha+1.) * factorial(N)) * ( 1. / lp**2. )
# Add x0 and w0
x = np.hstack((0.0,x))
w0 = (alpha+1.) * gammaF(alpha+1.)**2. * gammaF(N+1) / gammaF(N+alpha+2.)
w = np.hstack((w0,w))
return (x,w)
[docs] def GaussQuadrature(self,N,normed=False):
"""
GaussQuadrature(): Generates list of nodes and weights for the Gauss quadrature rule using the selected Polynomial basis
Syntax:
``(x,w) = GaussQuadrature(N,[normed=False])``
Input:
* ``N`` = (int) accuracy level required
Output:
* ``x`` = (1d-array,float) containing the nodes
* ``w`` = (1d-array,float) containing the weights
"""
if (self.poly == JACOBI):
(x,w) = self.__JacobiGQ(N)
elif (self.poly == HERMITEP):
(x,w) = self.__HermitePGQ(N)
elif (self.poly == HERMITEF):
(x,w) = self.__HermiteFGQ(N)
elif (self.poly == HERMITEP_PROB):
(x,w) = self.__HermiteP_Prob_GQ(N)
elif (self.poly == LAGUERREP):
(x,w) = self.__LaguerrePGQ(N)
elif (self.poly == LAGUERREF):
(x,w) = self.__LaguerreFGQ(N)
if normed:
w = w / self.Gamma(0)
return (x,w)
[docs] def GaussLobattoQuadrature(self,N):
"""
GaussLobattoQuadrature(): Generates list of nodes for the Gauss-Lobatto quadrature rule using selected Polynomial basis
Syntax:
x = GaussLobattoQuadrature(N)
Input:
* N = (int) accuracy level required
Output:
* x = (1d-array,float) containing the nodes
.. note:: Available only for Jacobi Polynomials
"""
if (self.poly == JACOBI):
return self.__JacobiGL(N)
else:
print "Gauss Lobatto quadrature does not apply to the selected Polynomials/Function."
[docs] def GaussRadauQuadrature(self,N):
"""
GaussRadauQuadrature(): Generates list of nodes for the Gauss-Radau quadrature rule using selected Polynomial basis
Syntax:
``x = GaussRadauQuadrature(N)''
Input:
* ``N'' = (int) accuracy level required
Output:
* ``x'' = (1d-array,float) containing the nodes
* ``w'' = (1d-array,float) weights
.. note:: Available only for Laguerre Polynomials/Functions
"""
if (self.poly == LAGUERREP):
return self.__LaguerrePGR(N)
elif (self.poly == LAGUERREF):
return self.__LaguerreFGR(N)
else:
print "Gauss Lobatto quadrature does not apply to the selected Polynomials/Function."
def __qLEvaluation(self,N,x):
"""
q,dq,LN = qLEvaluation(N,x)
Evaluate Legendre Polynomial LN and
q = L_N+1 - L_N-1
q' = L'_N+1 - L'_N-1
at point x.
Algorithm (24) taken from [1]_
"""
L = np.zeros((N+2))
DL = np.zeros((N+2))
L[0] = 1.
L[1] = x
DL[0] = 0.
DL[1] = 1.
for k in range(2,N+2):
L[k] = (2.*k-1.)/k * x * L[k-1] - (k-1.)/k * L[k-2]
DL[k] = DL[k-2] + (2.*k-1.) * L[k-1]
q = L[N+1] - L[N-1]
dq = DL[N+1] - DL[N-1]
return (q,dq,L[N])
def __JacobiP(self,x,N,alpha=None,beta=None):
"""
Evaluate Jacobi Polynomial of type (alpha,beta) > -1
(alpha+beta <> -1) at points x for order N and returns P[1:length(xp))]
If alpha=-0.5 and beta=-0.5, then this function refers to the generation
of Chebyshev Polynomials (routine ChebyshevP)
Note: normalized to be orthonormal
"""
if (self.poly != JACOBI):
print "The method __JacobiP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
if (alpha is None) and (beta is None):
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,-0.5) and self.__AlmostEqual(beta,-0.5) ):
return self.__ChebyshevP(x,N)
else:
xp = np.mat(x);
dims = np.shape(xp)
if (dims[1] == 1): xp = np.transpose(xp)
PL = np.mat(np.zeros((N+1,np.shape(xp)[1])))
gamma0 = 2**(alpha+beta+1)/(alpha+beta+1)*gammaF(alpha+1)*gammaF(beta+1)/gammaF(alpha+beta+1)
PL[0,:] = 1.0/np.sqrt(gamma0)
if (N == 0): return np.transpose(PL)
gamma1 = (alpha+1.)*(beta+1.)/(alpha+beta+3.)*gamma0
PL[1,:] = ((alpha+beta+2.)*xp/2. + (alpha-beta)/2.)/np.sqrt(gamma1)
if (N == 1): return np.transpose(PL[N,:])
# Recurrence
aold = 2./(2.+alpha+beta)*np.sqrt((alpha+1.)*(beta+1.)/(alpha+beta+3.))
for i in range(1,N):
h1 = 2.*i+alpha+beta
anew = 2./(h1+2.)*np.sqrt( (i+1.)*(i+1.+alpha+beta)*(i+1.+alpha)*(i+1.+beta)/(h1+1.)/(h1+3.) )
bnew = - (alpha**2. - beta**2.)/h1/(h1+2.)
PL[i+1,:] = 1.0 / anew*( -aold*PL[i-1,:] + np.multiply((xp-bnew),PL[i,:]) )
aold = anew
return np.transpose(PL[N,:])
def __LegendreP(self,r, N):
"""
P = LegendreP(r,N)
Returns an 1d-array with the Legendre polynomial of order N at points r.
This calls JacobiP(r,0.,0.,N)
Note: normalized to be orthonormal
"""
if (self.poly != JACOBI):
print "The method __LegendreP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,0.0) and self.__AlmostEqual(beta,0.0) ):
return self.__JacobiP(r,N,0.,0.)
def __ChebyshevP(self,r, N):
"""
T = ChebyshevP(r, N)
Returns an 1d-array with the Chebyshev (first type) polynomial
of order N at points r
Algorithm (21) taken from [1]_
"""
if (self.poly != JACOBI):
print "The method __ChebyshevP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,-0.5) and self.__AlmostEqual(beta,-0.5) ):
# shape r
r = np.reshape(np.array(r),(r.shape[0],1))
Ks = 50
rShape = r.shape
if ( N == 0 ): return np.ones(rShape)
if ( N == 1 ): return r
if ( N <= Ks):
T2 = np.ones(rShape)
T1 = r
for j in range(2,N+1):
T = 2 * r * T1 - T2
T2 = T1
T1 = T
else:
T = np.cos(N * np.arccos(r) )
return T
def __HermiteP(self,r, N):
"""
Returns the N-th Hermite Physicist polynomial using the recurrence relation
"""
if (self.poly != HERMITEP):
print "The method __HermiteP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
old2 = np.ones( (len(r),1) )
if (N == 0):
return old2
rs = np.reshape(r,(len(r),1))
old1 = 2. * rs
if (N == 1):
return old1
new = 2. * rs * old1 - 2. * old2
for i in xrange(3,N+1):
old2 = old1
old1 = new
new = 2. * rs * old1 - 2. * (i-1) * old2
return new
def __HermitePnorm(self,r, N):
"""
Returns the N-th Hermite Physicist normalized polynomial using the
recurrence relation in [3]
"""
if (self.poly != HERMITEP):
print "The method __HermitePnorm cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
old2 = np.zeros( (len(r),1) )
if (N == -1):
return old2
old1 = np.ones( (len(r),1) ) / (np.pi**(1./4.))
if (N == 0):
return old1
rs = np.reshape(r,(len(r),1))
new = np.sqrt(2) * rs * old1
for i in xrange(2,N+1):
old2 = old1
old1 = new
new = np.sqrt(2./i) * rs * old1 - np.sqrt((i-1.)/i) * old2
return new
def __HermiteF(self,r, N):
"""
Returns the N-th Hermite function using the recurrence relation
Reference: [2]
"""
if (self.poly != HERMITEF):
print "The method __HermiteF cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
r = np.reshape(r,(len(r),1))
old2 = np.exp(-r**2. / 2.)
if (N == 0):
return old2
old1 = np.sqrt(2.) * r * np.exp(-r**2. / 2.)
if (N == 1):
return old1
new = r * old1 - np.sqrt(1./2.) * old2
for i in xrange(3,N+1):
old2 = old1
old1 = new
new = r * np.sqrt(2./i) * old1 - np.sqrt((i-1)/i) * old2
return new
def __HermiteP_Prob(self,r, N):
"""
Returns the N-th Hermite Probabilistic polynomial using the recurrence relation
Use the Probabilistic Hermite Polynomial
"""
if (self.poly != HERMITEP_PROB):
print "The method __HermiteP_Prob cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
old2 = np.ones( (len(r),1) )
if (N == 0):
return old2
rs = np.reshape(r,(len(r),1))
old1 = rs
if (N == 1):
return old1
new = rs * old1 - old2
for i in xrange(3,N+1):
old2 = old1
old1 = new
new = rs * old1 - (i-1) * old2
return new
def __LaguerreP(self,r,N,alpha=None):
"""
__LaguerreP(): Generates the N-th Laguerre polynomial using the recurrence relation
Syntax:
``__LaguerreP(r,N,[alpha=None])``
Input:
* ``r`` = (1d-array,float) set of points
* ``N`` = (int) Polynomial order
* ``alpha`` = (optional,float) Laguerre parameter
Output:
* ``P`` = (1d-array,float) ``N``-th Laguerre polynomial evaluated on ``r``
"""
if (self.poly != LAGUERREP):
print "The method __LaguerreP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
if (alpha == None):
(alpha,) = self.params
# Recurrence relation
old2 = np.ones( (len(r),1) )
if (N == 0):
return old2
rs = np.reshape(r,(len(r),1))
old1 = alpha + 1. - rs
if (N == 1):
return old1
new = ( (3. + alpha - rs) * old1 - (1. + alpha) * old2 ) / 2.0
for i in xrange(3,N+1):
n = i - 1.
old2 = old1
old1 = new
new = ( (2*n + alpha + 1. - rs) * old1 - (n + alpha) * old2 ) / (n + 1.)
return new
def __LaguerreF(self,r,N,alpha=None):
"""
__LaguerreF(): Generates the N-th Laguerre function using the recurrence relation
Syntax:
``__LaguerreF(r,N,[alpha=None])``
Input:
* ``r`` = (1d-array,float) set of points
* ``N`` = (int) Polynomial order
* ``alpha`` = (optional,float) Laguerre parameter
Output:
* ``P`` = (1d-array,float) ``N``-th Laguerre function evaluated on ``r``
"""
if (self.poly != LAGUERREF):
print "The method __LaguerreF cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
if (alpha == None):
(alpha,) = self.params
# Recurrence relation
rs = np.reshape(r,(len(r),1))
old2 = np.exp(-rs/2.)
if (N == 0):
return old2
old1 = (alpha + 1. - rs) * np.exp(-rs/2.)
if (N == 1):
return old1
new = ( (3. + alpha - rs) * old1 - (1. + alpha) * old2 ) / 2.0
for i in xrange(3,N+1):
n = i - 1.
old2 = old1
old1 = new
new = ( (2*n + alpha + 1. - rs) * old1 - (n + alpha) * old2 ) / (n + 1.)
return new
def __GradLaguerreP(self,r,N,k,alpha=None):
"""
__GradLaguerreP(): Generates the k-th derivative of the N-th Laguerre polynomial using the recurrence relation
Syntax:
``__LaguerreP(r,N,k,[alpha=None])``
Input:
* ``r`` = (1d-array,float) set of points
* ``N`` = (int) Polynomial order
* ``k`` = (int) derivative order
* ``alpha`` = (optional,float) Laguerre parameter
Output:
* ``P`` = (1d-array,float) ``k``-th derivative of the ``N``-th Laguerre polynomial evaluated on ``r``
"""
if (self.poly != LAGUERREP):
print "The method __LaguerreP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
if (alpha == None):
(alpha,) = self.params
if (k == 0):
dP = self.__LaguerreP(r,N,alpha)
else:
if (N == 0):
dP = np.zeros(r.shape)
else:
dP = - self.__GradLaguerreP(r,N-1,k-1,alpha+1.)
return dP
def __GradLaguerreF(self,r,N,k,alpha=None):
"""
__GradLaguerreF(): Generates the k-th derivative of the N-th Laguerre function using the recurrence relation
Syntax:
``__LaguerreF(r,N,k,[alpha=None])``
Input:
* ``r`` = (1d-array,float) set of points
* ``N`` = (int) Polynomial order
* ``k`` = (int) derivative order
* ``alpha`` = (optional,float) Laguerre parameter
Output:
* ``P`` = (1d-array,float) ``k``-th derivative of the ``N``-th Laguerre function evaluated on ``r``
"""
if (self.poly != LAGUERREF):
print "The method __LaguerreF cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
if (alpha == None):
(alpha,) = self.params
if (k == 0):
dP = self.__LaguerreF(r,N,alpha)
else:
if (N == 0):
dP = -0.5 * self.__GradLaguerreF(r,N,k-1,alpha)
else:
dP = - self.__GradLaguerreF(r,N-1,k-1,alpha+1.) - 0.5 * self.__GradLaguerreF(r,N,k-1,alpha)
return dP
def __GradHermiteP(self,r,N,k):
"""
Compute the first derivative of the N-th Hermite Physicist Polynomial
using the recursion relation in [2]
"""
if (self.poly != HERMITEP):
print "The method __GradHermiteP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
if (k == 0):
dP = self.__HermiteP(r,N)
else:
if (N == 0):
dP = np.zeros(r.shape)
else:
dP = 2. * N * self.__GradHermiteP(r,N-1,k-1)
return dP
def __GradHermitePnorm(self,r,N,k):
"""
Compute the first derivative of the N-th Hermite Physicist Normalized Polynomial
using the recursion relation in [3]
"""
if (self.poly != HERMITEP):
print "The method __GradHermitePnorm cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
if (k == 0):
dP = self.__HermitePnorm(r,N)
else:
if (N == 0):
dP = np.zeros(r.shape)
else:
dP = np.sqrt(2.*N) * self.__GradHermitePnorm(r,N-1,k-1)
return dP
def __GradHermiteF(self,r,N,k):
"""
Compute the first derivative of the N-th Hermite Function using the
recursion relation in [2]
"""
if (self.poly != HERMITEF):
print "The method __GradHermiteF cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
if (k == 0):
dP = self.__HermiteF(r,N)
else:
if (N == 0):
# Compute the derivative of the first Hermite function using
# direct derivatives (See formula in the notes)
mat = np.zeros((k+1,k+1))
mat[0,0] = 1
for i in range(1,k+1):
for j in range(np.remainder(i,2),k+1,2):
if (j != 0):
mat[i,j] = mat[i,j] - mat[i-1,j-1]
if (j < i):
mat[i,j] = mat[i,j] + mat[i-1,j+1] * (j+1)
dP = np.zeros(r.shape)
for i in range(0,k+1):
dP = dP + mat[-1,i] * r**i * np.exp(-r**2./2.)
else:
dP = np.sqrt(N/2.) * self.__GradHermiteF(r,N-1,k-1) - np.sqrt((N+1.)/2.) * self.__GradHermiteF(r,N+1,k-1)
return dP
def __GradHermiteP_Prob(self,r,N,k):
"""
Compute the first derivative of the N-th Hermite Probabilistic Polynomial
"""
if (self.poly != HERMITEP_PROB):
print "The method __GradHermiteP_Prob cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
if (k == 0):
dP = self.__HermiteP_Prob(r,N)
else:
if (N == 0):
dP = np.zeros(r.shape)
else:
dP = N * self.__GradHermiteP_Prob(r,N-1,k-1)
return dP
def __GradJacobiP(self,r, N, k):
"""
dP = GradJacobiP(r, alpha, beta, N);
Purpose: Evaluate the kth-derivative of the Jacobi polynomial of type (alpha,beta)>-1,
at points r for order N and returns dP[1:length(r))]
"""
if (self.poly != JACOBI):
print "The method __GradJacobiP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,-0.5) and self.__AlmostEqual(beta,-0.5) ):
return self.__GradChebyshevP(r, N, k, 0)
else:
r = np.array(r)
dP = np.mat(np.zeros((r.shape[0]))).T
if (N >= k):
dP = gammaF(alpha+beta+N+1.+k)/(2.**k * gammaF(alpha+beta+N+1.)) * np.sqrt(self.__GammaJacobiP(N-k,alpha+k,beta+k)/self.__GammaJacobiP(N,alpha,beta)) * self.__JacobiP(r,N-k,alpha+k,beta+k)
return dP
def __GradChebyshevP(self,r, N, k, method=0):
"""
dP = GradChebyshevP(r,N,k,method)
Returns the k-th derivative of the Chebyshev polynomial of order N at
points r.
Method: 0 -> Matrix multiplication
1 -> Fast Chebyshev Transform
"""
if (self.poly != JACOBI):
print "The method __GradChebyshevP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,-0.5) and self.__AlmostEqual(beta,-0.5) ):
dP = np.zeros((N+1))
if ( k == 0 ):
dP = self.__ChebyshevP(r,N)
elif ( method == 0 ):
D = self.PolynomialDerivativeMatrix(r,k)
P = self.__ChebyshevP(r,N)
dP = np.dot(D,P)
return dP
def __GradORTHPOL(self,r,N,k):
"""
dP = GradORTHPOL(r,N,k)
Returns the 0-th derivative of the N-th orthogonal polynomial defined for the supplied measure
"""
# TODO: finish here.
[docs] def GradEvaluate(self,r,N,k):
"""
GradEvaluate(): evaluate the ``k``-th derivative of the ``N``-th order polynomial at points ``r``
Syntax:
``P = GradEvaluate(r,N,k)``
Input:
* ``r`` = (1d-array,float) set of points on which to evaluate the polynomial
* ``N`` = (int) order of the polynomial
* ``k`` = (int) order of the derivative
Output:
* ``P`` = Polynomial evaluated on ``r``
"""
if (self.poly == JACOBI):
return self.__GradJacobiP(r,N,k)
elif (self.poly == HERMITEP):
return self.__GradHermiteP(r,N,k)
elif (self.poly == HERMITEF):
return self.__GradHermiteF(r,N,k)
elif (self.poly == HERMITEP_PROB):
return self.__GradHermiteP_Prob(r,N,k)
elif (self.poly == LAGUERREP):
return self.__GradLaguerreP(r,N,k)
elif (self.poly == LAGUERREF):
return self.__GradLaguerreF(r,N,k)
elif (self.poly == ORTHPOL):
if k != 0:
print "DABISpectral1D: Error. Derivatives of Polynomials obtained using ORTHPOL package are not implemented"
return
return self.__GradORTHPOL(r,N)
def __GammaLaguerreF(self,N,alpha=None):
"""
__GammaLaguerreF(): evaluate the normalization constant for the ``N``-th order Laguerre function
Syntax:
``g = __GammaLaguerreF(N,[alpha=None])``
Input:
* ``N`` = (int) order of the polynomial
* ``alpha'' = (optional,float) Laguerre constant
Output:
* ``g`` = Normalization constant
"""
if (self.poly != LAGUERREF):
print "The method __GammaLaguerreF cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
if (alpha is None):
# Unpack parameters
(alpha,) = self.params
g = gammaF(N+alpha+1.)/gammaF(N+1)
return g
def __GammaLaguerreP(self,N,alpha=None):
"""
__GammaLaguerreP(): evaluate the normalization constant for the ``N``-th order Laguerre polynomial
Syntax:
``g = __GammaLaguerreP(N,[alpha=None])``
Input:
* ``N`` = (int) order of the polynomial
* ``alpha`` = (optional,float) Laguerre constant
Output:
* ``g`` = Normalization constant
"""
if (self.poly != LAGUERREP):
print "The method __GammaLaguerreP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
if (alpha is None):
# Unpack parameters
(alpha,) = self.params
g = gammaF(N+alpha+1.)/gammaF(N+1)
return g
def __GammaJacobiP(self,N,alpha=None,beta=None):
"""
gamma = GammaJacobiP(alpha,beta,N)
Generate the normalization constant for the
Jacobi Polynomial (alpha,beta) of order N.
"""
if (self.poly != JACOBI):
print "The method __GammaJacobiP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
if (alpha is None) and (beta is None):
# Unpack parameters
alpha,beta = self.params
g = 2**(alpha+beta+1.) * (gammaF(N+alpha+1.)*gammaF(N+beta+1.)) / (factorial(N,exact=True) * (2.*N + alpha + beta + 1.)*gammaF(N+alpha+beta+1.))
return g
def __GammaHermiteP(self,N):
"""
Returns the normalization contant for the Probabilistic Hermite Physicist
polynomial of order N
"""
if (self.poly != HERMITEP):
print "The method __GammaHermiteP cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
return np.sqrt(np.pi) * 2.**N * factorial(N,exact=True)
def __GammaHermiteF(self,N):
"""
Returns the normalization contant for the Hermite function of order N
"""
if (self.poly != HERMITEF):
print "The method __GammaHermiteF cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
return np.sqrt(2.*np.pi)
def __GammaHermiteP_Prob(self,N):
"""
Returns the normalization contant for the Probabilistic Hermite
Probabilistic polynomial of order N
"""
if (self.poly != HERMITEP_PROB):
print "The method __GammaHermiteP_Prob cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
return factorial(N,exact=True)
[docs] def Gamma(self,N):
"""
Gamma(): returns the normalization constant for the N-th polynomial
Syntax:
``g = Gamma(N)``
Input:
* ``N`` = polynomial order
Output:
* ``g`` = normalization constant
"""
if (self.poly == JACOBI):
return self.__GammaJacobiP(N)
elif (self.poly == HERMITEP):
return self.__GammaHermiteP(N)
elif (self.poly == HERMITEF):
return self.__GammaHermiteF(N)
elif (self.poly == HERMITEP_PROB):
return self.__GammaHermiteP_Prob(N)
elif (self.poly == LAGUERREP):
return self.__GammaLaguerreP(N)
elif (self.poly == LAGUERREF):
return self.__GammaLaguerreF(N)
[docs] def GradVandermonde1D(self,r,N,k,norm=True):
"""
GradVandermonde1D(): Initialize the ``k``-th gradient of the modal basis ``N`` at ``r``
Syntax:
``V = GradVandermonde1D(r,N,k,[norm])``
Input:
* ``r`` = (1d-array,float) set of ``M`` points on which to evaluate the polynomials
* ``N`` = (int) maximum order in the vanermonde matrix
* ``k`` = (int) derivative order
* ``norm`` = (optional,boolean) True -> orthonormal polynomials, False -> non orthonormal polynomials
Output:
* ``V`` = (2d-array(``MxN``),float) Generalized Vandermonde matrix
"""
DVr = np.mat(np.zeros((r.shape[0],N+1)))
for i in range(0,N+1):
DVr[:,i] = self.GradEvaluate(r,i,k)
DVr = np.asarray(DVr)
# Normalization
if (self.poly == JACOBI):
if (not norm):
# return non orthonormal basis
for i in range(0,N+1):
DVr[:,i] = DVr[:,i] * np.sqrt(self.__GammaJacobiP(i))
else:
if norm:
# return orthornormal basis
for i in range(0,N+1):
DVr[:,i] = DVr[:,i] / np.sqrt(self.Gamma(i))
return np.asarray(DVr)
[docs] def AssemblyDerivativeMatrix(self, x, N, k):
"""
AssemblyDerivativeMatrix(): Assemble the k-th derivative matrix using polynomials of order N.
Syntax:
``Dk = AssemblyDerivativeMatrix(x,N,k)``
Input:
* x = (1d-array,float) Set of points on which to evaluate the polynomials
* N = (int) maximum order in the vanermonde matrix
* k = (int) derivative order
Output:
* Dk = Derivative matrix
Description:
This function performs ``D = linalg.solve(V.T, Vx.T)`` where ``V`` and ``Vx`` are a Generalized Vandermonde Matrix and its derivative respectively.
Notes:
For Chebyshev Polynomial, this function refers to the recursion form implemented in ``PolynomialDerivativeMatrix``
"""
# Unpack parameters
if (self.poly == JACOBI):
alpha,beta = self.params
if (self.poly == JACOBI) and ( self.__AlmostEqual(alpha,-0.5) and self.__AlmostEqual(beta,-0.5) ):
return self.PolynomialDerivativeMatrix(x,k)
else:
V = self.GradVandermonde1D(x, N, 0)
Vx = self.GradVandermonde1D(x, N ,1)
D = LA.solve(V.T, Vx.T)
D = D.T
Dk = np.asarray(np.mat(D)**k)
return Dk
[docs] def FirstPolynomialDerivativeMatrix(self,x):
"""
PolynomialDerivativeMatrix(): Assemble the first Polynomial Derivative Matrix using matrix multiplication.
Syntax:
``D = FirstPolynomialDerivativeMatrix(x)``
Input:
* ``x`` = (1d-array,float) set of points on which to evaluate the derivative matrix
Output:
* ``D`` = derivative matrix
Notes:
Algorithm (37) from [1]_
"""
N = x.shape[0]
w = self.BarycentricWeights(x)
D = np.zeros((N,N))
for i in range(0,N):
for j in range(0,N):
if (j != i):
D[i,j] = w[j]/w[i] * 1./(x[i] - x[j])
D[i,i] = D[i,i] - D[i,j]
return D
[docs] def PolynomialDerivativeMatrix(self,x,k):
"""
PolynomialDerivativeMatrix(): Assemble the Polynomial ``k``-th Derivative Matrix using the matrix recursion. This algorithm is generic for every types of polynomials.
Syntax:
``D = PolynomialDerivativeMatrix(x,k)``
Input:
* ``x`` = (1d-array,float) set of points on which to evaluate the derivative matrix
* ``k`` = derivative order
Output:
* ``D`` = ``k``-th derivative matrix
Notes:
Algorithm (38) taken from [1]_
"""
N = x.shape[0]
w = self.BarycentricWeights(x)
D = np.zeros((N,N,k))
D[:,:,0] = self.FirstPolynomialDerivativeMatrix(x)
if ( k == 1 ): return D[:,:,k-1]
for m in range(2,k+1):
for i in range(0,N):
for j in range(0,N):
if ( j != i ):
D[i,j,m-1] = m/(x[i]-x[j]) * ( w[j]/w[i] * D[i,i,m-2] - D[i,j,m-2])
D[i,i,m-1] = D[i,i,m-1] - D[i,j,m-1]
return D[:,:,k-1]
[docs] def BarycentricWeights(self,x):
"""
BarycentricWeights(): Returns a 1-d array of weights for Lagrange Interpolation
Syntax:
``w = BarycentricWeights(x)``
Input:
* ``x`` = (1d-array,float) set of points
Output:
* ``w`` = (1d-array,float) set of barycentric weights
Notes:
Algorithm (30) from [1]_
"""
N = x.shape[0]
w = np.zeros((N))
for j in range(0,N):
w[j] = 1.
for j in range(1,N):
for k in range(0,j):
w[k] = w[k] * (x[k] - x[j])
w[j] = w[j] * (x[j] - x[k])
for j in range(0,N):
w[j] = 1. / w[j]
return w
[docs] def LagrangeInterpolationMatrix(self,x, w, xi):
"""
LagrangeInterpolationMatrix(): constructs the Lagrange Interpolation Matrix from points ``x`` to points ``xi``
Syntax:
``T = LagrangeInterpolationMatrix(x, w, xi)``
Input:
* ``x`` = (1d-array,float) set of ``N`` original points
* ``w`` = (1d-array,float) set of ``N`` barycentric weights
* ``xi`` = (1d-array,float) set of ``M`` interpolating points
Output:
* ``T`` = (2d-array(``MxN``),float) Lagrange Interpolation Matrix
Notes:
Algorithm (32) from [1]_
"""
N = x.shape[0]
M = xi.shape[0]
T = np.zeros((M,N))
for k in range(0,M):
rowHasMatch = False
for j in range(0,N):
T[k,j] = 0.
if self.__AlmostEqual(xi[k],x[j]):
rowHasMatch = True
T[k,j] = 1.
if (rowHasMatch == False):
s = 0.
for j in range(0,N):
t = w[j] / (xi[k] - x[j])
T[k,j] = t
s = s + t
for j in range(0,N):
T[k,j] = T[k,j] / s
return T
[docs] def LagrangeInterpolate(self, x, f, xi):
"""
LagrangeInterpolate(): Interpolate function values ``f`` from points ``x`` to points ``xi`` using Lagrange weights
Syntax:
``fi = LagrangeInterpolate(x, f, xi)``
Input:
* ``x`` = (1d-array,float) set of ``N`` original points where ``f`` is evaluated
* ``f`` = (1d-array,float) set of ``N`` function values
* ``xi`` = (1d-array,float) set of ``M`` points where the function is interpolated
Output:
* ``fi`` = (1d-array,float) set of ``M`` function values
Notes:
Modification of Algorithm (33) from [1]_
"""
# Obtain barycentric weights
w = self.BarycentricWeights(x)
# Generate the Interpolation matrix
T = self.LagrangeInterpolationMatrix(x, w, xi)
# Compute interpolated values
fi = np.dot(T,f)
return fi
[docs] def PolyInterp(self, x, f, xi, order):
"""
PolyInterp(): Interpolate function values ``f`` from points ``x`` to points ``xi`` using Forward and Backward Polynomial Transform
Syntax:
``fi = PolyInterp(x, f, xi)``
Input:
* ``x`` = (1d-array,float) set of ``N`` original points where ``f`` is evaluated
* ``f`` = (1d-array,float) set of ``N`` function values
* ``xi`` = (1d-array,float) set of ``M`` points where the function is interpolated
* ``order`` = (integer) order of polynomial interpolation
Output:
* ``fi`` = (1d-array,float) set of ``M`` function values
Notes:
"""
(fhat, residues, rank, s) = self.DiscretePolynomialTransform(x,f,order)
return self.InverseDiscretePolynomialTransform(xi,fhat,order)
[docs] def LegendreDerivativeCoefficients(self,fhat):
"""
LegendreDerivativeCoefficients(): computes the Legendre coefficients of the derivative of a function
Syntax:
``dfhat = LegendreDerivativeCoefficients(fhat)``
Input:
* ``fhat`` = (1d-array,float) list of Legendre coefficients of the original function
Output:
* ``dfhat`` = (1d-array,float) list of Legendre coefficients of the derivative of the original function
Notes:
Algorithm (4) from [1]_
"""
if (self.poly != JACOBI):
print "The method cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,0.) and self.__AlmostEqual(beta,0.) ):
N = fhat.shape[0]-1
dfhat = np.zeros((N+1))
dfhat[N-1] = (2.*N - 1.) * fhat[N]
for k in range(N-2, -1, -1):
dfhat[k] = (2.*k + 1.) * (fhat[k+1] + dfhat[k+2]/(2.*k + 5.) )
return dfhat
[docs] def ChebyshevDerivativeCoefficients(self,fhat):
"""
ChebyshevDerivativeCoefficients(): computes the Chebyshev coefficients of the derivative of a function
Syntax:
``dfhat = ChebyshevDerivativeCoefficients(fhat)``
Input:
* ``fhat`` = (1d-array,float) list of Chebyshev coefficients of the original function
Output:
* ``dfhat`` = (1d-array,float) list of Chebyshev coefficients of the derivative of the original function
Notes:
Algorithm (5) from [1]_
"""
if (self.poly != JACOBI):
print "The method cannot be called with the actual type of polynomials. Actual type: '%s'" % self.poly
else:
# Unpack parameters
alpha,beta = self.params
if ( self.__AlmostEqual(alpha,-0.5) and self.__AlmostEqual(beta,-0.5) ):
N = fhat.shape[0]-1
dfhat = np.zeros((N+1))
dfhat[N-1] = (2.*N) * fhat[N]
for k in range(N-2, 0, -1):
dfhat[k] = 2. * (k + 1.) * fhat[k+1] + dfhat[k+2]
dfhat[0] = fhat[1] + dfhat[2]/2.
return dfhat
[docs] def GramSchmidt(self, p, N, w):
"""
GramSchmidt(): creates a Generalized Vandermonde Matrix of orthonormal polynomials with respect to the weights ``w``
Syntax:
``V = GramSchmidt(p, N, w)``
Input:
* ``p`` = (1d-array,float) points at which to evaluate the new polynomials
* ``N`` = (int) the maximum order of the polynomials
* ``w`` = (1d-array,float) weights to be used for the orthogonoalization
Output:
* ``V`` = Generalized Vandermonde Matrix containing the new orthogonalized polynomials
Description:
Takes the points where the polynomials have to be evaluated and computes a Generalized Gram Schmidth procedure, where a weighted projection is used. If ``w==1`` then the usual inner product is used for the orthogonal projection.
"""
# Evaluate Vandermonde matrix
V = np.vander(p,N+1)
V = V[:,::-1]
# Evaluate Gram-Shmidt orthogonalization
gs = np.zeros(N+1) # Vector of gamma values for the new polynomials
for k in range(0,N+1):
for j in range(0,N+1):
for i in range(0,j):
# Use numerical quadrature to evaluate the projection
V[:,j] = V[:,j] - np.dot(V[:,j] * V[:,i], w) / np.sqrt(gs[i]) * V[:,i]
# Compute the gamma value for the new polynomial
gs[j] = np.dot(V[:,j]*V[:,j],w)
# Normalize the vector if required
V[:,j] = V[:,j]/np.sqrt(gs[j])
return V
def TestOrthogonality( V, w):
N = V.shape[1]
orth = np.zeros((N,N))
for i in range(0,N):
for j in range(0,N):
# Use numerical quadrature to compute the orthogonality constants
orth[i,j] = np.dot( V[:,i]* V[:,j], w)
return orth
def VandermondeAB(x,n,alpha,beta):
if (len(alpha) <= n or len(beta) <= n):
print "The number of arguments provided in alpha and beta are not consistent with n."
else:
V = np.zeros((len(x),n+1));
if n >= 0:
V[:,0] = 1;
if n >= 1:
V[:,1] = (x - alpha[0] );
for i in range(2,n+1):
V[:,i] = (x - alpha[i-1]) * V[:,i-1] - beta[i-1] * V[:,i-2];
return V