Representation of Grid Data¶
The grid_data
module provides representations of data on
uniform grids as well as for data on refined grid hirachies. Standard
arithmetic operations are supported for those data grids, further methods
to interpolate and resample. The number of dimensions is arbitrary.
Rudimentary vector and matrix oprations are also supported, using
Vectors of data grids (instead of grids of vectors).
- The important classes defined here are
- class postcactus.grid_data.RegGeom(shape, origin, delta=None, x1=None, reflevel=- 1, component=- 1, nghost=None, time=None, iteration=None)¶
Describes the geometry of a regular rectangular dataset, as well as information needed if part of refined grid hierachy, namely component number and refinement level. Also stores the number of ghost zones, which is however not used anywhere in this class.
- Parameters
shape (1d numpy arrary or list of int.) – Number of points in each dimension.
origin (1d numpy array or list of float.) – Position of cell center with lowest coordinate.
delta (1d numpy array or list of float.) – If not None, specifies grid spacing, else grid spacing is computed from x0, x1, and shape.
x1 (1d numpy array or list of float.) – If grid spacing is None, this specifies the position of the cell center with largest coordinates.
reflevel (int) – Refinement level if this belongs to a hierachy, else -1.
component (int) – Component number if this belongs to a hierachy, else -1.
nghost (1d numpy arrary or list of int.) – Number of ghost zones (default=0)
time (float or None) – Time if that makes sense, else None.
iteration (float or None) – Iteration if that makes sense, else None.
- component()¶
Component number if this grid belongs to a hierachy, else -1.
- Returns
Component number.
- Return type
int
- contains(pos)¶
Test if a coordinate is contained in the grid. The size of the grid cells is taken into account, resulting in a cube larger by dx/2 on each side compared to the one given by x0, x1.
- Parameters
pos (1d numpy array or list of float.) – Coordinate to test.
- Returns
If pos is contained.
- Return type
bool
- coords1d()¶
Get coordinates of the grid points as 1D arrays.
- Returns
The coordinate array of each dimension.
- Return type
list of 1d numpy arrays
- coords2d()¶
Get coordinates of the grid points as 2D arrays with the same shape as the grid. Useful for arithmetic computations involving both data and coordinates.
- Returns
The coordinate array of each dimension.
- Return type
list of numpy arrays with same shape as grid.
- dim_ext()¶
- Returns
Whether the extend in each dimension is larger than one gridpoint.
- Return type
array of bools
- dv()¶
- Returns
Volume of a grid cell.
- Return type
float
- dx()¶
- Returns
Grid spacing.
- Return type
1d numpy array of float
- flatten()¶
Remove dimensions which are only one gridpoint across
- ind2pos(i)¶
Compute coordinate corresponding to a grid point.
- Parameters
i (1d array or list of int.) – Grid index.
- Returns
The coordinate of grid point i.
- Return type
1d numpy array of float.
- nfdims()¶
- Returns
The number of dimensions with size larger than one gridpoint.
- Return type
int
- nghost()¶
- Returns
Number of ghost zones in each dimension.
- Return type
1d numpy array of int
- num_dims()¶
- Returns
the number of dimensions.
- Return type
int
- pos2ind(p)¶
Find the grid point nearest to a given position.
- Parameters
p (1d numpy array or list of float) – Coordinate.
- Returns
grid index of nearest point.
- Return type
tuple of int
- reflevel()¶
Refinement level if this grid belongs to a hierachy, else -1.
- Returns
Level.
- Return type
int
- scale_coords(scale)¶
Rescale all coordinates.
- Parameters
scale (float or 1d numpy array) – Factor to scale by.
- shape()¶
- Returns
Number of grid points in each dimension.
- Return type
1d numpy array of int
- shift_coords(shift)¶
Shift all coordinates.
- Parameters
shift (float or 1d numpy array) – vector added to the origin.
- volume()¶
- Returns
Volume of the whole grid.
- Return type
float
- x0()¶
Position of the corner with minimum coordinates (cell center).
- Returns
Corner coordinate.
- Return type
1d numpy array of float
- x1()¶
Position of the corner with maximum coordinates (cell center).
- Returns
Corner coordinate.
- Return type
1d numpy array of float
- class postcactus.grid_data.RegData(origin, delta, data, outside_value=0, reflevel=- 1, component=- 1, nghost=None, time=None, iteration=None)¶
Represents a rectangular data grid with coordinates, supporting common arithmetic operations. The standard mathematical unary functions, e.g. sin(), are implemented as member functions. Supports interpolation and resampling. The objects can also be used as a function, mapping coordinate to interpolated data. Supports numercial differentiation.
- Variables
data – The actual data (numpy array).
- Parameters
origin (1d numpy array or list of float.) – Position of cell center with lowest coordinate.
delta (1d numpy array or list of float.) – If not None, specifies grid spacing, else grid spacing is computed from x0, x1, and shape.
data (A numpy array.) – The data.
outside_value (float) – Value to use when interpolating outside the grid.
reflevel (int) – Refinement level if this belongs to a hierachy, else -1.
component (int) – Component number if this belongs to a hierachy, else -1.
nghost (1d numpy arrary or list of int.) – Number of ghost zones (default=0)
time (float or None) – Time if that makes sense, else None.
iteration (float or None) – Iteration if that makes sense, else None.
- apply_binary(a, b, f)¶
Apply a binary function to two regular data sets.
- apply_unary(op)¶
Apply a unary function to the data.
- Parameters
op (function operating on a numpy array) – unary function.
- Returns
result.
- Return type
RegData
instance.
- coords()¶
Get coordinates as regular datasets. Useful for arithmetic operations involving data and coordinates.
- Returns
list of coordinates for each dimension
- Return type
list of
RegData
instances
- diff(dim, order=2)¶
Computes the partial derivative along a given dimension. Uses either a 3-point central difference stencil for 2nd order accuracy, or a five point central stencil for 4th order. At the boundaries, 1 point (2 for 4th order) is computed to first order accuracy using one-sided derivatives. The array size in the dimension dim needs to be at least the stencil size.
- Parameters
dim (int) – Dimension of partial derivative.
order (int) – Order of accuracy (2 or 4).
- Returns
The derivative.
- Return type
RegData
instance
- flatten()¶
Remove dimensions which are only one gridpoint large.
- grad(order=2)¶
Compute the gradient. See diff for details.
- Parameters
order (int) – Order of accuracy (2 or 4).
- Returns
The gradient vector.
- Return type
Vec
instance.
- histogram(weights=None, vmin=None, vmax=None, nbins=400)¶
1D Histogram of the data. :param weights: the weight for each cell. Default is one. :type weights: RegData or numpy array of same shape or None. :param vmin: Lower bound of data to consider. Default is data range. :type vmin: float or None :param vmax: Upper bound of data to consider. Default is data range. :type vmax: float or None :param nbins: Number of bins to create. :type nbins: integer > 1
- Returns
the positions of the data bins and the distribution.
- Return type
tuple of two 1D numpy arrays.
- integral()¶
Compute the integral over the whole volume of the grid.
- Returns
The integral.
- Return type
float (or complex if data is complex)
- interp_mlinear(x)¶
Perform a multilinear interpolation to the coordinate x. If x is outside the grid, return the outside value.
- Parameters
x (1d numpy array or list of float) – Coordinate to interpolate to.
- Returns
Value at nearest grid point or outside_value.
- Return type
same as data.
- interp_nearest(x)¶
Zeroth order interpolation. For positions outside the grid, it uses the outside value specified when creating the object.
- Parameters
x (1d numpy array or list of float) – Coordinate to interpolate to.
- Returns
Value at nearest grid point or outside_value.
- Return type
same as data.
- interp_nearest_index(i)¶
Zeroth order interpolation. For positions outside the grid, it uses the outside value specified when creating the object.
- Parameters
i (1d numpy array or list of float) – Position in terms of real valued grid index (not a coordinate!)
- Returns
Value at nearest grid point or outside_value.
- Return type
same as data.
- max(axis=None)¶
- Parameters
axis (int or None) – optional, compute max along single axis
- Returns
Maximum value over all points on the grid if axis=None, else RegData containing max along given axis.
- Return type
float (or complex if data is complex)
- maximum(b)¶
- mean(axis=None)¶
- Parameters
axis (int or None) – optional, compute mean along single axis
- Returns
Mean value over all points on the grid if axis=None, else RegData containing mean along given axis.
- Return type
float (or complex if data is complex)
- min(axis=None)¶
- Parameters
axis (int or None) – optional, compute min along single axis
- Returns
Minimum value over all points on the grid if axis=None, else RegData containing min along given axis.
- Return type
float (or complex if data is complex)
- minimum(b)¶
- percentiles(fractions, weights=None, relative=True, vmin=None, vmax=None, nbins=400)¶
Find values for which a given fraction(s) of the data is smaller.
Optionally, the cells can have an optional weight, and absolute counts can be used insted of fraction.
- Parameters
fractions (list or array of floats) – list of fraction/absolute values
weights (RegData or numpy array of same shape or None.) – the weight for each cell. Default is one.
relative (bool) – whether fractions refer to relative or absolute count.
vmin (float or None) – Lower bound of data to consider. Default is data range.
vmax (float or None) – Upper bound of data to consider. Default is data range.
nbins (integer > 1) – Number of bins to create.
- Returns
data values corresponding to the given fractions.
- Return type
1D numpy array
- reflect(dim, parity=1)¶
Fill points assuming reflection symmetry along given axis.
- sample(geom, order=0, output=None, mode='constant')¶
Resample the data to a regular grid geometry specified. Coordinates outside the grid are treated according to the mode parameter, where ‘constant’ means use the outside value, ‘nearest’ means zero-order extrapolation. Can create new data for results or overwrite a given dataset (for efficiency).
- Parameters
order (int) – Interpolation order.
output (None or numpy array) – If not None, where to write results to.
mode (string) – How to treat outside positions.
- Returns
Interpolated data.
- Return type
RegData
or None
- sample_generic(coords, order=0, mode='constant', output=None)¶
Interpolate the data to arbitrary (irregular) coordinates. Coordinates outside the grid are treated according to the mode parameter, where ‘constant’ means use the outside value, ‘nearest’ means zero-order extrapolation. Do not use this method repeatedly, instead use the spline() method which is then more efficient.
- Parameters
coords (list of numpy arrays) – List of coordinate arrays for each dimension.
order (Int) – Order of the spline, must be in [2,5]. Default is 3.
mode ('constant', 'nearest', 'reflect' or 'wrap') – How to treat values outside the data range.
- Returns
Interpolated data.
- Return type
numpy array with same dimension as the coordinate arrays
- sample_intersect(other, order=0)¶
Replace data in the intersection with another object by resampled data from the other object. Ghost zones in the other object are excluded from the intersection, but used during the interpolation.
- spline(order=3, mode='constant')¶
Returns interpolating spline that represents the data as a function. Always use this instead of calling sample_generic() repeatedly, because computing the spline coefficients is expensive.
- Parameters
order (Int) – Order of the spline, must be in [2,5]. Default is 3.
mode ('constant', 'nearest', 'reflect' or 'wrap') – How to treat values outside the data range.
- class postcactus.grid_data.CompData(alldat, outside_value=0)¶
Composite data consisting of one or more regular datasets with different grid spacings, i.e. a mesh refinement hirachy. The grid spacings should differ by powers of two. Origins of the components are shifted relative to each other only by multiples of the finest spacing. Basic arithmetic operations are defined for this class, as well as interpolation and resampling. This class can be iterated over to get all the regular datasets, ordered by refinement level and componen number.
- Parameters
alldat (list of
RegData
instances.) – list of regular datasets.outside_value (float) – Value to use when interpolating outside the grid.
- apply_binary(a, b, op)¶
Apply a binary function to two data sets.
- apply_unary(f)¶
Apply a unary function to the data.
- Parameters
op (function operating on a numpy array) – unary function.
- Returns
result.
- Return type
- coords()¶
Get coordinates as composite datasets with the same structure. Useful for arithmetic operations involving data and coordinates.
- Returns
list of coordinates for each dimension
- Return type
list of
CompData
instances
- diff(dim, order=2)¶
Computes the partial derivative along a given dimension. Uses either a 3-point central difference stencil for 2nd order accuracy, or a five point central stencil for 4th order. At the boundaries, 1 point (2 for 4th order) is computed to first order accuracy using one-sided derivatives. The array size in the dimension dim needs to be at least the stencil size.
Note
The derivative is computed for each refinement level and component independently, if there are not enough ghost zones the order of accuracy drops to one at each component boundary.
- Parameters
dim (int) – Dimension of partial derivative.
order (int) – Order of accuracy (2 or 4).
- Returns
The derivative.
- Return type
CompData
instance with same structure.
- dim_ext()¶
- Returns
Whether the extend in each dimension is larger than one gridpoint.
- Return type
array of bools
- dx_coarse()¶
Grid spacing of coarsest level
- dx_finest()¶
Grid spacing of finest level
- grad(order=2)¶
Compute the gradient. See diff for details.
- Parameters
order (int) – Order of accuracy (2 or 4).
- Returns
The gradient vector.
- Return type
Vec
instance.
- integral()¶
Compute the integral over the whole volume of the grid. Note this is currently only implemented for the case of one refinement level with only one component, otherwise raises an exception.
- Returns
The integral.
- Return type
float (or complex if data is complex)
- interp_mlinear(x)¶
Perform a multilinear interpolation to the coordinate x, using the finest grid containing x. If x is outside the grid, return the outside value.
- Parameters
x (1d numpy array or list of float) – Coordinate to interpolate to.
- Returns
Value at nearest grid point or outside_value.
- Return type
same as data.
- interp_nearest(x)¶
Zeroth order interpolation, using the finest available grid covering the given coordinate. For positions outside the grid, it uses the outside value specified when creating the object.
- Parameters
x (1d numpy array or list of float) – Coordinate to interpolate to.
- Returns
Value at nearest grid point or outside_value.
- Return type
same as data.
- max()¶
- Returns
The maximum of the data.
- Return type
same as data.
- maximum(b)¶
- min()¶
- Returns
The minimum of the data.
- Return type
same as data.
- minimum(b)¶
- num_dims()¶
- Returns
Number of dimensions.
- Return type
int
- sample(geom, order=0, dest=None, adjust_spacing=True)¶
Resamples to a regular grid using.
- sample_alt(geom, order=0, dest=None, adjust_spacing=True)¶
Resamples to a regular grid with alternative method.
- scale_coords(scale)¶
Rescale all coordinates.
- Parameters
scale (float or 1d numpy array) – Factor to scale by.
- x0()¶
The bounding box corner with minimum coordinates.
- Returns
Corner coordinate.
- Return type
1d numpy array of float
- x1()¶
The bounding box corner with maximum coordinates.
- Returns
Corner coordinate.
- Return type
1d numpy array of float
- postcactus.grid_data.sample(func, x0, x1, shape)¶
Create a regular dataset by sampling a scalar function on a regular grid.
- Parameters
func (Anything accepting 1d numpy array as input and returns scalar.) – The function to sample.
x0 (1d numpy array or list of float) – Minimum corner of regular sample grid.
x0 – Maximum corner of regular sample grid.
shape (1d numpy array or list of int) – Number of sample points in each dimension.
- Returns
Sampled data.
- Return type
- postcactus.grid_data.merge_comp_data_1d(comps)¶
Convert 1D refined grid hierachy into irregularly (piecewise regularly) spaced 1D data, using the finest level available at each interval.
- Parameters
comps (
CompData
) – Grid hierachy.- Returns
coordinates and data values.
- Return type
two 1d numpy arrays