opendrift.readers.basereader.structured
Attributes
Classes
A structured reader. Data is gridded on a regular grid. Used by e.g.: |
Module Contents
- opendrift.readers.basereader.structured.logger
- class opendrift.readers.basereader.structured.StructuredReader[source]
Bases:
opendrift.readers.basereader.variables.Variables
A structured reader. Data is gridded on a regular grid. Used by e.g.:
opendrift.readers.reader_netCDF_CF_generic.Reader
.Attributes:
projected: is True if
fakeproj.fakeproj
is used because of missing projection information. The data points are assumed to be approximately equidistant on the surface (i.e. in meters).clipped: pixels to to remove along boundary (e.g. in case of bad data).
See also
- clipped = 0
- x = None
- y = None
- interpolation = 'linearNDFast'
- convolve = None
- save_interpolator = None
- interpolator_filename = None
- __lonlat2xy_parallel__ = None
- __disable_parallel__ = False
- var_block_before
- var_block_after
- abstract get_variables(variables, time=None, x=None, y=None, z=None)[source]
Obtain a _block_ of values of the requested variables at all positions (x, y, z) closest to given time. These will be stored in
opendrift.readers.interpolation.structured.ReaderBlock
and accessed from there.- Arguments:
variables: list of variables.
time: datetime or None, time at which data are requested.
x, y: float or ndarrays; coordinates of requested points.
z: float or ndarray; vertical position (in meters, positive up)
- Returns:
Dictionary
keywords: variables (string) values: 2D ndarray bounding x and y.
- prepare(extent, start_time, end_time, max_speed)[source]
Prepare reader for given simulation coverage in time and space.
- set_convolution_kernel(convolve)[source]
Set a convolution kernel or kernel size (of array of ones) used by get_variables on read variables.
- _get_variables_interpolated_(variables, profiles, profiles_depth, time, reader_x, reader_y, z)[source]
This method _must_ be implemented by every reader. Usually by subclassing one of the reader types (e.g.
structured.StructuredReader
).Arguments are in _native projection_ of reader.
- __check_env_arrays__(env)[source]
For the StructuredReader the variables are checked before entered into the ReaderBlock interpolator. This methods makes the second check a no-op.
- get_ocean_depth_area_volume(lonmin, lonmax, latmin, latmax)[source]
Get depth, area and volume of ocean basin within given coordinates
- _make_projected_grid_(lon, lat, eq_eps=0.1)[source]
Make the projected grid in cases where lon and lat are present as 2D variables, but not x and y and assert that it is approximately equidistant.
Args:
eq_eps: tolerance for equidistance checks.
- __validate_projected_grid__(eq_eps=0.1)[source]
Validate that the projected grid is approximately equidistant.
Args:
eq_eps: tolerance for equidistance checks.
Raises:
AssertionError if not equidistant within eq_eps.
- _slice_variable_(var, indxTime=None, indy=None, indx=None, indz=None, indrealization=None)[source]
Slice variable depending on number of dimensions available.
Args:
All arguments can be slice objects or index.
Returns:
var sliced using the slices or indexes necessary to use depending on number of dimensions available.
Raises:
Unsupported number of dimensions (outside 2..5) raises an exception.