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().

openquake.hazardlib.calc.disagg.build_disagg_matrix(bdata, bin_edges, sid, mon=<class 'openquake.baselib.performance.Monitor'>)[source]
Parameters:
  • bdata – a dictionary of probabilities of no exceedence
  • bin_edges – bin edges for each site
  • sid – site index
  • mon – a Monitor instance
Returns:

a dictionary key -> matrix|pmf for each key in bdata

openquake.hazardlib.calc.disagg.collect_bin_data(sources, sitecol, cmaker, iml4, truncation_level, n_epsilons, monitor=<Monitor dummy>)[source]
Parameters:
  • sources – a list of sources
  • sitecol – a SiteCollection instance
  • cmaker – a ContextMaker instance
  • iml4 – an ArrayWrapper of intensities of shape (N, R, M, P)
  • truncation_level – the truncation level
  • n_epsilons – the number of epsilons
  • monitor – a Monitor instance
Returns:

a dictionary (poe, imt, rlzi) -> probabilities of shape (N, E)

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>)[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.dist_pmf(matrix)[source]

Fold full disaggregation matrix to distance PMF.

Returns:1d array, a histogram representing distance PMF.
openquake.hazardlib.calc.disagg.lon_lat_bins(bb, coord_bin_width)[source]

Define bin edges for disaggregation histograms.

Given bins data as provided by collect_bin_data(), this function finds edges of histograms, taking into account maximum and minimum values of magnitude, distance and coordinates as well as requested sizes/numbers of bins.

openquake.hazardlib.calc.disagg.lon_lat_pmf(matrix)[source]

Fold full disaggregation matrix to longitude / latitude PMF.

Returns:2d array. First dimension represents longitude histogram bins, second one – latitude histogram bins.
openquake.hazardlib.calc.disagg.mag_dist_eps_pmf(matrix)[source]

Fold full disaggregation matrix to magnitude / distance / epsilon PMF.

Returns:3d array. First dimension represents magnitude histogram bins, second one – distance histogram bins, third one – epsilon histogram bins.
openquake.hazardlib.calc.disagg.mag_dist_pmf(matrix)[source]

Fold full disaggregation matrix to magnitude / distance PMF.

Returns:2d array. First dimension represents magnitude histogram bins, second one – distance histogram bins.
openquake.hazardlib.calc.disagg.mag_lon_lat_pmf(matrix)[source]

Fold full disaggregation matrix to magnitude / longitude / latitude PMF.

Returns:3d array. First dimension represents magnitude histogram bins, second one – longitude histogram bins, third one – latitude histogram bins.
openquake.hazardlib.calc.disagg.mag_pmf(matrix)[source]

Fold full disaggregation matrix to magnitude PMF.

Returns:1d array, a histogram representing magnitude PMF.
openquake.hazardlib.calc.disagg.make_iml4(R, imtls, iml_disagg, poes_disagg=(None, ), curves=())[source]
Returns:an ArrayWrapper over a 4D array of shape (N, R, M, P)
openquake.hazardlib.calc.disagg.trt_pmf(matrix)[source]

Fold full disaggregation matrix to tectonic region type PMF.

Returns:a scalar

Filters (filters)

Module filters contain filter functions for calculators.

Filters are functions (or other callable objects) that should take generators and return generators. There are two different kinds of filter functions:

  1. Source-site filters. Those functions take a generator of two-item tuples, each pair consists of seismic source object (that is, an instance of a subclass of BaseSeismicSource) and a site collection (instance of SiteCollection).
  2. Rupture-site filters. Those also take a generator of pairs, but in this case the first item in the pair is a rupture object (instance of Rupture). The second element in generator items is still site collection.

The purpose of both kinds of filters is to limit the amount of calculation to be done based on some criteria, like the distance between the source and the site. So common design feature of all the filters is the loop over pairs of the provided generator, filtering the sites collection, and if there are no items left in it, skipping the pair and continuing to the next one. If some sites need to be considered together with that source / rupture, the pair gets generated out, with a (possibly) limited site collection.

Consistency of filters’ input and output stream format allows several filters (obviously, of the same kind) to be chained together.

Filter functions should not make assumptions about the ordering of items in the original generator or draw more than one pair at once. Ideally, they should also perform reasonably fast (filtering stage that takes longer than the actual calculation on unfiltered collection only decreases performance).

Module openquake.hazardlib.calc.filters exports one distance-based filter function (see filter_sites_by_distance_to_rupture()) as well as a “no operation” filter (source_site_noop_filter). There is a class SourceFilter to determine the sites affected by a given source: the second one uses an R-tree index and it is faster if there are a lot of sources, i.e. if the initial time to prepare the index can be compensed. Finally, there is a function filter_sites_by_distance_to_rupture based on the Joyner-Boore distance.

exception openquake.hazardlib.calc.filters.FarAwayRupture[source]

Bases: exceptions.Exception

Raised if the rupture is outside the maximum distance for all sites

class openquake.hazardlib.calc.filters.IntegrationDistance(dic)[source]

Bases: _abcoll.Mapping

Pickleable object wrapping a dictionary of integration distances per tectonic region type. The integration distances can be scalars or list of pairs (magnitude, distance). Here is an example using ‘default’ as tectonic region type, so that the same values will be used for all tectonic region types:

>>> maxdist = IntegrationDistance({'default': [
...          (3, 30), (4, 40), (5, 100), (6, 200), (7, 300), (8, 400)]})
>>> maxdist('Some TRT', mag=2.5)
30
>>> maxdist('Some TRT', mag=3)
30
>>> maxdist('Some TRT', mag=3.1)
40
>>> maxdist('Some TRT', mag=8)
400
>>> maxdist('Some TRT', mag=8.5)  # 2000 km are used above the maximum
2000

It has also a method .get_closest(sites, rupture) returning the closest sites to the rupture and their distances. The integration distance can be missing if the sites have been already filtered (empty dictionary): in that case the method returns all the sites and all the distances.

get_bounding_box(lon, lat, trt=None, mag=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
  • mag – magnitude, possibly None
Returns:

min_lon, min_lat, max_lon, max_lat

get_closest(sites, rupture, distance_type='rrup', filter=True)[source]
Parameters:
  • sites – a (Filtered)SiteColletion
  • rupture – a rupture
  • distance_type – default ‘rrup’
Returns:

(close sites, close distances)

Raises:

a FarAwayRupture exception if the rupture is far away

class openquake.hazardlib.calc.filters.Piecewise(x, y)[source]

Bases: object

Given two arrays x and y of non-decreasing values, build a piecewise function associating to each x the corresponding y. If x is smaller then the minimum x, the minimum y is returned; if x is larger than the maximum x, the maximum y is returned.

class openquake.hazardlib.calc.filters.SourceFilter(sitecol, integration_distance, use_rtree=True)[source]

Bases: object

The SourceFilter uses the rtree library if available. The index is generated at instantiation time and kept in memory. The filter should be instantiated only once per calculation, after the site collection is known. It should be used as follows:

ss_filter = SourceFilter(sitecol, integration_distance)
for src, sites in ss_filter(sources):
   do_something(...)

As a side effect, sets the .nsites attribute of the source, i.e. the number of sites within the integration distance. Notice that SourceFilter instances can be pickled, but when unpickled the use_rtree flag is set to false and the index is lost: the reason is that libspatialindex indices cannot be properly pickled (https://github.com/Toblerity/rtree/issues/65).

Parameters:
get_affected_box(src)[source]

Get the enlarged bounding box of a source.

Parameters:src – a source object
Returns:a bounding box (min_lon, min_lat, max_lon, max_lat)
get_bounding_boxes(trt=None, mag=None)[source]
Parameters:
  • trt – a tectonic region type (used for the integration distance)
  • mag – a magnitude (used for the integration distance)
Returns:

a list of bounding boxes, one per site

get_close_sites(source)[source]

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

get_rectangle(src)[source]
Parameters:src – a source object
Returns:((min_lon, min_lat), width, height), useful for plotting
openquake.hazardlib.calc.filters.angular_distance(km, lat)[source]

Return the angular distance of two points at the given latitude.

openquake.hazardlib.calc.filters.context(*args, **kwds)[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.filter_sites_by_distance_to_rupture(rupture, integration_distance, sites)[source]

Filter out sites from the collection that are further from the rupture than some arbitrary threshold.

Parameters:
  • rupture – Instance of Rupture that was generated by :meth: openquake.hazardlib.source.base.BaseSeismicSource.iter_ruptures of an instance of this class.
  • integration_distance – Threshold distance in km.
  • sites – Instance of openquake.hazardlib.site.SiteCollection to filter.
Returns:

Filtered SiteCollection.

This function is similar to openquake.hazardlib.source.base.BaseSeismicSource.filter_sites_by_distance_to_source(). The same notes about filtering criteria apply. Site should not be filtered out if it is not further than the integration distance from the rupture’s surface projection along the great circle arc (this is known as Joyner-Boore distance, :meth:` openquake.hazardlib.geo.surface.base.BaseQuadrilateralSurface.get_joyner_boore_distance`).

openquake.hazardlib.calc.filters.get_distances(rupture, mesh, param)[source]
Parameters:
  • rupture – a rupture
  • mesh – a mesh of points
  • param – the kind of distance to compute (default rjb)
Returns:

an array of distances from the given mesh

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’

Ground Motion Fields (gmf)

Module gmf exports ground_motion_fields().

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

Bases: exceptions.Exception

class openquake.hazardlib.calc.gmf.GmfComputer(rupture, sites, imts, gsims, truncation_level=None, correlation_model=None, maximum_distance={})[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 sites:
Sites of interest to calculate GMFs.
Parameters:
  • imts – a sorted list of Intensity Measure Type strings
  • gsims – a set of GSIM instances
  • 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.
  • maximum_distance – Instance of openquake.hazardlib.calc.filters.IntegrationDistance By default there is no maximum_distance.
compute(gsim, num_events, seed=None)[source]
Parameters:
  • gsim – a GSIM instance
  • num_events – the number of seismic events
  • seed – a random seed or None
Returns:

a 32 bit array of shape (num_imts, num_sites, num_events)

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.

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
import logging
from openquake.baselib import parallel
from openquake.hazardlib.calc.filters import SourceFilter
from openquake.hazardlib.calc.hazard_curve import calc_hazard_curves
from openquake.commonlib import readinput

def main(job_ini):
    logging.basicConfig(level=logging.INFO)
    oq = readinput.get_oqparam(job_ini)
    sitecol = readinput.get_site_collection(oq)
    src_filter = SourceFilter(sitecol, oq.maximum_distance)
    csm = readinput.get_composite_source_model(oq).filter(src_filter)
    rlzs_assoc = csm.info.get_rlzs_assoc()
    for i, sm in enumerate(csm.source_models):
        for rlz in rlzs_assoc.rlzs_by_smodel[i]:
            gsim_by_trt = rlzs_assoc.gsim_by_trt[rlz.ordinal]
            hcurves = calc_hazard_curves(
                sm.src_groups, src_filter, oq.imtls,
                gsim_by_trt, oq.truncation_level,
                parallel.Starmap.apply)
        print('rlz=%s, hcurves=%s' % (rlz, hcurves))

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_curves(groups, ss_filter, imtls, gsim_by_trt, truncation_level=None, apply=<bound method type.apply of <class 'openquake.baselib.parallel.Sequential'>>)[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).
  • ss_filter – 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.
  • maximum_distance – The integration distance, if any
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.get_probability_no_exceedance(rupture, sctx, rctx, dctx, imtls, gsims, trunclevel, pne_mons)[source]
Parameters:
  • rupture – a Rupture instance
  • sctx – the corresponding SiteContext instance
  • rctx – the corresponding RuptureContext instance
  • dctx – the corresponding DistanceContext instance
  • imtls – a dictionary-like object providing the intensity levels
  • gsims – the list of GSIMs to use
  • trunclevel – the truncation level
  • pne_mons – monitors for the probability of no exceedance
Returns:

an array of shape (num_sites, num_levels, num_gsims)

openquake.hazardlib.calc.hazard_curve.pmap_from_grp(group, src_filter, gsims, param, monitor=<Monitor dummy>)[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 {grp_id: ProbabilityMap instance}
openquake.hazardlib.calc.hazard_curve.pmap_from_trt(sources, src_filter, gsims, param, monitor=<Monitor dummy>)[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.

Returns:a dictionary {grp_id: pmap} with attributes .grp_ids, .calc_times, .eff_ruptures
openquake.hazardlib.calc.hazard_curve.poe_map(src, s_sites, imtls, cmaker, trunclevel, ctx_mon, pne_mons, rup_indep=True)[source]

Compute the ProbabilityMap generated by the given source. Also, store some information in the monitors and optionally in the bounding boxes.

openquake.hazardlib.calc.hazard_curve.rupture_weight_pairs(src)[source]

Generator yielding (rupture, weight) for each rupture in the source

Stochastic Event Set (stochastic)

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

openquake.hazardlib.calc.stochastic.stochastic_event_set(sources, sites=None, source_site_filter=<openquake.hazardlib.calc.filters.SourceFilter object>)[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).
  • sites – A list of sites to consider (or None)
  • source_site_filter – The source filter to use (only meaningful is sites is not None)
  • source_site_filter – The rupture filter to use (only meaningful is sites is not None)
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.