List of metrics¶
|
TropD Eddy Driven Jet (EDJ) Metric |
|
TropD Outgoing Longwave Radiation (OLR) Metric |
|
TropD Precipitation Minus Evaporation (PE) Metric |
|
TropD Mass Streamfunction (PSI) Metric |
|
TropD Sea-level Pressure (PSL) Metric |
|
TropD Subtropical Jet (STJ) Metric |
|
TropD Tropopause Break (TPB) Metric |
|
TropD Near-Surface Zonal Wind (UAS) Metric |
Description of metrics¶
- skytropd.metrics.Shah_2020_1sigma(tracer: ndarray, lat: ndarray, zonal_mean_tracer=False) ndarray[source]¶
Computes the one-sigma width from 3-D tracer data Reference: Shah et al., JGR-A, 2020 https://doi.org/10.1029/2020JD033081
- Parameters:
tracer (numpy.ndarray (...,lat[, lon])) – N-dimensional array of tracer data for computing 1:math:sigma-width. If
zonal_mean_tracer=False, the last dimension should correspond to the longitude axis, otherwise it should be the latitude axislatitude (numpy.ndarray (lat,)) – array of latitudes in degrees (1-D). If not increasing, it will be sorted
zonal_mean_tracer (bool, optional) – whether the input tracer data is zonally symmetric (True) or zonally-varying (False) (i.e., has a longitude dimension), by default False
- Returns:
tracer_sigma_lat – N-2(N-1 for
zonal_mean_tracer=True) dimenional array of 1:math:sigma-widths in the NH- Return type:
numpy.ndarray
- skytropd.metrics.Shah_2020_GWL(tracer: ndarray, lat: ndarray, zonal_mean_tracer=False) ndarray[source]¶
Computes the gradient-weighted latitude (GWL) from tracer data Reference: Shah et al., JGR-A, 2020 https://doi.org/10.1029/2020JD033081
- Parameters:
tracer (numpy.ndarray (..., lat [, lon])) – N-dimensional array of tracer data for computing gradient. If
zonal_mean_tracer=False, the last dimension should correspond to the longitude axis, otherwise it should be the latitude axislatitude (numpy.ndarray (lat,)) – latitude array in degrees
zonal_mean_tracer (bool, optional) – whether the input tracer data is zonally symmetric (True) or zonally-varying (False) (i.e., has a longitude dimension), by default False
- Returns:
tracer_steep_lat – N-2(N-1 for
zonal_mean_tracer=True) dimenional array of GWL latitudes in the SH- Return type:
numpy.ndarray
- skytropd.metrics.TropD_Metric_EDJ(U: ndarray, lat: ndarray, lev: Optional[ndarray] = None, method: str = 'peak', n_fit: int = 1, **maxlat_kwargs) ndarray[source]¶
TropD Eddy Driven Jet (EDJ) Metric
Latitude of maximum of the zonal wind at the level closest to 850 hPa
- Parameters:
U (numpy.ndarray (dim1, ..., lat[, lev])) – N-dimensional array of zonal-mean zonal wind data. Also accepts surface/ 850hPa wind
lat (numpy.ndarray (lat,)) – latitude array
lev (numpy.ndarray, optional (lev,)) – vertical level array in hPa, used to find wind closest to 850hPa. if not provided, last axis is assumed to be lat
method ({"peak", "max", "fit"}, optional) –
Method for determining latitude of maximum zonal wind, by default “peak”:
- ”peak”: Latitude of the maximum of the zonal wind at the level closest
to 850hPa (smoothing parameter
n=30)
- ”max”: Latitude of the maximum of the zonal wind at the level closest to
850hPa (smoothing parameter
n=6)
- ”fit”: Latitude of the maximum of the zonal wind at the level closest to
850hPa using a quadratic polynomial fit of data from grid points surrounding the grid point of the maximum
n_fit (int, optional) – used when
method="fit", determines the number of points around the max to use while fitting, by default 1**kwargs (optional) –
additional keyword arguments for
TropD_Calculate_MaxLat()(not used formethod="fit")- nint, optional
Rank of moment used to calculate the location of max, e.g.,
n=1,2,4,6,8,..., by default 6 ifmethod="max", 30 ifmethod="peak"
- Returns:
Phi (numpy.ndarray (dim1, …, dimN-2[, dimN-1])) – N-2(N-1) dimensional latitudes of the EDJ
Umax (numpy.ndarray (dim1, …, dimN-2), conditional) – N-2 dimensional EDJ strength, returned if
method="fit"
- skytropd.metrics.TropD_Metric_OLR(olr: ndarray, lat: ndarray, method: str = '250W', Cutoff: Optional[float] = None, **maxlat_kwargs) ndarray[source]¶
TropD Outgoing Longwave Radiation (OLR) Metric
Latitude of maximum OLR or first latitude poleward of maximum where OLR reaches crosses a predefined cutoff
- Parameters:
olr (numpy.ndarray (dim1, ..., lat)) – zonal mean TOA olr (positive)
lat (numpy.ndarray (lat,)) – latitude array
method ({"250W", "20W", "cutoff", "10Perc", "max", "peak"}, optional) –
Method for determining the OLR maximum/threshold, by default “250W”:
- ”250W”: the 1st latitude poleward of the tropical OLR max in each hemisphere
where OLR crosses \(250W/m^2\)
- ”20W”: the 1st latitude poleward of the tropical OLR max in each hemisphere
where OLR crosses [the tropical OLR max minus \(20W/m^2\)]
- ”cutoff”: the 1st latitude poleward of the tropical OLR max in each
hemisphere where OLR crosses a specified cutoff value
- ”10Perc”: the 1st latitude poleward of the tropical OLR max in each
hemisphere where OLR is 10% smaller than the tropical OLR max
- ”max”: the latitude of the tropical olr max in each hemisphere with
smoothing parameter
n=6
- ”peak”: the latitude of maximum of tropical olr in each hemisphere with
smoothing parameter
n=30
Cutoff (float, optional) – if
method="cutoff", specifies the OLR cutoff value in \(W/m^2\)**kwargs (optional) –
additional keyword arguments for
TropD_Calculate_MaxLat()- nint, optional
Rank of moment used to calculate the location of max, e.g.,
n=1,2,4,6,8,..., by default 6 ifmethod="max", 30 ifmethod="peak"
- Returns:
Phi – N-1 dimensional latitudes of the near-equator OLR max/threshold crossing
- Return type:
numpy.ndarray (dim1, …, dimN-1)
- skytropd.metrics.TropD_Metric_PE(pe: ndarray, lat: ndarray, method: str = 'zero_crossing', lat_uncertainty: float = 0.0) ndarray[source]¶
TropD Precipitation Minus Evaporation (PE) Metric
Find the first zero crossing of zonal-mean precipitation minus evaporation poleward of the subtropical minimum
- Parameters:
pe (numpy.ndarray (dim1, ..., lat,)) – zonal-mean precipitation minus evaporation
lat (numpy.ndarray (lat,)) – latitude array
method ({"zero_crossing"}, optional) –
Method to compute the zero crossing for precipitation minus evaporation, by default “zero_crossing”.
- ”zero_crossing”: the first latitude poleward of the subtropical P-E min
where P-E changes from negative to positive.
lat_uncertainty (float, optional) – The minimal distance allowed between adjacent zero crossings in degrees, by default 0.0
- Returns:
Phi – N-1 dimensional latitudes of the 1st subtropical P-E zero crossing
- Return type:
numpy.ndarray (dim1, …, dimN-1)
- skytropd.metrics.TropD_Metric_PSI(V: ndarray, lat: ndarray, lev: ndarray, method: str = 'Psi_500', lat_uncertainty: float = 0.0, field_type: str = 'V', threshold: Optional[float] = 0.1) ndarray[source]¶
TropD Mass Streamfunction (PSI) Metric
Latitude of the subtropical zero crossing of the meridional mass streamfunction. The input field can be either zonal-mean meridional wind or a precomputed mass streamfunction.
- Parameters:
V (numpy.ndarray (dim1, ..., lat, lev) or (dim1, ..., lev, lat)) – N-dimensional zonal-mean meridional wind when
field_type="V"or a precomputed mass streamfunction whenfield_type="PSI". The final two axes may be ordered either as(..., lat, lev)or(..., lev, lat).lat (numpy.ndarray (lat,)) – latitude array
lev (numpy.ndarray (lev,)) – vertical level array in hPa
method ({"Psi_500", "Psi_500_10Perc", "Psi_<p1>_<p2>", "Psi_<p1>_<p2>_<x>Perc", "Psi_<p1>_<p2>_<x>Perc_center2d", "Psi_<p1>_<p2>_<x>Perc_profile", "Psi_500_Int", "Psi_Int"},) –
optional –
Method of determining which Psi zero crossing to return, by default “Psi_500”:
”Psi_500”: Zero crossing of the streamfunction (Psi) at 500hPa
- ”Psi_500_10Perc”: Crossing of
thresholdtimes the extremum value of Psi in each hemisphere at the 500hPa level. The historical method name is retained and defaults to
threshold=0.1.
- ”Psi_500_10Perc”: Crossing of
- ”Psi_<p1>_<p2>”: Zero crossing of Psi vertically averaged between the
p1andp2hPa levels, e.g."Psi_300_700"
- ”Psi_<p1>_<p2>_<x>Perc”: Latitude where the layer-mean streamfunction
between
p1andp2hPa decreases tox% of the Hadley-cell-center value from the full(lat, lev)field. Omitting the suffix defaults to thecenter2dvariant below, so"Psi_500_800_5Perc"remains the Hill-style default. In this repository, if the resulting latitude lies poleward of 60N or is not finite, the metric falls back to the same layer’s standard zero-crossing edge"Psi_<p1>_<p2>".
- ”Psi_<p1>_<p2>_<x>Perc_center2d”: As above, with the threshold taken
from the 2-D cell-center maximum.
- ”Psi_<p1>_<p2>_<x>Perc_profile”: As above, but with the threshold
taken from the layer-mean profile value at the cell-center latitude.
”Psi_500_Int”: Zero crossing of the vertically-integrated Psi at 500 hPa
”Psi_Int” : Zero crossing of the column-averaged Psi
The
Psi_<p1>_<p2>_<x>Perc*family follows the Hadley-cell descending-edge interpretation discussed by Hill, S. A., S. Bordoni, J. L. Mitchell, and J. M. Lora, 2025: Interpreting Seasonal and Interannual Hadley Cell Descending Edge Migrations via the Cell-Mean Rossby Number, Journal of Climate, 38, 5505-5520, https://doi.org/10.1175/JCLI-D-24-0678.1.lat_uncertainty (float, optional) – The minimal distance allowed between the adjacent zero crossings, same units as lat, by default 0.0. e.g., for
lat_uncertainty=10, this function will return NaN if another zero crossing is within 10 degrees of the most equatorward zero crossing.field_type ({"V", "PSI"}, optional) – Specifies whether the first argument is zonal-mean meridional wind or a precomputed mass streamfunction, by default
"V"threshold (float or None, optional) – Threshold fraction used for the threshold-crossing HC edge. For
method="Psi_500_10Perc", this is the threshold applied directly. For the other PSI methods, if the zero crossing cannot be found and the result would otherwise be NaN, the metric falls back to the latitude where Psi first crossesthresholdtimesPmaxbetween the subtropical maximum and minimum. Set toNoneto disable this fallback. By default 0.1.
- Returns:
Phi – N-2 dimensional latitudes of the Psi zero crossing
- Return type:
numpy.ndarray (dim1, …, dimN-2)
- skytropd.metrics.TropD_Metric_PSL(ps: ndarray, lat: ndarray, method: str = 'peak', **maxlat_kwargs) ndarray[source]¶
TropD Sea-level Pressure (PSL) Metric
Latitude of maximum of the subtropical sea-level pressure
- Parameters:
ps (np.ndarray (dim1, ..., lat)) – N-dimensional sea-level pressure
lat (np.ndarray (lat,)) – latitude array
method ({"peak", "max"}, optional) –
Method for determining latitude of max PSL, by default “peak”:
- ”peak”: latitude of the maximum of subtropical sea-level pressure (smoothing
parameter
n=30)
- ”max”: latitude of the maximum of subtropical sea-level pressure (smoothing
parameter
n=6)
**kwargs (optional) –
additional keyword arguments for
TropD_Calculate_MaxLat()(not used formethod="fit")- nint, optional
Rank of moment used to calculate the location of max, e.g.,
n=1,2,4,6,8,..., by default 6 ifmethod="max", 30 ifmethod="peak"
- Returns:
Phi – N-1 dimensional latitudes of the subtropical PSL maximum
- Return type:
numpy.ndarray (dim1, …, dimN-1)
- skytropd.metrics.TropD_Metric_STJ(U: ndarray, lat: ndarray, lev: ndarray, method: str = 'adjusted_peak', n_fit: int = 1, **maxlat_kwargs) ndarray[source]¶
TropD Subtropical Jet (STJ) Metric
Latitude of the (adjusted) maximum upper-level zonal-mean zonal wind
- Parameters:
U (numpy.ndarray (dim1,...,lat,lev,)) – N-dimensional zonal-mean zonal wind
lat (numpy.ndarray (lat,)) – latitude array
lev (numpy.ndarray (lev,)) – vertical level array in hPa
method ({"adjusted_peak", "adjusted_max", "core_peak", "core_max", "fit"}, optional) –
Method for determing the latitude of the STJ maximum, by default “adjusted_peak”:
- ”adjusted_peak”: Latitude of maximum (smoothing parameter``n=30``) of [the
zonal wind averaged between 100 and 400 hPa] minus [the zonal mean zonal wind at the level closest to 850hPa], poleward of 10 degrees and equatorward of the Eddy Driven Jet latitude
- ”adjusted_max”: Latitude of maximum (smoothing parameter
n=6) of [the zonal wind averaged between 100 and 400 hPa] minus [the zonal mean zonal wind at the level closest to 850hPa], poleward of 10 degrees and equatorward of the Eddy Driven Jet latitude
- ”adjusted_max”: Latitude of maximum (smoothing parameter
- ”core_peak”: Latitude of maximum (smoothing parameter
n=30) of the zonal wind averaged between 100 and 400 hPa, poleward of 10 degrees and equatorward of 60 degrees
- ”core_peak”: Latitude of maximum (smoothing parameter
- ”core_max”: Latitude of maximum (smoothing parameter
n=6) of the zonal wind averaged between 100 and 400 hPa, poleward of 10 degrees and equatorward of 60 degrees`r`n * “fit”: Latitude of the maximum of [the zonal wind averaged between 100 and 400`r`n hPa] minus [the zonal mean zonal wind at the level closest to 850 hPa]`r`n using a quadratic polynomial fit of data from grid points surrounding`r`n the grid point of the maximum
- ”core_max”: Latitude of maximum (smoothing parameter
n_fit (int, optional) – used when
method="fit", determines the number of points around the max to use while fitting, by default 1**kwargs (optional) –
additional keyword arguments for
TropD_Calculate_MaxLat()(not used formethod="fit")- nint, optional
Rank of moment used to calculate the location of max, e.g.,
n=1,2,4,6,8,..., by default 6 ifmethod="core_max"ormethod="adjusted_max", 30 ifmethod="core_peak"ormethod="adjusted_peak"
- Returns:
Phi (numpy.ndarray (dim1, …, dimN-2)) – N-2 dimensional latitudes of the STJ
Umax (numpy.ndarray (dim1, …, dimN-2), conditional) – N-2 dimensional STJ strength, returned if
method="fit"
- skytropd.metrics.TropD_Metric_TPB(T: ndarray, lat: ndarray, lev: ndarray, method: str = 'max_gradient', Z: Optional[ndarray] = None, Cutoff: float = 15000.0, **maxlat_kwargs) ndarray[source]¶
TropD Tropopause Break (TPB) Metric
Finds the latitude of the tropopause break
- Parameters:
T (numpy.ndarray (dim1, ..., lat, lev)) – N-dimensional temperature array (K)
lat (numpy.ndarray (lat,)) – latitude array
lev (numpy.ndarray (lev,)) – vertical levels array in hPa
method ({"max_gradient", "cutoff", "max_potemp"}, optional) –
Method to identify tropopause break, by default “max_gradient”:
- ”max_gradient”: The latitude of maximal poleward gradient of the tropopause
height
- ”cutoff”: The most equatorward latitude where the tropopause crosses
a prescribed cutoff value
- ”max_potemp”: The latitude of maximal difference between the potential
temperature at the tropopause and at the surface
Z (Optional[numpy.ndarray], optional) – N-dimensional geopotential height array (m), required by
method="cutoff"Cutoff (float, optional) – Geopotential height cutoff (m) that marks the location of the tropopause break, by default 1.5e4, used when
method="cutoff"**kwargs (optional) –
additional keyword arguments for
TropD_Calculate_MaxLat()(not used formethod="cutoff")- nint, optional
Rank of moment used to calculate the location of max, e.g.,
n=1,2,4,6,8,..., by default 6 ifmethod="max_gradient", 30 ifmethod="max_pottemp"
- Returns:
Phi – N-2 dimensional latitudes of the tropopause break
- Return type:
numpy.ndarray (dim1, …, dimN-2)
- skytropd.metrics.TropD_Metric_UAS(U: ndarray, lat: ndarray, lev: Optional[ndarray] = None, method: str = 'zero_crossing', lat_uncertainty: float = 0.0) ndarray[source]¶
TropD Near-Surface Zonal Wind (UAS) Metric
- Parameters:
U (numpy.ndarray (dim1, ..., lat[, lev])) – N-dimensional zonal mean zonal wind array. Also accepts surface wind
lat (numpy.ndarray (lat,)) – latitude array
lev (numpy.ndarray, optional (lev,)) – vertical level array in hPa, required if U has final dimension lev
method ({"zero_crossing"}, optional) –
Method for identifying the surface wind zero crossing, by default “zero_crossing”.
- ”zero_crossing”: the first subtropical latitude where near-surface zonal wind
changes from negative to positive
lat_uncertainty (float, optional) – the minimal distance allowed between adjacent zero crossings in degrees, by default 0.0
- Returns:
Phi – N-2(N-1) dimensional latitudes of the first subtropical zero crossing of near-surface zonal wind
- Return type:
numpy.ndarray (dim1, …, dimN-2[, dimN-1])
- skytropd.metrics.hemisphere_handler(metric_func: Callable[[...], Union[ndarray, Tuple[ndarray, ndarray]]]) Callable[source]¶
Wrapper for metrics to allow one or two hemispheres of data to be passed
- Parameters:
metric_func (Callable) – the skytropd metric function to wrap
- Returns:
wrapped function with same name, docstring, and call signature as the original function, but which operates over each hemisphere independently
- Return type:
Callable