tropy.analysis_tools¶
The tropy.analysis package is a collection of tools for dedicated analysis of scientific data. It has been especially used to explore atmospheric data, but the application to other types of data might be rather straight forward and meaningful.
tropy.analysis_tools.derived_variables¶
This module contains / will contain variables derived from some basic input fields. The name “basic” refers to the typical set of variables stored in climate or weather model output.
-
tropy.analysis_tools.derived_variables.calc_lwp(p, T, qv, qc, g=9.8, axis=-1)¶ Calculates the liquid water path.
- Parameters
p (numpy array) – atmospheric pressure [hPa]
T (numpy array) – temperature [K]
qv (numpy array) – water vapor mixing ratio [kg / kg]
qc (numpy array) – liquid water mixing ratio [kg / kg]
- Returns
lwp – liquid water path [g / m**2]
- Return type
numpy array
Notes
hydrostatic approximation is used to convert vertical integration from pressure to height levels
impact of condesate mass on air density is ignored.
tropy.analysis_tools.grid_and_interpolation¶
-
tropy.analysis_tools.grid_and_interpolation.COSMO_mean_vertical_levels(Nlev)¶ Outputs the mean vertical levels of COSMO-DE taken from the model description of COSMO (see reference below).
- Parameters
Nlev (int) –
number of vertical levels, i.e.
Nlev = 51 for intermediate half-levels (e.g. for w) Nlev = 50 for main levels (e.g. for u)
- Returns
lev – mean vertical COSMO-DE levels
- Return type
numpy array (1-dim)
References
M. Baldauf, J. Foerstner, S. Klink, T. Reinhardt, C. Schraff, A. Seifert und K. Stephan “Kurze Beschreibung des Lokal-Modells Kuerzestfrist COSMO-DE (LMK) und seiner Datenbanken auf dem Datenserver des DWD”, Version 1.6, Stand: 31. 03. 2011
Tabelle 1, Seite 29
-
tropy.analysis_tools.grid_and_interpolation.create_interpolation_index(lon1, lat1, lon2, lat2, xy=False)¶ Given georeference of two grids in lon and lat, an index is built to map a field on grid 1 to grid 2.
- Parameters
lon1 (numpy array (2-dim)) – longitude of input grid 1
lat1 (numpy array (2-dim)) – latitude of input grid 1
lon2 (numpy array (2-dim)) – longitude of target grid 2
lat2 (numpy array (2-dim)) – latitude of target grid 2
xy (bool, optional, default = False) – switch if georef field are interpreted as local Cartesian coordinates: Then no internal transformation is used.
- Returns
ir (numpy array, 2-dim, dtype = int) – index field for transforming rows
ic (numpy array, 2-dim, dtype = int) – index field for transforming columns
-
tropy.analysis_tools.grid_and_interpolation.curve_flow_filter(f, numberOfIterations=5)¶ Smoothing filter depending on isoline curvature. Interface for curvature flow filter from simpleITK toolkit.
- Parameters
f (numpy array (2-dim)) – 2d field to be filtered (smoothed)
numberOfIterations (int, optional, default = 5) – number of iterations, increases smooting effect
- Returns
f_sm – smoothed 2d field
- Return type
numpy array (2-dim)
Notes
Only works if SimpleITK is installed !!!
-
tropy.analysis_tools.grid_and_interpolation.cutout_cluster(c, nc, nedge=30)¶ Makes a cutout of a categorial field c for an object of class / number nc.
-
tropy.analysis_tools.grid_and_interpolation.cutout_field4box(f, ind, bsize, **kwargs)¶ Cuts out a field based on center pix index and box size.
-
tropy.analysis_tools.grid_and_interpolation.cutout_fields(fin, slices, vaxis=0)¶ Cuts out a field or a list of similar fields given row and column slices.
- Parameters
- Returns
fout – resulting cutout (possibly a list of cutted fields)
- Return type
2d or 3d numpy array or list of fields
-
tropy.analysis_tools.grid_and_interpolation.get_index(p, lon, lat)¶ Given a point p the nearest grid point is returned.
- Parameters
p (numpy array) – point in the 2d grid
lon (numpy array (2-dim)) – grid values of the 1st dimension (can be longitude, x, etc.)
lat (numpy array (2-dim)) – grid values of the 2nd dimension (can be latitude, y, etc.)
- Returns
ir (numpy array (2-dim)) – row index (the 1st)
ic (numpy array (2-dim)) – column index (the 2nd)
-
tropy.analysis_tools.grid_and_interpolation.i2iset(v, i)¶ Returns an index set for a n-dim numpy array given an index i which points to a position within the same but flattened array.
It is assumed that the array was flattened in C-mode.
- Parameters
v (numpy array (n-dim)) – field which determines the rank and shape
i (int) – input index which points to the position in v.flatten()
- Returns
iset – index set which points to the position in v
- Return type
Notes
This routine helps you to locate a maximum in a n-dimensional array using e.g. the
np.argmaxfunction.
-
tropy.analysis_tools.grid_and_interpolation.interpolate_field_to_hor_pos(pos, lon, lat, field, **kwargs)¶ Given a sequence of points (or just one) a field (2d or 3d) given a grid is interpolated to the new positions.
- Parameters
pos (numpy array (2-dim)) – sequence of points given as np.array((x,y)) where x and y can be n-dimensional
lon (numpy array (2-dim)) – grid values of the 1st dimension (can be longitude, x, etc.)
lat (numpy array (2-dim)) – grid values of the 2nd dimension (can be latitude, y, etc.)
field (numpy array) – 2d or 3d field to be interpolated to pos
**kwargs (dict content) – keyword arguments for the griddata routine of the scipy.interpolate package can be passed through
- Returns
f – field interpolated to the position sequence pos, f has the dimension n x Nlev, where n is the dimension of pos, Nlev the vertical dimension
- Return type
numpy array
-
tropy.analysis_tools.grid_and_interpolation.ldiff(var, axis=-1)¶ Calculates difference of a field at adjacent levels. Now, wrapper for
numpy.diffwith default on last axis.- Parameters
var (numpy array (n-dim where n = [1,2,3])) – field; 1d, 2d or 3d array
axis (int, optional, default = -1) – axis used for layer avering
- Returns
dvar – level difference
- Return type
numpy array (n-dim where n = [1,2,3])
-
tropy.analysis_tools.grid_and_interpolation.ll2xy(lon, lat, lon0=10.0, lat0=50.0, R=6380.0)¶ Transformation between longitude and latitude and local Cartesian coordinates in west-east direction (x) and south-north direction (y). The local projection is called “sinusoidal”.
- Parameters
- Returns
x (numpy array) – Cartesian coordinate in west-east direction (increasing towards the East)
y (numpy array) – Cartesian coordinate in south-north direction (increasing towards the North)
Notes
The unit of Earth’s radius determines the units of x and y. Default (km).
References
See https://en.wikipedia.org/wiki/Sinusoidal_projection for further info.
-
tropy.analysis_tools.grid_and_interpolation.ll2xyc(lon, lat, mlon=None, mlat=None, lon0=0, lat0=0.0)¶ Applies a centered version of ll2xy. The centering is around the mean lon/lat value and lon0, lat0 only define the x/y offset.
- Parameters
lon (numpy array) – longitude
lat (numpy array) – latitude
mlon (float, optional, default = None) –
longitude center point of projection
if
Nonethen average oflonis usedmlat (float, optional, defaulft = None) –
latitutde center point of projection
if
Nonethen average oflatis usedlon0 (float, optional, default = 0.) – zero longitude (for x-offset)
lat0 (float, optional, default = 0.) – zero latitude (for y-offset)
- Returns
x (numpy array) – x-coordinate
y (numpy array) – y-coordinate
-
tropy.analysis_tools.grid_and_interpolation.lmean(var, axis=-1)¶ Calculates the layer mean of a field.
- Parameters
var (numpy array (n-dim where n = [1,2,3])) – field; 1d, 2d or 3d array
axis (int, optional, default = -1) – axis used for layer avering
- Returns
varm – layer mean
- Return type
numpy array (n-dim where n = [1,2,3])
-
tropy.analysis_tools.grid_and_interpolation.low2hres_region(region)¶ Based on low-res. region of MSG channels an correspoding HRV region is calculated.
- Parameters
region (tuple of two tubles) – region slice for cuuting a field ((row1, row2), (col1, col2))
- Returns
hrv_region
- Return type
corresponding HRV region as ((row1, row2), (col1, col2))
Notes
Calculation
It is assumed that the hrv region is a dependent variable which depends on the low res region attribute That means: Given the low res region cutout hrv region is determined as corresponding cutout which
exactly fits with the low res region and
refers to the artifical hrv full disk of 11136 x 11136
low res pixel with the index (I,J) = (0,0) has a high res pixel with index (i_m, j_m) = (2,2) in the middle and (i_u, j_u) = (1,1) in the upper corner see doc “MSG Level 1.5 Image Data Format Description”, Fig.8 which then leads to the relation between the low res pixel (I,J) and its corresponding upper corner high res pixel (i_u, j_u): (i_u, j_u) = 3 * (I, J ) + 1
-
tropy.analysis_tools.grid_and_interpolation.make_add_edge(var)¶ Adds egde region to 2d data field. Values in added left and right column, and lower and upper row are linearily extrapolated from their neighbors.
- Parameters
var (numpy array (2-dim)) – 2d variable field (Nrow, Ncols)
- Returns
var_new – 2d variable field with added edge values (Nrows + 2, Ncols + 2)
- Return type
numpy array (2-dim)
-
tropy.analysis_tools.grid_and_interpolation.make_hrv_upscaling(v)¶ A 2d field is upscaled to 3-fold higher resolution using linear interpolation. Edge values are determined by linear extrapolation.
- Parameters
v (numpy array (2-dim)) – 2d variable field (Nrow, Ncols)
- Returns
vhigh – variable field with 3-times higher resolution (3*Nrows, 3*Ncols)
- Return type
numpy array (2-dim)
-
tropy.analysis_tools.grid_and_interpolation.make_index_set(nrows, ncols)¶ This simple routine makes an 2d index set for given row and column number of a 2d field.
-
tropy.analysis_tools.grid_and_interpolation.make_vert_cut(p1, p2, lon, lat, vg, **kwargs)¶ Makes a horizontal - vertical cut through a 3d field.
- Parameters
p1 (list or numpy array) – 1st end of the cutting line (lon, lat)
p2 (list or numpy array) – 2nd end of the cutting line (lon, lat)
lon (numpy array (2-dim)) – longitude
lat (numpy array (2-dim)) – latitude
vg (numpy array (3-dim)) – 3d field to be interpolated
**kwargs (dict content) – keyword arguments for the griddata routine of the scipy.interpolate package can be passed through
- Returns
s (numpy array (1-dim)) – distance of km from one to the other end of the cutting line
lo (numpy array (1-dim)) – longitude along the cut
la (numpy array (1-dim)) – latitude along the cut
v (numpy array (2-dim)) – field interpolated to the cutting line
Vertical levels are kept.
-
tropy.analysis_tools.grid_and_interpolation.mid2edge_gridvector(x_mid, x_edge0)¶ Converts an mid-point based grid vector into an edge-based grid vector.
- Parameters
x_mid (numpy array (1-dim)) – mid-point-based grid vector (length N)
x_edge0 (float) – left edge
- Returns
x_edge – edge-based grid vector (length N + 1)
- Return type
numpy array (1-dim)
-
tropy.analysis_tools.grid_and_interpolation.region2slice(lon, lat, region)¶ Given a region of ((lon1,lon2),(lat1, lat2)) the routine provides a slicing ((row1, row2), (col1, col2)) than accomodates the region.
- Parameters
lon (numpy array (2-dim)) – longitude
lat (numpy array (2-dim)) – latitude
region (tuple of two tubles) – ((lon1, lon2), (lat1, lat2))
- Returns
slice – region slice for cuuting a field ((row1, row2), (col1, col2))
- Return type
tuple of two tubles
-
tropy.analysis_tools.grid_and_interpolation.remap_field(lon, lat, f, dr=0.5)¶ Do nearest nearbor remapping on a simple equi-distant local cartesian coordiante grid.
Uses curvature flow filter for smoothing and pixel edge correction.
- Parameters
lon (numpy array, (2-dim)) – longitude
lat (numpy array, (2-dim)) – latitude
f (numpy array, (2-dim)) – 2d field to be interpolated
dr (float, optional, default = 0.5) – grid spacing in km
- Returns
xnew (numpy array, (2-dim)) – regridded, equi-distant x-coordinate in km
ynew (numpy array, (2-dim)) – regridded, equi-distant y-coordinate in km
fint (numpy array, (2-dim)) – remapped field
-
tropy.analysis_tools.grid_and_interpolation.simple_pixel_area(lon, lat, xy=False, uncertainty=False)¶ Approximates the grid box area of a lon - lat grid.
The grid is treated as if it is an equidistant grid in local Cartesian coordinates on Earth’s surface (i.e. in x and y). The grid point is placed at each of the for corners of the grid box and the corresponding average area ist calculated.
- Parameters
lon (numpy array (2-dim)) – longitude
lat (numpy array (2-dim)) – latitude
xy (bool, optional, default = False) – switch if
lon&latare interpreted as local Cartesian coordinatesuncertainty (bool, optional, default = False) –
switch if std of pixel area is output
Deviation are caused by effects of the not-rectangular grid.
- Returns
a (numpy array (2-dim)) – grid box area
da (numpy array (2-dim), optional) – standard deviation of grid box area, if
uncertainty == True
-
tropy.analysis_tools.grid_and_interpolation.spline_smooting(x, y, s=0.1, fixed_endpoints=True)¶ Smoothes a curve with smooting spline.
- Parameters
x (numpy array (1-dim)) – abscissa values
y (numpy array (1-dim)) – ordinate values
s (float, optional, default = 0.1) – smooting parameter
fixed_endpoints (nool, optional, default = True) – switch, if endpoinds should be hold fixed
- Returns
y_smooth – smoothed version of y
- Return type
numpy array (1-dim)
Notes
Uses the function
scipy.interpolate.UnivariateSpline.
-
tropy.analysis_tools.grid_and_interpolation.tube_cutout4box(v3d, i0, boxsize)¶ Performs tube cutout with given index set. The center index might not be constant.
- Parameters
- Returns
vtube – 3d tube, quadratic cutout in the last two dimensions
- Return type
numpy array (3-dim)
Notes
This function can be used top make vertical cutouts, but a transpose command has to be applied before and after.
-
tropy.analysis_tools.grid_and_interpolation.xy2ll(x, y, lon0=10.0, lat0=50.0, R=6380)¶ Transformation between local Cartesian coordinates in west-east direction (x) and south-north direction (y) and longitude and latitude. This assumes a sinusoidal projection. Inversion of the function
ll2xy.- Parameters
x (numpy array) – Cartesian coordinate in west-east direction (increasing towards the East)
y (numpy array) – Cartesian coordinate in south-north direction (increasing towards the North)
lon0 (float, optional, default = 10.) – arbitrary initial longitude for x = 0
lat0 (float, optional, default = 50.) – arbitrary initial latitude for y = 0
R (float, optional, default = 6380.) – approximation for the radius of Earth
- Returns
lon (numpy array) – longitude
lat (numpy array) – latitude
Notes
The unit of Earth’s radius determines the units of x and y. Default (km).
References
See https://en.wikipedia.org/wiki/Sinusoidal_projection for further info.
tropy.analysis_tools.optical_flow¶
tropy.analysis_tools.parallax¶
tropy.analysis_tools.segmentation¶
tropy.analysis_tools.statistics¶
tropy.analysis_tools.thermodynamic_variables¶
Package for the calculation of several thermodynamical properties of moist air (including condensate).
ToDo
Some final tests of the routines should be done.
There was an issue with mixing ratio vs. specific humidity, but might be resolved …
-
tropy.analysis_tools.thermodynamic_variables.H2r(hrel, p, T)¶ Converts relative humidity H into mixing ratio r.
- Parameters
hrel (numpy array) – relative humidity in [%]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
es – water vapor pressure in [Pa]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.absolute_humidity(r, p, T)¶ Calculates absolute humidity.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
rho_v – absolute humidity / m**3]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.both_mass_fractions2mixing_ratios(qv, qw)¶ Converts water and the vapor mass fraction to condensed water and vapor mixing ratio.
- Parameters
qv (numpy array) – water vapor mass fraction [kg / kg]
qw (numpy array) – condensed water mass fraction [kg / kg]
- Returns
r (numpy array) – water vapor mixing ratio [kg / kg]
rw (numpy array) – condensed water mixing ratio [kg / kg]
-
tropy.analysis_tools.thermodynamic_variables.dew_point(r, p, T)¶ Calculates dew point temperature after Markowski book p 13 eq (2.25)
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
Td – dew point temperature in [K]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.dry_air_density(r, p, T)¶ Calculates dry air density.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
rho_d – dry air density in [kg / m**3]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.dry_air_potential_temperature(r, p, T)¶ Calculates dry air potential temperature.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
theta_d – dry air potential temperature (K)
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.equivalent_potential_temperature(r, p, T, eq=15)¶ Calculates equivalent potential temperature. after Bolton (1980) MWR 108, p.1046, eq. (43)
Davies-Jones (2009) MWR137, eq. (6.3)
also using either one of (15), (21) or (22)
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
eq (int, optional, default = 15) – which eq. of Bolton (1980) to be chosen for lifting condensation level temperature eq in [15, 21, 22]
- Returns
theta_E – equivalent potential temperature
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.lifting_condensation_level_temperature(r, p, T, eq=15)¶ Calculates lifting condensation level temperature. after Bolton (1980) MWR 108, p.1046, eq. (15), (21) or (22)
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
eq (int, optional, default = 15) – which eq. of Bolton (1980) to be chosen eq in [15, 21, 22]
- Returns
TL – lifting condensation level temperature in [K]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.moist_gas_constant(r)¶ Calculates gas constant of moist air.
-
tropy.analysis_tools.thermodynamic_variables.moist_potential_temperature(r, p, T)¶ Calculates moist air potential temperature.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
theta_m – moist air potential temperature (K)
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.moist_specific_heat(r)¶ Calculates specific heat at constant pressure of moist air.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
- Returns
c_pm – specific heat at constant pressure of moist air
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.relative_humidity(r, p, T)¶ Calculates relative humidity.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
H – relative humidity in [%]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.saturation_over_ice(T, method='PK')¶ Calculates saturation water vapor pressure over ice.
- Parameters
T (numpy array) – temperature [K]
method (str, optional, default = 'PK') – use an equation which approximates the vapor pressure over ice method in [‘PK’, ‘Murphy’]
- Returns
es – water vapor pressure in [Pa]
- Return type
numpy array
Notes
‘PK’ refers to Pruppacher and Klett (1997) ‘Murphy’ refers to Murphy and Koop (2005) eq. (7)
References
Murphy, D. M., and T. Koop (2005), Review of the vapour pressures of ice and supercooled water for atmospheric applications, Quart. J. Roy. Meteor. Soc., 131(608), 1539-1565, doi:10.1256/qj.04.94.
-
tropy.analysis_tools.thermodynamic_variables.saturation_pressure(T)¶ Calculates saturation water vapor pressure after Bolton (1980) MWR 108, p.1046, eq. (10)
- Parameters
T (numpy array) – temperature [K]
- Returns
es – water vapor pressure in [Pa]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.specific_humidity(r)¶ Calculates specific humidity.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
- Returns
q – specific humidity i.e. vapor mass per moist air mass [kg / kg]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.thermodynamic_constants()¶ Sets some thermodynamic constants.
- Parameters
None –
- Returns
const – set of thermodynamic constans saved in dictionary
- Return type
-
tropy.analysis_tools.thermodynamic_variables.total_density(r, r_w, p, T)¶ Calculates total density of moist air with hydrometeors.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
r_w (numpy array) – mixing ratio of condensed water i.e. condensed water mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
rho_d – dry air density in [kg / m**3]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.water_vapor_pressure(r, p, T)¶ Calculates water vapor pressure.
- Parameters
r (numpy array) – mixing ratio i.e. vapor mass per dry air mass in [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
e – water vapor pressure in [Pa]
- Return type
numpy array
-
tropy.analysis_tools.thermodynamic_variables.watermass_fraction2rw(qw, r)¶ Converts water mass fraction to condensed water mixing ratio.
- Parameters
qw (numpy array) – condensed water mass fraction [kg / kg]
r (numpy array) – mixing ration of water vapor [kg / kg]
p (numpy array) – total gas pressure [Pa]
T (numpy array) – temperature [K]
- Returns
rw – condensed water mixing ratio [kg / kg]
- Return type
numpy array