climada.hazard.trop_cyclone package#

climada.hazard.trop_cyclone.trop_cyclone module#

class climada.hazard.trop_cyclone.trop_cyclone.TropCyclone(category: ndarray | None = None, basin: List | None = None, windfields: List[csr_matrix] | None = None, **kwargs)[source]#

Bases: Hazard

Contains tropical cyclone events.

category#

for every event, the TC category using the Saffir-Simpson scale:

  • -1 tropical depression

  • 0 tropical storm

  • 1 Hurrican category 1

  • 2 Hurrican category 2

  • 3 Hurrican category 3

  • 4 Hurrican category 4

  • 5 Hurrican category 5

Type:

np.ndarray of ints

basin#

Basin where every event starts:

  • ‘NA’ North Atlantic

  • ‘EP’ Eastern North Pacific

  • ‘WP’ Western North Pacific

  • ‘NI’ North Indian

  • ‘SI’ South Indian

  • ‘SP’ Southern Pacific

  • ‘SA’ South Atlantic

Type:

list of str

windfields#

For each event, the full velocity vectors at each centroid and track position in a sparse matrix of shape (npositions, ncentroids * 2) that can be reshaped to a full ndarray of shape (npositions, ncentroids, 2).

Type:

list of csr_matrix

intensity_thres = 17.5#

intensity threshold for storage in m/s

vars_opt = {'category'}#

Name of the variables that are not needed to compute the impact.

__init__(category: ndarray | None = None, basin: List | None = None, windfields: List[csr_matrix] | None = None, **kwargs)[source]#

Initialize values.

Parameters:
  • category (np.ndarray of int, optional) – For every event, the TC category using the Saffir-Simpson scale:

    • ‘-1 tropical depression

    • ‘0 tropical storm

    • ‘1 Hurrican category 1

    • ‘2 Hurrican category 2

    • ‘3 Hurrican category 3

    • ‘4 Hurrican category 4

    • ‘5 Hurrican category 5

  • basin (list of str, optional) – Basin where every event starts:

    • ‘NA’ North Atlantic

    • ‘EP’ Eastern North Pacific

    • ‘WP’ Western North Pacific

    • ‘NI’ North Indian

    • ‘SI’ South Indian

    • ‘SP’ Southern Pacific

    • ‘SA’ South Atlantic

  • windfields (list of csr_matrix, optional) – For each event, the full velocity vectors at each centroid and track position in a sparse matrix of shape (npositions, ncentroids * 2) that can be reshaped to a full ndarray of shape (npositions, ncentroids, 2).

  • **kwargs (Hazard properties, optional) – All other keyword arguments are passed to the Hazard constructor.

set_from_tracks(*args, **kwargs)[source]#

This function is deprecated, use TropCyclone.from_tracks instead.

classmethod from_tracks(tracks: TCTracks, centroids: Centroids, pool: ProcessPool | None = None, model: str = 'H08', model_kwargs: dict | None = None, ignore_distance_to_coast: bool = False, store_windfields: bool = False, metric: str = 'equirect', intensity_thres: float = 17.5, max_latitude: float = 61, max_dist_inland_km: float = 1000, max_dist_eye_km: float = 300, max_memory_gb: float = 8)[source]#

Create new TropCyclone instance that contains windfields from the specified tracks.

This function sets the intensity attribute to contain, for each centroid, the maximum wind speed (1-minute sustained winds at 10 meters above ground) experienced over the whole period of each TC event in m/s. The wind speed is set to 0 if it doesn’t exceed the threshold intensity_thres.

The category attribute is set to the value of the category-attribute of each of the given track data sets.

The basin attribute is set to the genesis basin for each event, which is the first value of the basin-variable in each of the given track data sets.

Optionally, the time dependent, vectorial winds can be stored using the store_windfields function parameter (see below).

Parameters:
  • tracks (climada.hazard.TCTracks) – Tracks of storm events.

  • centroids (Centroids, optional) – Centroids where to model TC. Default: global centroids at 360 arc-seconds resolution.

  • pool (pathos.pool, optional) – Pool that will be used for parallel computation of wind fields. Default: None

  • model (str, optional) – Parametric wind field model to use. Default: “H08”.

    • "H1980" (the prominent Holland 1980 model) from the paper:

      Holland, G.J. (1980): An Analytic Model of the Wind and Pressure Profiles in Hurricanes. Monthly Weather Review 108(8): 1212–1218. https://doi.org/10.1175/1520-0493(1980)108<1212:AAMOTW>2.0.CO;2

    • "H08" (Holland 1980 with b-value from Holland 2008) from the paper:

      Holland, G. (2008). A revised hurricane pressure-wind model. Monthly Weather Review, 136(9), 3432–3445. https://doi.org/10.1175/2008MWR2395.1

    • "H10" (Holland et al. 2010) from the paper:

      Holland et al. (2010): A Revised Model for Radial Profiles of Hurricane Winds. Monthly Weather Review 138(12): 4393–4401. https://doi.org/10.1175/2010MWR3317.1

    • "ER11" (Emanuel and Rotunno 2011) from the paper:

      Emanuel, K., Rotunno, R. (2011): Self-Stratification of Tropical Cyclone Outflow. Part I: Implications for Storm Structure. Journal of the Atmospheric Sciences 68(10): 2236–2249. https://dx.doi.org/10.1175/JAS-D-10-05024.1

  • model_kwargs (dict, optional) – If given, forward these kwargs to the selected wind model. None of the parameters is currently supported by the ER11 model. Default: None. The Holland models support the following parameters, in alphabetical order:

    gradient_to_surface_windsfloat, optional

    The gradient-to-surface wind reduction factor to use. In H1980, the wind profile is computed on the gradient level, and wind speeds are converted to the surface level using this factor. In H08 and H10, the wind profile is computed on the surface level, but the clipping interval of the B-value depends on this factor. Default: 0.9

    rho_air_constfloat or None, optional

    The constant value for air density (in kg/m³) to assume in the formulas from Holland 1980. By default, the constant value suggested in Holland 1980 is used. If set to None, the air density is computed from pressure following equation (9) in Holland et al. 2010. Default: 1.15

    vmax_from_cenboolean, optional

    Only used in H10. If True, replace the recorded value of vmax along the track by an estimate from pressure, following equation (8) in Holland et al. 2010. Default: True

    vmax_in_bracketsbool, optional

    Only used in H10. Specifies which of the two formulas in equation (6) of Holland et al. 2010 to use. If False, the formula with vmax outside of the brackets is used. Note that, a side-effect of the formula with vmax inside of the brackets is that the wind speed maximum is attained a bit farther away from the center than according to the recorded radius of maximum winds (RMW). Default: False

    cyclostrophicbool, optional

    If True, do not apply the influence of the Coriolis force (set the Coriolis terms to 0). Default: True for H10 model, False otherwise.

  • ignore_distance_to_coast (boolean, optional) – If True, centroids far from coast are not ignored. If False, the centroids’ distances to the coast are calculated with the Centroids.get_dist_coast() method (unless there is “dist_coast” column in the centroids’ GeoDataFrame) and centroids far from coast are ignored. Default: False.

  • store_windfields (boolean, optional) – If True, the Hazard object gets a list windfields of sparse matrices. For each track, the full velocity vectors at each centroid and track position are stored in a sparse matrix of shape (npositions, ncentroids * 2) that can be reshaped to a full ndarray of shape (npositions, ncentroids, 2). Default: False.

  • metric (str, optional) – Specify an approximation method to use for earth distances:

    • “equirect”: Distance according to sinusoidal projection. Fast, but inaccurate for large distances and high latitudes.

    • “geosphere”: Exact spherical distance. Much more accurate at all distances, but slow.

    Default: “equirect”.

  • intensity_thres (float, optional) – Wind speeds (in m/s) below this threshold are stored as 0. Default: 17.5

  • max_latitude (float, optional) – No wind speed calculation is done for centroids with latitude larger than this parameter. Default: 61

  • max_dist_inland_km (float, optional) – No wind speed calculation is done for centroids with a distance (in km) to the coast larger than this parameter. Default: 1000

  • max_dist_eye_km (float, optional) – No wind speed calculation is done for centroids with a distance (in km) to the TC center (“eye”) larger than this parameter. Default: 300

  • max_memory_gb (float, optional) – To avoid memory issues, the computation is done for chunks of the track sequentially. The chunk size is determined depending on the available memory (in GB). Note that this limit applies to each thread separately if a pool is used. Default: 8

Raises:

ValueError

Return type:

TropCyclone

apply_climate_scenario_knu(percentile: str = '50', scenario: str = '4.5', target_year: int = 2050, **kwargs)[source]#

From current TC hazard instance, return new hazard set with future events for a given RCP scenario and year based on the parametrized values derived by Jewson 2021 (https://doi.org/10.1175/JAMC-D-21-0102.1) based on those published by Knutson 2020 (https://doi.org/10.1175/BAMS-D-18-0194.1). The scaling for different years and RCP scenarios is obtained by linear interpolation.

Note: Only frequency changes are applied as suggested by Jewson 2022 (https://doi.org/10.1007/s00477-021-02142-6). Applying only frequency anyway changes mean intensities and most importantly avoids possible inconsistencies (including possible double-counting) that may arise from the application of both frequency and intensity changes, as the relationship between these two is non trivial to resolve.

Parameters:
  • percentile (str) – percentiles of Knutson et al. 2020 estimates, representing the mode uncertainty in future changes in TC activity. These estimates come from a review of state-of-the-art literature and models. For the ‘cat05’ variable (i.e. frequency of all tropical cyclones) the 5th, 25th, 50th, 75th and 95th percentiles are provided. For ‘cat45’ and ‘intensity’, the provided percentiles are the 10th, 25th, 50th, 75th and 90th. Please refer to the mentioned publications for more details. possible percentiles:

    • ‘5/10’ either the 5th or 10th percentile depending on variable (see text above)

    • ‘25’ for the 25th percentile

    • ‘50’ for the 50th percentile

    • ‘75’ for the 75th percentile

    • ‘90/95’ either the 90th or 95th percentile depending on variable (see text above)

    Default: ‘50’

  • scenario (str) – possible scenarios:

    • ‘2.6’ for RCP 2.6

    • ‘4.5’ for RCP 4.5

    • ‘6.0’ for RCP 6.0

    • ‘8.5’ for RCP 8.5

  • target_year (int) – future year to be simulated, between 2000 and 2100. Default: 2050.

Returns:

haz_cc – Tropical cyclone with frequencies and intensity scaled according to the Knutson criterion for the given year, RCP and percentile. Returns a new instance of climada.hazard.TropCyclone, self is not modified.

Return type:

climada.hazard.TropCyclone

set_climate_scenario_knu(*args, **kwargs)[source]#

This function is deprecated, use TropCyclone.apply_climate_scenario_knu instead.

classmethod video_intensity(track_name: str, tracks: ~climada.hazard.tc_tracks.TCTracks, centroids: ~climada.hazard.centroids.centr.Centroids, file_name: str | None = None, writer: <module 'matplotlib.animation' from '/home/docs/checkouts/readthedocs.org/user_builds/climada-python/conda/stable/lib/python3.11/site-packages/matplotlib/animation.py'> = <matplotlib.animation.PillowWriter object>, figsize: ~typing.Tuple[float, float] = (9, 13), adapt_fontsize: bool = True, **kwargs)[source]#

Generate video of TC wind fields node by node and returns its corresponding TropCyclone instances and track pieces.

Parameters:
  • track_name (str) – name of the track contained in tracks to record

  • tracks (climada.hazard.TCTracks) – tropical cyclone tracks

  • centroids (climada.hazard.Centroids) – centroids where wind fields are mapped

  • file_name (str, optional) – file name to save video (including full path and file extension)

  • writer (matplotlib.animation., optional*) – video writer. Default is pillow with bitrate=500

  • figsize (tuple, optional) – figure size for plt.subplots

  • adapt_fontsize (bool, optional) – If set to true, the size of the fonts will be adapted to the size of the figure. Otherwise the default matplotlib font size is used. Default is True.

  • kwargs (optional) – arguments for pcolormesh matplotlib function used in event plots

Returns:

tc_list, tc_coord

Return type:

list(TropCyclone), list(np.ndarray)

Raises:

ValueError

frequency_from_tracks(tracks: List)[source]#

Set hazard frequency from tracks data.

Parameters:

tracks (list of xarray.Dataset)

classmethod from_single_track(track: Dataset, centroids: Centroids, idx_centr_filter: ndarray, model: str = 'H08', model_kwargs: dict | None = None, store_windfields: bool = False, metric: str = 'equirect', intensity_thres: float = 17.5, max_dist_eye_km: float = 300, max_memory_gb: float = 8)[source]#

Generate windfield hazard from a single track dataset

Parameters:
  • track (xr.Dataset) – Single tropical cyclone track.

  • centroids (Centroids) – Centroids instance.

  • idx_centr_filter (np.ndarray) – Indices of centroids to restrict to (e.g. sufficiently close to coast).

  • model (str, optional) – Parametric wind field model, one of “H1980” (the prominent Holland 1980 model), “H08” (Holland 1980 with b-value from Holland 2008), “H10” (Holland et al. 2010), or “ER11” (Emanuel and Rotunno 2011). Default: “H08”.

  • model_kwargs (dict, optional) – If given, forward these kwargs to the selected model. Default: None

  • store_windfields (boolean, optional) – If True, store windfields. Default: False.

  • metric (str, optional) – Specify an approximation method to use for earth distances: “equirect” (faster) or “geosphere” (more accurate). See dist_approx function in climada.util.coordinates. Default: “equirect”.

  • intensity_thres (float, optional) – Wind speeds (in m/s) below this threshold are stored as 0. Default: 17.5

  • max_dist_eye_km (float, optional) – No wind speed calculation is done for centroids with a distance (in km) to the TC center (“eye”) larger than this parameter. Default: 300

  • max_memory_gb (float, optional) – To avoid memory issues, the computation is done for chunks of the track sequentially. The chunk size is determined depending on the available memory (in GB). Default: 8

Raises:

ValueError, KeyError

Returns:

haz

Return type:

TropCyclone

climada.hazard.trop_cyclone.trop_cyclone_windfields module#

climada.hazard.trop_cyclone.trop_cyclone_windfields.MBAR_TO_PA = 100.0#

Unit conversion factors for JIT functions that can’t use ureg

climada.hazard.trop_cyclone.trop_cyclone_windfields.DEF_MAX_DIST_EYE_KM = 300#

Default value for the maximum distance (in km) of a centroid to the TC center at which wind speed calculations are done.

climada.hazard.trop_cyclone.trop_cyclone_windfields.DEF_INTENSITY_THRES = 17.5#

Default value for the threshold below which wind speeds (in m/s) are stored as 0.

climada.hazard.trop_cyclone.trop_cyclone_windfields.DEF_MAX_MEMORY_GB = 8#

Default value of the memory limit (in GB) for windfield computations (in each thread).

climada.hazard.trop_cyclone.trop_cyclone_windfields.MODEL_VANG = {'ER11': 3, 'H08': 0, 'H10': 2, 'H1980': 1}#

Enumerate different symmetric wind field models.

climada.hazard.trop_cyclone.trop_cyclone_windfields.DEF_RHO_AIR = 1.15#

Default value for air density (in kg/m³), following Holland 1980.

climada.hazard.trop_cyclone.trop_cyclone_windfields.DEF_GRADIENT_TO_SURFACE_WINDS = 0.9#

Default gradient-to-surface wind reduction factor, following the 90%-rule mentioned in:

Franklin, J.L., Black, M.L., Valde, K. (2003): GPS Dropwindsonde Wind Profiles in Hurricanes and Their Operational Implications. Weather and Forecasting 18(1): 32–44. https://doi.org/10.1175/1520-0434(2003)018<0032:GDWPIH>2.0.CO;2

According to Table 2, this is a reasonable factor for the 750 hPa level in the eyewall region. For other regions and levels, values of 0.8 or even 0.75 might be justified.

climada.hazard.trop_cyclone.trop_cyclone_windfields.T_ICE_K = 273.16#

Freezing temperatur of water (in K), for conversion between K and °C

climada.hazard.trop_cyclone.trop_cyclone_windfields.compute_angular_windspeeds(si_track: Dataset, d_centr: ndarray, mask_centr_close: ndarray, model: int, cyclostrophic: bool | None = False, model_kwargs: dict | None = None)[source]#

Compute (absolute) angular wind speeds according to a parametric wind profile

Parameters:
  • si_track (xr.Dataset) – Output of tctrack_to_si. Which data variables are used in the computation of the wind profile depends on the selected model.

  • d_centr (np.ndarray of shape (npositions, ncentroids)) – Distance (in m) between centroids and track positions.

  • mask_centr_close (np.ndarray of shape (npositions, ncentroids)) – For each track position one row indicating which centroids are within reach.

  • model (int) – Wind profile model selection according to MODEL_VANG.

  • model_kwargs (dict, optional) – If given, forward these kwargs to the selected model. Default: None

  • cyclostrophic (bool, optional, deprecated) – This argument is deprecated and will be removed in a future release. Include cyclostrophic as model_kwargs instead.

Returns:

containing the magnitude of the angular windspeed per track position per centroid location

Return type:

ndarray of shape (npositions, ncentroids)

climada.hazard.trop_cyclone.trop_cyclone_windfields.V_ANG_EARTH = 7.29e-05#

Earth angular velocity (in radians per second)

climada.hazard.trop_cyclone.trop_cyclone_windfields.compute_windfields_sparse(track: Dataset, centroids: Centroids, idx_centr_filter: ndarray, model: str = 'H08', model_kwargs: dict | None = None, store_windfields: bool = False, metric: str = 'equirect', intensity_thres: float = 17.5, max_dist_eye_km: float = 300, max_memory_gb: float = 8) Tuple[csr_matrix, csr_matrix | None][source]#

Version of compute_windfields that returns sparse matrices and limits memory usage

Parameters:
  • track (xr.Dataset) – Single tropical cyclone track.

  • centroids (Centroids) – Centroids instance.

  • idx_centr_filter (np.ndarray) – Indices of centroids to restrict to (e.g. sufficiently close to coast).

  • model (str, optional) – Parametric wind field model, one of “H1980” (the prominent Holland 1980 model), “H08” (Holland 1980 with b-value from Holland 2008), “H10” (Holland et al. 2010), or “ER11” (Emanuel and Rotunno 2011). Default: “H08”.

  • model_kwargs (dict, optional) – If given, forward these kwargs to the selected model. Default: None

  • store_windfields (boolean, optional) – If True, store windfields. Default: False.

  • metric (str, optional) – Specify an approximation method to use for earth distances: “equirect” (faster) or “geosphere” (more accurate). See dist_approx function in climada.util.coordinates. Default: “equirect”.

  • intensity_thres (float, optional) – Wind speeds (in m/s) below this threshold are stored as 0. Default: 17.5

  • max_dist_eye_km (float, optional) – No wind speed calculation is done for centroids with a distance (in km) to the TC center (“eye”) larger than this parameter. Default: 300

  • max_memory_gb (float, optional) – To avoid memory issues, the computation is done for chunks of the track sequentially. The chunk size is determined depending on the available memory (in GB). Default: 8

Raises:

ValueError

Returns:

  • intensity (csr_matrix) – Maximum wind speed in each centroid over the whole storm life time.

  • windfields (csr_matrix or None) – If store_windfields is True, the full velocity vectors at each centroid and track position are stored in a sparse matrix of shape (npositions, ncentroids * 2) that can be reshaped to a full ndarray of shape (npositions, ncentroids, 2). If store_windfields is False, None is returned.

climada.hazard.trop_cyclone.trop_cyclone_windfields.tctrack_to_si(track: Dataset, metric: str = 'equirect') Dataset[source]#

Convert track variables to SI units and prepare for wind field computation

In addition to unit conversion, the variable names are shortened, the longitudinal coordinates are normalized and additional variables are defined:

  • cp (coriolis parameter)

  • pdelta (difference between environmental and central pressure, always strictly positive)

  • vtrans (translational velocity vectors)

  • vtrans_norm (absolute value of translational speed)

Furthermore, some scalar variables are stored as attributes:

  • latsign (1.0 if the track is located on the northern and -1.0 if on southern hemisphere)

  • mid_lon (the central longitude that was used to normalize the longitudinal coordinates)

Finally, some corrections are applied to variables:

  • clip central pressure values so that environmental pressure values are never exceeded

  • extrapolate radius of max wind from pressure if missing

Parameters:
  • track (xr.Dataset) – Track information.

  • metric (str, optional) – Specify an approximation method to use for earth distances: “equirect” (faster) or “geosphere” (more accurate). See dist_approx function in climada.util.coordinates. Default: “equirect”.

Return type:

xr.Dataset

climada.hazard.trop_cyclone.trop_cyclone_windfields.get_close_centroids(si_track: Dataset, centroids: ndarray, buffer_km: float, metric: str = 'equirect') ndarray[source]#

Check whether centroids lay within a buffer around track positions

Note that, hypothetically, a problem occurs when a TC track is travelling so far in longitude that adding a buffer exceeds 360 degrees (i.e. crosses the antimeridian), which is physically impossible, but might happen with synthetical or test data.

Parameters:
  • si_track (xr.Dataset with dimension “time”) – Track information as returned by tctrack_to_si. Hence, longitudinal coordinates are normalized around the central longitude stored in the “mid_lon” attribute. This makes sure that the buffered bounding box around the track does not cross the antimeridian. The data variables used by this function are “lat”, and “lon”.

  • centroids (np.ndarray of shape (ncentroids, 2)) – Coordinates of centroids, each row is a pair [lat, lon]. The longitudinal coordinates are normalized within this function to be consistent with the track coordinates.

  • buffer_km (float) – Size of the buffer (in km). The buffer is converted to a lat/lon buffer, rescaled in longitudinal direction according to the t_lat coordinates.

  • metric (str, optional) – Specify an approximation method to use for earth distances: “equirect” (faster) or “geosphere” (more accurate). See dist_approx function in climada.util.coordinates. Default: “equirect”.

Returns:

  • centroids_close_normalized (np.ndarray of shape (nclose, 2)) – Coordinates of close centroids, each row is a pair [lat, lon]. The normalization of longitudinal coordinates is consistent with the track coordinates.

  • mask_centr (np.ndarray of shape (ncentroids,)) – Mask that is True for close centroids and False for other centroids.

  • mask_centr_alongtrack (np.ndarray of shape (npositions, nclose)) – Each row is a mask that indicates the centroids within reach for one track position. Note that these masks refer only to the “close centroids” to reduce memory requirements. The number of positions npositions corresponds to the size of the “time” dimension of si_track.