This module provides a number of objects (mostly functions) useful for
dealing with Hermite series, including a Hermite
class that
encapsulates the usual arithmetic operations. (General information
on how this module represents and works with such polynomials is in the
docstring for its "parent" sub-package, numpy.polynomial
).
Class | Hermite |
An Hermite series class. |
Function | herm2poly |
Convert a Hermite series to a polynomial. |
Function | hermadd |
Add one Hermite series to another. |
Function | hermcompanion |
Return the scaled companion matrix of c. |
Function | hermder |
Differentiate a Hermite series. |
Function | hermdiv |
Divide one Hermite series by another. |
Function | hermfit |
Least squares fit of Hermite series to data. |
Function | hermfromroots |
Generate a Hermite series with given roots. |
Function | hermgauss |
Gauss-Hermite quadrature. |
Function | hermgrid2d |
Evaluate a 2-D Hermite series on the Cartesian product of x and y. |
Function | hermgrid3d |
Evaluate a 3-D Hermite series on the Cartesian product of x, y, and z. |
Function | hermint |
Integrate a Hermite series. |
Function | hermline |
Hermite series whose graph is a straight line. |
Function | hermmul |
Multiply one Hermite series by another. |
Function | hermmulx |
Multiply a Hermite series by x. |
Function | hermpow |
Raise a Hermite series to a power. |
Function | hermroots |
Compute the roots of a Hermite series. |
Function | hermsub |
Subtract one Hermite series from another. |
Function | hermval |
Evaluate an Hermite series at points x. |
Function | hermval2d |
Evaluate a 2-D Hermite series at points (x, y). |
Function | hermval3d |
Evaluate a 3-D Hermite series at points (x, y, z). |
Function | hermvander |
Pseudo-Vandermonde matrix of given degree. |
Function | hermvander2d |
Pseudo-Vandermonde matrix of given degrees. |
Function | hermvander3d |
Pseudo-Vandermonde matrix of given degrees. |
Function | hermweight |
Weight function of the Hermite polynomials. |
Function | poly2herm |
poly2herm(pol) |
Variable | hermdomain |
Undocumented |
Variable | hermone |
Undocumented |
Variable | hermx |
Undocumented |
Variable | hermzero |
Undocumented |
Function | _normed_hermite_n |
Evaluate a normalized Hermite polynomial. |
Convert a Hermite series to a polynomial.
Convert an array representing the coefficients of a Hermite series, ordered from lowest degree to highest, to an array of the coefficients of the equivalent polynomial (relative to the "standard" basis) ordered from lowest to highest degree.
poly2herm
The easy way to do conversions between polynomial basis sets is to use the convert method of a class instance.
>>> from numpy.polynomial.hermite import herm2poly >>> herm2poly([ 1. , 2.75 , 0.5 , 0.375]) array([0., 1., 2., 3.])
Add one Hermite series to another.
Returns the sum of two Hermite series c1
+ c2
. The arguments
are sequences of coefficients ordered from lowest order term to
highest, i.e., [1,2,3] represents the series P_0 + 2*P_1 + 3*P_2.
hermsub, hermmulx, hermmul, hermdiv, hermpow
Unlike multiplication, division, etc., the sum of two Hermite series is a Hermite series (without having to "reproject" the result onto the basis set) so addition, just like that of "standard" polynomials, is simply "component-wise."
>>> from numpy.polynomial.hermite import hermadd >>> hermadd([1, 2, 3], [1, 2, 3, 4]) array([2., 4., 6., 4.])
Return the scaled companion matrix of c.
The basis polynomials are scaled so that the companion matrix is
symmetric when c
is an Hermite basis polynomial. This provides
better eigenvalue estimates than the unscaled case and for basis
polynomials the eigenvalues are guaranteed to be real if
numpy.linalg.eigvalsh
is used to obtain them.
Differentiate a Hermite series.
Returns the Hermite series coefficients c
differentiated m
times
along axis
. At each iteration the result is multiplied by scl
(the
scaling factor is for use in a linear change of variable). The argument
c
is an array of coefficients from low to high degree along each
axis, e.g., [1,2,3] represents the series 1*H_0 + 2*H_1 + 3*H_2
while [[1,2],[1,2]] represents 1*H_0(x)*H_0(y) + 1*H_1(x)*H_0(y) +
2*H_0(x)*H_1(y) + 2*H_1(x)*H_1(y) if axis=0 is x and axis=1 is
y.
c
is multidimensional the
different axis correspond to different variables with the degree in
each axis given by the corresponding index.scl
. The end result is
multiplication by scl**m. This is for use in a linear change of
variable. (Default: 1)Axis over which the derivative is taken. (Default: 0).
hermint
In general, the result of differentiating a Hermite series does not resemble the same operation on a power series. Thus the result of this function may be "unintuitive," albeit correct; see Examples section below.
>>> from numpy.polynomial.hermite import hermder >>> hermder([ 1. , 0.5, 0.5, 0.5]) array([1., 2., 3.]) >>> hermder([-0.5, 1./2., 1./8., 1./12., 1./16.], m=2) array([1., 2., 3.])
Divide one Hermite series by another.
Returns the quotient-with-remainder of two Hermite series
c1
/ c2
. The arguments are sequences of coefficients from lowest
order "term" to highest, e.g., [1,2,3] represents the series
P_0 + 2*P_1 + 3*P_2.
hermadd, hermsub, hermmulx, hermmul, hermpow
In general, the (polynomial) division of one Hermite series by another results in quotient and remainder terms that are not in the Hermite polynomial basis set. Thus, to express these results as a Hermite series, it is necessary to "reproject" the results onto the Hermite basis set, which may produce "unintuitive" (but correct) results; see Examples section below.
>>> from numpy.polynomial.hermite import hermdiv >>> hermdiv([ 52., 29., 52., 7., 6.], [0, 1, 2]) (array([1., 2., 3.]), array([0.])) >>> hermdiv([ 54., 31., 52., 7., 6.], [0, 1, 2]) (array([1., 2., 3.]), array([2., 2.])) >>> hermdiv([ 53., 30., 52., 7., 6.], [0, 1, 2]) (array([1., 2., 3.]), array([1., 1.]))
Least squares fit of Hermite series to data.
Return the coefficients of a Hermite series of degree deg
that is the
least squares fit to the data values y
given at points x
. If y
is
1-D the returned coefficients will also be 1-D. If y
is 2-D multiple
fits are done, one for each column of y
, and the resulting
coefficients are stored in the corresponding columns of a 2-D return.
The fitted polynomial(s) are in the form
where n
is deg
.
deg
is a single integer
all terms up to and including the deg
'th term are included in the
fit. For NumPy versions >= 1.11.0 a list of integers specifying the
degrees of the terms to include may be used instead.M
,), optionaly
was 2-D,
the coefficients for the data in column k of y
are in column
k
.These values are only returned if full == True
rcond
.For more details, see numpy.linalg.lstsq
.
The rank of the coefficient matrix in the least-squares fit is deficient. The warning is only raised if full == False. The warnings can be turned off by
>>> import warnings >>> warnings.simplefilter('ignore', np.RankWarning)
numpy.polynomial.chebyshev.chebfit numpy.polynomial.legendre.legfit numpy.polynomial.laguerre.lagfit numpy.polynomial.polynomial.polyfit numpy.polynomial.hermite_e.hermefit hermval : Evaluates a Hermite series. hermvander : Vandermonde matrix of Hermite series. hermweight : Hermite weight function numpy.linalg.lstsq : Computes a least-squares fit from the matrix. scipy.interpolate.UnivariateSpline : Computes spline fits.
The solution is the coefficients of the Hermite series p
that
minimizes the sum of the weighted squared errors
where the wj are the weights. This problem is solved by setting up the (typically) overdetermined matrix equation
where V
is the weighted pseudo Vandermonde matrix of x
, c
are the
coefficients to be solved for, w
are the weights, y
are the
observed values. This equation is then solved using the singular value
decomposition of V
.
If some of the singular values of V
are so small that they are
neglected, then a RankWarning
will be issued. This means that the
coefficient values may be poorly determined. Using a lower order fit
will usually get rid of the warning. The rcond
parameter can also be
set to a value smaller than its default, but the resulting fit may be
spurious and have large contributions from roundoff error.
Fits using Hermite series are probably most useful when the data can be
approximated by sqrt(w(x)) * p(x), where w(x)
is the Hermite
weight. In that case the weight sqrt(w(x[i])) should be used
together with data values y[i]/sqrt(w(x[i])). The weight function is
available as hermweight
.
[1] | Wikipedia, "Curve fitting", https://en.wikipedia.org/wiki/Curve_fitting |
>>> from numpy.polynomial.hermite import hermfit, hermval >>> x = np.linspace(-10, 10) >>> err = np.random.randn(len(x))/10 >>> y = hermval(x, [1, 2, 3]) + err >>> hermfit(x, y, 2) array([1.0218, 1.9986, 2.9999]) # may vary
Generate a Hermite series with given roots.
The function returns the coefficients of the polynomial
in Hermite form, where the r_n
are the roots specified in roots
.
If a zero has multiplicity n, then it must appear in roots
n times.
For instance, if 2 is a root of multiplicity three and 3 is a root of
multiplicity 2, then roots
looks something like [2, 2, 2, 3, 3]. The
roots can appear in any order.
If the returned coefficients are c
, then
The coefficient of the last term is not generally 1 for monic polynomials in Hermite form.
out
is a
real array, if some of the roots are complex, then out
is complex
even if all the coefficients in the result are real (see Examples
below).numpy.polynomial.polynomial.polyfromroots numpy.polynomial.legendre.legfromroots numpy.polynomial.laguerre.lagfromroots numpy.polynomial.chebyshev.chebfromroots numpy.polynomial.hermite_e.hermefromroots
>>> from numpy.polynomial.hermite import hermfromroots, hermval >>> coef = hermfromroots((-1, 0, 1)) >>> hermval((-1, 0, 1), coef) array([0., 0., 0.]) >>> coef = hermfromroots((-1j, 1j)) >>> hermval((-1j, 1j), coef) array([0.+0.j, 0.+0.j])
Gauss-Hermite quadrature.
Computes the sample points and weights for Gauss-Hermite quadrature. These sample points and weights will correctly integrate polynomials of degree 2*deg − 1 or less over the interval [ − inf, inf] with the weight function f(x) = exp( − x2).
The results have only been tested up to degree 100, higher degrees may be problematic. The weights are determined by using the fact that
where c is a constant independent of k and xk is the k'th root of Hn, and then scaling the results to get the right value when integrating 1.
Evaluate a 2-D Hermite series on the Cartesian product of x and y.
This function returns the values:
where the points (a, b)
consist of all pairs formed by taking
a
from x
and b
from y
. The resulting points form a grid with
x
in the first dimension and y
in the second.
The parameters x
and y
are converted to arrays only if they are
tuples or a lists, otherwise they are treated as a scalars. In either
case, either x
and y
or their elements must support multiplication
and addition both with themselves and with the elements of c
.
If c
has fewer than two dimensions, ones are implicitly appended to
its shape to make it 2-D. The shape of the result will be c.shape[2:] +
x.shape.
x
and y
. If x
or y
is a list or
tuple, it is first converted to an ndarray, otherwise it is left
unchanged and, if it isn't an ndarray, it is treated as a scalar.c
has dimension
greater than two the remaining indices enumerate multiple sets of
coefficients.x
and y
.hermval, hermval2d, hermval3d, hermgrid3d
Evaluate a 3-D Hermite series on the Cartesian product of x, y, and z.
This function returns the values:
where the points (a, b, c)
consist of all triples formed by taking
a
from x
, b
from y
, and c
from z
. The resulting points form
a grid with x
in the first dimension, y
in the second, and z
in
the third.
The parameters x
, y
, and z
are converted to arrays only if they
are tuples or a lists, otherwise they are treated as a scalars. In
either case, either x
, y
, and z
or their elements must support
multiplication and addition both with themselves and with the elements
of c
.
If c
has fewer than three dimensions, ones are implicitly appended to
its shape to make it 3-D. The shape of the result will be c.shape[3:] +
x.shape + y.shape + z.shape.
x
, y
, and z
. If x
,`y`, or z
is a
list or tuple, it is first converted to an ndarray, otherwise it is
left unchanged and, if it isn't an ndarray, it is treated as a
scalar.c
has dimension
greater than two the remaining indices enumerate multiple sets of
coefficients.x
and y
.hermval, hermval2d, hermgrid2d, hermval3d
Integrate a Hermite series.
Returns the Hermite series coefficients c
integrated m
times from
lbnd
along axis
. At each iteration the resulting series is
multiplied by scl
and an integration constant, k
, is added.
The scaling factor is for use in a linear change of variable. ("Buyer
beware": note that, depending on what one is doing, one may want scl
to be the reciprocal of what one might expect; for more information,
see the Notes section below.) The argument c
is an array of
coefficients from low to high degree along each axis, e.g., [1,2,3]
represents the series H_0 + 2*H_1 + 3*H_2 while [[1,2],[1,2]]
represents 1*H_0(x)*H_0(y) + 1*H_1(x)*H_0(y) + 2*H_0(x)*H_1(y) +
2*H_1(x)*H_1(y) if axis=0 is x and axis=1 is y.
scl
before the integration constant is added. (Default: 1)Axis over which the integral is taken. (Default: 0).
hermder
Note that the result of each integration is multiplied by scl
.
Why is this important to note? Say one is making a linear change of
variable u = ax + b in an integral relative to x
. Then
dx = du ⁄ a, so one will need to set scl
equal to
1 ⁄ a - perhaps not what one would have first thought.
Also note that, in general, the result of integrating a C-series needs to be "reprojected" onto the C-series basis set. Thus, typically, the result of this function is "unintuitive," albeit correct; see Examples section below.
>>> from numpy.polynomial.hermite import hermint >>> hermint([1,2,3]) # integrate once, value 0 at 0. array([1. , 0.5, 0.5, 0.5]) >>> hermint([1,2,3], m=2) # integrate twice, value & deriv 0 at 0 array([-0.5 , 0.5 , 0.125 , 0.08333333, 0.0625 ]) # may vary >>> hermint([1,2,3], k=1) # integrate once, value 1 at 0. array([2. , 0.5, 0.5, 0.5]) >>> hermint([1,2,3], lbnd=-1) # integrate once, value 0 at -1 array([-2. , 0.5, 0.5, 0.5]) >>> hermint([1,2,3], m=2, k=[1,2], lbnd=-1) array([ 1.66666667, -0.5 , 0.125 , 0.08333333, 0.0625 ]) # may vary
Hermite series whose graph is a straight line.
numpy.polynomial.polynomial.polyline numpy.polynomial.chebyshev.chebline numpy.polynomial.legendre.legline numpy.polynomial.laguerre.lagline numpy.polynomial.hermite_e.hermeline
>>> from numpy.polynomial.hermite import hermline, hermval >>> hermval(0,hermline(3, 2)) 3.0 >>> hermval(1,hermline(3, 2)) 5.0
Multiply one Hermite series by another.
Returns the product of two Hermite series c1
* c2
. The arguments
are sequences of coefficients, from lowest order "term" to highest,
e.g., [1,2,3] represents the series P_0 + 2*P_1 + 3*P_2.
hermadd, hermsub, hermmulx, hermdiv, hermpow
In general, the (polynomial) product of two C-series results in terms that are not in the Hermite polynomial basis set. Thus, to express the product as a Hermite series, it is necessary to "reproject" the product onto said basis set, which may produce "unintuitive" (but correct) results; see Examples section below.
>>> from numpy.polynomial.hermite import hermmul >>> hermmul([1, 2, 3], [0, 1, 2]) array([52., 29., 52., 7., 6.])
Multiply a Hermite series by x.
Multiply the Hermite series c
by x, where x is the independent
variable.
hermadd, hermsub, hermmul, hermdiv, hermpow
The multiplication uses the recursion relationship for Hermite polynomials in the form
>>> from numpy.polynomial.hermite import hermmulx >>> hermmulx([1, 2, 3]) array([2. , 6.5, 1. , 1.5])
Raise a Hermite series to a power.
Returns the Hermite series c
raised to the power pow
. The
argument c
is a sequence of coefficients ordered from low to high.
i.e., [1,2,3] is the series P_0 + 2*P_1 + 3*P_2.
hermadd, hermsub, hermmulx, hermmul, hermdiv
>>> from numpy.polynomial.hermite import hermpow >>> hermpow([1, 2, 3], 2) array([81., 52., 82., 12., 9.])
Compute the roots of a Hermite series.
Return the roots (a.k.a. "zeros") of the polynomial
out
is also real, otherwise it is complex.numpy.polynomial.polynomial.polyroots numpy.polynomial.legendre.legroots numpy.polynomial.laguerre.lagroots numpy.polynomial.chebyshev.chebroots numpy.polynomial.hermite_e.hermeroots
The root estimates are obtained as the eigenvalues of the companion matrix, Roots far from the origin of the complex plane may have large errors due to the numerical instability of the series for such values. Roots with multiplicity greater than 1 will also show larger errors as the value of the series near such points is relatively insensitive to errors in the roots. Isolated roots near the origin can be improved by a few iterations of Newton's method.
The Hermite series basis polynomials aren't powers of x
so the
results of this function may seem unintuitive.
>>> from numpy.polynomial.hermite import hermroots, hermfromroots >>> coef = hermfromroots([-1, 0, 1]) >>> coef array([0. , 0.25 , 0. , 0.125]) >>> hermroots(coef) array([-1.00000000e+00, -1.38777878e-17, 1.00000000e+00])
Subtract one Hermite series from another.
Returns the difference of two Hermite series c1
- c2
. The
sequences of coefficients are from lowest order term to highest, i.e.,
[1,2,3] represents the series P_0 + 2*P_1 + 3*P_2.
hermadd, hermmulx, hermmul, hermdiv, hermpow
Unlike multiplication, division, etc., the difference of two Hermite series is a Hermite series (without having to "reproject" the result onto the basis set) so subtraction, just like that of "standard" polynomials, is simply "component-wise."
>>> from numpy.polynomial.hermite import hermsub >>> hermsub([1, 2, 3, 4], [1, 2, 3]) array([0., 0., 0., 4.])
Evaluate an Hermite series at points x.
If c
is of length n + 1
, this function returns the value:
The parameter x
is converted to an array only if it is a tuple or a
list, otherwise it is treated as a scalar. In either case, either x
or its elements must support multiplication and addition both with
themselves and with the elements of c
.
If c
is a 1-D array, then p(x)
will have the same shape as x
. If
c
is multidimensional, then the shape of the result depends on the
value of tensor
. If tensor
is true the shape will be c.shape[1:] +
x.shape. If tensor
is false the shape will be c.shape[1:]. Note that
scalars have shape (,).
Trailing zeros in the coefficients will be used in the evaluation, so they should be avoided if efficiency is a concern.
x
is a list or tuple, it is converted to an ndarray, otherwise
it is left unchanged and treated as a scalar. In either case, x
or its elements must support addition and multiplication with
with themselves and with the elements of c
.c
is multidimensional the
remaining indices enumerate multiple polynomials. In the two
dimensional case the coefficients may be thought of as stored in
the columns of c
.If True, the shape of the coefficient array is extended with ones
on the right, one for each dimension of x
. Scalars have dimension 0
for this action. The result is that every column of coefficients in
c
is evaluated for every element of x
. If False, x
is broadcast
over the columns of c
for the evaluation. This keyword is useful
when c
is multidimensional. The default value is True.
hermval2d, hermgrid2d, hermval3d, hermgrid3d
The evaluation uses Clenshaw recursion, aka synthetic division.
>>> from numpy.polynomial.hermite import hermval >>> coef = [1,2,3] >>> hermval(1, coef) 11.0 >>> hermval([[1,2],[3,4]], coef) array([[ 11., 51.], [115., 203.]])
Evaluate a 2-D Hermite series at points (x, y).
This function returns the values:
The parameters x
and y
are converted to arrays only if they are
tuples or a lists, otherwise they are treated as a scalars and they
must have the same shape after conversion. In either case, either x
and y
or their elements must support multiplication and addition both
with themselves and with the elements of c
.
If c
is a 1-D array a one is implicitly appended to its shape to make
it 2-D. The shape of the result will be c.shape[2:] + x.shape.
(x, y)
,
where x
and y
must have the same shape. If x
or y
is a list
or tuple, it is first converted to an ndarray, otherwise it is left
unchanged and if it isn't an ndarray it is treated as a scalar.c
has
dimension greater than two the remaining indices enumerate multiple
sets of coefficients.x
and y
.hermval, hermgrid2d, hermval3d, hermgrid3d
Evaluate a 3-D Hermite series at points (x, y, z).
This function returns the values:
The parameters x
, y
, and z
are converted to arrays only if
they are tuples or a lists, otherwise they are treated as a scalars and
they must have the same shape after conversion. In either case, either
x
, y
, and z
or their elements must support multiplication and
addition both with themselves and with the elements of c
.
If c
has fewer than 3 dimensions, ones are implicitly appended to its
shape to make it 3-D. The shape of the result will be c.shape[3:] +
x.shape.
(x, y, z)
, where x
, y
, and z
must have the same shape. If
any of x
, y
, or z
is a list or tuple, it is first converted
to an ndarray, otherwise it is left unchanged and if it isn't an
ndarray it is treated as a scalar.c
has dimension
greater than 3 the remaining indices enumerate multiple sets of
coefficients.x
, y
, and z
.hermval, hermval2d, hermgrid2d, hermgrid3d
Pseudo-Vandermonde matrix of given degree.
Returns the pseudo-Vandermonde matrix of degree deg
and sample points
x
. The pseudo-Vandermonde matrix is defined by
where 0 <= i <= deg
. The leading indices of V
index the elements of
x
and the last index is the degree of the Hermite polynomial.
If c
is a 1-D array of coefficients of length n + 1
and V
is the
array V = hermvander(x, n), then np.dot(V, c) and
hermval(x, c) are the same up to roundoff. This equivalence is
useful both for least squares fitting and for the evaluation of a large
number of Hermite series of the same degree and sample points.
x
is
scalar it is converted to a 1-D array.x
.>>> from numpy.polynomial.hermite import hermvander >>> x = np.array([-1, 0, 1]) >>> hermvander(x, 3) array([[ 1., -2., 2., 4.], [ 1., 0., -2., -0.], [ 1., 2., 2., -4.]])
Pseudo-Vandermonde matrix of given degrees.
Returns the pseudo-Vandermonde matrix of degrees deg
and sample
points (x, y)
. The pseudo-Vandermonde matrix is defined by
where 0 <= i <= deg[0]
and 0 <= j <= deg[1]
. The leading indices of
V
index the points (x, y)
and the last index encodes the degrees of
the Hermite polynomials.
If V = hermvander2d(x, y, [xdeg, ydeg]), then the columns of V
correspond to the elements of a 2-D coefficient array c
of shape
(xdeg + 1, ydeg + 1) in the order
and np.dot(V, c.flat) and hermval2d(x, y, c) will be the same up to roundoff. This equivalence is useful both for least squares fitting and for the evaluation of a large number of 2-D Hermite series of the same degrees and sample points.
x
and y
.hermvander, hermvander3d, hermval2d, hermval3d
Pseudo-Vandermonde matrix of given degrees.
Returns the pseudo-Vandermonde matrix of degrees deg
and sample
points (x, y, z)
. If l, m, n
are the given degrees in x, y, z
,
then The pseudo-Vandermonde matrix is defined by
where 0 <= i <= l
, 0 <= j <= m
, and 0 <= j <= n
. The leading
indices of V
index the points (x, y, z)
and the last index encodes
the degrees of the Hermite polynomials.
If V = hermvander3d(x, y, z, [xdeg, ydeg, zdeg]), then the columns
of V
correspond to the elements of a 3-D coefficient array c
of
shape (xdeg + 1, ydeg + 1, zdeg + 1) in the order
and np.dot(V, c.flat) and hermval3d(x, y, z, c) will be the same up to roundoff. This equivalence is useful both for least squares fitting and for the evaluation of a large number of 3-D Hermite series of the same degrees and sample points.
x
, y
, and z
.hermvander, hermvander3d, hermval2d, hermval3d
Weight function of the Hermite polynomials.
The weight function is exp( − x2) and the interval of integration is [ − inf, inf]. the Hermite polynomials are orthogonal, but not normalized, with respect to this weight function.
x
.poly2herm(pol)
Convert a polynomial to a Hermite series.
Convert an array representing the coefficients of a polynomial (relative to the "standard" basis) ordered from lowest degree to highest, to an array of the coefficients of the equivalent Hermite series, ordered from lowest to highest degree.
herm2poly
The easy way to do conversions between polynomial basis sets is to use the convert method of a class instance.
>>> from numpy.polynomial.hermite import poly2herm >>> poly2herm(np.arange(4)) array([1. , 2.75 , 0.5 , 0.375])
Evaluate a normalized Hermite polynomial.
Compute the value of the normalized Hermite polynomial of degree n at the points x.
This function is needed for finding the Gauss points and integration weights for high degrees. The values of the standard Hermite functions overflow when n >= 207.