Wrapper Functions

wrapper.barotropic_eqlat_lwa(ylat, vort, area, dmu, n_points, planet_radius=falwa.constant.EARTH_RADIUS)[source]

Compute equivalent latitude and wave activity on a barotropic sphere.

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • vort (numpy.ndarray) – 2-d numpy array of vorticity values; dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • dmu (numpy.array) – 1-d numpy array of latitudinal differential length element (e.g. dmu = planet_radius * cos(lat) d(lat)). Size = nlat.

  • n_points (int, default None) – Analysis resolution to calculate equivalent latitude. If input as None, it will be initialized as nlat.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

Returns

  • qref (numpy.array) – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).

  • lwa_result (numpy.ndarray) –

    2-d numpy array of local wave activity values;

    dimension = [nlat_s x nlon]

Examples

Example Notebook

wrapper.barotropic_input_qref_to_compute_lwa(ylat, qref, vort, area, dmu, planet_radius=falwa.constant.EARTH_RADIUS)[source]

This function computes LWA based on a prescribed Qref instead of Qref obtained from the vorticity field on a barotropic sphere.

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • qref (numpy.array) – 1-d numpy array of prescribed reference value of vorticity at each latitude; dimension = nlat.

  • vort (numpy.ndarray) – 2-d numpy array of vorticity values; dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • dmu (numpy.array) – 1-d numpy array of latitudinal differential length element (e.g. dmu = planet_radius * cos(lat) d(lat)). Size = nlat.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

Returns

lwa_result – 2-d numpy array of local wave activity values; dimension = [nlat_s x nlon]

Return type

numpy.ndarray

wrapper.eqvlat_bracket_hemispheric(ylat, vort, area, nlat_s=None, n_points=None, planet_radius=falwa.constant.EARTH_RADIUS, vgrad=None)[source]

Compute equivalent latitude and <…>_Q in Nakamura and Zhu (2010) in a hemispheric domain.

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • vort (numpy.ndarray) – 2-d numpy array of vorticity values; dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • nlat_s (int, default None) – The index of grid point that defines the extent of hemispheric domain from the pole. If input as None, it will be initialize as nlat // 2.

  • n_points (int, default None) – Analysis resolution to calculate equivalent latitude. If input as None, it will be initialized as nlat_s.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

  • vgrad (numpy.ndarray, optional) – 2-d numpy array of laplacian (or higher-order laplacian) values; dimension = (nlat, nlon)

Returns

  • q_part (numpy.ndarray) – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).

  • brac (numpy.ndarray or None) – 1-d numpy array of averaged vgrad in the square bracket. If vgrad = None, brac = None.

wrapper.eqvlat_hemispheric(ylat, vort, area, nlat_s=None, n_points=None, planet_radius=falwa.constant.EARTH_RADIUS)[source]

Compute equivalent latitude in a hemispheric domain.

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • vort (numpy.ndarray) – 2-d numpy array of vorticity values; dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • nlat_s (int, default None) – The index of grid point that defines the extent of hemispheric domain from the pole. If input as None, it will be initialize as nlat // 2.

  • n_points (int, default None) – Analysis resolution to calculate equivalent latitude. If input as None, it will be initialized as nlat_s.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

Returns

q_part – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).

Return type

numpy.ndarray

wrapper.qgpv_eqlat_lwa(ylat, vort, area, dmu, nlat_s=None, n_points=None, planet_radius=falwa.constant.EARTH_RADIUS)[source]

Compute equivalent latitutde qref and local wave activity lwa_result based on Quasi-geostrophic potential vorticity field vort at a pressure level as outlined in Huang and Nakamura (2017).

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • vort (numpy.ndarray) – 2-d numpy array of Quasi-geostrophic potential vorticity field; dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • dmu (numpy.array) – 1-d numpy array of latitudinal differential length element (e.g. dmu = planet_radius * cos(lat) d(lat)). Size = nlat.

  • nlat_s (int, default None) – The index of grid point that defines the extent of hemispheric domain from the pole. If input as None, it will be initialize as nlat // 2.

  • n_points (int, default None) – Analysis resolution to calculate equivalent latitude. If input as None, it will be initialized as nlat_s.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

Returns

  • qref (numpy.ndarray) – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).

  • lwa_result (numpy.ndarray) –

    2-d numpy array of local wave activity values;

    dimension = [nlat_s x nlon]

Examples

Example Notebook

wrapper.qgpv_eqlat_lwa_ncforce(ylat, vort, ncforce, area, dmu, nlat_s=None, n_points=None, planet_radius=falwa.constant.EARTH_RADIUS)[source]

Compute equivalent latitutde qref, local wave activity lwa_result and non-conservative force on wave activity capsigma based on Quasi- geostrophic potential vorticity field vort at a pressure level as outlined in Huang and Nakamura (2017).

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • vort (numpy.ndarray) – 2-d numpy array of Quasi-geostrophic potential vorticity field; dimension = (nlat, nlon).

  • ncforce (numpy.ndarray) – 2-d numpy array of non-conservative force field (i.e. theta in NZ10(a) in equation (23a) and (23b)); dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • dmu (numpy.array) – 1-d numpy array of latitudinal differential length element (e.g. dmu = planet_radius * cos(lat) d(lat)). Size = nlat.

  • nlat_s (int, default None) – The index of grid point that defines the extent of hemispheric domain from the pole. If input as None, it will be initialize as nlat // 2.

  • n_points (int, default None) – Analysis resolution to calculate equivalent latitude. If input as None, it will be initialized as nlat_s.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

Returns

  • qref (numpy.ndarray) – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).

  • lwa_result (numpy.ndarray) – 2-d numpy array of local wave activity values; dimension = (nlat, nlon).

  • capsigma (numpy.ndarray) – 2-d numpy array of non-conservative force contribution value; dimension = (nlat, nlon).

wrapper.qgpv_eqlat_lwa_options(ylat, vort, area, dmu, nlat_s=None, n_points=None, vgrad=None, ncforce=None, planet_radius=falwa.constant.EARTH_RADIUS)[source]

Compute equivalent latitutde qref, local wave activity lwa_result and non-conservative force on wave activity capsigma based on Quasi- geostrophic potential vorticity field vort at a pressure level as outlined in Huang and Nakamura (2017).

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • vort (numpy.ndarray) – 2-d numpy array of Quasi-geostrophic potential vorticity field; dimension = (nlat, nlon).

  • ncforce (numpy.ndarray) – 2-d numpy array of non-conservative force field (i.e. theta in NZ10(a) in equation (23a) and (23b)); dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • dmu (numpy.array) – 1-d numpy array of latitudinal differential length element (e.g. dmu = planet_radius * cos(lat) d(lat)). Size = nlat.

  • nlat_s (int, default None) – The index of grid point that defines the extent of hemispheric domain from the pole. If input as None, it will be initialize as nlat // 2.

  • n_points (int, default None) – Analysis resolution to calculate equivalent latitude. If input as None, it will be initialized as nlat_s.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

Returns

  • return_dict (dictionary) – A dictionary that consist of the 4 computed outputs listed below.

  • (1) qref (numpy.ndarray) – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).

  • (2) brac_result (numpy.ndarray) – 1-d numpy array of <…>_Q(y) in NZ10 where latitude y is given by ylat; dimension = (nlat).

  • (3) lwa_result (numpy.ndarray) – 2-d numpy array of local wave activity values; dimension = (nlat, nlon).

  • (4) capsigma (numpy.ndarray) – 2-d numpy array of non-conservative force contribution value; dimension = (nlat, nlon).

wrapper.qgpv_input_qref_to_compute_lwa(ylat, qref, vort, area, dmu, nlat_s=None, planet_radius=falwa.constant.EARTH_RADIUS)[source]

Compute equivalent latitutde qref and local wave activity lwa_result based on Quasi-geostrophic potential vorticity field vort at a pressure level as outlined in Huang and Nakamura (2017). This function computes lwa based on a prescribed qref instead of qref obtained from the QGPV field.

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • qref (numpy.ndarray) – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).

  • vort (numpy.ndarray) – 2-d numpy array of Quasi-geostrophic potential vorticity field; dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • dmu (numpy.array) – 1-d numpy array of latitudinal differential length element (e.g. dmu = planet_radius * cos(lat) d(lat)). Size = nlat.

  • nlat_s (int, default None) – The index of grid point that defines the extent of hemispheric domain from the pole. If input as None, it will be initialize as nlat // 2.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

Returns

lwa_result – 2-d numpy array of local wave activity values; dimension = (nlat, nlon).

Return type

numpy.ndarray

wrapper.theta_lwa(ylat, theta, area, dmu, nlat_s=None, n_points=None, planet_radius=falwa.constant.EARTH_RADIUS)[source]

Compute the surface wave activity B based on surface potential temperature. See Nakamura and Solomon (2010a) for details.

Parameters
  • ylat (numpy.array) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.

  • theta (numpy.ndarray) – 2-d numpy array of surface potential temperature field; dimension = (nlat, nlon).

  • area (numpy.ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).

  • dmu (numpy.array) – 1-d numpy array of latitudinal differential length element (e.g. dmu = planet_radius * cos(lat) d(lat)). Size = nlat.

  • nlat_s (int, default None) – The index of grid point that defines the extent of hemispheric domain from the pole. If input as None, it will be initialize as nlat // 2.

  • planet_radius (float, default 6.378e+6) – radius of spherical planet of interest consistent with input ‘area’.

Returns

  • qref (numpy.array) – 1-d numpy array of value reference potential temperature Theta(y) approximated by box counting method, where latitude y is given by ylat; dimension = (nlat).

  • lwa_result (numpy.ndarray) – 2-d numpy array of local surface wave activity values; dimension = (nlat, nlon).