Package 'orthopolynom'

Title: Collection of Functions for Orthogonal and Orthonormal Polynomials
Description: A collection of functions to construct sets of orthogonal polynomials and their recurrence relations. Additional functions are provided to calculate the derivative, integral, value and roots of lists of polynomial objects.
Authors: Frederick Novomestky <[email protected]>
Maintainer: Frederick Novomestky <[email protected]>
License: GPL (>= 2)
Version: 1.0-6.1
Built: 2024-11-06 06:17:00 UTC
Source: CRAN

Help Index


Inner products of Chebyshev polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Chebyshev polynomial of the first kind, Ck(x)C_k \left( x\right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.c.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=CnCn={4πn08πn=0h_n = \left\langle {C_n |C_n } \right\rangle = \left\{ {\begin{array}{cc} {4\,\pi } & {n \ne 0} \\ {8\,\pi } & {n = 0} \\ \end{array} } \right.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

gegenbauer.inner.products

Examples

###
### generate the inner products vector for the
### C Chebyshev polynomials of orders 0 to 10
###
h <- chebyshev.c.inner.products( 10 )
print( h )

Create list of Chebyshev polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Chebyshev polynomials of the first kind, Ck(x)C_k \left( x\right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.c.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function chebyshev.c.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Chebyshev polynomial

2

order 1 Chebyshev polynomial

...

n+1

order nn Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

chebyshev.c.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized C Chebyshev polynomials of orders 0 to 10
###
normalized.p.list <- chebyshev.c.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized C Chebyshev polynomials of orders 0 to 10
###
unnormalized.p.list <- chebyshev.c.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Chebyshev polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Chebyshev polynomial of the first kind, Ck(x)C_k \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.c.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

chebyshev.c.inner.products

Examples

###
### generate the recurrences data frame for 
### the normalized Chebyshev C polynomials
### of orders 0 to 10.
###
normalized.r <- chebyshev.c.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the normalized Chebyshev C polynomials
### of orders 0 to 10.
###
unnormalized.r <- chebyshev.c.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the Chebyshev polynomial

Description

This function returns the value of the weight function for the order kk Chebyshev polynomial of the first kind, Ck(x)C_k \left( x \right).

Usage

chebyshev.c.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (2,2)\left( -2,2 \right). The formula used to compute the weight function is as follows.

w(x)=11x24w\left( x \right) = \frac{1} {{\sqrt {1 - \frac{{x^2 }} {4}} }}

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the C Chebyshev weight function for arguments between -3 and 3
###
x <- seq( -3, 3, .01 )
y <- chebyshev.c.weight( x )
plot( x, y )

Inner products of Chebyshev polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Chebyshev polynomial of the second kind, Sk(x)S_k \left( x\right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.s.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=SnSn=πh_n = \left\langle {S_n |S_n } \right\rangle = \pi.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the 
### S Chebyshev polynomials of orders 0 to 10
###
h <- chebyshev.s.inner.products( 10 )
print( h )

Create list of Chebyshev polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Chebyshev polynomials of the second kind, Sk(x)S_k \left( x\right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.s.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function chebyshev.s.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Chebyshev polynomial

2

order 1 Chebyshev polynomial

...

n+1

order nn Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

chebyshev.s.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized S Chebyshev polynomials of orders 0 to 10
###
normalized.p.list <- chebyshev.s.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized S Chebyshev polynomials of orders 0 to 10
###
unnormalized.p.list <- chebyshev.s.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Chebyshev polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Chebyshev polynomial of the second kind, Sk(x)S_k \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.s.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

chebyshev.s.inner.products

Examples

###
### generate the recurrences data frame for 
### the normalized Chebyshev S polynomials
### of orders 0 to 10.
###
normalized.r <- chebyshev.s.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the normalized Chebyshev S polynomials
### of orders 0 to 10.
###
unnormalized.r <- chebyshev.s.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the Chebyshev polynomial

Description

This function returns the value of the weight function for the order kk Chebyshev polynomial of the second kind, Sk(x)S_k \left( x \right).

Usage

chebyshev.s.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (2,2)\left( -2,2 \right). The formula used to compute the weight function is as follows.

w(x)=1x24w\left( x \right) = \sqrt {1 - \frac{{x^2 }}{4}}

Value

The value of the weight function.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the S Chebyshev weight function for arguments between -2 and 2
###
x <- seq( -2, 2, .01 )
y <- chebyshev.s.weight( x )
plot( x, y )

Inner products of Chebyshev polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Chebyshev polynomial of the first kind, Tk(x)T_k \left( x\right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.t.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=TnTn={π2n0πn=0h_n = \left\langle {T_n |T_n } \right\rangle = \left\{ {\begin{array}{cc} {\frac{\pi } {2}} & {n \ne 0} \\ \pi & {n = 0} \\ \end{array} } \right.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the
### T Chybyshev polynomials of orders 0 to 10
###
h <- chebyshev.t.inner.products( 10 )
print( h )

Create list of Chebyshev polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Chebyshev polynomials of the first kind, Tk(x)T_k \left( x\right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.t.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function chebyshev.t.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Chebyshev polynomial

2

order 1 Chebyshev polynomial

...

n+1

order nn Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

chebyshev.u.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized T Chebyshev polynomials of orders 0 to 10
###
normalized.p.list <- chebyshev.t.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized T Chebyshev polynomials of orders 0 to 10
###
unnormalized.p.list <- chebyshev.t.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Chebyshev polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Chebyshev polynomial of the first kind, Tk(x)T_k \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.t.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

chebyshev.t.inner.products

Examples

###
### generate the recurrence relations for 
### the normalized T Chebyshev polynomials
### of orders 0 to 10
###
normalized.r <- chebyshev.t.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrence relations for 
### the normalized T Chebyshev polynomials
### of orders 0 to 10
###
unnormalized.r <- chebyshev.t.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the Chebyshev polynomial

Description

This function returns the value of the weight function for the order kk Chebyshev polynomial of the first kind, Tk(x)T_k \left( x \right).

Usage

chebyshev.t.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (1,1)\left( -1,1 \right). The formula used to compute the weight function is as follows.

w(x)=11x2w\left( x \right) = \frac{1}{{\sqrt {1 - x^2 } }}

Value

The value of the weight function.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the T Chebyshev function for argument values between -2 and 2
x <- seq( -1, 1, .01 )
y <- chebyshev.t.weight( x )
plot( x, y )

Inner products of Chebyshev polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Chebyshev polynomial of the second kind, Uk(x)U_k \left( x\right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.u.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=UnUn=π2h_n = \left\langle {U_n |U_n } \right\rangle = \frac{\pi }{2}

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the 
### U Chebyshev polynomials of orders 0 to 10
###
h <- chebyshev.u.inner.products( 10 )
print( h )

Create list of Chebyshev polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Chebyshev polynomials of the second kind, Uk(x)U_k \left( x\right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.u.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function chebyshev.u.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Chebyshev polynomial

2

order 1 Chebyshev polynomial

...

n+1

order nn Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

chebyshev.u.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized U Chebyshev polynomials of orders 0 to 10
###
normalized.p.list <- chebyshev.u.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized T Chebyshev polynomials of orders 0 to 10
###
unnormalized.p.list <- chebyshev.u.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Chebyshev polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Chebyshev polynomial of the second kind, Uk(x)U_k \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

chebyshev.u.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

chebyshev.u.inner.products

Examples

###
### generate the recurrence relations for 
### the normalized U Chebyshev polynomials
### of orders 0 to 10
###
normalized.r <- chebyshev.u.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrence relations for 
### the unnormalized U Chebyshev polynomials
### of orders 0 to 10
###
unnormalized.r <- chebyshev.u.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the Chebyshev polynomial

Description

This function returns the value of the weight function for the order kk Chebyshev polynomial of the second kind, Uk(x)U_k \left( x \right).

Usage

chebyshev.u.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (1,1)\left( -1,1 \right). The formula used to compute the weight function is as follows.

w(x)=1x2w\left( x \right) = \sqrt {1 - x^2 }

Value

The value of the weight function.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the U Chebyshev function for argument values between -2 and 2
###
x <- seq( -1, 1, .01 )
y <- chebyshev.u.weight( x )
plot( x, y )

Inner products of Gegenbauer polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Gegenbauer polynomial, Ck(α)(x)C_k^{\left( \alpha \right)} \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

gegenbauer.inner.products(n,alpha)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for the polynomial parameter

Details

The formula used to compute the inner products is as follows.

hn=Cn(α)Cn(α)={π  212αΓ(n+2α)n!  (n+α)[Γ(α)]2α02  πn2α=0h_n = \left\langle {C_n^{\left( \alpha \right)} |C_n^{\left( \alpha \right)} } \right\rangle = \left\{ {\begin{array}{cc} {\frac{{\pi \;2^{1 - 2\,\alpha } \,\Gamma \left( {n + 2\,\alpha } \right)}} {{n!\;\left( {n + \alpha } \right)\,\left[ {\Gamma \left( \alpha \right)} \right]^2 }}} & {\alpha \ne 0} \\ {\frac{{2\;\pi }} {{n^2 }}} & {\alpha = 0} \\ \end{array} } \right..

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

ultraspherical.inner.products

Examples

###
### generate the inner products vector for the 
### Gegenbauer polynomials of orders 0 to 10
### the polynomial parameter is 1.0
###
h <- gegenbauer.inner.products( 10, 1 )
print( h )

Create list of Gegenbauer polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Gegenbauer polynomials, Ck(α)(x)C_k^{\left( \alpha \right)} \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

gegenbauer.polynomials(n, alpha, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

alpha

polynomial parameter

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function gegenbauer.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Gegenbauer polynomial

2

order 1 Gegenbauer polynomial

...

n+1

order nn Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

gegenbauer.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized Gegenbauer polynomials of orders 0 to 10
### polynomial parameter is 1.0
###
normalized.p.list <- gegenbauer.polynomials( 10, 1, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized Gegenbauer polynomials of orders 0 to 10
### polynomial parameter is 1.0
###
unnormalized.p.list <- gegenbauer.polynomials( 10, 1, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Gegenbauer polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Gegenbauer polynomial, Ck(α)(x)C_k^{\left( \alpha \right)} \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

gegenbauer.recurrences(n, alpha, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for polynomial parameter

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

gegenbauer.inner.products

Examples

###
### generate the recurrences data frame for 
### the normalized Gegenbauer polynomials
### of orders 0 to 10.
### polynomial parameter value is 1.0
###
normalized.r <- gegenbauer.recurrences( 10, 1, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the normalized Gegenbauer polynomials
### of orders 0 to 10.
### polynomial parameter value is 1.0
###
unnormalized.r <- gegenbauer.recurrences( 10, 1, normalized=FALSE )
print( unnormalized.r )

Weight function for the Gegenbauer polynomial

Description

This function returns the value of the weight function for the order kk Gegenbauer polynomial, Ck(α)(x)C_k^{\left( \alpha \right)} \left( x \right).

Usage

gegenbauer.weight(x,alpha)

Arguments

x

the function argument which can be a vector

alpha

polynomial parameter

Details

The function takes on non-zero values in the interval (1,1)\left( -1,1 \right). The formula used to compute the weight function is as follows.

w(x)=(1x2)α0.5w\left( x \right) = \left( {1 - x^2 } \right)^{\alpha - 0.5}

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the Gegenbauer weight function for argument values between -1 and 1
###
x <- seq( -1, 1, .01 )
y <- gegenbauer.weight( x, 1 )
plot( x, y )

Inner products of generalized Hermite polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk generalized Hermite polynomial, Hk(μ)(x)H_k^{\left( \mu \right)} \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

ghermite.h.inner.products(n, mu)

Arguments

n

n integer value for the highest polynomial order

mu

mu polynomial parameter

Details

The parameter μ\mu must be greater than -0.5. The formula used to compute the inner products is as follows.

hn(μ)=Hm(μ)Hn(μ)=22n[n2]!  Γ([n+12]+μ+12)h_n \left( \mu \right) = \left\langle {H_m^{\left( \mu \right)} |H_n^{\left( \mu \right)} } \right\rangle = 2^{2\,n} \,\left[ {\frac{n} {2}} \right]!\;\Gamma \left( {\left[ {\frac{{n + 1}} {2}} \right] + \mu + \frac{1} {2}} \right)

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the
### generalized Hermite polynomials of orders 0 to 10
### polynomial parameter is 1
###
h <- ghermite.h.inner.products( 10, 1 )
print( h )

Create list of generalized Hermite polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk generalized Hermite polynomials, Hk(μ)(x)H_k^{\left( \mu \right)} \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

ghermite.h.polynomials(n, mu, normalized = FALSE)

Arguments

n

integer value for the highest polynomial order

mu

numeric value for the polynomial parameter

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Details

The parameter μ\mu must be greater than -0.5. The function ghermite.h.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 generalized Hermite polynomial

2

order 1 generalized Hermite polynomial

...

n+1

order nn generalized Hermite polynomial

Author(s)

Frederick Novomestky [email protected]

References

Alvarez-Nordase, R., M. K. Atakishiyeva and N. M. Atakishiyeva, 2004. A q-extension of the generalized Hermite polynomials with continuous orthogonality property on R, International Journal of Pure and Applied Mathematics, 10(3), 335-347.

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

ghermite.h.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized generalized Hermite polynomials of orders 0 to 10
### polynomial parameter is 1.0
###
normalized.p.list <- ghermite.h.polynomials( 10, 1, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized generalized Hermite polynomials of orders 0 to 10
### polynomial parameter is 1.0
###
unnormalized.p.list <- ghermite.h.polynomials( 10, 1, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for generalized Hermite polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk generalized Hermite polynomial, Hk(μ)(x)H_k^{\left( \mu \right)} \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

ghermite.h.recurrences(n, mu, normalized = FALSE)

Arguments

n

integer value for the highest polynomial order

mu

numeric value for the polynomial parameter

normalized

normalized boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Details

The parameter μ\mu must be greater than -0.5.

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Alvarez-Nordase, R., M. K. Atakishiyeva and N. M. Atakishiyeva, 2004. A q-extension of the generalized Hermite polynomials with continuous orthogonality property on R, International Journal of Pure and Applied Mathematics, 10(3), 335-347.

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

ghermite.h.inner.products

Examples

###
### generate the recurrences data frame for 
### the normalized generalized Hermite polynomials
### of orders 0 to 10.
### polynomial parameter value is 1.0
###
normalized.r <- ghermite.h.recurrences( 10, 1, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the unnormalized generalized Hermite polynomials
### of orders 0 to 10.
### polynomial parameter value is 1.0
###
unnormalized.r <- ghermite.h.recurrences( 10, 1, normalized=FALSE )
print( unnormalized.r )

Weight function for the generalized Hermite polynomial

Description

This function returns the value of the weight function for the order kk generalized Hermite polynomial, Hk(μ)(x)H_k^{\left( \mu \right)} \left( x \right).

Usage

ghermite.h.weight(x, mu)

Arguments

x

a numeric vector function argument

mu

polynomial parameter

Details

The function takes on non-zero values in the interval (,)\left( -\infty,\infty \right). The parameter μ\mu must be greater than -0.5. The formula used to compute the generalized Hermite weight function is as follows.

w(x,μ)=x2  μ  exp(x2)w\left( {x,\mu } \right) = \left| x \right|^{2\;\mu } \;\exp \left( { - x^2 } \right)

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the generalized Hermite weight function for argument values 
### between -3 and 3
###
x <- seq( -3, 3, .01 )
y <- ghermite.h.weight( x, 1 )

Inner products of generalized Laguerre polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk generalized Laguerre polynomial, Ln(α)(x)L_n^{\left( \alpha \right)} \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

glaguerre.inner.products(n,alpha)

Arguments

n

integer highest polynomial order

alpha

polynomial parameter

Details

The formula used to compute the inner products is as follows.

hn=Ln(α)Ln(α)=Γ(α+n+1)n!h_n = \left\langle {L_n^{\left( \alpha \right)} |L_n^{\left( \alpha \right)} } \right\rangle = \frac{{\Gamma \left( {\alpha + n + 1} \right)}} {{n!}}.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the
### generalized Laguerre polynomial inner products of orders 0 to 10
### polynomial parameter is 1.
###
h <- glaguerre.inner.products( 10, 1 )
print( h )

Create list of generalized Laguerre polynomials

Description

This function returns a list with n+1n + 1 elements containing the order nn generalized Laguerre polynomials, Ln(α)(x)L_n^{\left( \alpha \right)} \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

glaguerre.polynomials(n, alpha, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for the polynomial parameter

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function glaguerre.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 generalized Laguerre polynomial

2

order 1 generalized Laguerre polynomial

...

n+1

order nn generalized Laguerre polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

glaguerre.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized generalized Laguerre polynomials of orders 0 to 10
### polynomial parameter is 1.0
###
normalized.p.list <- glaguerre.polynomials( 10, 1, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized generalized Laguerre polynomials of orders 0 to 10
### polynomial parameter is 1.0
###
unnormalized.p.list <- glaguerre.polynomials( 10, 1, normalized=FALSE )

Recurrence relations for generalized Laguerre polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk generalized Laguerre polynomial, Ln(α)(x)L_n^{\left( \alpha \right)} \left( x \right) and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

glaguerre.recurrences(n, alpha, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for the polynomial parameter

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

glaguerre.inner.products

Examples

###
### generate the recurrences data frame for 
### the normalized generalized Laguerre polynomials
### of orders 0 to 10.  the polynomial parameter value is 1.0.
###
normalized.r <- glaguerre.recurrences( 10, 1, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the unnormalized generalized Laguerre polynomials
### of orders 0 to 10.  the polynomial parameter value is 1.0.
###
unnormalized.r <- glaguerre.recurrences( 10, 1, normalized=FALSE )
print( unnormalized.r )

Weight function for the generalized Laguerre polynomial

Description

This function returns the value of the weight function for the order kk generalized Laguerre polynomial, Ln(α)(x)L_n^{\left( \alpha \right)} \left( x \right).

Usage

glaguerre.weight(x,alpha)

Arguments

x

the function argument which can be a vector

alpha

polynomial parameter

Details

The function takes on non-zero values in the interval (0,)\left( 0,\infty \right). The formula used to compute the weight function is as follows.

w(x)=exxαw\left( x \right) = e^{ - x} \,x^\alpha

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the generalized Laguerre weight function for argument values
### between -3 and 3
### polynomial parameter value is 1.0
###
x <- seq( -3, 3, .01 )
y <- glaguerre.weight( x, 1 )

Inner products of Hermite polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Hermite polynomial, Hk(x)H_k \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

hermite.h.inner.products(n)

Arguments

n

integer value for highest polynomial order

Details

The formula used to compute the innner product is as follows.

hn=HnHn=π  2n  n!h_n = \left\langle {H_n |H_n } \right\rangle = \sqrt \pi \;2^n \;n!.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the
### Hermite polynomials of orders 0 to 10
###
h <- hermite.h.inner.products( 10 )
print( h )

Create list of Hermite polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Hermite polynomials, Hk(x)H_k \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

hermite.h.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function hermite.h.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Hermite polynomial

2

order 1 Hermite polynomial

...

n+1

order nn Hermite polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

hermite.h.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized Hermite polynomials of orders 0 to 10
###
normalized.p.list <- hermite.h.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized Hermite polynomials of orders 0 to 10
###
unnormalized.p.list <- hermite.h.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Hermite polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Hermite polynomial, Hk(x)H_k \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

hermite.h.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

hermite.h.inner.products,

Examples

###
### generate the recurrences data frame for 
### the normalized Hermite H polynomials
### of orders 0 to 10.
###
normalized.r <- hermite.h.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the unnormalized Hermite H polynomials
### of orders 0 to 10.
###
unnormalized.r <- hermite.h.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the Hermite polynomial

Description

This function returns the value of the weight function for the order kk Hermite polynomial, Hk(x)H_k \left( x \right).

Usage

hermite.h.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (,)\left( -\infty,\infty \right). The formula used to compute the weight function.

w(x)=exp(x2)w\left( x \right) = \exp \left( { - x^2 } \right)

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the Hermite weight function for argument values
### between -3 and 3
x <- seq( -3, 3, .01 )
y <- hermite.h.weight( x )
plot( x, y )

Inner products of Hermite polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Hermite polynomial, Hek(x)He_k \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

hermite.he.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=HenHen=2π  n!h_n = \left\langle {He_n |He_n } \right\rangle = \sqrt {2\,\pi } \;n!.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the
### scaled Hermite polynomials of orders 0 to 10
###
h <- hermite.he.inner.products( 10 )
print( h )

Create list of Hermite polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Hermite polynomials, Hek(x)He_k \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

hermite.he.polynomials(n, normalized=FALSE)

Arguments

n

integer value for thehighest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function hermite.he.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Hermite polynomial

2

order 1 Hermite polynomial

...

n+1

order nn Hermite polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

hermite.he.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized Hermite polynomials of orders 0 to 10
###
normalized.p.list <- hermite.he.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized Hermite polynomials of orders 0 to 10
###
unnormalized.p.list <- hermite.he.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Hermite polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Hermite polynomial, Hek(x)He_k \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

hermite.he.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

hermite.he.inner.products

Examples

###
### generate the recurrences data frame for 
### the normalized Hermite H polynomials
### of orders 0 to 10.
###
normalized.r <- hermite.he.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the unnormalized Hermite H polynomials
### of orders 0 to 10.
###
unnormalized.r <- hermite.he.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the Hermite polynomial

Description

This function returns the value of the weight function for the order kk Hermite polynomial, Hek(x)He_k \left( x \right).

Usage

hermite.he.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (,)\left( -\infty,\infty \right). The formula used to compute the weight function is as follows.

w(x)=exp(x22)w\left( x \right) = \exp \left( { - \frac{{x^2 }}{2}} \right)

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the scaled Hermite weight function for argument values
### between -3 and 3
###
x <- seq( -3, 3, .01 )
y <- hermite.he.weight( x )

Inner products of Jacobi polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Jacobi polynomial, Gk(p,q,x)G_k \left( {p,q,x} \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

jacobi.g.inner.products(n,p,q)

Arguments

n

integer value for the highest polynomial order

p

numeric value for the first polynomial parameter

q

numeric value for the first polynomial parameter

Details

The formula used to compute the inner products is as follows.

hn=GnGn=n!  Γ(n+q)Γ(n+p)Γ(n+pq+1)(2n+p)[Γ(2n+p)]2h_n = \left\langle {G_n |G_n } \right\rangle = \frac{{n!\;\Gamma \left( {n + q} \right)\,\Gamma \left( {n + p} \right)\,\Gamma \left( {n + p - q + 1} \right)}} {{\left( {2\,n + p} \right)\,\left[ {\Gamma \left( {2\,n + p} \right)} \right]^2 }}.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the 
### G Jacobi polynomials of orders 0 to 10
### parameter p is 3 and parameter q is 2
###
h <- jacobi.g.inner.products( 10, 3, 2 )
print( h )

Create list of Jacobi polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Jacobi polynomials, Gk(p,q,x)G_k \left( {p,q,x} \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

jacobi.g.polynomials(n, p, q, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

p

numeic value for the first polynomial parameter

q

numeric value for the second polynomial parameter

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function jacobi.g.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Jacobi polynomial

2

order 1 Jacobi polynomial

...

n+1

order nn Jacobi polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

jacobi.g.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized Jacobi G polynomials of orders 0 to 10
### first parameter value p is 3 and second parameter value q is 2
###
normalized.p.list <- jacobi.g.polynomials( 10, 3, 2, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of normalized Jacobi G polynomials of orders 0 to 10
### first parameter value p is 3 and second parameter value q is 2
###
unnormalized.p.list <- jacobi.g.polynomials( 10, 3, 2, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Jacobi polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Jacobi polynomial, Gk(p,q,x)G_k \left( {p,q,x} \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

jacobi.g.recurrences(n, p, q, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

p

numeric value for the first polynomial parameter

q

numeric value for the second polynomial parameter

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

jacobi.g.inner.products, pochhammer

Examples

###
### generate the recurrences data frame for 
### the normalized Jacobi G polynomials
### of orders 0 to 10.
### parameter p is 3 and parameter q is 2
###
normalized.r <- jacobi.g.recurrences( 10, 3, 2, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the normalized Jacobi G polynomials
### of orders 0 to 10.
### parameter p is 3 and parameter q is 2
###
unnormalized.r <- jacobi.g.recurrences( 10, 3, 2, normalized=FALSE )
print( unnormalized.r )

Weight function for the Jacobi polynomial

Description

This function returns the value of the weight function for the order kk Jacobi polynomial, Gk(p,q,x)G_k \left( {p,q,x} \right).

Usage

jacobi.g.weight(x,p,q)

Arguments

x

the function argument which can be a vector

p

the first polynomial parameter

q

the second polynomial parameter

Details

The function takes on non-zero values in the interval (0,1)\left( 0,1 \right). The formula used to compute the weight function is as follows.

w(x)=(1x)pq  xq1w\left( x \right) = \left( {1 - x} \right)^{p - q} \;x^{q - 1}

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the Jacobi G weight function for argument values
### between 0 and 1
### parameter p is 3 and q is 2
###
x <- seq( 0, 1, .01 )
y <- jacobi.g.weight( x, 3, 2 )

Create list of Jacobi matrices from monic recurrence parameters

Description

Return a list of $n$ real symmetric, tri-diagonal matrices which are the principal minors of the n×nn \times n Jacobi matrix derived from the monic recurrence parameters, aa and bb, for orthogonal polynomials.

Usage

jacobi.matrices(r)

Arguments

r

a data frame containing the parameters aa and bb

Value

A list of symmetric, tri-diagnonal matrices

1

a 1×11 \times 1 matrix

2

a 2×22 \times 2 matrix

...

n

an n×nn \times n matrix

Author(s)

Frederick Novomestky [email protected]

References

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Examples

r <- chebyshev.t.recurrences( 5 )
m.r <- monic.polynomial.recurrences( r )
j.m <- jacobi.matrices( m.r )

Inner products of Jacobi polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Jacobi polynomial, Pk(α,β)(x)P_k^{\left( {\alpha ,\beta } \right)} \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

jacobi.p.inner.products(n,alpha,beta)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for the first polynomial parameter

beta

numeric value for the first polynomial parameter

Details

The formula used to compute the innser products is as follows.

hn=Pn(α,β)Pn(α,β)=2α+β+12n+α+β+1Γ(n+α+1)Γ(n+β+1)n!  Γ(n+α+β+1)h_n = \left\langle {P_n^{\left( {\alpha ,\beta } \right)} |P_n^{\left( {\alpha ,\beta } \right)} } \right\rangle = \frac{{2^{\alpha + \beta + 1} }} {{2\,n + \alpha + \beta + 1}}\frac{{\Gamma \left( {n + \alpha + 1} \right)\,\Gamma \left( {n + \beta + 1} \right)}} {{n!\;\Gamma \left( {n + \alpha + \beta + 1} \right)}}.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner product vector for the P Jacobi polynomials of orders 0 to 10
###
h <- jacobi.p.inner.products( 10, 2, 2 )
print( h )

Create list of Jacobi polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Jacobi polynomials, Pk(α,β)(x)P_k^{\left( {\alpha ,\beta } \right)} \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

jacobi.p.polynomials(n, alpha, beta, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for the first polynomial parameter

beta

numeric value for the second polynomial parameter

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function jacobi.p.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Jacobi polynomial

2

order 1 Jacobi polynomial

...

n+1

order nn Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

jacobi.p.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized Jacobi P polynomials of orders 0 to 10
### first parameter value a is 2 and second parameter value b is 2
###
normalized.p.list <- jacobi.p.polynomials( 10, 2, 2, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized Jacobi P polynomials of orders 0 to 10
### first parameter value a is 2 and second parameter value b is 2
###
unnormalized.p.list <- jacobi.p.polynomials( 10, 2, 2, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Jacobi polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Jacobi polynomial, Pk(α,β)(x)P_k^{\left( {\alpha ,\beta } \right)} \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

jacobi.p.recurrences(n, alpha, beta, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for the first polynomial parameter

beta

numeric value for the second polynomial parameter

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

jacobi.p.inner.products, pochhammer

Examples

###
### generate the recurrences data frame for 
### the normalized Jacobi P polynomials
### of orders 0 to 10.
### parameter a is 2 and parameter b is 2
###
normalized.r <- jacobi.p.recurrences( 10, 2, 2, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the unnormalized Jacobi P polynomials
### of orders 0 to 10.
### parameter a is 2 and parameter b is 2
###
unnormalized.r <- jacobi.p.recurrences( 10, 2, 2, normalized=FALSE )
print( unnormalized.r )

Weight function for the Jacobi polynomial

Description

This function returns the value of the weight function for the order kk Jacobi polynomial, Pk(α,β)(x)P_k^{\left( {\alpha ,\beta } \right)} \left( x \right).

Usage

jacobi.p.weight(x,alpha,beta)

Arguments

x

the function argument which can be a vector

alpha

the first polynomial parameter

beta

the second polynomial parameter

Details

The function takes on non-zero values in the interval (1,1)\left( -1,1 \right). The formula used to compute the weight function is as follows.

w(x)=(1x)α  (1+x)βw\left( x \right) = \left( {1 - x} \right)^\alpha \;\left( {1 + x} \right)^\beta

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the Jacobi P weight function for argument values
### between -1 and 1
###
x <- seq( -1, 1, .01 )
y <- jacobi.p.weight( x, 2, 2 )

Inner products of Laguerre polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Laguerre polynomial, Ln(x)L_n \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

laguerre.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=LnLn=1h_n = \left\langle {L_n |L_n } \right\rangle = 1.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the
### Laguerre polynomial inner products of orders 0 to 10
###
h <- laguerre.inner.products( 10 )
print( h )

Create list of Laguerre polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Laguerre polynomials, Ln(x)L_n \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

laguerre.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function laguerre.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Laguerre polynomial

2

order 1 Laguerre polynomial

...

n+1

order nn Laguerre polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

laguerre.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized Laguerre polynomials of orders 0 to 10
###
normalized.p.list <- laguerre.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized Laguerre polynomials of orders 0 to 10
###
unnormalized.p.list <- laguerre.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Laguerre polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Laguerre polynomial, Ln(x)L_n \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

laguerre.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

glaguerre.recurrences

Examples

###
### generate the recurrences data frame for 
### the normalized Laguerre polynomials
### of orders 0 to 10.
###
normalized.r <- laguerre.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the normalized Laguerre polynomials
### of orders 0 to 10.
###
unnormalized.r <- laguerre.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the Laguerre polynomial

Description

This function returns the value of the weight function for the order kk Laguerre polynomial, Ln(x)L_n \left( x \right).

Usage

laguerre.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (0,)\left( 0,\infty \right). The formula used to compute the weight function is as follows. w(x)=exw\left( x \right) = e^{ - x}

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the Laguerre weight function for argument values
### between 0 and 3
x <- seq( -0, 3, .01 )
y <- laguerre.weight( x )
plot( x, y )

Inner products of Legendre polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk Legendre polynomial, Pk(x)P_k \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

legendre.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used compute the inner products is as follows.

hn=PnPn=22n+1h_n = \left\langle {P_n |P_n } \right\rangle = \frac{2} {{2\,n + 1}}.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

spherical.inner.products

Examples

###
### compute the inner product for the
###  Legendre polynomials of orders 0 to 1
###
h <- legendre.inner.products( 10 )
print( h )

Create list of Legendre polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk Legendre polynomials, Pk(x)P_k \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

legendre.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function legendre.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 Legendre polynomial

2

order 1 Legendre polynomial

...

n+1

order nn Legendre polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

legendre.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized Laguerre polynomials of orders 0 to 10
###
normalized.p.list <- legendre.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized Laguerre polynomials of orders 0 to 10
###
unnormalized.p.list <- legendre.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for Legendre polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk Legendre polynomial, Pk(x)P_k \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

legendre.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

legendre.inner.products

Examples

###
### generate the recurrences data frame for 
### the normalized Legendre polynomials
### of orders 0 to 10.
###
normalized.r <- legendre.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrences data frame for 
### the normalized Legendre polynomials
### of orders 0 to 10.
###
unnormalized.r <- legendre.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the Legendre polynomial

Description

This function returns the value of the weight function for the order kk Legendre polynomial, Pk(x)P_k \left( x \right).

Usage

legendre.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (1,1)\left( -1,1 \right). The formula used to compute the weight function is as follows.

w(x)=1w\left( x \right) = 1

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the Legendre weight function for argument values
### between -1 and 1
###
x <- seq( -1, 1, .01 )
y <- legendre.weight( x )
plot( x, y )

Calculate the logarithm of Pochhammer's symbol

Description

lpochhammer returns the value of the natural logarithm of Pochhammer's symbol calculated as

ln[(z)n]=lnΓ(z+n)lnΓ(z)\ln \left[ {\left( z \right)_n } \right] = \ln \Gamma \left( {z + n} \right) - \ln \Gamma \left( z \right)

where Γ(z)\Gamma \left( z \right) is the Gamma function

Usage

lpochhammer(z, n)

Arguments

z

argument of the symbol

n

integer number of terms in the symbol

Value

The value of the logarithm of the symbol

Author(s)

Frederick Novomestky [email protected]

See Also

pochhammer

Examples

lpochhammer( pi, 5 )

Create data frame of monic recurrences

Description

This function returns a data frame with parameters required to construct monic orthogonal polynomials based on the standard recurrence relation for the non-monic polynomials. The recurrence relation for monic orthogonal polynomials is as follows.

qk+1(x)=(xak)  qk(x)bk  qk1(x)q_{k + 1} \left( x \right) = \left( {x - a_k } \right)\;q_k \left( x \right) - b_k \;q_{k - 1} \left( x \right)

We require that q1(x)=0q_{-1} \left( x \right) = 0 and q0(x)=1q_0 \left( x \right) = 1. The recurrence for non-monic orthogonal polynomials is given by

ck  pk+1(x)=(dk+ek  x)  pk(x)fk  pk1(x)c_k \;p_{k + 1} \left( x \right) = \left( {d_k + e_k \;x} \right)\;p_k \left( x \right) - f_k \;p_{k - 1} \left( x \right)

We require that p1(x)=0p_{-1} \left( x \right) = 0 and p0(x)=1p_0 \left( x \right) = 1. The monic polynomial recurrence parameters, a and b, are related to the non-monic polynomial parameter vectors c, d, e and f in the following manner.

ak=dkeka_k = - \frac{{d_k }}{{e_k }}

bk=ck1  fkek1  ekb_k = \frac{{c_{k - 1} \;f_k }}{{e_{k - 1} \;e_k }}

with b0=0b_0 = 0.

Usage

monic.polynomial.recurrences(recurrences)

Arguments

recurrences

the data frame of recurrence parameter vectors c, d, e and f

Value

A data frame with n+1n + 1 rows and two named columns, a and b.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

orthogonal.polynomials,

Examples

###
### construct a list of the recurrences for the T Chebyshev polynomials of
### orders 0 to 10
###
r <- chebyshev.t.recurrences( 10, normalized=TRUE )
###
### construct the monic polynomial recurrences from the above list
###
m.r <- monic.polynomial.recurrences( r )

Create list of monic orthogonal polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk monic polynomials for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

monic.polynomials(monic.recurrences)

Arguments

monic.recurrences

a data frame containing the parameters a and b

Value

A list with n+1n + 1 polynomial objects

1

order 0 monic orthogonal polynomial

2

order 1 monic orthogonal polynomial

...

n+1

order nn monic orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

See Also

monic.polynomial.recurrences

Examples

###
### generate the recurrences for the T Chebyshev polynomials
### of orders 0 to 10
###
r <- chebyshev.t.recurrences( 10, normalized=TRUE )
###
### get the corresponding monic polynomial recurrences
###
m.r <- monic.polynomial.recurrences( r )
###
### obtain the list of monic polynomials
###
p.list <- monic.polynomials( m.r )

Create orthogonal polynomials

Description

Create list of orthogonal polynomials from the following recurrence relations for k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

ckpk+1(x)=(dk+ekx)pk(x)fkpk1(x)c_k p_{k+1}\left( x \right) = \left( d_k + e_k x \right) p_k \left( x \right) - f_k p_{k-1} \left( x \right)

We require that p1(x)=0p_{-1} \left( x \right) = 0 and p0(x)=1p_0 \left( x \right) = 1. The coefficients are the column vectors c{\bf{c}}, d{\bf{d}}, e{\bf{e}} and f{\bf{f}}.

Usage

orthogonal.polynomials(recurrences)

Arguments

recurrences

a data frame containing the parameters of the orthogonal polynomial recurrence relations

Details

The argument is a data frame with n+1n + 1 rows and four named columns. The column names are c, d, e and f. These columns correspond to the column vectors described above.

Value

A list of n+1n + 1 polynomial objects

1

Order 0 orthogonal polynomial

2

Order 1 orthogonal polynomial

...

n+1

Order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the recurrence relations for T Chebyshev polynomials of orders 0 to 10
###
r <- chebyshev.t.recurrences( 10, normalized=FALSE )
print( r )
###
### generate the list of orthogonal polynomials
###
p.list <- orthogonal.polynomials( r )
print( p.list )

Create orthonormal polynomials

Description

Create list of orthonormal polynomials from the following recurrence relations for k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

ckpk+1(x)=(dk+ekx)pk(x)fkpk1(x)c_k p_{k+1}\left( x \right) = \left( d_k + e_k x \right) p_k \left( x \right) - f_k p_{k-1} \left( x \right)

We require that p1(x)=0p_{-1} \left( x \right) = 0 and p0(x)=1p_0 \left( x \right) = 1. The coefficients are the column vectors c{\bf{c}}, d{\bf{d}}, e{\bf{e}} and f{\bf{f}}.

Usage

orthonormal.polynomials(recurrences, p.0)

Arguments

recurrences

a data frame containing the parameters of the orthonormal polynomial recurrence relations

p.0

a polynomial object for the order 0 orthonormal polynomial

Details

The argument is a data frame with n+1n + 1 rows and four named columns. The column names are c, d, e and f. These columns correspond to the column vectors described above.

Value

A list of n+1n + 1polynomial objects

1

Order 0 orthonormal polynomial

2

Order 1 orthonormal polynomial

...

n+1

Order nn orthonormal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate a data frame with the recurrences parameters for normalized T Chebyshev
### polynomials of orders 0 to 10
###
r <- chebyshev.t.recurrences( 10, normalized=TRUE )
print( r )
norm <- sqrt( pi )
###
### create the order 0 orthonormal polynomial
###
library("polynom")
p.0 <- polynomial( c( 1 / norm ) )
###
### generate a list of orthonormal polynomial objects
###
p.list <- orthonormal.polynomials( r, p.0 )
print( p.list )

Calculate the value of Pochhammer's symbol

Description

pochhammer returns the value of Pochhammer's symbol calculated as

(z)n=z  (z+1)    (z+n1)=Γ(z+n)Γ(z)\left( z \right)_n = z\;\left( {z + 1} \right)\; \ldots \;\left( {z + n - 1} \right) = \frac{{\Gamma \left( {z + n} \right)}}{{\Gamma \left( z \right)}}

where Γ(z)\Gamma \left( z \right) is the Gamma function

Usage

pochhammer(z, n)

Arguments

z

numeric value for the argument of the symbol

n

integer value for the number of terms in the symbol

Value

The value of Pochhammer's symbol

Author(s)

Frederick Novomestky [email protected]

Examples

###
### compute the Pochhamer's symbol fo z equal to 1 and
### n equal to 5
###
pochhammer( 1, 5 )

Create list of polynomial coefficient vectors

Description

This function returns a list with n+1n + 1 elements containing the vector of coefficients of the order kk polynomials for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n. Each element in the list is a vector.

Usage

polynomial.coefficients(polynomials)

Arguments

polynomials

list of polynomial objects

Value

A list of n+1n + 1 polynomial objects where each element is a vector of coefficients.

1

Coefficient(s) of the order 0 polynomial

2

Coefficient(s) of the order 1 polynomial

...

n+1

Coefficient(s) of the order nn polynomial

Author(s)

Frederick Novomestky [email protected]

Examples

###
### generate a list of normalized T Chebyshev polynomials
### of orders 0 to 10
###
p.list <- chebyshev.t.polynomials( 10, normalized=TRUE )
###
### obtain the list of coefficients for these polynomials
###
p.coef <- polynomial.coefficients( p.list )

Create list of polynomial derivatives

Description

This function returns a list with n+1n + 1 elements containing polynomial objects which are the derivatives of the order kk polynomials for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

polynomial.derivatives(polynomials)

Arguments

polynomials

list of polynomial objects

Details

The polynomial objects in the argument polynomials are as follows

  • 1order 0 polynomial

  • 2order 1 polynomial ...

  • n+1order nn polynomial

Value

List of n+1n + 1 polynomial objects

1

derivative of polynomials[[1]]

2

derivative of polynomials[[2]]

...

n+1

derivative of polynomials[[n+1]]

Author(s)

Frederick Novomestky [email protected]

Examples

###
### generate a list of normalized T Chebyshev polynomials of
### orders 0 to 10
###
p.list <- chebyshev.t.polynomials( 10, normalized=TRUE )
###
### generate the corresponding list of polynomial derivatives
###
p.deriv <- polynomial.derivatives( p.list )

Coerce polynomials to functions

Description

This function returns a list with n+1n + 1 elements containing the functions of the order $k$ polynomials for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n and for the given argument xx.

Usage

polynomial.functions(polynomials, ...)

Arguments

polynomials

a list of polynomial objects

...

further arguments to be passed to or from methods

Details

The function uses the method as.function.polynomial to coerce each polynomial object to a function object.

Value

A list of n+1n + 1 polynomial objects where each element is the function for the polynomial.

1

Function for the order 0 polynomial

2

Function for the order 1 polynomial

...

n+1

Function for the order nn polynomial

Author(s)

Frederick Novomestky [email protected]

Examples

###
### generate a list of T Chebyshev polynomials of
### orders 0 to 10
###
p.list <- chebyshev.t.polynomials( 10, normalized=FALSE )
###
### create the list of functions for each polynomial
###
f.list <- polynomial.functions( p.list )

Create list of polynomial integrals

Description

This function returns a list with n+1n + 1 elements containing polynomial objects which are the indefinite integrals of the order kk polynomials for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

polynomial.integrals(polynomials)

Arguments

polynomials

list of polynomial objects

Details

The polynomial objects in the argument polynomials are as follows

  • 1order 0 polynomial

  • 2order 1 polynomial ...

  • n+1order n polynomial

Value

List of n+1n + 1 polynomial objects

1

integral of polynomials[[1]]

2

integral of polynomials[[2]]

...

n+1

integral of polynomials[[n+1]]

Author(s)

Frederick Novomestky [email protected]

Examples

###
### generate a list of normalized T Chebyshev polynomials
### of orders 0 to 10
###
p.list <- chebyshev.t.polynomials( 10, normalized=TRUE )
###
### generate the corresponding list of polynomial integrals
###
p.int <- polynomial.integrals( p.list )

Create vector of polynomial orders

Description

This function returns a vector with nn elements containing the orders of the polynomials

Usage

polynomial.orders(polynomials)

Arguments

polynomials

list of $n$ polynomial objects

Value

A vector of nn values

1

Order of polynomials[[1]]

2

Order of polynomials[[2]]

...

n

Order of polynomials[[n]]

Author(s)

Frederick Novomestky [email protected]

Examples

###
### generate a list of normalized T Chebyshev polynomials
### of orders 0 to 10
###
p.list <- chebyshev.t.polynomials( 10, normalized=TRUE )
###
### get the vector of polynomial orders
###
p.order <- polynomial.orders( p.list )

Create a list of polynomial linear combinations

Description

This function returns a list with n+1n + 1 elements containing the vector of linear combinations of the order kk polynomials for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n. Each element in the list is a vector.

Usage

polynomial.powers(polynomials)

Arguments

polynomials

A list of polynomials

Details

The jj-th component in the list is a vector of linear combinations of the order kk polynomials for orders k=0,  1,  ,  j1k = 0,\;1,\; \ldots ,\;j - 1 equal to the monomial x raised to the power j1j - 1.

Value

A list of n+1n + 1 elements where each element is a vector of linear combinations.

1

Linear combination(s) of polynomials up to order 0

2

Linear combination(s) of polynomials up to order 1

...

n+1

Linear combination(s) of polynomials up to order nn

Author(s)

Frederick Novomestky [email protected]

Examples

###
### generate Legendre polynomials of orders 0 to 10
###
polynomials <- legendre.polynomials( 10 )
###
### generate list of linear combinations of these polynomials
###
alphas <- polynomial.powers( polynomials )
print( alphas )

Create a list of polynomial roots

Description

This function returns a list with nn elements containing the roots of the order $k$ monic orthogonal polynomials for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n using a data frame with the monic polynomial recurrence parameter vectors a\bf{a} and b\bf{b}

Usage

polynomial.roots(m.r)

Arguments

m.r

monic recurrence data frame with parameters a and b

Details

The parameter m.r is a data frame with $n$+1 rows and two names columns. The columns which are names a and b correspond to the above referenced vectors. Function jacobi.matrices is used to create a list of symmetric, tridiagonal Jacobi matrices from these named columns. The eigenvalues of the k×kk \times k Jacobi matrix are the roots or zeros of the order $k$ monic orthogonal polynomial.

Value

A list with nn elements each of which is a vector of polynomial roots

1

roots of the order 1 monic polynomial

2

roots of the order 2 monic polynomial

...

n

roots of the order nn monic polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

monic.polynomial.recurrences, jacobi.matrices

Examples

###
### generate the recurrences data frame for
### the normalized Chebyshev polynomials of
### orders 0 to 10
###
r <- chebyshev.t.recurrences( 10, normalized=TRUE )
###
### construct the corresponding monic polynomial
### recurrences
###
m.r <- monic.polynomial.recurrences( r )
###
### obtain the polynomial roots from the monic polynomial
### recurrences
p.roots <- polynomial.roots( m.r )

Create vector of polynomial values

Description

This function returns a list with n+1n + 1 elements containing the values of the order kk polynomials for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n and for the given argument xx.

Usage

polynomial.values( polynomials, x )

Arguments

polynomials

list of polynomial objects

x

the argument which can be any numeric object

Value

A list of n+1n + 1 polynomial objects where each element is the value of the polynomial.

1

Value(s) for the order 0 polynomial

2

Value(s) for the order 1 polynomial

...

n+1

Value(s) for the order nn polynomial

Author(s)

Frederick Novomestky [email protected]

Examples

###
### generate a list of T Chebyshev polynomials of
### orders 0 to 10
###
p.list <- chebyshev.t.polynomials( 10, normalized=FALSE )
x <- seq( -2, 2, .01 )
###
### compute the value of the polynomials for the given range of values in x
###
y <- polynomial.values( p.list, x )
print( y )

Scale values from [a,b] to [u.v]

Description

This function returns a vector of values that have been mapped from the interval [a,b] to the interval [u.v].

Usage

scaleX(x, a = min(x, na.rm = TRUE), b = max(x, na.rm = TRUE), u, v)

Arguments

x

A numerical vector of values to be mapped into a target interval

a

A numerical lower bound for the domain interval with min(x) as the default value

b

A numerical upper bound for the domain interval with max(x) as the default value

u

A numerical lower bound for the target interval

v

A numerical upper bound for the target interval

Details

Target lower and/or upper bounds can be -\infty and \infty, respectively. This accomodates finite target intervals, semi-infinite target intervals and infinite target intervals.

Value

A vector of transformed values with four attributes. The first attribute is called "a" and it is the domain interval lower bound. The second attribute is called "b" and it is the domain interval upper bound. The third attribute is called "u" and it is the target interval lower bound. The fourth attribute is called "v" and it is the target interval upper bound.

Author(s)

Frederick Novomestky [email protected], Gregor Gorjanc [email protected]

References

Seber, G. A. F. (1997) Linear Regression Analysis, New York.

Examples

x <- rnorm( 1000, 0, 10 )
y0 <- scaleX( x, u=0 , v=1 )
y1 <- scaleX( x, u=-1, v=1 )
y2 <- scaleX( x, u=-Inf, v=0 )
y3 <- scaleX( x, u=0, v=Inf )
y4 <- scaleX( x, u=-Inf, v=Inf )

Inner products of shifted Chebyshev polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk shifted Chebyshev polynomial of the first kind, Tk(x)T_k^* \left( x\right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

schebyshev.t.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=TnTn={π2n0πn=0h_n = \left\langle {T_n^* |T_n^* } \right\rangle = \left\{ {\begin{array}{cc} {\frac{\pi } {2}} & {n \ne 0} \\ \pi & {n = 0} \\ \end{array} } \right.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., NY.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### generate the inner products vector for the
### shifted T Chebyshev polynomials of orders 0 to 10
###
h <- schebyshev.t.inner.products( 10 )
print( h )

Create list of shifted Chebyshev polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk shifted Chebyshev polynomials of the first kind, Tk(x)T_k^* \left( x\right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

schebyshev.t.polynomials(n, normalized)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function schebyshev.t.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 shifted Chebyshev polynomial

2

order 1 shifted Chebyshev polynomial

...

n+1

order nn shifted Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

schebyshev.u.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized shifted T Chebyshev polynomials of orders 0 to 10
###
normalized.p.list <- schebyshev.t.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized shifted T Chebyshev polynomials of orders 0 to 10
###
unnormalized.p.list <- schebyshev.t.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for shifted Chebyshev polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk shifted Chebyshev polynomial of the first kind, Tk(x)T_k^* \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

schebyshev.t.recurrences(n, normalized)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

schebyshev.t.inner.products

Examples

###
### generate the recurrence relations for 
### the normalized shifted T Chebyshev polynomials
### of orders 0 to 10
###
normalized.r <- schebyshev.t.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrence relations for 
### the unnormalized shifted T Chebyshev polynomials
### of orders 0 to 10
###
unnormalized.r <- schebyshev.t.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the shifted Chebyshev polynomial

Description

This function returns the value of the weight function for the order kk shifted Chebyshev polynomial of the first kind, Tk(x)T_k^* \left( x \right).

Usage

schebyshev.t.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (0,1)\left( 0,1 \right). The formula used to compute the weight function is as follows.

w(x)=1xx2w\left( x \right) = \frac{1}{{\sqrt {x - x^2 } }}

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the shifted T Chebyshev weight function for argument values
### between 0 and 1
x <- seq( 0, 1, .01 )
y <- schebyshev.t.weight( x )
plot( x, y )

Inner products of shifted Chebyshev polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk shifted Chebyshev polynomial of the second kind, Uk(x)U_k^* \left( x\right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

schebyshev.u.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=UnUn=π8h_n = \left\langle {U_n^* |U_n^* } \right\rangle = \frac{\pi }{8}.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

h <- schebyshev.u.inner.products( 10 )

Create list of shifted Chebyshev polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk shifted Chebyshev polynomials of the second kind, Uk(x)U_k^* \left( x\right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

schebyshev.u.polynomials(n, normalized)

Arguments

n

integer value for highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function schebyshev.u.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 shifted Chebyshev polynomial

2

order 1 shifted Chebyshev polynomial

...

n+1

order nn shifted Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

schebyshev.u.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized shifted U Chebyshev polynomials of orders 0 to 10
###
normalized.p.list <- schebyshev.u.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized shifted U Chebyshev polynomials of orders 0 to 10
###
unnormalized.p.list <- schebyshev.u.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for shifted Chebyshev polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk shifted Chebyshev polynomial of the second kind, Uk(x)U_k^* \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

schebyshev.u.recurrences(n, normalized)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

schebyshev.u.inner.products

Examples

###
### generate the recurrence relations for 
### the normalized shifted U Chebyshev polynomials
### of orders 0 to 10
###
normalized.r <- schebyshev.u.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrence relations for 
### the unnormalized shifted T Chebyshev polynomials
### of orders 0 to 10
unnormalized.r <- schebyshev.u.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the shifted Chebyshev polynomial

Description

This function returns the value of the weight function for the order kk shifted Chebyshev polynomial of the second kind, Uk(x)U_k^* \left( x \right).

Usage

schebyshev.u.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (0,1)\left( 0,1 \right). The formula used to compute the weight function is as follows.

w(x)=xx2w\left( x \right) = \sqrt {x - x^2 }

Value

The value of the weight function.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the shifted U Chebyshev weight function for argument values
### between 0 and 1
###
x <- seq( 0, 1, .01 )
y <- schebyshev.u.weight( x )
plot( x, y )

Inner products of shifted Legendre polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk shifted Legendre polynomial, Pk(x)P_k^* \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

slegendre.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products is as follows.

hn=PnPn=12n+1h_n = \left\langle {P_n^* |P_n^* } \right\rangle = \frac{1}{{2\,n + 1}}.

Value

A vector with $n$+1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the inner products vector for the
### shifted Legendre polynomials of orders 0 to 10
###
h <- slegendre.inner.products( 10 )
print( h )

Create list of shifted Legendre polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk shifted Legendre polynomials, Pk(x)P_k^* \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

slegendre.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function slegendre.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 shifted Legendre polynomial

2

order 1 shifted Legendre polynomial

...

n+1

order nn shifted Legendre polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

slegendre.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized shifted Legendre polynomials of orders 0 to 10
###
normalized.p.list <- slegendre.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized shifted Legendre polynomials of orders 0 to 10
###
unnormalized.p.list <- slegendre.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for shifted Legendre polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk shifted Legendre polynomial, Pk(x)P_k^* \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

slegendre.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

slegendre.inner.products,

Examples

###
### generate the recurrence relations for normalized shifted Legendre polynomials
### of orders 0 to 10
###
normalized.r <- slegendre.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrence relations for normalized shifted Legendre polynomials
### of orders 0 to 10
###
unnormalized.r <- slegendre.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the shifted Legendre polynomial

Description

This function returns the value of the weight function for the order kk shifted Legendre polynomial, Pk(x)P_k^* \left( x \right).

Usage

slegendre.weight(x)

Arguments

x

the function argument which can be a vector

Details

The function takes on non-zero values in the interval (0,1)\left( 0,1 \right). The formula used to compute the weight function is as follows.

w(x)=1w\left( x \right) = 1

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the shifted Legendre weight function for argument values
### between 0 and 1
###
x <- seq( 0, 1, .01 )
y <- slegendre.weight( x )

Inner products of spherical polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk spherical polynomial, Pk(x)P_k \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

spherical.inner.products(n)

Arguments

n

integer value for the highest polynomial order

Details

The formula used to compute the inner products of the spherical orthogonal polynomials is the same as that used for the Legendre orthogonal polynomials.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

legendre.inner.products

Examples

###
### generate the inner products vector for the spherical polynomals
### of orders 0 to 10.
###
h <- spherical.inner.products( 10 )
print( h )

Create list of spherical polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk spherical polynomials, Pk(x)P_k \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

spherical.polynomials(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function spherical.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 spherical polynomial

2

order 1 spherical polynomial

...

n+1

order nn Chebyshev polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

spherical.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### generate a list of spherical orthonormal polynomials of orders 0 to 10
###
normalized.p.list <- spherical.polynomials( 10, normalized=TRUE )
print( normalized.p.list )
###
### generate a list of spherical orthogonal polynomials of orders 0 to 10
###
unnormalized.p.list <- spherical.polynomials( 10, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for spherical polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk spherical polynomial, Pk(x)P_k \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

spherical.recurrences(n, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

spherical.inner.products

Examples

###
### generate the recurrence relations for 
### the normalized spherical polynomials
### of orders 0 to 10
###
normalized.r <- spherical.recurrences( 10, normalized=TRUE )
print( normalized.r )
###
### generate the recurrence relations for 
### the unnormalized spherical polynomials
### of orders 0 to 10
###
unnormalized.r <- spherical.recurrences( 10, normalized=FALSE )
print( unnormalized.r )

Weight function for the spherical polynomial

Description

This function returns the value of the weight function for the order kk spherical polynomial, Pk(x)P_k \left( x \right).

Usage

spherical.weight(x)

Arguments

x

the function argument which can be a vector or matrix

Details

The function takes on non-zero values in the interval (1,1)\left( -1,1 \right). The formula used to compute the weight function is as follows.

w(x)=1w\left( x \right) = 1

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the spherical weight function for a sequence of values between -2 and 2
###
x <- seq( -2, 2, .01 )
y <- spherical.weight( x )
plot( x, y )

Inner products of ultraspherical polynomials

Description

This function returns a vector with n+1n + 1 elements containing the inner product of an order kk ultraspherical polynomial, Ck(α)(x)C_k^{\left( \alpha \right)} \left( x \right), with itself (i.e. the norm squared) for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

ultraspherical.inner.products(n,alpha)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for the polynomial parameter

Details

This function uses the same formula as the function gegenbauer.inner.products.

Value

A vector with n+1n + 1 elements

1

inner product of order 0 orthogonal polynomial

2

inner product of order 1 orthogonal polynomial

...

n+1

inner product of order nn orthogonal polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., NY.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

gegenbauer.inner.products

Examples

###
### generate the inner products vector for the
### ultraspherical polynomials of orders 0 to 10.
### the polynomial parameter is 1.0
###
h <- ultraspherical.inner.products( 10, 1 )
print( h )

Create list of ultraspherical polynomials

Description

This function returns a list with n+1n + 1 elements containing the order kk ultraspherical polynomials, Ck(α)(x)C_k^{\left( \alpha \right)} \left( x \right), for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

ultraspherical.polynomials(n, alpha, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

alpha

polynomial parameter

normalized

a boolean value which, if TRUE, returns a list of normalized orthogonal polynomials

Details

The function ultraspherical.recurrences produces a data frame with the recurrence relation parameters for the polynomials. If the normalized argument is FALSE, the function orthogonal.polynomials is used to construct the list of orthogonal polynomial objects. Otherwise, the function orthonormal.polynomials is used to construct the list of orthonormal polynomial objects.

Value

A list of n+1n + 1 polynomial objects

1

order 0 ultraspherical polynomial

2

order 1 ultraspherical polynomial

...

n+1

order nn ultraspherical polynomial

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

gegenbauer.recurrences, orthogonal.polynomials, orthonormal.polynomials

Examples

###
### gemerate a list of normalized ultra spherical polynomials 
### of orders 0 to 10
###
normalized.p.list <- ultraspherical.polynomials( 10, 1, normalized=TRUE )
print( normalized.p.list )
###
### gemerate a list of unnormalized ultra spherical polynomials 
### of orders 0 to 10
###
unnormalized.p.list <- ultraspherical.polynomials( 10, 1, normalized=FALSE )
print( unnormalized.p.list )

Recurrence relations for ultraspherical polynomials

Description

This function returns a data frame with n+1n + 1 rows and four named columns containing the coefficient vectors c, d, e and f of the recurrence relations for the order kk ultraspherical polynomial, Ck(α)(x)C_k^{\left( \alpha \right)} \left( x \right), and for orders k=0,  1,  ,  nk = 0,\;1,\; \ldots ,\;n.

Usage

ultraspherical.recurrences(n, alpha, normalized=FALSE)

Arguments

n

integer value for the highest polynomial order

alpha

numeric value for the polynomial parameter

normalized

boolean value which, if TRUE, returns recurrence relations for normalized polynomials

Value

A data frame with the recurrence relation parameters.

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

See Also

ultraspherical.recurrences

Examples

###
### generate the recurrence relations for 
### the normalized ultraspherical polynomials
### of orders 0 to 10
### polynomial parameter value is 1.0
###
normalized.r <- ultraspherical.recurrences( 10, 1, normalized=TRUE )
print( normalized.r )
###
### generate the recurrence relations for 
### the normalized ultraspherical polynomials
### of orders 0 to 10
### polynomial parameter value is 1.0
###
unnormalized.r <- ultraspherical.recurrences( 10, 1, normalized=FALSE )
print( unnormalized.r )

Weight function for the ultraspherical polynomial

Description

This function returns the value of the weight function for the order kk ultraspherical polynomial, Ck(α)(x)C_k^{\left( \alpha \right)} \left( x \right).

Usage

ultraspherical.weight(x,alpha)

Arguments

x

the function argument which can be a vector

alpha

polynomial parameter

Details

The function takes on non-zero values in the interval (1,1)\left( -1,1 \right). The formula used to compute the weight function is as follows.

w(x)=(1x2)α0.5w\left( x \right) = \left( {1 - x^2 } \right)^{\alpha - 0.5}

Value

The value of the weight function

Author(s)

Frederick Novomestky [email protected]

References

Abramowitz, M. and I. A. Stegun, 1968. Handbook of Mathematical Functions with Formulas, Graphs, and Mathematical Tables, Dover Publications, Inc., New York.

Courant, R., and D. Hilbert, 1989. Methods of Mathematical Physics, John Wiley, New York, NY.

Press, W. H., S. A. Teukolsky, W. T. Vetterling, and B. P. Flannery, 1992. Numerical Recipes in C, Cambridge University Press, Cambridge, U.K.

Szego, G., 1939. Orthogonal Polynomials, 23, American Mathematical Society Colloquium Publications, Providence, RI.

Examples

###
### compute the ultraspherical weight function for arguments between -2 and 2
### polynomial parameter is 1.0
###
x <- seq( -2, 2, .01 )
y <- ultraspherical.weight( x, 1 )
plot( x, y )