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. More samples 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:

from gvar.ode import Integrator

def f(x, y):
    ...

odeint = Integrator(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 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.

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).
  • 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.

Table Of Contents

Previous topic

gvar.dataset - Random Data Sets