openquake.risklib package

openquake.risklib.riskinput module

class openquake.risklib.riskinput.CompositeRiskModel(oqparam, rmdict, retrodict)[source]

Bases: collections.abc.Mapping

A container (imt, taxonomy) -> riskmodel

Parameters:
gen_outputs(riskinput, monitor, hazard=None)[source]

Group the assets per taxonomy and compute the outputs by using the underlying riskmodels. Yield the outputs generated as dictionaries out_by_lr.

Parameters:
  • riskinput – a RiskInput instance
  • monitor – a monitor object used to measure the performance
get_assets_ratios(assets, gmvs, imts)[source]
Parameters:
  • assets – assets on the same site
  • imts – intensity measure types
Params gmvs:

hazard on the given site, shape (E, M)

Returns:

a list of (assets, loss_ratios) for each taxonomy on the site

get_extra_imts(imts)[source]

Returns the extra IMTs in the risk functions, i.e. the ones not in the imts set (the set of IMTs for which there is hazard).

get_loss_ratios()[source]
Returns:a 1-dimensional composite array with loss ratios by loss type
init(oqparam)[source]
make_curve_params(oqparam)[source]
reduce(taxonomies)[source]
Parameters:taxonomies – a set of taxonomies
Returns:a new CompositeRiskModel reduced to the given taxonomies
taxonomy_dict
Returns:a dict taxonomy string -> taxonomy index
class openquake.risklib.riskinput.EpsilonMatrix0(num_assets, seeds)[source]

Bases: object

Mock-up for a matrix of epsilons of size N x E, used when asset_correlation=0.

Parameters:
  • num_assets – N assets
  • seeds – E seeds, set before calling numpy.random.normal
make_eps()[source]

Builds a matrix of A x E epsilons

class openquake.risklib.riskinput.EpsilonMatrix1(num_assets, num_events, seed)[source]

Bases: object

Mock-up for a matrix of epsilons of size A x E, used when asset_correlation=1.

Parameters:
  • num_assets – number of assets
  • num_events – number of events
  • seed – seed used to generate E epsilons
class openquake.risklib.riskinput.RiskInput(hazard_getter, assets_by_site, eps_dict=None)[source]

Bases: object

Contains all the assets and hazard values associated to a given imt and site.

Parameters:
  • hazard_getter – a callable returning the hazard data for a given realization
  • assets_by_site – array of assets, one per site
  • eps_dict – dictionary of epsilons (can be None)
epsilon_getter(aid, eids)[source]
Parameters:
  • aid – asset ordinal
  • eids – an array of event indices
Returns:

an array of E epsilons

imt_taxonomies

Return a list of pairs (imt, taxonomies) with a single element

exception openquake.risklib.riskinput.ValidationError[source]

Bases: Exception

openquake.risklib.riskinput.make_eps(assetcol, num_samples, seed, correlation)[source]
Parameters:
  • assetcol – an AssetCollection instance
  • num_samples (int) – the number of ruptures
  • seed (int) – a random seed
  • correlation (float) – the correlation coefficient
Returns:

epsilons matrix of shape (num_assets, num_samples)

openquake.risklib.riskinput.make_epsilon_getter(n_assets, n_events, correlation, master_seed, no_eps)[source]
Returns:a function (start, stop) -> matrix of shape (n_assets, n_events)
openquake.risklib.riskinput.read_composite_risk_model(dstore)[source]
Parameters:dstore – a DataStore instance
Returns:a CompositeRiskModel instance
openquake.risklib.riskinput.rsi2str(rlzi, sid, imt)[source]

Convert a triple (XXXX, YYYY, ZZZ) into a string of the form ‘rlz-XXXX/sid-YYYY/ZZZ’

openquake.risklib.riskinput.str2rsi(key)[source]

Convert a string of the form ‘rlz-XXXX/sid-YYYY/ZZZ’ into a triple (XXXX, YYYY, ZZZ)

openquake.risklib.riskmodels module

class openquake.risklib.riskmodels.Classical(taxonomy, vulnerability_functions, hazard_imtls, lrem_steps_per_interval, conditional_loss_poes, poes_disagg, insured_losses=False)[source]

Bases: openquake.risklib.riskmodels.RiskModel

Classical PSHA-Based RiskModel. Computes loss curves and insured curves.

kind = 'vulnerability'
class openquake.risklib.riskmodels.ClassicalBCR(taxonomy, vulnerability_functions_orig, vulnerability_functions_retro, hazard_imtls, lrem_steps_per_interval, interest_rate, asset_life_expectancy)[source]

Bases: openquake.risklib.riskmodels.RiskModel

kind = 'vulnerability'
class openquake.risklib.riskmodels.ClassicalDamage(taxonomy, fragility_functions, hazard_imtls, investigation_time, risk_investigation_time)[source]

Bases: openquake.risklib.riskmodels.Damage

Implements the ClassicalDamage riskmodel. Computes the damages.

kind = 'fragility'
class openquake.risklib.riskmodels.Damage(taxonomy, fragility_functions)[source]

Bases: openquake.risklib.riskmodels.RiskModel

Implements the ScenarioDamage riskmodel. Computes the damages.

kind = 'fragility'
class openquake.risklib.riskmodels.ProbabilisticEventBased(taxonomy, vulnerability_functions, conditional_loss_poes)[source]

Bases: openquake.risklib.riskmodels.RiskModel

Implements the Probabilistic Event Based riskmodel. Computes loss ratios and event IDs.

get_loss_ratios(gmvs, imti)[source]
Parameters:
  • gmvs – an array of shape (E, M)
  • imti – a dictionary imt -> imt index
Returns:

loss_ratios of shape (L, E)

kind = 'vulnerability'
class openquake.risklib.riskmodels.RiskModel(taxonomy, risk_functions, insured_losses)[source]

Bases: object

Base class. Can be used in the tests as a mock.

compositemodel = None
get_loss_types(imt)[source]
Parameters:imt – Intensity Measure Type string
Returns:loss types with risk functions of the given imt
get_output(assets, data_by_lt, epsgetter)[source]
Parameters:
  • assets – a list of assets with the same taxonomy
  • data_by_lt – hazards for each loss type
  • epsgetter – an epsilon getter function
Returns:

an ArrayWrapper of shape (L, …)

kind = None
loss_types

The list of loss types in the underlying vulnerability functions, in lexicographic order

time_event = None
class openquake.risklib.riskmodels.Scenario(taxonomy, vulnerability_functions, insured_losses, time_event=None)[source]

Bases: openquake.risklib.riskmodels.RiskModel

Implements the Scenario riskmodel. Computes the loss matrix.

kind = 'vulnerability'
openquake.risklib.riskmodels.build_vf_node(vf)[source]

Convert a VulnerabilityFunction object into a Node suitable for XML conversion.

openquake.risklib.riskmodels.filter_vset(elem)[source]
openquake.risklib.riskmodels.get_risk_files(inputs)[source]
Parameters:inputs – a dictionary key -> path name
Returns:a pair (file_type, {risk_type: path})
openquake.risklib.riskmodels.get_risk_models(oqparam, kind=None)[source]
Parameters:
  • oqparam – an OqParam instance
  • kind – vulnerability|vulnerability_retrofitted|fragility|consequence; if None it is extracted from the oqparam.file_type attribute
Returns:

a dictionary taxonomy -> loss_type -> function

openquake.risklib.riskmodels.get_riskmodel(taxonomy, oqparam, **extra)[source]

Return an instance of the correct riskmodel class, depending on the attribute calculation_mode of the object oqparam.

Parameters:
  • taxonomy – a taxonomy string
  • oqparam – an object containing the parameters needed by the riskmodel class
  • extra – extra parameters to pass to the riskmodel class
openquake.risklib.riskmodels.get_values(loss_type, assets, time_event=None)[source]
Returns:a numpy array with the values for the given assets, depending on the loss_type.
openquake.risklib.riskmodels.rescale(curves, values)[source]

Multiply the losses in each curve of kind (losses, poes) by the corresponding value.

openquake.risklib.scientific module

This module includes the scientific API of the oq-risklib

class openquake.risklib.scientific.BetaDistribution[source]

Bases: openquake.risklib.scientific.Distribution

sample(means, _covs, stddevs, _idxs=None)[source]
survival(loss_ratio, mean, stddev)[source]
class openquake.risklib.scientific.ConsequenceFunction(id, dist, params)

Bases: tuple

dist

Alias for field number 1

id

Alias for field number 0

params

Alias for field number 2

class openquake.risklib.scientific.ConsequenceModel(id, assetCategory, lossCategory, description, limitStates)[source]

Bases: dict

Container for a set of consequence functions. You can access each function given its name with the square bracket notation.

Parameters:
  • id (str) – ID of the model
  • assetCategory (str) – asset category (i.e. buildings, population)
  • lossCategory (str) – loss type (i.e. structural, contents, …)
  • description (str) – description of the model
  • limitStates – a list of limit state strings
  • consequence_functions – a dictionary name -> ConsequenceFunction
class openquake.risklib.scientific.CurveParams(index, loss_type, curve_resolution, ratios, user_provided)

Bases: tuple

curve_resolution

Alias for field number 2

index

Alias for field number 0

loss_type

Alias for field number 1

ratios

Alias for field number 3

user_provided

Alias for field number 4

class openquake.risklib.scientific.DegenerateDistribution[source]

Bases: openquake.risklib.scientific.Distribution

The degenerate distribution. E.g. a distribution with a delta corresponding to the mean.

sample(means, _covs, _stddev, _idxs)[source]
survival(loss_ratio, mean, _stddev)[source]
class openquake.risklib.scientific.DiscreteDistribution[source]

Bases: openquake.risklib.scientific.Distribution

sample(loss_ratios, probs)[source]
seed = None
survival(loss_ratios, probs)[source]

Required for the Classical Risk and BCR Calculators. Currently left unimplemented as the PMF format is used only for the Scenario and Event Based Risk Calculators.

Parameters:steps (int) – number of steps between loss ratios.
class openquake.risklib.scientific.Distribution[source]

Bases: object

A Distribution class models continuous probability distribution of random variables used to sample losses of a set of assets. It is usually registered with a name (e.g. LN, BT, PM) by using openquake.baselib.general.CallableDict

sample(means, covs, stddevs, idxs)[source]
Returns:

sample a set of losses

Parameters:
  • means – an array of mean losses
  • covs – an array of covariances
  • stddevs – an array of stddevs
survival(loss_ratio, mean, stddev)[source]

Return the survival function of the distribution with mean and stddev applied to loss_ratio

class openquake.risklib.scientific.FragilityFunctionContinuous(limit_state, mean, stddev)[source]

Bases: object

class openquake.risklib.scientific.FragilityFunctionDiscrete(limit_state, imls, poes, no_damage_limit=None)[source]

Bases: object

interp
class openquake.risklib.scientific.FragilityFunctionList(array, **attrs)[source]

Bases: list

A list of fragility functions with common attributes; there is a function for each limit state.

build(limit_states, discretization, steps_per_interval)[source]
Parameters:
  • limit_states – a sequence of limit states
  • discretization – continouos fragility discretization parameter
  • steps_per_interval – steps_per_interval parameter
Returns:

a populated FragilityFunctionList instance

mean_loss_ratios_with_steps(steps)[source]

For compatibility with vulnerability functions

class openquake.risklib.scientific.FragilityModel(id, assetCategory, lossCategory, description, limitStates)[source]

Bases: dict

Container for a set of fragility functions. You can access each function given the IMT and taxonomy with the square bracket notation.

Parameters:
  • id (str) – ID of the model
  • assetCategory (str) – asset category (i.e. buildings, population)
  • lossCategory (str) – loss type (i.e. structural, contents, …)
  • description (str) – description of the model
  • limitStates – a list of limit state strings
build(continuous_fragility_discretization, steps_per_interval)[source]

Return a new FragilityModel instance, in which the values have been replaced with FragilityFunctionList instances.

Parameters:
  • continuous_fragility_discretization – configuration parameter
  • steps_per_interval – configuration parameter
class openquake.risklib.scientific.LogNormalDistribution(epsilons=None)[source]

Bases: openquake.risklib.scientific.Distribution

Model a distribution of a random variable whoose logarithm are normally distributed.

Attr epsilons:An array of random numbers generated with numpy.random.multivariate_normal() with size E
sample(means, covs, _stddevs, idxs)[source]
survival(loss_ratio, mean, stddev)[source]
class openquake.risklib.scientific.LossCurvesMapsBuilder(conditional_loss_poes, return_periods, loss_dt, weights, num_events, eff_time, risk_investigation_time)[source]

Bases: object

Build losses curves and maps for all loss types at the same time.

Parameters:
  • conditional_loss_poes – a list of PoEs, possibly empty
  • return_periods – ordered array of return periods
  • loss_dt – composite dtype for the loss types
  • weights – weights of the realizations
  • num_events – number of events for each realization
  • eff_time – ses_per_logic_tree_path * hazard investigation time
build(losses_by_event, stats=())[source]
Parameters:
  • losses_by_event – the aggregate loss table with shape R -> (E, L)
  • stats – list of pairs [(statname, statfunc), …]
Returns:

two arrays with shape (P, R, L) and (P, S, L)

build_curve(asset_value, loss_ratios, rlzi)[source]
build_curves_maps(loss_arrays, rlzi)[source]
build_loss_maps(losses, clp, stats=())[source]
Parameters:
  • losses – an array of shape R, E
  • clp – a list of C conditional loss poes
  • stats – list of pairs [(statname, statfunc), …]
Returns:

two arrays of shape (C, R) and (C, S)

build_maps(losses, clp, stats=())[source]
Parameters:
  • losses – an array of shape (A, R, P)
  • clp – a list of C conditional loss poes
  • stats – list of pairs [(statname, statfunc), …]
Returns:

an array of loss_maps of shape (A, R, C, LI)

build_pair(losses, stats)[source]
Parameters:losses – a list of lists with R elements
Returns:two arrays of shape (P, R) and (P, S) respectively
pair(array, stats)[source]

:return (array, array_stats) if stats, else (array, None)

class openquake.risklib.scientific.VulnerabilityFunction(vf_id, imt, imls, mean_loss_ratios, covs=None, distribution='LN')[source]

Bases: object

dtype = dtype([('iml', '<f4'), ('loss_ratio', '<f4'), ('cov', '<f4')])
init()[source]
interpolate(gmvs)[source]
Parameters:gmvs – array of intensity measure levels
Returns:(interpolated loss ratios, interpolated covs, indices > min)
loss_ratio_exceedance_matrix(steps)[source]

Compute the LREM (Loss Ratio Exceedance Matrix).

Parameters:steps (int) – Number of steps between loss ratios.
mean_imls()[source]

Compute the mean IMLs (Intensity Measure Level) for the given vulnerability function.

Parameters:vulnerability_function – the vulnerability function where the IMLs (Intensity Measure Level) are taken from.
mean_loss_ratios_with_steps(steps)[source]

Split the mean loss ratios, producing a new set of loss ratios. The new set of loss ratios always includes 0.0 and 1.0

Parameters:steps (int) –

the number of steps we make to go from one loss ratio to the next. For example, if we have [0.5, 0.7]:

steps = 1 produces [0.0,  0.5, 0.7, 1]
steps = 2 produces [0.0, 0.25, 0.5, 0.6, 0.7, 0.85, 1]
steps = 3 produces [0.0, 0.17, 0.33, 0.5, 0.57, 0.63,
                    0.7, 0.8, 0.9, 1]
sample(means, covs, idxs, epsilons)[source]

Sample the epsilons and apply the corrections to the means. This method is called only if there are nonzero covs.

Parameters:
  • means – array of E’ loss ratios
  • covs – array of E’ floats
  • idxs – array of E booleans with E >= E’
  • epsilons – array of E floats
Returns:

array of E’ loss ratios

set_distribution(epsilons=None)[source]
strictly_increasing()[source]
Returns:a new vulnerability function that is strictly increasing. It is built by removing piece of the function where the mean loss ratio is constant.
class openquake.risklib.scientific.VulnerabilityFunctionWithPMF(vf_id, imt, imls, loss_ratios, probs, seed=42)[source]

Bases: openquake.risklib.scientific.VulnerabilityFunction

Vulnerability function with an explicit distribution of probabilities

Parameters:
  • vf_id (str) – vulnerability function ID
  • imt (str) – Intensity Measure Type
  • imls – intensity measure levels (L)
  • ratios – an array of mean ratios (M)
  • probs – a matrix of probabilities of shape (M, L)
init()[source]
interpolate(gmvs)[source]
Parameters:gmvs – array of intensity measure levels
Returns:(interpolated probabilities, None, indices > min)
loss_ratio_exceedance_matrix(steps)[source]

Compute the LREM (Loss Ratio Exceedance Matrix). Required for the Classical Risk and BCR Calculators. Currently left unimplemented as the PMF format is used only for the Scenario and Event Based Risk Calculators.

Parameters:steps (int) – Number of steps between loss ratios.
sample(probs, _covs, idxs, epsilons)[source]

Sample the .loss_ratios with the given probabilities.

Parameters:
  • probs – array of E’ floats
  • _covs – ignored, it is there only for API consistency
  • idxs – array of E booleans with E >= E’
  • epsilons – array of E floats
Returns:

array of E’ probabilities

set_distribution(epsilons=None)[source]
class openquake.risklib.scientific.VulnerabilityModel(id=None, assetCategory=None, lossCategory=None)[source]

Bases: dict

Container for a set of vulnerability functions. You can access each function given the IMT and taxonomy with the square bracket notation.

Parameters:
  • id (str) – ID of the model
  • assetCategory (str) – asset category (i.e. buildings, population)
  • lossCategory (str) – loss type (i.e. structural, contents, …)

All such attributes are None for a vulnerability model coming from a NRML 0.4 file.

openquake.risklib.scientific.annual_frequency_of_exceedence(poe, t_haz)[source]
Parameters:
  • poe – array of probabilities of exceedence
  • t_haz – hazard investigation time
Returns:

array of frequencies (with +inf values where poe=1)

openquake.risklib.scientific.average_loss(losses_poes)[source]

Given a loss curve with poes over losses defined on a given time span it computes the average loss on this period of time.

Note:As the loss curve is supposed to be piecewise linear as it is a result of a linear interpolation, we compute an exact integral by using the trapeizodal rule with the width given by the loss bin width.
openquake.risklib.scientific.bcr(eal_original, eal_retrofitted, interest_rate, asset_life_expectancy, asset_value, retrofitting_cost)[source]

Compute the Benefit-Cost Ratio.

BCR = (EALo - EALr)(1-exp(-r*t))/(r*C)

Where:

  • BCR – Benefit cost ratio
  • EALo – Expected annual loss for original asset
  • EALr – Expected annual loss for retrofitted asset
  • r – Interest rate
  • t – Life expectancy of the asset
  • C – Retrofitting cost
openquake.risklib.scientific.broadcast(func, composite_array, *args)[source]

Broadcast an array function over a composite array

openquake.risklib.scientific.build_imls(ff, continuous_fragility_discretization, steps_per_interval=0)[source]

Build intensity measure levels from a fragility function. If the function is continuous, they are produced simply as a linear space between minIML and maxIML. If the function is discrete, they are generated with a complex logic depending on the noDamageLimit and the parameter steps per interval.

Parameters:
  • ff – a fragility function object
  • continuous_fragility_discretization – .ini file parameter
  • steps_per_interval – .ini file parameter
Returns:

generated imls

openquake.risklib.scientific.build_loss_curve_dt(curve_resolution, insured_losses=False)[source]
Parameters:
  • curve_resolution – dictionary loss_type -> curve_resolution
  • insured_losses – configuration parameter
Returns:

loss_curve_dt

openquake.risklib.scientific.classical(vulnerability_function, hazard_imls, hazard_poes, steps=10)[source]
Parameters:
  • vulnerability_function – an instance of openquake.risklib.scientific.VulnerabilityFunction representing the vulnerability function used to compute the curve.
  • hazard_imls – the hazard intensity measure type and levels
  • steps (int) – Number of steps between loss ratios.
openquake.risklib.scientific.classical_damage(fragility_functions, hazard_imls, hazard_poes, investigation_time, risk_investigation_time)[source]
Parameters:
  • fragility_functions – a list of fragility functions for each damage state
  • hazard_imls – Intensity Measure Levels
  • hazard_poes – hazard curve
  • investigation_time – hazard investigation time
  • risk_investigation_time – risk investigation time
Returns:

an array of M probabilities of occurrence where M is the numbers of damage states.

openquake.risklib.scientific.conditional_loss_ratio(loss_ratios, poes, probability)[source]

Return the loss ratio corresponding to the given PoE (Probability of Exceendance). We can have four cases:

  1. If probability is in poes it takes the bigger corresponding loss_ratios.
  2. If it is in (poe1, poe2) where both poe1 and poe2 are in poes, then we perform a linear interpolation on the corresponding losses
  3. if the given probability is smaller than the lowest PoE defined, it returns the max loss ratio .
  4. if the given probability is greater than the highest PoE defined it returns zero.
Parameters:
  • loss_ratios – an iterable over non-decreasing loss ratio values (float)
  • poes – an iterable over non-increasing probability of exceedance values (float)
  • probability (float) – the probability value used to interpolate the loss curve
openquake.risklib.scientific.fine_graining(points, steps)[source]
Parameters:
  • points – a list of floats
  • steps (int) – expansion steps (>= 2)
>>> fine_graining([0, 1], steps=0)
[0, 1]
>>> fine_graining([0, 1], steps=1)
[0, 1]
>>> fine_graining([0, 1], steps=2)
array([0. , 0.5, 1. ])
>>> fine_graining([0, 1], steps=3)
array([0.        , 0.33333333, 0.66666667, 1.        ])
>>> fine_graining([0, 0.5, 0.7, 1], steps=2)
array([0.  , 0.25, 0.5 , 0.6 , 0.7 , 0.85, 1.  ])

N points become S * (N - 1) + 1 points with S > 0

openquake.risklib.scientific.insured_loss_curve(curve, deductible, insured_limit)[source]

Compute an insured loss ratio curve given a loss ratio curve

Parameters:
  • curve – an array 2 x R (where R is the curve resolution)
  • deductible (float) – the deductible limit in fraction form
  • insured_limit (float) – the insured limit in fraction form
>>> losses = numpy.array([3, 20, 101])
>>> poes = numpy.array([0.9, 0.5, 0.1])
>>> insured_loss_curve(numpy.array([losses, poes]), 5, 100)
array([[ 3.        , 20.        ],
       [ 0.85294118,  0.5       ]])
openquake.risklib.scientific.insured_losses(losses, deductible, insured_limit)[source]
Parameters:
  • losses – an array of ground-up loss ratios
  • deductible (float) – the deductible limit in fraction form
  • insured_limit (float) – the insured limit in fraction form

Compute insured losses for the given asset and losses, from the point of view of the insurance company. For instance:

>>> insured_losses(numpy.array([3, 20, 101]), 5, 100)
array([ 0, 15, 95])
  • if the loss is 3 (< 5) the company does not pay anything
  • if the loss is 20 the company pays 20 - 5 = 15
  • if the loss is 101 the company pays 100 - 5 = 95
openquake.risklib.scientific.loss_maps(curves, conditional_loss_poes)[source]
Parameters:
  • curves – an array of loss curves
  • conditional_loss_poes – a list of conditional loss poes
Returns:

a composite array of loss maps with the same shape

openquake.risklib.scientific.losses_by_period(losses, return_periods, num_events=None, eff_time=None)[source]
Parameters:
  • losses – array of simulated losses
  • return_periods – return periods of interest
  • num_events – the number of events (>= to the number of losses)
  • eff_time – investigation_time * ses_per_logic_tree_path
Returns:

interpolated losses for the return periods, possibly with NaN

NB: the return periods must be ordered integers >= 1. The interpolated losses are defined inside the interval min_time < time < eff_time where min_time = eff_time /num_events. Outside the interval they have NaN values. Here is an example:

>>> losses = [3, 2, 3.5, 4, 3, 23, 11, 2, 1, 4, 5, 7, 8, 9, 13]
>>> losses_by_period(losses, [1, 2, 5, 10, 20, 50, 100], 20)
array([ nan,  nan,  0. ,  3.5,  8. , 13. , 23. ])

If num_events is not passed, it is inferred from the number of losses; if eff_time is not passed, it is inferred from the longest return period.

openquake.risklib.scientific.make_epsilons(matrix, seed, correlation)[source]

Given a matrix N * R returns a matrix of the same shape N * R obtained by applying the multivariate_normal distribution to N points and R samples, by starting from the given seed and correlation.

openquake.risklib.scientific.mean_std(fractions)[source]

Given an N x M matrix, returns mean and std computed on the rows, i.e. two M-dimensional vectors.

openquake.risklib.scientific.normalize_curves_eb(curves)[source]

A more sophisticated version of normalize_curves, used in the event based calculator.

Parameters:curves – a list of pairs (losses, poes)
Returns:first losses, all_poes
openquake.risklib.scientific.pairwise_diff(values)[source]

Differences between a value and the next value in a sequence

openquake.risklib.scientific.pairwise_mean(values)[source]

Averages between a value and the next value in a sequence

openquake.risklib.scientific.return_periods(eff_time, num_losses)[source]
Parameters:
  • eff_time – ses_per_logic_tree_path * investigation_time
  • num_losses – used to determine the minimum period
Returns:

an array of 32 bit periods

Here are a few examples:

>>> return_periods(1, 1)
Traceback (most recent call last):
   ...
AssertionError: eff_time too small: 1
>>> return_periods(2, 2)
array([1, 2], dtype=uint32)
>>> return_periods(2, 10)
array([1, 2], dtype=uint32)
>>> return_periods(100, 2)
array([ 50, 100], dtype=uint32)
>>> return_periods(1000, 1000)
array([   1,    2,    5,   10,   20,   50,  100,  200,  500, 1000],
      dtype=uint32)
openquake.risklib.scientific.scenario_damage(fragility_functions, gmvs)[source]
Parameters:
  • fragility_functions – a list of D - 1 fragility functions
  • gmvs – an array of E ground motion values
Returns:

an array of (D, E) damage fractions

openquake.risklib.utils module

openquake.risklib.utils.memoized(func)[source]

A simple memoize implementation. It works by adding a .cache dictionary to the decorated function. The cache will grow indefinitely, so it is your responsibility to clear it, if needed.

openquake.risklib.utils.numpy_map(f, *args)[source]
openquake.risklib.utils.pairwise(iterable)[source]

s -> (s0,s1), (s1,s2), (s2, s3), …

Module contents