Numerical Analysis Modules in gvar

gvar.GVars can be used in many numerical algorithms, to propagates errors through the algorithm. A code that is written in pure Python is likely to work well with gvar.GVars, perhaps with minor modifications. Here we describe some sample numerical codes, included in gvar, that have been adapted to work with gvar.GVars, as well as with floats. More examples will follow with time.

Cubic Splines

The module gvar.cspline implements a class for smoothing and/or interpolating one-dimensional data using cubic splines:

class gvar.cspline.CSpline(x, y, deriv=(None, None), warn=True)

Cubic spline resentation of a function.

Given N values of a function y[i] at N points x[i] for i=0..N-1 (the ‘knots’ of the spline), the code

from gvar.cspline import CSpline

f = CSpline(x, y)

defines a function f such that: a) f(x[i]) = y[i] for all i; and b) f(x) is continuous, as are its first and second derivatives. Function f(x) is a cubic polynomial between the knots x[i].

CSpline(x,y) creates a natural spline, which has zero second derivative at the end points, x[0] and x[-1]. More generally one can specify the derivatives of f(x) at one or both of the endpoints:

f = CSpline(x, y, deriv=[dydx_i, dydx_f])

where dydx_i is the derivative at x[0] and dydx_f is the derivative at x[-1]. Replacing either (or both) of these with None results in a derivative corresponding to zero second derivative at that boundary (i.e., a natural boundary).

Derivatives and integrals of the spline function can also be evaluated:

f.D(x) — first derivative at x;

f.D2(x) — second derivative at x;

f.integ(x) — integral from x[0] to x.

Parameters:
  • x (1-d sequence of numbers) – The knots of the spline, where the function values are specified.
  • y (1-d sequence) – Function values at the locations specified by x[i].
  • deriv (2-component sequence) – Derivatives at initial and final boundaries of the region specified by x[i]. Default value is None for each boundary.
  • warn – If True, warnings are generated when the spline function is called for x values that fall outside of the original range of xs used to define the spline. Default value is True; out-of-range warnings are suppressed if set to False.

Ordinary Differential Equations

The module gvar.ode implements two classes for integrating systems of first-order differential equations using an adaptive Runge-Kutta algorithm. One integrates scalar- or array-valued equations, while the other integrates dictionary-valued equations:

class gvar.ode.Integrator(deriv, tol=1e-05, h=None, hmin=None, analyzer=None)

Integrate dy/dx = deriv(x,y).

An Integrator object odeint integrates dy/dx = f(x,y) to obtain y(x1) from y0 = y(x0). y and f(x,y) can be scalars or numpy arrays. Typical usage is illustrated by the following code for integrating dy/dx = y:

from gvar.ode import Integrator

def f(x, y):
    return y

odeint = Integrator(deriv=f,  tol=1e-8)
y0 = 1.
y1 = odeint(y0, interval=(0, 1.))
y2 = odeint(y1, interval=(1., 2.))
...

Here the first call to odeint integrates the differential equation from x=0 to x=1 starting with y=y0 at x=0; the result is y1=exp(1), of course. Similarly the second call to odeint continues the integration from x=1 to x=2, giving y2=exp(2).

An alternative interface creates a new function which is the solution of the differential equation for specific initial conditions. The code above could be rewritten:

x0 = 0.         # initial conditions
y0 = 1.
y = Integrator(deriv=f, tol=1e-8).solution(x0, y0)
y1 = y(1)
y2 = y(2)
...

Here method Integrator.solution() returns a function y(x) where: a) y(x0) = y0; and b) y(x) uses the integator to integrate the differential equation to point x starting from the last point at which y was evaluated (or from x0 for the first call to y(x)).

The integrator uses an adaptive Runge-Kutta algorithm that adjusts the integrator’s step size to obtain relative accuracy tol in the solution. An initial step size can be set in the Integrator by specifying parameter h. A minimum step size hmin can also be specified; the Integrator raises an exception if the step size becomes smaller than hmin. The Integrator keeps track of the number of good steps, where h is increased, and the number of bad steps, where h is decreased and the step is repeated: odeint.ngood and odeint.nbad, respectively.

A custom criterion for step-size changes can be implemented by specifying a function for parameter delta. This is a function delta(yerr, y, delta_y) — of the estimated error yerr after a given step, the proposed value for y, and the proposed change delta_y in y — that returns a number to compare with tolerance tol. The step size is decreased and the step repeated if delta(yerr, y, delta_y) > tol; otherwise the step is accepted and the step size increased. The default definition of delta is roughly equivalent to:

import numpy as np
import gvar as gv

def delta(yerr, y, delta_y):
    return np.max(
        np.fabs(yerr) / (np.fabs(y) + np.fabs(delta_y) + gv.ode.TINY)
        )

A custom definition can be used to allow an Integrator to work with data types other than floats or numpy arrays of floats. All that is required of the data type is that it support ordinary arithmetic. Therefore, for example, defining delta(yerr, y, delta_y) with np.abs() instead of np.fabs() allows y to be complex valued. (Actually the default delta allows this as well.)

An analyzer analyzer(x,y) can be specified using parameter analyzer. This function is called after every full step of the integration, with the current values of x and y. Objects of type gvar.ode.Solution are examples of (simple) analyzers.

Parameters:
  • deriv – Function of x and y that returns dy/dx. The return value should have the same shape as y if arrays are used.
  • tol (float) – Relative accuracy in y relative to |y| + h|dy/dx| for each step in the integration. Any integration step that achieves less precision is repeated with a smaller step size. The step size is increased if precision is higher than needed. Default is 1e-5.
  • h (float or None) – Absolute value of initial step size. The default value equals the entire width of the integration interval.
  • hmin (float or None) – Smallest step size allowed. An exception is raised if a smaller step size is needed. This is mostly useful for preventing infinite loops caused by programming errors. The default value is zero (which does not prevent infinite loops).
  • delta – Function delta(yerr, y, delta_y) that returns a number to be compared with tol at each integration step: if it is larger than tol, the step is repeated with a smaller step size; if it is smaller the step is accepted and a larger step size used for the subsequent step. Here yerr is an estimate of the error in y on the last step; y is the proposed value; and delta_y is the change in y over the last step.
  • analyzer – Function of x and y that is called after each step of the integration. This can be used to analyze intermediate results.
class gvar.ode.DictIntegrator(deriv, tol=1e-05, h=None, hmin=None, analyzer=None)

Integrate dy/dx = deriv(x,y) where y is a dictionary.

An DictIntegrator object odeint integrates dy/dx = f(x,y) to obtain y(x1) from y0 = y(x0). y and f(x,y) are dictionary types having the same keys, and containing scalars and/or numpy arrays as values. Typical usage is:

from gvar.ode import DictIntegrator

def f(x, y):
    ...

odeint = DictIntegrator(deriv=f,  tol=1e-8)
y1 = odeint(y0, interval=(x0, x1))
y2 = odeint(y1, interval=(x1, x2))
...

The first call to odeint integrates from x=x0 to x=x1, returning y1=y(x1). The second call continues the integration to x=x2, returning y2=y(x2). Any of the initial parameters can be reset in the calls to odeint: for example,

y2 = odeint(y1, interval=(x1, x2), tol=1e-10)

The integrator uses an adaptive Runge-Kutta algorithm that adjusts the integrator’s step size to obtain relative accuracy tol in the solution. An initial step size can be set in the DictIntegrator by specifying parameter h. A minimum ste psize hmin can also be specified; the Integrator raises an exception if the step size becomes smaller than hmin. The DictIntegrator keeps track of the number of good steps, where h is increased, and the number of bad steps, where h is decreases and the step is repeated: odeint.ngood and odeint.nbad, respectively.

An analyzer analyzer(x,y) can be specified using parameter analyzer. This function is called after every full step of the integration with the current values of x and y. Objects of type gvar.ode.Solution are examples of (simple) analyzers.

Parameters:
  • deriv – Function of x and y that returns dy/dx. The return value should be a dictionary with the same keys as y, and values that have the same shape as the corresponding values in y.
  • tol (float) – Relative accuracy in y relative to |y| + h|dy/dx| for each step in the integration. Any integration step that achieves less precision is repeated with a smaller step size. The step size is increased if precision is higher than needed.
  • h (float) – Absolute value of initial step size. The default value equals the entire width of the integration interval.
  • hmin (float) – Smallest step size allowed. An exception is raised if a smaller step size is needed. This is mostly useful for preventing infinite loops caused by programming errors. The default value is zero (which does not prevent infinite loops).
  • analyzer – Function of x and y that is called after each step of the integration. This can be used to analyze intermediate results.

A simple analyzer class is:

class gvar.ode.Solution

ODE analyzer for storing intermediate values.

Usage: eg, given

odeint = Integrator(...)
soln = Solution()
y0 = ...
y = odeint(y0, interval=(x0, x), analyzer=soln)

then the soln.x[i] are the points at which the integrator evaluated the solution, and soln.y[i] is the solution of the differential equation at that point.

Power Series

This module provides tools for manipulating power series approximations of functions. A function’s power series is specified by the coefficients in its Taylor expansion with respect to an independent variable, say x:

f(x) = f(0) + f'(0)*x + (f''(0)/2)*x**2 + (f'''(0)/6)*x**3 + ...
     = f0 + f1*x + f2*x**2 + f3*x**3 + ...

In practice a power series is different from a polynomial because power series, while infinite order in principle, are truncated at some finite order in numerical applications. The order of a power series is the highest power of x that is retained in the approximation; coefficients for still higher-order terms are assumed to be unknown (as opposed to zero).

Taylor’s theorem can be used to generate power series for functions of power series:

g(f(x)) = g(f0) + g'(f0)*(f(x)-f0) + (g''(f0)/2)*(f(x)-f0)**2 + ...
        = g0 + g1*x + g2*x**2 + ...

This allows us to define a full calculus for power series, where arithmetic expressions and (sufficiently differentiable) functions of power series return new power series.

Power series arithmetic

Class PowerSeries provides a numerical implementation of the power series calculus. PowerSeries([f0,f1,f2,f3...]) is a numerical representation of a power series with coefficients f0, f1, f2, f3... (as in f(x) above). Thus, for example, we can define a 4th-order power series approximation f to exp(x)=1+x+x**2/2+... using

>>> from gvar.powerseries import *
>>> f = PowerSeries([1., 1., 1/2., 1/6., 1/24.])
>>> print f             # print the coefficients
[ 1.          1.          0.5         0.16666667  0.04166667]

Arithmetic expressions involving instances of class PowerSeries are themselves PowerSeries as in, for example,

>>> print 1/f           # power series for exp(-x)
[ 1.         -1.          0.5        -0.16666667  0.04166667]
>>> print log(f)        # power series for x
[ 0.  1.  0. -0.  0.]
>>> print f/f           # power series for 1
[ 1.  0.  0.  0.  0.]

The standard arithmetic operators (+,-,*,/,=,**) are supported, as are the usual elementary functions (exp, log, sin, cos, tan ...). Different PowerSeries can be combined arithmetically to create new PowerSeries; the order of the result is that of the operand with the lowest order.

PowerSeries can be differentiated and integrated:

>>> print f.deriv()     # derivative of exp(x)
[ 1.          1.          0.5         0.16666667]
>>> print f.integ()     # integral of exp(x) (from x=0)
[ 0.          1.          0.5         0.16666667  0.04166667  0.00833333]

Each PowerSeries represents a function. The PowerSeries for a function of a function is easily obtained. For example, assume f represents function f(x)=exp(x), as above, and g represents g(x)=log(1+x):

>>> g = PowerSeries([0, 1, -1/2., 1/3., -1/4.])

Then f(g) gives the PowerSeries for exp(log(1+x)) = 1 + x:

>>> print f(g)
[  1.0000e+00   1.0000e+00   0.0000e+00  -2.7755e-17 -7.6327e-17]

Individual coefficients from the powerseries can be accessed using array-element notation: for example,

>>> print f[0], f[1], f[2], f[3]
1.0 1.0 0.5 0.166666666667
>>> f[0] = f[0] - 1.
>>> print f             # f is now the power series for exp(x)-1
[ 0.          1.          0.5         0.16666667  0.04166667]

Numerical evaluation of power series

The power series can also be evaluated for a particular numerical value of x: continuing the example,

>>> x = 0.01
>>> print f(x)          # should be exp(0.01)-1 approximately
0.0100501670833
>>> print exp(x)-1      # verify that it is
0.0100501670842

The independent variable x could be of any arithmetic type (it need not be a float).

Taylor expansions of Python functions

PowerSeries can be used to compute Taylor series for more-or-less arbitrary pure-Python functions provided the functions are locally analytic (or at least sufficiently differentiable). To compute the N-th order expansion of a Python function g(x), first create a N-th order PowerSeries variable that represents the expansion parameter: say, x = PowerSeries([0.,1.],order=N). The Taylor series for function g is then given by g_taylor = g(x) which is a PowerSeries instance. For example, consider:

>>> from gvar.powerseries import *
>>> def g(x):           # an example of a Python function
...     return 0.5/sqrt(1+x) + 0.5/sqrt(1-x)
... 
>>> x = PowerSeries([0.,1.],order=5)    # Taylor series for x
>>> print x
[ 0.  1.  0.  0.  0.  0.]
>>> g_taylor = g(x)     # Taylor series for g(x) about x=0
>>> print g_taylor
[ 1.         0.         0.375      0.         0.2734375  0.       ]
>>> exp_taylor = exp(x) # Taylor series for exp(x) about x=0
>>> print exp_taylor
[ 1.          1.          0.5         0.16666667  0.04166667  0.00833333]
class gvar.powerseries.PowerSeries(c=None, order=None)

Power series representation of a function.

The power series created by PowerSeries(c) corresponds to:

c[0] + c[1]*x + c[2]*x**2 + ... .

The order of the power series is normally determined by the length of the input list c. This can be overridden by specifying the order of the power series using the order parameter. The list of c[i]s is then padded with zeros if c is too short, or truncated if it is too long. Omitting c altogether results in a power series all of whose coefficients are zero. Individual series coefficients are accessed using array/list notation: for example, the 3rd-order coefficient of PowerSeries p is p[3]. The order of p is p.order. PowerSeries should work for coefficients of any data type that supports ordinary arithmetic.

Arithmetic expressions of PowerSeries variables yield new PowerSeries results that represent the power series expansion of the expression. Expressions can include the standard mathematical functions (log, exp, sqrt, sin, cos, tan...). PowerSeries can also be differentiated (p.deriv()) and integrated (p.integ()).

Parameters:
  • c (list or array) – Power series coefficients (optional if parameter order specified).
  • order (integer) – Highest power in power series (optional if parameter c specified).
deriv(n=1)

Compute n-th derivative of self.

Parameters:n (positive integer) – Number of derivatives.
Returns:n-th derivative of self.
integ(n=1, x0=None)

Compute n-th indefinite integral of self.

If x0 is specified, then the definite integral, integrating from point x0, is returned.

Parameters:
  • n (integer) – Number of integrations.
  • x0 – Starting point for definite integral (optional).
Returns:

n-th integral of self.

order

Highest power in power series.