Other operations

Var.diff(axis=0, n=1)

Computes the forward difference along the given axis. Mimics the same behaviour of the np.diff() function.

Parameters
axisstring, Axis class, or int

Axis along which to compute differences.

nint (optional)

Number of times values are differenced.

Returns
dvarVar

New variable containing n-th differenced values.

Examples

>>> import pygeode as pyg
>>> v = pyg.yearlessn(5)
>>> v[:]
array([0., 1., 2., 3., 4.])
>>> v.diff('time')[:]
array([1., 1., 1., 1.])
>>> v.diff('time', 2)[:]
array([0., 0., 0.])
Var.deriv(daxis, dx=None, df='centre')

Computes derivative along the given axis.

Parameters
daxisstring, Axis class, or int

Axis along which to compute derivative.

dxVar, or None (optional)

Coordinate with respect to which to take the derivative (see notes). Must share axis along which the derivative is being taken. If None, the coordinate axis is used.

dfstring (optional)

Type of derivative to take. One of ‘left’, ‘right’, ‘centre’, or ‘2’. See notes.

Returns
dvarVar

Numerical derivative of var

Notes

The derivative is computed using a forward (df = ‘right’), backward (df = ‘left’), or centred (df = ‘centre’) difference approximation; for instance, the forward difference is computed as:

dvar[i] = (var[i+1] - var[i]) / (dx[i+1] - dx[i]).

One-sided differences are used at the axis boundaries so that dvar is defined on the same axis as var. The second derivative can also be computed (df = ‘2’)

Examples

>>> import pygeode as pyg, numpy as np
>>> from pygeode.tutorial import t1
>>> print(t1.Temp.deriv('lon')) # Compute simple derivative
<Var 'dTemp'>:
  Shape:  (lat,lon)  (31,60)
  Axes:
    lat <Lat>      :  90 S to 90 N (31 values)
    lon <Lon>      :  0 E to 354 E (60 values)
  Attributes:
    {}
  Type:  DerivativeVar (dtype="float64")
>>> x = 6317e3 * pyg.cosd(t1.lat) * np.pi / 180. * t1.lon
>>> print(t1.Temp.deriv('lon', dx=x, df='2')) # Compute 2nd derivative with respect to geometric length
<Var 'd2Temp'>:
  Shape:  (lat,lon)  (31,60)
  Axes:
    lat <Lat>      :  90 S to 90 N (31 values)
    lon <Lon>      :  0 E to 354 E (60 values)
  Attributes:
    {}
  Type:  SecondDerivativeVar (dtype="float64")
Var.integrate(iaxis, dx=None, v0=None, order=1, type='trapz')

Computes an indefinite integral along the given axis.

Parameters
iaxisstring, Axis class, or int

Axis along which to compute integral.

dxVar, or None (optional)

Coordinate with respect to which to integrate (see notes). Must share axis along which the derivative is being taken. If None, the coordinate axis is used.

v0float, Var, or None (optional)

Constant of integration. See notes.

order: int (1 or -1)

Direction along axis to integrate. 1 corresponds to an integration from the first to the last element, while -1 integrates in the other direction.

typestring (optional)

Type of numerical integral to take. One of ‘trapz’, ‘rectr’, or ‘rectl’; defaults to ‘trapz’. See notes.

Returns
ivarVar

Numerical integral of var

Notes

Possible integration methods are
  • ‘trapz’: trapezoidal rule

  • ‘rectl’: left rectangle method or Riemann sum

  • ‘rectr’: right rectangle method or Riemann sum

Examples

>>> import pygeode as pyg
>>> from pygeode.tutorial import t1
>>> print(t1.Temp.integrate('lon')) # Compute simple derivative
<Var 'iTemp'>:
  Shape:  (lat,lon)  (31,60)
  Axes:
    lat <Lat>      :  90 S to 90 N (31 values)
    lon <Lon>      :  0 E to 354 E (60 values)
  Attributes:
    {}
  Type:  IntegrateVar (dtype="float64")
Var.interpolate(inaxis, outaxis, inx=None, outx=None, interp_type='cspline', d_below=nan, d_above=nan, transpose=True, omit_nonmonotonic=False)

Interpolates a variable along a single dimension.

Parameters
invarVar

The input variable

inaxisAxis

The input axis being interpolated from

outaxisAxis

The output axis being interpolated to.

inxVar (optional)

The coordinates from which we are interpolating (must be conformable to the input var). If not provided, the values of inaxis are used. This can be either a one dimensional field defined on inaxis, or a multidimensional field.

outxVar (optional)

The coordinates to which we are interpolating (must be conformable to the output var). If not provided, the values of outaxis are used. This can be either a one dimensional field defined on outaxis, or a multidimensional field.

interp_typestring (optional)

The type of interpolation. One of ‘linear’, ‘polynomial’, ‘cspline’, ‘cspline_periodic’, ‘akima’, ‘akima_periodic’. Default is ‘cspline’ (cubic spline interpolation)

d_belowfloat (optional)

The slope for linearly extrapolating below the input data.

d_abovefloat (optional)

The slope for linearly extrapolating above the input data. By default, no extrapolation is done.

transposeboolean (optional)

If True, tranposes the output axes so that the axis to which we are interpolating replaces the input axis in the axes order. Otherwise the new axis becomes the last dimension. Default is True.

omit_nonmonotonicboolean (optional)

If True, the interpolation routine skips over any non-monotonic datapoints in the coordinate field from which we are interpolating. This is an experimental feature; do not use unless you know what you are doing! Default is False.

Returns
interpolatedVar

The interpolated variable

Notes

if inx or outx are not explicitly given, they will take on the values of inaxis and outaxis respectively. Note that inx and outx are assumed to have the same units, so if inaxis and outaxis have different units, you’ll need to override either inx or outx explicitly, and provide some kind of mapping field there.

Examples:

1) Simple interpolation along one axis (re-gridding):

To interpolate some data to 20 evenly spaced longitude values (starting at 0):
newvar = Interp (invar = data, inaxis = 'lon', outaxis = Lon(20))

2) Interpolation involving a change of coordinates:

Suppose we have some ozone data (o3) on a model vertical coordinate.
Suppose we also have pre-computed a pressure field (pfield) over these coordinates.
Suppose finally we have our desired pressure axis (paxis) that we want to interpolate to.

To interpolate onto pressure levels:
newvar = Interp (invar = o3, inaxis = 'eta', outaxis = paxis, inx = pfield)
levels).

Now, you may be asking yourself "shouldn't pressure be interpolated on a
log scale?".  This is true, but I tried to simplify things by assuming a
linear scale just to show the basics.  If you wanted to interpolate over
a log scale, then your interpolation coordinate would be log(pressure),
instead of pressure.  You'll have to explicitly provide this log coordinate
for both the inputs and outputs, for example:

newvar = Interp (invar = o3, inaxis = 'eta', outaxis = paxis,
                 inx = pfield.log(), outx = paxis.log()  )

Here, our input and output axes remain the same as before, but now we're
using log(pressure) internally as the coordinate over which to perform
the interpolation.
Var.smooth(saxis, kernel=15, fft=True)

Smooths this variable along saxis by convolving it with an averaging kernel. The returned variable is defined on the same axes.

Parameters
saxisany axis identifier (string, Axis, or int)

Axis over which the smoothing should be performed

kernelsequence or int (optional)

Averaging kernel with which to convolve this variable. Does not need to be normalized. If an integer is provided, a Hanning window is used of length kernel (numpy.hanning())

fftboolean (optional, True by default)

If True, scipy.signal.fftconvolve() is used to compute the convolution. Otherwise, scipy.signal.convolve() is used. In many cases the former is more efficient.

Returns
outVar

Smoothed variable.

Notes

When the convolution is performed, the source data is extended on either end of the axis being smoothed by reflecting the data by enough to ensure the returned variable is defined on the same grid as the source variable. That is, if the original data is t1, t2, .., tN, and the kernel is L items long, the convolved sequence is tj, t_j-1, t1, t1, t2, .. tN-1, tN, tN, tN-1, .. tN-l, where j = floor(L/2) and l = L - j - 1.

Var.fft_smooth(saxis, maxharm)

Smooths this variable along saxis by retaining leading Fourier components.

Parameters
saxisany axis identifier (string, Axis, or int)

Axis over which the smoothing should be performed

maxharmint

Maximum harmonic to retain.

Returns
outVar

Smoothed variable.

Notes

The variable data is Fourier transformed using np.fft.rfft() (if real) or np.fft.fft() (if complex). The coefficients for harmonics equal to and greater than maxharm are set to zero, then an inverse transform is applied.

Examples

>>> import pygeode as pyg, numpy as np
>>> tm = pyg.modeltime365n('1 Jan 2000', 365)
>>> v = pyg.cos(3 * 2 * np.pi * tm / 365.)
>>> np.std(v[:])
0.7071067811865476
>>> np.std(v.fft_smooth('time', 3)[:]) # This retains only the annual and semi-annual cycle
1.0983314772510289e-16
>>> np.std(v.fft_smooth('time', 4)[:]) # This retains up to the third harmonic
0.7071067811865476
Var.composite(**kwargs)

Creates a composite based on this variable.

Parameters
<axis selection>string

A single axis selection string (similar to Var.__call__()) that specifies the central ‘date’ of each event to composite (although composites need not be constructed along time axes). See Notes.

evlenint

Length of segement around each central ‘date’ to extract.

evoffint

Offset from the central ‘dates’ to include. A positive value will lead to dates prior to the central ‘date’ being included.

Returns
cvarVar

Composite variable. The axis along which composites are to be constructed is replaced by an Event axis and an Offset axis.

Notes

The axis matched is used as the composite axis and the returned indices are the key dates. evlen is required; it can either be an integer or a list of integers specifying the length of each event. evoff is a single integer. If the requested composite extends past the ends of the underlying axis, the variable will contain NaNs.

Examples

>>> import pygeode as pyg
>>> from pygeode.tutorial import t2
>>> dts = ['12 May 2012', '15 Aug 2015', '1 Feb 2018', '28 Dec 2020']
>>> cT = t2.Temp.composite(l_time = dts, evlen = 15, evoff = 10)
>>> print(cT)
<Var 'Temp'>:
  Shape:  (event,time,pres,lat,lon)  (4,15,20,31,60)
  Axes:
    event <Event>  :  1  to 4  (4 values)
    time <Yearless>:  day -10, 00:00:00 to day 4, 00:00:00 (15 values)
    pres <Pres>    :  1000 hPa to 50 hPa (20 values)
    lat <Lat>      :  90 S to 90 N (31 values)
    lon <Lon>      :  0 E to 354 E (60 values)
  Attributes:
    {}
  Type:  CompositeVar (dtype="float64")
>>> cT(s_event=1, s_pres=50, s_lat=0, s_lon = 180)[:]
array([199.66766628, 199.69594407, 199.72479591, 199.75418884,
       199.78408917, 199.81446257, 199.84527408, 199.87648814,
       199.90806866, 199.939979  , 199.97218208, 200.00464036,
       200.03731592, 200.07017048, 200.10316548])
>>> cT(s_event=4, s_pres=50, s_lat=0, s_lon = 180)[:]
array([201.02847025, 201.04367089, 201.05782423, 201.07091237,
       201.08291874, 201.09382812, 201.10362669, 201.112302  ,
       201.11984303, 201.12624021, 201.1314854 , 201.13557193,
       201.1384946 , 201.14024969,          nan])
Var.flatten(inner=- 1, naxis=None)
Var.lag(iaxis, lags, reverse=False)

Adds a lag axis with offset values.

See Also:

Var class overview