5. gridData.core — Core functionality for storing n-D grids

The core module contains classes and functions that are independent of the grid data format. In particular this module contains the Grid class that acts as a universal constructor for specific formats:

g = Grid(**kwargs)           # construct
g.export(filename, format)   # export to the desired format

Some formats can also be read:

g = Grid()                   # make an empty Grid
g.load(filename)             # populate with data from filename

5.1. Classes and functions

class gridData.core.Grid(grid=None, edges=None, origin=None, delta=None, metadata={}, interpolation_spline_order=3, file_format=None)[source]

Class to manage a multidimensional grid object.

The export(format=’dx’) method always exports a 3D object, the rest should work for an array of any dimension.

The grid (Grid.grid) can be manipulated as a standard numpy array.

The attribute Grid.metadata holds a user-defined dictionary that can be used to annotate the data. It is saved with save().

Create a Grid object from data.

From a numpy.histogramdd():

grid,edges = numpy.histogramdd(...)
g = Grid(grid,edges=edges)

From an arbitrary grid:

g = Grid(grid,origin=origin,delta=delta)

From a saved file:

g = Grid(filename)

or

g = Grid()
g.load(filename)
Arguments
grid

histogram or density, defined on numpy nD array, or filename

edges

list of arrays, the lower and upper bin edges along the axes (both are output by numpy.histogramdd())

origin

cartesian coordinates of the center of grid[0,0,…,0]

delta

Either n x n array containing the cell lengths in each dimension, or n x 1 array for rectangular arrays.

metadata

a user defined dictionary of arbitrary values associated with the density; the class does not touch metadata[] but stores it with save()

interpolation_spline_order

order of interpolation function for resampling; cubic splines = 3 [3]

file_format

file format; only necessary when grid is a filename (see Grid.load()); default is None and the file format is autodetected.

Changed in version 0.5.0: New file_format keyword argument.

centers()[source]

Returns the coordinates of the centers of all grid cells as an iterator.

check_compatible(other)[source]

Check if other can be used in an arithmetic operation.

  1. other is a scalar

  2. other is a grid defined on the same edges

Raises

TypeError if not compatible.

export(filename, file_format=None, type=None, typequote='"')[source]

export density to file using the given format.

The format can also be deduced from the suffix of the filename though the format keyword takes precedence.

The default format for export() is ‘dx’. Use ‘dx’ for visualization.

Implemented formats:

dx

OpenDX

pickle

pickle (use Grid.load() to restore); Grid.save() is simpler than export(format='python').

Parameters
  • filename (str) – name of the output file

  • file_format ({'dx', 'pickle', None} (optional)) – output file format, the default is “dx”

  • type (str (optional)) –

    for DX, set the output DX array type, e.g., “double” or “float”. By default (None), the DX type is determined from the numpy dtype of the array of the grid (and this will typically result in “double”).

    New in version 0.4.0.

  • typequote (str (optional)) –

    For DX, set the character used to quote the type string; by default this is a double-quote character, ‘”’. Custom parsers like the one from NAMD-GridForces (backend for MDFF) expect no quotes, and typequote=’’ may be used to appease them.

    New in version 0.5.0.

interpolated

B-spline function over the data grid(x,y,z).

The interpolated() function allows one to obtain data values for any values of the coordinates:

interpolated([x1,x2,...],[y1,y2,...],[z1,z2,...]) -> F[x1,y1,z1],F[x2,y2,z2],...

The interpolation order is set in Grid.interpolation_spline_order.

The interpolated function is computed once and is cached for better performance. Whenever interpolation_spline_order is modified, Grid.interpolated() is recomputed.

The value for unknown data is set in Grid.interpolation_cval (TODO: also recompute when interpolation_cval value is changed.)

Example

Example usage for resampling:

XX, YY, ZZ = numpy.mgrid[40:75:0.5, 96:150:0.5, 20:50:0.5]
FF = interpolated(XX, YY, ZZ)

Note

Values are interpolated with a spline function. It is possible that the spline will generate values that would not normally appear in the data. For example, a density is non-negative but a cubic spline interpolation can generate negative values, especially at the boundary between 0 and high values.

interpolation_spline_order

Order of the B-spline interpolation of the data.

3 = cubic; 4 & 5 are also supported

Only choose values that are acceptable to scipy.ndimage.spline_filter()!

See also

interpolated

load(filename, file_format=None)[source]

Load saved (pickled or dx) grid and edges from <filename>.pickle

Grid.load(<filename>.pickle) Grid.load(<filename>.dx)

The load() method calls the class’s constructor method and completely resets all values, based on the loaded data.

resample(edges)[source]

Resample data to a new grid with edges edges.

This method creates a new grid with the data from the current grid resampled to a regular grid specified by edges. The order of the interpolation is set by Grid.interpolation_spline_order: change the value before calling resample().

Parameters

edges (tuple of arrays or Grid) – edges of the new grid or a Grid instance that provides Grid.edges

Returns

a new Grid with the data interpolated over the new grid cells

Return type

Grid

Examples

Providing edges (a tuple of three arrays, indicating the boundaries of each grid cell):

g = grid.resample(edges)

As a convenience, one can also supply another Grid as the argument for this method

g = grid.resample(othergrid)

and the edges are taken from Grid.edges.

resample_factor(factor)[source]

Resample to a new regular grid.

Parameters

factor (float) – The number of grid cells are scaled with factor in each dimension, i.e., factor * N_i cells along each dimension i.

Returns

Return type

Grid

See also

resample()

save(filename)[source]

Save a grid object to <filename>.pickle

Internally, this calls Grid.export(filename, format="python"). A grid can be regenerated from the saved data with

g = Grid(filename="grid.pickle")

Note

The pickle format depends on the Python version and therefore it is not guaranteed that a grid saved with, say, Python 2.7 can also be read with Python 3.5. The OpenDX format is a better alternative for portability.

gridData.core.ndmeshgrid(*arrs)[source]

Return a mesh grid for N dimensions.

The input are N arrays, each of which contains the values along one axis of the coordinate system. The arrays do not have to have the same number of entries. The function returns arrays that can be fed into numpy functions so that they produce values for all points spanned by the axes arrs.

Original from http://stackoverflow.com/questions/1827489/numpy-meshgrid-in-3d and fixed.