Additional functions¶
|
Find the index of the item in the array nearest to the value |
|
Find latitude of absolute maximum value for a given interval |
|
Calculate unweighted seasonal means from monthly time series. |
|
Calculates the meridional mass streamfunction by integrating meridional wind from top of the atmosphere to surface |
Calculate the tropopause height in isobaric coordinates |
|
|
Find the first (with increasing index) zero crossing of the function F. |
Description of additional functions¶
- skytropd.functions.TropD_Calculate_MaxLat(F: ndarray, lat: ndarray, n: int = 6, axis: int = -1) ndarray[source]¶
Find latitude of absolute maximum value for a given interval
Note: assumes a smoothly varying function
- Parameters:
F (numpy.ndarray (dim1, ..., dimN-1, lat)) – N-dimensional array w/ lat as specified axis (default last), data assumed contingous with invalid data only on ends. interior nans are untested and will prompt warning
lat (numpy.ndarray) – latitude array
n (int, optional) –
- rank of moment used to calculate the position of max value.
n = 1,2,4,6,8,… , by default 6
axis (int, optional) – axis corresponding to latitude, by default -1
- Returns:
N-1 dimensional array of latitude(s) of max values of F
- Return type:
numpy.ndarray (dim1, …, dimN-1)
- skytropd.functions.TropD_Calculate_Mon2Season(Fm: ndarray, season: List[int], m: Optional[int] = None, first_jan_idx: int = 0, axis: int = 0) ndarray[source]¶
Calculate unweighted seasonal means from monthly time series. Dec is averaged with Jan and Feb of the same year, not the following year
- Parameters:
Fm (numpy.ndarray (time, dim2, ..., dimN)) – N-dimensional array, only utilizes full years-worth of data will be utilized
season (list of int) – 0-based list of month indices, e.g., [11,0,1] for DJF
m (int, optional) – deprecated. use first_jan_idx instead.
first_jan_idx (int, optional) – index of first occurence of Jan in Fm, by default 0
axis (int, optional) – time axis to resample along, by default 0
- Returns:
the annual time series of the seasonal means (not weighted for unequal month lengths)
- Return type:
numpy.ndarray (time, dim2, …, dimN)
- skytropd.functions.TropD_Calculate_StreamFunction(V: ndarray, lat: ndarray, lev: ndarray) ndarray[source]¶
Calculates the meridional mass streamfunction by integrating meridional wind from top of the atmosphere to surface
- Parameters:
V (numpy.ndarray (dim1, ..., dimN-2, lat, lev) or (dim1, ..., dimN-2, lev, lat)) – N-dimensional array of zonal-mean meridional wind with final 2 axes corresponding to latitude and vertical level, respectively, or the swapped trailing order
(..., lev, lat)lat (numpy.ndarray) – equally-spaced latitude array
lev (numpy.ndarray) – vertical level array in hPa
- Returns:
the meridional mass overturning streamfunction (psi)
- Return type:
numpy.ndarray (same shape as V)
- skytropd.functions.TropD_Calculate_TropopauseHeight(T: ndarray, P: ndarray, Z: None) ndarray[source]¶
- skytropd.functions.TropD_Calculate_TropopauseHeight(T: ndarray, P: ndarray, Z: ndarray) Tuple[ndarray, ndarray]
Calculate the tropopause height in isobaric coordinates
Based on the method in Birner (2010), according to the WMO definition: first level where the lapse rate <= 2K/km AND where the lapse rate <= 2K/km for at least 2km above that level
- Parameters:
T (numpy.ndarray (dim1, ..., dimN-1, lev)) – N-dimensional temperature array
P (numpy.ndarray (lev,)) – pressure levels in hPa
Z (numpy.ndarray (same shape as T), optional) – N-dimensional array to be interpolated at tropopause, such as geopotential height (m) to yield tropopause height
- Returns:
numpy.ndarray (dim1, …, dimN-1) – the tropopause pressure in hPa
numpy.ndarray (dim1, …, dimN-1), optional – returned if Z is provided. Corresponds to Z evaluated at the tropopause. If Z is geopotential height, it is tropopause altitude (m)
- skytropd.functions.TropD_Calculate_ZeroCrossing(F: ndarray, lat: ndarray, lat_uncertainty: float = 0.0, axis: int = -1) ndarray[source]¶
Find the first (with increasing index) zero crossing of the function F.
This function uses the compiled backend for the per-latitude-band search.
- Parameters:
F (numpy.ndarray (dim1, ..., dimN-1, lat)) – N dimensional array to search
lat (numpy.ndarray (lat,)) – latitude array
lat_uncertainty (float, optional) – same unit as lat. The minimum distance allowed between adjacent zero crossings of identical sign change, by default 0.0. For example, for lat_uncertainty = 10, if the most equatorward zero crossing is from positive to negative, NaN is returned if an additional zero crossing from positive to negative is found within 10 degrees of the first one.
axis (int, optional) – axis corresponding to latitude, by default -1 (last)
- Returns:
latitude(s) of zero crossing by linear interpolation
- Return type:
np.ndarray
- skytropd.functions.find_nearest(array: ndarray, value: float, axis: Optional[int] = None, skipna: bool = False) ndarray[source]¶
Find the index of the item in the array nearest to the value
- Parameters:
array (numpy.ndarray) – array to search
value (float) – value to be found
axis (int, optional) – the axis of the array to search. if not given, the array is flattened prior to being searched. If given, output will be (n-1)-dimensional, returning the position of the value within the specified axis
skipna (bool, optional) – whether to skip over NaN in the array, by default False
- Returns:
index(es) nearest to value in array
- Return type:
numpy.ndarray