Hazard: RiverFlood

A river flood hazard is generated by the class RiverFlood() that extracts flood data simulated within the Inter-Sectoral Impact Model Intercomparison Project (ISIMIP, https://www.isimip.org/). The method from_nc() generates a data set with flood depth in m and the flooded fraction in each centroid. The data is derived from global hydrological models driven by various climate forcings. A link to the ISIMIP data repository will be provided soon. In this tutorial we show how flood depth and fractions can be translated into socio-economic impacts.

Besides, all other general Hazard Attributes, the class RiverFlood() has further Attributes related to the flooded area and flood volume:

additional Attributes (always calculated)

Name

Data Type

Description

fla_ann_av

float

average flooded area per year

fla_ev_av

float

average flooded area per event

fla_event

1d array(n_events)

average flooded area per year

fla_annual

1d array(n_years)

average flooded area per event

additional Attributes (only calculated if ‘save_centr’ = True in ‘set_flooded_area()’)

Name

Data Type

Description

fla_ev_centr

2d array(n_events x n_centroids)

flooded area in every centroid for every event

fla_ann_centr

2d array(n_years x n_centroids)

flooded area in every centroid for every year

fv_ann_centr

2d array(n_years x n_centroids)

flood volume in every centroid for every year

How is this tutorial structured?

Part 1: Data availability and use

Part 2: Generating a RiverFlood Hazard

Part 3: Calculating Flooded Area

Part 4: Generating ISIMIP Exposure

Part 5: Setting JRC damage functions

Part 6: Deriving flood impact with LitPop exposure

Part 1: Data availability and use

To work with the CLIMADA RiverFlood-Module input data containing spatially explicit flood depth and flooded fraction is required.

The input data can be found at https://files.isimip.org/cama-flood/ *.

On this page, data from the ISIMIP2a and ISIMIP2b simulation rounds can be accessed. The simulations contain the output of the river routing model CaMa-Flood for runoff input data generated by various combinations of global hydrological models (GHMs) and climate forcings.

ISIMIP2a

In the ISIMIP2a simulation round, 12 GHMs were driven by 4 climate reanalysis data sets and covers the time period 1971-2010. The runoff was used as input for CaMa-Flood to derive spatially explicit flood depth (flddph) and flooded fraction (fldfrc) of the maximum flood event of each year on 150 arcsec (~ 5 km) and 300 arcsec (~ 10 km) resolution. Data are provided for different protection standards including ‘0’- no protection, ‘100’- protection against all events smaller than 100 year return period, and’Flopros’- merged layer in the Flopros data base on global protection standards. File naming conventions follow the scheme:

<indicator_resolution_GHM_ClimateForcingDataset_ProtectionStandard.nc>

ISIMIP2b

In the ISIMIP2b simulation round, 6 GHMs were driven by 4 global circulation models (GCMs) and covers the time period 2005-2100 for RCP 2.6, 6.0 and RCP8.5 (only a smaller ensemble). Additionally, historical and preindustrial control runs are provided. Resolution and protection assumptions are the same as under ISIMIP2a. File naming conventions follow the scheme:

<indicator_resolution_GHM_GCM_ProtectionStandard.nc>

For futher information on flood data generation see also:

Willner, S. N. et al. (2018) ‘Adaptation required to preserve future high-end river flood risk at present levels’, Science Advances. American Association for the Advancement of Science, 4(1), p. eaao1914. doi:10.1126/sciadv.aao1914.

Willner, S. N., Otto, C. and Levermann, A. (2018) ‘Global economic response to river floods’, Nature Climate Change. Nature Publishing Group, 8(7), pp. 594–598. doi:10.1038/s41558-018-0173-2.

Sauer, I. et al. (2020) ‘Climate Signals in River Flood Damages Emerge under Sound Regional Disaggregation’. doi:10.21203/rs.3.rs-37259/v1.

*Currently, log-in data are required, please contact inga.sauer@pik-potsdam.de to obtain access.

Part 2: Generating a RiverFlood Hazard

A river flood is generated with the method from_nc(). There are different options for choosing centroids. You can set centroids for: - countries - regions - global hazards - with random coordinates - with random shape files (under development)

Countries or regions can either be set with corresponding ISIMIPNatID centroids (ISINatIDGrid = True) or with Natural Earth Multipolygons (default). It is obligatory to set paths for flood depth and flood fraction, here we present example files from floods for the year 2000.

Setting floods for countries with Natural Earth Multipolygons:

[8]:
import numpy as np
import matplotlib.pyplot as plt
from climada_petals.hazard.river_flood import RiverFlood
from climada.hazard.centroids import Centroids
from climada_petals.util.constants import HAZ_DEMO_FLDDPH, HAZ_DEMO_FLDFRC

years = [2000]
# generating RiverFlood hazard from netCDF file
# uses centroids from Natural Earth Multipolygon for Germany and Switzerland
rf = RiverFlood.from_nc(countries = ['DEU','CHE'], years=years, dph_path=HAZ_DEMO_FLDDPH, frc_path=HAZ_DEMO_FLDFRC)
rf.event_name
2021-12-06 17:09:44,170 - climada.hazard.base - WARNING - The use of Hazard.set_raster is deprecated.Use Hazard.from_raster instead.
2021-12-06 17:09:44,188 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
2021-12-06 17:09:44,203 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
2021-12-06 17:09:44,219 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/fldfrc_2000_DEMO.nc
[8]:
['2000']
[4]:
# Note: Points outside the selected countries are masked in further analysis.
# plot centroids:
rf.centroids.plot()
# get resolution
print('resolution:')
rf.centroids.meta['transform'][0]
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/pyproj/crs/crs.py:1256: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems
  return self._crs.to_proj4(version=version)
resolution:
[4]:
0.04166666666666662
../_images/tutorial_climada_hazard_RiverFlood_9_3.png
[5]:
# plotting intensity (Flood depth in m)
rf.plot_intensity(event=0, smooth = False)
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/cartopy/mpl/geoaxes.py:1797: MatplotlibDeprecationWarning: shading='flat' when X and Y have the same dimensions as C is deprecated since 3.3.  Either specify the corners of the quadrilaterals with X and Y, or pass shading='auto', 'nearest' or 'gouraud', or set rcParams['pcolor.shading'].  This will become an error two minor releases later.
  result = matplotlib.axes.Axes.pcolormesh(self, *args, **kwargs)
[5]:
<GeoAxesSubplot:title={'center':'RF max intensity at each point'}>
../_images/tutorial_climada_hazard_RiverFlood_10_2.png

Setting flood with ISIMIP NatIDGrid:

[7]:
# generating RiverFlood hazard from netCDF file, using the ISIMIP NatIDGrid (according to ISIMIP standards) with a resolution of 150as (aprox 5km)
# setting centroids for a region
rf_isi = RiverFlood.from_nc(countries = ['DEU','CHE'], years=years, dph_path=HAZ_DEMO_FLDDPH, frc_path=HAZ_DEMO_FLDFRC, ISINatIDGrid=True)
rf_isi.centroids.plot()
rf_isi.plot_intensity(event=0, smooth = False)
2021-12-06 16:50:52,281 - climada.hazard.base - WARNING - The use of Hazard.set_raster is deprecated.Use Hazard.from_raster instead.
2021-12-06 16:50:52,303 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
2021-12-06 16:50:52,319 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
2021-12-06 16:50:52,329 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/fldfrc_2000_DEMO.nc
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/pyproj/crs/crs.py:1256: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems
  return self._crs.to_proj4(version=version)
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/cartopy/mpl/geoaxes.py:1797: MatplotlibDeprecationWarning: shading='flat' when X and Y have the same dimensions as C is deprecated since 3.3.  Either specify the corners of the quadrilaterals with X and Y, or pass shading='auto', 'nearest' or 'gouraud', or set rcParams['pcolor.shading'].  This will become an error two minor releases later.
  result = matplotlib.axes.Axes.pcolormesh(self, *args, **kwargs)
[7]:
<GeoAxesSubplot:title={'center':'RF max intensity at each point'}>
../_images/tutorial_climada_hazard_RiverFlood_12_3.png
../_images/tutorial_climada_hazard_RiverFlood_12_4.png

Setting flood with random points as coordinates:

[9]:
lat = np.arange(47, 56, 0.2)
lon = np.arange(5.8, 15, 0.2)
lon, lat = np.meshgrid(lon, lat)
rand_centroids = Centroids.from_lat_lon(lat.flatten(), lon.flatten())
rf_rand = RiverFlood.from_nc(dph_path=HAZ_DEMO_FLDDPH, frc_path=HAZ_DEMO_FLDFRC,
                    centroids=rand_centroids, ISINatIDGrid=False)
rf_rand.centroids.plot()
rf_rand.plot_intensity(event = 0)

2021-12-06 16:51:28,940 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/fldfrc_2000_DEMO.nc
2021-12-06 16:51:28,946 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/pyproj/crs/crs.py:1256: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems
  return self._crs.to_proj4(version=version)
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/cartopy/mpl/geoaxes.py:1797: MatplotlibDeprecationWarning: shading='flat' when X and Y have the same dimensions as C is deprecated since 3.3.  Either specify the corners of the quadrilaterals with X and Y, or pass shading='auto', 'nearest' or 'gouraud', or set rcParams['pcolor.shading'].  This will become an error two minor releases later.
  result = matplotlib.axes.Axes.pcolormesh(self, *args, **kwargs)
[9]:
<GeoAxesSubplot:title={'center':'RF max intensity at each point'}>
../_images/tutorial_climada_hazard_RiverFlood_14_3.png
../_images/tutorial_climada_hazard_RiverFlood_14_4.png
[10]:
# setting random poits using raster
min_lat, max_lat, min_lon, max_lon = 45.6 , 49., 9., 14.
cent = Centroids.from_pnt_bounds((min_lon, min_lat, max_lon, max_lat), res=0.1)
rf_rast = RiverFlood.from_nc(dph_path=HAZ_DEMO_FLDDPH, frc_path=HAZ_DEMO_FLDFRC,
                    centroids=cent, ISINatIDGrid=False)
rf_rast.plot_intensity(event=0)
2021-12-06 16:51:48,944 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/fldfrc_2000_DEMO.nc
2021-12-06 16:51:48,951 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/pyproj/crs/crs.py:1256: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems
  return self._crs.to_proj4(version=version)
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/cartopy/mpl/geoaxes.py:1797: MatplotlibDeprecationWarning: shading='flat' when X and Y have the same dimensions as C is deprecated since 3.3.  Either specify the corners of the quadrilaterals with X and Y, or pass shading='auto', 'nearest' or 'gouraud', or set rcParams['pcolor.shading'].  This will become an error two minor releases later.
  result = matplotlib.axes.Axes.pcolormesh(self, *args, **kwargs)
[10]:
<GeoAxesSubplot:title={'center':'RF max intensity at each point'}>
../_images/tutorial_climada_hazard_RiverFlood_15_3.png

Part 3: Calculating Flooded Area

The fraction indicates the flooded part of a grid cell. It is possible to calculate the flooded area for each grid cell and for the whole area under consideration

As ISIMIP simulations currently provide yearly data with the maximum event, event and yearly flooded area are the same.

[11]:
#setting river flood
rf_DEU = RiverFlood.from_nc(countries = ['DEU'], years=years, dph_path=HAZ_DEMO_FLDDPH, frc_path=HAZ_DEMO_FLDFRC)
rf_DEU.plot_fraction(event=0, smooth = False)
# calculating flooded area
rf_DEU.set_flooded_area()
print("Total flooded area for year " + str(years[0]) + " in Germany:")
print(str(rf_DEU.fla_annual[0]) + " m2")

print("Total flooded area at first event in Germany:")
print(str(rf_DEU.fla_event[0]) + " m2")

2021-12-06 16:51:57,602 - climada.hazard.base - WARNING - The use of Hazard.set_raster is deprecated.Use Hazard.from_raster instead.
2021-12-06 16:51:57,621 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
2021-12-06 16:51:57,634 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
2021-12-06 16:51:57,649 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/fldfrc_2000_DEMO.nc
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/cartopy/mpl/geoaxes.py:1797: MatplotlibDeprecationWarning: shading='flat' when X and Y have the same dimensions as C is deprecated since 3.3.  Either specify the corners of the quadrilaterals with X and Y, or pass shading='auto', 'nearest' or 'gouraud', or set rcParams['pcolor.shading'].  This will become an error two minor releases later.
  result = matplotlib.axes.Axes.pcolormesh(self, *args, **kwargs)
2021-12-06 16:52:01,262 - climada.hazard.centroids.centr - INFO - Convert centroids to GeoSeries of Point shapes.
/Users/evelynm/climada_python/climada/hazard/centroids/centr.py:822: UserWarning: Geometry is in a geographic CRS. Results from 'buffer' are likely incorrect. Use 'GeoSeries.to_crs()' to re-project geometries to a projected CRS before this operation.

  xy_pixels = self.geometry.buffer(res / 2).envelope
Total flooded area for year 2000 in Germany:
2437074832.038374 m2
Total flooded area at first event in Germany:
2437074832.038374 m2
../_images/tutorial_climada_hazard_RiverFlood_17_5.png
[12]:
#calculate flooded area
rf_DEU.set_flooded_area(save_centr = True)
print("affected area in each affected centroid and each event:")
rf_DEU.fla_ev_centr.data
affected area in each affected centroid and each event:
[12]:
array([  584715.81718056,  5053615.42987047,   584715.81718053,
         772660.21477222,  4635961.19139062,  4844788.31063079,
         877073.72577079,    62648.12847878,   793542.93642048,
         710012.09844878,   104512.31458989,  7942935.27615602,
       12980428.9975574 ,  8862643.59587924, 11015597.79960522,
        8151960.31900815,  8026545.7605037 ,  8110155.46617279,
        3428003.77254628,  2550100.30565756,  5748176.99827292,
       10743865.08816197, 10116791.51696033,   167219.69117698,
         167377.62651458,  2426975.61490724,  7866748.53143359,
       12448711.57332617, 10691246.31833798,  3807841.25590698,
         795043.7594347 ,  9958969.25866094, 10858623.66475107,
       10586635.08712338,  9938046.70065339,  9414992.10340605,
        4100752.00183631,    20922.20331432,  4020851.5658597 ,
        8628077.07461108, 12271973.68427333, 11036399.71240725,
        6135986.84213667,  5863741.76936025,  5319251.23373377,
        7936993.24829797, 12104437.71477011, 12041612.40883577,
        9444812.29294032,   628258.03278596,  1907510.5122277 ,
        3731174.16397122,  9034472.37433647, 13017186.72931901,
        8300814.5901365 ,  2242896.97806137,  4003675.58205206,
        8636201.25119162, 10795251.36876784,  7839658.2240182 ,
        7986389.9370346 ,   859427.78209491,  1573601.99913756,
        6755998.31492206, 10910307.66299093, 12127226.06553792,
       12630778.64273264, 12295077.44567993,  8329600.41176178,
         524533.99971255,   629440.8191955 ,   314720.40959772,
         881217.10779283,   251776.3374484 ,    84004.22287586,
        5649284.31387142,  1722086.62090805,  2142107.81474468,
        7413373.17151569, 11676586.99808064, 12075607.56251821,
        9051455.12336371, 11823594.53327545,  7224363.59517164,
          42041.479908  ,  4666604.29120064,  4099044.47915377,
        3951899.15417685,  8975856.20966038, 12255091.63636668,
        8198088.95830755,  9333208.58240128,  9627499.23235424,
         588580.71259418,    84161.65049093,  2566930.35528228,
        4060799.89490784,   841616.56614496,  1283465.17764114,
          63121.24246087,   778495.28694269,   862656.94355715,
         252484.96984349,    68119.64171433,   340598.1953547 ,
         476837.45234946,  1977089.19581193,  9340041.72726154,
       11998885.92035427,  5272237.71440184,    45450.32455399,
         704480.06861677,  3431499.68736199,  2727019.56583407,
        5022260.92439447,  2727019.56583394,  3317873.85448186,
        4526009.59398457, 11417371.18678022,  2001451.49577753,
        2092426.62153581,    90975.06618462,   545850.43682358,
        1182675.89349668,  3047664.70063641,  4662472.38853107,
       10484877.04798255, 10507620.35283201, 11758527.53769297,
        7914831.06916955,  5344785.09035629,  3411564.99185211,
         159206.36416826,  3824076.16877868,  5963738.09206057,
          45524.71676818,   751157.84158056,  2412809.97877638,
        6669371.29305688,  9833339.32540515, 11290130.36798235,
       10334111.21979271, 11358416.90487568,  8035113.05446658,
        4666283.71550905,  1616127.4966119 ,   637346.02812978,
         159466.56281482,   888456.57183127,  2505903.20571988,
        4715654.29155094,  4077787.82812764,    91123.75112706,
        1412418.21871583,  4419502.17166185, 11322126.22505713,
        7289900.62057462,  3576607.18698346,  2323655.79628766,
         797332.8671151 ,   432837.83608639,  1117175.77323329,
        2439547.17032033,  4559900.9592355 ,  2941136.28432971,
        1117175.77323329,    68398.51916611,    45599.00945964,
        2553544.50319797,  3807517.28822123,  7843030.09154589,
        3579522.19779226,   569987.61990444,   136908.36527585,
         114090.30439655,    22818.05921908,  1300629.43027509,
        4152887.05347047,   912722.43517239,   136908.36527586,
         410725.06926384,  1460355.79002089,  1551628.08666564,
        2624076.8815839 ,  7507142.09570224,  2760985.22029591,
         250998.65639055,   159726.41287331,   319452.82574662,
        2329333.4694195 ,  3950732.208999  ,  4750013.19093546,
        1164666.73470981,   137019.61741205,   936300.67910442,
        1210339.86075792,  1027647.03754179,   137019.61741205,
         685098.06047493,  1415869.35333926,  4932706.12049274,
       10459163.73034007, 12423111.55332744,  4658666.93883924,
         411058.82565086,  1278849.65617129,   913464.11608036,
         274261.58936941,  4731012.2436776 ,   137130.7946847 ,
        6925104.95863252,  4822432.61382616,  4548171.0776705 ,
        4593881.26274455,  3771096.65427761,  4868142.79890019,
        4616736.35528191,  9530590.19067662,  8776370.85982062,
        8959211.60011772,   297116.70851342,  1898512.78921824,
        7502556.71861722,   343104.72928449,   548967.5881579 ,
        7273820.62297752,   160115.53323185,  3133689.83594515,
        6427495.36489242,  8874975.84878219,  7479682.68299801,
        4048635.70969465, 10430383.85545952,  7777040.88549384,
        3682657.53061663,    45784.30480911, 10805096.61452431,
        3845881.56399029,  5471224.64954781,    91568.60961822,
         114460.77035089,  5448332.53545232,  1304852.74202514,
        3731420.9935142 , 11995488.25307395, 10667743.07715331,
        4715783.60520665,  2541028.92856497,  4051911.01724683,
          68676.46221054,  4375933.02392118,  7514691.61028689,
         618587.43214996,  8568581.198868  ,  2084868.78441238,
          45821.28892426,  2703456.21656234,  6254606.13986513,
        3642792.77453402,    91642.57784853,  3001294.52289039,
        4673771.75699323,    91642.57784852,  6094231.47026817,
         824783.20730461,    68731.9383873 ,   458212.92258201,
        2909651.88503098,  4765414.18148062,  5819303.77006197,
          68787.37702238,   275149.50808951,  2476345.5194196 ,
        8850642.43013321,  1467463.93637202,  5525919.03388077,
       10547397.32962378,  1467463.93637202,   710802.88255137,
        1903117.31083391,  4838045.18357794,  7612469.24333563,
        1329889.20902028,   321007.73274474,   802519.38524791,
        6213792.8841836 ,   550299.016179  ,  1444534.85073734,
        4539966.64323985,    45858.24801163,  2292912.40725453,
        3714518.2791294 ,  4677541.58413571,  5021478.29574301,
        4654612.49850087,  5663493.97477634,  1811400.75475137,
         389795.12311364,   343936.87176537,  4058454.9907367 ,
        4933732.36289857,  3878142.80501333,  2616025.48722424,
        3786352.59451557,  8743032.51003858,   688427.75417171,
        1032641.57782859,  6677749.3543814 ,  6471220.95332924,
         206528.32090862,  3740457.48926669,  7526810.08378189,
        3327400.68716238,  6035216.59860048,  4222356.94924433,
        1468645.82569923,  4061723.65344085,  3373296.00612736,
        1973492.94515943,   114737.96348104,   160633.1355162 ,
         481899.40654861,  5507422.03337424,  7848076.67538843,
        4773099.06709529,  4658361.30397309,  1858754.96832113,
        2340654.42829876,   688427.75417175,  3465086.43034104,
        4314147.15974187,  2180021.13249548,  7389124.76803563,
        5094365.23126939,  2273638.55111259,  7027610.14485272,
        2296604.5581947 ,  8727097.23558467,   206694.41130696,
         620083.26065688,  6384560.66322575,   987540.01563459,
        1377962.77769441,  4409481.0597325 ,  5833375.63770317,
         160762.31693474,  2664061.31317241,  2916687.81885159,
         183728.36412085,  2503299.04970953,  3031518.06815012,
        1056438.14382492,  2939653.82593355,  6476425.11933014,
        4478379.08097884,  1377962.77769448,   826777.64522785,
        2687027.32025452,  5029564.10650116,    45932.09103021,
         160891.41056458, 10848678.76457667,   436705.2801817 ,
        2206510.95693681,  9354686.58283404,  1333100.28954366,
        3585580.03744784,  1103255.47846846,   436705.2801817 ,
        2183526.29387852,  3585580.03744784,   344767.32364266,
        4527944.15886665,  2987983.50724635,  2022634.93682888,
        3884378.30254839,    22984.48746242,   137906.93480856,
        3930347.62866516,  4665851.28097746,  1195193.38149251,
         551627.73923423,   344767.32364266,  3792440.50655434,
        3861393.85355002,  3585580.03744784,  1103255.47846846,
        1356084.84557202,  4117522.08266347,  6026764.41462543,
        4485568.99342198,  2484315.14800192,  3841487.11382578,
         368046.66974845,   782099.20334174,  7199913.35353265,
         713090.46113229,  2001253.84542028,  4761603.96225989,
         230029.18532961,  4025510.56920494,    23002.91685928,
        2553323.78309566,  2070262.48051402,  4462566.11505762,
         644081.66536506,   552070.04479106,   966896.00520992,
        6353888.18738154,    46042.66739333,  7597040.27065209,
         368341.33914669,  1473365.35658675,  1427322.76624445,
        6468995.09204306,  4028733.5158435 ,   184170.66957333,
          23021.33369667,  1588472.04684531,  3107880.20817732,
        2394218.77145437,   713661.38312219,  1151066.68818345,
        4028733.5158435 ,  3430178.76937919,   851789.36855201,
          69064.00611507,   138238.43785095,  2027496.98119351,
        2741728.95032279,    92158.95186185,  9123736.41537035,
          46079.47593093,  1105907.50280763,  1797099.55795347,
       10367882.75835618,  5183941.37917809,  4792265.63092521,
          46079.47593093,  4630987.30591261,    69119.21892548,
        1727980.35243889,  1105907.50280763,  4907464.12797093,
        7142318.83299559,    46079.47593093,  2119655.99340428,
        4354510.48385451,    23039.73796546,   990708.79118756,
        6220729.56918368,   484220.71775385,   322813.81183591,
        6686857.57690364,   184465.03725218,  1867708.59780709,
        1567952.8770407 ,  2582510.49468729,  8854322.43234125,
       11690472.35828338,   622569.53260236,  1775476.01878391,
        6940497.00815814,   299755.69392329,  4196579.87598521,
         484220.71775387, 10883437.88237945,   391988.21926018,
        1867708.59780728,  4035172.80900787,    46116.25931304,
        1292284.47717146,  7915242.58386276,    69229.53131946,
        1615355.7039227 ,  3023022.74672555,   207688.58052607,
        6230657.68442802,  9969052.29508483,   138459.06263891,
       11261336.34242281,  2099962.43659108,    23076.50876078,
        5561438.76078262,  1476896.56069001,  2076885.75153151,
       11953632.16604402,   369224.1401725 ,  1823044.15012581,
        8007548.46443489,   276918.12527781,   346423.14080041,
        5658244.43590979,  3210187.66745817,  3210187.66745817,
         323328.25040788,  1662831.03282448,  4087792.9646555 ,
         484992.37561182,  3810654.49503272,  4572785.39403927,
       14388107.73464195,   577371.88341003,  1570451.57879818,
        4318741.76103709,  3048523.59602618,   415707.7582061 ,
        5450390.9063244 ,  6674418.96791818,  3325662.06564879,
        2632815.89159201, 11131730.17885479, 14226443.44812288,
        8175585.9293106 ,  2218870.16200874,  3258965.45627431,
         901415.93604765,   184905.83338038,  3582550.55369736,
         254245.52594315,   138679.38512554,  2704247.80814296,
         346698.44936018, 13428786.49400053,  1479246.66704298,
        3490097.79172438,  1433020.28605656,  1317454.11833161,
        3883022.46062701,  8806140.64603687,  7350007.4923756 ,
        6818402.8194786 ,  9846235.50978571,  8783027.02502618,
        4691985.84996017, 12481144.39222511, 13729258.40093017,
       12527370.77321223, 10239161.03972273,   138679.38512554,
         670815.55181327,  3053367.36157378,  3539130.39008407,
         300710.42442124,   208184.13584106,  2775788.54969063,
        1734867.78969926,    69394.71643513,  3770446.13846321,
         532026.14587166,  3284683.10995292,  1919920.36685962,
        2752657.01793862,  2382551.86361778,   693947.13742262,
        1387894.27484531,  1434157.4460641 ,  3354077.70520894,
       12699232.90566386,  1272236.40065568,  2081841.3045532 ,
        2428814.92712191,  7656550.36654593,  6800681.96828465,
         717078.72303205,    23149.8991283 ,  3403035.29818859,
        3310435.64103781,   347248.49871514,  4375330.94367087,
         902846.06431938,  4954078.74696315,   810246.51496873,
         578747.47989191,  4143831.90859385,  4467930.81642188,
        8171914.94645002,  2500189.23386921,  5486526.39827995,
         925996.03250709,    46299.7982566 ,   208349.09383909,
        3266718.48119637,  2919195.1880162 ,   417027.87629623,
         671878.24814075,  4957998.37854305,  4934829.77061048,
        1876625.52424711,  3336223.01036971,   810887.5762013 ,
        8386894.52601908,    46336.43032498,  4401960.85053028,
        3730082.87210306,    46373.03711109,   231865.20242587,
        2921501.41560279,  1159325.93115142,  3918521.55659626,
          69559.56072776,   347797.79014247,  1669429.34949569,
        2156346.28808637,  6469038.64831783, 12798958.41811345,
        1808548.44395856,  3895335.07684266,    69559.56072776,
        1669429.34949561,  1738988.89672704,  2133159.80833278,
          23186.51855555,    23204.80929843,   116024.05493402,
         835373.14149704,  6195684.08774647,    69614.43296041,
        3202263.7270813 ,  5105058.14695738,    46409.61859686,
        3619950.24380173,  3411107.09349742,   139228.86592083,
        3643155.01426779,  7843225.38416266,   348072.15129509,
          46409.61859686,  9142695.12359154,   696144.30259014,
         348346.32255877,  1254046.78283975,  6618580.07454613,
        1045038.91360584,   441238.6788458 ,   464461.78143516,
        9266012.36390282,  2090077.82721168,  1394481.2151967 ,
        2788962.43039341,  5229304.71932662,  7111854.17585799,
         418344.35373641,  1301515.74312973,  2672755.64442255,
        2091721.71456896,  3114341.2579497 ,  3230548.04392058,
        1882549.67298334,  5024075.05294284,  4748680.79119445,
         116480.37319711,  4286477.62517271,   163200.01148239,
         536228.64017519])

Part 4: Generating ISIMIP Exposure

The exposed assets are calculated by means of national GDP converted to total national wealth as a proxy for asset distribution, downscaled by means of data from spatially explicit GDP distribution. Data for past (1971-2010) and future (2005-2100) periods can be accessed at ISIMIP, https://www.isimip.org/ .

More information on spatially explicit GDP time series:

Geiger, T. (2018) ‘Continuous national gross domestic product (GDP) time series for 195 countries: Past observations (1850-2005) harmonized with future projections according to the Shared Socio-economic Pathways (2006-2100)’, Earth System Science Data. Copernicus GmbH, pp. 847–856. doi: 10.5194/essd-10-847-2018.

Murakami, D. and Yamagata, Y. (2019) ‘Estimation of gridded population and GDP scenarios with spatially explicit statistical downscaling’, Sustainability (Switzerland). Multidisciplinary Digital Publishing Institute, 11(7), p. 2106. doi: 10.3390/su11072106.

[1]:
# set exposure for damage calculation
from climada_petals.entity.exposures.gdp_asset import GDP2Asset
from climada_petals.util.constants import DEMO_GDP2ASSET
gdpa = GDP2Asset()
gdpa.set_countries(countries = ['CHE'], ref_year = 2000, path=DEMO_GDP2ASSET)
gdpa
[1]:
<climada_petals.entity.exposures.gdp_asset.GDP2Asset at 0x7fe7ba7278e0>
[2]:
from matplotlib import colors
norm=colors.LogNorm(vmin=1.0e2, vmax=1.0e10)
gdpa.plot_scatter()
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/pyproj/crs/crs.py:1256: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems
  return self._crs.to_proj4(version=version)
[2]:
<GeoAxesSubplot:title={'center':'CHE GDP2Asset  GDP2Asset 2018'}>
../_images/tutorial_climada_hazard_RiverFlood_21_2.png

Part 5: Setting JRC damage functions

In CLIMADA we currently calculate damage by translating flood-depth into a damage factors. Damage assessments implemented in CLIMADA base on the residential damage functions basing on an empirical estimate published in the JRC report. Individual damage functions are available for six continents:

RF1: Africa RF2: Asia RF3: Europe RF4: North America RF5: Oceania RF6: South America

For further information on depth-damage functions, see also:

Huizinga, J., Moel, H. de and Szewczyk, W. (2017) Global flood depth-damage functions : Methodology and the Database with Guidelines, Joint Research Centre (JRC). doi: 10.2760/16510.

[3]:
# import impact function set for RiverFlood using JRC damage functions () for 6 regions
from climada_petals.entity.impact_funcs.river_flood import ImpfRiverFlood,flood_imp_func_set
impf_set = flood_imp_func_set()
impf_AFR = impf_set.get_func(fun_id=1)
impf_AFR[0].plot()
impf_EUR = impf_set.get_func(fun_id=3)
impf_EUR[0].plot()
impf_OCE = impf_set.get_func(fun_id=6)
impf_OCE[0].plot()
[3]:
<AxesSubplot:title={'center':'RF 6: Flood South America JRC Residential noPAA'}, xlabel='Intensity (m)', ylabel='Impact (%)'>
../_images/tutorial_climada_hazard_RiverFlood_23_1.png
../_images/tutorial_climada_hazard_RiverFlood_23_2.png
../_images/tutorial_climada_hazard_RiverFlood_23_3.png

The plots illustrate how flood-depth is translated into a damage factor (0%-100%). The damage factor is then multiplied with the exposed asset in each grid cell to derive a local damage.

Linking exposures to the correct impact function

If the ISIMIP exposure presented above is used, the correct impact function ID is automatically provided in the GeoDataFrame:

[4]:
gdpa
[4]:
<climada_petals.entity.exposures.gdp_asset.GDP2Asset at 0x7fe7ba7278e0>

The column ‘impf_RF’ indicates the ID of the impact function (in this case 3 for Europe). If other Exposure data is used the impact function needs to be set manually.

Part 6: Deriving flood impact with LitPop exposure

[5]:
from climada.entity import LitPop
lp_exp = LitPop.from_countries(['DEU'], fin_mode='pc')
lp_exp
2021-12-06 17:08:09,150 - climada.entity.exposures.litpop.litpop - INFO -
 LitPop: Init Exposure for country: DEU (276)...

2021-12-06 17:08:09,152 - climada.entity.exposures.litpop.gpw_population - WARNING - Reference year: 2018. Using nearest available year for GPW data: 2020
2021-12-06 17:08:09,153 - climada.entity.exposures.litpop.gpw_population - INFO - GPW Version v4.11
2021-12-06 17:08:13,072 - climada.util.finance - INFO - GDP DEU 2014: 3.884e+12.
2021-12-06 17:08:13,547 - climada.util.finance - INFO - GDP DEU 2018: 3.975e+12.
2021-12-06 17:08:13,900 - climada.entity.exposures.base - INFO - Hazard type not set in impf_
2021-12-06 17:08:13,900 - climada.entity.exposures.base - INFO - category_id not set.
2021-12-06 17:08:13,901 - climada.entity.exposures.base - INFO - cover not set.
2021-12-06 17:08:13,901 - climada.entity.exposures.base - INFO - deductible not set.
2021-12-06 17:08:13,902 - climada.entity.exposures.base - INFO - centr_ not set.
[5]:
<climada.entity.exposures.litpop.litpop.LitPop at 0x7fe6f54ecc70>
[6]:
# In the LitPop exposure the damage function for river floods needs
# to be specified manually.
import pandas as pd
from climada_petals.util.constants import RIVER_FLOOD_REGIONS_CSV

info = pd.read_csv(RIVER_FLOOD_REGIONS_CSV)
lp_exp.gdf['impf_RF'] = info.loc[info['ISO']=='DEU','if_RF'].values[0]
lp_exp
lp_exp.plot_hexbin(pop_name=True)
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/pyproj/crs/crs.py:1256: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems
  return self._crs.to_proj4(version=version)
[6]:
<GeoAxesSubplot:title={'center':"LitPop Exposure for ['DEU'] at 30 as, year: 2018, financial mode: pc,\nexp: (1, 1), admin1_calc: False"}>
../_images/tutorial_climada_hazard_RiverFlood_29_2.png
[9]:
from climada.engine import Impact

rf = RiverFlood.from_nc(countries = ['DEU'], years=years, dph_path=HAZ_DEMO_FLDDPH, frc_path=HAZ_DEMO_FLDFRC)
imp=Impact()
imp.calc(lp_exp, impf_set,rf,save_mat=True)
rf.plot_intensity(0)
imp.plot_scatter_eai_exposure()
2021-12-06 17:09:50,333 - climada.hazard.base - WARNING - The use of Hazard.set_raster is deprecated.Use Hazard.from_raster instead.
2021-12-06 17:09:50,352 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
2021-12-06 17:09:50,366 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/flddph_2000_DEMO.nc
2021-12-06 17:09:50,381 - climada.util.coordinates - INFO - Reading /Users/evelynm/climada_petals/doc/tutorial/data/demo/fldfrc_2000_DEMO.nc
2021-12-06 17:09:50,407 - climada.entity.exposures.base - INFO - Matching 661396 exposures with 41548 centroids.
2021-12-06 17:09:50,449 - climada.engine.impact - INFO - Calculating damage for 656383 assets (>0) and 1 events.
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/pyproj/crs/crs.py:1256: UserWarning: You will likely lose important projection information when converting to a PROJ string from another format. See: https://proj.org/faq.html#what-is-the-best-format-for-describing-coordinate-reference-systems
  return self._crs.to_proj4(version=version)
/Users/evelynm/opt/anaconda3/envs/climada_env/lib/python3.8/site-packages/cartopy/mpl/geoaxes.py:1797: MatplotlibDeprecationWarning: shading='flat' when X and Y have the same dimensions as C is deprecated since 3.3.  Either specify the corners of the quadrilaterals with X and Y, or pass shading='auto', 'nearest' or 'gouraud', or set rcParams['pcolor.shading'].  This will become an error two minor releases later.
  result = matplotlib.axes.Axes.pcolormesh(self, *args, **kwargs)
[9]:
<GeoAxesSubplot:title={'center':'Expected annual impact'}>
../_images/tutorial_climada_hazard_RiverFlood_30_3.png
../_images/tutorial_climada_hazard_RiverFlood_30_4.png
[ ]: