Science calculations¶
Calculus¶
-
aeolus.calc.
d_dx
(cube, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate a derivative w.r.t. x-coordinate.
-
aeolus.calc.
d_dy
(cube, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate a derivative w.r.t. y-coordinate.
-
aeolus.calc.
d_dz
(cube, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate a derivative w.r.t. z-coordinate.
-
aeolus.calc.
deriv
(cube, coord)[source]¶ Calculate a derivative w.r.t. the given coordinate.
Uses iris.analysis.calculus.differentiate and then interpolates the result to the grid points of the original cube.
- Parameters
cube (iris.cube.Cube) – Input cube containing the given coordinate.
coord (str or iris.coords.Coord) – Coordinate for differentiation.
- Returns
d(cube)/d(coord).
- Return type
See also
aeolus.calc.calculus.d_dx()
,aeolus.calc.calculus.d_dy()
,aeolus.calc.calculus.d_dz()
-
aeolus.calc.
div_h
(i_cube, j_cube, r_planet=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate horizontal divergence.
Note: currently works only in spherical coordinates.
- Parameters
i_cube – i-component (e.g. u-wind)
j_cube – j-component (e.g. v-wind)
r_planet (float, optional) – Radius of the planet (m). If not given, an attempt is made to get it from the cube metadata.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of horizontal divergence.
- Return type
Notes
Divergence in spherical coordinates is defined as
\[\nabla\cdot \vec A = \frac{1}{r cos \theta} ( \frac{\partial \vec A_\lambda}{\partial \lambda} + \frac{\partial}{\partial \theta} (\vec A_\theta cos \theta))\]where lambda is longitude, theta is latitude.
-
aeolus.calc.
integrate
(cube, coord)[source]¶ Integrate the cube along a 1D coordinate using the trapezoidal rule.
Note: coord must be one of the dimensional coordinates of the cube.
- Parameters
cube (iris.cube.Cube) – Input cube containing the given coordinate.
coord (str or iris.coords.Coord) – Coordinate for integration.
- Returns
Integrated cube.
- Return type
Statistics¶
-
aeolus.calc.
cumsum
(cube, axis, axis_weights=False, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Cumulative sum of a cube.
- Parameters
cube (iris.cube.Cube) – Input cube.
axis (str) – Coordinate axis of operation (t|z|y|x).
axis_weights (bool, optional) – If True, multiply data by the coordinate spacings.
model (aeolus.model.Model, optional) – Model class with a relevant coordinate name.
- Returns
Cube of cumulative sums with the same dimensions as the input cube.
- Return type
-
aeolus.calc.
spatial
(cube, aggr, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate spatial statistic with geographic grid weights.
- Parameters
cube (iris.cube.Cube) – Cube with longitude and latitude coordinates.
aggr (str) – Statistical aggregator (see iris.analysis for available aggregators).
model (aeolus.model.Model, optional) – Model class with relevant coordinate names.
- Returns
Collapsed cube.
- Return type
Examples
>>> spatial(my_data_cube, "mean")
-
aeolus.calc.
spatial_quartiles
(cube, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate quartiles over horizontal coordinates.
-
aeolus.calc.
minmaxdiff
(cubelist, name)[source]¶ Spatial maximum minus spatial minimum for a given cube.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
name (str) – Cube name.
- Returns
Difference between the extrema with collapsed spatial dimensions.
- Return type
-
aeolus.calc.
region_mean_diff
(cubelist, name, region_a, region_b)[source]¶ Difference between averages over two regions for a given cube.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
name (str) – Cube name.
region_a (aeolus.region.Region) – First region.
region_b (aeolus.region.Region) – Second region.
- Returns
Difference between the region averages with collapsed spatial dimensions.
- Return type
-
aeolus.calc.
meridional_mean
(cube, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate cube’s meridional average.
- Parameters
cube (iris.cube.Cube) – Cube with a latitude coordinate.
model (aeolus.model.Model, optional) – Model class with a relevant coordinate name.
- Returns
Collapsed cube.
- Return type
-
aeolus.calc.
normalize_cube
(cube)[source]¶ Normalize cube data, i.e. make the values range from 0 to 1.
\[z_i = (x_i - min(x)) / (max(x) - min(x))\]- Parameters
cube (iris.cube.Cube) – The input cube.
-
aeolus.calc.
zonal_mean
(cube, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate cube’s zonal average.
- Parameters
cube (iris.cube.Cube) – Cube with a latitude coordinate.
model (aeolus.model.Model, optional) – Model class with a relevant coordinate name.
- Returns
Collapsed cube.
- Return type
-
aeolus.calc.
last_n_day_mean
(cube, days=365, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Average the cube over the last n days of its time dimension.
-
aeolus.calc.
vertical_mean
(cube, weight_by=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Vertical mean of a cube with optional weighting.
- Parameters
cube (iris.cube.Cube) – Cube to average.
weight_by (str or iris.coords.Coord or iris.cube.Cube, optional) – Coordinate of the given cube or another cube used for weighting.
model (aeolus.model.Model, optional) – Model class with a relevant coordinate name.
- Returns
Collapsed cube.
- Return type
Diagnostics¶
-
aeolus.calc.
air_density
(cubelist, const=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Get air density from the given cube list.
If not present, it is attempted to calculate it from pressure and temperature.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes containing temperature.
const (aeolus.const.const.ConstContainer, optional) – Must have dry_air_gas_constant as an attribute. If not given, an attempt is made to retrieve it from cube attributes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of air density.
- Return type
-
aeolus.calc.
air_temperature
(cubelist, const=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Get the real temperature from the given cube list.
If not present, it is attempted to calculate it from other variables, Exner pressure or air pressure, and potential temperature.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes containing temperature.
const (aeolus.const.const.ConstContainer, optional) – Must have reference_surface_pressure and dry_air_gas_constant as attributes. If not given, an attempt is made to retrieve it from cube attributes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of air temperature.
- Return type
-
aeolus.calc.
bond_albedo
(cubelist, const=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Bold albedo.
\[4 \frac{OSR_{TOA}}{S_{0}}\]- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
const (aeolus.const.const.ConstContainer, optional) – Must have a ScalarCube of condensible_density as an attribute. If not given, attempt to retrieve it from cube attributes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of bond albedo.
- Return type
-
aeolus.calc.
dry_lapse_rate
(cubelist, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Dry lapse rate, or the change of air temperature with altitude.
\[\gamma = \partial T_{air} / \partial z\]- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes containing temperature.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of dT/dz.
- Return type
-
aeolus.calc.
flux
(cubelist, quantity, axis, weight_by_density=True, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate horizontal or vertical flux of some quantity.
\[F_{x} = u (\rho q), F_{y} = v (\rho q), F_{z} = w (\rho q)\]- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
quantity (str or iris.Constraint) – Quantity (present in the cube list).
axis (str) – Axis of the flux component (x|y|z)
weight_by_density (bool, optional) – Multiply by a cube of air density (must be present in the input cube list).
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of a flux component with the same dimensions as input cubes.
- Return type
-
aeolus.calc.
geopotential_height
(cubelist, const=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Get geopotential height from the given cube list.
If not present, the altitude coordinate is transformed into a cube.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes containing temperature.
const (aeolus.const.const.ConstContainer, optional) – Must have gravity as an attribute. If not given, an attempt is made to retrieve it from cube attributes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of geopotential height.
- Return type
-
aeolus.calc.
ghe_norm
(cubelist, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Normalised greenhouse effect parameter.
\[GHE = 1 - \left(\frac{T_{eff}}{T_{sfc}}\right)^{1/4}\]- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of greenhouse effect parameter with collapsed spatial dimensions.
- Return type
-
aeolus.calc.
heat_redist_eff
(cubelist, region_a, region_b, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Heat redistribution efficiency (Leconte et al. 2013).
\[\eta=\frac{OLR_{TOA,night}}{OLR_{TOA,day}}\]- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
region_a (aeolus.region.Region) – First region (usually nightside).
region_b (aeolus.region.Region) – Second region (usually dayside).
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of eta parameter with collapsed spatial dimensions.
- Return type
-
aeolus.calc.
horiz_wind_cmpnts
(cubelist, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Extract u and v wind components from a cube list and interpolate v on u’s grid if necessary.
- Parameters
cubelist (iris.cube.CubeList) – List of cubes with horizontal wind components
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
u, v – Cubes of wind components.
- Return type
-
aeolus.calc.
precip_sum
(cubelist, ptype='total', const=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate a sum of different types of precipitation [\(mm~day^{-1}\)].
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
ptype (str, optional) – Precipitation type (total|stra|conv|rain|snow).
const (aeolus.const.const.ConstContainer, optional) – Must have a ScalarCube of condensible_density as an attribute. If not given, attempt to retrieve it from cube attributes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Sum of cubes of precipitation with units converted to mm per day.
- Return type
-
aeolus.calc.
sfc_net_energy
(cubelist, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate domain-average surface energy flux.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes with net LW and SW radiation, sensible and latent surface fluxes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of total surface downward energy flux.
- Return type
-
aeolus.calc.
sfc_water_balance
(cubelist, const=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate domain-average precipitation minus evaporation.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
const (aeolus.const.const.ConstContainer, optional) – Must have a ScalarCube of condensible_density as an attribute. If not given, attempt is made to retrieve it from cube attributes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of total surface downward water flux (P-E).
- Return type
-
aeolus.calc.
superrotation_index
(cubelist, const=None, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Local superrotation index.
\[ \begin{align}\begin{aligned}s = \frac{m}{\Omega a^2} - 1,\\m = a cos\phi(\Omega a cos\phi + u)\end{aligned}\end{align} \]- Parameters
cubelist (iris.cube.CubeList) – List of cubes containing a cube of zonal velocity (u).
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
s_idx – Cubes of superrotation index.
- Return type
References
Read (1986), https://doi.org/10.1002/qj.49711247114
-
aeolus.calc.
toa_cloud_radiative_effect
(cubelist, kind, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate domain-average TOA cloud radiative effect (CRE).
\[CRE_{TOA} = F_{up,clear-sky} - F_{up,all-sky}\]- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes
kind (str) – Shortwave (‘sw’), longwave (‘lw’), or ‘total’ CRE
- Returns
Cube of CRE.
- Return type
-
aeolus.calc.
toa_eff_temp
(cubelist, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate effective temperature from TOA OLR.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of \(T_{eff}\).
- Return type
-
aeolus.calc.
toa_net_energy
(cubelist, model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Calculate domain-average TOA energy flux.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
model (aeolus.model.Model, optional) – Model class with relevant variable names.
- Returns
Cube of total TOA downward energy flux.
- Return type
-
aeolus.calc.
water_path
(cubelist, kind='water_vapour', model=<class 'aeolus.model.base.Model'>(53 fields))[source]¶ Water vapour or condensate path, i.e. a vertical integral of a water phase.
\[WP = \int_{z_{sfc}}^{z_{top}} \rho q dz\]- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes containing appropriate mixing ratio and air density.
kind (str, optional) – Short name of the water phase to be integrated. Options are water_vapour (default) | liquid_water | ice_water | cloud_water cloud_water is the sum of liquid and ice phases.
model (aeolus.model.Model, optional) – Model class with relevant coordinate names. model.z is used as a vertical coordinate for integration.
- Returns
Cube of water path with collapsed vertical dimension.
- Return type