gridded.pysgrid package

Submodules

gridded.pysgrid.lookup module

Created on Mar 24, 2015

@author: ayan

gridded.pysgrid.processing_2d module

Created on Apr 2, 2015

@author: ayan

gridded.pysgrid.processing_2d.avg_to_cell_center(data_array, avg_dim)[source]

Given a two-dimensional numpy.array, average adjacent row values (avg_dim=1) or adjacent column values (avg_dim=0) to the grid cell center.

Parameters
  • data_array (numpy.array) – 2-dimensional data

  • avg_dim (int) – integer specify array axis to be averaged

Returns

averages

Return type

numpy.array

gridded.pysgrid.processing_2d.rotate_vectors(x_arr, y_arr, angle_arr)[source]

Given x and y vectors in a projected coordinate system, rotate them by angles into a different coordinate system.

All arrays must have the same dimensions.

Parameters
  • x_arr (numpy.array) – array of x-directed vectors

  • y_arr (numpy.array) – array of y-directed vectors

  • angle_arr (numpy.array) – array of angles in radians

Returns

x and y arrays of rotated vectors

Return type

tuple

gridded.pysgrid.processing_2d.vector_sum(x_arr, y_arr)[source]

Calculate the vector sum of arrays of x and y vectors.

Parameters
  • x_arr (numpy.array) – array of x-directed vectors

  • y_arr (numpy.array) – array of y-directed vectors

Returns

array of vector sums

Return type

numpy.array

gridded.pysgrid.read_netcdf module

Created on Mar 19, 2015

@author: ayan

class gridded.pysgrid.read_netcdf.NetCDFDataset(nc)[source]

Bases: object

find_coordinates_by_location(location_str, topology_dim)[source]

Find a grid coordinates variables with a location attribute equal to location_str. This method can be used to infer edge, face, or volume coordinates from the location attribute of a variable.

Location is a required attribute per SGRID conventions.

Parameters
  • location_str (str) – the location value to search for

  • topology_dim (int) – the topology dimension of the grid

find_node_coordinates(node_dimensions)[source]

Find the variables for the grid cell vertices.

find_variables_by_attr(**kwargs)[source]
sgrid_compliant_file()[source]

Determine whether a dataset is SGRID compliant.

Returns

True if dataset is compliant, raise an exception if it is not

Return type

bool

gridded.pysgrid.read_netcdf.find_grid_topology_var(nc)[source]

Get the variable from a netCDF dataset that have a cf_role `attribute of ‘grid_topology’ and `topology_dimension of 2.

Params nc

netCDF dataset

Returns

variable name that contain grid topologies

Return type

string

gridded.pysgrid.read_netcdf.parse_axes(axes_attr)[source]
gridded.pysgrid.read_netcdf.parse_padding(padding_str, mesh_topology_var)[source]

Use regex expressions to break apart an attribute string containing padding types for each variable with a cf_role of ‘grid_topology’.

Padding information is returned within a named tuple for each node dimension of an edge, face, or vertical dimension. The named tuples have the following attributes: mesh_topology_var, dim_name, dim_var, and padding. Padding information is returned as a list of these named tuples.

Parameters

padding_str (str) – string containing padding types from a netCDF attribute.

Returns

named tuples with padding information.

Return type

list.

gridded.pysgrid.read_netcdf.parse_vector_axis(variable_standard_name)[source]

gridded.pysgrid.sgrid module

Created on Apr 20, 2015

@author: ayan

class gridded.pysgrid.sgrid.SGrid(node_lon=None, node_lat=None, node_mask=None, center_lon=None, center_lat=None, center_mask=None, edge1_lon=None, edge1_lat=None, edge1_mask=None, edge2_lon=None, edge2_lat=None, edge2_mask=None, edges=None, node_padding=None, edge1_padding=None, edge2_padding=None, grid_topology_var=None, variables=None, grid_variables=None, dimensions=None, node_dimensions=None, node_coordinates=None, edge1_coordinates=None, edge2_coordinates=None, angles=None, edge1_dimensions=None, edge2_dimensions=None, faces=None, face_padding=None, face_coordinates=None, face_dimensions=None, vertical_padding=None, vertical_dimensions=None, tree=None, use_masked_boundary=False, grid_topology=None, *args, **kwargs)[source]

Bases: object

all_padding()[source]
apply_padding_to_idxs(idxs, padding=('none', 'none'))[source]

Given a list of indexes, increment each dimension to compensate for padding. Input indexes are assumed to be cell indexes

build_celltree(use_mask=True)[source]

Builds the celltree across the grid defined by nodes (self.node_lon, self.node_lat) If center masking is provided in self.center_mask, it will remove masked cells, and take precedence over any node masking for celltree insertion.

If node masking is provided in self.node_mask and self.center_mask is not provided, it will remove masked nodes from the grid, which also removes all adjacent cells

Parameters

use_mask – If False, ignores all masks and builds the celltree over the raw arrays. Does nothing if self.node_mask or self.center_mask are not present

build_kdtree(grid='node')[source]

Builds the kdtree for the specified grid

property center_padding
property centers
property edge1_padding
property edge2_padding
geo_to_logical(points, indices=None, _memo=False, _copy=False, _hash=None)[source]

Given a list of lon/lat points, converts them to l/m coordinates in logical cell space.

get_all_edge_padding()[source]
get_all_face_padding()[source]
get_efficient_slice(points=None, indices=None, location=None, _memo=False, _copy=False, _hash=None)[source]

Computes the minimum 2D slice that captures all the provided points/indices within. :param points: Nx2 array of longitude/latitude. (Optional) :param indices: Nx2 array of logical cell indices (Optional, but required if points omitted) :param location: ‘center’, ‘edge1’, ‘edge2’,’node’

get_padding_by_location(location)[source]
get_padding_slices(padding=('none', 'none'))[source]

Given a pair of padding types, return a numpy slice object you can use directly on data or lon/lat variables

get_variable_at_index(var, index)[source]

Given a list of 2D indices, return the value of var at each index

get_variable_by_index(var, index)[source]

index = index arr of quads (maskedarray only) var = ndarray/ma.array returns ndarray/ma.array

ordering is idx, idx+[0,1], idx+[1,1], idx+[1,0] masked values from var remain masked

Function to get the node values of a given face index. Emulates the ‘self.grid.nodes[self.grid.nodes.faces[index]]’ paradigm of unstructured grids.

infer_location(variable)[source]

Assuming default is psi grid, check variable dimensions to determine which grid it is on.

property info

Summary of information about the grid

This needs to be implimented – see UGrid for example

interpolate(points, variable, location=None, fill_value=0, indices=None, alphas=None, padding=None, slices=None, unmask=False, _memo=False, _hash=None, _copy=False)

Interpolates a variable on one of the grids to an array of points. :param points: Nx2 Array of lon/lat coordinates to be interpolated to.

Parameters
  • variable – Array-like of values to associate at location on grid (node, center, edge1, edge2). This may be more than a 2 dimensional array, but you must pass ‘slices’ kwarg with appropriate slice collection to reduce it to 2 dimensions.

  • location – One of (‘node’, ‘center’, ‘edge1’, ‘edge2’, ‘face’). ‘edge1’ is conventionally associated with the ‘vertical’ edges and likewise ‘edge2’ with the ‘horizontal’. Determines type of interpolation, see below for details

  • fill_value – If masked values are encountered in interpolation, this value takes the place of the masked value

  • indices – If computed already, array of Nx2 cell indices can be passed in to increase speed.

  • alphas – If computed already, array of alphas can be passed in to increase speed.

  • unmask – If true, unmask results, using fill value for masked values.

Depending on the location specified, different interpolation will be used.

For ‘center’, no interpolation

For ‘edge1’ or ‘edge2’, interpolation is linear, edge to edge across the cell

For ‘node’, interpolation is bilinear from the four nodes of each cell

The variable specified may be any array-like. - With a numpy array:

sgrid.interpolate_var_to_points(points, sgrid.u[time_idx, depth_idx]) - With a raw netCDF Variable:

sgrid.interpolate_var_to_points(points, nc.variables[‘u’], slices=[time_idx, depth_idx])

If you have pre-computed information, you can pass it in to avoid unnecessary computation and increase performance.

  • ind = # precomputed indices of points. This may be a masked array

  • alphas = # precomputed alphas (useful if interpolating to the same points frequently)

sgrid.interpolate_var_to_points(points, sgrid.u, indices=ind, alphas=alphas, slices=[time_idx, depth_idx])

interpolate_var_to_points(points, variable, location=None, fill_value=0, indices=None, alphas=None, padding=None, slices=None, unmask=False, _memo=False, _hash=None, _copy=False)[source]

Interpolates a variable on one of the grids to an array of points. :param points: Nx2 Array of lon/lat coordinates to be interpolated to.

Parameters
  • variable – Array-like of values to associate at location on grid (node, center, edge1, edge2). This may be more than a 2 dimensional array, but you must pass ‘slices’ kwarg with appropriate slice collection to reduce it to 2 dimensions.

  • location – One of (‘node’, ‘center’, ‘edge1’, ‘edge2’, ‘face’). ‘edge1’ is conventionally associated with the ‘vertical’ edges and likewise ‘edge2’ with the ‘horizontal’. Determines type of interpolation, see below for details

  • fill_value – If masked values are encountered in interpolation, this value takes the place of the masked value

  • indices – If computed already, array of Nx2 cell indices can be passed in to increase speed.

  • alphas – If computed already, array of alphas can be passed in to increase speed.

  • unmask – If true, unmask results, using fill value for masked values.

Depending on the location specified, different interpolation will be used.

For ‘center’, no interpolation

For ‘edge1’ or ‘edge2’, interpolation is linear, edge to edge across the cell

For ‘node’, interpolation is bilinear from the four nodes of each cell

The variable specified may be any array-like. - With a numpy array:

sgrid.interpolate_var_to_points(points, sgrid.u[time_idx, depth_idx]) - With a raw netCDF Variable:

sgrid.interpolate_var_to_points(points, nc.variables[‘u’], slices=[time_idx, depth_idx])

If you have pre-computed information, you can pass it in to avoid unnecessary computation and increase performance.

  • ind = # precomputed indices of points. This may be a masked array

  • alphas = # precomputed alphas (useful if interpolating to the same points frequently)

sgrid.interpolate_var_to_points(points, sgrid.u, indices=ind, alphas=alphas, slices=[time_idx, depth_idx])

classmethod load_grid(nc)[source]
locate_faces(points, _memo=False, _copy=False, _hash=None, use_mask=True)[source]

Given a list of points, returns a list of x, y indices of the cell that contains each respective point

Points that are not on the node grid will have an index of -1

If a single point is passed in, a single index will be returned. If a sequence of points is passed in an array of indexes will be returned.

Parameters
  • points (array-like containing one or more points: shape (2,) for one point, shape (N, 2) for more than one point.) – The points that you want to locate – (lon, lat). If the shape of point is 1D, function will return a scalar index. If it is 2D, it will return a 1D array of indices.

  • grid (Name of the grid ('node', 'center', 'edge1', 'edge2)) – The grid on which you want to locate the points

This version utilizes the CellTree data structure.

locate_nearest(points, grid, _memo=False, _copy=False, _hash=None)[source]
nearest_var_to_points(points, variable, indices=None, grid=None, alphas=None, mask=None, slices=None, _memo=False, slice_grid=True, _hash=None, _copy=False)[source]
property node_padding
property nodes
property non_grid_variables
padding_slices = {'both': (1, -1), 'high': (None, 1), 'low': (1, None), 'none': (None, None)}
save_as_netcdf(filepath)[source]

save the grid as a netcdf file

Parameters

filepath – path to the file to be created and saved to

topology_dimension = 2
static x_to_l(x, y, a, b)[source]

Params: x: x coordinate of point y: y coordinate of point a: x coefficients b: y coefficients

Returns: (l,m) - coordinate in logical space to use for interpolation

Eqns: m = (-bb +- sqrt(bb^2 - 4*aa*cc))/(2*aa) l = (l-a1 - a3*m)/(a2 + a4*m)

class gridded.pysgrid.sgrid.SGridAttributes(nc, topology_dim, topology_variable)[source]

Bases: object

Class containing methods to help with getting the attributes for either SGrid.

get_angles()[source]
get_attr_coordinates(attr_name)[source]
get_attr_dimension(attr_name)[source]
get_cell_center_lat_lon()[source]
get_cell_edge1_lat_lon()[source]
get_cell_edge2_lat_lon()[source]
get_cell_node_lat_lon()[source]
get_dimensions()[source]
get_masks(node, center, edge1, edge2)[source]
get_node_coordinates()[source]
get_topology_var()[source]
get_variable_attributes(sgrid)[source]
gridded.pysgrid.sgrid.load_grid(nc)[source]

Get a SGrid object from a netCDF4.Dataset or file/URL.

Parameters

nc (str or netCDF4.Dataset) – a netCDF4 Dataset or URL/filepath to the netCDF file

Returns

SGrid object

Return type

sgrid.SGrid

gridded.pysgrid.utils module

Created on Mar 23, 2015

@author: ayan

class gridded.pysgrid.utils.GridPadding(mesh_topology_var, face_dim, node_dim, padding)

Bases: tuple

face_dim

Alias for field number 1

mesh_topology_var

Alias for field number 0

node_dim

Alias for field number 2

padding

Alias for field number 3

gridded.pysgrid.utils.calculate_angle_from_true_east(lon_lat_1, lon_lat_2)[source]

Return the angle from true east in radians.

gridded.pysgrid.utils.calculate_bearing(lon_lat_1, lon_lat_2)[source]

return bearing from true north in degrees

gridded.pysgrid.utils.check_element_equal(lst)[source]

Check that all elements in an iterable are the same.

Params lst

iterable object to be checked

Returns

result of element equality check

Return type

bool

gridded.pysgrid.utils.determine_variable_slicing(sgrid_obj, nc_variable, method='center')[source]

Figure out how to slice a variable. This function only knows who to figure out slices that would be used to trim data before averaging to grid cell centers; grid cell nodes will be supported later.

Parameters
  • sgrid_obj (sgrid.SGrid) – an SGrid object derived from a netCDF file or netCDF4.Dataset object.

  • nc_dataset (netCDF4.Dataset) – a netCDF4.Dataset object from which the sgrid_obj was derived.

  • variable (str) – the name of a variable to be sliced.

  • method (str) – slice method for analysis at grid cell centers or grid cell nodes; accepts either ‘center’ or ‘node’.

Returns

the slice for the variable for the given method.

Return type

tuple

gridded.pysgrid.utils.does_intersection_exist(a, b)[source]
gridded.pysgrid.utils.infer_avg_axes(sgrid_obj, nc_var_obj)[source]

Infer which numpy axis to average over given the a variable defined on the grid. Works well for 2D. Not so sure about 3D.

gridded.pysgrid.utils.infer_variable_location(sgrid, variable)[source]
gridded.pysgrid.utils.pair_arrays(x_array, y_array)[source]

Given two arrays to equal dimensions, pair their values element-wise.

For example given arrays [[1, 2], [3, 4]] and [[-1, -2], [-3, -4]], this function will return [[[1, -1], [2, -2]], [[3, -3], [4, -4]]].

Parameters
  • x_array (np.array) – a numpy array containing “x” coordinates

  • y_array (np.array) – a numpy array containing “y” coordinates

Returns

array containing (x, y) arrays

Return type

np.array

gridded.pysgrid.utils.points_in_polys(points, polys, polyy=None)[source]
Parameters
  • points – Numpy array of Nx2 points

  • polys – Numpy array of N polygons of degree M represented by Mx2 points (NxMx2) for each point, see if respective poly contains it. Returns array of True/False

gridded.pysgrid.variables module

Created on Apr 15, 2015

@author: ayan

class gridded.pysgrid.variables.SGridVariable(center_axis=None, center_slicing=None, coordinates=None, dimensions=None, dtype=None, grid=None, location=None, node_axis=None, node_slicing=None, standard_name=None, variable=None, vector_axis=None, x_axis=None, y_axis=None, z_axis=None, data=None)[source]

Bases: object

Object of variables found and inferred from an SGRID compliant dataset.

classmethod create_var(nc_var_obj)[source]
classmethod create_variable(nc_var_obj, sgrid_obj)[source]
property data
property max
property min
property ndim
property shape

Module contents