1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302
|
# (C) British Crown Copyright 2015 - 2018, Met Office
#
# This file is part of cartopy.
#
# cartopy is free software: you can redistribute it and/or modify it under
# the terms of the GNU Lesser General Public License as published by the
# Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# cartopy 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 Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with cartopy. If not, see <https://www.gnu.org/licenses/>.
#
# cython: embedsignature=True
"""
This module defines the Geodesic class which can interface with the Proj
geodesic functions.
"""
from cpython.mem cimport PyMem_Malloc
from cpython.mem cimport PyMem_Free
import numpy as np
cimport numpy as np
cimport cython
from cython.parallel cimport prange
import shapely.geometry as sgeom
cdef extern from "geodesic.h":
# External imports of Proj4.9 functions
cdef struct geod_geodesic:
pass
ctypedef geod_geodesic* geodesic_t
void geod_init(geodesic_t, double, double)
void geod_direct(geodesic_t, double, double, double, double,
double*, double*, double*) nogil
void geod_inverse(geodesic_t, double, double, double, double,
double*, double*, double*) nogil
cdef class Geodesic:
"""
Define an ellipsoid on which to solve geodesic problems.
"""
cdef geod_geodesic* geod
cdef double radius
cdef double flattening
def __init__(self, radius=6378137.0, flattening=1/298.257223563):
"""
Parameters
----------
radius: float, optional
Equatorial radius (metres). Defaults to the WGS84 semimajor axis
(6378137.0 metres).
flattening: float, optional
Flattening of ellipsoid. Setting flattening = 0 gives a sphere.
Negative flattening gives a prolate ellipsoid. If flattening > 1,
set flattening to 1/flattening.
Defaults to the WGS84 flattening (1/298.257223563).
"""
# This method exists solely for docstrings. The __cinit__ method is
# where the real work happens. Make sure to keep the signatures synced.
def __cinit__(self, radius=6378137.0, flattening=1/298.257223563):
# allocate some memory (filled with random data)
self.geod = <geod_geodesic*> PyMem_Malloc(sizeof(geod_geodesic))
if not self.geod:
raise MemoryError()
geod_init(self.geod, radius, flattening)
self.radius = radius
self.flattening = flattening
def __dealloc__(self):
# Free allocated memory.
PyMem_Free(self.geod)
def __str__(self):
fmt = self.radius, 1/self.flattening
return '<Geodesic: radius=%0.3f, flattening=1/%0.3f>' % fmt
@cython.boundscheck(False)
def direct(self, points, azimuths, distances):
"""
Solve the direct geodesic problem where the length of the geodesic is
specified in terms of distance.
Can accept and broadcast length 1 arguments. For example, given a
single start point and distance, an array of different azimuths can be
supplied to locate multiple endpoints.
Parameters
----------
points: array_like, shape=(n *or* 1, 2)
The starting longitude-latitude point(s) from which to travel.
azimuths: float or array_like with shape=(n, )
List of azimuth values (degrees).
distances : float or array_like with shape(n, )
List of distances values (metres).
Returns
-------
`numpy.ndarray`, shape=(n, 3)
The longitudes, latitudes, and forward azimuths of the located
endpoint(s).
"""
cdef int n_points, i
cdef double[:, :] pts, orig_pts
cdef double[:] azims, dists
# Create numpy arrays from inputs, and ensure correct shape. Note:
# reshape(-1) returns a 1D array from a 0 dimensional array as required
# for broadcasting.
pts = np.array(points, dtype=np.float64).reshape((-1, 2))
azims = np.array(azimuths, dtype=np.float64).reshape(-1)
dists = np.array(distances, dtype=np.float64).reshape(-1)
sizes = [pts.shape[0], azims.size, dists.size]
n_points = max(sizes)
if not all(size in [1, n_points] for size in sizes):
raise ValueError("Inputs must have common length n or length one.")
# Broadcast any length 1 arrays to the correct size.
if pts.shape[0] == 1:
orig_pts = pts
pts = np.empty([n_points, 2], dtype=np.float64)
pts[:, :] = orig_pts
if azims.size == 1:
azims = np.repeat(azims, n_points)
if dists.size == 1:
dists = np.repeat(dists, n_points)
cdef double[:, :] return_pts = np.empty((n_points, 3),
dtype=np.float64)
with nogil:
for i in prange(n_points):
geod_direct(self.geod,
pts[i, 1], pts[i, 0], azims[i], dists[i],
&return_pts[i, 1], &return_pts[i, 0],
&return_pts[i, 2])
return return_pts
@cython.boundscheck(False)
def inverse(self, points, endpoints):
"""
Solve the inverse geodesic problem.
Can accept and broadcast length 1 arguments. For example, given a
single start point, an array of different endpoints can be supplied to
find multiple distances.
Parameters
----------
points: array_like, shape=(n *or* 1, 2)
The starting longitude-latitude point(s) from which to travel.
endpoints: array_like, shape=(n *or* 1, 2)
The longitude-latitude point(s) to travel to.
Returns
-------
`numpy.ndarray`, shape=(n, 3)
The distances, and the (forward) azimuths of the start and end
points.
"""
cdef int n_points, i
cdef double[:, :] pts, epts, orig_pts
# Create numpy arrays from inputs, and ensure correct shape.
points = np.array(points, dtype=np.float64)
endpoints = np.array(endpoints, dtype=np.float64)
if points.ndim > 2 or (points.ndim == 2 and points.shape[1] != 2):
raise ValueError(
'Expecting input points to be (N, 2), got {}'
''.format(points.shape))
pts = points.reshape((-1, 2))
epts = endpoints.reshape((-1, 2))
sizes = [pts.shape[0], epts.shape[0]]
n_points = max(sizes)
if not all(size in [1, n_points] for size in sizes):
raise ValueError("Inputs must have common length n or length one.")
# Broadcast any length 1 arrays to the correct size.
if pts.shape[0] == 1:
orig_pts = pts
pts = np.empty([n_points, 2], dtype=np.float64)
pts[:, :] = orig_pts
if epts.shape[0] == 1:
orig_pts = epts
epts = np.empty([n_points, 2], dtype=np.float64)
epts[:, :] = orig_pts
cdef double[:, :] results = np.empty((n_points, 3), dtype=np.float64)
with nogil:
for i in prange(n_points):
geod_inverse(self.geod, pts[i, 1], pts[i, 0], epts[i, 1],
epts[i, 0], &results[i, 0], &results[i, 1],
&results[i, 2])
return results
def circle(self, double lon, double lat, double radius, int n_samples=180,
endpoint=False):
"""
Find a geodesic circle of given radius at a given point.
Parameters
----------
lon : float
Longitude coordinate of the centre.
lat : float
Latitude coordinate of the centre.
radius : float
The radius of the circle (metres).
n_samples: int, optional
Integer number of sample points of circle.
endpoint: bool, optional
Whether to repeat endpoint at the end of returned array.
Returns
-------
`numpy.ndarray`, shape=(n_samples, 2)
The evenly spaced longitude-latitude points on the circle.
"""
cdef int i
# Put the input arguments into c-typed values.
cdef double[:, :] center = np.array([lon, lat]).reshape((1, 2))
cdef double[:] radius_m = np.asarray(radius).reshape(1)
azimuths = np.linspace(360., 0., n_samples,
endpoint=endpoint).astype(np.double)
return self.direct(center, azimuths, radius_m)[:, 0:2]
def geometry_length(self, geometry):
"""
Return the distance (in physical meters) of the given Shapely geometry.
The geometry is assumed to be in spherical (lon, lat) coordinates.
Parameters
----------
geometry : `shapely.geometry.BaseGeometry`
The Shapely geometry to compute the length of. For polygons, the
exterior length will be calculated. For multi-part geometries, the
sum of the parts will be computed.
"""
result = None
if hasattr(geometry, 'geoms'):
# Multi-geometry.
result = sum(self.geometry_length(geom) for geom in geometry.geoms)
elif hasattr(geometry, 'exterior'):
# Polygon.
result = self.geometry_length(geometry.exterior)
elif hasattr(geometry, 'coords'):
coords = np.array(geometry.coords)
# LinearRings are (N, 2), whereas LineStrings are (2, N).
if not isinstance(geometry, sgeom.LinearRing):
coords = coords.T
result = self.geometry_length(coords)
elif isinstance(geometry, np.ndarray):
coords = geometry
distances, _, _ = np.array(
self.inverse(coords[:-1, :], coords[1:, :]).T)
result = distances.sum()
else:
raise TypeError('Unhandled type {}'.format(geometry.__class__))
return result
|