diff --git a/CHANGELOG.md b/CHANGELOG.md index f2b121151..b88436381 100644 --- a/CHANGELOG.md +++ b/CHANGELOG.md @@ -70,7 +70,11 @@ Removed: ### Added - `climada.hazard.tc_tracks.TCTracks.subset_years` function [#1023](https://github.com/CLIMADA-project/climada_python/pull/1023) -- `climada.hazard.tc_tracks.TCTracks.from_FAST` function, add Australia basin (AU) [#993](https://github.com/CLIMADA-project/climada_python/pull/993) +-`climada.hazard.tc_tracks.compute_track_density` function, `climada.hazard.tc_tracks.compute_genesis_density` function, + `climada.util.coordinates.compute_grid_cell_area` function, `climada.util.coordinates.compute_grid_cell_area_validation` function, + `climada.hazard.tc_tracks.normalize_hist` function, `climada.hazard.plot.plot_track_density` function + [#1003](https://github.com/CLIMADA-project/climada_python/pull/1003) +-`climada.hazard.tc_tracks.TCTracks.from_FAST` function, add Australia basin (AU) [#993](https://github.com/CLIMADA-project/climada_python/pull/993) - Add `osm-flex` package to CLIMADA core [#981](https://github.com/CLIMADA-project/climada_python/pull/981) - `doc.tutorial.climada_entity_Exposures_osm.ipynb` tutorial explaining how to use `osm-flex`with CLIMADA - `climada.util.coordinates.bounding_box_global` function [#980](https://github.com/CLIMADA-project/climada_python/pull/980) diff --git a/climada/hazard/plot.py b/climada/hazard/plot.py index 3ab1cec8b..1591b74cf 100644 --- a/climada/hazard/plot.py +++ b/climada/hazard/plot.py @@ -21,6 +21,9 @@ import logging +import cartopy.crs as ccrs +import cartopy.feature as cfeature +import matplotlib.colors as mcolors import matplotlib.pyplot as plt import numpy as np from deprecation import deprecated @@ -158,6 +161,136 @@ def plot_intensity( raise ValueError("Provide one event id or one centroid id.") + @staticmethod + def plot_track_density( + hist: np.ndarray, + axis=None, + projection=ccrs.Mollweide(), + add_features: dict = None, + title: str = None, + figsize=(12, 6), + div_cmap=False, + cbar=True, + cbar_kwargs: dict = { + "orientation": "horizontal", + "pad": 0.05, + "shrink": 0.8, + "label": "n° tracks per 1° x 1° grid cell", + }, + **kwargs, + ): + """ + Plot the track density of tropical cyclone tracks on a customizable world map. + + Parameters: + ---------- + hist: np.ndarray + 2D histogram of track density. + axis: GeoAxes, optional + Existing Cartopy axis. + projection: cartopy.crs, optional + Projection for the map. + add_features: dict + Dictionary of map features to add. Keys can be 'land', 'coastline', 'borders', and + 'lakes'. Values are Booleans indicating whether to include each feature. + title: str + Title of the plot. + figsize: tuple + Figure size when creating a new figure. + div_cmap: bool, default = False + If True, the colormap will be centered to 0. + cbar: bool, Default = True + If True, the color bar is added + cbar_kwargs: dict + dictionary containing keyword arguments passed to cbar + kwargs: + Additional keyword arguments passed to `ax.contourf`. + + Returns: + ------- + axis: GeoAxes + The plot axis. + + + Example: + -------- + >>> axis = plot_track_density( + ... hist=hist, + ... cmap='Spectral_r', + ... cbar_kwargs={'shrink': 0.8, + 'label': 'Cyclone Density [n° tracks / km²]', + 'pad': 0.1}, + ... add_features={ + ... 'land': True, + ... 'coastline': True, + ... 'borders': False, + ... 'lakes': False + ... }, + ... title='My Tropical Cyclone Track Density Map', + ... figsize=(10, 5), + ... levels=20 + ... ) + + """ + + # Default features + default_features = { + "land": True, + "coastline": True, + "borders": False, + "lakes": False, + } + add_features = add_features or default_features + + # Sample data + lon = np.linspace(-180, 180, hist.shape[1]) + lat = np.linspace(-90, 90, hist.shape[0]) + + # Create figure and axis if not provided + if axis is None: + _, axis = plt.subplots( + figsize=figsize, subplot_kw={"projection": projection} + ) + + # Add requested features + if add_features.get("land", False): + land = cfeature.NaturalEarthFeature( + category="physical", + name="land", + scale="50m", + facecolor="lightgrey", + alpha=0.6, + ) + axis.add_feature(land) + if add_features.get("coastline", False): + axis.add_feature(cfeature.COASTLINE, linewidth=0.5) + if add_features.get("borders", False): + axis.add_feature(cfeature.BORDERS, linestyle=":") + if add_features.get("lakes", False): + axis.add_feature(cfeature.LAKES, alpha=0.4, edgecolor="black") + + if div_cmap: + norm = mcolors.TwoSlopeNorm( + vmin=np.nanmin(hist), vcenter=0, vmax=np.nanmax(hist) + ) + kwargs["norm"] = norm + + # contourf = axis.contourf(lon, lat, hist, transform=ccrs.PlateCarree(), **kwargs) + contourf = axis.imshow( + hist, + extent=[lon.min(), lon.max(), lat.min(), lat.max()], + transform=ccrs.PlateCarree(), + origin="lower", + **kwargs, + ) + + if cbar: + plt.colorbar(contourf, ax=axis, **cbar_kwargs) + if title: + axis.set_title(title, fontsize=16) + + return axis + def plot_fraction(self, event=None, centr=None, smooth=True, axis=None, **kwargs): """Plot fraction values for a selected event or centroid. diff --git a/climada/hazard/tc_tracks.py b/climada/hazard/tc_tracks.py index 2a7441e70..e02332705 100644 --- a/climada/hazard/tc_tracks.py +++ b/climada/hazard/tc_tracks.py @@ -39,7 +39,6 @@ import matplotlib.cm as cm_mp import matplotlib.pyplot as plt import netCDF4 as nc -import numba import numpy as np import pandas as pd import pathos @@ -52,6 +51,7 @@ from matplotlib.lines import Line2D from shapely.geometry import LineString, MultiLineString, Point from sklearn.metrics import DistanceMetric +from tqdm import tqdm import climada.hazard.tc_tracks_synth import climada.util.coordinates as u_coord @@ -3094,3 +3094,196 @@ def _zlib_from_dataarray(data_var: xr.DataArray) -> bool: if np.issubdtype(data_var.dtype, float) or np.issubdtype(data_var.dtype, int): return True return False + + +def compute_track_density( + tc_track: TCTracks, + res: int = 5, + bounds: tuple = None, + genesis: bool = False, + norm: str = None, + filter_tracks: bool = True, + wind_min: float = None, + wind_max: float = None, +) -> tuple[np.ndarray, np.ndarray, np.ndarray]: + """Compute absolute and normalized tropical cyclone track density. Before using this function, + apply the same temporal resolution to all tracks by calling :py:meth:`equal_timestep` on the + TCTrack object. Due to the computational cost of the this function, it is not recommended to + use a grid resolution higher tha 0.1°. First, this function creates 2D bins of the specified + resolution (e.g. 1° x 1°). Second, since tracks are not lines but a series of points, it counts + the number of points per bin. Lastly, it returns the absolute or normalized count per bin. + To plot the output of this function use :py:meth:`plot_track_density`. + + Parameters: + ---------- + tc_track: TCTracks object + track object containing a list of all tracks + res: int (optional), default: 5° + resolution in degrees of the grid bins in which the density will be computed + bounds: tuple, (optional) dafault: None + (lon_min, lat_min, lon_max, lat_max) latitude and longitude bounds. + genesis: bool, (optional) default = False + If true the function computes the track density of only the genesis location of tracks + norm: str (optional), default: None + If None the function returns the number of samples in each bin. If True, it normalize the + bin count as specified: if norm = area -> normalize by gird cell area. If norm = sum -> + normalize by the total sum across all bins. + filter_tracks: bool (optional) default: True + If True the track density is computed as the number of different tracks crossing a grid + cell. If False, the track density takes into account how long the track stayed in each + grid cell. Hence slower tracks increase the density if the parameter is set to False. + wind_min: float (optional), default: None + Minimum wind speed above which to select tracks (inclusive). + wind_max: float (optional), default: None + Maximal wind speed below which to select tracks (exclusive if wind_min is also provided, + otherwise inclusive). + Returns: + ------- + hist_count: np.ndarray + 2D matrix containing the the absolute count per grid cell of track point or the normalized + number of track points, depending on the norm parameter. + lat_bins: np.ndarray + latitude bins in which the point were counted + lon_bins: np.ndarray + longitude bins in which the point were counted + + Example: + -------- + >>> tc_tracks = TCTrack.from_ibtracs_netcdf("path_to_file") + >>> tc_tracks.equal_timestep(time_step_h = 1) + >>> hist_count, *_ = compute_track_density(tc_track = tc_tracks, res = 1) + >>> ax = plot_track_density(hist_count) + + """ + + limit_ratio: float = 1.12 * 1.1 # record tc speed 112km/h -> 1.12°/h + 10% margin + time_value: float = tc_track.data[0].time_step[0].values.astype(float) + + if time_value > (res / limit_ratio): + warnings.warn( + "The time step is too big for the current resolution. For the desired resolution, \n" + f"apply a time step of {res/limit_ratio}h." + ) + elif res < 0.1: + warnings.warn( + "The resolution is too high. The computation might take several minutes \n" + "to hours. Consider using a resolution below 0.1°." + ) + + # define grid resolution and bounds for density computation + if not bounds: + lon_min, lat_min, lon_max, lat_max = -180, -90, 180, 90 + else: + lon_min, lat_min, lon_max, lat_max = bounds[0], bounds[1], bounds[2], bounds[3] + + lat_bins: np.ndarray = np.linspace(lat_min, lat_max, int(180 / res)) + lon_bins: np.ndarray = np.linspace(lon_min, lon_max, int(360 / res)) + + # compute 2D density + if genesis: + hist_count = compute_genesis_density( + tc_track=tc_track, lat_bins=lat_bins, lon_bins=lon_bins + ) + else: + hist_count = np.zeros((len(lat_bins) - 1, len(lon_bins) - 1)) + for track in tqdm(tc_track.data, desc="Processing Tracks"): + + # select according to wind speed + wind_speed = track.max_sustained_wind.values + if wind_min and wind_max: + index = np.where((wind_speed >= wind_min) & (wind_speed < wind_max))[0] + elif wind_min and not wind_max: + index = np.where(wind_speed >= wind_min)[0] + elif wind_max and not wind_min: + index = np.where(wind_speed <= wind_max)[0] + else: + index = slice(None) # select all the track + + # compute 2D density + hist_new, _, _ = np.histogram2d( + track.lat.values[index], + track.lon.values[index], + bins=[lat_bins, lon_bins], + density=False, + ) + if filter_tracks: + hist_new[hist_new > 1] = 1 + + hist_count += hist_new + + if norm: + hist_count = normalize_hist(res=res, hist_count=hist_count, norm=norm) + + return hist_count, lat_bins, lon_bins + + +def compute_genesis_density( + tc_track: TCTracks, lat_bins: np.ndarray, lon_bins: np.ndarray +) -> np.ndarray: + """Compute the density of track genesis locations. This function works under the hood + of :py:meth:`compute_track_density`. If it is called with the parameter genesis = True, + the function return the number of genesis points per grid cell. + + Parameters: + ----------- + + tc_track: TCT track object + TC track object containing a list of all tracks. + lat_bins: 1D np.array + array containg the latitude bins. + lon_bins: 1D np.array + array containg the longitude bins. + + Returns: + -------- + hist_count: 2D np.array + array containing the number of genesis points per grid cell + """ + # Extract the first lat and lon from each dataset + first_lats = np.array([ds.lat.values[0] for ds in tc_track.data]) + first_lons = np.array([ds.lon.values[0] for ds in tc_track.data]) + + # compute 2D density of genesis points + hist_count, _, _ = np.histogram2d( + first_lats, + first_lons, + bins=[lat_bins, lon_bins], + density=False, + ) + + return hist_count + + +def normalize_hist( + res: int, + hist_count: np.ndarray, + norm: str, +) -> np.ndarray: + """Normalize the number of points per grid cell by the grid cell area or by the total sum of + the histogram count. + + Parameters: + ----------- + res: float + resolution of grid cells + hist_count: 2D np.array + Array containing the count of tracks or genesis locations + norm: str + if norm = area normalize by gird cell area, if norm = sum normalize by the total sum + Returns: + --------- + + """ + + if norm == "area": + grid_area = u_coord.compute_grid_cell_area(res=res) + norm_hist: np.ndarray = hist_count / grid_area + elif norm == "sum": + norm_hist: np.ndarray = hist_count / hist_count.sum() + else: + raise ValueError( + "Invalid value for input parameter 'norm':\n" + "it should be either 'area' or 'sum'" + ) + + return norm_hist diff --git a/climada/hazard/test/test_tc_tracks.py b/climada/hazard/test/test_tc_tracks.py index f1943c7f3..41525adac 100644 --- a/climada/hazard/test/test_tc_tracks.py +++ b/climada/hazard/test/test_tc_tracks.py @@ -1359,6 +1359,111 @@ def test_track_land_params(self): on_land, ) + def test_compute_density_tracks(self): + """Test `compute_track_density` to ensure proper density count.""" + # create track + track = xr.Dataset( + { + "time_step": ("time", np.timedelta64(1, "h") * np.arange(4)), + "max_sustained_wind": ("time", [10, 20, 30, 20]), + "central_pressure": ("time", [1, 1, 1, 1]), + "radius_max_wind": ("time", [1, 1, 1, 1]), + "environnmental_pressure": ("time", [1, 1, 1, 1]), + "basin": ("time", ["NA", "NA", "NA", "NA"]), + }, + coords={ + "time": ("time", pd.date_range("2025-01-01", periods=4, freq="12H")), + "lat": ("time", [-90, -90, -90, -90]), + "lon": ("time", [-179, -169, -159, -149]), + }, + attrs={ + "max_sustained_wind_unit": "m/s", + "central_pressure_unit": "hPa", + "name": "storm_0", + "sid": "0", + "orig_event_flag": True, + "data_provider": "FAST", + "id_no": "0", + "category": "1", + }, + ) + + tc_tracks = tc.TCTracks([track]) + + hist_abs, *_ = tc.compute_track_density( + tc_tracks, + res=10, + norm=None, + ) + # hist_norm, *_ = tc.compute_track_density(tc_tracks, res=10, density=True) + hist_wind_min, *_ = tc.compute_track_density( + tc_tracks, res=10, norm=None, wind_min=11, wind_max=None + ) + hist_wind_max, *_ = tc.compute_track_density( + tc_tracks, res=10, norm=None, wind_min=None, wind_max=30 + ) + hist_wind_max, *_ = tc.compute_track_density( + tc_tracks, res=10, norm=None, wind_min=None, wind_max=30 + ) + hist_wind_both, *_ = tc.compute_track_density( + tc_tracks, res=10, norm=None, wind_min=11, wind_max=29 + ) + self.assertEqual(hist_abs.shape, (17, 35)) + self.assertEqual(hist_abs.sum(), 4) + self.assertEqual(hist_wind_min.sum(), 3) + self.assertEqual(hist_wind_max.sum(), 4) + self.assertEqual(hist_wind_both.sum(), 2) + # the track defined above occupy positions 0 to 4 + np.testing.assert_array_equal(hist_abs[0, 0:4], [1, 1, 1, 1]) + + def test_normalize_hist(self): + """test the correct normalization of a 2D matrix by grid cell area and sum of + values of the matrix.""" + + M = np.ones((10, 10)) + M_norm = tc.normalize_hist(res=10, hist_count=M, norm="sum") + np.testing.assert_array_equal(M, M_norm * 100) + + def test_compute_genesis_density(self): + """Check that the correct number of grid point is computed per grid cell for the starting + location of cyclone tracks""" + # create track + track = xr.Dataset( + { + "time_step": ("time", np.timedelta64(1, "h") * np.arange(4)), + "max_sustained_wind": ("time", [10, 20, 30, 20]), + "central_pressure": ("time", [1, 1, 1, 1]), + "radius_max_wind": ("time", [1, 1, 1, 1]), + "environnmental_pressure": ("time", [1, 1, 1, 1]), + "basin": ("time", ["NA", "NA", "NA", "NA"]), + }, + coords={ + "time": ("time", pd.date_range("2025-01-01", periods=4, freq="12H")), + "lat": ("time", [-90, -89, -88, -87]), + "lon": ("time", [-179, -169, -159, -149]), + }, + attrs={ + "max_sustained_wind_unit": "m/s", + "central_pressure_unit": "hPa", + "name": "storm_0", + "sid": "0", + "orig_event_flag": True, + "data_provider": "FAST", + "id_no": "0", + "category": "1", + }, + ) + res = 10 + tc_tracks = tc.TCTracks([track]) + lat_bins = np.linspace(-90, 90, int(180 / res)) + lon_bins = np.linspace(-180, 180, int(360 / res)) + hist = tc.compute_genesis_density( + tc_track=tc_tracks, lat_bins=lat_bins, lon_bins=lon_bins + ) + self.assertEqual(hist.shape, (17, 35)) + self.assertEqual(hist.sum(), 1) + self.assertEqual(hist[0, 0], 1) + # Execute Tests if __name__ == "__main__": diff --git a/climada/util/coordinates.py b/climada/util/coordinates.py index 5b2a72cf6..78feec287 100644 --- a/climada/util/coordinates.py +++ b/climada/util/coordinates.py @@ -46,7 +46,8 @@ import shapely.vectorized import shapely.wkt from cartopy.io import shapereader -from shapely.geometry import MultiPolygon, Point, box +from pyproj import Geod +from shapely.geometry import MultiPolygon, Point, Polygon, box from sklearn.neighbors import BallTree import climada.util.hdf5_handler as u_hdf5 @@ -459,6 +460,136 @@ def get_gridcellarea(lat, resolution=0.5, unit="ha"): return area +def compute_grid_cell_area_validation(res: float) -> tuple[np.ndarray]: + """ + This function computes the area of each grid cell on a sphere (Earth), using latitude and + longitude bins based on the given resolution. The area is computed using the spherical cap + approximation for each grid cell. The function return a 2D matrix with the corresponding area. + + The formula used to compute the area of each grid cell is derived from the integral of the + surface area of a spherical cap between squared latitude and longitude bins: + + A = R**2 * (Δλ) * (sin(ϕ2) - sin(ϕ1)) + + Where: + - R is the radius of the Earth (in km). + - Δλ is the width of the grid cell in terms of longitude (in radians). + - sin(ϕ2) - sin(ϕ1) is the difference in the sine of the latitudes, which + accounts for the varying horizontal distance between longitudinal lines at different latitudes. + + This formula is the direct integration over λ and ϕ2 of: + + A = R**2 * Δλ * Δϕ * cos(ϕ1) + + which approximate the grid cell area as a square computed by the multiplication of two + arc legths, with the logitudal arc length ajdusted by latitude. + + Parameters: + ---------- + res: int + The resolution of the grid in degrees. The grid will have cells of size `res x res` in + terms of latitude and longitude. + + Returns: + ------- + grid_area: np.ndarray + A 2D array of shape `(len(lat_bins)-1, len(lon_bins)-1)` containing the area of each grid + cell, in square kilometers. Each entry in the array corresponds to the area of the + corresponding grid cell on Earth. + + Example: + -------- + >>> res = 10 #10° resolution + >>> area = compute_grid_cell_area(res = res) + >>> print(area.shape) # (180/10, 360/10) + (18, 36) + """ + + lat_bins = np.linspace(-90, 90, int(180 / res)) # lat bins + lon_bins = np.linspace(-180, 180, int(360 / res)) + + r_earth = 6371 # Earth's radius [km] + # Convert to radians + lat_bin_edges = np.radians(lat_bins) + lon_res_rad = np.radians(res) + + # Compute area + areas = ( + r_earth**2 + * lon_res_rad + * np.abs(np.sin(lat_bin_edges[1:]) - np.sin(lat_bin_edges[:-1])) + ) + # Broadcast to create a full 2D grid of areas + grid_area = np.tile( + areas[:, np.newaxis], (1, len(lon_bins) - 1) + ) # each row as same area + + return grid_area, lat_bins, lon_bins + + +def compute_grid_cell_area( + res: float = None, + bounds: tuple = None, + projection: str = "WGS84", + units: str = "km^2", +) -> np.ndarray: + """ + Compute the area of each grid cell in a latitude-longitude grid. + + Parameters: + ----------- + res: float + Grid resolution in degrees + bounds: tuple, dafault: None + (lon_min, lat_min, lon_max, lat_max) latitude and longitude bounds to compute grid cell area + projection: str + Ellipsoid or spherical projection to approximate Earth. To get the complete list of + projections call :py:meth:`pyproj.get_ellps_map()`. Widely used projections: + - "WGS84": Uses the WGS84 ellipsoid (default) + - "GRS80" + - "IAU76" + - "sphere": Uses a perfect sphere with Earth's mean radius (6371 km) + units: str (optional) Default "km^2" + units of the area. Either km^2 or m^2. + + Returns: + -------- + area: np.ndarray + A 2D numpy array of grid cell areas in km² + Example: + -------- + >>> area = compute_grid_areas(res = 1, projection ="sphere", units = "m^2") + """ + geod = Geod(ellps=projection) + + if not bounds: + lon_min, lat_min, lon_max, lat_max = -180, -90, 180, 90 + else: + lon_min, lat_min, lon_max, lat_max = bounds[0], bounds[1], bounds[2], bounds[3] + + lat_edges: np.ndarray = np.linspace(lat_min, lat_max, int(180 / res)) + lon_edges: np.ndarray = np.linspace(lon_min, lon_max, int(360 / res)) + + area = np.zeros((len(lat_edges) - 1, len(lon_edges) - 1)) + + # Iterate over consecutive latitude and longitude edges + for i, (lat1, lat2) in enumerate(zip(lat_edges[:-1], lat_edges[1:])): + for j, (lon1, lon2) in enumerate(zip(lon_edges[:-1], lon_edges[1:])): + + # 5th point to close the loop + poly_lons = [lon1, lon2, lon2, lon1, lon1] + poly_lats = [lat1, lat1, lat2, lat2, lat1] + + # Compute the area of the grid cell + poly_area, _ = geod.polygon_area_perimeter(poly_lons, poly_lats) + area[i, j] = abs(poly_area) # Convert from m² to km² + + if units == "km^2": + area = area / 1e6 + + return area + + def grid_is_regular(coord): """Return True if grid is regular. If True, returns height and width. diff --git a/climada/util/test/test_coordinates.py b/climada/util/test/test_coordinates.py index 56ba3af91..3e335d46e 100644 --- a/climada/util/test/test_coordinates.py +++ b/climada/util/test/test_coordinates.py @@ -565,6 +565,31 @@ def test_get_gridcellarea(self): self.assertAlmostEqual(area2[0], 1781.5973363005) self.assertTrue(area2[0] <= 2500) + def test_compute_grid_area(self): + """Test that the two twin functions calculate the area of a gridcell correctly. Using + an absolute reference and mutual validation by comparison.""" + res = 1 + area = u_coord.compute_grid_cell_area( + res=res, projection="sphere", units="km^2" + ) + area_test, *_ = u_coord.compute_grid_cell_area_validation(res=res) + + self.assertEqual(area.shape, (179, 359)) + self.assertEqual(area_test.shape, (179, 359)) + # check that all rows have equal area with 1e-5 tolerance for relative and absolute precision + for i in range(area.shape[0]): + self.assertTrue( + np.all(np.isclose(area[i, :], area[i, 0], rtol=1e-5, atol=1e-5)) + ) + self.assertTrue( + np.all( + np.isclose(area_test[i, :], area_test[i, 0], rtol=1e-5, atol=1e-5) + ) + ) + + # check that both methods give similar results with 0.01% tolerance in relative difference + self.assertTrue(np.allclose(area, area_test, rtol=1e-2, atol=1e-5)) + def test_read_vector_pass(self): """Test one columns data""" shp_file = shapereader.natural_earth(