Function Reference¶
-
ocgis.
format_return
(ret_path, ops, with_auxiliary_files=False)[source]¶ Format an OpenClimateGIS path returning an absolute path to a zip file or simply passing ret_path through.
>>> import ocgis >>> ops = ocgis.OcgOperations(...) >>> ret = ops.execute() >>> new_path = ocgis.format_return(ret,ops)
Parameters: - ret_path (str) – The path returned by
ocgis.OcgOperations.execute()
. - ops (
ocgis.OcgOperations
) – Instance ofocgis.OcgOperations
. - with_auxiliary_files (bool) – If True, zip additional contents of CSV and netCDF returns.
Raises: NotImplementedError
Returns: str
- ret_path (str) – The path returned by
-
ocgis.util.helpers.
add_shapefile_unique_identifier
(in_path, out_path, name=None, template=None)[source]¶ >>> add_shapefile_unique_identifier('/path/to/foo.shp', '/path/to/new_foo.shp') '/path/to/new_foo.shp'
Parameters: - in_path (str) – Full path to the input shapefile.
- out_path (str) – Full path to the output shapefile.
- name (str) – The name of the unique identifer. If
None
, defaults toocgis.constants.OCGIS_UNIQUE_GEOMETRY_IDENTIFIER
. - template (str) – The integer attribute to copy as the unique identifier.
Returns: Path to the copied shapefile with the addition of a unique integer attribute called
name
.Return type:
-
ocgis.util.helpers.
get_bounds_from_1d
(centroids)[source]¶ Parameters: centroids ( numpy.ndarray
) – Vector representing center coordinates from which to interpolate bounds.Returns: A n-by-2 array with n equal to the shape of centroids
.>>> import numpy as np >>> centroids = np.array([1,2,3]) >>> get_bounds_from_1d(centroids) np.array([[0, 1],[1, 2],[2, 3]])
Return type: numpy.ndarray
Raises: NotImplementedError, ValueError
-
ocgis.util.helpers.
get_sorted_uris_by_time_dimension
(uris, variable=None)[source]¶ Sort a sequence of NetCDF URIs by the maximum time extent in ascending order.
Parameters: uris (list[str]) – The sequence of NetCDF URIs to sort. >>> uris = ['/path/to/file2.nc', 'path/to/file1.nc']
Parameters: variable (str) – The target variable for sorting. If None
is provided, then the variable will be autodiscovered.Returns: A sequence of sorted URIs. Return type: list[str]
-
ocgis.util.large_array.
compute
(ops, tile_dimension, verbose=False, use_optimizations=True)[source]¶ Used for computations on large arrays where memory limitations are a consideration. It is is also useful for extracting data from a server that has limitations on the size of requested data arrays. This function creates an empty destination NetCDF file that is then filled by executing the operations on chunks of the requested target dataset(s) and filling the destination NetCDF file.
Parameters: - ops (
ocgis.OcgOperations
) – The target operations to tile. There must be a calculation associated with the operations. - tile_dimension (int) – The target tile/chunk dimension. This integer value must be greater than zero.
- verbose (bool) – If
True
, print more verbose information to terminal. - use_optimizations (bool) – If
True
, cacheField
andTemporalGroupVariable
objects for reuse during tile iteration.
Raises: AssertionError, ValueError
Returns: Path to the output NetCDF file.
Return type: >>> from ocgis import RequestDataset, OcgOperations >>> from ocgis.util.large_array import compute >>> rd = RequestDataset(uri='/path/to/file', variable='tas') >>> ops = OcgOperations(dataset=rd, calc=[{'func':'mean','name':'mean'}],output_format='nc') >>> ret = compute(ops, 25)
- ops (
-
ocgis.regrid.base.
iter_esmf_fields
(ofield, regrid_method='auto', value_mask=None, split=True)[source]¶ For all data or a single time coordinate, yield an ESMF
Field
from the input OCGIS fieldofield
. Only one realization and level are allowed.Parameters: - ofield (
ocgis.Field
) – The input OCGIS Field object. - regrid_method – See
create_esmf_grid()
. - value_mask – See
create_esmf_grid()
. - split (bool) – If
True
, yield a single extra dimension slice from the source field. IfFalse
, yield all data. WhenFalse
, OCGIS uses ESMF’sndbounds
argument to field creation. UseTrue
if there are memory limitations. UseFalse
for faster performance.
Return type: tuple(str,
ESMF.driver.field.Field
, int)The returned tuple elements are:
Index Type Description 0 str The variable name currently being converted. 1 Field
The ESMF Field object. 2 int The current time index of the yielded ESMF field. If split=False
, This will beslice(None)
.Raises: AssertionError - ofield (
-
ocgis.regrid.base.
regrid_field
(source, destination, regrid_method='auto', value_mask=None, split=True, fill_value=None, weights_in=None, weights_out=None)[source]¶ Regrid
source
data to match the grid ofdestination
.Parameters: - source (
ocgis.Field
) – The source field. - destination (
ocgis.Field
) – The destination field. - regrid_method – See
create_esmf_grid()
. - value_mask (
numpy.ndarray
) – Seeiter_esmf_fields()
. - split (bool) – See
iter_esmf_fields()
. - fill_value (int | float) – Destination fill value used to fill the destination field before regridding. If
None
, then the default fill value for the destination field data type will be used. - weights_in (str) – Optional path to an input weights file. The route handle will be created from weights in this file. Assumes a SCRIP-like structure for the input weight file.
- weights_out (str) – Optional path to an output weight file. Does NOT do any regridding - just writes the weights.
Return type: - source (
-
ocgis.regrid.base.
create_esmf_grid
(ogrid, regrid_method='auto', value_mask=None)[source]¶ Create an ESMF
Grid
object from an OCGISocgis.Grid
object.Parameters: - ogrid (
Grid
) – The target OCGIS grid to convert to an ESMF grid. - regrid_method – If
'auto'
orESMF.api.constants.RegridMethod.CONSERVE
, use corners/bounds from the grid object. IfESMF.api.constants.RegridMethod.CONSERVE
, the corners/bounds must be present on the grid. - value_mask (
numpy.array
) – If anbool
numpy.array
with same shape asogrid
, use a logical or operation with the OCGIS field mask when creating the input grid’s mask. Values ofTrue
invalue_mask
are assumed to be masked. IfNone
is provided (the default) then the mask will be set using the spatial mask onogrid
.
Return type: ESMF.driver.grid.Grid
- ogrid (
-
ocgis.variable.
stack
(targets, stack_dim)[source]¶ Stack targets vertically using the stack dimension. For example, this function may be used to concatenate variables along the time dimension.
- For variables, the stack dimension object, value, and mask on the new variable will be a deep copy.
- For variables, the collection hierarchy is not traversed. Use the parent collections directly for the stack method.
- For collections, the returned value is a copy of the first field with stacked variables as deep copies. If a variable’s dimensions does not contain the stacked dimension, it is a returned as a copy.
Parameters: - targets ([
Variable
|VariableCollection
|Field
, …]) – List of variables, variable collections, or fields to stack vertically along the stack dimension. - stack_dim (
str
|Dimension
) – Dimension to use for vertical stacking.
Return type: Same as input type to
targets
.