Hazard forecasts and impacts forecasts#

In this tutorial, we showcase CLIMADA functionality to handle forecast data using the classes HazardForecast and ImpactForecast. Exemplary forecast data will be downloaded from the Open Government Data (OGD) portal of MeteoSwiss.

Table of Contents#

  1. Quickstart

  2. Reading & processing forecast data

  3. Impact Forecast calculation

  4. Products based on ImpactForecast

# import necessary packages
from pathlib import Path

import cartopy.crs as ccrs
import cartopy.feature as cfeature
import matplotlib.colors as mcolors

# import numbers
import matplotlib.pyplot as plt
import numpy as np
import pandas as pd
import xarray as xr

# Import CLIMADA
from climada.engine import ImpactCalc, ImpactForecast
from climada.entity import ImpactFunc, ImpactFuncSet
from climada.entity.exposures.litpop import LitPop
from climada.hazard.forecast import HazardForecast
from climada.util.api_client import Client

# Features for maps
admin1 = cfeature.NaturalEarthFeature(
    category="cultural", name="admin_1_states_provinces_lines", scale="10m"
)
lakes = cfeature.LAKES.with_scale("10m")
borders = cfeature.BORDERS.with_scale("10m")

0. Quickstart#

The next cell downloads some exemplary forecast data and loads it as a HazardForecast object. It then demonstrates the capabilities of this class:

  • Computing reductions (in this case maximum) over the lead time or member dimensions.

  • Regular impact calculation with ImpactCalc yielding an ImpactForecast object.

  • More reductions and the use of other methods from Impact.

def fix_lead_time(ds: xr.Dataset) -> xr.Dataset:
    """Make lead_time a valid timedelta

    Only needed for the pre-computed demonstrator data from the data API
    """
    ds = ds.copy()  # Shallow copy
    ds["lead_time"] = ds["time"] - ds["forecast_reference_time"]
    return ds
# Load example hazard data
client = Client()
ds_info_icon_fcst = client.get_dataset_info_by_uuid(
    "60229fa0-4d5d-4f39-876c-b4344579d6eb"
)
_, path_to_forecast_nc = client.download_dataset(ds_info_icon_fcst)

# Read data as HazardForecast object
ds = xr.open_dataset(path_to_forecast_nc[0])
ds = fix_lead_time(ds)
hazard_forecast = HazardForecast.from_xarray_raster(
    ds,
    hazard_type="PR",
    intensity_unit="mm/h",
    coordinate_vars={
        "longitude": "x",
        "latitude": "y",
        "lead_time": "lead_time",
        "member": "realization",
    },
    intensity="precipitation_amount_1hsum",
    crs="EPSG:2056",
    open_dataset_kws={"decode_coords": "all"},
)
# transform to EPSG:4326
hazard_forecast.centroids.gdf = hazard_forecast.centroids.gdf.to_crs(epsg=4326)

# For each member, use the maximum intensity over the lead time
hazard_forecast = hazard_forecast.max(dim="lead_time")

# Read exposure for population
exp_che = LitPop.from_population(countries=["CHE", "LIE"], res_arcsec=30)
# Assign centroids because hazard not on regular grid
exp_che.assign_centroids(hazard_forecast, threshold=100)

# Define impact functions
threshold_for_affected = 50
imp_fun_precip = ImpactFunc.from_step_impf(
    (0, threshold_for_affected, 200),
    haz_type="PR",
    impf_id=1,
    intensity_unit="mm/h",
    name="Severe Precipitation Step Function",
)
impf_set = ImpactFuncSet([imp_fun_precip])

# Compute impact (forecast)
impact_forecast = ImpactCalc(exp_che, impf_set, hazard_forecast).impact(
    assign_centroids=False, save_mat=True
)

# Compute quantile over remaining ensemble members and plot
impact_forecast_q075 = impact_forecast.quantile(0.92)
impact_forecast_q075.plot_hexbin_impact_exposure()
2026-06-08 13:58:51,220 - climada.hazard.xarray - WARNING - Failed to read values of 'event' as dates. Hazard.event_name will be empty strings
2026-06-08 13:58:51,220 - climada.hazard.xarray - WARNING - Failed to read values of 'event' as dates or ordinals. Hazard.date will be ones only
2026-06-08 13:58:51,894 - climada.entity.exposures.litpop.litpop - INFO - 
 LitPop: Init Exposure for country: CHE (756)...

2026-06-08 13:58:52,027 - climada.entity.exposures.litpop.litpop - INFO - 
 LitPop: Init Exposure for country: LIE (438)...

2026-06-08 13:58:52,063 - climada.entity.exposures.base - INFO - Hazard type not set in impf_
2026-06-08 13:58:52,063 - climada.entity.exposures.base - INFO - category_id not set.
2026-06-08 13:58:52,063 - climada.entity.exposures.base - INFO - cover not set.
2026-06-08 13:58:52,063 - climada.entity.exposures.base - INFO - deductible not set.
2026-06-08 13:58:52,064 - climada.entity.exposures.base - INFO - centr_ not set.
2026-06-08 13:58:52,064 - climada.entity.exposures.base - INFO - Matching 70525 exposures with 27000 centroids.
2026-06-08 13:58:52,102 - climada.util.coordinates - INFO - No exact centroid match found. Reprojecting coordinates to nearest neighbor closer than the threshold = 100 degree
2026-06-08 13:58:52,120 - climada.entity.exposures.base - INFO - No specific impact function column found for hazard PR. Using the anonymous 'impf_' column.
2026-06-08 13:58:52,123 - climada.engine.impact_calc - INFO - Calculating impact for 69394 assets (>0) and 21 events.
2026-06-08 13:58:52,131 - climada.engine.impact_calc - WARNING - eai_exp and aai_agg are undefined with forecasts. Setting them to NaN arrays.
2026-06-08 13:58:52,144 - climada.engine.impact_forecast - INFO - at_event gives the total impact for one specific combination of member and lead_time.
<GeoAxes: >
../_images/4d96abee552ba17fdd96233a4ccda9e1db0e0879a512e52a32e1a0d88b3063a1.png

1. Reading & processing forecast data#

1.1 Using MeteoSwiss Open Government Data weather forecasts#

This example retrieves the latest OGD forecast data from MeteoSwiss. You can also skip this step and go to section 1.2 for using stored demo data from the CLIMADA API.

Requires the package meteodata-lab, available via pip.

You can find out all about data openly published here (info and notebooks on how to retrieve):

Check out:

  • https://github.com/MeteoSwiss/meteodata-lab and

  • https://github.com/MeteoSwiss/opendata-nwp-demos/blob/main/01_retrieve_process_precip.ipynb

  • https://www.meteoschweiz.admin.ch/service-und-publikationen/service/open-data.html

!pip install meteodata-lab
Requirement already satisfied: meteodata-lab in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (0.7.2)
Requirement already satisfied: click<9.0.0,>=8.1.7 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (8.3.3)
Requirement already satisfied: earthkit-data<1,>=0.11 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.20.0)
Requirement already satisfied: earthkit-meteo>=0.4.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (0.6.2)
Requirement already satisfied: eccodes<2.40,>=2.38 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (2.39.2)
Requirement already satisfied: eccodes-cosmo-resources-python<2.39,>=2.38 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (2.38.3.1)
Requirement already satisfied: numpy<2.4.0,>=1.26.4 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (2.2.6)
Requirement already satisfied: pydantic<3.0,>=2.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (2.13.4)
Requirement already satisfied: pyyaml<7.0.0,>=6.0.1 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (6.0.3)
Requirement already satisfied: setuptools in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (82.0.1)
Requirement already satisfied: xarray>=2024 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from meteodata-lab) (2026.4.0)
Requirement already satisfied: cfgrib>=0.9.10.1 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.9.15.1)
Requirement already satisfied: dask in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2026.3.0)
Requirement already satisfied: deprecation in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2.1.0)
Requirement already satisfied: earthkit-utils<0.99,>=0.2 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.3.0)
Requirement already satisfied: entrypoints in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.4)
Requirement already satisfied: filelock in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (3.29.0)
Requirement already satisfied: jinja2 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (3.1.6)
Requirement already satisfied: jsonschema in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (4.26.0)
Requirement already satisfied: lru-dict in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.4.1)
Requirement already satisfied: markdown in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (3.10.2)
Requirement already satisfied: multiurl>=0.3.3 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.3.7)
Requirement already satisfied: netcdf4 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.7.4)
Requirement already satisfied: pandas in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2.2.2)
Requirement already satisfied: pdbufr>=0.11 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.14.2)
Requirement already satisfied: tqdm>=4.63 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (4.67.3)
Requirement already satisfied: covjsonkit>=0.2.2 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.2.20)
Requirement already satisfied: array-api-compat in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-utils<0.99,>=0.2->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.15.0)
Requirement already satisfied: pint in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from earthkit-utils<0.99,>=0.2->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.25.3)
Requirement already satisfied: attrs in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from eccodes<2.40,>=2.38->meteodata-lab) (26.1.0)
Requirement already satisfied: cffi in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from eccodes<2.40,>=2.38->meteodata-lab) (2.0.0)
Requirement already satisfied: findlibs in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from eccodes<2.40,>=2.38->meteodata-lab) (0.1.2)
Requirement already satisfied: annotated-types>=0.6.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pydantic<3.0,>=2.0->meteodata-lab) (0.7.0)
Requirement already satisfied: pydantic-core==2.46.4 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pydantic<3.0,>=2.0->meteodata-lab) (2.46.4)
Requirement already satisfied: typing-extensions>=4.14.1 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pydantic<3.0,>=2.0->meteodata-lab) (4.15.0)
Requirement already satisfied: typing-inspection>=0.4.2 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pydantic<3.0,>=2.0->meteodata-lab) (0.4.2)
Requirement already satisfied: orjson in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (3.11.9)
Requirement already satisfied: covjson-pydantic in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.7.0)
Requirement already satisfied: conflator in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.1.8)
Requirement already satisfied: scipy in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.17.1)
Requirement already satisfied: python-dateutil>=2.8.2 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pandas->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2.9.0.post0)
Requirement already satisfied: pytz>=2020.1 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pandas->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2026.2)
Requirement already satisfied: tzdata>=2022.7 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pandas->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2026.2)
Requirement already satisfied: requests in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from multiurl>=0.3.3->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2.33.1)
Requirement already satisfied: six>=1.5 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from python-dateutil>=2.8.2->pandas->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.17.0)
Requirement already satisfied: packaging>=24.2 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from xarray>=2024->meteodata-lab) (26.2)
Requirement already satisfied: pycparser in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from cffi->eccodes<2.40,>=2.38->meteodata-lab) (2.22)
Requirement already satisfied: rich-argparse>1.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from conflator->covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.8.0)
Requirement already satisfied: rich>=11.0.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from rich-argparse>1.0->conflator->covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (15.0.0)
Requirement already satisfied: markdown-it-py>=2.2.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from rich>=11.0.0->rich-argparse>1.0->conflator->covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (4.2.0)
Requirement already satisfied: pygments<3.0.0,>=2.13.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from rich>=11.0.0->rich-argparse>1.0->conflator->covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2.20.0)
Requirement already satisfied: mdurl~=0.1 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from markdown-it-py>=2.2.0->rich>=11.0.0->rich-argparse>1.0->conflator->covjsonkit>=0.2.2->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.1.2)
Requirement already satisfied: cloudpickle>=3.0.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from dask->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (3.1.2)
Requirement already satisfied: fsspec>=2021.09.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from dask->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2026.4.0)
Requirement already satisfied: partd>=1.4.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from dask->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.4.2)
Requirement already satisfied: toolz>=0.12.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from dask->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.1.0)
Requirement already satisfied: locket in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from partd>=1.4.0->dask->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.0.0)
Requirement already satisfied: MarkupSafe>=2.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from jinja2->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (3.0.3)
Requirement already satisfied: jsonschema-specifications>=2023.03.6 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from jsonschema->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2025.9.1)
Requirement already satisfied: referencing>=0.28.4 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from jsonschema->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.37.0)
Requirement already satisfied: rpds-py>=0.25.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from jsonschema->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2026.5.1)
Requirement already satisfied: cftime in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from netcdf4->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (1.6.5)
Requirement already satisfied: certifi in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from netcdf4->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2026.4.22)
Requirement already satisfied: flexcache>=0.3 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pint->earthkit-utils<0.99,>=0.2->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.3)
Requirement already satisfied: flexparser>=0.4 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pint->earthkit-utils<0.99,>=0.2->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (0.4)
Requirement already satisfied: platformdirs>=2.1.0 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from pint->earthkit-utils<0.99,>=0.2->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (4.9.6)
Requirement already satisfied: charset_normalizer<4,>=2 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from requests->multiurl>=0.3.3->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (3.4.7)
Requirement already satisfied: idna<4,>=2.5 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from requests->multiurl>=0.3.3->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (3.13)
Requirement already satisfied: urllib3<3,>=1.26 in /Users/ldr.riedel/miniforge3/envs/climada_env_3.12/lib/python3.12/site-packages (from requests->multiurl>=0.3.3->earthkit-data<1,>=0.11->earthkit-data[covjsonkit]<1,>=0.11->meteodata-lab) (2.6.3)
# optional packages for this section only
from earthkit.data import config
from meteodatalab import ogd_api
from meteodatalab.operators import time_operators as time_ops

config.set("cache-policy", "user")
# Define target lead times for ICON-CH2-EPS: up 5d, here +0h, to +12h
lead_times = [f"P0DT{i}H" for i in np.arange(0, 13)]
model = "ogd-forecasting-icon-ch2"
nwp_vars = ["TOT_PREC", "VMAX_10M"]  # precip. and wind; temperature: "T_2M"
reftime = "latest"  # also try: 'latest'
dict_ens_data = {}
for var in nwp_vars:
    req = ogd_api.Request(
        collection=model,
        variable=var,
        ref_time=reftime,
        perturbed=True,  # ensembles
        lead_time=lead_times,
    )
    dict_ens_data[var] = ogd_api.get_from_ogd(req)

# We extract the hourly delta in total precipitation to receive the hourly precipitation
dict_ens_data["TOT_PREC"] = time_ops.delta(
    dict_ens_data["TOT_PREC"], np.timedelta64(1, "h")
)
                                                                                                                                                                     
ds_icon_latest = xr.Dataset(dict_ens_data)
init_time_str = np.datetime_as_string(ds_icon_latest.ref_time.values[0], unit="h")

After loading the data as a netcdf file, we read it into a HazardForecast object (see next subsection for more details on the structuer of this class).

haz_fc_icon_latest = HazardForecast.from_xarray_raster(
    ds_icon_latest,
    hazard_type="PR",
    intensity_unit="mm/h",
    coordinate_vars={
        "longitude": "lon",
        "latitude": "lat",
        "lead_time": "lead_time",
        "member": "eps",
    },
    intensity="TOT_PREC",
)
2026-06-08 13:59:21,773 - climada.hazard.xarray - WARNING - Failed to read values of 'event' as dates. Hazard.event_name will be empty strings
2026-06-08 13:59:21,774 - climada.hazard.xarray - WARNING - Failed to read values of 'event' as dates or ordinals. Hazard.date will be ones only
# Plot the maximum precipitation intensity for each pixel in the entire dataset
ax = haz_fc_icon_latest.plot_intensity(event=0)
lead_times = np.unique(haz_fc_icon_latest.lead_time).astype("timedelta64[h]")
ax.set_title(
    f"Forecasted hourly precipication (model: {model}) for init time {init_time_str} "
    f"(UTC),\nmaximized over {len(np.unique(haz_fc_icon_latest.member))} ensemble "
    f"members and {len(lead_times)} lead times ({lead_times[0]} - {lead_times[-1]})"
)
Text(0.5, 1.0, 'Forecasted hourly precipication (model: ogd-forecasting-icon-ch2) for init time 2026-06-08T06 (UTC),\nmaximized over 20 ensemble members and 13 lead times (0 hours - 12 hours)')
../_images/0296f81b892da09bbb8f653a06dccaaa82717924a8bea958c1906a17d9710f84.png

1.2 Using CLIMADA API demo data#

As an example, we load existing forecast data for an event including heavy precipitation and flooding on June 21, 2024 in in Misox, Switzerland.

After loading the data as a netcdf file, we read it into a HazardForecast object.

# get the demo forecast netCDF file fcst_icon_seamless_2024-06-18_12_00.nc via CLIMADA API
client = Client()
ds_info_icon_fcst = client.get_dataset_info_by_uuid(
    "60229fa0-4d5d-4f39-876c-b4344579d6eb"
)
_, path_to_forecast_nc = client.download_dataset(ds_info_icon_fcst)

Load into HazardForecast file from NetCDF on disk.

ds = xr.open_dataset(path_to_forecast_nc[0])
ds = fix_lead_time(ds)
hazard_forecast = HazardForecast.from_xarray_raster(
    ds,
    hazard_type="PR",
    intensity_unit="mm/h",
    coordinate_vars={
        "longitude": "x",
        "latitude": "y",
        "lead_time": "lead_time",
        "member": "realization",
    },
    intensity="precipitation_amount_1hsum",
    crs="EPSG:2056",
)
# transform to EPSG:4326
hazard_forecast.centroids.gdf = hazard_forecast.centroids.gdf.to_crs(epsg=4326)
2026-06-08 13:59:40,792 - climada.hazard.xarray - WARNING - Failed to read values of 'event' as dates. Hazard.event_name will be empty strings
2026-06-08 13:59:40,793 - climada.hazard.xarray - WARNING - Failed to read values of 'event' as dates or ordinals. Hazard.date will be ones only

HazardForecast objects include forecasts produced at a specific initializatio (or reference) time, but for several several lead times and several forecast ensemble members. For this, it includes the two array-like attributes:

  • member to indicate the forecast ensemble member, and

  • lead_time to indicate the lead time of the data,

both of which can be used to select specific forecasts.

pd.DataFrame.from_records(
    {
        "member": hazard_forecast.member,
        "lead_time": hazard_forecast.lead_time,
        "event_name": hazard_forecast.event_name,
        "event_id": hazard_forecast.event_id,
    }
).head(25)
event_id event_name lead_time member
0 1 lt_0h_m_0 0 days 00:00:00 0
1 2 lt_0h_m_1 0 days 00:00:00 1
2 3 lt_0h_m_2 0 days 00:00:00 2
3 4 lt_0h_m_3 0 days 00:00:00 3
4 5 lt_0h_m_4 0 days 00:00:00 4
5 6 lt_0h_m_5 0 days 00:00:00 5
6 7 lt_0h_m_6 0 days 00:00:00 6
7 8 lt_0h_m_7 0 days 00:00:00 7
8 9 lt_0h_m_8 0 days 00:00:00 8
9 10 lt_0h_m_9 0 days 00:00:00 9
10 11 lt_0h_m_10 0 days 00:00:00 10
11 12 lt_0h_m_11 0 days 00:00:00 11
12 13 lt_0h_m_12 0 days 00:00:00 12
13 14 lt_0h_m_13 0 days 00:00:00 13
14 15 lt_0h_m_14 0 days 00:00:00 14
15 16 lt_0h_m_15 0 days 00:00:00 15
16 17 lt_0h_m_16 0 days 00:00:00 16
17 18 lt_0h_m_17 0 days 00:00:00 17
18 19 lt_0h_m_18 0 days 00:00:00 18
19 20 lt_0h_m_19 0 days 00:00:00 19
20 21 lt_0h_m_20 0 days 00:00:00 20
21 22 lt_1h_m_0 0 days 01:00:00 0
22 23 lt_1h_m_1 0 days 01:00:00 1
23 24 lt_1h_m_2 0 days 01:00:00 2
24 25 lt_1h_m_3 0 days 01:00:00 3
# select hazard forecast for a specific member and a specific lead time
test_member = 10
test_lead_time = np.timedelta64(78, "h")

# plot
ax = hazard_forecast.select(
    member=test_member,
    lead_time=test_lead_time,
).plot_intensity(0)
ax.set_title(
    f"Precipation Forecast for member {test_member} and lead time "
    f"{test_lead_time.astype('timedelta64[h]')}"
)
Text(0.5, 1.0, 'Precipation Forecast for member 10 and lead time 78 hours')
../_images/bb883b4f1295ec002186d74a9b9e60af03f8eb93ab846628960189c82bc1ac6f.png

The two forecast attributes member and lead_time can be used to generate (centroid-wise) statistics by means of methods like, e.g.,

  • .max() for the centroid-wise maximum of all lead times and ensemble members

  • .quantile(q=0.3, dim="member") for the 0.3-quantile of the hazard intensity over all ensemble members, for each lead time and centroid

# Reduce over dimension 'lead_time'
haz_max = hazard_forecast.max(dim="lead_time")
haz_max.member
array([ 0,  1,  2,  3,  4,  5,  6,  7,  8,  9, 10, 11, 12, 13, 14, 15, 16,
       17, 18, 19, 20])

This operation reduces over the dimension ‘lead_time’. The associated values are therefore set to placeholders; in this case, Not-a-Time (NaT).

haz_max.lead_time
array(['NaT', 'NaT', 'NaT', 'NaT', 'NaT', 'NaT', 'NaT', 'NaT', 'NaT',
       'NaT', 'NaT', 'NaT', 'NaT', 'NaT', 'NaT', 'NaT', 'NaT', 'NaT',
       'NaT', 'NaT', 'NaT'], dtype=timedelta64)
# Plot the centroid-wise maximum over all lead times for 4 select members
fig, axes = plt.subplots(2, 2, subplot_kw={"projection": ccrs.Mercator()})
norm = mcolors.Normalize(vmin=0)
for idx, member in enumerate((0, 6, 12, 18)):
    ax = axes.flat[idx]
    haz_max.select(member=member).plot_intensity(0, axis=ax, norm=norm)
    ax.set_title(f"Member {member}")
../_images/d8c572329b11450115d2684d92169833928053d58e90680e22db7b265cf48d26.png

We can also reduce over the ‘member’, leaving us with the ensemble median for each lead time.

haz_med = hazard_forecast.median(dim="member")
haz_med.member, haz_med.lead_time
(array([-1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
        -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
        -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
        -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
        -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
        -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1,
        -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1, -1]),
 array([              0,   3600000000000,   7200000000000,  10800000000000,
         14400000000000,  18000000000000,  21600000000000,  25200000000000,
         28800000000000,  32400000000000,  36000000000000,  39600000000000,
         43200000000000,  46800000000000,  50400000000000,  54000000000000,
         57600000000000,  61200000000000,  64800000000000,  68400000000000,
         72000000000000,  75600000000000,  79200000000000,  82800000000000,
         86400000000000,  90000000000000,  93600000000000,  97200000000000,
        100800000000000, 104400000000000, 108000000000000, 111600000000000,
        115200000000000, 118800000000000, 122400000000000, 126000000000000,
        129600000000000, 133200000000000, 136800000000000, 140400000000000,
        144000000000000, 147600000000000, 151200000000000, 154800000000000,
        158400000000000, 162000000000000, 165600000000000, 169200000000000,
        172800000000000, 176400000000000, 180000000000000, 183600000000000,
        187200000000000, 190800000000000, 194400000000000, 198000000000000,
        201600000000000, 205200000000000, 208800000000000, 212400000000000,
        216000000000000, 219600000000000, 223200000000000, 226800000000000,
        230400000000000, 234000000000000, 237600000000000, 241200000000000,
        244800000000000, 248400000000000, 252000000000000, 255600000000000,
        259200000000000, 262800000000000, 266400000000000, 270000000000000,
        273600000000000, 277200000000000, 280800000000000, 284400000000000,
        288000000000000, 291600000000000, 295200000000000, 298800000000000,
        302400000000000, 306000000000000, 309600000000000, 313200000000000,
        316800000000000, 320400000000000, 324000000000000, 327600000000000,
        331200000000000, 334800000000000, 338400000000000, 342000000000000,
        345600000000000, 349200000000000, 352800000000000, 356400000000000,
        360000000000000, 363600000000000, 367200000000000, 370800000000000,
        374400000000000, 378000000000000, 381600000000000, 385200000000000,
        388800000000000, 392400000000000, 396000000000000, 399600000000000,
        403200000000000, 406800000000000, 410400000000000],
       dtype='timedelta64[ns]'))
# plot median over member dimension for given a lead time
ax = haz_med.select(lead_time=test_lead_time).plot_intensity(0)
ax.set_title(
    f"Median of members for Precipitation Forecast at lead time {test_lead_time}"
)
Text(0.5, 1.0, 'Median of members for Precipitation Forecast at lead time 78 hours')
../_images/cdb6273383f41f73e74e0cca2b6ea543d014eae9ac85acc395e630979a22e32a.png

While most meteorological data come as netCDF (.nc), HazardForecast objects can be written to and read from HDF5 (.h5) files using the methods:

  • write_hdf5()

  • from_hdf5()

from tempfile import TemporaryDirectory

# write hazard forecast object to HDF5 file
with TemporaryDirectory() as tmpdir:
    hazard_forecast.write_hdf5(Path(tmpdir) / "demo_hazard_forecast.h5")

    # read hazard forecast from HDF5 file
    hazard_forecast_read = HazardForecast.from_hdf5(
        Path(tmpdir) / "demo_hazard_forecast.h5"
    )

np.array_equal(
    hazard_forecast_read.intensity.toarray(), hazard_forecast.intensity.toarray()
)
2026-06-08 13:59:55,005 - climada.hazard.io - INFO - Writing /var/folders/9m/w030_m9j2j5gyjdjczm2tf680000gq/T/tmp87icnb4i/demo_hazard_forecast.h5
2026-06-08 13:59:55,145 - climada.hazard.centroids.centr - INFO - Writing /var/folders/9m/w030_m9j2j5gyjdjczm2tf680000gq/T/tmp87icnb4i/demo_hazard_forecast.h5
2026-06-08 13:59:55,270 - climada.hazard.io - INFO - Reading /var/folders/9m/w030_m9j2j5gyjdjczm2tf680000gq/T/tmp87icnb4i/demo_hazard_forecast.h5
True

2. Impact Forecast calculation#

HazardForecast objects can be handed in to the standard CLIMADA ImpactCalc engine to produce ImpactForecast objects, consisting of one impact forecast per hazard forecast.

We exemplify this below by loading a population exposure layer for Switzerland and Liechtenstein, and by defining a step impact function. This impact function results in impacts describing how much population (how many people) are exposed to a precipitation of more than 50mm/h for the severe precipitation event in Misox loaded above.

# Read exposure for population
exp_che = LitPop.from_population(countries=["CHE", "LIE"], res_arcsec=30)
# Assign centroids because hazard not on regular grid
exp_che.assign_centroids(hazard_forecast, threshold=100)

# define impact functions
threshold_for_affected = 50
imp_fun_precip = ImpactFunc.from_step_impf(
    (0, threshold_for_affected, 200),
    haz_type="PR",
    impf_id=1,
    intensity_unit="mm/h",
    name="Severe Precipitation Step Function",
)
impf_set = ImpactFuncSet([imp_fun_precip])
impf_set.plot()
2026-06-08 13:59:55,450 - climada.entity.exposures.litpop.litpop - INFO - 
 LitPop: Init Exposure for country: CHE (756)...

2026-06-08 13:59:55,586 - climada.entity.exposures.litpop.litpop - INFO - 
 LitPop: Init Exposure for country: LIE (438)...

2026-06-08 13:59:55,622 - climada.entity.exposures.base - INFO - Hazard type not set in impf_
2026-06-08 13:59:55,623 - climada.entity.exposures.base - INFO - category_id not set.
2026-06-08 13:59:55,623 - climada.entity.exposures.base - INFO - cover not set.
2026-06-08 13:59:55,623 - climada.entity.exposures.base - INFO - deductible not set.
2026-06-08 13:59:55,623 - climada.entity.exposures.base - INFO - centr_ not set.
2026-06-08 13:59:55,630 - climada.entity.exposures.base - INFO - Matching 70525 exposures with 27000 centroids.
2026-06-08 13:59:55,667 - climada.util.coordinates - INFO - No exact centroid match found. Reprojecting coordinates to nearest neighbor closer than the threshold = 100 degree
<Axes: title={'center': 'PR 1: Severe Precipitation Step Function'}, xlabel='Intensity (mm/h)', ylabel='Impact (%)'>
../_images/a5a14608214c7f6881134fee31c62305647cd78baa7bed1423d85b7b78abe857.png
# calculate impact
impact_forecast = ImpactCalc(exp_che, impf_set, hazard_forecast).impact(
    assign_centroids=False,  # assign_centroid=False because we do not want to reassign centroids
)  # impact matrix is saved by default
2026-06-08 13:59:55,763 - climada.entity.exposures.base - INFO - No specific impact function column found for hazard PR. Using the anonymous 'impf_' column.
2026-06-08 13:59:55,766 - climada.engine.impact_calc - INFO - Calculating impact for 69394 assets (>0) and 2415 events.
2026-06-08 13:59:55,930 - climada.engine.impact_calc - WARNING - eai_exp and aai_agg are undefined with forecasts. Setting them to NaN arrays.
pd.DataFrame.from_records(
    {
        "member": impact_forecast.member,
        "lead_time": impact_forecast.lead_time,
        "event_name": impact_forecast.event_name,
        "event_id": impact_forecast.event_id,
    }
).head(25)
event_id event_name lead_time member
0 1 lt_0h_m_0 0 days 00:00:00 0
1 2 lt_0h_m_1 0 days 00:00:00 1
2 3 lt_0h_m_2 0 days 00:00:00 2
3 4 lt_0h_m_3 0 days 00:00:00 3
4 5 lt_0h_m_4 0 days 00:00:00 4
5 6 lt_0h_m_5 0 days 00:00:00 5
6 7 lt_0h_m_6 0 days 00:00:00 6
7 8 lt_0h_m_7 0 days 00:00:00 7
8 9 lt_0h_m_8 0 days 00:00:00 8
9 10 lt_0h_m_9 0 days 00:00:00 9
10 11 lt_0h_m_10 0 days 00:00:00 10
11 12 lt_0h_m_11 0 days 00:00:00 11
12 13 lt_0h_m_12 0 days 00:00:00 12
13 14 lt_0h_m_13 0 days 00:00:00 13
14 15 lt_0h_m_14 0 days 00:00:00 14
15 16 lt_0h_m_15 0 days 00:00:00 15
16 17 lt_0h_m_16 0 days 00:00:00 16
17 18 lt_0h_m_17 0 days 00:00:00 17
18 19 lt_0h_m_18 0 days 00:00:00 18
19 20 lt_0h_m_19 0 days 00:00:00 19
20 21 lt_0h_m_20 0 days 00:00:00 20
21 22 lt_1h_m_0 0 days 01:00:00 0
22 23 lt_1h_m_1 0 days 01:00:00 1
23 24 lt_1h_m_2 0 days 01:00:00 2
24 25 lt_1h_m_3 0 days 01:00:00 3

3. Products based on ImpactForecast#

We can now analyse the impact forecasts by selecting specific ensemble member or lead times (or other impact attributes), or by computing statistics (min, max, quantile, mean) over the forecast dimensions member and lead_time.

For the specific case of Misox, we will zoom in to the canton of Ticino which was most heavily affected.

3.1 Per ensemble member, per lead time#

# select for a specific member and a specific lead time and plot Ticino extent

impact_forecast_1m_1lt = impact_forecast.select(
    member=test_member,
    lead_time=test_lead_time,
)
extent_ticino = [8.2, 9.3, 45.58, 47]

fig = plt.figure()
ax = plt.axes(projection=ccrs.Mercator())
ax.add_feature(borders, linewidth=1)
ax.add_feature(admin1, edgecolor="gray", linewidth=1.0, facecolor="none")
ax.add_feature(lakes, facecolor="lightblue")

# plot impact forecast
impact_forecast_1m_1lt.plot_hexbin_impact_exposure(
    event_id=impact_forecast_1m_1lt.event_id, pop_name=False, axis=ax
)
ax.set_extent(extent_ticino)
ax.set_title(
    f"Affected people from precip>{threshold_for_affected}mm/h\nfor member "
    f"{test_member} and lead time {test_lead_time.astype('timedelta64[h]')}"
)
2026-06-08 13:59:55,979 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
Text(0.5, 1.0, 'Affected people from precip>50mm/h\nfor member 10 and lead time 78 hours')
../_images/f437c1ee08e7a4dc845b27da190f74dde1dabb4fb710d77df175cd0fb01cf214.png

3.2 All ensemble members, per lead time#

We can now easily compute the centroid-wise maximum impact over all lead times for each ensemble member.

imp_fct_all_mem_1lt = impact_forecast.max(dim="lead_time")
fig, ax = plt.subplots(
    3, 7, figsize=(20, 10), subplot_kw={"projection": ccrs.Mercator()}
)
fig.suptitle(f"Affected people predicted by each member, at lead time {test_lead_time}")
ax = ax.flatten()
for j, member in enumerate(imp_fct_all_mem_1lt.member):
    ax[j].add_feature(borders, linewidth=1)
    ax[j].add_feature(admin1, edgecolor="gray", linewidth=1.0, facecolor="none")
    ax[j].add_feature(lakes, facecolor="lightblue")

    # plot impact forecast
    imp_fct_all_mem_1lt.plot_hexbin_impact_exposure(
        event_id=imp_fct_all_mem_1lt.event_id[j],
        pop_name=False,
        axis=ax[j],
        vmin=0,
        vmax=110,
    )
    ax[j].set_extent(extent_ticino)
    ax[j].set_title(f"Member {j}")
../_images/6297c91d47bf58ab3217cf22e562f8cf91601d17ba94c0eddb1ddcea1f24c259.png

Next, we can analyse the evolution of predicted number of affected people over Switzerland (member distribution shown as a histogram) with increasing forecast lead time.

3.3 All ensemble members and lead times#

lead_time_window = np.timedelta64(1, "h") * np.arange(69, 87, 2)
n_lead_times = len(lead_time_window)

fig, ax = plt.subplots(
    nrows=n_lead_times,
    ncols=1,
    figsize=(4, 1.2 * n_lead_times),
    sharex=True,
    sharey=True,
)
ax = ax.flatten()
log_histogram_bins = np.arange(7)

for j, lead_time in enumerate(lead_time_window):
    imp_fct_all_mem_1lt = impact_forecast.select(lead_time=lead_time)
    imp_per_member = imp_fct_all_mem_1lt.imp_mat.sum(axis=1)
    ax[j].hist(np.log10(imp_per_member + 1), bins=log_histogram_bins, alpha=1)
    ax[j].set_yticks([4, 8, 12, 16])
    ax[j].text(
        0.98,
        0.90,
        f"lead time: {lead_time.astype('timedelta64[h]')}",
        ha="right",
        va="top",
        transform=ax[j].transAxes,
    )

ax[-1].set_xlabel("")

bin_centers = 0.5 * (log_histogram_bins[:-1] + log_histogram_bins[1:])
ax[-1].set_xticks(bin_centers)
ax[-1].set_xticklabels(["0", "<10", "<100", "<1000", "<10k", "<100k"])

plt.suptitle(
    "Ensemble distribution of predicted affected population over CH per lead time"
)
plt.tight_layout()
plt.subplots_adjust(hspace=0)
2026-06-08 14:00:32,395 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
2026-06-08 14:00:32,399 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
2026-06-08 14:00:32,403 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
2026-06-08 14:00:32,407 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
2026-06-08 14:00:32,411 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
2026-06-08 14:00:32,422 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
2026-06-08 14:00:32,430 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
2026-06-08 14:00:32,437 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
2026-06-08 14:00:32,445 - climada.engine.impact - INFO - The eai_exp and aai_agg are computed for the selected subset of events WITHOUT modification of the frequencies.
../_images/b0a5fb0d33f3d7dcb6fad6249245a2556e588d98d5b73f0d3fbd09e8cc00c9d4.png

Attributes of an ImpactForecast object that do not make sense if impacts correspond to different ensemble members and lead times (instead of impacts of different events in a probabilistic event set as usual in Impact objects) are filled with NaNs.

print("average annual imact: ", impact_forecast.aai_agg)
print("average annual impact per centroid: ", impact_forecast.eai_exp)
try:
    impact_forecast.calc_freq_curve()
except NotImplementedError as e:
    print(e)
average annual imact:  nan
average annual impact per centroid:  [nan nan nan ... nan nan nan]
2026-06-08 14:00:32,579 - climada.engine.impact_forecast - ERROR - calc_freq_curve is not defined for ImpactForecast
calc_freq_curve is not defined for ImpactForecast

ImpactForecast objects can be written to and read from hdf5 files:

# write hazard forecast object to h5 file
with TemporaryDirectory() as tmpdir:
    impact_forecast.write_hdf5(Path(tmpdir) / "demo_impact_forecast.h5")

    # read hazard forecast from h5 file
    imp_fct = ImpactForecast.from_hdf5(Path(tmpdir) / "demo_impact_forecast.h5")

np.array_equal(impact_forecast.imp_mat.toarray(), imp_fct.imp_mat.toarray())
2026-06-08 14:00:32,604 - climada.engine.impact_forecast - INFO - at_event gives the total impact for one specific combination of member and lead_time.
True