# -*- coding: utf-8 -*-# vim: tabstop=4 shiftwidth=4 softtabstop=4## Copyright (C) 2012-2025 GEM Foundation## OpenQuake is free software: you can redistribute it and/or modify it# under the terms of the GNU Affero General Public License as published# by the Free Software Foundation, either version 3 of the License, or# (at your option) any later version.## OpenQuake is distributed in the hope that it will be useful,# but WITHOUT ANY WARRANTY; without even the implied warranty of# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the# GNU Affero General Public License for more details.## You should have received a copy of the GNU Affero General Public License# along with OpenQuake. If not, see <http://www.gnu.org/licenses/>."""Module :mod:`openquake.hazardlib.geo.mesh` defines classes :class:`Mesh` andits subclass :class:`RectangularMesh`."""importwarningsimportnumpyfromscipy.spatial.distanceimportcdistimportshapely.geometryimportshapely.opsfromalpha_shapesimportAlpha_Shaperfromopenquake.baselib.generalimportcached_property,gen_slicesfromopenquake.hazardlib.geo.pointimportPointfromopenquake.hazardlib.geoimportgeodeticfromopenquake.hazardlib.geoimportutilsasgeo_utilsF32=numpy.float32
[docs]defsqrt(array):# due to numerical errors an array of positive values can become negative;# for instance: 1 - array([[ 0.99999989, 1.00000001, 1. ]]) =# array([[ 1.08272703e-07, -5.19256105e-09, -3.94126065e-10]])# here we replace the small negative values with zerosarray[array<0]=0returnnumpy.sqrt(array)
[docs]defsurface_to_arrays(surface):""" :param surface: a (Multi)Surface object :returns: a list of S arrays of shape (3, N, M) """ifhasattr(surface,'surfaces'):# multiplanar surfaceslst=[]forsurfinsurface.surfaces:arr=surf.mesh.arrayiflen(arr.shape)==2:# PlanarSurfacearr=arr.reshape(3,1,4)lst.append(arr)returnlstmesh=surface.meshiflen(mesh.lons.shape)==1:# 1D meshshp=(3,1)+mesh.lons.shapeelse:# 2D meshshp=(3,)+mesh.lons.shapereturn[mesh.array.reshape(shp)]
[docs]defcalc_inclination(earth_surface_tangent_normal,tl_normal,tl_area,br_normal,br_area):# inclination calculation# top-left trianglesen=earth_surface_tangent_normal[:-1,:-1]# cosine of inclination of the triangle is scalar product# of vector normal to triangle plane and (normalized) vector# pointing to top left corner of a triangle from earth centerincl_cos=numpy.sum(en*tl_normal,axis=-1).clip(-1.0,1.0)# we calculate average angle using mean of circular quantities# formula: define 2d vector for each triangle where length# of the vector corresponds to triangle's weight (we use triangle# area) and angle is equal to inclination angle. then we calculate# the angle of vector sum of all those vectors and that angle# is the weighted average.xx=numpy.sum(tl_area*incl_cos)# express sine via cosine using Pythagorean trigonometric identity,# this is a bit faster than sin(arccos(incl_cos))yy=numpy.sum(tl_area*sqrt(1-incl_cos*incl_cos))# bottom-right trianglesen=earth_surface_tangent_normal[1:,1:]# we need to clip scalar product values because in some cases# they might exceed range where arccos is defined ([-1, 1])# because of floating point imprecisionincl_cos=numpy.sum(en*br_normal,axis=-1).clip(-1.0,1.0)# weighted angle vectors are calculated independently for top-left# and bottom-right triangles of each cell in a mesh. here we# combine both and finally get the weighted mean anglexx+=numpy.sum(br_area*incl_cos)yy+=numpy.sum(br_area*sqrt(1-incl_cos*incl_cos))inclination=numpy.degrees(numpy.arctan2(yy,xx))returninclination
[docs]defcalc_azimuth(points,along_azimuth,tl_area,br_area):# azimuth calculation is done similar to one for inclination. we also# do separate calculations for top-left and bottom-right triangles# and also combine results using mean of circular quantities approach# unit vector along z axisz_unit=numpy.array([0.0,0.0,1.0])# unit vectors pointing west from each point of the mesh, they define# planes that contain meridian of respective pointnorms_west=geo_utils.normalized(numpy.cross(points+z_unit,points))# unit vectors parallel to planes defined by previous ones. they are# directed from each point to a point lying on z axis on the same# distance from earth centernorms_north=geo_utils.normalized(numpy.cross(points,norms_west))# need to normalize triangles' azimuthal edges because we will project# them on other normals and thus calculate an angle in betweenalong_azimuth=geo_utils.normalized(along_azimuth)# process top-left triangles# here we identify the sign of direction of the triangles' azimuthal# edges: is edge pointing west or east? for finding that we project# those edges to vectors directing to west by calculating scalar# product and get the sign of resulting value: if it is negative# than the resulting azimuth should be negative as top edge is pointing# west.sign=numpy.sign(numpy.sign(numpy.sum(along_azimuth[:-1]*norms_west[:-1,:-1],axis=-1))# we run numpy.sign(numpy.sign(...) + 0.1) to make resulting values# be only either -1 or 1 with zero values (when edge is pointing# strictly north or south) expressed as 1 (which means "don't# change the sign")+0.1)# the length of projection of azimuthal edge on norms_north is cosine# of edge's azimuthaz_cos=numpy.sum(along_azimuth[:-1]*norms_north[:-1,:-1],axis=-1)# use the same approach for finding the weighted mean# as for inclination (see above)xx=numpy.sum(tl_area*az_cos)# the only difference is that azimuth is defined in a range# [0, 360), so we need to have two reference planes and change# sign of projection on one normal to sign of projection to another oneyy=numpy.sum(tl_area*sqrt(1-az_cos*az_cos)*sign)# bottom-right trianglessign=numpy.sign(numpy.sign(numpy.sum(along_azimuth[1:]*norms_west[1:,1:],axis=-1))+0.1)az_cos=numpy.sum(along_azimuth[1:]*norms_north[1:,1:],axis=-1)xx+=numpy.sum(br_area*az_cos)yy+=numpy.sum(br_area*sqrt(1-az_cos*az_cos)*sign)azimuth=numpy.degrees(numpy.arctan2(yy,xx))returnazimuth
[docs]classMesh(object):""" Mesh object represent a collection of points and provides the most efficient way of keeping those collections in memory. :param lons: A numpy array of longitudes. Can be 1D or 2D. :param lats: Numpy array of latitudes. The array must be of the same shape as ``lons``. :param 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 :meth:`from_points_list`. """# Tolerance level to be used in various spatial operations when# approximation is required -- set to 5 meters.# NB: it affects the rjb distance and therefore nearly every calculationDIST_TOLERANCE=0.005@propertydeflons(self):returnself.array[0]@propertydeflats(self):returnself.array[1]@propertydefdepths(self):try:returnself.array[2]exceptIndexError:returnnumpy.zeros(self.shape)def__init__(self,lons,lats,depths=None,round=None):assert((lons.shape==lats.shape)andlen(lons.shape)in(1,2)and(depthsisNoneordepths.shape==lats.shape)),(lons.shape,lats.shape)assertlons.size>0ifdepthsisNone:self.array=numpy.array([lons,lats])else:self.array=numpy.array([lons,lats,depths])ifroundisnotNone:numpy.round(self.array,round,self.array)
[docs]@classmethoddeffrom_coords(cls,coords,sort=True):""" Create a mesh object from a list of 3D coordinates (by sorting them) :params coords: list of coordinates :param sort: flag (default True) :returns: a :class:`Mesh` instance """coords=list(coords)ifsort:coords.sort()iflen(coords[0])==2:# 2D coordinateslons,lats=zip(*coords)depths=Noneelse:# 3D coordinateslons,lats,depths=zip(*coords)depths=numpy.array(depths)returncls(numpy.array(lons),numpy.array(lats),depths)
[docs]@classmethoddeffrom_points_list(cls,points,round=None):""" Create a mesh object from a collection of points. :param point: List of :class:`~openquake.hazardlib.geo.point.Point` objects. :returns: An instance of :class:`Mesh` with one-dimensional arrays of coordinates from ``points``. """lons=numpy.zeros(len(points),dtype=float)lats=lons.copy()depths=lons.copy()foriinrange(len(points)):lons[i]=points[i].longitudelats[i]=points[i].latitudedepths[i]=points[i].depthifnotdepths.any():# all points have zero depth, no need to waste memorydepths=Nonereturncls(lons,lats,depths,round)
@propertydefshape(self):""" Return the shape of this mesh. :returns tuple: The shape of this mesh as (rows, columns) """returnself.array.shape[1:]@cached_propertydefxyz(self):""" :returns: an array of shape (N, 3) with the cartesian coordinates """returngeo_utils.spherical_to_cartesian(self.lons.flat,self.lats.flat,self.depths.flat)def__iter__(self):""" Generate :class:`~openquake.hazardlib.geo.point.Point` objects the mesh is composed of. Coordinates arrays are processed sequentially (as if they were flattened). """lons=self.lons.flatlats=self.lats.flatdepths=self.depths.flatforiinrange(self.lons.size):yieldPoint(lons[i],lats[i],depths[i])def__getitem__(self,item):""" Get a submesh of this mesh. :param item: Indexing is only supported by slices. Those slices are used to cut the portion of coordinates (and depths if it is available) arrays. These arrays are then used for creating a new mesh. :returns: A new object of the same type that borrows a portion of geometry from this mesh (doesn't copy the array, just references it). """ifisinstance(item,int):raiseValueError('You must pass a slice, not an index: %s'%item)lons=self.lons[item]lats=self.lats[item]depths=self.depths[item]returntype(self)(lons,lats,depths)def__len__(self):""" Return the number of points in the mesh. """returnself.lons.sizedef__eq__(self,mesh,tol=1.0E-7):""" Compares the mesh with another returning True if all elements are equal to within the specific tolerance, False otherwise :param mesh: Mesh for comparison as instance of :class: openquake.hazardlib.geo.mesh.Mesh :param float tol: Numerical precision for equality """ifself.shape!=mesh.shape:returnFalseeliflen(self.array)!=len(mesh.array):# 3D vs 2D arraysok=(numpy.allclose(self.array[0],mesh.array[0],atol=tol)andnumpy.allclose(self.array[1],mesh.array[1],atol=tol))iflen(self.array)==2:returnokand(mesh.array[2]==0).all()eliflen(mesh.array)==2:returnokand(self.array[2]==0).all()returnnumpy.allclose(self.array,mesh.array,atol=tol)
[docs]defget_min_distance(self,mesh):""" 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. """# mesh.xyz has shape (N, 3); we split in slices to avoid running out of memory# in the large array of shape len(self)*Ndists=[]forslcingen_slices(0,len(mesh),10_000):dists.append(cdist(self.xyz,mesh.xyz[slc]).min(axis=0))returnnumpy.concatenate(dists)
[docs]defget_closest_points(self,mesh):""" Find closest point of this mesh for each point in the other mesh :returns: :class:`Mesh` object of the same shape as `mesh` with closest points from this one at respective indices. """min_idx=cdist(self.xyz,mesh.xyz).argmin(axis=0)# lose shapeifhasattr(mesh,'shape'):min_idx=min_idx.reshape(mesh.shape)lons=self.lons.take(min_idx)lats=self.lats.take(min_idx)deps=self.depths.take(min_idx)returnMesh(lons,lats,deps)
[docs]defget_distance_matrix(self):""" 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). Uses :func:`openquake.hazardlib.geo.geodetic.geodetic_distance`. """assertself.lons.ndim==1distances=geodetic.geodetic_distance(self.lons.reshape(self.lons.shape+(1,)),self.lats.reshape(self.lats.shape+(1,)),self.lons,self.lats)returndistances
def_get_proj_convex_hull(self):""" Create a projection centered in the center of this mesh and define a convex polygon in that projection, enveloping all the points of the mesh. :returns: Tuple of two items: projection function and shapely 2d polygon. Note that the result geometry can be line or point depending on number of points in the mesh and their arrangement. """# create a projection centered in the center of points collectionproj=geo_utils.OrthographicProjection.from_(self.lons.flatten(),self.lats.flatten())# project all the points and create a shapely multipoint object.# need to copy an array because otherwise shapely misinterprets itcoords=proj(self.lons.flatten(),self.lats.flatten()).Tmultipoint=shapely.geometry.MultiPoint(coords)# create a 2d polygon from a convex hull around that multipointreturnproj,multipoint.convex_hull
[docs]defget_joyner_boore_distance(self,mesh,unstructured=False):""" Compute and return Joyner-Boore distance to each point of ``mesh``. Point's depth is ignored. See :meth:`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. """# we perform a hybrid calculation (geodetic mesh-to-mesh distance# and distance on the projection plane for close points). first,# we find the closest geodetic distance for each point of target# mesh to this one. in general that distance is greater than# the exact distance to enclosing polygon of this mesh and it# depends on mesh spacing. but the difference can be neglected# if calculated geodetic distance is over some threshold.# get the highest slice from the 3D meshdistances=geodetic.min_geodetic_distance((self.lons,self.lats),(mesh.lons,mesh.lats))# here we find the points for which calculated mesh-to-mesh# distance is below a threshold. this threshold is arbitrary:# lower values increase the maximum possible error, higher# values reduce the efficiency of that filtering. the maximum# error is equal to the maximum difference between a distance# from site to two adjacent points of the mesh and distance# from site to the line connecting them. thus the error is# a function of distance threshold and mesh spacing. the error# is maximum when the site lies on a perpendicular to the line# connecting points of the mesh and that passes the middle# point between them. the error then can be calculated as# ``err = trsh - d = trsh - \sqrt(trsh^2 - (ms/2)^2)``, where# ``trsh`` and ``d`` are distance to mesh points (the one# we found on the previous step) and distance to the line# connecting them (the actual distance) and ``ms`` is mesh# spacing. the threshold of 40 km gives maximum error of 314# meters for meshes with spacing of 10 km and 5.36 km for# meshes with spacing of 40 km. if mesh spacing is over# ``(trsh / \sqrt(2)) * 2`` then points lying in the middle# of mesh cells (that is inside the polygon) will be filtered# out by the threshold and have positive distance instead of 0.# so for threshold of 40 km mesh spacing should not be more# than 56 km (typical values are 5 to 10 km).idxs=(distances<40).nonzero()[0]# indices on the first dimensionifnotlen(idxs):# no point is close enough, return distances as they arereturndistances# for all the points that are closer than the threshold we need# to recalculate the distance and set it to zero, if point falls# inside the enclosing polygon of the mesh. for doing that we# project both this mesh and the points of the second mesh--selected# by distance threshold--to the same Cartesian space, define# minimum shapely polygon enclosing the mesh and calculate point# to polygon distance, which gives the most accurate value# of distance in km (and that value is zero for points inside# the polygon).ifunstructured:proj=geo_utils.OrthographicProjection.from_(self.lons,self.lats)# Points at distances lower than 40 kmmesh_xx,mesh_yy=proj(mesh.lons[idxs],mesh.lats[idxs])# Points representing the surface of the rupturesfc_xx,sfc_yy=proj(self.lons.flatten(),self.lats.flatten())points=[(lo,la)forlo,lainzip(sfc_xx,sfc_yy)]shaper=Alpha_Shaper(points)_alpha_opt,polygon=shaper.optimize()else:proj,polygon=self._get_proj_enclosing_polygon()ifnotisinstance(polygon,shapely.geometry.Polygon):# either line or point is our enclosing polygon. draw# a square with side of 10 m around in order to have# a proper polygon instead.polygon=polygon.buffer(self.DIST_TOLERANCE,1)mesh_xx,mesh_yy=proj(mesh.lons[idxs],mesh.lats[idxs])# replace geodetic distance values for points-closer-than-the-threshold# by more accurate point-to-polygon distance values.distances[idxs]=geo_utils.point_to_polygon_distance(polygon,mesh_xx,mesh_yy)returndistances
def_get_proj_enclosing_polygon(self):""" See :meth:`Mesh._get_proj_enclosing_polygon`. :class:`RectangularMesh` contains an information about relative positions of points, so it allows to define the minimum polygon, containing the projection of the mesh, which doesn't necessarily have to be convex (in contrast to :class:`Mesh` implementation). :returns: Same structure as :meth:`Mesh._get_proj_convex_hull`. """ifself.lons.size<4:# the mesh doesn't contain even a single cellreturnself._get_proj_convex_hull()iflen(self.lons.shape)==1:# 1D meshlons=self.lons.reshape(len(self.lons),1)lats=self.lats.reshape(len(self.lats),1)else:# 2D meshlons=self.lons.Tlats=self.lats.Tproj=geo_utils.OrthographicProjection.from_(lons,lats)mesh2d=proj(lons,lats).Tlines=iter(mesh2d)# we iterate over horizontal stripes, keeping the "previous"# line of points. we keep it reversed, such that together# with the current line they define the sequence of points# around the stripe.prev_line=next(lines)[::-1]polygons=[]fori,lineinenumerate(lines):coords=numpy.concatenate((prev_line,line,prev_line[0:1]))# create the shapely polygon object from the stripe# coordinates and simplify it (remove redundant points,# if there are any lying on the straight line).stripe=shapely.geometry.LineString(coords).simplify(self.DIST_TOLERANCE)# ignore RuntimeWarning: divide by zero in .buffer# since it is raised by shapely and there is nothing we can do# see https://github.com/gem/oq-engine/issues/10009withwarnings.catch_warnings():warnings.filterwarnings('ignore',category=RuntimeWarning)stripe=stripe.buffer(self.DIST_TOLERANCE,2)polygons.append(shapely.geometry.Polygon(stripe.exterior))prev_line=line[::-1]# ignore RuntimeWarning: divide by zero encountered in unary_union# since it is raised by shapely and there is nothing we can do# see https://github.com/gem/oq-engine/issues/10009withwarnings.catch_warnings():warnings.filterwarnings('ignore',category=RuntimeWarning)polygon=shapely.ops.unary_union(polygons).simplify(self.DIST_TOLERANCE)# debug_plot(polygons)returnproj,polygon
[docs]defget_convex_hull(self):""" Get a convex polygon object that contains projections of all the points of the mesh. :returns: Instance of :class:`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. """proj,polygon2d=self._get_proj_convex_hull()# if mesh had only one point, the convex hull is a point. if there# were two, it is a line string. we need to return a convex polygon# object, so extend that area-less geometries by some arbitrarily# small distance.ifisinstance(polygon2d,(shapely.geometry.LineString,shapely.geometry.Point)):polygon2d=polygon2d.buffer(self.DIST_TOLERANCE,1)# avoid circular importsfromopenquake.hazardlib.geo.polygonimportPolygonreturnPolygon._from_2d(polygon2d,proj)
[docs]defreduce(self,n):""" Reduce the mesh by `n` times """returnMesh(reduce1d(self.lons,n),reduce1d(self.lats,n),reduce1d(self.depths,n))
[docs]classRectangularMesh(Mesh):""" A specification of :class:`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. """def__init__(self,lons,lats,depths=None,round=None):super().__init__(lons,lats,depths,round)assertlons.ndim==2
[docs]@classmethoddeffrom_points_list(cls,points,round=None):""" Create a rectangular mesh object from a list of lists of points. Lists in a list are supposed to have the same length. :param point: List of lists of :class:`~openquake.hazardlib.geo.point.Point` objects. """assertpointsisnotNoneandlen(points)>0andlen(points[0])>0, \
'list of at least one non-empty list of points is required'lons=numpy.zeros((len(points),len(points[0])),dtype=float)lats=lons.copy()depths=lons.copy()num_cols=len(points[0])fori,rowinenumerate(points):assertlen(row)==num_cols, \
'lists of points are not of uniform length'forj,pointinenumerate(row):lons[i,j]=point.longitudelats[i,j]=point.latitudedepths[i,j]=point.depthifnotdepths.any():depths=Nonereturncls(lons,lats,depths,round)
[docs]defget_middle_point(self):""" Return the middle point of the mesh. :returns: An instance of :class:`~openquake.hazardlib.geo.point.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. """num_rows,num_cols=self.lons.shapemid_row=num_rows//2depth=0ifnum_rows&1==1:# there are odd number of rowsmid_col=num_cols//2ifnum_cols&1==1:# odd number of columns, we can easily take# the middle pointdepth=self.depths[mid_row,mid_col]returnPoint(self.lons[mid_row,mid_col],self.lats[mid_row,mid_col],depth)else:# even number of columns, need to take two middle# points on the middle rowlon1,lon2=self.lons[mid_row,mid_col-1:mid_col+1]lat1,lat2=self.lats[mid_row,mid_col-1:mid_col+1]depth1=self.depths[mid_row,mid_col-1]depth2=self.depths[mid_row,mid_col]else:# there are even number of rows. take the row just above# and the one just below the middle and find middle point# of eachsubmesh1=self[mid_row-1:mid_row]submesh2=self[mid_row:mid_row+1]p1,p2=submesh1.get_middle_point(),submesh2.get_middle_point()lon1,lat1,depth1=p1.longitude,p1.latitude,p1.depthlon2,lat2,depth2=p2.longitude,p2.latitude,p2.depth# we need to find the middle between two pointsdepth=(depth1+depth2)/2.0lon,lat=geo_utils.get_middle_point(lon1,lat1,lon2,lat2)returnPoint(lon,lat,depth)
[docs]defget_mean_inclination_and_azimuth(self):""" 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. """assert1notinself.lons.shape,("inclination and azimuth are only defined for mesh of more than ""one row and more than one column of points")assert((self.depths[1:]-self.depths[:-1])>=0).all(),("get_mean_inclination_and_azimuth() requires next mesh row ""to be not shallower than the previous one")points,along_azimuth,updip,diag=self.triangulate()# define planes that are perpendicular to each point's vector# as normals to those planesearth_surface_tangent_normal=geo_utils.normalized(points)# calculating triangles' area and normals for top-left trianglese1=along_azimuth[:-1]e2=updip[:,:-1]tl_area=geo_utils.triangle_area(e1,e2,diag)tl_normal=geo_utils.normalized(numpy.cross(e1,e2))# ... and bottom-right trianglese1=along_azimuth[1:]e2=updip[:,1:]br_area=geo_utils.triangle_area(e1,e2,diag)br_normal=geo_utils.normalized(numpy.cross(e1,e2))if(self.depths==0).all():# mesh is on earth surface, inclination is zeroinclination=0else:inclination=calc_inclination(earth_surface_tangent_normal,tl_normal,tl_area,br_normal,br_area)azimuth=calc_azimuth(points,along_azimuth,tl_area,br_area)ifazimuth<0:azimuth+=360ifinclination>90:# average inclination is over 90 degree, that means that we need# to reverse azimuthal direction in order for inclination to be# in range [0, 90]inclination=180-inclinationazimuth=(azimuth+180)%360returninclination,azimuth
[docs]defget_cell_dimensions(self):""" 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. """points,along_azimuth,updip,diag=self.triangulate()top=along_azimuth[:-1]left=updip[:,:-1]tl_area=geo_utils.triangle_area(top,left,diag)top_length=numpy.sqrt(numpy.sum(top*top,axis=-1))left_length=numpy.sqrt(numpy.sum(left*left,axis=-1))bottom=along_azimuth[1:]right=updip[:,1:]br_area=geo_utils.triangle_area(bottom,right,diag)bottom_length=numpy.sqrt(numpy.sum(bottom*bottom,axis=-1))right_length=numpy.sqrt(numpy.sum(right*right,axis=-1))cell_area=tl_area+br_areatl_center=(points[:-1,:-1]+points[:-1,1:]+points[1:,:-1])/3br_center=(points[:-1,1:]+points[1:,:-1]+points[1:,1:])/3cell_center=((tl_center*tl_area.reshape(tl_area.shape+(1,))+br_center*br_area.reshape(br_area.shape+(1,)))/cell_area.reshape(cell_area.shape+(1,)))cell_length=((top_length*tl_area+bottom_length*br_area)/cell_area)cell_width=((left_length*tl_area+right_length*br_area)/cell_area)returncell_center,cell_length,cell_width,cell_area
[docs]deftriangulate(self):""" 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. """points=geo_utils.spherical_to_cartesian(self.lons,self.lats,self.depths)# triangulate the mesh by defining vectors of triangles edges:# →along_azimuth=points[:,1:]-points[:,:-1]# ↑updip=points[:-1]-points[1:]# ↗diag=points[:-1,1:]-points[1:,:-1]returnpoints,along_azimuth,updip,diag
[docs]defget_mean_width(self):""" 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. """assert1notinself.lons.shape,("mean width is only defined for mesh of more than ""one row and more than one column of points")_,cell_length,cell_width,cell_area=self.get_cell_dimensions()# compute widths along each mesh columnwidths=numpy.sum(cell_width,axis=0)# compute (weighted) mean cell length along each mesh columncolumn_areas=numpy.sum(cell_area,axis=0)mean_cell_lengths=numpy.sum(cell_length*cell_area,axis=0)/ \
column_areas# compute and return weighted meanreturnnumpy.sum(widths*mean_cell_lengths)/ \
numpy.sum(mean_cell_lengths)
[docs]defreduce(self,n):""" Reduce the mesh by `n^2` times """returnRectangularMesh(reduce2d(self.lons,n),reduce2d(self.lats,n),reduce2d(self.depths,n))