Model Grid

The Grid class encapsulates all functionality related to the spatial discretization of the sphere. It provides access to coordinates, transformations between grid and spectral space, derivatives, quadrature, differential equation solvers, filtering, interpolation and a few other routines. Creating a Grid is usually the first thing to do after importing barotropic.

Grid

class barotropic.Grid(resolution=2.5, rsphere=6371200.0, omega=7.292e-05, ntrunc=None, legfunc='stored')

Bases: object

Regular lon-lat grid and operations for a spherical planet.

Note

Latitudes start at the North Pole, as is convention in spharm. This results in horizontally flipped images when plotting fields in matplotlib without explicitly specifying the latitude grid.

Parameters:
  • resolution (number) – Grid resolution (uniform) in degrees.

  • rsphere (number) – Radius of the planet in m.

  • omega (number) – Angular velocity of the planet in 1/s.

  • ntrunc (int) – Threshold for triangular truncation.

  • legfunc (str) – Parameter given to spharm.Spharmt.

The default values of rsphere and omega correspond to those of planet Earth. Consult the spharm documentation for further information on the ntrunc and legfunc parameters.

nlon, nlat

Number of longitudes/latitudes.

dlon, dlat

Longitude/latitude spacing in degrees.

lon, lat

Longitude/latitude coordinates in degrees (1D).

lon2, lat2

Longitude/latitude coordinates in degrees (2D).

dlam, dphi

Longitude/latitude spacing in radians.

lam, phi

Longitude/latitude coordinates in radians (1D).

lam2, phi2

Longitude/latitude coordinates in radians (2D).

fcor

Coriolis parameter in 1/s (1D).

fcor2

Coriolis parameter in 1/s (2D).

circumference(lat)

Circumference (in m) of the sphere at (a) given latitude(s).

Parameters:

lat (number | array) – Input latitude(s) in degrees.

Returns:

Circumference or array of circumferences.

coriolis(lat)

Coriolis parameter (in m) for a given latitude (in degrees).

Parameters:

lat (number | array) – Input latitude(s) in degrees.

Returns:

Coriolis parameter or array of Coriolis parameters.

derivative_meridional(f, order=4)

Finite difference first derivative in meridional direction (∂f/∂φ).

Parameters:
  • f (array) – 1D or 2D input field.

  • order (2 | 4) – Order of the finite-difference approximation.

Returns:

Array containing the derivative.

Accepts both 2D f = f(φ,λ) and 1D f = f(φ) fields. For 1D input, the derivatives at the poles are always set to zero, as the input is assumed to represent a zonally symmetric profile of some quantity (e.g. zonal-mean PV).

The 2nd order approximation uses a 3-point stencil, the 4th order approximation a 5-point stencil (centered in inner domain, offset at the edges). Stencil coefficients were obtained from http://web.media.mit.edu/~crtaylor/calculator.html.

divergence(u, v)

Gridded divergence from vector components.

Parameters:
  • u (2D array) – Zonal component of input vector field.

  • v (2D array) – Meridional component of input vector field.

Returns:

2D divergence field.

divergence_spectral(u, v)

divergence() with spectral output field.

equivalent_latitude(area)

Latitude such that surface up to North Pole is area-sized.

Parameters:

area (number | array) – input area (in m²).

Returns:

Latitude(s) in degrees.

Solve formula for area A, north of latitude φ, for φ:

A = 2 * pi * (1 - sin(φ)) * r²
φ = arcsin(1 - A / 2 / pi / r²)
filter_meridional(field, window, width=None)

Filter the input in meridional direction with the given window.

Parameters:
  • field (array) – 1 or 2-dimensional input signal.

  • window (str | array) – Window used in the convolution.

  • width (number) – Width of the window.

Returns:

Filtered field.

If width is None, window must be a gridded window given as a 1D array. Otherwise, window and width are given to get_filter_window() to obtain a window function.

filter_zonal(field, window, width=None)

Filter the input in zonal direction with the given window.

Parameters:
  • field (array) – 1 or 2-dimensional input signal.

  • window (str | array) – Window used in the convolution.

  • width (number) – Width of the window.

Returns:

Filtered field.

If width is None, window must be a gridded window given as a 1D array. Otherwise, window and width are given to get_filter_window() to obtain a window function.

get_filter_window(window, width)

Wraps scipy.signal.get_window() for degree input.

Parameters:
  • window (str) – The type of window to create.

  • width (number) – Width of the window in degrees.

Returns:

Gridded window function.

Window widths are restricted to odd numbers of gridpoints so windows can properly be centered on a gridpoint during convolution. The returned window array is normalized such that it sums to 1.

gradient(f)

Gridded vector gradient of a horizontal (2D) field.

Parameters:

f (array) – 2D input field.

Returns:

A tuple with the two components of the horizontal gradient.

Horizontal gradient on the sphere:

∇f(φ,λ) = (1/r ∂f/∂φ, 1/(r sin(φ)) ∂f/∂λ)
property gridpoint_area

Surface area of each gridpoint as a function of latitude.

The area associated with a gridpoint (λ, φ) in a regular lon-lat grid:

r² * dλ * ( sin(φ + dφ/2) - sin(φ - dφ/2) )
interp1d_meridional(field, lat, pole_values, kind='linear')

Meridional interpolation function.

Parameters:
  • field (array) – 1D meridional profile or 2D field on the input latitudes.

  • lat (array) – Input latitudes in degrees.

  • pole_values ((number, number)) – Values at the poles (North, South), required to complete the interpolation.

  • kind (string) – Given to scipy.interpolate.interp1d().

Returns:

Configured scipy.interpolate.interp1d() interpolation function.

See also: Grid.interpolate_meridional().

interpolate_meridional(*interp1d_args, **interp1d_kwargs)

Interpolation onto the regular latitude grid.

Constructs and immediately evaluates the interpolation function from Grid.interp1d_meridional() for Grid.lat.

laplace(f)

Laplacian of an input field.

Parameters:

f (array) – 2D input field.

Returns:

Gridded Laplacian of f.

laplace_spectral(f)

laplace() with spectral in- and output fields.

mean(field, axis=None, region=None)

Area-weighted mean of the input field.

Parameters:
  • field (array) – 2D input field.

  • axis (None | int | "meridional" | "zonal") – Axis over which to compute the mean.

  • region (None | grid.GridRegion) – Region to which the mean is restricted. Construct regions with region or specify any object that implements a compatible extract method.

Returns:

Area-weighted mean.

quad(y, z=None, method='sptrapz')

Surface integral with an optional custom domain of integration.

Parameters:
  • y (array) – 2D input field.

  • z (array) – Custom domain of integration given by the region z >= 0.

  • method ("boxcount" | "sptrapz") – quadrature scheme.

Returns:

Value of the integral.

See also quad_meridional().

quad_meridional(y, z=None, method='sptrapz')

Line integral along meridians.

Parameters:
  • y (array) – 2D input field.

  • z (array) – Custom domain of integration given by the region z >= 0.

  • method ("boxcount" | "sptrapz") – quadrature scheme.

Returns:

nlon-sized array of integral values.

The sptrapz quadrature scheme is based on the trapezoidal rule adapted for sperical surface domains using the antiderivate of r * (a*φ + b) * cos(φ) to integrate over the piecewise-linear, interpolated segments between gridpoints in the meridional direction. If given, the field that defines the custom domain of integration is also linearly interpolated. This implementation is accurate but slow compared to the much simpler boxcounting. Note that no interpolation is carried out in the zonal direction (since lines of constant latitude are not great-circles, linear interpolation is non-trivial). The boundaries of the domain of integration are therefore not continuous in the zonal direction even for the sptrapz scheme.

See also quad().

property region

Region extractor with indexing syntax based on coordinates:

grid.region[lon_slice,lat_slice]

with the slices’ start and stop values in the range of 0° to 360°. Slices must not have the step parameter specified but start and stop can be left unspecified. If a single numeric value x is given instead of a slice, the value is treated like slice(x, x) but additionally verified for existence on the Grid.

Returns:

Configured barotropic.grid.GridRegion instance for a valid selection.

Note

Coordinate order is [lon,lat] (plotting convention) and inclusive on all sides (like xarray selections). This is different from indexing into the arrays directly on purpose.

Selections can cross the grid boundary in the (periodic) longitude dimension and the order of the selection ranges does not matter:

>>> grid = bt.Grid(10.)
>>> grid.region[300:50,-5:-50]
Region -50.00 to -10.00 °N (nlat = 5)
       300.00 to  50.00 °E (nlon = 12)
    ···································· 90°N
    ····································
    ····································
    ····································
    XXXXXX························XXXXXX EQ
    XXXXXX························XXXXXX
    XXXXXX························XXXXXX
    ····································
    ···································· 90°S
    0°E              180°E
of <barotropic.Grid> ...

The latitude dimension accepts "N" and "S" as special keywords for easy reference to the northern and southern hemisphere, respectively, e.g.:

>>> grid.region[:,"N"]
Region   0.00 to  90.00 °N (nlat = 10)
         0.00 to 350.00 °E (nlon = 36)
    XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX 90°N
    XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
    XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
    XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX
    XXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX EQ
    ····································
    ····································
    ····································
    ···································· 90°S
    0°E              180°E
of <barotropic.Grid> ...
property shape

Tuple of grid dimensions (nlat, nlon).

solve_diffusion(field_grid, coeff, dt, order=1)

Advance diffusion equations of various order with an implicit step.

Parameters:
  • field_grid (array) – 2D input field.

  • coeff (number) – Diffusion coefficient.

  • dt (number) – Time step in s.

  • order (int) – Order of diffusion. 1: regular diffusion, 2: hyperdiffusion, …

Returns:

Solution after one dt-sized step.

Solves:

∂f/∂t = κ·∇²f

(here order = 1) with an implicit Euler step.

solve_diffusion_spectral(field_spectral, coeff, dt, order=1)

solve_diffusion() with spectral in- and output fields.

solve_poisson(rhs_grid, op_add=0.0)

Solve Poisson(-like) equations.

Parameters:
  • rhs_grid (array) – Input field for right hand side of Poisson equation.

  • op_add (number) – Optional term with constant coefficient.

Returns:

Solution to the Poisson(-like) equation.

Discretized equation:

(∆ - op_add) f = rhs
solve_poisson_spectral(rhs_spec, op_add=0.0)

solve_poisson() with spectral in- and output fields.

to_grid(field_spectral)

Transform a spectral field into grid space.

Parameters:

field_spectral (array) – Spectral representation of input field.

Returns:

Gridded representation of input field.

to_spectral(field_grid)

Transform a gridded field into spectral space.

Parameters:

field_grid (array) – Gridded representation of input field.

Returns:

Spectral representation of input field.

vorticity(u, v)

Vorticity from vector components.

Parameters:
  • u (array) – Zonal component of input field (gridded).

  • v (array) – Meridional component of input field (gridded).

Returns:

Gridded vorticity field.

vorticity_divergence(u, v)

Vorticity and divergence from vector components.

Parameters:
  • u (array) – Zonal component of input field (gridded).

  • v (array) – Meridional component of input field (gridded).

Returns:

Tuple of gridded (vorticity, divergence) fields.

vorticity_divergence_spectral(u, v)

vorticity_divergence() with spectral output fields.

vorticity_spectral(u, v)

vorticity() with spectral output field.

wind(vorticity, divergence)

Wind components from vorticity and divergence fields.

Parameters:
  • vorticity (array) – Vorticity field (gridded).

  • divergence (array) – Divergence field (gridded).

Returns:

Tuple of (zonal, meridional) wind compontents (gridded).

zonalize(field, levels=None, interpolate=True, quad='sptrapz')

Zonalize the field with equivalent latitude coordinates.

Parameters:
  • field (array) – 2D input field.

  • levels (None | array | int) – Contours generated for the zonalization. If an array of contour-values is given, these are used directly. If value is an integer, this number of contours is sampled between the maximum and minimum values in the input field automatically. By default, the number of levels is set equal to a multiple of the number of latitudes resolved by the grid (2 for boxcounting quadrature, 1 otherwise).

  • interpolate (bool) – Interpolate output onto regular latitudes of grid. Without interpolation, values are returned on the equivalent latitudes that arise in the computation. These may be irregular and unordered.

  • quad (str) – Quadrature rule used in the surface integrals of the computation. See quad().

Returns:

If interpolate is True, return the contour values interpolated onto the regular latitudes as an array. Otherwise, return a tuple of contour values and associated equivalent latitudes (both arrays).

Implements the zonalization procedure of Nakamura and Zhu (2010).

Region Selection

class barotropic.grid.GridRegion(grid, lat_indices, lon_indices)

Bases: object

Region selection for convenient data extraction.

Note

Use barotropic.Grid.region to create instances of this class.

Parameters:
  • grid (Grid) – Lon-lat grid instance that the selection is tied to.

  • lat_indices (sequence) – Latitude dimension indices included in the selection.

  • lon_indices (sequence) – Longitude dimension indices included in the selection.

extract(*fields)

Extract data from the selected region.

Parameters:

fields (array) – Field or fields from which to extract the selected region. The shape of each field must match that of the Grid associated with the GridRegion instance.

Returns:

A single array (one-argument call) or a tuple of arrays (multi-argument call, order preserved) containing the data selection.

>>> grid = bt.Grid(10.)
>>> sel = grid.region[45:90,-20:30]
>>> sel.extract(grid.lon2)
array([[50., 60., 70., 80., 90.],
       [50., 60., 70., 80., 90.],
       [50., 60., 70., 80., 90.],
       [50., 60., 70., 80., 90.],
       [50., 60., 70., 80., 90.],
       [50., 60., 70., 80., 90.]])

If the region crosses the periodic boundary in the longitude dimension, the two parts of the fields from the original array are reassembled into a field connected across the boundary:

>>> grid.region[340:10,60:70].extract(grid.lon2)
array([[340., 350.,   0.,  10.],
       [340., 350.,   0.,  10.]])

Use GridRegion.lon_mono to obtain continuous longitude coordinates in such situations.

Always returns fields of full dimensionality, even for single-value selections:

>>> grid.region[0,50].extract(grid.fcor2)
array([[0.00011172]])
property gridpoint_area

See barotropic.Grid.gridpoint_area.

property lat

Latitude coordinate array for the region.

property lon

Longitude coordinate array for the region.

Note

The array will not be monotonic if the region crosses the 0° meridian. In case a monotonic longitude coordinate is required, e.g. for plotting, use lon_mono.

property lon_mono

Longitude coordinate array for the region, monotonic across 0°.

Longitudes west of 0° (the western edge of a Grid) are reduced by 360° to obtain monotonic values:

>>> sel = bt.Grid(10.).region[320:30,:]
>>> sel.lon
array([320., 330., 340., 350.,   0.,  10.,  20.,  30.])
>>> sel.lon_mono
array([-40., -30., -20., -10.,   0.,  10.,  20.,  30.])
property mask

Boolean array that is True in the selected region.

mean(field)

Area-weighted mean in the selected region.

Wraps barotropic.Grid.mean().

property shape

Shape of extracted 2D fields.