API Reference#

This page provides an auto-generated summary of pop-tool’s API. For more details and examples, refer to the relevant chapters in the main part of the documentation.

Grid#

get_grid(grid_name[, scrip])

Return a xarray.Dataset() with POP grid variables.

Equation of State#

eos(salt, temp[, return_coefs])

Compute density as a function of salinity, temperature, and depth (or pressure).

compute_pressure

Convert depth in meters to pressure in bars.

Calculations#

cfc11sol(SALT, TEMP)

Compute CFC11 solubility in seawater. Reference: Warner & Weiss (1985) , Deep Sea Research, vol32 doi:10.1016/0198-0149(85)90099-8.

cfc12sol(SALT, TEMP)

Compute CFC12 solubility in seawater. Reference: Warner & Weiss (1985) , Deep Sea Research, vol32 doi:10.1016/0198-0149(85)90099-8.

Region masks#

region_mask_3d(grid_name[, mask_name, ...])

Generate as 3-D region mask with dimensions ('region', 'nlat', 'nlon').

list_region_masks(grid_name)

Get a list of the pre-defined region masks.

Utilities#

lateral_fill(da_in, isvalid_mask[, ...])

Perform lateral fill on xarray.DataArray

xgcm utilities#

to_xgcm_grid_dataset(ds[, metrics])

Modify POP model output to be compatible with xgcm.

pop_tools.get_grid(grid_name, scrip=False)[source]#

Return a xarray.Dataset() with POP grid variables.

Parameters:
  • grid_name (str) – Name of grid (i.e., POP_gx3v7, POP_gx1v7, POP_tx0.1v3)

  • scrip (boolean, optional) – Return grid in SCRIP format

Returns:

dso (xarray.Dataset) – Dataset containing POP grid variables.

pop_tools.eos(salt, temp, return_coefs=False, **kwargs)[source]#

Compute density as a function of salinity, temperature, and depth (or pressure).

McDougall, T.J., D.R. Jackett, D.G. Wright, and R. Feistel, 2003: Accurate and Computationally Efficient Algorithms for Potential Temperature and Density of Seawater. J. Atmos. Oceanic Technol., 20, 730–741, https://doi.org/10.1175/1520-0426(2003)20<730:AACEAF>2.0.CO;2.

test value:

rho = 1033.213387 kg/m^3; S = 35.0 PSU, theta = 20.0 C, pressure = 2000.0 dbars

Parameters:
  • salt (float) – salinity, psu

  • temp (float) – potential temperature, degree C

  • return_coefs (boolean, optional [default=False]) – Logical, if true function returns 2 additional arguments: dRHOdS and dRHOdT

  • depth (float, optional) – depth in meters, if not provided, pressure must be provided.

  • pressure (float, optional) – depth in dbar

Returns:

  • rho (float) – density kg/m^3

  • dRHOdS (float, optional) – Derivative of rho with respect to salinity.

  • dRHOdT (float, optional) – Derivative of rho with respect to temperature.

pop_tools.compute_pressure(depth)[source]#

Convert depth in meters to pressure in bars.

Parameters:

depth (float) – Depth in meters

Returns:

pressure (float) – Pressure in dbar

pop_tools.region_mask_3d(grid_name, mask_name=None, region_defs=None)[source]#

Generate as 3-D region mask with dimensions (‘region’, ‘nlat’, ‘nlon’).

Parameters:
  • grid_name (str) – Name of grid (i.e., POP_gx3v7, POP_gx1v7, POP_tx0.1v3)

  • mask_name (str, optional) – Name of the pre-defined mask. If not provided and region_defs is not provided, the default REGION_MASK is used.

  • region_defs (dict, optional) – Dictionary containing region definitions. The dictionary has the following form:

    region_defs = {region1_name: list_of_criteria_dicts_1, region2_name: list_of_criteria_dicts_2,…}

    list_of_criteria_dicts is a list of dictionaries; each must include the keys ‘match’ or ‘bounds’. For instance:

    list_of_criteria_dicts = [{‘match’: {‘REGION_MASK’: [1, 2, 3, 6]}, ‘bounds’: {‘TLAT’: [-90., -30.]}}]

    will identify a region where the default POP grid REGION_MASK values match the specified list and TLAT is between 90°S and 30°S. Multiple entries in the list_of_criteria_dicts are applied with an “or” condition.

Returns:

mask3d (xarray.DataArray) – Region mask DataArray with dimensions (‘region’, ‘nlat’, ‘nlon’). mask3d has “ones” within a region and “zeroes” outside it. A region coordinate is included that contains the region names.

pop_tools.list_region_masks(grid_name)[source]#

Get a list of the pre-defined region masks.

Parameters:

grid_name (str) – Name of grid (i.e., POP_gx3v7, POP_gx1v7, POP_tx0.1v3)

Returns:

region_mask_list (list) – List of defined region mask names.

pop_tools.lateral_fill(da_in, isvalid_mask, ltripole=False, tol=0.0001, use_sor=False, rc=1.8, max_iter=10000)[source]#

Perform lateral fill on xarray.DataArray

Parameters:
  • da_in (xarray.DataArray) – DataArray on which to fill NaNs. Fill is performed on the two rightmost dimenions. Grid is assumed periodic in x direction (last dimension).

  • isvalid_mask (xarray.DataArray, boolean) – Valid values mask: True where data should be filled. Must have the same rightmost dimenions as da_in.

  • ltripole (boolean, optional [default=False]) – Logical flag; if True then treat the top row of the grid as periodic

    in the sense of a tripole grid.

    tolfloat, optional [default=1.0e-4]

    Convergence criteria: stop filling when values change is less or equal to tol * var; i.e. delta <= tol * np.abs(var[j, i]).

  • use_sor (boolean, optional [default=False]) – switch to select SOR fill algorithm over progressive fill algorithm

  • rc (float, optional [default=1.8, valid bounds=(1.0,2.0)]) – over-relaxation coefficient to use in SOR fill algorithm. Larger arrrays (or extent of region to be filled if not global) typically converge faster with larger coefficients. For completely land-filling a 1 deg. grid (360x180) a coefficient in the range 1.85-1.9 is near optimal.

  • max_iter (integer, optional, [default=10000]) – maximum number of iterations to do before giving up if tol is not reached.

Returns:

da_out (xarray.DataArray) – DataArray with NaNs filled by iterative smoothing.

pop_tools.to_xgcm_grid_dataset(ds, metrics='detect', **kwargs)[source]#

Modify POP model output to be compatible with xgcm.

Parameters:
  • ds (xarray.Dataset) – An xarray Dataset

  • metrics ({"detect"} or dict, optional) – Dictionary providing metrics to the xgcm.Grid contructor. If "detect", will autodetect metrics that are present by searching for variables named DXU, DXT, DYU, DYT, DZU, DZT, UAREA, TAREA. If None, no metrics will be assigned.

  • kwargs – Additional keyword arguments are passed through to xgcm.Grid class.

Returns:

  • grid (xgcm.Grid) – An xgcm Grid object

  • ds_new (xarray.Dataset) – Xarray dataset with distinct dimensions for variables at different grid points.

Examples

>>> from pop_tools import DATASETS
>>> import xarray as xr
>>> filepath = DATASETS.fetch("g.e20.G.TL319_t13.control.001_hfreq-coarsen.nc")
>>> ds = xr.open_dataset(filepath)
>>> ds
<xarray.Dataset>
Dimensions:          (d2: 2, nlat: 240, nlon: 360, time: 1, z_t: 62, z_t_150m: 15, z_w: 62, z_w_top: 62)
Coordinates:
* time             (time) object 0050-01-01 00:00:00
* z_t              (z_t) float32 500.0 1500.0 2500.0 ... 562499.06 587499.06
* z_t_150m         (z_t_150m) float32 500.0 1500.0 2500.0 ... 13500.0 14500.0
* z_w              (z_w) float32 0.0 1000.0 2000.0 ... 549999.06 574999.06
* z_w_top          (z_w_top) float32 0.0 1000.0 2000.0 ... 549999.06 574999.06
    TLONG            (nlat, nlon) float64 ...
    TLAT             (nlat, nlon) float64 ...
Dimensions without coordinates: d2, nlat, nlon
Data variables:
    dz               (z_t) float32 ...
    dzw              (z_w) float32 ...
    KMT              (nlat, nlon) float64 ...
    REGION_MASK      (nlat, nlon) float64 ...
    TAREA            (nlat, nlon) float64 ...
    DXU              (nlat, nlon) float64 ...
    DYU              (nlat, nlon) float64 ...
    DYT              (nlat, nlon) float64 ...
    time_bound       (time, d2) object ...
    UEU              (time, z_t, nlat, nlon) float32 ...
    VNU              (time, z_t, nlat, nlon) float32 ...
    S_FLUX_ROFF_VSF  (time, z_w_top, nlat, nlon) float32 ...
>>> grid, ds_new = to_xgcm_grid_dataset(ds)
>>> grid
<xgcm.Grid>
X Axis (periodic):
    * center   nlon_t --> right
    * right    nlon_u --> center
Z Axis (periodic):
    * center   z_t --> left
    * left     z_w_top --> center
Y Axis (periodic):
    * center   nlat_t --> right
    * right    nlat_u --> center
>>> ds_new
<xarray.Dataset>
Dimensions:          (d2: 2, nlat_t: 240, nlat_u: 240, nlon_t: 360, nlon_u: 360, time: 1, z_t: 62, z_t_150m: 15, z_w: 62, z_w_top: 62)
Coordinates:
* time             (time) object 0050-01-01 00:00:00
* z_t              (z_t) float32 500.0 1500.0 2500.0 ... 562499.06 587499.06
* z_t_150m         (z_t_150m) float32 500.0 1500.0 2500.0 ... 13500.0 14500.0
* z_w              (z_w) float32 0.0 1000.0 2000.0 ... 549999.06 574999.06
* z_w_top          (z_w_top) float32 0.0 1000.0 2000.0 ... 549999.06 574999.06
  TLONG            (nlat_t, nlon_t) float64 nan nan nan nan ... nan nan nan
  TLAT             (nlat_t, nlon_t) float64 nan nan nan nan ... nan nan nan
* nlon_u           (nlon_u) int64 1 2 3 4 5 6 7 ... 355 356 357 358 359 360
* nlat_u           (nlat_u) int64 1 2 3 4 5 6 7 ... 235 236 237 238 239 240
* nlon_t           (nlon_t) float64 0.5 1.5 2.5 3.5 ... 357.5 358.5 359.5
* nlat_t           (nlat_t) float64 0.5 1.5 2.5 3.5 ... 237.5 238.5 239.5
Dimensions without coordinates: d2
Data variables:
    dz               (z_t) float32 ...
    dzw              (z_w) float32 ...
    KMT              (nlat_t, nlon_t) float64 nan nan nan nan ... nan nan nan
    REGION_MASK      (nlat_t, nlon_t) float64 nan nan nan nan ... nan nan nan
    TAREA            (nlat_t, nlon_t) float64 nan nan nan nan ... nan nan nan
    DXU              (nlat_u, nlon_u) float64 nan nan nan nan ... nan nan nan
    DYU              (nlat_u, nlon_u) float64 nan nan nan nan ... nan nan nan
    DYT              (nlat_t, nlon_t) float64 nan nan nan nan ... nan nan nan
    time_bound       (time, d2) object ...
    UEU              (time, z_t, nlat_u, nlon_u) float32 inf inf inf ... inf inf
    VNU              (time, z_t, nlat_t, nlon_u) float32 inf inf inf ... inf inf
    S_FLUX_ROFF_VSF  (time, z_w, nlat_t, nlon_t) float32 inf inf inf ... inf inf
pop_tools.cfc11sol(SALT, TEMP)[source]#

Compute CFC11 solubility in seawater. Reference:

Warner & Weiss (1985) , Deep Sea Research, vol32 doi:10.1016/0198-0149(85)90099-8

Parameters:
  • SALT (float) – salinity, psu

  • TEMP (float) – potential temperature, degree C

Returns:

SOLUBILITY_CFC11 (returned value in mol/m3/patm)

pop_tools.cfc12sol(SALT, TEMP)[source]#

Compute CFC12 solubility in seawater. Reference:

Warner & Weiss (1985) , Deep Sea Research, vol32 doi:10.1016/0198-0149(85)90099-8

Parameters:
  • SALT (float) – salinity, psu

  • TEMP (float) – potential temperature, degree C

Returns:

SOLUBILITY_CFC12 (returned value in mol/m3/patm)