Science calculations¶
Calculus¶
-
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.
spatial
(cube, aggr, coords=['latitude', 'longitude'])[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).
coords (list, optional) – List of names of spatial coordinates.
- Returns
Collapsed cube.
- Return type
Examples
>>> spatial(my_data_cube, "mean")
-
aeolus.calc.
meridional_mean
(cube, lat_name='latitude')[source]¶ Calculate cube’s meridional average.
- Parameters
cube (iris.cube.Cube) – Cube with a latitude coordinate.
lat_name (str, optional) – Name of the latitude coordinate.
- Returns
Collapsed cube.
- Return type
Diagnostics¶
-
aeolus.calc.
bond_albedo
(cubelist)[source]¶ Bold albedo.
\[1 - 4 \frac{OSR_{TOA}}{S_{0}}\]- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes. Cubes must have “planet_conf” attribute to get solar constant.
- Returns
Difference between the region averages.
- Return type
-
aeolus.calc.
ghe_norm
(cubelist)[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.
- Returns
Difference between the region averages.
- Return type
-
aeolus.calc.
heat_redist_eff
(cubelist, region_a, region_b)[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).
- Returns
Eta parameter.
- Return type
-
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.
- Return type
-
aeolus.calc.
sfc_net_energy
(cubelist)[source]¶ Calculate domain-average surface energy flux.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
- Returns
Cube with reduced dimensions.
- Return type
-
aeolus.calc.
sfc_water_balance
(cubelist)[source]¶ Calculate domain-average precipitation minus evaporation.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
- Returns
Cube with reduced 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.
- Return type
-
aeolus.calc.
toa_cloud_radiative_effect
(cubelist, kind)[source]¶ Calculate domain-average TOA cloud radiative effect (CRE).
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes
kind (str) – Shortwave (‘sw’), longwave (‘lw’), or ‘total’ CRE
- Returns
Cube of CRE with reduced dimensions.
- Return type
-
aeolus.calc.
toa_eff_temp
(cubelist)[source]¶ Calculate effective temperature from TOA OLR.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
- Returns
- Return type
-
aeolus.calc.
toa_net_energy
(cubelist)[source]¶ Calculate domain-average TOA energy flux.
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
- Returns
Cube with reduced dimensions.
- Return type
-
aeolus.calc.
total_precip
(cubelist, ptype=None)[source]¶ Calculate total precipitation flux [\(mm~day^{-1}\)].
- Parameters
cubelist (iris.cube.CubeList) – Input list of cubes.
ptype (str, optional) – Precipitation type (stra|conv).
- Returns
Sum of cubes of precipitation with units converted to mm per day.
- Return type
-
aeolus.calc.
water_path
(cubelist, kind='water_vapour', coord_name='level_height')[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.
coord_name (str or iris.coords.Coord, optional) – Vertical coordinate for integration.
- Returns
Difference between the region averages.
- Return type