Fieldsets and Fields#

parcels.fieldset module#

class parcels.fieldset.FieldSet(U, V, fields=None)[source]#

Bases: object

FieldSet class that holds hydrodynamic data needed to execute particles.

Parameters
  • U (parcels.field.Field) – Field object for zonal velocity component

  • V (parcels.field.Field) – Field object for meridional velocity component

  • fields (dict mapping str to Field) – Additional fields to include in the FieldSet. These fields can be used in custom kernels.

Methods

add_constant(name, value)

Add a constant to the FieldSet.

add_constant_field(name, value[, mesh])

Wrapper function to add a Field that is constant in space,

add_field(field[, name])

Add a parcels.field.Field object to the FieldSet.

add_periodic_halo([zonal, meridional, halosize])

Add a 'halo' to all parcels.field.Field objects in a FieldSet, through extending the Field (and lon/lat) by copying a small portion of the field on one side of the domain to the other.

add_vector_field(vfield)

Add a parcels.field.VectorField object to the FieldSet.

computeTimeChunk([time, dt])

Load a chunk of three data time steps into the FieldSet.

from_b_grid_dataset(filenames, variables, ...)

Initialises FieldSet object from NetCDF files of Bgrid fields.

from_c_grid_dataset(filenames, variables, ...)

Initialises FieldSet object from NetCDF files of Curvilinear NEMO fields.

from_data(data, dimensions[, transpose, ...])

Initialise FieldSet object from raw data.

from_mitgcm(filenames, variables, dimensions)

Initialises FieldSet object from NetCDF files of MITgcm fields. All parameters and keywords are exactly the same as for FieldSet.from_nemo(), except that gridindexing is set to 'mitgcm' for grids that have the shape::.

from_mom5(filenames, variables, dimensions)

Initialises FieldSet object from NetCDF files of MOM5 fields.

from_nemo(filenames, variables, dimensions)

Initialises FieldSet object from NetCDF files of Curvilinear NEMO fields.

from_netcdf(filenames, variables, dimensions)

Initialises FieldSet object from NetCDF files.

from_parcels(basename[, uvar, vvar, ...])

Initialises FieldSet data from NetCDF files using the Parcels FieldSet.write() conventions.

from_pop(filenames, variables, dimensions[, ...])

Initialises FieldSet object from NetCDF files of POP fields.

from_xarray_dataset(ds, variables, dimensions)

Initialises FieldSet data from xarray Datasets.

get_fields()

Returns a list of all the parcels.field.Field and parcels.field.VectorField objects associated with this FieldSet.

write(filename)

Write FieldSet to NetCDF file using NEMO convention.

add_UVfield

check_complete

checkvaliddimensionsdict

parse_wildcards

add_constant(name, value)[source]#

Add a constant to the FieldSet. Note that all constants are stored as 32-bit floats. While constants can be updated during execution in SciPy mode, they can not be updated in JIT mode.

Parameters
  • name (str) – Name of the constant

  • value – Value of the constant (stored as 32-bit float)

Examples

Tutorials using fieldset.add_constant: Analytical advection Diffusion Periodic boundaries

add_constant_field(name, value, mesh='flat')[source]#
Wrapper function to add a Field that is constant in space,

useful e.g. when using constant horizontal diffusivity

Parameters
  • name (str) – Name of the parcels.field.Field object to be added

  • value (float) – Value of the constant field (stored as 32-bit float)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also this tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

add_field(field, name=None)[source]#

Add a parcels.field.Field object to the FieldSet.

Parameters

Examples

For usage examples see the following tutorials:

add_periodic_halo(zonal=False, meridional=False, halosize=5)[source]#

Add a ‘halo’ to all parcels.field.Field objects in a FieldSet, through extending the Field (and lon/lat) by copying a small portion of the field on one side of the domain to the other.

Parameters
  • zonal (bool) – Create a halo in zonal direction (Default value = False)

  • meridional (bool) – Create a halo in meridional direction (Default value = False)

  • halosize (int) – size of the halo (in grid points). Default is 5 grid points

add_vector_field(vfield)[source]#

Add a parcels.field.VectorField object to the FieldSet.

Parameters

vfield (parcels.VectorField) – class:parcels.field.VectorField object to be added

computeTimeChunk(time=0.0, dt=1)[source]#

Load a chunk of three data time steps into the FieldSet. This is used when FieldSet uses data imported from netcdf, with default option deferred_load. The loaded time steps are at or immediatly before time and the two time steps immediately following time if dt is positive (and inversely for negative dt)

Parameters
  • time – Time around which the FieldSet chunks are to be loaded. Time is provided as a double, relatively to Fieldset.time_origin. Default is 0.

  • dt – time step of the integration scheme, needed to set the direction of time chunk loading. Default is 1.

classmethod from_b_grid_dataset(filenames, variables, dimensions, indices=None, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, tracer_interp_method='bgrid_tracer', chunksize=None, **kwargs)[source]#

Initialises FieldSet object from NetCDF files of Bgrid fields.

Parameters
  • filenames – Dictionary mapping variables to file(s). The filepath may contain wildcards to indicate multiple files, or be a list of file. filenames can be a list [files], a dictionary {var:[files]}, a dictionary {dim:[files]} (if lon, lat, depth and/or data not stored in same files as data), or a dictionary of dictionaries {var:{dim:[files]}} time values are in filenames[data]

  • variables (dict) – Dictionary mapping variables to variable names in the netCDF file(s).

  • dimensions (dict) –

    Dictionary mapping data dimensions (lon, lat, depth, time, data) to dimensions in the netCF file(s). Note that dimensions can also be a dictionary of dictionaries if dimension names are different for each variable. U and V velocity nodes are not located as W velocity and T tracer nodes (see http://www2.cesm.ucar.edu/models/cesm1.0/pop2/doc/sci/POPRefManual.pdf ).

    +-----------------------------+-----------------------------+-----------------------------+
    |U[k,j+1,i],V[k,j+1,i]        |                             |U[k,j+1,i+1],V[k,j+1,i+1]    |
    +-----------------------------+-----------------------------+-----------------------------+
    |                             |W[k:k+2,j+1,i+1],T[k,j+1,i+1]|                             |
    +-----------------------------+-----------------------------+-----------------------------+
    |U[k,j,i],V[k,j,i]            |                             |U[k,j,i+1],V[k,j,i+1]        |
    +-----------------------------+-----------------------------+-----------------------------+
    

    In 2D: U and V nodes are on the cell vertices and interpolated bilinearly as a A-grid. T node is at the cell centre and interpolated constant per cell as a C-grid. In 3D: U and V nodes are at the midlle of the cell vertical edges, They are interpolated bilinearly (independently of z) in the cell. W nodes are at the centre of the horizontal interfaces. They are interpolated linearly (as a function of z) in the cell. T node is at the cell centre, and constant per cell.

  • indices – Optional dictionary of indices for each dimension to read from file(s), to allow for reading of subset of data. Default is to read the full extent of each dimension. Note that negative indices are not allowed.

  • fieldtype – Optional dictionary mapping fields to fieldtypes to be used for UnitConverter. (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also this tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • tracer_interp_method (str) – Method for interpolation of tracer fields. It is recommended to use ‘bgrid_tracer’ (default) Note that in the case of from_pop() and from_b_grid_dataset(), the velocity fields are default to ‘bgrid_velocity’

  • chunksize – size of the chunks in dask loading (Default value = None)

  • **kwargs – Keyword arguments passed to the Fieldset.from_netcdf() constructor.

classmethod from_c_grid_dataset(filenames, variables, dimensions, indices=None, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, tracer_interp_method='cgrid_tracer', gridindexingtype='nemo', chunksize=None, **kwargs)[source]#

Initialises FieldSet object from NetCDF files of Curvilinear NEMO fields.

See here for a more detailed explanation of the different methods that can be used for c-grid datasets.

Parameters
  • filenames – Dictionary mapping variables to file(s). The filepath may contain wildcards to indicate multiple files, or be a list of file. filenames can be a list [files], a dictionary {var:[files]}, a dictionary {dim:[files]} (if lon, lat, depth and/or data not stored in same files as data), or a dictionary of dictionaries {var:{dim:[files]}} time values are in filenames[data]

  • variables (dict) – Dictionary mapping variables to variable names in the netCDF file(s).

  • dimensions (dict) –

    Dictionary mapping data dimensions (lon, lat, depth, time, data) to dimensions in the netCF file(s). Note that dimensions can also be a dictionary of dictionaries if dimension names are different for each variable. Watch out: NEMO is discretised on a C-grid: U and V velocities are not located on the same nodes (see https://www.nemo-ocean.eu/doc/node19.html ).

    +-----------------------------+-----------------------------+-----------------------------+
    |                             |         V[k,j+1,i+1]        |                             |
    +-----------------------------+-----------------------------+-----------------------------+
    |U[k,j+1,i]                   |W[k:k+2,j+1,i+1],T[k,j+1,i+1]|U[k,j+1,i+1]                 |
    +-----------------------------+-----------------------------+-----------------------------+
    |                             |         V[k,j,i+1]          |                             |
    +-----------------------------+-----------------------------+-----------------------------+
    

    To interpolate U, V velocities on the C-grid, Parcels needs to read the f-nodes, which are located on the corners of the cells. (for indexing details: https://www.nemo-ocean.eu/doc/img360.png ) In 3D, the depth is the one corresponding to W nodes.

  • indices – Optional dictionary of indices for each dimension to read from file(s), to allow for reading of subset of data. Default is to read the full extent of each dimension. Note that negative indices are not allowed.

  • fieldtype – Optional dictionary mapping fields to fieldtypes to be used for UnitConverter. (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also this tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • tracer_interp_method (str) – Method for interpolation of tracer fields. It is recommended to use ‘cgrid_tracer’ (default) Note that in the case of from_nemo() and from_cgrid(), the velocity fields are default to ‘cgrid_velocity’

  • gridindexingtype (str) – The type of gridindexing. Set to ‘nemo’ in FieldSet.from_nemo() See also the Grid indexing documentation on oceanparcels.org (Default value = ‘nemo’)

  • chunksize – size of the chunks in dask loading. (Default value = None)

  • **kwargs – Keyword arguments passed to the Fieldset.from_netcdf() constructor.

classmethod from_data(data, dimensions, transpose=False, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, **kwargs)[source]#

Initialise FieldSet object from raw data.

Parameters
  • data

    Dictionary mapping field names to numpy arrays. Note that at least a ‘U’ and ‘V’ numpy array need to be given, and that the built-in Advection kernels assume that U and V are in m/s

    1. If data shape is [xdim, ydim], [xdim, ydim, zdim], [xdim, ydim, tdim] or [xdim, ydim, zdim, tdim], whichever is relevant for the dataset, use the flag transpose=True

    2. If data shape is [ydim, xdim], [zdim, ydim, xdim], [tdim, ydim, xdim] or [tdim, zdim, ydim, xdim], use the flag transpose=False (default value)

    3. If data has any other shape, you first need to reorder it

  • dimensions (dict) – Dictionary mapping field dimensions (lon, lat, depth, time) to numpy arrays. Note that dimensions can also be a dictionary of dictionaries if dimension names are different for each variable (e.g. dimensions[‘U’], dimensions[‘V’], etc).

  • transpose (bool) – Whether to transpose data on read-in (Default value = False)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also this tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • **kwargs – Keyword arguments passed to the Field constructor.

Examples

For usage examples see the following tutorials:

classmethod from_mitgcm(filenames, variables, dimensions, indices=None, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, tracer_interp_method='cgrid_tracer', chunksize=None, **kwargs)[source]#

Initialises FieldSet object from NetCDF files of MITgcm fields. All parameters and keywords are exactly the same as for FieldSet.from_nemo(), except that gridindexing is set to ‘mitgcm’ for grids that have the shape:

+-----------------------------+-----------------------------+-----------------------------+
|                             |         V[k,j+1,i]          |                             |
+-----------------------------+-----------------------------+-----------------------------+
|U[k,j,i]                     |    W[k-1:k,j,i], T[k,j,i]   |U[k,j,i+1]                   |
+-----------------------------+-----------------------------+-----------------------------+
|                             |         V[k,j,i]            |                             |
+-----------------------------+-----------------------------+-----------------------------+

For indexing details: https://mitgcm.readthedocs.io/en/latest/algorithm/algorithm.html#spatial-discretization-of-the-dynamical-equations Note that vertical velocity (W) is assumed positive in the positive z direction (which is upward in MITgcm)

classmethod from_mom5(filenames, variables, dimensions, indices=None, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, tracer_interp_method='bgrid_tracer', chunksize=None, **kwargs)[source]#

Initialises FieldSet object from NetCDF files of MOM5 fields.

Parameters
  • filenames – Dictionary mapping variables to file(s). The filepath may contain wildcards to indicate multiple files, or be a list of file. filenames can be a list [files], a dictionary {var:[files]}, a dictionary {dim:[files]} (if lon, lat, depth and/or data not stored in same files as data), or a dictionary of dictionaries {var:{dim:[files]}} time values are in filenames[data]

  • variables (dict) – Dictionary mapping variables to variable names in the netCDF file(s). Note that the built-in Advection kernels assume that U and V are in m/s

  • dimensions (dict) –

    Dictionary mapping data dimensions (lon, lat, depth, time, data) to dimensions in the netCF file(s). Note that dimensions can also be a dictionary of dictionaries if dimension names are different for each variable.

    +-------------------------------+-------------------------------+-------------------------------+
    |U[k,j+1,i],V[k,j+1,i]          |                               |U[k,j+1,i+1],V[k,j+1,i+1]      |
    +-------------------------------+-------------------------------+-------------------------------+
    |                               |W[k-1:k+1,j+1,i+1],T[k,j+1,i+1]|                               |
    +-------------------------------+-------------------------------+-------------------------------+
    |U[k,j,i],V[k,j,i]              |                               |U[k,j,i+1],V[k,j,i+1]          |
    +-------------------------------+-------------------------------+-------------------------------+
    

    In 2D: U and V nodes are on the cell vertices and interpolated bilinearly as a A-grid. T node is at the cell centre and interpolated constant per cell as a C-grid. In 3D: U and V nodes are at the middle of the cell vertical edges, They are interpolated bilinearly (independently of z) in the cell. W nodes are at the centre of the horizontal interfaces, but below the U and V. They are interpolated linearly (as a function of z) in the cell. Note that W is normally directed upward in MOM5, but Parcels requires W in the positive z-direction (downward) so W is multiplied by -1. T node is at the cell centre, and constant per cell.

  • indices – Optional dictionary of indices for each dimension to read from file(s), to allow for reading of subset of data. Default is to read the full extent of each dimension. Note that negative indices are not allowed.

  • fieldtype – Optional dictionary mapping fields to fieldtypes to be used for UnitConverter. (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also the Unit converters tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • tracer_interp_method (str) – Method for interpolation of tracer fields. It is recommended to use ‘bgrid_tracer’ (default) Note that in the case of from_mom5() and from_b_grid_dataset(), the velocity fields are default to ‘bgrid_velocity’

  • chunksize – size of the chunks in dask loading (Default value = None)

  • **kwargs – Keyword arguments passed to the Fieldset.from_b_grid_dataset() constructor.

classmethod from_nemo(filenames, variables, dimensions, indices=None, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, tracer_interp_method='cgrid_tracer', chunksize=None, **kwargs)[source]#

Initialises FieldSet object from NetCDF files of Curvilinear NEMO fields.

See here for a detailed tutorial on the setup for 2D NEMO fields and here for the tutorial on the setup for 3D NEMO fields.

See here for a more detailed explanation of the different methods that can be used for c-grid datasets.

Parameters
  • filenames – Dictionary mapping variables to file(s). The filepath may contain wildcards to indicate multiple files, or be a list of file. filenames can be a list [files], a dictionary {var:[files]}, a dictionary {dim:[files]} (if lon, lat, depth and/or data not stored in same files as data), or a dictionary of dictionaries {var:{dim:[files]}} time values are in filenames[data]

  • variables (dict) – Dictionary mapping variables to variable names in the netCDF file(s). Note that the built-in Advection kernels assume that U and V are in m/s

  • dimensions (dict) –

    Dictionary mapping data dimensions (lon, lat, depth, time, data) to dimensions in the netCF file(s). Note that dimensions can also be a dictionary of dictionaries if dimension names are different for each variable. Watch out: NEMO is discretised on a C-grid: U and V velocities are not located on the same nodes (see https://www.nemo-ocean.eu/doc/node19.html).

    +-----------------------------+-----------------------------+-----------------------------+
    |                             |         V[k,j+1,i+1]        |                             |
    +-----------------------------+-----------------------------+-----------------------------+
    |U[k,j+1,i]                   |W[k:k+2,j+1,i+1],T[k,j+1,i+1]|U[k,j+1,i+1]                 |
    +-----------------------------+-----------------------------+-----------------------------+
    |                             |         V[k,j,i+1]          |                             |
    +-----------------------------+-----------------------------+-----------------------------+
    

    To interpolate U, V velocities on the C-grid, Parcels needs to read the f-nodes, which are located on the corners of the cells. (for indexing details: https://www.nemo-ocean.eu/doc/img360.png ) In 3D, the depth is the one corresponding to W nodes The gridindexingtype is set to ‘nemo’. See also the Grid indexing documentation on oceanparcels.org

  • indices – Optional dictionary of indices for each dimension to read from file(s), to allow for reading of subset of data. Default is to read the full extent of each dimension. Note that negative indices are not allowed.

  • fieldtype – Optional dictionary mapping fields to fieldtypes to be used for UnitConverter. (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also this tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • tracer_interp_method (str) – Method for interpolation of tracer fields. It is recommended to use ‘cgrid_tracer’ (default) Note that in the case of from_nemo() and from_cgrid(), the velocity fields are default to ‘cgrid_velocity’

  • chunksize – size of the chunks in dask loading. Default is None (no chunking)

  • **kwargs – Keyword arguments passed to the Fieldset.from_c_grid_dataset() constructor.

classmethod from_netcdf(filenames, variables, dimensions, indices=None, fieldtype=None, mesh='spherical', timestamps=None, allow_time_extrapolation=None, time_periodic=False, deferred_load=True, chunksize=None, **kwargs)[source]#

Initialises FieldSet object from NetCDF files.

Parameters
  • filenames – Dictionary mapping variables to file(s). The filepath may contain wildcards to indicate multiple files or be a list of file. filenames can be a list [files], a dictionary {var:[files]}, a dictionary {dim:[files]} (if lon, lat, depth and/or data not stored in same files as data), or a dictionary of dictionaries {var:{dim:[files]}}. time values are in filenames[data]

  • variables (dict) – Dictionary mapping variables to variable names in the netCDF file(s). Note that the built-in Advection kernels assume that U and V are in m/s

  • dimensions (dict) – Dictionary mapping data dimensions (lon, lat, depth, time, data) to dimensions in the netCF file(s). Note that dimensions can also be a dictionary of dictionaries if dimension names are different for each variable (e.g. dimensions[‘U’], dimensions[‘V’], etc).

  • indices – Optional dictionary of indices for each dimension to read from file(s), to allow for reading of subset of data. Default is to read the full extent of each dimension. Note that negative indices are not allowed.

  • fieldtype – Optional dictionary mapping fields to fieldtypes to be used for UnitConverter. (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None) (Default value = None)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also this tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • timestamps – list of lists or array of arrays containing the timestamps for each of the files in filenames. Outer list/array corresponds to files, inner array corresponds to indices within files. Default is None if dimensions includes time.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • deferred_load (bool) – boolean whether to only pre-load data (in deferred mode) or fully load them (default: True). It is advised to deferred load the data, since in that case Parcels deals with a better memory management during particle set execution. deferred_load=False is however sometimes necessary for plotting the fields.

  • interp_method (str) – Method for interpolation. Options are ‘linear’ (default), ‘nearest’, ‘linear_invdist_land_tracer’, ‘cgrid_velocity’, ‘cgrid_tracer’ and ‘bgrid_velocity’

  • gridindexingtype (str) – The type of gridindexing. Either ‘nemo’ (default) or ‘mitgcm’ are supported. See also the Grid indexing documentation on oceanparcels.org

  • chunksize – size of the chunks in dask loading. Default is None (no chunking). Can be None or False (no chunking), ‘auto’ (chunking is done in the background, but results in one grid per field individually), or a dict in the format {parcels_varname: {netcdf_dimname : (parcels_dimname, chunksize_as_int)}, ...}, where parcels_dimname is one of (‘time’, ‘depth’, ‘lat’, ‘lon’)

  • netcdf_engine – engine to use for netcdf reading in xarray. Default is ‘netcdf’, but in cases where this doesn’t work, setting netcdf_engine=’scipy’ could help

  • **kwargs – Keyword arguments passed to the parcels.Field constructor.

Examples

For usage examples see the following tutorials:

classmethod from_parcels(basename, uvar='vozocrtx', vvar='vomecrty', indices=None, extra_fields=None, allow_time_extrapolation=None, time_periodic=False, deferred_load=True, chunksize=None, **kwargs)[source]#

Initialises FieldSet data from NetCDF files using the Parcels FieldSet.write() conventions.

Parameters
  • basename (str) – Base name of the file(s); may contain wildcards to indicate multiple files.

  • indices – Optional dictionary of indices for each dimension to read from file(s), to allow for reading of subset of data. Default is to read the full extent of each dimension. Note that negative indices are not allowed.

  • fieldtype – Optional dictionary mapping fields to fieldtypes to be used for UnitConverter. (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None)

  • extra_fields – Extra fields to read beyond U and V (Default value = None)

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • deferred_load (bool) – boolean whether to only pre-load data (in deferred mode) or fully load them (default: True). It is advised to deferred load the data, since in that case Parcels deals with a better memory management during particle set execution. deferred_load=False is however sometimes necessary for plotting the fields.

  • chunksize – size of the chunks in dask loading (Default value = None)

  • uvar – (Default value = ‘vozocrtx’)

  • vvar – (Default value = ‘vomecrty’)

  • **kwargs – Keyword arguments passed to the Fieldset.from_netcdf() constructor.

classmethod from_pop(filenames, variables, dimensions, indices=None, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, tracer_interp_method='bgrid_tracer', chunksize=None, depth_units='m', **kwargs)[source]#
Initialises FieldSet object from NetCDF files of POP fields.

It is assumed that the velocities in the POP fields is in cm/s.

Parameters
  • filenames – Dictionary mapping variables to file(s). The filepath may contain wildcards to indicate multiple files, or be a list of file. filenames can be a list [files], a dictionary {var:[files]}, a dictionary {dim:[files]} (if lon, lat, depth and/or data not stored in same files as data), or a dictionary of dictionaries {var:{dim:[files]}} time values are in filenames[data]

  • variables (dict) – Dictionary mapping variables to variable names in the netCDF file(s). Note that the built-in Advection kernels assume that U and V are in m/s

  • dimensions (dict) –

    Dictionary mapping data dimensions (lon, lat, depth, time, data) to dimensions in the netCF file(s). Note that dimensions can also be a dictionary of dictionaries if dimension names are different for each variable. Watch out: POP is discretised on a B-grid: U and V velocity nodes are not located as W velocity and T tracer nodes (see http://www2.cesm.ucar.edu/models/cesm1.0/pop2/doc/sci/POPRefManual.pdf ).

    +-----------------------------+-----------------------------+-----------------------------+
    |U[k,j+1,i],V[k,j+1,i]        |                             |U[k,j+1,i+1],V[k,j+1,i+1]    |
    +-----------------------------+-----------------------------+-----------------------------+
    |                             |W[k:k+2,j+1,i+1],T[k,j+1,i+1]|                             |
    +-----------------------------+-----------------------------+-----------------------------+
    |U[k,j,i],V[k,j,i]            |                             |U[k,j,i+1],V[k,j,i+1]        |
    +-----------------------------+-----------------------------+-----------------------------+
    

    In 2D: U and V nodes are on the cell vertices and interpolated bilinearly as a A-grid. T node is at the cell centre and interpolated constant per cell as a C-grid. In 3D: U and V nodes are at the middle of the cell vertical edges, They are interpolated bilinearly (independently of z) in the cell. W nodes are at the centre of the horizontal interfaces. They are interpolated linearly (as a function of z) in the cell. T node is at the cell centre, and constant per cell. Note that Parcels assumes that the length of the depth dimension (at the W-points) is one larger than the size of the velocity and tracer fields in the depth dimension.

  • indices – Optional dictionary of indices for each dimension to read from file(s), to allow for reading of subset of data. Default is to read the full extent of each dimension. Note that negative indices are not allowed.

  • fieldtype – Optional dictionary mapping fields to fieldtypes to be used for UnitConverter. (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also this tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • tracer_interp_method (str) – Method for interpolation of tracer fields. It is recommended to use ‘bgrid_tracer’ (default) Note that in the case of from_pop() and from_bgrid(), the velocity fields are default to ‘bgrid_velocity’

  • chunksize – size of the chunks in dask loading (Default value = None)

  • depth_units – The units of the vertical dimension. Default in Parcels is ‘m’, but many POP outputs are in ‘cm’

  • **kwargs – Keyword arguments passed to the Fieldset.from_b_grid_dataset() constructor.

classmethod from_xarray_dataset(ds, variables, dimensions, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, **kwargs)[source]#

Initialises FieldSet data from xarray Datasets.

Parameters
  • ds (xr.Dataset) – xarray Dataset. Note that the built-in Advection kernels assume that U and V are in m/s

  • variables (dict) – Dictionary mapping parcels variable names to data variables in the xarray Dataset.

  • dimensions (dict) – Dictionary mapping data dimensions (lon, lat, depth, time, data) to dimensions in the xarray Dataset. Note that dimensions can also be a dictionary of dictionaries if dimension names are different for each variable (e.g. dimensions[‘U’], dimensions[‘V’], etc).

  • fieldtype – Optional dictionary mapping fields to fieldtypes to be used for UnitConverter. (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation, see also this tutorial:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • **kwargs – Keyword arguments passed to the Field.from_xarray() constructor.

get_fields()[source]#

Returns a list of all the parcels.field.Field and parcels.field.VectorField objects associated with this FieldSet.

write(filename)[source]#

Write FieldSet to NetCDF file using NEMO convention.

Parameters

filename (str) – Basename of the output fileset.

parcels.field module#

class parcels.field.Field(name, data, lon=None, lat=None, depth=None, time=None, grid=None, mesh='flat', timestamps=None, fieldtype=None, transpose=False, vmin=None, vmax=None, cast_data_dtype='float32', time_origin=None, interp_method='linear', allow_time_extrapolation=None, time_periodic=False, gridindexingtype='nemo', to_write=False, **kwargs)[source]#

Bases: object

Class that encapsulates access to field data.

Parameters
  • name (str) – Name of the field

  • data (np.ndarray) –

    2D, 3D or 4D numpy array of field data.

    1. If data shape is [xdim, ydim], [xdim, ydim, zdim], [xdim, ydim, tdim] or [xdim, ydim, zdim, tdim], whichever is relevant for the dataset, use the flag transpose=True

    2. If data shape is [ydim, xdim], [zdim, ydim, xdim], [tdim, ydim, xdim] or [tdim, zdim, ydim, xdim], use the flag transpose=False

    3. If data has any other shape, you first need to reorder it

  • lon (np.ndarray or list) – Longitude coordinates (numpy vector or array) of the field (only if grid is None)

  • lat (np.ndarray or list) – Latitude coordinates (numpy vector or array) of the field (only if grid is None)

  • depth (np.ndarray or list) – Depth coordinates (numpy vector or array) of the field (only if grid is None)

  • time (np.ndarray) – Time coordinates (numpy vector) of the field (only if grid is None)

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation: (only if grid is None)

    1. spherical: Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat (default): No conversion, lat/lon are assumed to be in m.

  • timestamps (np.ndarray) – A numpy array containing the timestamps for each of the files in filenames, for loading from netCDF files only. Default is None if the netCDF dimensions dictionary includes time.

  • grid (parcels.grid.Grid) – parcels.grid.Grid object containing all the lon, lat depth, time mesh and time_origin information. Can be constructed from any of the Grid objects

  • fieldtype (str) – Type of Field to be used for UnitConverter (either ‘U’, ‘V’, ‘Kh_zonal’, ‘Kh_meridional’ or None)

  • transpose (bool) – Transpose data to required (lon, lat) layout

  • vmin (float) – Minimum allowed value on the field. Data below this value are set to zero

  • vmax (float) – Maximum allowed value on the field. Data above this value are set to zero

  • cast_data_dtype (str) – Cast Field data to dtype. Supported dtypes are “float32” (np.float32 (default)) and “float64 (np.float64). Note that dtype can only be “float32” in JIT mode

  • time_origin (parcels.tools.converters.TimeConverter) – Time origin of the time axis (only if grid is None)

  • interp_method (str) – Method for interpolation. Options are ‘linear’ (default), ‘nearest’, ‘linear_invdist_land_tracer’, ‘cgrid_velocity’, ‘cgrid_tracer’ and ‘bgrid_velocity’

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation in time (i.e. beyond the last available time snapshot)

  • time_periodic (bool, float or datetime.timedelta) – To loop periodically over the time component of the Field. It is set to either False or the length of the period (either float in seconds or datetime.timedelta object). The last value of the time series can be provided (which is the same as the initial one) or not (Default: False) This flag overrides the allow_time_extrapolation and sets it to False

  • chunkdims_name_map (str, optional) – Gives a name map to the FieldFileBuffer that declared a mapping between chunksize name, NetCDF dimension and Parcels dimension; required only if currently incompatible OCM field is loaded and chunking is used by ‘chunksize’ (which is the default)

  • to_write (bool) – Write the Field in NetCDF format at the same frequency as the ParticleFile outputdt, using a filenaming scheme based on the ParticleFile name

Examples

For usage examples see the following tutorials:

Attributes
ctypes_struct

Returns a ctypes struct object containing all relevant pointers and sizes for this field.

Methods

add_periodic_halo(zonal, meridional[, ...])

Add a 'halo' to all Fields in a FieldSet.

calc_cell_edge_sizes()

Method to calculate cell sizes based on numpy.gradient method.

cell_areas()

Method to calculate cell sizes based on cell_edge_sizes.

eval(time, z, y, x[, particle, applyConversion])

Interpolate field values in space and time.

from_netcdf(filenames, variable, dimensions)

Create field from netCDF file.

from_xarray(da, name, dimensions[, mesh, ...])

Create field from xarray Variable.

set_depth_from_field(field)

Define the depth dimensions from another (time-varying) field.

set_scaling_factor(factor)

Scales the field data by some constant factor.

spatial_interpolation(ti, z, y, x, time[, ...])

Interpolate horizontal field values using a SciPy interpolator.

temporal_interpolate_fullfield(ti, time)

Calculate the data of a field between two snapshots using linear interpolation.

time_index(time)

Find the index in the time array associated with a given time.

write(filename[, varname])

Write a Field to a netcdf file.

ccode_convert

ccode_eval

chunk_data

chunk_setup

collect_timeslices

computeTimeChunk

data_concatenate

get_block

get_block_id

get_dim_filenames

interpolator2D

interpolator3D

reconnect_bnd_indices

rescale_and_set_minmax

reshape

search_indices

search_indices_curvilinear

search_indices_rectilinear

search_indices_vertical_s

search_indices_vertical_z

add_periodic_halo(zonal, meridional, halosize=5, data=None)[source]#

Add a ‘halo’ to all Fields in a FieldSet.

Add a ‘halo’ to all Fields in a FieldSet, through extending the Field (and lon/lat) by copying a small portion of the field on one side of the domain to the other. Before adding a periodic halo to the Field, it has to be added to the Grid on which the Field depends

See this tutorial for a detailed explanation on how to set up periodic boundaries

Parameters
  • zonal (bool) – Create a halo in zonal direction.

  • meridional (bool) – Create a halo in meridional direction.

  • halosize (int) – Size of the halo (in grid points). Default is 5 grid points

  • data – if data is not None, the periodic halo will be achieved on data instead of self.data and data will be returned (Default value = None)

calc_cell_edge_sizes()[source]#

Method to calculate cell sizes based on numpy.gradient method.

Currently only works for Rectilinear Grids

cell_areas()[source]#

Method to calculate cell sizes based on cell_edge_sizes.

Currently only works for Rectilinear Grids

property ctypes_struct[source]#

Returns a ctypes struct object containing all relevant pointers and sizes for this field.

eval(time, z, y, x, particle=None, applyConversion=True)[source]#

Interpolate field values in space and time.

We interpolate linearly in time and apply implicit unit conversion to the result. Note that we defer to scipy.interpolate to perform spatial interpolation.

classmethod from_netcdf(filenames, variable, dimensions, indices=None, grid=None, mesh='spherical', timestamps=None, allow_time_extrapolation=None, time_periodic=False, deferred_load=True, **kwargs)[source]#

Create field from netCDF file.

Parameters
  • filenames (list of str or dict) – list of filenames to read for the field. filenames can be a list [files] or a dictionary {dim:[files]} (if lon, lat, depth and/or data not stored in same files as data) In the latter case, time values are in filenames[data]

  • variable (dict, tuple of str or str) – Dict or tuple mapping field name to variable name in the NetCDF file.

  • dimensions (dict) – Dictionary mapping variable names for the relevant dimensions in the NetCDF file

  • indices – dictionary mapping indices for each dimension to read from file. This can be used for reading in only a subregion of the NetCDF file. Note that negative indices are not allowed. (Default value = None)

  • mesh

    String indicating the type of mesh coordinates and units used during velocity interpolation:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • timestamps – A numpy array of datetime64 objects containing the timestamps for each of the files in filenames. Default is None if dimensions includes time.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation in time (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – boolean whether to loop periodically over the time component of the FieldSet This flag overrides the allow_time_extrapolation and sets it to False (Default value = False)

  • deferred_load (bool) – boolean whether to only pre-load data (in deferred mode) or fully load them (default: True). It is advised to deferred load the data, since in that case Parcels deals with a better memory management during particle set execution. deferred_load=False is however sometimes necessary for plotting the fields.

  • gridindexingtype (str) – The type of gridindexing. Either ‘nemo’ (default) or ‘mitgcm’ are supported. See also the Grid indexing documentation on oceanparcels.org

  • chunksize – size of the chunks in dask loading

  • netcdf_decodewarning (bool) – Whether to show a warning id there is a problem decoding the netcdf files. Default is True, but in some cases where these warnings are expected, it may be useful to silence them by setting netcdf_decodewarning=False.

  • grid – (Default value = None)

  • **kwargs – Keyword arguments passed to the Field constructor.

Examples

For usage examples see the following tutorial:

classmethod from_xarray(da, name, dimensions, mesh='spherical', allow_time_extrapolation=None, time_periodic=False, **kwargs)[source]#

Create field from xarray Variable.

Parameters
  • da (xr.DataArray) – Xarray DataArray

  • name (str) – Name of the Field

  • dimensions (dict) – Dictionary mapping variable names for the relevant dimensions in the DataArray

  • mesh (str) –

    String indicating the type of mesh coordinates and units used during velocity interpolation:

    1. spherical (default): Lat and lon in degree, with a correction for zonal velocity U near the poles.

    2. flat: No conversion, lat/lon are assumed to be in m.

  • allow_time_extrapolation (bool) – boolean whether to allow for extrapolation in time (i.e. beyond the last available time snapshot) Default is False if dimensions includes time, else True

  • time_periodic (bool, float or datetime.timedelta) – boolean whether to loop periodically over the time component of the FieldSet This flag overrides the allow_time_extrapolation and sets it to False (Default value = False)

  • **kwargs – Keyword arguments passed to the Field constructor.

set_depth_from_field(field)[source]#

Define the depth dimensions from another (time-varying) field.

Notes

See this tutorial for a detailed explanation on how to set up time-evolving depth dimensions.

set_scaling_factor(factor)[source]#

Scales the field data by some constant factor.

Parameters

factor – scaling factor

Examples

For usage examples see the following tutorial:

spatial_interpolation(ti, z, y, x, time, particle=None)[source]#

Interpolate horizontal field values using a SciPy interpolator.

temporal_interpolate_fullfield(ti, time)[source]#

Calculate the data of a field between two snapshots using linear interpolation.

Parameters
  • ti – Index in time array associated with time (via time_index())

  • time – Time to interpolate to

time_index(time)[source]#

Find the index in the time array associated with a given time.

Note that we normalize to either the first or the last index if the sampled value is outside the time value range.

write(filename, varname=None)[source]#

Write a Field to a netcdf file.

Parameters
  • filename (str) – Basename of the file (i.e. ‘{filename}{Field.name}.nc’)

  • varname (str) – Name of the field, to be appended to the filename. (Default value = None)

class parcels.field.NestedField(name, F, V=None, W=None)[source]#

Bases: list

NestedField is a class that allows for interpolation of fields on different grids of potentially varying resolution.

The NestedField class is a list of Fields where the first Field that contains the particle within the domain is then used for interpolation. This induces that the order of the fields in the list matters. Each one it its turn, a field is interpolated: if the interpolation succeeds or if an error other than ErrorOutOfBounds is thrown, the function is stopped. Otherwise, next field is interpolated. NestedField returns an ErrorOutOfBounds only if last field is as well out of boundaries. NestedField is composed of either Fields or VectorFields.

Parameters
  • name (str) – Name of the NestedField

  • F (list of Field) – List of fields (order matters). F can be a scalar Field, a VectorField, or the zonal component (U) of the VectorField

  • V (list of Field) – List of fields defining the meridional component of a VectorField, if F is the zonal component. (default: None)

  • W (list of Field) – List of fields defining the vertical component of a VectorField, if F and V are the zonal and meridional components (default: None)

Examples

See here for a detailed tutorial

Methods

append(object, /)

Append object to the end of the list.

clear(/)

Remove all items from list.

copy(/)

Return a shallow copy of the list.

count(value, /)

Return number of occurrences of value.

extend(iterable, /)

Extend list by appending elements from the iterable.

index(value[, start, stop])

Return first index of value.

insert(index, object, /)

Insert object before index.

pop([index])

Remove and return item at index (default last).

remove(value, /)

Remove first occurrence of value.

reverse(/)

Reverse IN PLACE.

sort(*[, key, reverse])

Sort the list in ascending order and return None.

class parcels.field.VectorField(name, U, V, W=None)[source]#

Bases: object

Class VectorField stores 2 or 3 fields which defines together a vector field. This enables to interpolate them as one single vector field in the kernels.

Parameters

Methods

spatial_c_grid_interpolation3D(ti, z, y, x, time)

Perform C grid interpolation in 3D. ::.

ccode_eval

dist

eval

jacobian

spatial_c_grid_interpolation2D

spatial_c_grid_interpolation3D_full

spatial_slip_interpolation

spatial_c_grid_interpolation3D(ti, z, y, x, time, particle=None, applyConversion=True)[source]#

Perform C grid interpolation in 3D.

+---+---+---+
|   |V1 |   |
+---+---+---+
|U0 |   |U1 |
+---+---+---+
|   |V0 |   |
+---+---+---+

The interpolation is done in the following by interpolating linearly U depending on the longitude coordinate and interpolating linearly V depending on the latitude coordinate. Curvilinear grids are treated properly, since the element is projected to a rectilinear parent element.