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, bins, mon=<class 'openquake.baselib.performance.Monitor'>)[source]¶
- Parameters: - bdata – a dictionary of probabilities of no exceedence
- bins – bin edges
- mon – a Monitor instance
 - Returns: - a dictionary key -> matrix|pmf for each key in bdata 
- 
openquake.hazardlib.calc.disagg.collect_bin_data(rupdata, sitecol, cmaker, iml3, truncation_level, n_epsilons, monitor=<Monitor [jenkins]>)[source]¶
- Parameters: - rupdata – array of ruptures
- sitecol – a SiteCollection instance with a single site
- cmaker – a ContextMaker instance
- iml3 – an ArrayWrapper of intensities of shape (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.disaggregate(cmaker, sitecol, rupdata, iml2, truncnorm, epsilons, monitor=<Monitor [jenkins]>)[source]¶
- Disaggregate (separate) PoE in different contributions. - Parameters: - cmaker – a ContextMaker instance
- sitecol – a SiteCollection with N=1 site
- ruptures – an iterator over ruptures with the same TRT
- iml2 – a 2D array of IMLs of shape (M, P)
- truncnorm – an instance of scipy.stats.truncnorm
- epsilons – the epsilon bins
- monitor – a Monitor instance
 - Returns: - an AccumDict with keys (poe, imt, rlzi) and mags, dists, lons, lats 
- 
openquake.hazardlib.calc.disagg.disaggregate_pne(gsim, rupture, sctx, dctx, imt, iml, truncnorm, epsilons, eps_bands)[source]¶
- Disaggregate (separate) PoE of - imlin different contributions each coming from- epsilonsdistribution bins.- Other parameters are the same as for gsim.get_poes, with differences that - truncation_levelis required to be positive.- Returns: - Contribution to probability of exceedance of - imlcoming from different sigma bands in the form of a 2D numpy array of probabilities with shape (n_sites, n_epsilons)
- 
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>, filter_distance='rjb')[source]¶
- Compute “Disaggregation” matrix representing conditional probability of an intensity mesaure type - imtexceeding, at least once, an intensity measure level- imlat 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
PSHAcalculator it should be an iterator of seismic sources.
- site – Siteof interest to calculate disaggregation matrix for.
- imt – Instance of intensity measure typeclass.
- 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.get_bins(bin_edges, sid)[source]¶
- Returns: - mags, dists, lons, lats, eps for the given sid 
- 
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.lon_lat_trt_pmf(matrices)[source]¶
- Fold full disaggregation matrices to lon / lat / TRT PMF. - Parameters: - matrices – a matrix with T submatrices - Returns: - 3d array. First dimension represents longitude histogram bins, second one latitude histogram bins, third one trt 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. 
Filters (filters)¶
- 
class openquake.hazardlib.calc.filters.IntegrationDistance(dic)[source]¶
- Bases: - collections.abc.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 - 
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_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 
 
- 
- 
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, filename=None)[source]¶
- Bases: - object- Filter objects have a .filter method yielding filtered sources, i.e. sources with an attribute .indices, containg the IDs of the sites within the given maximum distance. There is also a .new method that filters the sources in parallel and returns a dictionary src_group_id -> filtered sources. Filter the sources by using self.sitecol.within_bbox which is based on numpy. - 
close_sids(rec, trt, mag)[source]¶
- Parameters: - rec – a record with fields minlon, minlat, maxlon, maxlat
- trt – tectonic region type string
- mag – magnitude
 - Returns: - the site indices within the bounding box enlarged by the integration distance for the given TRT and magnitude 
 - 
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 
 - 
sitecol¶
- Read the site collection from .filename and cache it 
 
- 
- 
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. 
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, imts, cmaker, truncation_level=None, correlation_model=None)[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.SiteCollectionsitecol:
- a complete SiteCollection
 - Parameters: - imts – a sorted list of Intensity Measure Type strings
- cmaker – a openquake.hazardlib.gsim.base.ContextMakerinstance
- 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 beNone, in which case non-correlated ground motion fields are calculated. Correlation model is not used iftruncation_levelis zero.
 - 
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) and two arrays with shape (num_imts, num_events): sig for stddev_inter and eps for the random part 
 
- :param 
- 
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
GMPEorIPE.
- 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 beNone, in which case non-correlated ground motion fields are calculated. Correlation model is not used iftruncation_levelis 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, srcfilter=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, srcfilter, imtls, gsim_by_trt, truncation_level=None, apply=<function sequential_apply>, filter_distance='rjb', reqv=None)[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) toGMPEorIPEobjects.
- truncation_level – Float, number of standard deviations for truncation of the intensity distribution.
- apply – apply function to use (default sequential_apply)
- filter_distance – The distance used to filter the ruptures (default rjb)
- 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.
- groups – A sequence of groups of seismic sources objects (instances of
of 
- 
openquake.hazardlib.calc.hazard_curve.classical(group, src_filter, gsims, param, monitor=<Monitor [jenkins]>)[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: pmap} with attributes .grp_ids, .calc_times, .eff_ruptures 
Stochastic Event Set (stochastic)¶
openquake.hazardlib.calc.stochastic contains
stochastic_event_set().
- 
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, calc_times, eff_ruptures 
- 
openquake.hazardlib.calc.stochastic.sample_ruptures(sources, srcfilter, param, monitor=<Monitor [jenkins]>)[source]¶
- Parameters: - sources – a sequence of sources of the same group
- srcfilter – SourceFilter instance used also for bounding box post filtering
- param – a dictionary of additional parameters including ses_per_logic_tree_path
- monitor – monitor instance
 - Yields: - dictionaries with keys rup_array, calc_times, eff_ruptures 
- 
openquake.hazardlib.calc.stochastic.stochastic_event_set(sources, 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).
- source_site_filter – The source filter to use (default noop filter)
 - Returns: - Generator of - Ruptureobjects that are contained in an event set. Some ruptures can be missing from it, others can appear one or more times in a row.
- sources – An iterator of seismic sources objects (instances of subclasses
of