openquake.hazardlib.calc package

Hazardlib calculators

Disaggregation (disagg)

openquake.hazardlib.calc.disagg contains disaggregation() as well as several aggregation functions for extracting a specific PMF from the result of disaggregation().

class openquake.hazardlib.calc.disagg.BinData(dists, lons, lats, pnes)

Bases: tuple

dists

Alias for field number 0

lats

Alias for field number 2

lons

Alias for field number 1

pnes

Alias for field number 3

openquake.hazardlib.calc.disagg.assert_same_shape(arrays)[source]

Raises an AssertionError if the shapes are not consistent

openquake.hazardlib.calc.disagg.disaggregate(ctxs, tom, g_by_z, iml2dict, eps3, sid=0, bin_edges=())[source]
Parameters
  • ctxs – a list of U RuptureContexts

  • tom – a temporal occurrence model

  • g_by_z – an array of gsim indices

  • iml2dict – a dictionary of arrays imt -> (P, Z)

  • eps3 – a triplet (truncnorm, epsilons, eps_bands)

openquake.hazardlib.calc.disagg.disaggregation(sources, site, imt, iml, gsim_by_trt, truncation_level, n_epsilons, mag_bin_width, dist_bin_width, coord_bin_width, source_filter=<openquake.hazardlib.calc.filters.SourceFilter object>, **kwargs)[source]

Compute “Disaggregation” matrix representing conditional probability of an intensity mesaure type imt exceeding, at least once, an intensity measure level iml at a geographical location site, given rupture scenarios classified in terms of:

  • rupture magnitude

  • Joyner-Boore distance from rupture surface to site

  • longitude and latitude of the surface projection of a rupture’s point closest to site

  • epsilon: number of standard deviations by which an intensity measure level deviates from the median value predicted by a GSIM, given the rupture parameters

  • rupture tectonic region type

In other words, the disaggregation matrix allows to compute the probability of each scenario with the specified properties (e.g., magnitude, or the magnitude and distance) to cause one or more exceedences of a given hazard level.

For more detailed information about the disaggregation, see for instance “Disaggregation of Seismic Hazard”, Paolo Bazzurro, C. Allin Cornell, Bulletin of the Seismological Society of America, Vol. 89, pp. 501-520, April 1999.

Parameters
  • sources – Seismic source model, as for PSHA calculator it should be an iterator of seismic sources.

  • siteSite of interest to calculate disaggregation matrix for.

  • imt – Instance of intensity measure type class.

  • iml – Intensity measure level. A float value in units of imt.

  • gsim_by_trt – Tectonic region type to GSIM objects mapping.

  • truncation_level – Float, number of standard deviations for truncation of the intensity distribution.

  • n_epsilons – Integer number of epsilon histogram bins in the result matrix.

  • mag_bin_width – Magnitude discretization step, width of one magnitude histogram bin.

  • dist_bin_width – Distance histogram discretization step, in km.

  • coord_bin_width – Longitude and latitude histograms discretization step, in decimal degrees.

  • source_filter – Optional source-site filter function. See openquake.hazardlib.calc.filters.

Returns

A tuple of two items. First is itself a tuple of bin edges information for (in specified order) magnitude, distance, longitude, latitude, epsilon and tectonic region types.

Second item is 6d-array representing the full disaggregation matrix. Dimensions are in the same order as bin edges in the first item of the result tuple. The matrix can be used directly by pmf-extractor functions.

openquake.hazardlib.calc.disagg.get_edges_shapedic(oq, sitecol, mags_by_trt)[source]
Returns

(mag dist lon lat eps trt) edges and shape dictionary

openquake.hazardlib.calc.disagg.lon_lat_bins(lon, lat, size_km, coord_bin_width)[source]

Define lon, lat bin edges for disaggregation histograms.

Parameters
  • lon – longitude of the site

  • lat – latitude of the site

  • size_km – total size of the bins in km

  • coord_bin_width – bin width in degrees

Returns

two arrays lon bins, lat bins

openquake.hazardlib.calc.disagg.lon_lat_trt_pmf(matrices)[source]

Fold full disaggregation matrices to lon / lat / TRT PMF.

Parameters

matrices – a matrix with T submatrices

Returns

4d array. First dimension represents longitude histogram bins, second one latitude histogram bins, third one trt histogram bins, last dimension is the z index, associatd to the realization.

openquake.hazardlib.calc.disagg.set_mean_std(ctxs, cmaker)[source]

Filters (filters)

class openquake.hazardlib.calc.filters.IntegrationDistance[source]

Bases: dict

A dictionary trt -> [(mag, dist), …]

get_bounding_box(lon, lat, trt=None)[source]

Build a bounding box around the given lon, lat by computing the maximum_distance at the given tectonic region type and magnitude.

Parameters
  • lon – longitude

  • lat – latitude

  • trt – tectonic region type, possibly None

Returns

min_lon, min_lat, max_lon, max_lat

get_dist_bins(trt, nbins=51)[source]
Returns

an array of distance bins, from 10m to maxdist

max()[source]
Returns

a dictionary trt -> maxdist

classmethod new(value)[source]
Parameters

value – string to be converted

Returns

IntegrationDistance dictionary

>>> md = IntegrationDistance.new('50')
>>> md
{'default': [(2.5, 50), (10.2, 50)]}
>>> md.max()
{'default': 50}
class openquake.hazardlib.calc.filters.SourceFilter(sitecol, integration_distance={'default': [(2.5, 1000), (10.2, 1000)]})[source]

Bases: object

Filter objects have a .filter method yielding filtered sources and the IDs of the sites within the given maximum distance. Filter the sources by using self.sitecol.within_bbox which is based on numpy.

close_sids(src_or_rec, trt=None, maxdist=None)[source]
Parameters
  • src_or_rec – a source or a rupture record

  • trt – passed only if src_or_rec is a rupture record

Returns

the site indices within the maximum_distance of the hypocenter, plus the maximum size of the bounding box

filter(sources)[source]
Parameters

sources – a sequence of sources

Yields

pairs (sources, sites)

get_close_sites(source)[source]

Returns the sites within the integration distance from the source, or None.

get_enlarged_box(src, maxdist=None)[source]

Get the enlarged bounding box of a source.

Parameters
  • src – a source object

  • maxdist – a scalar maximum distance (or None)

Returns

a bounding box (min_lon, min_lat, max_lon, max_lat)

get_rectangle(src)[source]
Parameters

src – a source object

Returns

((min_lon, min_lat), width, height), useful for plotting

reduce(factor=100)[source]

Reduce the SourceFilter to a subset of sites

split(sources)[source]
Yields

pairs (split, sites)

openquake.hazardlib.calc.filters.context(src)[source]

Used to add the source_id to the error message. To be used as

with context(src):

operation_with(src)

Typically the operation is filtering a source, that can fail for tricky geometries.

openquake.hazardlib.calc.filters.floatdict(value)[source]
Parameters

value – input string corresponding to a literal Python number or dictionary

Returns

a Python dictionary key -> number

>>> floatdict("200")
{'default': 200}
>>> floatdict("{'active shallow crust': 250., 'default': 200}")
{'active shallow crust': 250.0, 'default': 200}
openquake.hazardlib.calc.filters.get_distances(rupture, sites, param)[source]
Parameters
  • rupture – a rupture

  • sites – a mesh of points or a site collection

  • param – the kind of distance to compute (default rjb)

Returns

an array of distances from the given sites

openquake.hazardlib.calc.filters.getdefault(dic_with_default, key)[source]
Parameters
  • dic_with_default – a dictionary with a ‘default’ key

  • key – a key that may be present in the dictionary or not

Returns

the value associated to the key, or to ‘default’

openquake.hazardlib.calc.filters.magdepdist(pairs)[source]
Parameters

pairs – a list of pairs [(mag, dist), …]

Returns

a scipy.interpolate.interp1d function

openquake.hazardlib.calc.filters.magstr(mag)[source]
Returns

a string representation of the magnitude

openquake.hazardlib.calc.filters.split_source(src)[source]
Parameters

src – a splittable (or not splittable) source

Returns

the underlying sources (or the source itself)

openquake.hazardlib.calc.filters.unique_sorted(items)[source]

Check that the items are unique and sorted

Ground Motion Fields (gmf)

Module gmf exports ground_motion_fields().

exception openquake.hazardlib.calc.gmf.CorrelationButNoInterIntraStdDevs(corr, gsim)[source]

Bases: Exception

class openquake.hazardlib.calc.gmf.GmfComputer(rupture, sitecol, cmaker, correlation_model=None, cross_correl=None, amplifier=None, sec_perils=())[source]

Bases: object

Given an earthquake rupture, the ground motion field computer computes ground shaking over a set of sites, by randomly sampling a ground shaking intensity model.

Parameters

rupture – Rupture to calculate ground motion fields radiated from.

:param openquake.hazardlib.site.SiteCollection sitecol:

a complete SiteCollection

Parameters
  • imts – a sorted list of Intensity Measure Type strings

  • cmaker – a openquake.hazardlib.gsim.base.ContextMaker instance

  • truncation_level – Float, number of standard deviations for truncation of the intensity distribution, or None.

  • correlation_model – Instance of correlation model object. See openquake.hazardlib.correlation. Can be None, in which case non-correlated ground motion fields are calculated. Correlation model is not used if truncation_level is zero.

  • amplifier – None or an instance of Amplifier

compute(gsim, num_events, mean_stds)[source]
Parameters
  • gsim – GSIM used to compute mean_stds

  • num_events – the number of seismic events

  • mean_stds – array of shape (4, M, N)

Returns

a 32 bit array of shape (num_imts, num_sites, num_events) and two arrays with shape (num_imts, num_events): sig for tau and eps for the random part

compute_all(sig_eps=None)[source]
Returns

(dict with fields eid, sid, gmv_X, …), dt

openquake.hazardlib.calc.gmf.exp(vals, imt)[source]

Exponentiate the values unless the IMT is MMI

openquake.hazardlib.calc.gmf.ground_motion_fields(rupture, sites, imts, gsim, truncation_level, realizations, correlation_model=None, seed=None)[source]

Given an earthquake rupture, the ground motion field calculator computes ground shaking over a set of sites, by randomly sampling a ground shaking intensity model. A ground motion field represents a possible ‘realization’ of the ground shaking due to an earthquake rupture.

Note

This calculator is using random numbers. In order to reproduce the same results numpy random numbers generator needs to be seeded, see http://docs.scipy.org/doc/numpy/reference/generated/numpy.random.seed.html

Parameters
  • rupture (openquake.hazardlib.source.rupture.Rupture) – Rupture to calculate ground motion fields radiated from.

  • sites (openquake.hazardlib.site.SiteCollection) – Sites of interest to calculate GMFs.

  • imts – List of intensity measure type objects (see openquake.hazardlib.imt).

  • gsim – Ground-shaking intensity model, instance of subclass of either GMPE or IPE.

  • truncation_level – Float, number of standard deviations for truncation of the intensity distribution, or None.

  • realizations – Integer number of GMF realizations to compute.

  • correlation_model – Instance of correlation model object. See openquake.hazardlib.correlation. Can be None, in which case non-correlated ground motion fields are calculated. Correlation model is not used if truncation_level is zero.

  • seed (int) – The seed used in the numpy random number generator

Returns

Dictionary mapping intensity measure type objects (same as in parameter imts) to 2d numpy arrays of floats, representing different realizations of ground shaking intensity for all sites in the collection. First dimension represents sites and second one is for realizations.

openquake.hazardlib.calc.gmf.rvs(distribution, *size)[source]

Hazard Curves (hazard_curve)

openquake.hazardlib.calc.hazard_curve implements calc_hazard_curves(). Here is an example of a classical PSHA parallel calculator computing the hazard curves per each realization in less than 20 lines of code:

import sys
from openquake.commonlib import logs
from openquake.calculators.base import calculators

def main(job_ini):
    with logs.init('job', job_ini) as log:
        calc = calculators(log.get_oqparam(), log.calc_id)
        calc.run(individual_rlzs='true', shutdown=True)
        print('The hazard curves are in %s::/hcurves-rlzs'
             % calc.datastore.filename)

if __name__ == '__main__':
    main(sys.argv[1])  # path to a job.ini file

NB: the implementation in the engine is smarter and more efficient. Here we start a parallel computation per each realization, the engine manages all the realizations at once.

openquake.hazardlib.calc.hazard_curve.calc_hazard_curve(site1, src, gsims, oqparam, monitor=<Monitor [runner]>)[source]
Parameters
  • site1 – site collection with a single site

  • src – a seismic source object

  • gsims – a list of GSIM objects

  • oqparam – an object with attributes .maximum_distance, .imtls

  • monitor – a Monitor instance (optional)

Returns

a ProbabilityCurve object

openquake.hazardlib.calc.hazard_curve.calc_hazard_curves(groups, srcfilter, imtls, gsim_by_trt, truncation_level=None, apply=<function sequential_apply>, reqv=None, **kwargs)[source]

Compute hazard curves on a list of sites, given a set of seismic source groups and a dictionary of ground shaking intensity models (one per tectonic region type).

Probability of ground motion exceedance is computed in different ways depending if the sources are independent or mutually exclusive.

Parameters
  • groups – A sequence of groups of seismic sources objects (instances of of BaseSeismicSource).

  • srcfilter – A source filter over the site collection or the site collection itself

  • imtls – Dictionary mapping intensity measure type strings to lists of intensity measure levels.

  • gsim_by_trt – Dictionary mapping tectonic region types (members of openquake.hazardlib.const.TRT) to GMPE or IPE objects.

  • truncation_level – Float, number of standard deviations for truncation of the intensity distribution.

  • apply – apply function to use (default sequential_apply)

  • reqv – If not None, an instance of RjbEquivalent

Returns

An array of size N, where N is the number of sites, which elements are records with fields given by the intensity measure types; the size of each field is given by the number of levels in imtls.

openquake.hazardlib.calc.hazard_curve.classical(group, sitecol, cmaker)[source]

Compute the hazard curves for a set of sources belonging to the same tectonic region type for all the GSIMs associated to that TRT. The arguments are the same as in calc_hazard_curves(), except for gsims, which is a list of GSIM instances.

Returns

a dictionary with keys pmap, source_data, rup_data, extra

Stochastic Event Set (stochastic)

openquake.hazardlib.calc.stochastic contains stochastic_event_set().

openquake.hazardlib.calc.stochastic.get_ebr_df(ebruptures, cmakerdict)[source]
Parameters
  • ebruptures – the output of sample_ebruptures

  • rlzs_by_gsim_trt – a double dictionary trt -> gsim -> rlzs

Returns

a DataFrame with fields eid, rlz indexed by rupture ordinal

openquake.hazardlib.calc.stochastic.get_rup_array(ebruptures, srcfilter=<openquake.hazardlib.calc.filters.SourceFilter object>)[source]

Convert a list of EBRuptures into a numpy composite array, by filtering out the ruptures far away from every site

openquake.hazardlib.calc.stochastic.sample_cluster(sources, srcfilter, num_ses, param)[source]

Yields ruptures generated by a cluster of sources.

Parameters
  • sources – A sequence of sources of the same group

  • num_ses – Number of stochastic event sets

  • param – a dictionary of additional parameters including ses_per_logic_tree_path

Yields

dictionaries with keys rup_array, source_data, eff_ruptures

openquake.hazardlib.calc.stochastic.sample_ebruptures(src_groups, cmakerdict)[source]

Sample independent sources without filtering.

Parameters
  • src_groups – a list of source groups

  • cmakerdict – a dictionary TRT -> cmaker

Returns

a list of EBRuptures

openquake.hazardlib.calc.stochastic.sample_ruptures(sources, cmaker, sitecol=None, monitor=<Monitor [runner]>)[source]
Parameters
  • sources – a sequence of sources of the same group

  • cmaker – a ContextMaker instance with ses_per_logic_tree_path, ses_seed

  • sitecol – SiteCollection instance used for filtering (None for no filtering)

  • monitor – monitor instance

Yields

dictionaries with keys rup_array, source_data

openquake.hazardlib.calc.stochastic.stochastic_event_set(sources, source_site_filter=<openquake.hazardlib.calc.filters.SourceFilter object>, **kwargs)[source]

Generates a ‘Stochastic Event Set’ (that is a collection of earthquake ruptures) representing a possible realization of the seismicity as described by a source model.

The calculator loops over sources. For each source, it loops over ruptures. For each rupture, the number of occurrence is randomly sampled by calling openquake.hazardlib.source.rupture.BaseProbabilisticRupture.sample_number_of_occurrences()

Note

This calculator is using random numbers. In order to reproduce the same results numpy random numbers generator needs to be seeded, see http://docs.scipy.org/doc/numpy/reference/generated/numpy.random.seed.html

Parameters
  • sources – An iterator of seismic sources objects (instances of subclasses of BaseSeismicSource).

  • source_site_filter – The source filter to use (default noop filter)

Returns

Generator of Rupture objects that are contained in an event set. Some ruptures can be missing from it, others can appear one or more times in a row.