STJ_PV package

Submodules

Runtime Configuration: run_stj

Run STJ: Main module “glue” that connects Subtropical Jet Metric calc, plot and diags.

python run_stj.py –help

Usage

usage: run_stj.py [-h] [–sample] [–sens] [–warn] [–file FILE] [–ys YS] [–ye YE]

Find the sub-tropical jet

optional arguments:
-h, --help

show this help message and exit

--sample

Perform a sample run

--sens

Perform a parameter sensitivity run

--warn

Show all warning messages

--file FILE

Configuration file path

--ys YS

Start Year

--ye YE

End Year

Authors: Penelope Maher, Michael Kelleher

class STJ_PV.run_stj.JetFindRun(config_file=None)[source]

Bases: object

Class containing properties about an individual attempt to find the subtropical jet.

Parameters
configstring, optional

Location of YAML-formatted configuration file, default None

Attributes
data_sourcestring

Path to input data configuration file

configdict

Dictionary of properties of the run

freqTuple

Output data frequency (time, spatial)

methodSTJMetric

Jet finder type

loglogging.Logger

Debug log

log_setup(self)[source]

Create a logger object with file location from self.config.

run(self, date_s=None, date_e=None, save=True)[source]

Find the jet, save location to a file.

Parameters
date_s, date_edatetime.datetime

Beginning and end dates, optional. If not included, use (Jan 1, self.year_s) and/or (Dec 31, self.year_e)

run_sensitivity(self, sens_param, sens_range, date_s=None, date_e=None)[source]

Perform a parameter sweep on a particular parameter of the JetFindRun.

Parameters
sens_paramstring

Configuration parameter of JetFindRun()

sens_rangeiterable

Range of values of sens_param over which to iterate

date_s, date_edatetime.datetime, optional

Start and end dates, respectively. Optional, defualts to config file defaults

STJ_PV.run_stj.check_config_req(cfg_file, required_keys_all, id_file=True)[source]

Check that required keys exist within a configuration file.

Parameters
cfg_filestring

Path to configuration file

required_keys_alllist

Required keys that must exist in configuration file

Returns
configdict

Dictionary of loaded configuration file

mkeysbool

True if required keys are missing

STJ_PV.run_stj.check_data_config(cfg_file)[source]

Check the settings in a data configuration file.

Parameters
cfg_filestring

Path to configuration file

STJ_PV.run_stj.check_run_config(cfg_file)[source]

Check the settings in a run configuration file.

Parameters
cfg_filestring

Path to configuration file

STJ_PV.run_stj.main(sample_run=True, sens_run=False, cfg_file=None, year_s=1979, year_e=2018)[source]

Run the STJ Metric given a configuration file.

STJ_PV.run_stj.make_parse()[source]

Make command line argument parser with argparse.

Metric Definitions: stj_metric

Calculate the position of the subtropical jet in both hemispheres.

class STJ_PV.stj_metric.STJDavisBirner(props, data)[source]

Bases: STJ_PV.stj_metric.STJMetric

Subtropical jet position metric using the method of Davis and Birner 2016.

Parameters
propsJetFindRun()

Class containing properties about the current search for the STJ

dataInputData()

Input data class containing a year (or more) of required data

Notes

Based on “Climate Model Biases in the Width of the Tropical Belt Parameters”, Davis and Birner (2016)

The logic for the method is to subtract the surface wind and then find the max in the upper level wind.

find_jet(self, shemis=True)[source]

Find the subtropical jet using method from Davis and Birner (2016).

doi:10.1175/JCLI-D-15-0336.1

Parameters
shemislogical, optional

If True, find jet position in Southern Hemisphere, if False, find N.H. jet

find_max_wind_surface(self, uzonal, lat, test_plot=False)[source]

Find most equatorward maximum wind on column maximum wind surface.

Parameters
uzonalarray_like

Zonal mean zonal wind

shemisbool

Flag to indicate hemisphere, true for SH

test_plot(self, lat, max_wind_surface, lat_idx, stj_lat, stj_intens)[source]

Make a test plot for this method.

class STJ_PV.stj_metric.STJKangPolvani(props, data)[source]

Bases: STJ_PV.stj_metric.STJMetric

Subtropical jet position metric: Kang and Polvani 2010.

Parameters
propsJetFindRun()

Class containing properties about the current search for the STJ

dataInputData()

Input data class containing a year (or more) of required data

find_jet(self, shemis=True)[source]

Find the subtropical jet using input parameters.

Parameters
shemislogical, optional

If True, find jet position in Southern Hemisphere, if False, find N.H. jet

find_single_jet(self, u_high, u_low, lat, sign_change)[source]

Get streamfunction inflection point with the largest zonal-mean zonal wind shear.

get_flux_div(self, lats)[source]

Calculate the meridional eddy momentum flux divergence.

get_jet_lat(self, del_f, hem_s)[source]

Find the 200hpa zero crossing of meridional eddy momentum flux divergence.

get_jet_loc(self, data, expected_lat, lat)[source]

Get jet location based on sign changes of Del(f).

loop_jet_lat(self, data, expected_lat, lat)[source]

Get jet location at multiple times.

class STJ_PV.stj_metric.STJMaxWind(props, data)[source]

Bases: STJ_PV.stj_metric.STJMetric

Subtropical jet position metric: maximum zonal mean zonal wind on a pressure level.

Parameters
propsJetFindRun()

Class containing properties about the current search for the STJ

dataInputData()

Input data class containing a year (or more) of required data

find_jet(self, shemis=True)[source]

Find the subtropical jet using input parameters.

Parameters
shemislogical, optional

If True, find jet position in Southern Hemisphere, if False, find N.H. jet

class STJ_PV.stj_metric.STJMetric(name=None, data=None, props=None)[source]

Bases: object

Generic Class containing Sub Tropical Jet metric methods and attributes.

Parameters
namestring

Name of this type of method

dataInputData

Object containing required input data

propsJetFindRun

Properties about the current jet finding attempt, including log file

append(self, other)[source]

Append another metric’s intensity, latitude, and theta positon to this one.

compute(self)[source]

Compute all dask arrays in self.out_data.

save_jet(self)[source]

Save jet position to file.

set_hemis(self, shemis)[source]

Select hemisphere data.

This function sets self.hemis to be an length N list of slices such that only the desired hemisphere is selected with N-D data (e.g. uwind and ipv) along all other axes. It also returns the latitude for the selected hemisphere, an index to select the hemisphere in output arrays, and the extrema function to find min/max of PV derivative in a particular hemisphere.

Parameters
shemisboolean

If true - use southern hemisphere data, if false, use NH data

Returns
extremafunction

Function used to identify extrema in meridional PV gradient, either scipy.signal.argrelmax() if SH, or scipy.signal.argrelmin() for NH

latstuple

(start, end) of latitude, for use with xarray.DataArray.sel(lat=slice(*lats))

hem_sstring

Abbreviation for hemisphere (NH or SH)

class STJ_PV.stj_metric.STJPV(props, data)[source]

Bases: STJ_PV.stj_metric.STJMetric

Subtropical jet position metric using dynamic tropopause on isentropic levels.

Parameters
propsJetFindRun

Class containing properties about the current search for the STJ

dataInputData

Input data class containing a year (or more) of required data

find_jet(self, shemis=True, debug=False)[source]

Find the subtropical jet using input parameters.

Parameters
shemislogical, optional

If True, find jet position in Southern Hemisphere, If False, find N.H. jet

debuglogical, optional

Enter debug mode if true, returns d(theta) / d(lat) values, polynomial fit, and jet latitude

find_single_jet(self, theta_xpv, lat, ushear, extrema, debug=False)[source]

Find jet location for a 1D array of theta on latitude.

Parameters
theta_xpvarray_like

Theta on PV level as a function of latitude

latarray_like

1D array of latitude same shape as theta_xpv from isolate_pv()

usheararray_like

1D array along latitude axis of maximum surface - troposphere u-wind shear

debugboolean

If True, returns debugging information about how jet position is found, if False (default) returns only jet location

Returns
jet_locint

If debug is False, Index of jet location on latitude axis

jet_loc, jet_loc_all, dtheta, theta_fit, lat, y_s, y_etuple

If debug is True, return lots of stuff TODO: document this better!!

isolate_pv(self, pv_lev)[source]

Get the potential temperature, zonal wind and zonal wind shear for a PV level.

Parameters
pv_levfloat

PV value (for a particular hemisphere, >0 for NH, <0 for SH) on which to interpolate potential temperature and wind

theta_bndstuple, optional

Start and end theta levels to use for interpolation. Default is None, if None, use all theta levels, otherwise restrict so theta_bnds[0] <= theta <= theta_bnds[1]

Returns
theta_xpvarray_like

N-1 dimensional array (where self.data.ipv is N-D) of potential temperature on pv_lev PVU

uwnd_xpvarray_like

N-1 dimensional array (where self.data.uwnd is N-D) of zonal wind on pv_lev PVU

usheararray_like

Wind shear between uwnd_xpv and “surface”, meaning the lowest valid level

select_jet(self, locs, ushear)[source]

Select correct jet latitude given list of possible jet locations.

Parameters
locslist

List of indicies of jet locations

usheararray_like

1D array along latitude axis of maximum surface - troposphere u-wind shear

Returns
jet_locint

Index of the jet location. Between [0, lat.shape[0] - 1]

Notes

  • If the list of locations is empty, return 0 as the location, this is interpreted by find_jet() as missing.

  • If the list of locations is exactly one return that location.

  • Otherwise use the location with maximum zonal wind shear between lowest provided level and the dynamical tropopause.

STJ_PV.stj_metric.get_season(month)[source]

Map month index to index of season [DJF -> 0, MAM -> 1, JJA -> 2, SON -> 3].

STJ_PV.stj_metric.lowest_valid(col)[source]

Given 1-D array find lowest (along axis) valid data.

Parameters
colarray_like

1D array of data

Returns
validcol.dtype

Lowest (in index) valid value in col

Input Data: input_data

Generate or load input data for STJ Metric.

class STJ_PV.input_data.InputData(props, date_s=None, date_e=None)[source]

Bases: object

Contains the relevant input data and routines for an JetFindRun.

Parameters
jet_findJetFindRun()

Object containing properties about the metric calculation to be performed. Used to locate correct files, and variables within those files.

yearint, optional

Year of data to load, not used when all years are in a single file

get_data(self)[source]

Get a single xarray.Dataset of required components for metric.

load_vars = []
write_data(self, out_file=None)[source]

Write netCDF file of the out_data.

class STJ_PV.input_data.InputDataSTJPV(props, date_s=None, date_e=None)[source]

Bases: STJ_PV.input_data.InputData

Contains the relevant input data and routines for an STJPV jet find.

Parameters
jet_findJetFindRun()

Object containing properties about the metric calculation to be performed. Used to locate correct files, and variables within those files.

yearint, optional

Year of data to load, not used when all years are in a single file

get_data(self)[source]

Load and compute required data, return xarray.Dataset.

load_vars = ['uwnd', 'vwnd', 'tair', 'epv', 'ipv']
class STJ_PV.input_data.InputDataUWind(props, date_s=None, date_e=None, vwnd=False)[source]

Bases: STJ_PV.input_data.InputData

Contains the relevant input data and routines for an STJPV jet find.

Parameters
jet_findJetFindRun()

Object containing properties about the metric calculation to be performed. Used to locate correct files, and variables within those files.

yearint, optional

Year of data to load, not used when all years are in a single file

get_data(self)[source]

Load and compute (if needed) U-Wind on selected pressure level.

STJ_PV.input_data.package_data(relpath, file_name)[source]

Get data relative to this installed package. Generally used for the sample data.

General Utilities: utils

Utility functions not specific to subtropical jet finding.

class STJ_PV.utils.NDSlicer(axis, ndim, start=None, stop=None, step=None)[source]

Bases: object

Create an n-dimensional slice list.

Parameters
axisinteger

Axis on which to apply the slice

ndiminteger

Total number of dimensions of array to be sliced

start, stop, stepinteger, optional

Index of beginning, stop and step width of the slice [start:stop:step] default for each is None.

Examples

Create random array, slice it:

x = np.random.randn(5, 3)

# Create slicer equivalent to [1:-1, :]
slc = NDSlicer(0, x.ndim)
print(x)
[[ 0.68470539  0.87880216 -0.45086367]
 [ 1.06804045  0.63094676 -0.76633033]
 [-1.69841915  0.35207064 -0.4582049 ]
 [-0.56431067  0.62833728 -0.04101542]
 [-0.02760744  2.02814338  0.13195714]]
print(x[slc[1:-1]])
[[ 1.06804045  0.63094676 -0.76633033]
 [-1.69841915  0.35207064 -0.4582049 ]
 [-0.56431067  0.62833728 -0.04101542]]
slice(self, start=None, stop=None, step=None)[source]

Create an n-dimensional slice list.

This is a legacy compatibility method, calls self.__getitem__.

Parameters
axisinteger

Axis on which to apply the slice

ndiminteger

Total number of dimensions of array to be sliced

start, stop, stepinteger, optional

Index of beginning, stop and step width of the slice [start:stop:step] default for each is None.

Returns
slicerlist

list of slices such that all data at other axes are kept, one axis is sliced

STJ_PV.utils.convert_radians_latlon(lat, lon)[source]

Convert input lat/lon array to radians if input is degrees, do nothing if radians.

Parameters
latarray_like

ND array of latitude

lonarray_like

ND array of longitude

Returns
latarray_like

ND array of latitude in radians

lonarray_like

ND array of longitude in radians

STJ_PV.utils.diff_cfd(data, axis=-1, cyclic=False)[source]

Calculate centered finite difference on a field along an axis with even spacing.

Parameters
dataarray_like

ND array of data of which to calculate the differences

axisinteger

Axis of data on which differences are calculated

cyclicbool

Flag to indicate whether data is cyclic on axis

Returns
diffarray_like

ND array of central finite differences of data along axis

STJ_PV.utils.diff_cfd_xr(data, dim='lon', cyclic=False)[source]

Calculate centered finite difference along an axis with even spacing.

Parameters
dataxarray.DataArray

N-D array of data of which to calculate the differences

dimstring

Dimension name of data along which differences are calculated

cyclicbool

Flag to indicate whether data is cyclic on axis

Returns
diffxarray.DataArray

N-D array of central finite differences of data along axis

STJ_PV.utils.diffz(data, vcoord, axis=None)[source]

Calculate vertical derivative for data on uneven vertical levels.

Parameters
dataarray_like

N-D array of input data to be differentiated, where data.shape[axis] == vcoord.shape[0]

vcoordarray_like

Vertical coordinate, 1D

axisinteger, optional

Axis where data.shape[axis] == vcoord.shape[0]

Returns
dxdzarray_like

N-D array of d(data)/d(vcoord), same shape as input data

STJ_PV.utils.dlon_dlat(lon, lat, cyclic=True)[source]

Calculate distance along lat/lon axes on spherical grid.

Parameters
latarray_like

ND array of latitude

lonarray_like

ND array of longitude

Returns
dlongarray_like

NLat x Nlon array of horizontal distnaces along longitude axis

dlatgarray_like

NLat x Nlon array of horizontal distnaces along latitude axis

STJ_PV.utils.dth_dp(theta_in, data_in)[source]

Calculate vertical derivative on even (theta) levels.

Parameters
theta_inarray_like, 1D

Vertical coordinate (potential temperature)

data_inarray_like, ND, where data_in.shape[1] == theta_in.shape[0]

Data on theta levels

Returns
out_dataarray_like ND, same as data

Centered finite difference for data_in[:, 1:-1, …], backward/forward for bottom/top boundaries

STJ_PV.utils.find_tropopause_mask(dtdz, d_z, thr=2.0)[source]

Use Reichler et al. 2003 method to calculate tropopause level.

Parameters
dtdzarray_like

1-D array of lapse rate in K/km

d_zarray_like

1-D array (same shape as dtdz) of difference in height between levels

thrfloat

Lapse rate threshold in K/km (WMO definition is 2.0 K/km)

Returns
trop_levelinteger

Index of tropopause level on the vertical axis

STJ_PV.utils.get_tropopause(t_air, pres, thr=2.0, vaxis=1)[source]

Return the tropopause temperature and pressure for WMO tropopause.

Parameters
t_airarray_like

ND array of temperature, where axis 1 is vertical axis

presarray_like

ND array of pressure levels, shape is same as t_air

thrfloat

Lapse rate threshold, default/WMO definition is 2.0 K km^-1

Returns
trop_temp, trop_presarray_like

Temperature and pressure at tropopause level, in (N-1)-D arrays, where dimension dropped is vertical axis, same units as input t_air and pres respectively

STJ_PV.utils.get_tropopause_pres(t_air, pres, thr=2.0)[source]

Return the tropopause temperature and pressure for WMO tropopause.

Parameters
t_airarray_like

ND array of temperature, where axis 1 is vertical axis

presarray_like

1D array of pressure levels, shape is same as t_air.shape[1]

thrfloat

Lapse rate threshold, default/WMO definition is 2.0 K km^-1

Returns
trop_temp, trop_presarray_like

Temperature and pressure at tropopause level, in (N-1)-D arrays, where dimension dropped is vertical axis, same units as input t_air and pres respectively

STJ_PV.utils.get_tropopause_theta(theta_in, pres, thr=2.0)[source]

Return the tropopause temperature and pressure for WMO tropopause.

Parameters
theta_inarray_like

1D array of potential temperature

presarray_like

ND array of pressure levels, with one axis matching shape of theta_in

thrfloat

Lapse rate threshold, default/WMO definition is 2.0 K km^-1

Returns
trop_temp, trop_presarray_like

Temperature and pressure at tropopause level, in (N-1)-D arrays, where dimension dropped is vertical axis, same units as input t_air and pres respectively

STJ_PV.utils.inc_with_z(vcoord, levname)[source]

Find what proportion of vcoord is increasing along the levname axis.

Parameters
vcoordxarray.DataArray

array of vertical coordinate to test

levnamestr

String name of vertical coordinate variable along which to test

Returns
pctincfloat

Percentage of valid points in vcoord which are increasing with increasing index

STJ_PV.utils.interp_nd(lat, theta_in, data, lat_hr, theta_hr)[source]

Perform interpolation on 2-dimensions on up to 4-dimensional numpy array.

Parameters
latarray_like

One dimensional latitude coordinate array, matches a dimension of data

theta_inarray_like

One dimensional theta (vertical) coordinate array, matches a dimension of data

dataarray_like

Data to be interpolated to high-resolution grid

lat_hrarray_like

1-D latitude coordinate array that data is interpolated to

theta_hrarray_like

1-D vertical coordinate array that data is interpolated to

Returns
data_interparray_like

Interpolated data where lat/theta dimensions are interpolated to lat_hr and theta_hr

STJ_PV.utils.inv_theta(thta, pres)[source]

Calculate potential temperature from temperature and pressure coordinate.

Parameters
thtaarray_like

Either ND array of potential temperature same shape as pres, or 1D array of potential temperature where its shape is same as one dimension of pres [in K]

presarray_like

ND array of pressure [in Pa]

Returns
tairarray_like

ND array of air temperature in K, same shape as tair input

STJ_PV.utils.ipv(uwnd, vwnd, tair, pres, lat, lon, th_levels=None)[source]

Calculate isentropic PV on theta surfaces.

Parameters
uwndarray_like

3 or 4-D zonal wind component (t, p, y, x) or (p, y, x)

vwndarray_like

3 or 4-D meridional wind component (t, p, y, x) or (p, y, x)

tairarray_like

3 or 4-D air temperature (t, p, y, x) or (p, y, x)

presarray_like

1D pressure in Pa

latarray_like

1D latitude in degrees

lonarray_like

1D longitude in degrees

th_levelsarray_like, optional

1D Theta levels on which to calculate PV. Defaults to 300K - 500K by 5K.

Returns
ipvarray_like

3 or 4-D isentropic potential vorticity in units of m-2 s-1 K kg-1 (e.g. 10^6 PVU)

p_tharray_like

Pressure on isentropic levels [Pa]

u_tharray_like

Zonal wind on isentropic levels [m/s]

Notes

Interpolation assumes pressure is monotonically increasing.

STJ_PV.utils.ipv_theta(uwnd, vwnd, pres, lat, lon, th_levels)[source]

Calculate isentropic PV on theta surfaces from data on theta levels.

Parameters
uwndarray_like

3 or 4-D zonal wind component (t, p, y, x) or (p, y, x)

vwndarray_like

3 or 4-D meridional wind component (t, p, y, x) or (p, y, x)

presarray_like

3 or 4-D pressure in Pa (t, p, y, x) or (p, y, x)

latarray_like

1D latitude in degrees

lonarray_like

1D longitude in degrees

th_levelsarray_like

1D Theta levels on which to calculate PV

Returns
ipvarray_like

3 or 4-D isentropic potential vorticity in units of m-2 s-1 K kg-1 (e.g. 10^6 PVU)

STJ_PV.utils.lapse_rate(t_air, pres, vaxis=None)[source]

Calculate the lapse rate of temperature in K/km from isobaric temperature data.

Parameters
t_airarray_like

N-D array of temperature on isobaric surfaces in K

presarray_like

1-D array of pressure levels, matching one dimension of t_air, in hPa

Returns
dtdzarray_like

N-D array of lapse rate in K/km matching dimensionality of t_air

d_zarray_like

N-D array of height differences in km between levels, same shape as t_air

STJ_PV.utils.rel_vort(uwnd, vwnd, lat, lon, cyclic=True)[source]

Calculate the relative vorticity given zonal (uwnd) and meridional (vwnd) winds.

Parameters
uwndarray_like

Zonal wind >=2 dimensional (time, level, lat, lon)

vwndarray_like

Meridional wind >=2 dimensional, same dimensions as uwnd

latarray_like

Latitude array, 1 dimensional with lat.shape[0] == uwnd.shape[-2]

lonarray_like

Longitude array, 1 dimensional with lon.shape[0] == uwnd.shape[-1]

cyclicboolean

Flag to indicate if data is cyclic in longitude direction

Returns
rel_vortarray_like

Relative vorticity V = U x V (cross product of U and V wind) If cyclic, returns same dimensions as uwnd & vwnd, if not it is (*vwnd.shape[0:-1], vwnd.shape[-1] - 2)

STJ_PV.utils.theta(tair, pres)[source]

Calculate potential temperature from temperature and pressure coordinate.

Parameters
tairarray_like

ND array of air temperature [in K]

presarray_like

Either ND array of pressure same shape as tair, or 1D array of pressure where its shape is same as one dimension of tair [in Pa]

Returns
thetaarray_like

ND array of potential temperature in K, same shape as tair input

STJ_PV.utils.trop_lev_1d(dtdz, d_z, thr=2.0, return_idx=False)[source]

Given 1D arrays for lapse rate and change in height, and a threshold, find tropopause.

Parameters
dtdzarray_like

1D array of lapse rate (in same units as thr, usually K km^-1

d_zarray_like

1D array of change of height from previous level, same shape as dtdz, in same units as dtdz denominator units

thrfloat

Lapse rate threshold for definition of tropopause, WMO/default is 2.0 K km^-1

return_idxbool

Flag to return index of tropopause level. Default is False, if True output is out_mask, idx if False, output is out_mask

Returns
out_maskarray_like

1D array, same shape as dtdz and d_z, of booleans, True everywhere except the tropopause level

STJ_PV.utils.vinterp(data, vcoord, vlevels)[source]

Perform linear vertical interpolation.

Parameters
dataarray_like (>= 2D)

Array of data to be interpolated

vcoordarray_like (>= 2D, where data.shape[1] == vcoord.shape[1])

Array representing the vertical structure (height/pressure/PV/theta/etc.) of data

vlevelsarray_like (1D)

Levels, in same units as vcoord, to interpolate to

Returns
out_dataarray_like, (data.shape[0], vlevels.shape[0], *data.shape[2:])

Data on vlevels

Examples

data and vcoord of the same shape

Given some u-wind data uwnd that is on (time, pressure, lat, lon), a potential temperature field theta on the same dimensions (t, p, lat, lon), and a few selected potential temperature levels to which u-wind will be interpolated

>>> th_levs = [275, 300, 325, 350, 375]
>>> print(uwnd.shape)
(365, 17, 73, 144)
>>> print(theta.shape)
(365, 17, 73, 144)
uwnd_theta = vinterp(uwnd, theta, th_levs)
>>> print(uwnd_theta.shape)
(365, 5, 73, 144)

This results in u-wind on the requested theta surfaces

data is 1D vcoord is 4D

Given some potential temperature levels th_levs = [275, 300, 325, 350, 375] and a 4D array of potential vorticity ipv on (time, theta, lat, lon), and several selected PV levels

>>> pv_levs = [1.0, 2.0, 3.0]
>>> print(th_levs.shape)
(5, )
>>> print(ipv.shape)
(365, 5, 73, 144)
>>> theta_pv = vinterp(th_levs, ipv, pv_levs)
>>> print(theta_pv.shape)
(365, 3, 73, 144)

This gives potential temperature on potential vorticity surfaces

data is 4D vcoord is 1D

Given some pressure levels levs = [1000., 900., 850., 500.] and a 4D array of potential vorticity epv on (time, pressure, lat, lon), and several selected pressure levels

>>> pres_new = [950.0, 800.0, 700.0]
>>> print(levs.shape)
(4, )
>>> print(epv.shape)
(365, 4, 73, 144)
>>> epv_intrp = vinterp(epv, levs, pres_new)
>>> print(epv_interp.shape)
(365, 3, 73, 144)

This gives potential vorticity on new pressure surfaces.

STJ_PV.utils.xr_dlon_dlat(data, vlon='lon', vlat='lat', cyclic=True)[source]

Calculate distance on lat/lon axes on spherical grid for xarray.DataArray.

Parameters
dataxarray.DataArray

DataArray with latitude and longitude coordinates

vlon, vlatstring

Variable names of latitude and longitude

Returns
dlongxarray.DataArray

NLat x Nlon array of horizontal distnaces along longitude axis in m

dlatgxarray.DataArray

NLat x Nlon array of horizontal distances along latitude axis in m

STJ_PV.utils.xr_inv_theta(thta, pvar='level')[source]

Calculate potential temperature from temperature and pressure coordinate.

Parameters
thtaxarray.DataArray

ND array of potential temperature [in K]

pvarstring

Pressure level coordinate name for tair

Returns
tairxarray.DataArray

ND array of air temperature in K, same shape as thta input

STJ_PV.utils.xr_rel_vort(uwnd, vwnd, dimvars, cyclic=True)[source]

Calculate the relative vorticity given zonal (u) and meridional (v) winds.

Parameters
uwndxarray.DataArray

Array of Zonal wind (with at least lat / lon dimensions)

vwndxarray.DataArray

Array of Meridional wind with same dimensions as uwnd

cyclicboolean

Flag to indicate if data is cyclic in longitude direction

Returns
rel_vortxarray.DataArray

Relative vorticity V = U x V (cross product of U and V wind) with same dimensions as uwnd & vwnd

STJ_PV.utils.xrdiffz(data, vcoord, dim='lev')[source]

Calculate vertical derivative for data on uneven vertical levels.

Parameters
dataxarray.DataArray

N-D array of input data to be differentiated, where data.shape[axis] == vcoord.shape[0]

vcoordxarray.DataArray

Vertical coordinate, 1D

dimstring

Vertical dimension name, must exist in both data and vcoord

Returns
dxdzxarray.DataArray

N-D array of d(data)/d(vcoord), same shape as input data

STJ_PV.utils.xripv(uwnd, vwnd, tair, dimvars=None, th_levels=None)[source]

Calculate isentropic PV on theta surfaces from xarray.DataArray.

Parameters
uwndxarray.DataArray

3 or 4-D zonal wind component (t, p, y, x) or (p, y, x)

vwndxarray.DataArray

3 or 4-D meridional wind component (t, p, y, x) or (p, y, x)

tairxarray.DataArray

3 or 4-D air temperature (t, p, y, x) or (p, y, x)

dimvarsdict

Mapping of variable names for standard coordinates. This will default to ‘lev’ -> ‘level’, ‘lat’ -> ‘lat’, ‘lon’ -> ‘lon’

th_levelsarray_like, optional

1D array of Theta levels on which to calculate PV. Defaults to 300K - 500K by 5K.

Returns
ipvxarray.DataArray

3 or 4-D isentropic potential vorticity in units of m-2 s-1 K kg-1 (e.g. 10^6 PVU)

p_thxarray.DataArray

Pressure on isentropic levels [Pa]

u_thxarray.DataArray

Zonal wind on isentropic levels [m/s]

Notes

Interpolation assumes pressure is monotonically increasing.

STJ_PV.utils.xripv_theta(uwnd, vwnd, pres, dimvars)[source]

Calculate isentropic PV on theta surfaces from data on theta levels.

Parameters
uwndxarray.DataArray

3 or 4-D zonal wind component (t, theta, y, x) or (theta, y, x)

vwndxarray.DataArray

3 or 4-D meridional wind component (t, theta, y, x) or (theta, y, x)

presxarray.DataArray

3 or 4-D pressure in Pa (t, theta, y, x) or (theta, y, x)

dimvarsdict

Mapping of variable names for standard coordinates. This will default to ‘theta’ -> ‘theta’, ‘lat’ -> ‘lat’, ‘lon’ -> ‘lon’

Returns
ipvxarray.DataArray

3 or 4-D isentropic potential vorticity in units of m-2 s-1 K kg-1 (e.g. 10^6 PVU)

STJ_PV.utils.xrtheta(tair, pvar='level')[source]

Calculate potential temperature from temperature and pressure coordinate.

Parameters
tairxarray.DataArray

ND array of air temperature [in K]

pvarstring

Pressure level coordinate name for tair

Returns
thetaxarray.DataArray

ND array of potential temperature in K, same shape as tair input

STJ_PV.utils.xrvinterp(data, vcoord, vlevs, levname, newlevname)[source]

Perform vertical interpolation for several levels for an xarray.DataArray.

Parameters
dataxarray.DataArray (>= 2D)

array of data to be interpolated

vcoordxarray.DataArray

array representing the vertical structure (height/pressure/PV/theta/etc.) of data. This should be >= 2D, and data[levname].shape == vcoord[levname].shape

vlevsarray_like (1D)

Levels, in same units as vcoord, to interpolate to

levnamestring

Name of the vertical level coordinate variable upon which to interpolate

newlevnamestring

Name of new vertical level coordinate variable

Returns
out_dataarray_like, (data.shape[0], vlevs.shape[0], *data.shape[2:])

Data on vlevels

Notes

If the input vertical coordinate data is increasing / decreasing with height in different places (e.g. potential vorticity across hemispheres), mask data and compute each part separately then add them together (if appropriate) or combine them, or whatever you’d like to do with them. We’re not the boss of you :)

Eddy Kinetic Energy Terms: eddy_terms

Compute terms related to Eddy Kinetic Energy.

class STJ_PV.eddy_terms.Kinetic_Eddy_Energies(uwnd, vwnd, cfg)[source]

Bases: object

calc_momentum_flux(self, integration_top=None, rh=None, tau=None)[source]

Calculate meridional change in eddy momentum flux.

Notes

\(S = \nabla\cdot\overline{F} = \frac{-1}{\cos(\phi)}\frac{\partial}{\partial phi}(\cos^2(\phi) [\overline{u'v'}])\)

get_components(self, zonal=True, time=True)[source]

Get time and zonal anomalies and means for each component.

Comparison Utilities

Compare Two Runs

Script to compare two or more runs of STJ Find.

class STJ_PV.compare_two_runs.FileDiag(info, opts_hem=None, file_path=None)[source]

Bases: object

Contains information about an STJ metric output file in a DataFrame.

append_metric(self, other)[source]

Append the DataFrame attribute (self.lats) to another FileDiag’s DataFrame.

make_dframe(self)[source]

Creates dataframe from input netCDF / xarray.

time_subset(self, other)[source]

Make two fds objects have matching times.

STJ_PV.compare_two_runs.main()[source]

Selects two files to compare, loads and plots them.

Compare Runs Map

Compare two STJ metrics. Plot limited timeseries and a map.

STJ_PV.compare_runs_map.draw_map_lines(pmap, axis)[source]

Draw lat/lon and coast lines on a Basemap.

Parameters
pmapmpl_toolkits.basemap.Basemap

A Basemap onto which lines will be drawn

axismatplotlib.axes.Axes

Axis reference for pmap

Returns
None
STJ_PV.compare_runs_map.get_pvgrad_pos(date)[source]

Get STJ position at all longitudes for PVGrad method.

Parameters
datedatetime.datetime

Selected date to compute STJ metric

Returns
jet_latnumpy.ndarray

An array of jet latitude locations (hemisphere, time, lon)

STJ_PV.compare_runs_map.main(width=174, figscale=1.0, extn='png')[source]

Load data, make plots.

STJ_PV.compare_runs_map.plot_annotations(fig, axes, cfill)[source]

Annotate the map and line plots.

STJ_PV.compare_runs_map.plot_labels(fig, figscale)[source]

Put a, b, c, … on plots.

Module contents