END-TO-END IMPACT CALCULATION

Goal of this tutorial

The goal of this tutorial is to show a full end-to-end impact computation. Note that this tutorial examplifies the work flow, but does not explore all possible features.

What is an Impact?

The impact is the combined effect of hazard events on a set of exposures mediated by a set of impact functions. By computing the impact for each event (historical and synthetic) and for each exposure value at each geographical location, the Impact class provides different risk measures, such as the expected annual impact per exposure, the probable maximum impact for different return periods, and the total average annual impact.

Impact class data structure

The impact class does not require any attributes to be defined by the user. All attributes are set by the method impact.calc(). This method requires three attributes: an Exposure, a Hazard, and an ImpactFuncSet. After calling impact.calc(Exposure, ImpactFuncSet, Hazard, save_mat=False), the Impact object has the following attributes:

Attributes from input

Data Type

Description

tag

(dict)

dictionary storing the tags of the inputs (Exposure.tag, ImpactFuncSet.tag Hazard.tag)

even_id

list(int)

id (>0) of each hazard event (Hazard.event_id)

event_name

(list(str))

name of each event (Hazard.event_name)

date

np.array

date of events (Hazard.date)

coord_exp

np.array

exposures coordinates [lat, lon] (in degrees) (Exposure.gdf.latidues, Exposure.gdf.longitude)

frequency

np.array

annual frequency of events (Hazard.frequency)

unit

str

value unit used (Exposure.value_unit)

csr

str

unit system for Exposure and Hazard geographical data (Exposure.csr)

Computed attributes

Data Type

Description

at_event

np.array

impact for each hazard event summed over all locations

eai_exp

np.array

expected annual impact for each locations, summed over all events weigthed by frequency

aai_agg

float

total annual average aggregated impact value (summed over events and locations)

impt_mat

sparse.csr_matrix

matrix (num_events X num_exp) with impact values (only filled if save_mat is True).

tot_value

float

total exposure value affected (sum of value all exposures locations affected by at least one hazard event)

All other methods compute values from the attributes set by Impact.cacl(). For example, one can compute the frequency exceedance curve, plot impact data, or compute traditional risk transfer over impact.

How do I compute an impact in CLIMADA?

In CLIMADA, impacts are computed using the Impact class. To computation of the impact requires an Exposure , an ImpactFuncSet, and a Hazard object. For details about how to define Exposures , Hazard, Impact Functions see the respective tutorials.

The steps of an impact caculations are typically:

  • Set exposure

  • Set hazard and hazard centroids

  • Set impact functions in impact function set

  • Compute impact

  • Visualize, save, use impact output

Hints: Before computing the impact of a given Exposure and Hazard, it is important to correctly match the Exposures' coordinates with the Hazard Centroids. Try to have similar resolutions in Exposures and Hazard. By the impact calculation the nearest neighbor for each Exposure to the Hazard's Centroids is searched.

Hint: Set first the Exposures and use its coordinates information to set a matching Hazard.

Hint: The configurable parameter max_matrix_size defined in the configuration file (located at /climada/conf/defaults.conf) controls the maximum matrix size contained in a chunk. You can decrease its value if you are having memory issues when using the Impact.calc() method. A high value will make the computation fast, but increase the memory use.

Structure of the tutorial

We begin with one very detailled example, and later present in quick and dirty examples.

Part1: Detailed impact calculation with Litpop and TropCyclone

Part2: Quick examples: raster and point exposures/hazards

Part3: Visualization methods

Detailed Impact calculation - LitPop + TropCyclone

We present a detailed example for the hazard Tropical Cyclones and the exposures from LitPop .

Define the exposure

Reminder: The exposures must be defined according to your problem either using CLIMADA exposures such as BlackMarble, LitPop,OSM, extracted from external sources (imported via csv, excel, api, …) or directly user defined. As a reminder, exposures are geopandas dataframes with at least columns ‘latitude’, ‘longitude’ and ‘value’ of exposures. For impact calculations, for each exposure value the corresponding impact function to use (defined by the column if_) and the associated hazard centroids must be defined. This is done after defining the impact function(s) and the hazard(s). See tutorials on Exposures , Hazard, ImpactFuncSet for more details.

Exposures are either defined as a series of (latitude/longitude) points or as a raster of (latitude/longitude) points. Fundamentally, this changes nothing for the impact computations. Note that for larger number of points, consider using a raster which might be more efficient (computationally). For a low number of points, avoid using a raster if this adds a lot of exposures values equal to 0.

We shall here use a raster example.

[1]:
# Exposure from the module Litpop
# Note that the file gpw_v4_population_count_rev11_2015_30_sec.tif must be downloaded (do not forget to unzip) if
# you want to execute this cell on your computer.

%matplotlib inline
import numpy as np
from climada.entity import LitPop

# Cuba with resolution 10km and financial_mode = income group.
exp_lp = LitPop()
exp_lp.set_country(countries=['CUB'], res_km=10, fin_mode='income_group')
exp_lp.check()
2020-10-20 17:08:42,434 - climada - DEBUG - Loading default config file: /Users/ckropf/Documents/Climada/climada_python/climada/conf/defaults.conf
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/pandas_datareader/compat/__init__.py:7: FutureWarning: pandas.util.testing is deprecated. Use the functions in the public API at pandas.testing instead.
  from pandas.util.testing import assert_frame_equal
2020-10-20 17:08:47,579 - climada.entity.exposures.litpop - INFO - Generating LitPop data at a resolution of 300.0 arcsec.
2020-10-20 17:08:50,282 - climada.entity.exposures.gpw_import - INFO - Reference year: 2016. Using nearest available year for GWP population data: 2015
2020-10-20 17:08:50,284 - climada.entity.exposures.gpw_import - INFO - GPW Version v4.11
2020-10-20 17:09:05,219 - climada.util.finance - INFO - GDP CUB 2016: 9.137e+10.
2020-10-20 17:09:05,386 - climada.util.finance - INFO - Income group CUB 2016: 3.
2020-10-20 17:09:05,566 - climada.entity.exposures.litpop - INFO - Creating the LitPop exposure took 19 s
2020-10-20 17:09:05,567 - climada.entity.exposures.base - INFO - Hazard type not set in if_
2020-10-20 17:09:05,567 - climada.entity.exposures.base - INFO - centr_ not set.
2020-10-20 17:09:05,568 - climada.entity.exposures.base - INFO - deductible not set.
2020-10-20 17:09:05,569 - climada.entity.exposures.base - INFO - cover not set.
2020-10-20 17:09:05,569 - climada.entity.exposures.base - INFO - category_id not set.
2020-10-20 17:09:05,570 - climada.entity.exposures.base - INFO - geometry not set.
2020-10-20 17:09:05,586 - climada.entity.exposures.base - INFO - Hazard type not set in if_
2020-10-20 17:09:05,588 - climada.entity.exposures.base - INFO - centr_ not set.
2020-10-20 17:09:05,588 - climada.entity.exposures.base - INFO - deductible not set.
2020-10-20 17:09:05,589 - climada.entity.exposures.base - INFO - cover not set.
2020-10-20 17:09:05,590 - climada.entity.exposures.base - INFO - category_id not set.
2020-10-20 17:09:05,590 - climada.entity.exposures.base - INFO - geometry not set.
[2]:
exp_lp.gdf.head()
[2]:
value latitude longitude region_id if_
0 9.818338e+05 21.875000 -84.875000 192 1
1 1.110009e+06 21.875000 -84.791667 192 1
2 0.000000e+00 21.958333 -84.708333 192 1
3 1.026952e+06 21.958333 -84.625000 192 1
4 1.109359e+06 21.958333 -84.541667 192 1
[3]:
# not needed for impact calculations
# visualize the define exposure
exp_lp.plot_raster()
print('\n Raster properties exposures:', exp_lp.meta)
2020-10-20 17:09:05,637 - climada.util.coordinates - INFO - Raster from resolution 0.0833333333333286 to 0.0833333333333286.
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/cartopy/mpl/gridliner.py:307: UserWarning: The .xlabels_top attribute is deprecated. Please use .top_labels to toggle visibility instead.
  warnings.warn('The .xlabels_top attribute is deprecated. Please '
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/cartopy/mpl/gridliner.py:343: UserWarning: The .ylabels_right attribute is deprecated. Please use .right_labels to toggle visibility instead.
  warnings.warn('The .ylabels_right attribute is deprecated. Please '
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/cartopy/mpl/feature_artist.py:215: MatplotlibDeprecationWarning: Using a string of single character colors as a color sequence is deprecated. Use an explicit list instead.
  **style)

 Raster properties exposures: {'width': 129, 'height': 41, 'crs': {'init': 'epsg:4326', 'no_defs': True}, 'transform': Affine(0.0833333333333286, 0.0, -84.91666666666669,
       0.0, 0.0833333333333286, 19.833333333333336)}
../_images/tutorial_climada_engine_Impact_23_3.png

Define the hazard

Let us define a tropical cyclone hazard using the TropCyclone and TCTracks modules.

[4]:
from climada.hazard import TCTracks, TropCyclone, Centroids

# Load histrocial tropical cyclone tracks from ibtracs over the North Atlantic basin between 2010-2012
ibtracks_na = TCTracks()
ibtracks_na.read_ibtracs_netcdf(provider='usa', basin='NA', year_range=(2010, 2012), correct_pres=True)
print('num tracks hist:', ibtracks_na.size)
#Add randomly generated tracks using the calc_random_walk method (2 per historical track)
ibtracks_na.calc_random_walk(ens_size=1)
print('num tracks hist+syn:', ibtracks_na.size)
ibtracks_na.equal_timestep(0.5)  # Interpolation to make the track smooth c.f.
2020-10-20 17:09:07,502 - climada.hazard.tc_tracks - WARNING - `correct_pres` is deprecated. Use `estimate_missing` instead.
2020-10-20 17:09:08,050 - climada.hazard.tc_tracks - INFO - Progress: 10%
2020-10-20 17:09:08,175 - climada.hazard.tc_tracks - INFO - Progress: 20%
2020-10-20 17:09:08,310 - climada.hazard.tc_tracks - INFO - Progress: 30%
2020-10-20 17:09:08,441 - climada.hazard.tc_tracks - INFO - Progress: 40%
2020-10-20 17:09:08,568 - climada.hazard.tc_tracks - INFO - Progress: 50%
2020-10-20 17:09:08,692 - climada.hazard.tc_tracks - INFO - Progress: 60%
2020-10-20 17:09:08,846 - climada.hazard.tc_tracks - INFO - Progress: 70%
2020-10-20 17:09:08,981 - climada.hazard.tc_tracks - INFO - Progress: 80%
2020-10-20 17:09:09,105 - climada.hazard.tc_tracks - INFO - Progress: 90%
num tracks hist: 60
2020-10-20 17:09:09,256 - climada.hazard.tc_tracks_synth - INFO - Computing 60 synthetic tracks.
/Users/ckropf/Documents/Climada/climada_python/climada/hazard/tc_tracks_synth.py:128: NumbaWarning:
Compilation is falling back to object mode WITH looplifting enabled because Function "_one_rnd_walk" failed type inference due to: non-precise type pyobject
During: typing of argument at /Users/ckropf/Documents/Climada/climada_python/climada/hazard/tc_tracks_synth.py (138)

File "../../climada/hazard/tc_tracks_synth.py", line 138:
def _one_rnd_walk(track, ens_size, ens_amp0, ens_amp, max_angle, rnd_vec):
    <source elided>
    """
    ens_track = list()
    ^

  @jit(parallel=True)
/Users/ckropf/Documents/Climada/climada_python/climada/hazard/tc_tracks_synth.py:128: NumbaWarning:
Compilation is falling back to object mode WITHOUT looplifting enabled because Function "_one_rnd_walk" failed type inference due to: cannot determine Numba type of <class 'numba.core.dispatcher.LiftedLoop'>

File "../../climada/hazard/tc_tracks_synth.py", line 147:
def _one_rnd_walk(track, ens_size, ens_amp0, ens_amp, max_angle, rnd_vec):
    <source elided>
    ens_track.append(track)
    for i_ens in range(ens_size):
    ^

  @jit(parallel=True)
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/numba/core/object_mode_passes.py:178: NumbaWarning: Function "_one_rnd_walk" was compiled in object mode without forceobj=True, but has lifted loops.

File "../../climada/hazard/tc_tracks_synth.py", line 138:
def _one_rnd_walk(track, ens_size, ens_amp0, ens_amp, max_angle, rnd_vec):
    <source elided>
    """
    ens_track = list()
    ^

  state.func_ir.loc))
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/numba/core/object_mode_passes.py:188: NumbaDeprecationWarning:
Fall-back from the nopython compilation path to the object mode compilation path has been detected, this is deprecated behaviour.

For more information visit https://numba.pydata.org/numba-doc/latest/reference/deprecation.html#deprecation-of-object-mode-fall-back-behaviour-when-using-jit

File "../../climada/hazard/tc_tracks_synth.py", line 138:
def _one_rnd_walk(track, ens_size, ens_amp0, ens_amp, max_angle, rnd_vec):
    <source elided>
    """
    ens_track = list()
    ^

  state.func_ir.loc))
/Users/ckropf/Documents/Climada/climada_python/climada/hazard/tc_tracks_synth.py:128: NumbaWarning:
Compilation is falling back to object mode WITHOUT looplifting enabled because Function "_one_rnd_walk" failed type inference due to: non-precise type pyobject
During: typing of argument at /Users/ckropf/Documents/Climada/climada_python/climada/hazard/tc_tracks_synth.py (147)

File "../../climada/hazard/tc_tracks_synth.py", line 147:
def _one_rnd_walk(track, ens_size, ens_amp0, ens_amp, max_angle, rnd_vec):
    <source elided>
    ens_track.append(track)
    for i_ens in range(ens_size):
    ^

  @jit(parallel=True)
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/numba/core/object_mode_passes.py:178: NumbaWarning: Function "_one_rnd_walk" was compiled in object mode without forceobj=True.

File "../../climada/hazard/tc_tracks_synth.py", line 147:
def _one_rnd_walk(track, ens_size, ens_amp0, ens_amp, max_angle, rnd_vec):
    <source elided>
    ens_track.append(track)
    for i_ens in range(ens_size):
    ^

  state.func_ir.loc))
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/numba/core/object_mode_passes.py:188: NumbaDeprecationWarning:
Fall-back from the nopython compilation path to the object mode compilation path has been detected, this is deprecated behaviour.

For more information visit https://numba.pydata.org/numba-doc/latest/reference/deprecation.html#deprecation-of-object-mode-fall-back-behaviour-when-using-jit

File "../../climada/hazard/tc_tracks_synth.py", line 147:
def _one_rnd_walk(track, ens_size, ens_amp0, ens_amp, max_angle, rnd_vec):
    <source elided>
    ens_track.append(track)
    for i_ens in range(ens_size):
    ^

  state.func_ir.loc))
num tracks hist+syn: 120
2020-10-20 17:09:17,313 - climada.hazard.tc_tracks - INFO - Interpolating 120 tracks to 0.5h time steps.
[5]:
# not needed for calculations
# visualize tracks
ax = ibtracks_na.plot()
ax.get_legend()._loc = 2
../_images/tutorial_climada_engine_Impact_27_0.png

From the tracks, we generate the hazards (the tracks are only the coordinates of the center of the cyclones, the full cyclones however affects a region around the tracks).

First thing we define the set of centroids which are geographical points where the hazard has a defined value. In our case, we want to define windspeeds from the tracks.

Remember: In the impact computations, for each exposure geographical point, one must assign a centroid from the hazard. By default, each exposure is assigned to the closest centroid from the hazard. But one can also define manually which centroid is assigned to which exposure point.

Examples: - Define the exposures from a given source (e.g., raster of asset values from LitPop). Define the hazard centroids from the exposures’ geolocations (e.g. compute Tropical Cyclone windspeed at each raster point and assign centroid to each raster point). - Define the exposures from a given source (e.g. houses position and value). Define the hazard from a given source (e.g. where lanslides occur). Use a metric to assign to each exposures point a hazard centroid (all houses in a radius of 5km around the lanslide are assigned to this centroid, if a house is within 5km of two landslides, choose the closest one). - Define a geographical raster. Define the exposures value on this raster. Define the hazard centroids on the geographical raster.

We shall pursue with the first case (Litpop + TropicalCyclone)

Hint: computing the wind speeds in many locations for many tc tracks is a computationally costly operation. Thus, we should define centroids only where we also have an exposure.

[6]:
# Define the centroids from the exposures position
centrs = Centroids()
lat = exp_lp.gdf['latitude'].values
lon = exp_lp.gdf['longitude'].values
centrs.set_lat_lon(lat, lon)
centrs.check()
[7]:
# Using the tracks, compute the windspeed at the location of the centroids
tc = TropCyclone()
tc.set_from_tracks(ibtracks_na, centrs)
tc.check()
2020-10-20 17:09:42,736 - climada.hazard.centroids.centr - INFO - Convert centroids to GeoSeries of Point shapes.
2020-10-20 17:09:43,110 - climada.util.coordinates - INFO - dist_to_coast: UTM 32616 (1/3)
2020-10-20 17:09:43,429 - climada.util.coordinates - INFO - dist_to_coast: UTM 32617 (2/3)
2020-10-20 17:09:44,115 - climada.util.coordinates - INFO - dist_to_coast: UTM 32618 (3/3)
2020-10-20 17:09:44,569 - climada.hazard.trop_cyclone - INFO - Mapping 120 tracks to 1387 centroids.
2020-10-20 17:09:45,232 - climada.hazard.trop_cyclone - INFO - Progress: 10%
2020-10-20 17:09:45,889 - climada.hazard.trop_cyclone - INFO - Progress: 20%
2020-10-20 17:09:46,626 - climada.hazard.trop_cyclone - INFO - Progress: 30%
2020-10-20 17:09:47,521 - climada.hazard.trop_cyclone - INFO - Progress: 40%
2020-10-20 17:09:48,289 - climada.hazard.trop_cyclone - INFO - Progress: 50%
2020-10-20 17:09:48,634 - climada.hazard.trop_cyclone - INFO - Progress: 60%
2020-10-20 17:09:49,142 - climada.hazard.trop_cyclone - INFO - Progress: 70%
2020-10-20 17:09:50,001 - climada.hazard.trop_cyclone - INFO - Progress: 80%
2020-10-20 17:09:50,660 - climada.hazard.trop_cyclone - INFO - Progress: 90%

Hint: The operation of computing the windspeed in different location is in general computationally expensive. Hence, if you have a lot of tropical cyclone tracks, you should first make sure that all your tropical cyclones actually affect your exposure (remove those that don’t). Then, be careful when defining the centroids. For a large country like China, there is no need for centroids 500km inlans (no tropical cyclones gets so far).

Impact function

For Tropical Cyclones, some calibrated default impact functions exist. Here we will use the one from Emanuel (2011).

[8]:
from climada.entity import ImpactFuncSet, IFTropCyclone
# impact function TC
impf_tc= IFTropCyclone()
impf_tc.set_emanuel_usa()

# add the impact function to an Impact function set
impf_set = ImpactFuncSet()
impf_set.append(impf_tc)
impf_set.check()
2020-10-20 17:09:51,501 - climada.entity.impact_funcs.base - WARNING - For intensity = 0, mdd != 0 or paa != 0. Consider shifting the origin of the intensity scale. In impact.calc the impact is always null at intensity = 0.

Recall that the exposures, hazards and impact functions must be matched in the impact calculations. Here it is simple, since there is a single impact function for all the hazards. We must simply make sure that the exposure is assigned this impact function through renaming the if_ column from the hazard type of the impact function in the impact function set and set the values of the column to the id of the impact function.

[39]:
# Get the hazard type and hazard id
[haz_type] = impf_set.get_hazard_types()
[haz_id] = impf_set.get_ids()[haz_type]
print(f"hazard type: {haz_type}, hazard id: {haz_id}")
hazard type: TC, hazard id: 1
[40]:
# Exposures: rename column and assign id
exp_lp.gdf.rename(columns={"if_": "if_" + haz_type}, inplace=True)
exp_lp.gdf['if_' + haz_type] = haz_id
exp_lp.check()
exp_lp.gdf.head()
2020-10-20 11:59:49,326 - climada.entity.exposures.base - INFO - centr_ not set.
2020-10-20 11:59:49,327 - climada.entity.exposures.base - INFO - deductible not set.
2020-10-20 11:59:49,328 - climada.entity.exposures.base - INFO - cover not set.
2020-10-20 11:59:49,329 - climada.entity.exposures.base - INFO - category_id not set.
2020-10-20 11:59:49,329 - climada.entity.exposures.base - INFO - geometry not set.
[40]:
value latitude longitude region_id if_TC
0 9.818338e+05 21.875000 -84.875000 192 1
1 1.110009e+06 21.875000 -84.791667 192 1
2 0.000000e+00 21.958333 -84.708333 192 1
3 1.026952e+06 21.958333 -84.625000 192 1
4 1.109359e+06 21.958333 -84.541667 192 1

Impact computation

We are finally ready for the impact computation. This is the simplest step. Just give the exposure, impact function and hazard to the Impact.calc() method.

Note: we did not specifically assign centroids to the exposures. Hence, the default is used - each exposure is associated with the closest centroids. Since we defined the centroids from the exposures, this is a one-to-one mapping.

Note: we did not define an Entity in this impact calculations. Recall that Entity is a container class for Exposures, Impact Functions, Discount Rates and Measures. Since we had only one Exposure and one Impact Function, the container would not have added any value, but for more complex projects, the Entity class is very useful.

[41]:
# Compute impact
from climada.engine import Impact
imp = Impact()
imp.calc(exp_lp, impf_set, tc, save_mat=False) #Do not save the results geographically resolved (only aggregate values)
2020-10-20 11:59:49,342 - climada.entity.exposures.base - INFO - Matching 1387 exposures with 1387 centroids.
2020-10-20 11:59:49,345 - climada.engine.impact - INFO - Calculating damage for 1373 assets (>0) and 120 events.
[9]:
exp_lp.gdf
[9]:
value latitude longitude region_id if_
0 9.818338e+05 21.875000 -84.875000 192 1
1 1.110009e+06 21.875000 -84.791667 192 1
2 0.000000e+00 21.958333 -84.708333 192 1
3 1.026952e+06 21.958333 -84.625000 192 1
4 1.109359e+06 21.958333 -84.541667 192 1
... ... ... ... ... ...
1382 4.015236e+05 20.291667 -74.291667 192 1
1383 2.292229e+06 20.208333 -74.291667 192 1
1384 2.293454e+06 20.125000 -74.291667 192 1
1385 2.290998e+06 20.291667 -74.208333 192 1
1386 2.292228e+06 20.208333 -74.208333 192 1

1387 rows × 5 columns

For example we can now obtain the aggregated average annual impact or plot the average annual impact in each exposure location.

[42]:
print(f"Aggregated average annual impact: {round(imp.aai_agg,0)} $")
Aggregated average annual impact: 346729533.0 $
[43]:
imp.plot_hexbin_eai_exposure(buffer=1)
[43]:
<cartopy.mpl.geoaxes.GeoAxesSubplot at 0x7f9301412990>
../_images/tutorial_climada_engine_Impact_46_1.png
[53]:
# Compute exceedance frequency curve
freq_curve = imp.calc_freq_curve()
freq_curve.plot()
[53]:
<matplotlib.axes._subplots.AxesSubplot at 0x7f8f7608dcd0>
../_images/tutorial_climada_engine_Impact_47_1.png

Quick examples - points, raster, custom

User defined point exposure and Tropical Cyclone hazard

[23]:
%matplotlib inline
# EXAMPLE: POINT EXPOSURES WITH POINT HAZARD
import numpy as np
from climada.entity import Exposures, ImpactFuncSet, IFTropCyclone
from climada.hazard import Centroids, TCTracks, TropCyclone
from climada.engine import Impact

# Set Exposures in points
exp_pnt = Exposures(crs='epsg:4326') #set coordinate system
exp_pnt['latitude'] = np.array([21.899326, 21.960728, 22.220574, 22.298390, 21.787977, 21.787977, 21.981732])
exp_pnt['longitude'] = np.array([88.307422, 88.565362, 88.378337, 87.806356, 88.348835, 88.348835, 89.246521])
exp_pnt['value'] = np.array([1.0e5, 1.2e5, 1.1e5, 1.1e5, 2.0e5, 2.5e5, 0.5e5])
exp_pnt.check()
exp_pnt.plot_scatter(buffer=0.05)

# Set Hazard in Exposures points
# set centroids from exposures coordinates
centr_pnt = Centroids()
centr_pnt.set_lat_lon(exp_pnt.latitude.values, exp_pnt.longitude.values, exp_pnt.crs)
# compute Hazard in that centroids
tr_pnt = TCTracks()
tr_pnt.read_ibtracs_netcdf(storm_id='2007314N10093')
tc_pnt = TropCyclone()
tc_pnt.set_from_tracks(tr_pnt, centroids=centr_pnt)
tc_pnt.check()
ax_pnt = tc_pnt.centroids.plot(c=np.array(tc_pnt.intensity[0,:].todense()).squeeze()) # plot intensity per point
ax_pnt.get_figure().colorbar(ax_pnt.collections[0], fraction=0.0175, pad=0.02).set_label('Intensity (m/s)') # add colorbar

# Set impact function
impf_pnt = ImpactFuncSet()
impf_tc = IFTropCyclone()
impf_tc.set_emanuel_usa()
impf_pnt.append(impf_tc)
impf_pnt.check()

# Get the hazard type and hazard id
[haz_type] = impf_set.get_hazard_types()
[haz_id] = impf_set.get_ids()[haz_type]
# Exposures: rename column and assign id
exp_lp.gdf.rename(columns={"if_": "if_" + haz_type}, inplace=True)
exp_lp.gdf['if_' + haz_type] = haz_id
exp_lp.gdf.head()

# Compute Impact
imp_pnt = Impact()
imp_pnt.calc(exp_pnt, impf_pnt, tc_pnt)
# nearest neighbor of exposures to centroids gives identity
print('Nearest neighbor hazard.centroids indexes for each exposure:', exp_pnt.centr_TC.values)
imp_pnt.plot_scatter_eai_exposure(ignore_zero=False, buffer=0.05)
2020-10-19 11:58:14,526 - climada.entity.exposures.base - INFO - tag metadata set to default value:  File:
 Description:
2020-10-19 11:58:14,527 - climada.entity.exposures.base - INFO - ref_year metadata set to default value: 2018
2020-10-19 11:58:14,527 - climada.entity.exposures.base - INFO - value_unit metadata set to default value: USD
2020-10-19 11:58:14,528 - climada.entity.exposures.base - INFO - meta metadata set to default value: None
2020-10-19 11:58:14,529 - climada.entity.exposures.base - INFO - Setting if_ to default impact functions ids 1.
2020-10-19 11:58:14,530 - climada.entity.exposures.base - INFO - centr_ not set.
2020-10-19 11:58:14,530 - climada.entity.exposures.base - INFO - deductible not set.
2020-10-19 11:58:14,531 - climada.entity.exposures.base - INFO - cover not set.
2020-10-19 11:58:14,531 - climada.entity.exposures.base - INFO - category_id not set.
2020-10-19 11:58:14,532 - climada.entity.exposures.base - INFO - region_id not set.
2020-10-19 11:58:14,533 - climada.entity.exposures.base - INFO - geometry not set.
2020-10-19 11:58:15,640 - climada.hazard.centroids.centr - INFO - Convert centroids to GeoSeries of Point shapes.
2020-10-19 11:58:15,903 - climada.util.coordinates - INFO - dist_to_coast: UTM 32645 (1/1)
2020-10-19 11:58:16,138 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 7 centroids.
2020-10-19 11:58:16,398 - climada.entity.impact_funcs.base - WARNING - For intensity = 0, mdd != 0 or paa != 0. Consider shifting the origin of the intensity scale. In impact.calc the impact is always null at intensity = 0.
2020-10-19 11:58:16,401 - climada.entity.exposures.base - INFO - Matching 7 exposures with 7 centroids.
2020-10-19 11:58:16,403 - climada.engine.impact - INFO - Calculating damage for 7 assets (>0) and 1 events.
2020-10-19 11:58:16,404 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
Nearest neighbor hazard.centroids indexes for each exposure: [0 1 2 3 4 5 6]
[23]:
<cartopy.mpl.geoaxes.GeoAxesSubplot at 0x7f92f65c0c50>
../_images/tutorial_climada_engine_Impact_50_2.png
../_images/tutorial_climada_engine_Impact_50_3.png
../_images/tutorial_climada_engine_Impact_50_4.png

Raster from file

[25]:
# EXAMPLE: RASTER EXPOSURES WITH RASTER HAZARD
from rasterio.warp import Resampling
from climada.entity import LitPop, ImpactFuncSet, ImpactFunc
from climada.hazard import Hazard
from climada.engine import Impact
from climada.util.constants import HAZ_DEMO_FL

# Exposures belonging to a raster (the raser information is contained in the meta attribute)
exp_ras = LitPop()
exp_ras.set_country(countries=['VEN'], res_km=10, fin_mode='income_group')
exp_ras.reset_index()
exp_ras.check()
exp_ras.plot_raster()
print('\n Raster properties exposures:', exp_ras.meta)

# Initialize hazard object with haz_type = 'FL' (for Flood)
hazard_type='FL'
haz_ras = Hazard(haz_type=hazard_type)
# Load a previously generated (either with CLIMADA or other means) hazard
# from file (HAZ_DEMO_FL) and resample the hazard raster to the exposures' ones
# Hint: check how other resampling methods affect to final impact
haz_ras.set_raster([HAZ_DEMO_FL], dst_crs=exp_ras.meta['crs'], transform=exp_ras.meta['transform'],
                   width=exp_ras.meta['width'], height=exp_ras.meta['height'],
                   resampling=Resampling.nearest)
haz_ras.intensity[haz_ras.intensity==-9999] = 0 # correct no data values
haz_ras.check()
haz_ras.plot_intensity(1)
print('Raster properties centroids:', haz_ras.centroids.meta)

# Set dummy impact function
impf_dum = ImpactFunc()
impf_dum.id = haz_id
impf_dum.name = 'dummy'
impf_dum.intensity_unit = 'm'
impf_dum.haz_type = hazard_type
impf_dum.intensity = np.linspace(0, 10, 100)
impf_dum.mdd = np.linspace(0, 10, 100)
impf_dum.paa = np.ones(impf_dum.intensity.size)
# Add the impact function to the impact function set
impf_ras = ImpactFuncSet()
impf_ras.append(impf_dum)
impf_ras.check()

# Exposures: rename column and assign id
exp_lp.gdf.rename(columns={"if_": "if_" + hazard_type}, inplace=True)
exp_lp.gdf['if_' + haz_type] = haz_id
exp_lp.gdf.head()

# Compute impact
imp_ras = Impact()
imp_ras.calc(exp_ras, impf_ras, haz_ras, save_mat=False)
# nearest neighbor of exposures to centroids is not identity because litpop does not contain data outside the country polygon
print('\n Nearest neighbor hazard.centroids indexes for each exposure:', exp_ras.centr_FL.values)
imp_ras.plot_raster_eai_exposure()
2020-10-19 12:07:38,904 - climada.entity.exposures.litpop - INFO - Generating LitPop data at a resolution of 300.0 arcsec.
2020-10-19 12:07:41,501 - climada.entity.exposures.gpw_import - INFO - Reference year: 2016. Using nearest available year for GWP population data: 2015
2020-10-19 12:07:41,501 - climada.entity.exposures.gpw_import - INFO - GPW Version v4.11
2020-10-19 12:07:53,519 - climada.util.finance - INFO - GDP VEN 2014: 4.824e+11.
2020-10-19 12:07:53,587 - climada.util.finance - INFO - Income group VEN 2016: 3.
2020-10-19 12:07:53,934 - climada.entity.exposures.litpop - INFO - Creating the LitPop exposure took 16 s
2020-10-19 12:07:53,935 - climada.entity.exposures.base - INFO - Hazard type not set in if_
2020-10-19 12:07:53,936 - climada.entity.exposures.base - INFO - centr_ not set.
2020-10-19 12:07:53,936 - climada.entity.exposures.base - INFO - deductible not set.
2020-10-19 12:07:53,937 - climada.entity.exposures.base - INFO - cover not set.
2020-10-19 12:07:53,937 - climada.entity.exposures.base - INFO - category_id not set.
2020-10-19 12:07:53,938 - climada.entity.exposures.base - INFO - geometry not set.
2020-10-19 12:07:53,956 - climada.entity.exposures.base - INFO - Hazard type not set in if_
2020-10-19 12:07:53,957 - climada.entity.exposures.base - INFO - centr_ not set.
2020-10-19 12:07:53,958 - climada.entity.exposures.base - INFO - deductible not set.
2020-10-19 12:07:53,959 - climada.entity.exposures.base - INFO - cover not set.
2020-10-19 12:07:53,959 - climada.entity.exposures.base - INFO - category_id not set.
2020-10-19 12:07:53,960 - climada.entity.exposures.base - INFO - geometry not set.
2020-10-19 12:07:53,966 - climada.util.coordinates - INFO - Raster from resolution 0.0833333333333286 to 0.0833333333333286.

 Raster properties exposures: {'width': 163, 'height': 138, 'crs': {'init': 'epsg:4326', 'no_defs': True}, 'transform': Affine(0.0833333333333286, 0.0, -73.41666666666669,
       0.0, -0.0833333333333286, 12.166666666666664)}
2020-10-19 12:07:56,272 - climada.util.coordinates - INFO - Reading /Users/ckropf/Documents/Climada/climada_python/data/demo/SC22000_VE__M1.grd.gz
Raster properties centroids: {'driver': 'GSBG', 'dtype': 'float32', 'nodata': 1.701410009187828e+38, 'width': 163, 'height': 138, 'count': 1, 'crs': {'init': 'epsg:4326', 'no_defs': True}, 'transform': Affine(0.0833333333333286, 0.0, -73.41666666666669,
       0.0, -0.0833333333333286, 12.166666666666664)}
2020-10-19 12:07:59,125 - climada.entity.impact_funcs.base - WARNING - For intensity = 0, mdd != 0 or paa != 0. Consider shifting the origin of the intensity scale. In impact.calc the impact is always null at intensity = 0.
2020-10-19 12:07:59,127 - climada.entity.exposures.base - INFO - Matching 10770 exposures with 22494 centroids.
2020-10-19 12:07:59,131 - climada.engine.impact - INFO - Calculating damage for 10717 assets (>0) and 1 events.
2020-10-19 12:07:59,132 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_FL. Using impact functions in if_.

 Nearest neighbor hazard.centroids indexes for each exposure: [5705 5543 5706 ... 7659 7822 7660]
2020-10-19 12:07:59,140 - climada.util.coordinates - INFO - Raster from resolution 0.0833333333333286 to 0.0833333333333286.
[25]:
<cartopy.mpl.geoaxes.GeoAxesSubplot at 0x7f92fdd4d3d0>
../_images/tutorial_climada_engine_Impact_52_2.png
../_images/tutorial_climada_engine_Impact_52_3.png
../_images/tutorial_climada_engine_Impact_52_4.png

## VISUALIZATION

Making plots

The expected annual impact per exposure can be visualized through different methods: plot_hexbin_eai_exposure(), plot_scatter_eai_exposur(), plot_raster_eai_exposure() and plot_basemap_eai_exposure() (similarly as with Exposures).

[4]:
imp_pnt.plot_basemap_eai_exposure(buffer=5000)
2020-10-16 11:22:36,585 - climada.util.coordinates - INFO - Setting geometry points.
2020-10-16 11:22:36,686 - climada.entity.exposures.base - INFO - Setting latitude and longitude attributes.
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/pyproj/crs/crs.py:53: FutureWarning: '+init=<authority>:<code>' syntax is deprecated. '<authority>:<code>' is the preferred initialization method. When making the change, be mindful of axis order changes: https://pyproj4.github.io/pyproj/stable/gotchas.html#axis-order-changes-in-proj-6
  return _prepare_from_string(" ".join(pjargs))
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/cartopy/mpl/gridliner.py:307: UserWarning: The .xlabels_top attribute is deprecated. Please use .top_labels to toggle visibility instead.
  warnings.warn('The .xlabels_top attribute is deprecated. Please '
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/cartopy/mpl/gridliner.py:343: UserWarning: The .ylabels_right attribute is deprecated. Please use .right_labels to toggle visibility instead.
  warnings.warn('The .ylabels_right attribute is deprecated. Please '
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/contextily/tile.py:268: FutureWarning: The url format using 'tileX', 'tileY', 'tileZ' as placeholders is deprecated. Please use '{x}', '{y}', '{z}' instead.
  FutureWarning,
2020-10-16 11:22:38,961 - climada.entity.exposures.base - INFO - Setting latitude and longitude attributes.
/Users/ckropf/opt/anaconda3/envs/climada_env2/lib/python3.7/site-packages/pyproj/crs/crs.py:53: FutureWarning: '+init=<authority>:<code>' syntax is deprecated. '<authority>:<code>' is the preferred initialization method. When making the change, be mindful of axis order changes: https://pyproj4.github.io/pyproj/stable/gotchas.html#axis-order-changes-in-proj-6
  return _prepare_from_string(" ".join(pjargs))
[4]:
<cartopy.mpl.geoaxes.GeoAxesSubplot at 0x7f80657fe150>
../_images/tutorial_climada_engine_Impact_56_5.png

Making videos

Given a fixed exposure and impact functions, a sequence of hazards can be visualized hitting the exposures.

[4]:
# exposure
from climada.entity import BlackMarble, add_sea
country_name = {'United States Of America': ['Florida']}
exp_video = BlackMarble()
exp_video.set_countries(country_name, 2016, res_km=2.5)
exp_video.check()

# impact function
if_def = IFTropCyclone()
if_def.set_emanuel_usa()
ifs_video = ImpactFuncSet()
ifs_video.append(if_def)
ifs_video.check()

# compute sequence of hazards using TropCyclone video_intensity method
exp_sea = add_sea(exp_video, (100, 5))
centr_video = Centroids()
centr_video.set_lat_lon(exp_sea.latitude.values, exp_sea.longitude.values)
centr_video.check()

track_name = '2017242N16333'
tr_irma = TCTracks()
tr_irma.read_ibtracs_netcdf(provider='usa', storm_id=track_name) # IRMA 2017

tc_video = TropCyclone()
tc_list, _ = tc_video.video_intensity(track_name, tr_irma, centr_video) # empty file name to not to write the video

# generate video of impacts
file_name='./results/irma_imp_fl.gif'
imp_video = Impact()
imp_list = imp_video.video_direct_impact(exp_video, ifs_video, tc_list, file_name)

2019-10-29 18:17:55,127 - climada.util.finance - INFO - GDP USA 2016: 1.871e+13.
2019-10-29 18:17:55,180 - climada.util.finance - INFO - Income group USA 2016: 4.
2019-10-29 18:17:55,180 - climada.entity.exposures.black_marble - INFO - Nightlights from NASA's earth observatory for year 2016.
2019-10-29 18:17:55,796 - climada.entity.exposures.nightlight - DEBUG - Found all required satellite data (4 files) in folder /Users/aznarsig/Documents/Python/climada_python/data/system
2019-10-29 18:17:55,797 - climada.entity.exposures.nightlight - DEBUG - All required files already exist. No downloads necessary.
2019-10-29 18:21:45,895 - climada.entity.exposures.black_marble - INFO - Processing country United States Of America.
2019-10-29 18:29:41,450 - climada.entity.exposures.black_marble - INFO - Generating resolution of approx 2.5 km.
2019-10-29 18:29:42,812 - climada.entity.exposures.base - INFO - Hazard type not set in if_
2019-10-29 18:29:42,815 - climada.entity.exposures.base - INFO - centr_ not set.
2019-10-29 18:29:42,820 - climada.entity.exposures.base - INFO - deductible not set.
2019-10-29 18:29:42,820 - climada.entity.exposures.base - INFO - cover not set.
2019-10-29 18:29:42,821 - climada.entity.exposures.base - INFO - category_id not set.
2019-10-29 18:29:42,822 - climada.entity.exposures.base - INFO - geometry not set.
2019-10-29 18:29:42,830 - climada.entity.exposures.base - INFO - Adding sea at 5 km resolution and 100 km distance from coast.
2019-10-29 18:29:48,407 - climada.hazard.tc_tracks - INFO - Reading 2017242N16333: IRMA
2019-10-29 18:29:48,653 - climada.hazard.centroids.centr - INFO - Setting geometry points.
2019-10-29 18:29:51,700 - climada.hazard.centroids.centr - DEBUG - Setting dist_coast 58555 points.
2019-10-29 18:30:01,391 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:01,481 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:01,492 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:01,498 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:01,578 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:01,589 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:01,593 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:01,671 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:01,684 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:01,689 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:01,774 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:01,788 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:01,793 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:01,879 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:01,894 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:01,899 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:01,987 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,003 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,008 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:02,091 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,108 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,113 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:02,193 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,210 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,214 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:02,293 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,310 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,315 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:02,397 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,416 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,421 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:02,507 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,526 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,535 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:02,664 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,684 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,688 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:02,778 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,802 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,807 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:02,907 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:02,928 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:02,933 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:03,030 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:03,051 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:03,055 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:03,139 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:03,159 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:03,163 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:03,245 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:03,270 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:03,277 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:03,418 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:03,443 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:03,449 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:03,573 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:03,594 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:03,729 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:04,110 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:04,130 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:04,137 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:04,253 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:04,268 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:04,275 - climada.hazard.trop_cyclone - INFO - Mapping 1 tracks to 58555 centroids.
2019-10-29 18:30:04,395 - climada.hazard.trop_cyclone - DEBUG - Append events.
2019-10-29 18:30:04,406 - climada.hazard.trop_cyclone - DEBUG - Compute frequency.
2019-10-29 18:30:04,414 - climada.entity.exposures.base - INFO - Matching 30986 exposures with 58555 centroids.
2019-10-29 18:30:06,830 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,831 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,835 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,837 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,838 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,843 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,845 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,845 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,851 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,853 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,854 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,859 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,861 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,862 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,867 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,869 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,869 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,874 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,876 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,877 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,882 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,884 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,884 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,889 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,893 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,894 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,899 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,902 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,902 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,916 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,918 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,919 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,930 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,932 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,933 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,942 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,946 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,947 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:06,966 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:06,981 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:06,985 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,143 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:07,147 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:07,150 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,182 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:07,197 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:07,211 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,246 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:07,261 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:07,262 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,276 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:07,280 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:07,282 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,292 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:07,295 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:07,296 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,305 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:07,309 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:07,310 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,318 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:07,321 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:07,323 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,334 - climada.engine.impact - INFO - Exposures matching centroids found in centr_TC
2019-10-29 18:30:07,337 - climada.engine.impact - INFO - Calculating damage for 24024 assets (>0) and 1 events.
2019-10-29 18:30:07,338 - climada.engine.impact - INFO - Missing exposures impact functions for hazard if_TC. Using impact functions in if_.
2019-10-29 18:30:07,353 - climada.engine.impact - INFO - Generating video ./results/irma_imp_fl.gif
23it [07:27, 28.51s/it]
../_images/tutorial_climada_engine_Impact_58_2.png