easyclimate.core.windspharm¶
Easy climate top interface for the windspharm
This is the top layer of packaging for the windspharm package.
See also
Spherical harmonic vector wind computations (xarray interface)
Dawson, A. (2016). Windspharm: A High-Level Library for Global Wind Field Computations Using Spherical Harmonics. Journal of Open Research Software, 4(1), e31.DOI: https://doi.org/10.5334/jors.129
Functions¶
|
Calculate the wind speed (magnitude of vector wind). |
Calculate relative vorticity and horizontal divergence. |
|
|
Calculate relative vorticity. |
|
Calculate horizontal divergence. |
|
Calculate planetary vorticity (Coriolis parameter). |
|
Calculate absolute vorticity (sum of relative and planetary vorticity). |
Calculate stream function and velocity potential. |
|
|
Calculate stream function. |
|
Calculate velocity potential. |
|
Calculate irrotational and non-divergent components of the vector wind. |
|
Calculate irrotational (divergent) component of the vector wind. |
|
Calculate non-divergent (rotational) component of the vector wind. |
|
Calculate Rossby wave sources (RWS). |
|
Computes the vector gradient of a scalar field on the sphere. |
Module Contents¶
- easyclimate.core.windspharm.calc_wind_speed(u_data: xarray.DataArray, v_data: xarray.DataArray, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.DataArray¶
Calculate the wind speed (magnitude of vector wind).
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
The wind speed (
xarray.DataArray).- u_data
- easyclimate.core.windspharm.calc_relative_vorticity_and_horizontal_divergence(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.Dataset¶
Calculate relative vorticity and horizontal divergence.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Relative vorticity and horizontal divergence (
xarray.Dataset).- u_data
- easyclimate.core.windspharm.calc_relative_vorticity(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.DataArray¶
Calculate relative vorticity.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Relative vorticity (
xarray.DataArray).- u_data
- easyclimate.core.windspharm.calc_divergence(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.DataArray¶
Calculate horizontal divergence.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Horizontal divergence (
xarray.DataArray).- u_data
- easyclimate.core.windspharm.calc_planetary_vorticity(u_data: xarray.DataArray, v_data: xarray.DataArray, omega: float = 7.292e-05, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.DataArray¶
Calculate planetary vorticity (Coriolis parameter).
Parameters¶
- u_data:
xarray.DataArray. The zonal component of the wind.
- v_data:
xarray.DataArray. The meridional component of vector wind.
- omega:
float. Earth’s angular velocity. The default value if not specified is \(7.292 \times 10^{-5} \mathrm{s^{-1}}\).
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Planetary vorticity (
xarray.DataArray).- u_data:
- easyclimate.core.windspharm.calc_absolute_vorticity(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, omega: float = 7.292e-05, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.DataArray¶
Calculate absolute vorticity (sum of relative and planetary vorticity).
Parameters¶
- u_data:
xarray.DataArray. The zonal component of the wind.
- v_data:
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- omega:
float. Earth’s angular velocity. The default value if not specified is \(7.292 \times 10^{-5} \mathrm{s^{-1}}\).
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Absolute vorticity (
xarray.DataArray).- u_data:
- easyclimate.core.windspharm.calc_streamfunction_and_velocity_potential(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.Dataset¶
Calculate stream function and velocity potential.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Stream function and velocity potential (
xarray.Dataset).- u_data
- easyclimate.core.windspharm.calc_streamfunction(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.DataArray¶
Calculate stream function.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
stream function (
xarray.DataArray).- u_data
- easyclimate.core.windspharm.calc_velocity_potential(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.DataArray¶
Calculate velocity potential.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Velocity potential (
xarray.DataArray).- u_data
- easyclimate.core.windspharm.calc_helmholtz(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.Dataset¶
Calculate irrotational and non-divergent components of the vector wind.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Irrotational and non-divergent components of the vector wind (
xarray.Dataset).- u_data
- easyclimate.core.windspharm.calc_irrotational_component(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.Dataset¶
Calculate irrotational (divergent) component of the vector wind.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Irrotational (divergent) component of the vector wind (
xarray.Dataset).- u_data
- easyclimate.core.windspharm.calc_nondivergent_component(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.Dataset¶
Calculate non-divergent (rotational) component of the vector wind.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Non-divergent (rotational) component of the vector wind (
xarray.Dataset).- u_data
- easyclimate.core.windspharm.calc_rossby_wave_source(u_data: xarray.DataArray, v_data: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.DataArray¶
Calculate Rossby wave sources (RWS).
\[RWS=-\nabla \cdot \left({v}_{x}\zeta \right)=-\left(\zeta \nabla \cdot {v}_{x}+{v}_{x}\cdot \nabla \zeta \right)\]with \(\zeta\) being the absolute vorticity.
Parameters¶
- u_data
xarray.DataArray. The zonal component of the wind.
- v_data
xarray.DataArray. The meridional component of vector wind.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
Rossby wave sources (
xarray.DataArray).Reference¶
Sardeshmukh, P. D., & Hoskins, B. J. (1988). The Generation of Global Rotational Flow by Steady Idealized Tropical Divergence. Journal of Atmospheric Sciences, 45(7), 1228-1251. https://doi.org/10.1175/1520-0469(1988)045<1228:TGOGRF>2.0.CO;2
James IN (1994) Low frequency variability of the circulation. Introduction to Circulating Atmospheres. Cambridge University Press, Cambridge, UK, pp 255–301
Trenberth, K. E., Branstator, G. W., Karoly, D., Kumar, A., Lau, N.-C., and Ropelewski, C. (1998), Progress during TOGA in understanding and modeling global teleconnections associated with tropical sea surface temperatures, J. Geophys. Res., 103(C7), 14291–14324, doi: https://doi.org/10.1029/97JC01444.
Nie, Y., Zhang, Y., Yang, X.-Q., & Ren, H.-L. (2019). Winter and summer Rossby wave sources in the CMIP5 models. Earth and Space Science, 6, 1831–1846. https://doi.org/10.1029/2019EA000674
Fuentes-Franco, R., Koenigk, T., Docquier, D. et al. Exploring the influence of the North Pacific Rossby wave sources on the variability of summer atmospheric circulation and precipitation over the Northern Hemisphere. Clim Dyn 59, 2025–2039 (2022). https://doi.org/10.1007/s00382-022-06194-4
- u_data
- easyclimate.core.windspharm.calc_gradient(data_input: xarray.DataArray, truncation: int = None, R: float = 6371200.0, legfunc: str = 'stored', lon_dim: str = 'lon', lat_dim: str = 'lat') xarray.Dataset¶
Computes the vector gradient of a scalar field on the sphere.
Parameters¶
- data_input
xarray.DataArray The spatio-temporal data to be calculated.
- truncation:
int. Truncation limit (triangular truncation) for the spherical harmonic computation.
- R:
float. The radius in metres of the sphere used in the spherical harmonic computations. Default is 6371200 m, the approximate mean spherical Earth radius.
- legfunc:
str, ‘stored’ (default) or ‘computed’. If ‘stored’, associated legendre functions are precomputed and stored when the class instance is created. This uses \(O(\mathrm{nlat}^3)\) memory, but speeds up the spectral transforms. If ‘computed’, associated legendre functions are computed on the fly when transforms are requested. This uses \(O(\mathrm{nlat}^2)\) memory, but slows down the spectral transforms a bit.
- lon_dim:
str, default: lon. Longitude coordinate dimension name. By default extracting is applied over the lon dimension.
- lat_dim:
str, default: lat. Latitude coordinate dimension name. By default extracting is applied over the lat dimension.
Returns¶
The zonal and meridional components of the vector gradient respectively (
xarray.Dataset).- data_input