xwmt.wmt module#

class xwmt.wmt.WaterMassTransformations(grid, xbudget_dict, mask=None, teos10=True, cp=3992.0, rho_ref=1035.0, t_var='conservative', s_var='absolute', method='default', rebin=False)#

Bases: WaterMass

An extension of the WaterMass class that includes methods for WaterMass transformation analysis.

available_processes(available=True)#

Get a list of all tendency processes that are both specified by xbudget_dict and available in the dataset.

Parameters:

available (bool) – Default True. If False, include processes that are not available in the dataset.

Returns:

  • processes (list of str)

  • >>> wmt = WaterMassTransformations(grid, xbudget_dict)

  • >>> processes = wmt.available_processes()

  • [‘diffusion’]

calc_hlamdot_and_lambda(lambda_name, term)#

Get layer-integrated extensive tracer tendencies (* kg/m^2/s) and corresponding scalar field of lambda

Parameters:
  • lambda_name (str) – Specifies lambda

  • term (str) – key for tendency variable in the xbudget_dict

Returns:

hlamdot, lam

Return type:

xr.DataArray, xr.DataArray

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> (hlamdot, lam) = wmt.calc_hlamdot_and_lambda("heat", "diffusion")
datadict(tracer, term)#

Get a dictionary that organizes the variables and metadata necessary to evaluate tendency terms.

Parameters:
  • component (str) – Supported options: [“heat”, “salt”]

  • term (str) – key for tendency variable in the xbudget_dict

Returns:

ddict

Return type:

dict

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> ddict = wmt.datadict("heat", "diffusion")
get_lambda_key(lambda_var)#

Search lambda dictionary for lambda corresponding to variable name

Parameters:

lambda_var (str) – Variable name of lambda.

Return type:

str

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> wmt.get_lambda_key("thetao")
'heat'
>>> wmt.get_lambda_var("sigma2")
'density'

See also

self.lambdas, self.get_lambda_var

get_lambda_var(lambda_name=None)#

Search lambda dictionary for variable name corresponding to specified lambda

Parameters:

lambda_name (str (default: None)) – Name of lambda variable. Supported options: [“heat”, “salt”, “density”] If None, return None

Return type:

str or list

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> wmt.get_lambda_var("heat")
'thetao'
>>> wmt.get_lambda_var("density")
['sigma0', 'sigma1', 'sigma2', 'sigma3', 'sigma4']

See also

self.lambdas, self.get_lambda_key

integrate_transformations(lambda_name, term=None, group_processes=False, sum_components=True, bins=None, mask=None)#

Lazily compute horizontally-integrated water mass transformations.

Parameters:
  • lambda_name (str) – Specifies lambda (e.g., ‘heat’, ‘salt’, ‘sigma0’, etc.). Use lambdas() for a list of available lambdas.

  • term (str or list of str (default: None)) – Specifies process term (e.g., ‘boundary_forcing’, ‘vertical_diffusion’, etc.). Use self.available_processes() to list all terms. If None, defaults to list of all terms.

  • group_processes (bool (default: False)) – Specify whether process terms are summed to categories forcing and diffusion.

  • sum_components (bool (default: True)) – Specify whether heat and salt tendencies are summed together (True) or kept separated (False).

  • bins (array like (default None)) – np.array with lambda values specifying the edges for each bin. If not specified, array will be automatically derived from the scalar field of lambda (e.g., temperature).

  • mask (xr.DataArray (default: None)) – Boolean region mask (with same X and Y grid dimensions as grid._ds variables). If None, generate an all-True mask for domain-wide calculations.

Returns:

transformations – Dataset containing components of water mass transformations, possibly grouped as specified by the arguments.

Return type:

xarray.Dataset

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> wmt_rates = wmt.map_transformations("heat", term=["diffusion"])

See also

self.transformations_from_hlamdot, self.map_transformations

lambdas(lambda_key=None)#

Return dictionary of desired lambdas.

Parameters:

lambda_key (str or list of str (default: None)) – Lambda(s) of interest. If None, return a list of all available lambdas.

Return type:

list

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> wmt.lambdas()
['thetao', 'so', 'sigma0', 'sigma1', 'sigma2', 'sigma3', 'sigma4']
map_transformations(lambda_name, term=None, group_processes=False, sum_components=True, bins=None, mask=None)#

Lazily compute column-wise water mass transformations.

Parameters:
  • lambda_name (str) – Specifies lambda (e.g., ‘heat’, ‘salt’, ‘sigma0’, etc.). Use lambdas() for a list of available lambdas.

  • term (str or list of str (default: None)) – Specifies process term (e.g., ‘boundary_forcing’, ‘vertical_diffusion’, etc.). Use self.available_processes() to list all terms. If None, defaults to list of all terms.

  • group_processes (bool (default: False)) – Specify whether process terms are summed to categories forcing and diffusion.

  • sum_components (bool (default: True)) – Specify whether heat and salt tendencies are summed together (True) or kept separated (False).

  • bins (array like (default None)) – np.array with lambda values specifying the edges for each bin. If not specified, array will be automatically derived from the scalar field of lambda (e.g., temperature).

  • mask (xr.DataArray (default: None)) – Boolean region mask (with same X and Y grid dimensions as grid._ds variables). If None, generate an all-True mask for domain-wide calculations.

Returns:

transformations – Dataset containing components of water mass transformations, possibly grouped as specified by the arguments.

Return type:

xarray.Dataset

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> wmt_rates = wmt.map_transformations("heat", term=["diffusion"])

See also

self.transformations_from_hlamdot, self.integrate_transformations

process_names(tracer, term)#

Get a tuple containing the names of variables for the ‘tracer’ and the ‘process’-specific tendency corresponding to a general ‘term’ in the tendency equation.

Parameters:
  • tracer (str) – Supported options: [“heat”, “salt”]

  • term (str) – key for tendency variable in the xbudget_dict

Returns:

names(tracer_name, process)

Return type:

tuple

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> wmt.process_names("heat", "diffusion")
('thetao', 'opottempdiff')
rho_tend(term)#

Get density tendency ‘term’ from underlying heat and salt tendencies.

Parameters:

term (str) – key for tendency variable in the xbudget_dict

Returns:

The two distinct components contributing to the overall density tendency.

Return type:

rho_tend_heat, rho_tend_salt

Example

>>> wmt = WaterMassTransformations(grid, xbudget_dict)
>>> wmt.get_density()
>>> (rho_tend_heat, rho_tend_salt) = wmt.rho_tend("diffusion")
transform_hlamdot_term(lambda_name, term, bins=None, mask=None, integrate=True)#

Lazily compute extensive tendencies and transform them into lambda space along the vertical (“Z”) dimension.

Parameters:
  • lambda_name (str) – Specifies lambda (e.g., ‘heat’, ‘salt’, ‘sigma0’, etc.). Use lambdas() for a list of available lambdas.

  • term (str) – Specifies process term (e.g., ‘boundary_forcing’, ‘vertical_diffusion’, etc.). Use self.available_processes() to see a list of all available terms.

  • bins (array like (default None)) – np.array with lambda values specifying the edges for each bin. If not specified, array will be automatically derived from the scalar field of lambda (e.g., temperature).

  • mask (xr.DataArray (default: None)) – Boolean region mask (with same X and Y grid dimensions as grid._ds variables). If None, generate an all-True mask for domain-wide calculations.

  • integrate (bool (default: True)) – Whether to integrate the result over (“X”, “Y”) dimensions.

Returns:

transformed_term – Data array containing transformed tendencies.

Return type:

xarray.DataArray

See also

self.transformations_from_hlamdot, self.integrate_transformations

transformations_from_hlamdot(lambda_name, term=None, bins=None, mask=None, integrate=True)#

Lazily compute extensive tendencies, transform them into lambda space along the vertical (“Z”) dimension, and optionally integrate along the horizontal dimensions (“X”, “Y”).

Parameters:
  • lambda_name (str) – Specifies lambda (e.g., ‘heat’, ‘salt’, ‘sigma0’, etc.). Use lambdas() for a list of available lambdas.

  • term (str or list of str (default: None)) – Specifies process term (e.g., ‘boundary_forcing’, ‘vertical_diffusion’, etc.). Use self.available_processes() to list all terms. If None, defaults to list of all terms.

  • bins (array like (default None)) – np.array with lambda values specifying the edges for each bin. If not specified, array will be automatically derived from the scalar field of lambda (e.g., temperature).

  • mask (xr.DataArray (default: None)) – Boolean region mask (with same X and Y grid dimensions as grid._ds variables). If None, generate an all-True mask for domain-wide calculations.

  • integrate (bool (default: True)) – Whether to integrate the result over (“X”, “Y”) dimensions.

Returns:

transformations – Dataset containing components of water mass transformations, possibly grouped as specified by the arguments.

Return type:

xarray.Dataset

See also

self.integrate_transformations, self.map_transformations