openquake.hazardlib package

Subpackages

const

Module openquake.hazardlib.const defines various constants.

class openquake.hazardlib.const.ConstantContainer[source]

Bases: object

Class that doesn’t support instantiation.

>>> ConstantContainer()
Traceback (most recent call last):
    ...
AssertionError: do not create objects ConstantContainer, use class properties instead
class openquake.hazardlib.const.IMC[source]

Bases: openquake.hazardlib.const.ConstantContainer

The intensity measure component is the component of interest of ground shaking for an intensity measure.

AVERAGE_HORIZONTAL = 'Average horizontal'

Usually defined as the geometric average of the maximum of the two horizontal components (which may not occur at the same time).

GMRotD100 = 'Average Horizontal (GMRotD100)'

The geometric mean of the records rotated into the most adverse direction for the structure.

GMRotI50 = 'Average Horizontal (GMRotI50)'

An orientation-independent alternative to AVERAGE_HORIZONTAL. Defined at Boore et al. (2006, Bull. Seism. Soc. Am. 96, 1502-1511) and is used for all the NGA GMPEs.

GREATER_OF_TWO_HORIZONTAL = 'Greater of two horizontal'

The largest value obtained from two perpendicular horizontal components.

HORIZONTAL = 'Horizontal'

The horizontal component.

MEDIAN_HORIZONTAL = 'Median horizontal'

The median horizontal component.

PEAK_SRSS_HORIZONTAL = 'Peak square root of sum of squares of horizontals'

“the peak square root of the sum of squares of two orthogonal horizontal components in the time domain” p. 880 of Kanno et al. (2006, Bull. Seism. Soc. Am. 96, 879-897)

RANDOM_HORIZONTAL = 'Random horizontal'

A randomly chosen horizontal component.

RotD100 = 'Horizontal Maximum Direction (RotD100)'
RotD50 = 'Average Horizontal (RotD50)'

An orientation-independent alternative to AVERAGE_HORIZONTAL. Defined at Boore et al. (2006, Bull. Seism. Soc. Am. 96, 1502-1511) and is used for all the NGA GMPEs.

VECTORIAL = 'Square root of sum of squares of peak horizontals'

“Vectorial addition: a_V = sqrt(max|a_1(t)|^2 + max|a_2(t)|^2)). This means that the maximum ground amplitudes occur simultaneously on the two horizontal components; this is a conservative assumption.” p. 53 of Douglas (2003, Earth-Sci. Rev. 61, 43-104)

VERTICAL = 'Vertical'

The vertical component.

class openquake.hazardlib.const.StdDev[source]

Bases: openquake.hazardlib.const.ConstantContainer

GSIM standard deviation represents ground shaking variability at a site.

INTER_EVENT = 'Inter event'

Standard deviation representing ground shaking variability within different events.

INTRA_EVENT = 'Intra event'

Standard deviation representing ground shaking variability within a single event.

TOTAL = 'Total'

Total standard deviation, defined as the square root of the sum of inter- and intra-event squared standard deviations, represents the total ground shaking variability, and is the only one that is used for calculating a probability of intensity exceedance (see openquake.hazardlib.gsim.base.GroundShakingIntensityModel.get_poes()).

class openquake.hazardlib.const.TRT[source]

Bases: openquake.hazardlib.const.ConstantContainer

Container for constants that define some of the common Tectonic Region Types.

ACTIVE_SHALLOW_CRUST = 'Active Shallow Crust'
GEOTHERMAL = 'Geothermal'
INDUCED = 'Induced'
STABLE_CONTINENTAL = 'Stable Shallow Crust'
SUBDUCTION_INTERFACE = 'Subduction Interface'
SUBDUCTION_INTRASLAB = 'Subduction IntraSlab'
VOLCANIC = 'Volcanic'

correlation

Module openquake.hazardlib.correlation defines correlation models for spatially-distributed ground-shaking intensities.

class openquake.hazardlib.correlation.BaseCorrelationModel[source]

Bases: object

Base class for correlation models for spatially-distributed ground-shaking intensities.

apply_correlation(sites, imt, residuals)[source]

Apply correlation to randomly sampled residuals.

Parameters:
  • sitesSiteCollection residuals were sampled for.
  • imt – Intensity measure type object, see openquake.hazardlib.imt.
  • residuals – 2d numpy array of sampled residuals, where first dimension represents sites (the length as sites parameter) and second one represents different realizations (samples).
Returns:

Array of the same structure and semantics as residuals but with correlations applied.

NB: the correlation matrix is cached. It is computed only once per IMT for the complete site collection and then the portion corresponding to the sites is multiplied by the residuals.

get_lower_triangle_correlation_matrix(sites, imt)[source]

Get lower-triangle matrix as a result of Cholesky-decomposition of correlation matrix.

The resulting matrix should have zeros on values above the main diagonal.

The actual implementations of BaseCorrelationModel interface might calculate the matrix considering site collection and IMT (like JB2009CorrelationModel does) or might have it pre-constructed for a specific site collection and IMT, in which case they will need to make sure that parameters to this function match parameters that were used to pre-calculate decomposed correlation matrix.

Parameters:
class openquake.hazardlib.correlation.JB2009CorrelationModel(vs30_clustering)[source]

Bases: openquake.hazardlib.correlation.BaseCorrelationModel

“Correlation model for spatially distributed ground-motion intensities” by Nirmal Jayaram and Jack W. Baker. Published in Earthquake Engineering and Structural Dynamics 2009; 38, pages 1687-1708.

Parameters:vs30_clustering – Boolean value to indicate whether “Case 1” or “Case 2” from page 1700 should be applied. True value means that Vs 30 values show or are expected to show clustering (“Case 2”), False means otherwise.
get_lower_triangle_correlation_matrix(sites, imt)[source]

See BaseCorrelationModel.get_lower_triangle_correlation_matrix().

imt

Module openquake.hazardlib.imt defines different intensity measure types.

class openquake.hazardlib.imt.PGA[source]

Bases: openquake.hazardlib.imt._IMT

Peak ground acceleration during an earthquake measured in units of g, times of gravitational acceleration.

class openquake.hazardlib.imt.PGV[source]

Bases: openquake.hazardlib.imt._IMT

Peak ground velocity during an earthquake measured in units of cm/sec.

class openquake.hazardlib.imt.PGD[source]

Bases: openquake.hazardlib.imt._IMT

Peak ground displacement during an earthquake measured in units of cm.

class openquake.hazardlib.imt.SA[source]

Bases: openquake.hazardlib.imt._IMT

Spectral acceleration, defined as the maximum acceleration of a damped, single-degree-of-freedom harmonic oscillator. Units are g, times of gravitational acceleration.

Parameters:
  • period – The natural period of the oscillator in seconds.
  • damping – The degree of damping for the oscillator in percents.
Raises:

ValueError – if period or damping is not positive.

damping

itemgetter(item, ...) –> itemgetter object

Return a callable object that fetches the given item(s) from its operand. After f = itemgetter(2), the call f(r) returns r[2]. After g = itemgetter(2, 5, 3), the call g(r) returns (r[2], r[5], r[3])

period

itemgetter(item, ...) –> itemgetter object

Return a callable object that fetches the given item(s) from its operand. After f = itemgetter(2), the call f(r) returns r[2]. After g = itemgetter(2, 5, 3), the call g(r) returns (r[2], r[5], r[3])

class openquake.hazardlib.imt.IA[source]

Bases: openquake.hazardlib.imt._IMT

Arias intensity. Determines the intensity of shaking by measuring the acceleration of transient seismic waves. Units are m/s.

class openquake.hazardlib.imt.CAV[source]

Bases: openquake.hazardlib.imt._IMT

Cumulative Absolute Velocity. Defins the integral of the absolute acceleration time series. Units are “g-sec”

class openquake.hazardlib.imt.RSD[source]

Bases: openquake.hazardlib.imt._IMT

Relative significant duration, 5-95% of Arias intensity, in seconds.

class openquake.hazardlib.imt.MMI[source]

Bases: openquake.hazardlib.imt._IMT

Modified Mercalli intensity, a Roman numeral describing the severity of an earthquake in terms of its effects on the earth’s surface and on humans and their structures.

near_fault

Module openquake.hazardlib.nearfault provides methods for near fault PSHA calculation.

openquake.hazardlib.near_fault.average_s_rad(site, hypocenter, reference, pp, normal, dist_to_plane, e, p0, p1, delta_slip)[source]

Gets the average S-wave radiation pattern given an e-path as described in: Spudich et al. (2013) “Final report of the NGA-West2 directivity working group”, PEER report, page 90- 92 and computes: the site to the direct point distance, rd, and the hypocentral distance, r_hyp.

Parameters:
  • sitePoint object representing the location of the target site
  • hypocenterPoint object representing the location of hypocenter
  • referencePoint object representing the location of the reference point for coordinate projection within the calculation. The suggested reference point is Epicentre.
  • pp – the projection point pp on the patch plane, a numpy array
  • normal – normal of the plane, describe by a normal vector[a, b, c]
  • dist_to_plane – d is the constant term in the plane equation, e.g., ax + by + cz = d
  • e – a float defining the E-path length, which is the distance from Pd(direction) point to hypocentre. In km.
  • p0Point object representing the location of the starting point on fault segment
  • p1Point object representing the location of the ending point on fault segment.
  • delta_slip – slip direction away from the strike direction, in decimal degrees. A positive angle is generated by a counter-clockwise rotation.
Returns:

fs, float value of the average S-wave radiation pattern. rd, float value of the distance from site to the direct point. r_hyp, float value of the hypocetre distance.

openquake.hazardlib.near_fault.directp(node0, node1, node2, node3, hypocenter, reference, pp)[source]

Get the Direct Point and the corresponding E-path as described in Spudich et al. (2013). This method also provides a logical variable stating if the DPP calculation must consider the neighbouring patch. To define the intersection point(Pd) of PpPh line segment and fault plane, we obtain the intersection points(Pd) with each side of fault plan, and check which intersection point(Pd) is the one fitting the definition in the Chiou and Spudich(2014) directivity model. Two possible locations for Pd, the first case, Pd locates on the side of the fault patch when Pp is not inside the fault patch. The second case is when Pp is inside the fault patch, then Pd=Pp.

For the first case, it follows three conditions: 1. the PpPh and PdPh line vector are the same, 2. PpPh >= PdPh, 3. Pd is not inside the fault patch.

If we can not find solution for all the four possible intersection points for the first case, we check if the intersection point fit the second case by checking if Pp is inside the fault patch.

Because of the coordinate system mapping(from geographic system to Catestian system), we allow an error when we check the location. The allow error will keep increasing after each loop when no solution in the two cases are found, until the solution get obtained.

Parameters:
  • node0Point object representing the location of one vertices on the target fault segment.
  • node1Point object representing the location of one vertices on the target fault segment. Note, the order should be clockwise.
  • node2Point object representing the location of one vertices on the target fault segment. Note, the order should be clockwise.
  • node3Point object representing the location of one vertices on the target fault segment. Note, the order should be clockwise.
  • hypocenterPoint object representing the location of floating hypocenter on each segment calculation. In the method, we take the direction point of the previous fault patch as hypocentre for the current fault patch.
  • referencePoint object representing the location of reference point for projection
  • pp – the projection of the site onto the plane containing the fault slipped area. A numpy array.
Returns:

Pd, a numpy array, representing the location of direction point E, the distance from direction point to hypocentre. go_next_patch, flag indicates if the calculation goes on the next fault patch. 1: yes, 0: no.

openquake.hazardlib.near_fault.get_plane_equation(p0, p1, p2, reference)[source]

Define the equation of target fault plane passing through 3 given points which includes two points on the fault trace and one point on the fault plane but away from the fault trace. Note: in order to remain the consistency of the fault normal vector direction definition, the order of the three given points is strickly defined.

Parameters:
  • p0 – The fault trace and is the closer points from the starting point of fault trace. Point object representing the location of the one vertex of the fault patch.
  • p1 – The fault trace and is the further points from the starting point of fault trace. Point object representing the location of the one vertex of the fault patch.
  • p2 – The point on the fault plane but away from the fault trace. Point object representing the location of the one vertex of the fault patch.
  • referencePoint object representing the origin of the cartesian system used the represent objects in a projected reference
Returns:

normal: normal vector of the plane (a,b,c) dist_to_plane: d in the plane equation, ax + by + cz = d

openquake.hazardlib.near_fault.get_xyz_from_ll(projected, reference)[source]

This method computes the x, y and z coordinates of a set of points provided a reference point

Parameters:
  • projectedPoint object representing the coordinates of target point to be projected
  • referencePoint object representing the coordinates of the reference point.
Returns:

x y z

openquake.hazardlib.near_fault.isochone_ratio(e, rd, r_hyp)[source]

Get the isochone ratio as described in Spudich et al. (2013) PEER report, page 88.

Parameters:
  • e – a float defining the E-path length, which is the distance from Pd(direction) point to hypocentre. In km.
  • rd – float, distance from the site to the direct point.
  • r_hyp – float, the hypocentre distance.
Returns:

c_prime, a float defining the isochone ratio

openquake.hazardlib.near_fault.projection_pp(site, normal, dist_to_plane, reference)[source]

This method finds the projection of the site onto the plane containing the slipped area, defined as the Pp(i.e. ‘perpendicular projection of site location onto the fault plane’ Spudich et al. (2013) - page 88) given a site.

Parameters:
  • site – Location of the site, [lon, lat, dep]
  • normal – Normal to the plane including the fault patch, describe by a normal vector[a, b, c]
  • dist_to_plane – D in the plane equation, ax + by + cz = d
  • referencePoint object representing the location of project reference point
Returns:

pp, the projection point, [ppx, ppy, ppz], in xyz domain , a numpy array.

openquake.hazardlib.near_fault.vectors2angle(v1, v2)[source]

Returns the angle in radians between vectors ‘v1’ and ‘v2’.

Parameters:
  • v1 – vector, a numpy array
  • v2 – vector, a numpy array
Returns:

the angle in radians between the two vetors

nrml

It is possible to save a Node object into a NRML file by using the function write(nodes, output) where output is a file object. If you want to make sure that the generated file is valid according to the NRML schema just open it in ‘w+’ mode: immediately after writing it will be read and validated. It is also possible to convert a NRML file into a Node object with the routine read(node, input) where input is the path name of the NRML file or a file object opened for reading. The file will be validated as soon as opened.

For instance an exposure file like the following:

<?xml version='1.0' encoding='utf-8'?>
<nrml xmlns="http://openquake.org/xmlns/nrml/0.4"
      xmlns:gml="http://www.opengis.net/gml">
  <exposureModel
      id="my_exposure_model_for_population"
      category="population"
      taxonomySource="fake population datasource">

    <description>
      Sample population
    </description>

    <assets>
      <asset id="asset_01" number="7" taxonomy="IT-PV">
          <location lon="9.15000" lat="45.16667" />
      </asset>

      <asset id="asset_02" number="7" taxonomy="IT-CE">
          <location lon="9.15333" lat="45.12200" />
      </asset>
    </assets>
  </exposureModel>
</nrml>

can be converted as follows:

>> nrml = read(<path_to_the_exposure_file.xml>)

Then subnodes and attributes can be conveniently accessed:

>> nrml.exposureModel.assets[0][‘taxonomy’] ‘IT-PV’ >> nrml.exposureModel.assets[0][‘id’] ‘asset_01’ >> nrml.exposureModel.assets[0].location[‘lon’] ‘9.15000’ >> nrml.exposureModel.assets[0].location[‘lat’] ‘45.16667’

The Node class provides no facility to cast strings into Python types; this is a job for the Node class which can be subclassed and supplemented by a dictionary of validators.

exception openquake.hazardlib.nrml.DuplicatedID[source]

Bases: exceptions.Exception

Raised when two sources with the same ID are found in a source model

class openquake.hazardlib.nrml.SourceModelParser(converter)[source]

Bases: object

A source model parser featuring a cache.

Parameters:converteropenquake.commonlib.source.SourceConverter instance
parse_groups(fname)[source]

Parse all the groups and return them ordered by number of sources. It does not count the ruptures, so it is relatively fast.

Parameters:fname – the full pathname of the source model file
parse_src_groups(fname, apply_uncertainties=None)[source]
Parameters:
  • fname – the full pathname of the source model file
  • apply_uncertainties – a function modifying the sources (or None)
openquake.hazardlib.nrml.convert(node)[source]

Convert a node into a string in NRML format

openquake.hazardlib.nrml.get_rupture_collection(node, fname, converter)[source]
openquake.hazardlib.nrml.get_source_model_04(node, fname, converter)[source]
openquake.hazardlib.nrml.get_source_model_05(node, fname, converter)[source]
openquake.hazardlib.nrml.get_tag_version(nrml_node)[source]

Extract from a node of kind NRML the tag and the version. For instance from ‘{http://openquake.org/xmlns/nrml/0.4}fragilityModel’ one gets the pair (‘fragilityModel’, ‘nrml/0.4’).

openquake.hazardlib.nrml.parse(fname, *args)[source]

Parse a NRML file and return an associated Python object. It works by calling nrml.read() and node_to_obj() in sequence.

openquake.hazardlib.nrml.read(source, chatty=True, stop=None)[source]

Convert a NRML file into a validated Node object. Keeps the entire tree in memory.

Parameters:source – a file name or file object open for reading
openquake.hazardlib.nrml.write(nodes, output=<open file '<stdout>', mode 'w'>, fmt='%.7E', gml=True, xmlns=None)[source]

Convert nodes into a NRML file. output must be a file object open in write mode. If you want to perform a consistency check, open it in read-write mode, then it will be read after creation and validated.

Params nodes:

an iterable over Node objects

Params output:

a file-like object in write or read-write mode

Parameters:

pmf

Module openquake.hazardlib.pmf implements PMF.

class openquake.hazardlib.pmf.PMF(data, epsilon=1e-15)[source]

Bases: object

Probability mass function is a function that gives the probability that a discrete random variable is exactly equal to some value.

Parameters:
  • data

    The PMF data in a form of list of tuples. Each tuple must contain two items with the first item being the probability of the implied random variable to take value of the second item.

    There must be at least one tuple in the list. All probabilities must sum up to 1 within the given precision.

    The type of values (second items in tuples) is not strictly defined, those can be objects of any (mixed or homogeneous) type.

  • epsilon – the tolerance for the sum of the probabilities (default 1E-15)
Raises:

ValueError – If probabilities do not sum up to 1 or there is zero or negative probability.

NB: in the engine the sum is checked to be exactly one by using Decimal numbers.

assert_equal(other)
sample(number_samples)[source]

Produces a list of samples from the probability mass function.

Parameters:data (int) – Number of samples
Returns:Samples from PMF as a list
sample_pairs(number_samples)[source]

Produces a list of samples from the probability mass function.

Parameters:data (int) – Number of samples
Returns:Samples from PMF as a list of pairs

probability_map

exception openquake.hazardlib.probability_map.AllEmptyProbabilityMaps[source]

Bases: exceptions.ValueError

Raised by get_shape(pmaps) if all passed probability maps are empty

class openquake.hazardlib.probability_map.ProbabilityCurve(array)[source]

Bases: object

This class is a small wrapper over an array of PoEs associated to a set of intensity measure types and levels. It provides a few operators, including the complement operator ~

~p = 1 - p

and the inclusive or operator |

p = p1 | p2 = ~(~p1 * ~p2)

Such operators are implemented efficiently at the numpy level, by dispatching on the underlying array.

Here is an example of use:

>>> poe = ProbabilityCurve(numpy.array([0.1, 0.2, 0.3, 0, 0]))
>>> ~(poe | poe) * .5
<ProbabilityCurve
[ 0.405  0.32   0.245  0.5    0.5  ]>
convert(imtls, idx=0)[source]

Convert a probability curve into a record of dtype imtls.dt.

Parameters:
  • imtls – DictArray instance
  • idx – extract the data corresponding to the given inner index
class openquake.hazardlib.probability_map.ProbabilityMap(shape_y, shape_z=1)[source]

Bases: dict

A dictionary site_id -> ProbabilityCurve. It defines the complement operator ~, performing the complement on each curve

~p = 1 - p

and the “inclusive or” operator |:

m = m1 | m2 = {sid: m1[sid] | m2[sid] for sid in all_sids}

Such operators are implemented efficiently at the numpy level, by dispatching on the underlying array. Moreover there is a classmethod .build(L, I, sids, initvalue) to build initialized instances of ProbabilityMap. The map can be represented as 3D array of shape (shape_x, shape_y, shape_z) = (N, L, I), where N is the number of site IDs, L the total number of hazard levels and I the number of GSIMs.

array

The underlying array of shape (N, L, I)

classmethod build(shape_y, shape_z, sids, initvalue=0.0)[source]
Parameters:
  • shape_y – the total number of intensity measure levels
  • shape_z – the number of inner levels
  • sids – a set of site indices
  • initvalue – the initial value of the probability (default 0)
Returns:

a ProbabilityMap dictionary

convert(imtls, nsites, idx=0)[source]

Convert a probability map into a composite array of length nsites and dtype imtls.dt.

Parameters:
  • imtls – DictArray instance
  • nsites – the total number of sites
  • idx – index on the z-axis (default 0)
convert2(imtls, sids)[source]

Convert a probability map into a composite array of shape (N, Z) and dtype imtls.dt.

Parameters:
  • imtls – DictArray instance
  • sids – the IDs of the sites we are interested in
Returns:

an array of curves of shape (N, Z)

extract(inner_idx)[source]

Extracts a component of the underlying ProbabilityCurves, specified by the index inner_idx.

filter(sids)[source]

Extracs a submap of self for the given sids.

classmethod from_array(array, sids)[source]
Parameters:
  • array – array of shape (N, L, I)
  • sids – array of N site IDs
nbytes

The size of the underlying array

setdefault(sid, value)[source]

Works like dict.setdefault: if the sid key is missing, it fills it with an array and returns it.

Parameters:
  • sid – site ID
  • value – value used to fill the returned array
sids

The ordered keys of the map as a numpy.uint32 array

openquake.hazardlib.probability_map.get_shape(pmaps)[source]
Parameters:pmaps – a set of homogenous ProbabilityMaps
Returns:the common shape (N, L, I)

site

Module openquake.hazardlib.site defines Site.

class openquake.hazardlib.site.FilteredSiteCollection(indices, complete)[source]

Bases: object

A class meant to store proper subsets of a complete collection of sites in a memory-efficient way.

Parameters:
  • indices – an array of indices referring to the complete site collection
  • complete – the complete site collection the filtered collection was derived from

Notice that if you filter a FilteredSiteCollection fsc, you will get a different FilteredSiteCollection referring to the complete SiteCollection fsc.complete, not to the filtered collection fsc.

assert_equal(other)
at_sea_level()[source]

True if all depths are zero

backarc

Extract backarc array from FilteredSiteCollection

depths

Extract depths array from FilteredSiteCollection

filter(mask)[source]

Create a FilteredSiteCollection with only a subset of sites from this one.

Parameters:mask – Numpy array of boolean values of the same length as this filtered sites collection. True values should indicate that site with that index should be included into the filtered collection.
Returns:A new FilteredSiteCollection instance, unless all the values in mask are True, in which case this site collection is returned, or if all the values in mask are False, in which case method returns None. New collection has data of only those sites that were marked for inclusion in mask.
lats

Extract lats array from FilteredSiteCollection

lons

Extract lons array from FilteredSiteCollection

mesh

Return a mesh with the given lons, lats, and depths

sids

Extract sids array from FilteredSiteCollection

total_sites

The total number of the original sites, without filtering

vs30

Extract vs30 array from FilteredSiteCollection

vs30measured

Extract vs30measured array from FilteredSiteCollection

z1pt0

Extract z1pt0 array from FilteredSiteCollection

z2pt5

Extract z2pt5 array from FilteredSiteCollection

class openquake.hazardlib.site.Site(location, vs30, vs30measured, z1pt0, z2pt5, backarc=False)[source]

Bases: object

Site object represents a geographical location defined by its position as well as its soil characteristics.

Parameters:
  • location – Instance of Point representing where the site is located.
  • vs30 – Average shear wave velocity in the top 30 m, in m/s.
  • vs30measured – Boolean value, True if vs30 was measured on that location and False if it was inferred.
  • z1pt0 – Vertical distance from earth surface to the layer where seismic waves start to propagate with a speed above 1.0 km/sec, in meters.
  • z2pt5 – Vertical distance from earth surface to the layer where seismic waves start to propagate with a speed above 2.5 km/sec, in km.
  • backarc" – Boolean value, True if the site is in the subduction backarc and False if it is in the subduction forearc or is unknown
Raises:

ValueError – If any of vs30, z1pt0 or z2pt5 is zero or negative.

Note

Sites are pickleable

assert_equal(other)
class openquake.hazardlib.site.SiteCollection(sites)[source]

Bases: object

A collection of sites.

Instances of this class are intended to represent a large collection of sites in a most efficient way in terms of memory usage.

Note

If a SiteCollection is created from sites containing only lon and lat, iterating over the collection will yield Sites with a reference depth of 0.0 (the sea level). Otherwise, it is possible to model the sites on a realistic topographic surface by specifying the depth of each site.

Parameters:sites – A list of instances of Site class.
assert_equal(other)
at_sea_level()[source]

True if all depths are zero

backarc

backarc array

dtype = dtype([('sids', '<u4'), ('lons', '<f8'), ('lats', '<f8'), ('depths', '<f8'), ('_vs30', '<f8'), ('_vs30measured', '?'), ('_z1pt0', '<f8'), ('_z2pt5', '<f8'), ('_backarc', '?')])
filter(mask)[source]

Create a FilteredSiteCollection with only a subset of sites from this one.

Parameters:mask – Numpy array of boolean values of the same length as this sites collection. True values should indicate that site with that index should be included into the filtered collection.
Returns:A new FilteredSiteCollection instance, unless all the values in mask are True, in which case this site collection is returned, or if all the values in mask are False, in which case method returns None. New collection has data of only those sites that were marked for inclusion in mask.
classmethod from_points(lons, lats, depths, sitemodel)[source]

Build the site collection from

Parameters:
  • lons – a sequence of longitudes
  • lats – a sequence of latitudes
  • depths – a sequence of depths
  • sitemodel – an object containing the attributes reference_vs30_value, reference_vs30_type, reference_depth_to_1pt0km_per_sec, reference_depth_to_2pt5km_per_sec, reference_backarc
indices

The full set of indices from 0 to total_sites - 1

mesh

Return a mesh with the given lons, lats, and depths

split_in_tiles(hint)[source]

Split a SiteCollection into a set of tiles (SiteCollection instances).

Parameters:hint – hint for how many tiles to generate
vs30

vs30 array

vs30measured

vs30measured array

z1pt0

z1pt0 array

z2pt5

z2pt5 array

openquake.hazardlib.site.getarray(sc, name='backarc')[source]
openquake.hazardlib.site.prop

Extract sids array from FilteredSiteCollection

sourceconverter

class openquake.hazardlib.sourceconverter.RuptureConverter(rupture_mesh_spacing, complex_fault_mesh_spacing=None)[source]

Bases: object

Convert ruptures from nodes into Hazardlib ruptures.

convert_complexFaultRupture(node)[source]

Convert a complexFaultRupture node.

Parameters:node – the rupture node
convert_griddedRupture(node)[source]

Convert a griddedRupture node.

Parameters:node – the rupture node
convert_multiPlanesRupture(node)[source]

Convert a multiPlanesRupture node.

Parameters:node – the rupture node
convert_node(node)[source]

Convert the given rupture node into a hazardlib rupture, depending on the node tag.

Parameters:node – a node representing a rupture
convert_ruptureCollection(node)[source]
Parameters:node – a ruptureCollection node
Returns:a dictionary grp_id -> EBRuptures
convert_simpleFaultRupture(node)[source]

Convert a simpleFaultRupture node.

Parameters:node – the rupture node
convert_singlePlaneRupture(node)[source]

Convert a singlePlaneRupture node.

Parameters:node – the rupture node
convert_surfaces(surface_nodes)[source]

Utility to convert a list of surface nodes into a single hazardlib surface. There are three possibilities:

  1. there is a single simpleFaultGeometry node; returns a openquake.hazardlib.geo.simpleFaultSurface instance
  2. there is a single complexFaultGeometry node; returns a openquake.hazardlib.geo.complexFaultSurface instance
  3. there is a list of PlanarSurface nodes; returns a openquake.hazardlib.geo.MultiSurface instance
Parameters:surface_nodes – surface nodes as just described
fname = None
geo_line(edge)[source]

Utility function to convert a node of kind edge into a openquake.hazardlib.geo.Line instance.

Parameters:edge – a node describing an edge
geo_lines(edges)[source]

Utility function to convert a list of edges into a list of openquake.hazardlib.geo.Line instances.

Parameters:edge – a node describing an edge
geo_planar(surface)[source]

Utility to convert a PlanarSurface node with subnodes topLeft, topRight, bottomLeft, bottomRight into a openquake.hazardlib.geo.PlanarSurface instance.

Parameters:surface – PlanarSurface node
get_mag_rake_hypo(node)[source]
class openquake.hazardlib.sourceconverter.SourceConverter(investigation_time=50.0, rupture_mesh_spacing=10.0, complex_fault_mesh_spacing=None, width_of_mfd_bin=1.0, area_source_discretization=None)[source]

Bases: openquake.hazardlib.sourceconverter.RuptureConverter

Convert sources from valid nodes into Hazardlib objects.

convert_UCERFSource(node)

Converts the Ucerf Source node into an SES Control object

convert_areaSource(node)[source]

Convert the given node into an area source object.

Parameters:node – a node with tag areaGeometry
Returns:a openquake.hazardlib.source.AreaSource instance
convert_characteristicFaultSource(node)[source]

Convert the given node into a characteristic fault object.

Parameters:node – a node with tag areaGeometry
Returns:a openquake.hazardlib.source.CharacteristicFaultSource instance
convert_complexFaultSource(node)[source]

Convert the given node into a complex fault object.

Parameters:node – a node with tag areaGeometry
Returns:a openquake.hazardlib.source.ComplexFaultSource instance
convert_hpdist(node)[source]

Convert the given node into a probability mass function for the hypo depth distribution.

Parameters:node – a hypoDepthDist node
Returns:a openquake.hazardlib.pmf.PMF instance
convert_mfdist(node)[source]

Convert the given node into a Magnitude-Frequency Distribution object.

Parameters:node – a node of kind incrementalMFD or truncGutenbergRichterMFD
Returns:a openquake.hazardlib.mfd.EvenlyDiscretizedMFD. or openquake.hazardlib.mfd.TruncatedGRMFD instance
convert_multiPointSource(node)[source]

Convert the given node into a MultiPointSource object.

Parameters:node – a node with tag multiPointGeometry
Returns:a openquake.hazardlib.source.MultiPointSource
convert_nonParametricSeismicSource(node)[source]

Convert the given node into a non parametric source object.

Parameters:node – a node with tag areaGeometry
Returns:a openquake.hazardlib.source.NonParametricSeismicSource instance
convert_npdist(node)[source]

Convert the given node into a Nodal Plane Distribution.

Parameters:node – a nodalPlaneDist node
Returns:a openquake.hazardlib.geo.NodalPlane instance
convert_pointSource(node)[source]

Convert the given node into a point source object.

Parameters:node – a node with tag pointGeometry
Returns:a openquake.hazardlib.source.PointSource instance
convert_simpleFaultSource(node)[source]

Convert the given node into a simple fault object.

Parameters:node – a node with tag areaGeometry
Returns:a openquake.hazardlib.source.SimpleFaultSource instance
convert_sourceGroup(node)[source]

Convert the given node into a SourceGroup object.

Parameters:node – a node with tag sourceGroup
Returns:a SourceGroup instance
convert_sourceModel(node)[source]
class openquake.hazardlib.sourceconverter.SourceGroup(trt, sources=None, name=None, src_interdep='indep', rup_interdep='indep', srcs_weights=None, grp_probability=None, min_mag=None, max_mag=None, id=0, eff_ruptures=-1)[source]

Bases: _abcoll.Sequence

A container for the following parameters:

Parameters:
  • trt (str) – the tectonic region type all the sources belong to
  • sources (list) – a list of hazardlib source objects
  • name – The name of the group
  • src_interdep – A string specifying if the sources in this cluster are independent or mutually exclusive
  • rup_indep – A string specifying if the ruptures within each source of the cluster are independent or mutually exclusive
  • weights – A dictionary whose keys are the source IDs of the cluster and the values are the weights associated with each source
  • min_mag – the minimum magnitude among the given sources
  • max_mag – the maximum magnitude among the given sources
  • id – an optional numeric ID (default None) useful to associate the model to a database object
  • eff_ruptures – the number of ruptures contained in the group; if -1, the number is unknown and has to be computed by using get_set_num_ruptures
classmethod collect(sources)[source]
Parameters:sources – dictionaries with a key ‘tectonicRegion’
Returns:an ordered list of SourceGroup instances
srcs_weights

The weights of the underlying sources. If not specified, returns an array of 1s.

tot_ruptures()[source]
update(src)[source]

Update the attributes sources, min_mag, max_mag according to the given source.

Parameters:src – an instance of :class: openquake.hazardlib.source.base.BaseSeismicSource
openquake.hazardlib.sourceconverter.area_to_point_sources(area_src)[source]

Split an area source into a generator of point sources.

MFDs will be rescaled appropriately for the number of points in the area mesh.

Parameters:area_srcopenquake.hazardlib.source.AreaSource
openquake.hazardlib.sourceconverter.collapse(array)[source]

Collapse a homogeneous array into a scalar; do nothing if the array is not homogenous

openquake.hazardlib.sourceconverter.get_key(node)[source]

Convert the given pointSource node into a tuple

openquake.hazardlib.sourceconverter.get_set_num_ruptures(src)[source]

Extract the number of ruptures and set it

openquake.hazardlib.sourceconverter.mfds2multimfd(mfds)[source]

Convert a list of MFD nodes into a single MultiMFD node

openquake.hazardlib.sourceconverter.split_coords_2d(seq)[source]
Parameters:seq – a flat list with lons and lats
Returns:a validated list of pairs (lon, lat)
>>> split_coords_2d([1.1, 2.1, 2.2, 2.3])
[(1.1, 2.1), (2.2, 2.3)]
openquake.hazardlib.sourceconverter.split_coords_3d(seq)[source]
Parameters:seq – a flat list with lons, lats and depths
Returns:a validated list of (lon, lat, depths) triplets
>>> split_coords_3d([1.1, 2.1, 0.1, 2.3, 2.4, 0.1])
[(1.1, 2.1, 0.1), (2.3, 2.4, 0.1)]
openquake.hazardlib.sourceconverter.split_fault_source(src)[source]

Generator splitting a fault source into several fault sources.

Parameters:src – an instance of openquake.hazardlib.source.base.SeismicSource
openquake.hazardlib.sourceconverter.split_source(src)[source]

Split an area source into point sources and a fault sources into smaller fault sources.

Parameters:src – an instance of openquake.hazardlib.source.base.SeismicSource
openquake.hazardlib.sourceconverter.update_source_model(sm_node)[source]
Parameters:sm_node – a sourceModel Node object containing sourceGroups

sourcewriter

Source model XML Writer

openquake.hazardlib.sourcewriter.build_arbitrary_mfd(mfd)[source]

Parses the arbitrary MFD as a Node

Parameters:mfd – MFD as instance of :class: openquake.hazardlib.mfd.arbitrary.ArbitraryMFD
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_area_source_geometry(area_source)[source]

Returns the area source geometry as a Node

Parameters:area_source – Area source model as an instance of the :class: openquake.hazardlib.source.area.AreaSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_area_source_node(area_source)[source]

Parses an area source to a Node class

Parameters:area_source – Area source as instance of :class: openquake.hazardlib.source.area.AreaSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_characteristic_fault_source_node(source)[source]
openquake.hazardlib.sourcewriter.build_complex_fault_geometry(fault_source)[source]

Returns the complex fault source geometry as a Node

Parameters:fault_source – Complex fault source model as an instance of the :class: openquake.hazardlib.source.complex_fault.ComplexFaultSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_complex_fault_source_node(fault_source)[source]

Parses a complex fault source to a Node class

Parameters:fault_source – Simple fault source as instance of :class: openquake.hazardlib.source.complex_fault.ComplexFaultSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_evenly_discretised_mfd(mfd)[source]

Returns the evenly discretized MFD as a Node

Parameters:mfd – MFD as instance of :class: openquake.hazardlib.mfd.evenly_discretized.EvenlyDiscretizedMFD
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_hypo_depth_dist(hdd)[source]

Returns the hypocentral depth distribution as a Node instance

Parameters:hdd – Hypocentral depth distribution as an instance of :class: openquake.hzardlib.pmf.PMF
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_hypo_list_node(hypo_list)[source]
Parameters:hypo_list – an array of shape (N, 3) with columns (alongStrike, downDip, weight)
Returns:a hypoList node containing N hypo nodes
openquake.hazardlib.sourcewriter.build_linestring_node(line, with_depth=False)[source]

Parses a line to a Node class

Parameters:
Returns:

Instance of openquake.baselib.node.Node

openquake.hazardlib.sourcewriter.build_multi_mfd(mfd)[source]

Parses the MultiMFD as a Node

Parameters:mfd – MFD as instance of :class: openquake.hazardlib.mfd.multi_mfd.MultiMFD
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_multi_point_source_node(multi_point_source)[source]

Parses a point source to a Node class

Parameters:point_source – MultiPoint source as instance of :class: openquake.hazardlib.source.point.MultiPointSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_nodal_plane_dist(npd)[source]

Returns the nodal plane distribution as a Node instance

Parameters:npd – Nodal plane distribution as instance of :class: openquake.hazardlib.pmf.PMF
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_nonparametric_source_node(source)[source]
openquake.hazardlib.sourcewriter.build_point_source_geometry(point_source)[source]

Returns the poing source geometry as a Node

Parameters:point_source – Point source model as an instance of the :class: openquake.hazardlib.source.point.PointSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_point_source_node(point_source)[source]

Parses a point source to a Node class

Parameters:point_source – Point source as instance of :class: openquake.hazardlib.source.point.PointSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_rupture_node(rupt, probs_occur)[source]
Parameters:
  • rupt – a hazardlib rupture object
  • probs_occur – a list of floats with sum 1
openquake.hazardlib.sourcewriter.build_simple_fault_geometry(fault_source)[source]

Returns the simple fault source geometry as a Node

Parameters:fault_source – Simple fault source model as an instance of the :class: openquake.hazardlib.source.simple_fault.SimpleFaultSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_simple_fault_source_node(fault_source)[source]

Parses a simple fault source to a Node class

Parameters:fault_source – Simple fault source as instance of :class: openquake.hazardlib.source.simple_fault.SimpleFaultSource
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_slip_list_node(slip_list)[source]
Parameters:slip_list – an array of shape (N, 2) with columns (slip, weight)
Returns:a hypoList node containing N slip nodes
openquake.hazardlib.sourcewriter.build_source_group(source_group)[source]
openquake.hazardlib.sourcewriter.build_source_model(csm)[source]
openquake.hazardlib.sourcewriter.build_truncated_gr_mfd(mfd)[source]

Parses the truncated Gutenberg Richter MFD as a Node

Parameters:mfd – MFD as instance of :class: openquake.hazardlib.mfd.truncated_gr.TruncatedGRMFD
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.build_youngs_coppersmith_mfd(mfd)[source]

Parses the Youngs & Coppersmith MFD as a node. Note that the MFD does not hold the total moment rate, but only the characteristic rate. Therefore the node is written to the characteristic rate version regardless of whether or not it was originally created from total moment rate

Parameters:mfd – MFD as instance of :class: openquake.hazardlib.mfd.youngs_coppersmith_1985. YoungsCoppersmith1985MFD
Returns:Instance of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.get_distributed_seismicity_source_nodes(source)[source]

Returns list of nodes of attributes common to all distributed seismicity source classes

Parameters:source – Seismic source as instance of :class: openquake.hazardlib.source.area.AreaSource or :class: openquake.hazardlib.source.point.PointSource
Returns:List of instances of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.get_fault_source_nodes(source)[source]

Returns list of nodes of attributes common to all fault source classes

Parameters:source – Fault source as instance of :class: openquake.hazardlib.source.simple_fault.SimpleFaultSource or :class: openquake.hazardlib.source.complex_fault.ComplexFaultSource
Returns:List of instances of openquake.baselib.node.Node
openquake.hazardlib.sourcewriter.get_source_attributes(source)[source]

Retreives a dictionary of source attributes from the source class

Parameters:source – Seismic source as instance of :class: openquake.hazardlib.source.base.BaseSeismicSource
Returns:Dictionary of source attributes
openquake.hazardlib.sourcewriter.hdf5write(h5file, obj, root='')[source]

Write a generic object serializable to a Node-like object into a :class: openquake.baselib.hdf5.File

openquake.hazardlib.sourcewriter.write_source_model(dest, sources_or_groups, name=None)[source]

Writes a source model to XML.

Parameters:
  • dest (str) – Destination path
  • sources_or_groups (list) – Source model as list of sources or a list of SourceGroups
  • name (str) – Name of the source model (if missing, extracted from the filename)

stats

Utilities to compute mean and quantile curves

openquake.hazardlib.stats.apply_stat(f, arraylist, *extra, **kw)[source]
Parameters:
  • f – a callable arraylist -> array (of the same shape and dtype)
  • arraylist – a list of arrays of the same shape and dtype
  • extra – additional positional arguments
  • kw – keyword arguments
Returns:

an array of the same shape and dtype

Broadcast statistical functions to composite arrays. Here is an example:

>>> dt = numpy.dtype([('a', (float, 2)), ('b', float)])
>>> a1 = numpy.array([([1, 2], 3)], dt)
>>> a2 = numpy.array([([4, 5], 6)], dt)
>>> apply_stat(mean_curve, [a1, a2])
array([([2.5, 3.5], 4.5)], 
      dtype=[('a', '<f8', (2,)), ('b', '<f8')])
openquake.hazardlib.stats.compute_pmap_stats(pmaps, stats, weights)[source]
Parameters:
  • pmaps – a list of R probability maps
  • stats – a sequence of S statistic functions
  • weights – a list of R weights
Returns:

a probability map with S internal values

openquake.hazardlib.stats.compute_stats(array, stats, weights)[source]
Parameters:
  • array – an array of R elements (which can be arrays)
  • stats – a sequence of S statistic functions
  • weights – a list of R weights
Returns:

an array of S elements (which can be arrays)

openquake.hazardlib.stats.compute_stats2(arrayNR, stats, weights)[source]
Parameters:
  • arrayNR – an array of (N, R) elements
  • stats – a sequence of S statistic functions
  • weights – a list of R weights
Returns:

an array of (N, S) elements

openquake.hazardlib.stats.max_curve(values, weights=None)[source]

Compute the maximum curve by taking the upper limits of the values; the weights are ignored and present only for API compatibility. The values can be arrays and then the maximum is taken pointwise:

>>> max_curve([numpy.array([.3, .2]), numpy.array([.1, .4])])
array([ 0.3,  0.4])
openquake.hazardlib.stats.mean_curve(values, weights=None)[source]

Compute the mean by using numpy.average on the first axis.

openquake.hazardlib.stats.quantile_curve(quantile, curves, weights=None)[source]

Compute the weighted quantile aggregate of a set of curves.

Parameters:
  • quantile – Quantile value to calculate. Should be in the range [0.0, 1.0].
  • curves – Array of R PoEs (possibly arrays)
  • weights – Array-like of weights, 1 for each input curve, or None
Returns:

A numpy array representing the quantile aggregate

tom

Module openquake.hazardlib.tom contains implementations of probability density functions for earthquake temporal occurrence modeling.

class openquake.hazardlib.tom.PoissonTOM(time_span)[source]

Bases: object

Poissonian temporal occurrence model.

Parameters:time_span – The time interval of interest, in years.
Raises:ValueError – If time_span is not positive.
assert_equal(other)
get_probability_no_exceedance(occurrence_rate, poes)[source]

Compute and return, for a number of ground motion levels and sites, the probability that a rupture with annual occurrence rate given by occurrence_rate and able to cause ground motion values higher than a given level at a site with probability poes, does not cause any exceedance in the time window specified by the time_span parameter given in the constructor.

The probability is computed using the following formula

(1 - e ** (-occurrence_rate * time_span)) ** poes
Parameters:
  • occurrence_rate – The average number of events per year.
  • poes – 2D numpy array containing conditional probabilities the the a rupture occurrence causes a ground shaking value exceeding a ground motion level at a site. First dimension represent sites, second dimension intensity measure levels. poes can be obtained calling the method.
Returns:

2D numpy array containing probabilities of no exceedance. First dimension represents sites, second dimensions intensity measure levels.

get_probability_one_occurrence(occurrence_rate)[source]

Calculate and return the probability of event to occur once within the time range defined by the constructor’s time_span parameter value.

get_probability_one_or_more_occurrences(occurrence_rate)[source]

Calculate and return the probability of event to happen one or more times within the time range defined by constructor’s time_span parameter value.

Calculates probability as 1 - e ** (-occurrence_rate*time_span).

Parameters:occurrence_rate – The average number of events per year.
Returns:Float value between 0 and 1 inclusive.
sample_number_of_occurrences(occurrence_rate, seeds=None)[source]

Draw a random sample from the distribution and return a number of events to occur.

The method uses the numpy random generator, which needs a seed in order to get reproducible results. If the seed is None, it should be set outside of this method.

Parameters:
  • occurrence_rate – The average number of events per year.
  • seeds – Random number generator seeds, one per each occurrence_rate
Returns:

Sampled integer number of events to occur within model’s time span.

valid

Validation library for the engine, the desktop tools, and anything else

class openquake.hazardlib.valid.Choice(*choices)[source]

Bases: object

Check if the choice is valid (case sensitive).

class openquake.hazardlib.valid.ChoiceCI(*choices)[source]

Bases: object

Check if the choice is valid (case insensitive version).

class openquake.hazardlib.valid.Choices(*choices)[source]

Bases: openquake.hazardlib.valid.Choice

Convert the choices, passed as a comma separated string, into a tuple of validated strings. For instance

>>> Choices('xml', 'csv')('xml,csv')
('xml', 'csv')
class openquake.hazardlib.valid.FloatRange(minrange, maxrange)[source]

Bases: object

class openquake.hazardlib.valid.MetaParamSet(name, bases, dic)[source]

Bases: type

Set the .name attribute of every Param instance defined inside any subclass of ParamSet.

class openquake.hazardlib.valid.NoneOr(cast)[source]

Bases: object

Accept the empty string (casted to None) or something else validated by the underlying cast validator.

class openquake.hazardlib.valid.Param(validator, default=<object object>, name=None)[source]

Bases: object

A descriptor for validated parameters with a default, to be used as attributes in ParamSet objects.

Parameters:
  • validator – the validator
  • default – the default value
NODEFAULT = <object object>
class openquake.hazardlib.valid.ParamSet(**names_vals)[source]

Bases: openquake.baselib.hdf5.LiteralAttrs

A set of valid interrelated parameters. Here is an example of usage:

>>> class MyParams(ParamSet):
...     a = Param(positiveint)
...     b = Param(positivefloat)
...
...     def is_valid_not_too_big(self):
...         "The sum of a and b must be under 10: a={a} and b={b}"
...         return self.a + self.b < 10
>>> mp = MyParams(a='1', b='7.2')
>>> mp
<MyParams a=1, b=7.2>
>>> MyParams(a='1', b='9.2').validate()
Traceback (most recent call last):
...
ValueError: The sum of a and b must be under 10: a=1 and b=9.2

The constrains are applied in lexicographic order. The attribute corresponding to a Param descriptor can be set as usual:

>>> mp.a = '2'
>>> mp.a
'2'

A list with the literal strings can be extracted as follows:

>>> mp.to_params()
[('a', "'2'"), ('b', '7.2')]

It is possible to build a new object from a dictionary of parameters which are assumed to be already validated:

>>> MyParams.from_(dict(a="'2'", b='7.2'))
<MyParams a='2', b=7.2>
classmethod check(dic)[source]

Convert a dictionary name->string into a dictionary name->value by converting the string. If the name does not correspond to a known parameter, just ignore it and print a warning.

classmethod from_(dic)[source]

Build a new ParamSet from a dictionary of string-valued parameters which are assumed to be already valid.

params = {}
to_params()[source]

Convert the instance dictionary into a sorted list of pairs (name, valrepr) where valrepr is the string representation of the underlying value.

validate()[source]

Apply the is_valid methods to self and possibly raise a ValueError.

class openquake.hazardlib.valid.Regex(regex)[source]

Bases: object

Compare the value with the given regex

class openquake.hazardlib.valid.SimpleId(length, regex='^[\w_\-]+$')[source]

Bases: object

Check if the given value is a valid ID.

Parameters:
  • length – maximum length of the ID
  • regex – accepted characters
class openquake.hazardlib.valid.SiteParam(lon, lat, depth, z1pt0, z2pt5, measured, vs30, backarc)

Bases: tuple

backarc

Alias for field number 7

depth

Alias for field number 2

lat

Alias for field number 1

lon

Alias for field number 0

measured

Alias for field number 5

vs30

Alias for field number 6

z1pt0

Alias for field number 3

z2pt5

Alias for field number 4

openquake.hazardlib.valid.ab_values(value)[source]

a and b values of the GR magniture-scaling relation. a is a positive float, b is just a float.

openquake.hazardlib.valid.boolean(value)[source]
Parameters:value – input string such as ‘0’, ‘1’, ‘true’, ‘false’
Returns:boolean
>>> boolean('')
False
>>> boolean('True')
True
>>> boolean('false')
False
>>> boolean('t')
Traceback (most recent call last):
    ...
ValueError: Not a boolean: t
openquake.hazardlib.valid.check_levels(imls, imt)[source]

Raise a ValueError if the given levels are invalid.

Parameters:
  • imls – a list of intensity measure and levels
  • imt – the intensity measure type
>>> check_levels([0.1, 0.2], 'PGA')  # ok
>>> check_levels([0.1], 'PGA')
Traceback (most recent call last):
   ...
ValueError: Not enough imls for PGA: [0.1]
>>> check_levels([0.2, 0.1], 'PGA')
Traceback (most recent call last):
   ...
ValueError: The imls for PGA are not sorted: [0.2, 0.1]
>>> check_levels([0.2, 0.2], 'PGA')
Traceback (most recent call last):
   ...
ValueError: Found duplicated levels for PGA: [0.2, 0.2]
openquake.hazardlib.valid.check_weights(nodes_with_a_weight)[source]

Ensure that the sum of the values is 1

Parameters:nodes_with_a_weight – a list of Node objects with a weight attribute
openquake.hazardlib.valid.compose(*validators)[source]

Implement composition of validators. For instance

>>> utf8_not_empty = compose(utf8, not_empty)
openquake.hazardlib.valid.coordinates(value)[source]

Convert a non-empty string into a list of lon-lat coordinates.

>>> coordinates('')
Traceback (most recent call last):
...
ValueError: Empty list of coordinates: ''
>>> coordinates('1.1 1.2')
[(1.1, 1.2, 0.0)]
>>> coordinates('1.1 1.2, 2.2 2.3')
[(1.1, 1.2, 0.0), (2.2, 2.3, 0.0)]
>>> coordinates('1.1 1.2 -0.4, 2.2 2.3 -0.5')
[(1.1, 1.2, -0.4), (2.2, 2.3, -0.5)]
>>> coordinates('0 0 0, 0 0 -1')
Traceback (most recent call last):
...
ValueError: There are overlapping points in 0 0 0, 0 0 -1
openquake.hazardlib.valid.decreasing_probabilities(value)[source]
Parameters:value – input string, comma separated or space separated
Returns:a list of decreasing probabilities
>>> decreasing_probabilities('1')
Traceback (most recent call last):
...
ValueError: Not enough probabilities, found '1'
>>> decreasing_probabilities('0.2 0.1')
[0.2, 0.1]
>>> decreasing_probabilities('0.1 0.2')
Traceback (most recent call last):
...
ValueError: The probabilities 0.1 0.2 are not in decreasing order
openquake.hazardlib.valid.depth(value)
Parameters:value – input string
Returns:a floating point number
openquake.hazardlib.valid.dictionary(value)[source]
Parameters:value – input string corresponding to a literal Python object
Returns:the Python object
>>> dictionary('')
{}
>>> dictionary('{}')
{}
>>> dictionary('{"a": 1}')
{'a': 1}
>>> dictionary('"vs30_clustering: true"')  # an error really done by a user
Traceback (most recent call last):
   ...
ValueError: '"vs30_clustering: true"' is not a valid Python dictionary
>>> dictionary('{"ls": logscale(0.01, 2, 5)}')
{'ls': [0.01, 0.037606030930863933, 0.14142135623730948, 0.53182958969449856, 1.9999999999999991]}
openquake.hazardlib.valid.disagg_outputs(value)[source]

Validate disaggregation outputs. For instance

>>> disagg_outputs('TRT Mag_Dist')
['TRT', 'Mag_Dist']
>>> disagg_outputs('TRT, Mag_Dist')
['TRT', 'Mag_Dist']
openquake.hazardlib.valid.float_(value)[source]
Parameters:value – input string
Returns:a floating point number
openquake.hazardlib.valid.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}
>>> text = "{'active shallow crust': 250., 'default': 200}"
>>> sorted(floatdict(text).items())
[('active shallow crust', 250.0), ('default', 200)]
openquake.hazardlib.valid.floats32(value)[source]
Parameters:value – string of whitespace separated floats
Returns:an array of 32 bit floats
openquake.hazardlib.valid.gsim(value, **kwargs)[source]

Make sure the given value is the name of an available GSIM class.

>>> gsim('BooreAtkinson2011')
'BooreAtkinson2011()'
openquake.hazardlib.valid.hazard_id(value)[source]
>>> hazard_id('')
()
>>> hazard_id('-1')
(-1,)
>>> hazard_id('42')
(42,)
>>> hazard_id('42,3')
(42, 3)
>>> hazard_id('42,3,4')
(42, 3, 4)
>>> hazard_id('42:3')
Traceback (most recent call last):
   ...
ValueError: Invalid hazard_id '42:3'
openquake.hazardlib.valid.hypo_list(nodes)[source]
Parameters:nodes – a hypoList node with N hypocenter nodes
Returns:a numpy array of shape (N, 3) with strike, dip and weight
openquake.hazardlib.valid.integers(value)[source]
Parameters:value – input string
Returns:non-empty list of integers
>>> integers('1, 2')
[1, 2]
>>> integers(' ')
Traceback (most recent call last):
   ...
ValueError: Not a list of integers: ' '
openquake.hazardlib.valid.intensity_measure_type(value)[source]

Make sure value is a valid intensity measure type and return it in a normalized form

>>> intensity_measure_type('SA(0.10)')  # NB: strips the trailing 0
'SA(0.1)'
>>> intensity_measure_type('SA')  # this is invalid
Traceback (most recent call last):
  ...
ValueError: Invalid IMT: 'SA'
openquake.hazardlib.valid.intensity_measure_types(value)[source]
Parameters:value – input string
Returns:non-empty list of Intensity Measure Type objects
>>> intensity_measure_types('PGA')
['PGA']
>>> intensity_measure_types('PGA, SA(1.00)')
['PGA', 'SA(1.0)']
>>> intensity_measure_types('SA(0.1), SA(0.10)')
Traceback (most recent call last):
  ...
ValueError: Duplicated IMTs in SA(0.1), SA(0.10)
openquake.hazardlib.valid.intensity_measure_types_and_levels(value)[source]
Parameters:value – input string
Returns:Intensity Measure Type and Levels dictionary
>>> intensity_measure_types_and_levels('{"SA(0.10)": [0.1, 0.2]}')
{'SA(0.1)': [0.1, 0.2]}
openquake.hazardlib.valid.latitude(value)[source]
Parameters:value – input string
Returns:latitude float, rounded to 5 digits, i.e. 1 meter maximum
>>> latitude('-0.123456')
-0.12346
openquake.hazardlib.valid.latitudes(value)[source]
Parameters:value – a comma separated string of latitudes
Returns:a list of latitudes
openquake.hazardlib.valid.logscale(x_min, x_max, n)[source]
Parameters:
  • x_min – minumum value
  • x_max – maximum value
  • n – number of steps
Returns:

an array of n values from x_min to x_max

openquake.hazardlib.valid.lon_lat(value)[source]
Parameters:value – a pair of coordinates
Returns:a tuple (longitude, latitude)
>>> lon_lat('12 14')
(12.0, 14.0)
openquake.hazardlib.valid.longitude(value)[source]
Parameters:value – input string
Returns:longitude float, rounded to 5 digits, i.e. 1 meter maximum
>>> longitude('0.123456')
0.12346
openquake.hazardlib.valid.longitudes(value)[source]
Parameters:value – a comma separated string of longitudes
Returns:a list of longitudes
openquake.hazardlib.valid.loss_ratios(value)[source]
Parameters:value – input string
Returns:dictionary loss_type -> loss ratios
>>> loss_ratios('{"structural": [0.1, 0.2]}')
{'structural': [0.1, 0.2]}
openquake.hazardlib.valid.mag_scale_rel(value)[source]
Parameters:value – name of a Magnitude-Scale relationship in hazardlib
Returns:the corresponding hazardlib object
openquake.hazardlib.valid.maximum_distance(value)[source]
Parameters:value – input string corresponding to a valid maximum distance
Returns:a IntegrationDistance mapping
openquake.hazardlib.valid.namelist(value)[source]
Parameters:value – input string
Returns:list of identifiers separated by whitespace or commas
>>> namelist('a,b')
['a', 'b']
>>> namelist('a1  b_2       _c')
['a1', 'b_2', '_c']
>>> namelist('a1 b_2 1c')
Traceback (most recent call last):
    ...
ValueError: List of names containing an invalid name: 1c
openquake.hazardlib.valid.nonzero(value)[source]
Parameters:value – input string
Returns:the value unchanged
>>> nonzero('1')
'1'
>>> nonzero('0')
Traceback (most recent call last):
  ...
ValueError: '0' is zero
openquake.hazardlib.valid.not_empty(value)[source]

Check that the string is not all blanks

openquake.hazardlib.valid.pmf(value)[source]

Comvert a string into a Probability Mass Function.

Parameters:value – a sequence of probabilities summing up to 1 (no commas)
Returns:a list of pairs [(probability, index), ...] with index starting from 0
>>> pmf("0.157 0.843")
[(0.157, 0), (0.843, 1)]
openquake.hazardlib.valid.point(value)[source]
Parameters:value – a tuple of coordinates as a string (2D or 3D)
Returns:a tuple of coordinates as a string (2D or 3D)
openquake.hazardlib.valid.point3d(value, lon, lat, depth)[source]

This is used to convert nodes of the form <hypocenter lon=”LON” lat=”LAT” depth=”DEPTH”/>

Parameters:
  • value – None
  • lon – longitude string
  • lat – latitude string
Returns:

a validated triple (lon, lat, depth)

openquake.hazardlib.valid.posList(value)[source]
Parameters:value – a string with the form lon1 lat1 [depth1] ... lonN latN [depthN] without commas, where the depts are optional.
Returns:a list of floats without other validations
openquake.hazardlib.valid.positivefloat(value)[source]
Parameters:value – input string
Returns:positive float
openquake.hazardlib.valid.positivefloats(value)[source]
Parameters:value – string of whitespace separated floats
Returns:a list of positive floats
openquake.hazardlib.valid.positiveint(value)[source]
Parameters:value – input string
Returns:positive integer
openquake.hazardlib.valid.positiveints(value)[source]
>>> positiveints('1, -1')
Traceback (most recent call last):
   ...
ValueError: -1 is negative in '1, -1'
openquake.hazardlib.valid.probabilities(value, rows=0, cols=0)[source]
Parameters:
  • value – input string, comma separated or space separated
  • rows – the number of rows if the floats are in a matrix (0 otherwise)
  • cols – the number of columns if the floats are in a matrix (or 0
Returns:

a list of probabilities

>>> probabilities('')
[]
>>> probabilities('1')
[1.0]
>>> probabilities('0.1 0.2')
[0.1, 0.2]
>>> probabilities('0.1, 0.2')  # commas are ignored
[0.1, 0.2]
openquake.hazardlib.valid.simple_slice(value)[source]
>>> simple_slice('2:5')
(2, 5)
>>> simple_slice('0:None')
(0, None)
openquake.hazardlib.valid.site_param(z1pt0, z2pt5, vs30Type, vs30, lon, lat, depth=0, backarc='false')[source]

Used to convert a node like

<site lon=”24.7125” lat=”42.779167” vs30=”462” vs30Type=”inferred” z1pt0=”100” z2pt5=”5” backarc=”False”/>

into a 7-tuple (z1pt0, z2pt5, measured, vs30, backarc, lon, lat)

openquake.hazardlib.valid.slip_list(nodes)[source]
Parameters:nodes – a slipList node with N slip nodes
Returns:a numpy array of shape (N, 2) with slip angle and weight
openquake.hazardlib.valid.utf8(value)[source]

Check that the string is UTF-8. Returns an encode bytestring.

>>> utf8(b'\xe0')  
Traceback (most recent call last):
...
ValueError: Not UTF-8: ...
openquake.hazardlib.valid.utf8_not_empty(value)[source]

Check that the string is UTF-8 and not empty

openquake.hazardlib.valid.weights(value)[source]

Space-separated list of weights:

>>> weights('0.1 0.2 0.7')
[0.1, 0.2, 0.7]
>>> weights('0.1 0.2 0.8')
Traceback (most recent call last):
  ...
ValueError: The weights do not sum up to 1: [0.1, 0.2, 0.8]
openquake.hazardlib.valid.wkt_polygon(value)[source]

Convert a string with a comma separated list of coordinates into a WKT polygon, by closing the ring.