openquake.hazardlib.geo package¶
Surface classes¶
- openquake.hazardlib.geo.surface package
- Submodules
- openquake.hazardlib.geo.surface.base module
- openquake.hazardlib.geo.surface.complex_fault module
- openquake.hazardlib.geo.surface.gridded module
- openquake.hazardlib.geo.surface.multi module
- openquake.hazardlib.geo.surface.planar module
- openquake.hazardlib.geo.surface.simple_fault module
- Module contents
Geographic primitives and utilities¶
geodetic¶
Module openquake.hazardlib.geo.geodetic
contains functions for geodetic
transformations, optimized for massive calculations.
- openquake.hazardlib.geo.geodetic.EARTH_ELEVATION = -8.848¶
Maximum elevation on Earth in km.
- openquake.hazardlib.geo.geodetic.EARTH_RADIUS = 6371.0¶
Earth radius in km.
- openquake.hazardlib.geo.geodetic.azimuth(lons1, lats1, lons2, lats2)[source]¶
Calculate the azimuth between two points or two collections of points.
Parameters are the same as for
geodetic_distance()
.Implements an “alternative formula” from http://williams.best.vwh.net/avform.htm#Crs
- Returns
Azimuth as an angle between direction to north from first point and direction to the second point measured clockwise in decimal degrees.
- openquake.hazardlib.geo.geodetic.distance(lons1, lats1, depths1, lons2, lats2, depths2)[source]¶
Calculate a distance between two points (or collections of points) considering points’ depth.
Calls
geodetic_distance()
, finds the “vertical” distance between points by subtracting one depth from another and combine both using Pythagoras theorem.- Returns
Distance in km, a square root of sum of squares of
geodetic
distance and vertical distance, which is just a difference between depths.
- openquake.hazardlib.geo.geodetic.distance_matrix(lons, lats, diameter=12742.0)[source]¶
- Parameters
lons – array of m longitudes
lats – array of m latitudes
- Returns
matrix of (m, m) distances
- openquake.hazardlib.geo.geodetic.distance_to_arc(alon, alat, aazimuth, plons, plats)[source]¶
Calculate a closest distance between a great circle arc and a point (or a collection of points).
- Parameters
alat (float alon,) – Arc reference point longitude and latitude, in decimal degrees.
azimuth – Arc azimuth (an angle between direction to a north and arc in clockwise direction), measured in a reference point, in decimal degrees.
plats (float plons,) – Longitudes and latitudes of points to measure distance. Either scalar values or numpy arrays of decimal degrees.
- Returns
Distance in km, a scalar value or numpy array depending on
plons
andplats
. A distance is negative if the target point lies on the right hand side of the arc.
Solves a spherical triangle formed by reference point, target point and a projection of target point to a reference great circle arc.
- openquake.hazardlib.geo.geodetic.distance_to_semi_arc(alon, alat, aazimuth, plons, plats)[source]¶
In this method we use a reference system centerd on (alon, alat) and with the y-axis corresponding to aazimuth direction to calculate the minimum distance from a semiarc with generates in (alon, alat).
Parameters are the same as for
distance_to_arc()
.
- openquake.hazardlib.geo.geodetic.distances(lon, lat, lons, lats)[source]¶
Calculate the geodetic distances between one point and a collection of points by using numba. Assume the parameters are in radians.
- openquake.hazardlib.geo.geodetic.distances_to_arc(lon, lat, azi, lons, lats)[source]¶
Calculate the distances between one arc and a collection of points by using numba. Assume the parameters are in degrees.
- openquake.hazardlib.geo.geodetic.fast_azimuth(lon, lat, lons, lats)[source]¶
Calculate the azimuths of a collection of points with respect to a reference point.
- openquake.hazardlib.geo.geodetic.fast_point_at(lon, lat, azimuth, distance)[source]¶
Perform a forward geodetic transformation: find points lying at a given distances from a given point on a great circle arc defined by azimuth.
- Parameters
lat (lon,) – Coordinates of the reference point, in radians
azimuth – An azimuth of a great circle arc of interest measured in a reference point in decimal degrees.
distance – Distance to target point in km.
- Returns
Array of shape (2, N) with longitudes and latitudes
Implements the same approach as
npoints_towards()
.
- openquake.hazardlib.geo.geodetic.fast_spherical_to_cartesian(lons, lats, deps)[source]¶
Return the position vectors (in Cartesian coordinates) of list of spherical coordinates.
For equations see: http://mathworld.wolfram.com/SphericalCoordinates.html.
Parameters are components of spherical coordinates in a form of scalars, lists or numpy arrays.
depths
can beNone
in which case it’s considered zero for all points.- Returns
np.array
of 3d vectors representing points’ coordinates in Cartesian space in km. The array has shape lons.shape + (3,). In particular, iflons
andlats
are scalars the result is a 3D vector and if they are vectors the result is a matrix of shape (N, 3).
See also
cartesian_to_spherical()
.
- openquake.hazardlib.geo.geodetic.geodetic_distance(lons1, lats1, lons2, lats2, diameter=12742.0)[source]¶
Calculate the geodetic distance between two points or two collections of points.
Parameters are coordinates in decimal degrees. They could be scalar float numbers or numpy arrays, in which case they should “broadcast together”.
Implements http://williams.best.vwh.net/avform.htm#Dist
- Returns
Distance in km, floating point scalar or numpy array of such.
- openquake.hazardlib.geo.geodetic.intervals_between(lon1, lat1, depth1, lon2, lat2, depth2, length)[source]¶
Find a list of points between two given ones that lie on the same great circle arc and are equally spaced by
length
km.- Parameters
depth1 (float lon1, lat1,) – Coordinates of a point to start placing intervals from. The first point in the resulting list has these coordinates.
depth2 (float lon2, lat2,) – Coordinates of the other end of the great circle arc segment to put intervals on. The last resulting point might be closer to the first reference point than the second one or further, since the number of segments is taken as rounded division of length between two reference points and
length
.length – Required distance between two subsequent resulting points, in km.
- Returns
Tuple of three 1d numpy arrays: longitudes, latitudes and depths of resulting points respectively.
Rounds the distance between two reference points with respect to
length
and callsnpoints_towards()
.
- openquake.hazardlib.geo.geodetic.min_distance_to_segment(seglons, seglats, lons, lats)[source]¶
This function computes the shortest distance to a segment in a 2D reference system.
- Parameters
seglons – A list or an array of floats specifying the longitude values of the two vertexes delimiting the segment.
seglats – A list or an array of floats specifying the latitude values of the two vertexes delimiting the segment.
lons – A list or a 1D array of floats specifying the longitude values of the points for which the calculation of the shortest distance is requested.
lats – A list or a 1D array of floats specifying the latitude values of the points for which the calculation of the shortest distance is requested.
- Returns
An array of the same shape as lons which contains for each point defined by (lons, lats) the shortest distance to the segment. Distances are negative for those points that stay on the ‘left side’ of the segment direction and whose projection lies within the segment edges. For all other points, distance is positive.
- openquake.hazardlib.geo.geodetic.min_geodetic_distance(a, b)[source]¶
Compute the minimum distance between first mesh and each point of the second mesh when both are defined on the earth surface.
- Parameters
a – a pair of (lons, lats) or an array of cartesian coordinates
b – a pair of (lons, lats) or an array of cartesian coordinates
- openquake.hazardlib.geo.geodetic.npoints_between(lon1, lat1, depth1, lon2, lat2, depth2, npoints)[source]¶
Find a list of specified number of points between two given ones that are equally spaced along the great circle arc connecting given points.
- Parameters
depth1 (float lon1, lat1,) – Coordinates of a point to start from. The first point in a resulting list has these coordinates.
depth2 (float lon2, lat2,) – Coordinates of a point to finish at. The last point in a resulting list has these coordinates.
npoints – Integer number of points to return. First and last points count, so if there have to be two intervals,
npoints
should be 3.
- Returns
Tuple of three 1d numpy arrays: longitudes, latitudes and depths of resulting points respectively.
Finds distance between two reference points and calls
npoints_towards()
.
- openquake.hazardlib.geo.geodetic.npoints_towards(lon, lat, depth, azimuth, hdist, vdist, npoints)[source]¶
Find a list of specified number of points starting from a given one along a great circle arc with a given azimuth measured in a given point.
- Parameters
depth (float lon, lat,) – Coordinates of a point to start from. The first point in a resulting list has these coordinates.
azimuth – A direction representing a great circle arc together with a reference point.
hdist – Horizontal (geodetic) distance from reference point to the last point of the resulting list, in km.
vdist – Vertical (depth) distance between reference and the last point, in km.
npoints – Integer number of points to return. First and last points count, so if there have to be two intervals,
npoints
should be 3.
- Returns
Tuple of three 1d numpy arrays: longitudes, latitudes and depths of resulting points respectively.
Implements “completely general but more complicated algorithm” from http://williams.best.vwh.net/avform.htm#LL
- openquake.hazardlib.geo.geodetic.point_at(lon, lat, azimuth, distance)[source]¶
Perform a forward geodetic transformation: find a point lying at a given distance from a given one on a great circle arc defined by azimuth.
- Parameters
lat (float lon,) – Coordinates of a reference point, in decimal degrees.
azimuth – An azimuth of a great circle arc of interest measured in a reference point in decimal degrees.
distance – Distance to target point in km.
- Returns
Tuple of two float numbers: longitude and latitude of a target point in decimal degrees respectively.
Implements the same approach as
npoints_towards()
.
- openquake.hazardlib.geo.geodetic.spherical_to_cartesian(lons, lats, depths=None)[source]¶
Return the position vectors (in Cartesian coordinates) of list of spherical coordinates.
For equations see: http://mathworld.wolfram.com/SphericalCoordinates.html.
Parameters are components of spherical coordinates in a form of scalars, lists or numpy arrays.
depths
can beNone
in which case it’s considered zero for all points.- Returns
np.array
of 3d vectors representing points’ coordinates in Cartesian space in km. The array has shape lons.shape + (3,). In particular, iflons
andlats
are scalars the result is a 3D vector and if they are vectors the result is a matrix of shape (N, 3).
See also
cartesian_to_spherical()
.
line¶
Module openquake.hazardlib.geo.line
defines Line
.
- class openquake.hazardlib.geo.line.Line(points)[source]¶
Bases:
object
This class represents a geographical line, which is basically a sequence of geographical points.
A line is defined by at least two points.
- Parameters
points (list of
Point
instances) – The sequence of points defining this line.
- average_azimuth()[source]¶
Calculate and return weighted average azimuth of all line’s segments in decimal degrees. Uses formula from http://en.wikipedia.org/wiki/Mean_of_circular_quantities >>> from openquake.hazardlib.geo.point import Point as P >>> ‘%.1f’ % Line([P(0, 0), P(1e-5, 1e-5)]).average_azimuth() ‘45.0’ >>> ‘%.1f’ % Line([P(0, 0), P(0, 1e-5), P(1e-5, 1e-5)]).average_azimuth() ‘45.0’ >>> line = Line([P(0, 0), P(-2e-5, 0), P(-2e-5, 1.154e-5)]) >>> ‘%.1f’ % line.average_azimuth() ‘300.0’
- classmethod from_vectors(lons, lats, deps=None)[source]¶
Creates a line from three numpy.ndarray instances containing longitude, latitude and depths values
- get_length() float [source]¶
Calculate the length of the line as a sum of lengths of all its segments.
- Returns
Total length in km.
- get_lengths() numpy.ndarray [source]¶
Calculate a numpy.ndarray instance with the length of the segments composing the polyline.
- Returns
Segments length in km.
- get_tu(mesh)[source]¶
Computes the U and T coordinates of the GC2 method for a mesh of points.
- Parameters
mesh – An instance of
openquake.hazardlib.geo.mesh.Mesh
- get_tu_hat()[source]¶
Return the unit vectors defining the local origin for each segment of the trace.
- Parameters
sx – The vector with the x coordinates of the trace
sy – The vector with the y coordinates of the trace
- Returns
Two arrays of size n x 3 (when n is the number of segments composing the trace
- get_ui_ti(mesh, uhat, that)[source]¶
Compute the t and u coordinates. ti and ui have shape (num_segments x num_sites).
- horizontal()[source]¶
Check if this line is horizontal (i.e. all depths of points are equal).
- Returns bool
True if this line is horizontal, false otherwise.
- keep_corners(delta)[source]¶
Removes the points where the change in direction is lower than a tolerance value.
- Parameters
delta – An angle in decimal degrees
- on_surface()[source]¶
Check if this line is defined on the surface (i.e. all points are on the surfance, depth=0.0).
- Returns bool
True if this line is on the surface, false otherwise.
- resample(section_length)[source]¶
Resample this line into sections. The first point in the resampled line corresponds to the first point in the original line. Starting from the first point in the original line, a line segment is defined as the line connecting the last point in the resampled line and the next point in the original line. The line segment is then split into sections of length equal to
section_length
. The resampled line is obtained by concatenating all sections. The number of sections in a line segment is calculated as follows:round(segment_length / section_length)
. Note that the resulting line has a length that is an exact multiple ofsection_length
, therefore its length is in general smaller or greater (depending on the rounding) than the length of the original line. For a straight line, the difference between the resulting length and the original length is at maximum half of thesection_length
. For a curved line, the difference my be larger, because of corners getting cut.- Parameters
section_length (float) – The length of the section, in km.
- Returns
A new line resampled into sections based on the given length.
- Return type
An instance of
Line
- openquake.hazardlib.geo.line.get_average_azimuth(azimuths, distances) float [source]¶
Computes the average azimuth.
- Parameters
azimuths – A
numpy.ndarray
instancedistances – A
numpy.ndarray
instance
- Returns
A float with the mean azimuth in decimal degrees
- openquake.hazardlib.geo.line.get_tu(ui, ti, sl, weights)[source]¶
Compute the T and U quantitities.
- Parameters
ui – A
numpy.ndarray
instance of cardinality (num segments x num sites)ti – A
numpy.ndarray
instance of cardinality (num segments x num sites)sl – A
numpy.ndarray
instance with the segments’ lengthweights – A
numpy.ndarray
instance of cardinality (num segments x num sites)
mesh¶
Module openquake.hazardlib.geo.mesh
defines classes Mesh
and
its subclass RectangularMesh
.
- class openquake.hazardlib.geo.mesh.Mesh(lons, lats, depths=None)[source]¶
Bases:
object
Mesh object represent a collection of points and provides the most efficient way of keeping those collections in memory.
- Parameters
lons – A numpy array of longitudes. Can be 1D or 2D.
lats – Numpy array of latitudes. The array must be of the same shape as
lons
.depths – Either
None
, which means that all points the mesh consists of are lying on the earth surface (have zero depth) or numpy array of the same shape as previous two.
Mesh object can also be created from a collection of points, see
from_points_list()
.- DIST_TOLERANCE = 0.005¶
Tolerance level to be used in various spatial operations when approximation is required – set to 5 meters.
- property depths¶
- classmethod from_coords(coords, sort=True)[source]¶
Create a mesh object from a list of 3D coordinates (by sorting them)
- Params coords
list of coordinates
- Parameters
sort – flag (default True)
- Returns
a
Mesh
instance
- get_closest_points(mesh)[source]¶
Find closest point of this mesh for each point in the other mesh
- Returns
Mesh
object of the same shape as mesh with closest points from this one at respective indices.
- get_convex_hull()[source]¶
Get a convex polygon object that contains projections of all the points of the mesh.
- Returns
Instance of
openquake.hazardlib.geo.polygon.Polygon
that is a convex hull around all the points in this mesh. If the original mesh had only one point, the resulting polygon has a square shape with a side length of 10 meters. If there were only two points, resulting polygon is a stripe 10 meters wide.
- get_distance_matrix()[source]¶
Compute and return distances between each pairs of points in the mesh.
This method requires that the coordinate arrays are one-dimensional. NB: the depth of the points is ignored
Warning
Because of its quadratic space and time complexity this method is safe to use for meshes of up to several thousand points. For mesh of 10k points it needs ~800 Mb for just the resulting matrix and four times that much for intermediate storage.
- Returns
Two-dimensional numpy array, square matrix of distances. The matrix has zeros on main diagonal and positive distances in kilometers on all other cells. That is, value in cell (3, 5) is the distance between mesh’s points 3 and 5 in km, and it is equal to value in cell (5, 3).
- get_joyner_boore_distance(mesh)[source]¶
Compute and return Joyner-Boore distance to each point of
mesh
. Point’s depth is ignored.See
openquake.hazardlib.geo.surface.base.BaseSurface.get_joyner_boore_distance()
for definition of this distance.- Returns
numpy array of distances in km of the same shape as
mesh
. Distance value is considered to be zero if a point lies inside the polygon enveloping the projection of the mesh or on one of its edges.
- get_min_distance(mesh)[source]¶
Compute and return the minimum distance from the mesh to each point in another mesh.
- Returns
numpy array of distances in km of shape (self.size, mesh.size)
Method doesn’t make any assumptions on arrangement of the points in either mesh and instead calculates the distance from each point of this mesh to each point of the target mesh and returns the lowest found for each.
- property lats¶
- property lons¶
- property shape¶
Return the shape of this mesh.
- Returns tuple
The shape of this mesh as (rows, columns)
- property xyz¶
- Returns
an array of shape (N, 3) with the cartesian coordinates
- class openquake.hazardlib.geo.mesh.RectangularMesh(lons, lats, depths=None)[source]¶
Bases:
openquake.hazardlib.geo.mesh.Mesh
A specification of
Mesh
that requires coordinate numpy-arrays to be two-dimensional.Rectangular mesh is meant to represent not just an unordered collection of points but rather a sort of table of points, where index of the point in a mesh is related to it’s position with respect to neighbouring points.
- classmethod from_points_list(points)[source]¶
Create a rectangular mesh object from a list of lists of points. Lists in a list are supposed to have the same length.
- Parameters
point – List of lists of
Point
objects.
- get_cell_dimensions()[source]¶
Calculate centroid, width, length and area of each mesh cell.
- Returns
Tuple of four elements, each being 2d numpy array. Each array has both dimensions less by one the dimensions of the mesh, since they represent cells, not vertices. Arrays contain the following cell information:
centroids, 3d vectors in a Cartesian space,
length (size along row of points) in km,
width (size along column of points) in km,
area in square km.
- get_mean_inclination_and_azimuth()[source]¶
Calculate weighted average inclination and azimuth of the mesh surface.
- Returns
Tuple of two float numbers: inclination angle in a range [0, 90] and azimuth in range [0, 360) (in decimal degrees).
The mesh is triangulated, the inclination and azimuth for each triangle is computed and average values weighted on each triangle’s area are calculated. Azimuth is always defined in a way that inclination angle doesn’t exceed 90 degree.
- get_mean_width()[source]¶
Calculate and return (weighted) mean width (km) of a mesh surface.
The length of each mesh column is computed (summing up the cell widths in a same column), and the mean value (weighted by the mean cell length in each column) is returned.
- get_middle_point()[source]¶
Return the middle point of the mesh.
- Returns
An instance of
Point
.
The middle point is taken from the middle row and a middle column of the mesh if there are odd number of both. Otherwise the geometric mean point of two or four middle points.
- triangulate()[source]¶
Convert mesh points to vectors in Cartesian space.
- Returns
Tuple of four elements, each being 2d numpy array of 3d vectors (the same structure and shape as the mesh itself). Those arrays are:
points vectors,
vectors directed from each point (excluding the last column) to the next one in a same row →,
vectors directed from each point (excluding the first row) to the previous one in a same column ↑,
vectors pointing from a bottom left point of each mesh cell to top right one ↗.
So the last three arrays of vectors allow to construct triangles covering the whole mesh.
nodalplane¶
Module openquake.hazardlib.geo.nodalplane
implements
NodalPlane
.
- class openquake.hazardlib.geo.nodalplane.NP(strike, dip, rake)¶
Bases:
tuple
- dip¶
Alias for field number 1
- rake¶
Alias for field number 2
- strike¶
Alias for field number 0
- class openquake.hazardlib.geo.nodalplane.NodalPlane(strike, dip, rake)[source]¶
Bases:
object
Nodal plane represents earthquake rupture orientation and propagation direction.
- Parameters
strike – Angle between line created by the intersection of rupture plane and the North direction (defined between 0 and 360 degrees).
dip – Angle between earth surface and fault plane (defined between 0 and 90 degrees).
rake – Angle describing rupture propagation direction (defined between -180 and +180 degrees).
- Raises
ValueError – If any of parameters exceeds the definition range.
- classmethod check_dip(dip)[source]¶
Check if
dip
is in range(0, 90]
and raiseValueError
otherwise.
point¶
Module openquake.hazardlib.geo.point
defines Point
.
- class openquake.hazardlib.geo.point.Point(longitude, latitude, depth=0.0)[source]¶
Bases:
object
This class represents a geographical point in terms of longitude, latitude, and depth (with respect to the Earth surface).
- Parameters
longitude (float) – Point longitude, in decimal degrees.
latitude (float) – Point latitude, in decimal degrees.
depth (float) – Point depth (default to 0.0), in km. Depth > 0 indicates a point below the earth surface, and depth < 0 above the earth surface.
- EQUALITY_DISTANCE = 0.001¶
The distance between two points for them to be considered equal, in km.
- azimuth(point)[source]¶
Compute the azimuth (in decimal degrees) between this point and the given point.
- Parameters
point (Instance of
Point
) – Destination point.- Returns
The azimuth, value in a range
[0, 360)
.- Return type
float
- closer_than(mesh, radius)[source]¶
Check for proximity of points in the
mesh
.- Parameters
mesh –
openquake.hazardlib.geo.mesh.Mesh
instance.radius – Proximity measure in km.
- Returns
Numpy array of boolean values in the same shape as the mesh coordinate arrays with
True
on indexes of points that are not further thanradius
km from this point. Functiondistance()
is used to calculate distances to points of the mesh. Points of the mesh that lie exactlyradius
km away from this point also haveTrue
in their indices.
- distance(point)[source]¶
Compute the distance (in km) between this point and the given point.
Distance is calculated using pythagoras theorem, where the hypotenuse is the distance and the other two sides are the horizontal distance (great circle distance) and vertical distance (depth difference between the two locations).
- Parameters
point (Instance of
Point
) – Destination point.- Returns
The distance.
- Return type
float
- distance_to_mesh(mesh, with_depths=True)[source]¶
Compute distance (in km) between this point and each point of
mesh
.- Parameters
mesh –
Mesh
of points to calculate distance to.with_depths – If
True
(by default), distance is calculated between actual point and the mesh, geodetic distance of projections is combined with vertical distance (difference of depths). If this is set toFalse
, only geodetic distance between projections is calculated.
- Returns
Numpy array of floats of the same shape as
mesh
with distance values in km in respective indices.
- equally_spaced_points(point, distance)[source]¶
Compute the set of points equally spaced between this point and the given point.
- classmethod from_vector(vector)[source]¶
Create a point object from a 3d vector in Cartesian space.
- Parameters
vector – Tuple, list or numpy array of three float numbers representing point coordinates in Cartesian 3d space.
- Returns
A
Point
object created from those coordinates.
- on_surface()[source]¶
Check if this point is defined on the surface (depth is 0.0).
- Returns bool
True if this point is on the surface, false otherwise.
- point_at(horizontal_distance, vertical_increment, azimuth)[source]¶
Compute the point with given horizontal, vertical distances and azimuth from this point.
- Parameters
horizontal_distance (float) – Horizontal distance, in km.
vertical_increment (float) – Vertical increment, in km. When positive, the new point has a greater depth. When negative, the new point has a smaller depth.
- Returns
The point at the given distances.
- Return type
Instance of
Point
- to_polygon(radius)[source]¶
Create a circular polygon with specified radius centered in the point.
- Parameters
radius – Required radius of a new polygon, in km.
- Returns
Instance of
Polygon
that approximates a circle around the point with specified radius.
- property wkt2d¶
Generate WKT (Well-Known Text) to represent this point in 2 dimensions (ignoring depth).
- property x¶
Alias for .longitude
- property y¶
Alias for .latitude
- property z¶
Alias for .depth
polygon¶
Module openquake.hazardlib.geo.polygon
defines Polygon
.
- class openquake.hazardlib.geo.polygon.Polygon(points)[source]¶
Bases:
object
Polygon objects represent an area on the Earth surface.
- Parameters
points – The list of
Point
objects defining the polygon vertices. The points are connected by great circle arcs in order of appearance. Polygon segment should not cross another polygon segment. At least three points must be defined.- Raises
ValueError – If
points
contains less than three unique points or if polygon perimeter intersects itself.
- property centroid¶
- Returns
the centroid of the underlying 2D polygon
- property coords¶
Coordinates of the polygon as a linear ring, rounded to 5 digits
- dilate(dilation)[source]¶
Extend the polygon to a specified buffer distance.
Note
In extreme cases where dilation of a polygon creates holes, thus resulting in a multi-polygon, we discard the holes and simply return the ‘exterior’ of the shape.
- Parameters
dilation – Distance in km to extend polygon borders to.
- Returns
New
Polygon
object with (in general) more vertices and border that is approximatelydilation
km far (measured perpendicularly to edges and circularly to vertices) from the border of original polygon.
- discretize(mesh_spacing)[source]¶
Get a mesh of uniformly spaced points inside the polygon area with distance of
mesh_spacing
km between.- Returns
An instance of
Mesh
that holds the points data. Mesh is created with no depth information (all the points are on the Earth surface).
- classmethod from_wkt(wkt_string)[source]¶
Create a polygon object from a WKT (Well-Known Text) string.
- Parameters
wkt_string – A standard WKT polygon string.
- Returns
New
Polygon
object.
- intersects(mesh)[source]¶
Check for intersection with each point of the
mesh
.Mesh coordinate values are in decimal degrees.
- Parameters
mesh –
openquake.hazardlib.geo.mesh.Mesh
instance.- Returns
Numpy array of bool values in the same shapes in the input coordinate arrays with
True
on indexes of points that lie inside the polygon or on one of its edges andFalse
for points that neither lie inside nor touch the boundary.
- property wkt¶
Generate WKT (Well-Known Text) to represent this polygon.
- openquake.hazardlib.geo.polygon.UPSAMPLING_STEP_KM = 100¶
Polygon upsampling step for long edges, in kilometers. See
get_resampled_coordinates()
.
- openquake.hazardlib.geo.polygon.get_resampled_coordinates(lons, lats)[source]¶
Resample polygon line segments and return the coordinates of the new vertices. This limits distortions when projecting a polygon onto a spherical surface.
Parameters define longitudes and latitudes of a point collection in the form of lists or numpy arrays.
- Returns
A tuple of two numpy arrays: longitudes and latitudes of resampled vertices.
utils¶
Module openquake.hazardlib.geo.utils
contains functions that are common
to several geographical primitives and some other low-level spatial operations.
- class openquake.hazardlib.geo.utils.OrthographicProjection(west, east, north, south)[source]¶
Bases:
object
Callable OrthographicProjection object that can perform both forward and reverse projection (converting from longitudes and latitudes to x and y values on 2d-space and vice versa). The call takes three arguments: first two are numpy arrays of longitudes and latitudes or abscissae and ordinates of points to project and the third one is a boolean that allows to choose what operation is requested – is it forward or reverse one.
True
value given to third positional argument (or keyword argument “reverse”) indicates that the projection of points in 2d space back to earth surface is needed. The default value for “reverse” argument isFalse
, which means forward projection (degrees to kilometers).Raises
ValueError
in forward projection mode if any of the target points is further than 90 degree (along the great circle arc) from the projection center.Parameters are given as floats, representing decimal degrees (first two are longitudes and last two are latitudes). They define a bounding box in a spherical coordinates of the collection of points that is about to be projected. The center point of the projection (coordinates (0, 0) in Cartesian space) is set to the middle point of that bounding box. The resulting projection is defined for spherical coordinates that are not further from the bounding box center than 90 degree on the great circle arc.
The result projection is of type Orthographic. This projection is prone to distance, area and angle distortions everywhere outside of the center point, but still can be used for checking shapes: verifying if line intersects itself (like in
line_intersects_itself()
) or if point is inside of a polygon (like inopenquake.hazardlib.geo.polygon.Polygon.discretize()
). It can be also used for measuring distance to an extent of around 700 kilometers (error doesn’t exceed 1 km up until then).
- class openquake.hazardlib.geo.utils.PolygonPlotter(ax)[source]¶
Bases:
object
Add polygons to a given axis object
- exception openquake.hazardlib.geo.utils.SiteAssociationError[source]¶
Bases:
Exception
Raised when there are no sites close enough
- class openquake.hazardlib.geo.utils.SphericalBB(west, east, north, south)¶
Bases:
tuple
- east¶
Alias for field number 1
- north¶
Alias for field number 2
- south¶
Alias for field number 3
- west¶
Alias for field number 0
- openquake.hazardlib.geo.utils.angular_distance(km, lat=0, lat2=None)[source]¶
Return the angular distance of two points at the given latitude.
>>> '%.3f' % angular_distance(100, lat=40) '1.174' >>> '%.3f' % angular_distance(100, lat=80) '5.179'
- openquake.hazardlib.geo.utils.assoc(objects, sitecol, assoc_dist, mode)[source]¶
Associate geographic objects to a site collection.
- Parameters
objects – something with .lons, .lats or [‘lon’] [‘lat’], or a list of lists of objects with a .location attribute (i.e. assets_by_site)
assoc_dist – the maximum distance for association
mode – if ‘strict’ fail if at least one site is not associated if ‘error’ fail if all sites are not associated
- Returns
(filtered site collection, filtered objects)
- openquake.hazardlib.geo.utils.assoc_to_polygons(polygons, data, sitecol, mode)[source]¶
Associate data from a shapefile with polygons to a site collection :param polygons: polygon shape data :param data: rest of the data belonging to the shapes :param sitecol: a (filtered) site collection :param mode: ‘strict’, ‘warn’ or ‘filter’ :returns: filtered site collection, filtered objects, discarded
- openquake.hazardlib.geo.utils.bbox2poly(bbox)[source]¶
- Parameters
bbox – a geographic bounding box West-East-North-South
- Returns
a list of pairs corrisponding to the bbox polygon
- openquake.hazardlib.geo.utils.cartesian_to_spherical(arrayN3)[source]¶
Return the spherical coordinates for coordinates in Cartesian space.
This function does an opposite to
spherical_to_cartesian()
.- Parameters
arrayN3 – Array of cartesian coordinates of shape (N, 3)
- Returns
Array of shape (3, N) representing longitude (decimal degrees), latitude (decimal degrees) and depth (km) in specified order.
- openquake.hazardlib.geo.utils.check_extent(lons, lats, msg='')[source]¶
- Parameters
lons – an array of longitudes (more than one)
lats – an array of latitudes (more than one)
- Params msg
message to display in case of too large extent
- Returns
(dx, dy, dz) in km (rounded)
- openquake.hazardlib.geo.utils.clean_points(points)[source]¶
Given a list of
Point
objects, return a new list with adjacent duplicate points removed.
- openquake.hazardlib.geo.utils.cross_idl(lon1, lon2, *lons)[source]¶
Return True if two longitude values define line crossing international date line.
>>> cross_idl(-45, 45) False >>> cross_idl(-180, -179) False >>> cross_idl(180, 179) False >>> cross_idl(45, -45) False >>> cross_idl(0, 0) False >>> cross_idl(-170, 170) True >>> cross_idl(170, -170) True >>> cross_idl(-180, 180) True
- openquake.hazardlib.geo.utils.fix_lon(lon)[source]¶
- Returns
a valid longitude in the range -180 <= lon < 180
>>> fix_lon(11) 11 >>> fix_lon(181) -179 >>> fix_lon(-182) 178
- openquake.hazardlib.geo.utils.geohash(lon, lat, length)[source]¶
Encode a position given in lon, lat into a geohash of the given lenght
>>> geohash(lon=10, lat=45, length=5) b'spzpg'
- openquake.hazardlib.geo.utils.get_bounding_box(obj, maxdist)[source]¶
Return the dilated bounding box of a geometric object.
- Parameters
obj – an object with method .get_bounding_box, or with an attribute .polygon or a list of locations
maxdist – maximum distance in km
- openquake.hazardlib.geo.utils.get_longitudinal_extent(lon1, lon2)[source]¶
Return the distance between two longitude values as an angular measure. Parameters represent two longitude values in degrees.
- Returns
Float, the angle between
lon1
andlon2
in degrees. Value is positive iflon2
is on the east fromlon1
and negative otherwise. Absolute value of the result doesn’t exceed 180 for valid parameters values.
- openquake.hazardlib.geo.utils.get_middle_point(lon1, lat1, lon2, lat2)[source]¶
Given two points return the point exactly in the middle lying on the same great circle arc.
Parameters are point coordinates in degrees.
- Returns
Tuple of longitude and latitude of the point in the middle.
- openquake.hazardlib.geo.utils.get_spherical_bounding_box(lons, lats)[source]¶
Given a collection of points find and return the bounding box, as a pair of longitudes and a pair of latitudes.
Parameters define longitudes and latitudes of a point collection respectively in a form of lists or numpy arrays.
- Returns
A tuple of four items. These items represent western, eastern, northern and southern borders of the bounding box respectively. Values are floats in decimal degrees.
- Raises
ValueError – If points collection has the longitudinal extent of more than 180 degrees (it is impossible to define a single hemisphere bound to poles that would contain the whole collection).
- openquake.hazardlib.geo.utils.line_intersects_itself(lons, lats, closed_shape=False)[source]¶
Return
True
if line of points intersects itself. Line with the last point repeating the first one considered intersecting itself.The line is defined by lists (or numpy arrays) of points’ longitudes and latitudes (depth is not taken into account).
- Parameters
closed_shape – If
True
the line will be checked twice: first time with its original shape and second time with the points sequence being shifted by one point (the last point becomes first, the first turns second and so on). This is useful for checking that the sequence of points defines a validPolygon
.
- openquake.hazardlib.geo.utils.min_distance(xyz, xyzs)[source]¶
- Parameters
xyz – an array of shape (3,)
xyzs – an array of shape (N, 3)
- Returns
the minimum euclidean distance between the point and the points
- openquake.hazardlib.geo.utils.normalized(vector)[source]¶
Get unit vector for a given one.
- Parameters
vector – Numpy vector as coordinates in Cartesian space, or an array of such.
- Returns
Numpy array of the same shape and structure where all vectors are normalized. That is, each coordinate component is divided by its vector’s length.
- openquake.hazardlib.geo.utils.plane_fit(points)[source]¶
This fits an n-dimensional plane to a set of points. See http://stackoverflow.com/questions/12299540/plane-fitting-to-4-or-more-xyz-points
- Parameters
points – An instance of :class:~numpy.ndarray. The number of columns must be equal to three.
- Returns
A point on the plane and the normal to the plane.
- openquake.hazardlib.geo.utils.point_to_polygon_distance(polygon, pxx, pyy)[source]¶
Calculate the distance to polygon for each point of the collection on the 2d Cartesian plane.
- Parameters
polygon – Shapely “Polygon” geometry object.
pxx – List or numpy array of abscissae values of points to calculate the distance from.
pyy – Same structure as
pxx
, but with ordinate values.
- Returns
Numpy array of distances in units of coordinate system. Points that lie inside the polygon have zero distance.
- openquake.hazardlib.geo.utils.triangle_area(e1, e2, e3)[source]¶
Get the area of triangle formed by three vectors.
Parameters are three three-dimensional numpy arrays representing vectors of triangle’s edges in Cartesian space.
- Returns
Float number, the area of the triangle in squared units of coordinates, or numpy array of shape of edges with one dimension less.
Uses Heron formula, see http://mathworld.wolfram.com/HeronsFormula.html.
Module contents¶
Package openquake.hazardlib.geo
contains implementations of different
geographical primitives, such as Point
,
Line
Polygon
and
Mesh
, as well as different
surface
implementations and utility
class NodalPlane
.