1 """
2 Cython wrapper to provide python interfaces to
3 PROJ.4 (https://github.com/OSGeo/proj.4/wiki) functions.
4
5 Performs cartographic transformations and geodetic computations.
6
7 The Proj class can convert from geographic (longitude,latitude)
8 to native map projection (x,y) coordinates and vice versa, or
9 from one map projection coordinate system directly to another.
10 The module variable pj_list is a dictionary containing all the
11 available projections and their descriptions.
12
13 The Geod class can perform forward and inverse geodetic, or
14 Great Circle, computations. The forward computation involves
15 determining latitude, longitude and back azimuth of a terminus
16 point given the latitude and longitude of an initial point, plus
17 azimuth and distance. The inverse computation involves
18 determining the forward and back azimuths and distance given the
19 latitudes and longitudes of an initial and terminus point.
20
21 Input coordinates can be given as python arrays, lists/tuples,
22 scalars or numpy/Numeric/numarray arrays. Optimized for objects
23 that support the Python buffer protocol (regular python and
24 numpy array objects).
25
26 Download: http://python.org/pypi/pyproj
27
28 Requirements: python 2.4 or higher.
29
30 Example scripts are in 'test' subdirectory of source distribution.
31 The 'test()' function will run the examples in the docstrings.
32
33 Contact: Jeffrey Whitaker <jeffrey.s.whitaker@noaa.gov
34
35 copyright (c) 2006 by Jeffrey Whitaker.
36
37 Permission to use, copy, modify, and distribute this software
38 and its documentation for any purpose and without fee is hereby
39 granted, provided that the above copyright notice appear in all
40 copies and that both the copyright notice and this permission
41 notice appear in supporting documentation. THE AUTHOR DISCLAIMS
42 ALL WARRANTIES WITH REGARD TO THIS SOFTWARE, INCLUDING ALL
43 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS. IN NO EVENT
44 SHALL THE AUTHOR BE LIABLE FOR ANY SPECIAL, INDIRECT OR
45 CONSEQUENTIAL DAMAGES OR ANY DAMAGES WHATSOEVER RESULTING FROM
46 LOSS OF USE, DATA OR PROFITS, WHETHER IN AN ACTION OF CONTRACT,
47 NEGLIGENCE OR OTHER TORTIOUS ACTION, ARISING OUT OF OR IN
48 CONNECTION WITH THE USE OR PERFORMANCE OF THIS SOFTWARE. """
49
50 import sys
51 from pyproj import _proj
52 from pyproj.datadir import pyproj_datadir
53 __version__ = _proj.__version__
54 set_datapath = _proj.set_datapath
55 from array import array
56 import os, math
57
58
59
60 if sys.version_info[0] == 2:
61 string_types = (basestring,)
62 else:
63 string_types = (str,)
64
65 pj_list={
66 'aea': "Albers Equal Area",
67 'aeqd': "Azimuthal Equidistant",
68 'airy': "Airy",
69 'aitoff': "Aitoff",
70 'alsk': "Mod. Stererographics of Alaska",
71 'apian': "Apian Globular I",
72 'august': "August Epicycloidal",
73 'bacon': "Bacon Globular",
74 'bipc': "Bipolar conic of western hemisphere",
75 'boggs': "Boggs Eumorphic",
76 'bonne': "Bonne (Werner lat_1=90)",
77 'cass': "Cassini",
78 'cc': "Central Cylindrical",
79 'cea': "Equal Area Cylindrical",
80 'chamb': "Chamberlin Trimetric",
81 'collg': "Collignon",
82 'crast': "Craster Parabolic (Putnins P4)",
83 'denoy': "Denoyer Semi-Elliptical",
84 'eck1': "Eckert I",
85 'eck2': "Eckert II",
86 'eck3': "Eckert III",
87 'eck4': "Eckert IV",
88 'eck5': "Eckert V",
89 'eck6': "Eckert VI",
90 'eqc': "Equidistant Cylindrical (Plate Caree)",
91 'eqdc': "Equidistant Conic",
92 'etmerc': "Extended Transverse Mercator" ,
93 'euler': "Euler",
94 'fahey': "Fahey",
95 'fouc': "Foucaut",
96 'fouc_s': "Foucaut Sinusoidal",
97 'gall': "Gall (Gall Stereographic)",
98 'geocent': "Geocentric",
99 'geos': "Geostationary Satellite View",
100 'gins8': "Ginsburg VIII (TsNIIGAiK)",
101 'gn_sinu': "General Sinusoidal Series",
102 'gnom': "Gnomonic",
103 'goode': "Goode Homolosine",
104 'gs48': "Mod. Stererographics of 48 U.S.",
105 'gs50': "Mod. Stererographics of 50 U.S.",
106 'hammer': "Hammer & Eckert-Greifendorff",
107 'hatano': "Hatano Asymmetrical Equal Area",
108 'healpix': "HEALPix",
109 'rhealpix': "rHEALPix",
110 'igh': "Interrupted Goode Homolosine",
111 'imw_p': "International Map of the World Polyconic",
112 'isea': "Icosahedral Snyder Equal Area",
113 'kav5': "Kavraisky V",
114 'kav7': "Kavraisky VII",
115 'krovak': "Krovak",
116 'labrd': "Laborde",
117 'laea': "Lambert Azimuthal Equal Area",
118 'lagrng': "Lagrange",
119 'larr': "Larrivee",
120 'lask': "Laskowski",
121 'lonlat': "Lat/long (Geodetic)",
122 'latlon': "Lat/long (Geodetic alias)",
123 'latlong': "Lat/long (Geodetic alias)",
124 'longlat': "Lat/long (Geodetic alias)",
125 'lcc': "Lambert Conformal Conic",
126 'lcca': "Lambert Conformal Conic Alternative",
127 'leac': "Lambert Equal Area Conic",
128 'lee_os': "Lee Oblated Stereographic",
129 'loxim': "Loximuthal",
130 'lsat': "Space oblique for LANDSAT",
131 'mbt_s': "McBryde-Thomas Flat-Polar Sine",
132 'mbt_fps': "McBryde-Thomas Flat-Pole Sine (No. 2)",
133 'mbtfpp': "McBride-Thomas Flat-Polar Parabolic",
134 'mbtfpq': "McBryde-Thomas Flat-Polar Quartic",
135 'mbtfps': "McBryde-Thomas Flat-Polar Sinusoidal",
136 'merc': "Mercator",
137 'mil_os': "Miller Oblated Stereographic",
138 'mill': "Miller Cylindrical",
139 'moll': "Mollweide",
140 'murd1': "Murdoch I",
141 'murd2': "Murdoch II",
142 'murd3': "Murdoch III",
143 'natearth': "Natural Earth",
144 'nell': "Nell",
145 'nell_h': "Nell-Hammer",
146 'nicol': "Nicolosi Globular",
147 'nsper': "Near-sided perspective",
148 'nzmg': "New Zealand Map Grid",
149 'ob_tran': "General Oblique Transformation",
150 'ocea': "Oblique Cylindrical Equal Area",
151 'oea': "Oblated Equal Area",
152 'omerc': "Oblique Mercator",
153 'ortel': "Ortelius Oval",
154 'ortho': "Orthographic",
155 'pconic': "Perspective Conic",
156 'poly': "Polyconic (American)",
157 'putp1': "Putnins P1",
158 'putp2': "Putnins P2",
159 'putp3': "Putnins P3",
160 'putp3p': "Putnins P3'",
161 'putp4p': "Putnins P4'",
162 'putp5': "Putnins P5",
163 'putp5p': "Putnins P5'",
164 'putp6': "Putnins P6",
165 'putp6p': "Putnins P6'",
166 'qua_aut': "Quartic Authalic",
167 'robin': "Robinson",
168 'rouss': "Roussilhe Stereographic",
169 'rpoly': "Rectangular Polyconic",
170 'sinu': "Sinusoidal (Sanson-Flamsteed)",
171 'somerc': "Swiss. Obl. Mercator",
172 'stere': "Stereographic",
173 'sterea': "Oblique Stereographic Alternative",
174 'gstmerc': "Gauss-Schreiber Transverse Mercator (aka Gauss-Laborde Reunion)",
175 'tcc': "Transverse Central Cylindrical",
176 'tcea': "Transverse Cylindrical Equal Area",
177 'tissot': "Tissot Conic",
178 'tmerc': "Transverse Mercator",
179 'tpeqd': "Two Point Equidistant",
180 'tpers': "Tilted perspective",
181 'ups': "Universal Polar Stereographic",
182 'urm5': "Urmaev V",
183 'urmfps': "Urmaev Flat-Polar Sinusoidal",
184 'utm': "Universal Transverse Mercator (UTM)",
185 'vandg': "van der Grinten (I)",
186 'vandg2': "van der Grinten II",
187 'vandg3': "van der Grinten III",
188 'vandg4': "van der Grinten IV",
189 'vitk1': "Vitkovsky I",
190 'wag1': "Wagner I (Kavraisky VI)",
191 'wag2': "Wagner II",
192 'wag3': "Wagner III",
193 'wag4': "Wagner IV",
194 'wag5': "Wagner V",
195 'wag6': "Wagner VI",
196 'wag7': "Wagner VII",
197 'weren': "Werenskiold I",
198 'wink1': "Winkel I",
199 'wink2': "Winkel II",
200 'wintri': "Winkel Tripel"}
201
202 pj_ellps={
203 "MERIT": {'a':6378137.0,'rf':298.257,'description':"MERIT 1983"},
204 "SGS85": {'a':6378136.0,'rf':298.257,'description':"Soviet Geodetic System 85"},
205 "GRS80": {'a':6378137.0,'rf':298.257222101,'description':"GRS 1980(IUGG, 1980)"},
206 "IAU76": {'a':6378140.0,'rf':298.257,'description':"IAU 1976"},
207 "airy": {'a':6377563.396,'b':6356256.910,'description':"Airy 1830"},
208 "APL4.9": {'a':6378137.0,'rf':298.25,'description':"Appl. Physics. 1965"},
209 "NWL9D": {'a':6378145.0,'rf':298.25,'description':" Naval Weapons Lab., 1965"},
210 "mod_airy": {'a':6377340.189,'b':6356034.446,'description':"Modified Airy"},
211 "andrae": {'a':6377104.43,'rf':300.0,'description':"Andrae 1876 (Den., Iclnd.)"},
212 "aust_SA": {'a':6378160.0,'rf':298.25,'description':"Australian Natl & S. Amer. 1969"},
213 "GRS67": {'a':6378160.0,'rf':298.2471674270,'description':"GRS 67(IUGG 1967)"},
214 "bessel": {'a':6377397.155,'rf':299.1528128,'description':"Bessel 1841"},
215 "bess_nam": {'a':6377483.865,'rf':299.1528128,'description':"Bessel 1841 (Namibia)"},
216 "clrk66": {'a':6378206.4,'b':6356583.8,'description':"Clarke 1866"},
217 "clrk80": {'a':6378249.145,'rf':293.4663,'description':"Clarke 1880 mod."},
218 "CPM": {'a':6375738.7,'rf':334.29,'description':"Comm. des Poids et Mesures 1799"},
219 "delmbr": {'a':6376428.,'rf':311.5,'description':"Delambre 1810 (Belgium)"},
220 "engelis": {'a':6378136.05,'rf':298.2566,'description':"Engelis 1985"},
221 "evrst30": {'a':6377276.345,'rf':300.8017,'description':"Everest 1830"},
222 "evrst48": {'a':6377304.063,'rf':300.8017,'description':"Everest 1948"},
223 "evrst56": {'a':6377301.243,'rf':300.8017,'description':"Everest 1956"},
224 "evrst69": {'a':6377295.664,'rf':300.8017,'description':"Everest 1969"},
225 "evrstSS": {'a':6377298.556,'rf':300.8017,'description':"Everest (Sabah & Sarawak)"},
226 "fschr60": {'a':6378166.,'rf':298.3,'description':"Fischer (Mercury Datum) 1960"},
227 "fschr60m": {'a':6378155.,'rf':298.3,'description':"Modified Fischer 1960"},
228 "fschr68": {'a':6378150.,'rf':298.3,'description':"Fischer 1968"},
229 "helmert": {'a':6378200.,'rf':298.3,'description':"Helmert 1906"},
230 "hough": {'a':6378270.0,'rf':297.,'description':"Hough"},
231 "intl": {'a':6378388.0,'rf':297.,'description':"International 1909 (Hayford)"},
232 "krass": {'a':6378245.0,'rf':298.3,'description':"Krassovsky, 1942"},
233 "kaula": {'a':6378163.,'rf':298.24,'description':"Kaula 1961"},
234 "lerch": {'a':6378139.,'rf':298.257,'description':"Lerch 1979"},
235 "mprts": {'a':6397300.,'rf':191.,'description':"Maupertius 1738"},
236 "new_intl": {'a':6378157.5,'b':6356772.2,'description':"New International 1967"},
237 "plessis": {'a':6376523.,'b':6355863.,'description':"Plessis 1817 (France)"},
238 "SEasia": {'a':6378155.0,'b':6356773.3205,'description':"Southeast Asia"},
239 "walbeck": {'a':6376896.0,'b':6355834.8467,'description':"Walbeck"},
240 "WGS60": {'a':6378165.0,'rf':298.3,'description':"WGS 60"},
241 "WGS66": {'a':6378145.0,'rf':298.25,'description':"WGS 66"},
242 "WGS72": {'a':6378135.0,'rf':298.26,'description':"WGS 72"},
243 "WGS84": {'a':6378137.0,'rf':298.257223563,'description':"WGS 84"},
244 "sphere": {'a':6370997.0,'b':6370997.0,'description':"Normal Sphere"},
245 }
246
247
248
249
250
251 set_datapath(pyproj_datadir)
252
253 -class Proj(_proj.Proj):
254 """
255 performs cartographic transformations (converts from
256 longitude,latitude to native map projection x,y coordinates and
257 vice versa) using proj (https://github.com/OSGeo/proj.4/wiki).
258
259 A Proj class instance is initialized with proj map projection
260 control parameter key/value pairs. The key/value pairs can
261 either be passed in a dictionary, or as keyword arguments,
262 or as a proj4 string (compatible with the proj command). See
263 http://www.remotesensing.org/geotiff/proj_list for examples of
264 key/value pairs defining different map projections.
265
266 Calling a Proj class instance with the arguments lon, lat will
267 convert lon/lat (in degrees) to x/y native map projection
268 coordinates (in meters). If optional keyword 'inverse' is True
269 (default is False), the inverse transformation from x/y to
270 lon/lat is performed. If optional keyword 'radians' is True
271 (default is False) lon/lat are interpreted as radians instead of
272 degrees. If optional keyword 'errcheck' is True (default is
273 False) an exception is raised if the transformation is invalid.
274 If errcheck=False and the transformation is invalid, no
275 exception is raised and 1.e30 is returned. If the optional keyword
276 'preserve_units' is True, the units in map projection coordinates
277 are not forced to be meters.
278
279 Works with numpy and regular python array objects, python
280 sequences and scalars.
281 """
282
283 - def __new__(self, projparams=None, preserve_units=False, **kwargs):
284 """
285 initialize a Proj class instance.
286
287 Proj4 projection control parameters must either be given in a
288 dictionary 'projparams' or as keyword arguments. See the proj
289 documentation (https://github.com/OSGeo/proj.4/wiki) for more information
290 about specifying projection parameters.
291
292 Example usage:
293
294 >>> from pyproj import Proj
295 >>> p = Proj(proj='utm',zone=10,ellps='WGS84') # use kwargs
296 >>> x,y = p(-120.108, 34.36116666)
297 >>> 'x=%9.3f y=%11.3f' % (x,y)
298 'x=765975.641 y=3805993.134'
299 >>> 'lon=%8.3f lat=%5.3f' % p(x,y,inverse=True)
300 'lon=-120.108 lat=34.361'
301 >>> # do 3 cities at a time in a tuple (Fresno, LA, SF)
302 >>> lons = (-119.72,-118.40,-122.38)
303 >>> lats = (36.77, 33.93, 37.62 )
304 >>> x,y = p(lons, lats)
305 >>> 'x: %9.3f %9.3f %9.3f' % x
306 'x: 792763.863 925321.537 554714.301'
307 >>> 'y: %9.3f %9.3f %9.3f' % y
308 'y: 4074377.617 3763936.941 4163835.303'
309 >>> lons, lats = p(x, y, inverse=True) # inverse transform
310 >>> 'lons: %8.3f %8.3f %8.3f' % lons
311 'lons: -119.720 -118.400 -122.380'
312 >>> 'lats: %8.3f %8.3f %8.3f' % lats
313 'lats: 36.770 33.930 37.620'
314 >>> p2 = Proj('+proj=utm +zone=10 +ellps=WGS84') # use proj4 string
315 >>> x,y = p2(-120.108, 34.36116666)
316 >>> 'x=%9.3f y=%11.3f' % (x,y)
317 'x=765975.641 y=3805993.134'
318 >>> p = Proj(init="epsg:32667")
319 >>> 'x=%12.3f y=%12.3f (meters)' % p(-114.057222, 51.045)
320 'x=-1783486.760 y= 6193833.196 (meters)'
321 >>> p = Proj("+init=epsg:32667",preserve_units=True)
322 >>> 'x=%12.3f y=%12.3f (feet)' % p(-114.057222, 51.045)
323 'x=-5851322.810 y=20320934.409 (feet)'
324 >>> p = Proj(proj='hammer') # hammer proj and inverse
325 >>> x,y = p(-30,40)
326 >>> 'x=%12.3f y=%12.3f' % (x,y)
327 'x=-2711575.083 y= 4395506.619'
328 >>> lon,lat = p(x,y,inverse=True)
329 >>> 'lon=%9.3f lat=%9.3f (degrees)' % (lon,lat)
330 'lon= -30.000 lat= 40.000 (degrees)'
331 """
332
333 if projparams is None:
334 if len(kwargs) == 0:
335 raise RuntimeError('no projection control parameters specified')
336 else:
337 projstring = _dict2string(kwargs)
338 elif isinstance(projparams, string_types):
339
340 projstring = projparams
341 else:
342 projstring = _dict2string(projparams)
343
344 if not projstring.count('+units=') and not preserve_units:
345 projstring = '+units=m '+projstring
346 else:
347 kvpairs = []
348 for kvpair in projstring.split():
349 if kvpair.startswith('+units') and not preserve_units:
350 k,v = kvpair.split('=')
351 kvpairs.append(k+'=m ')
352 else:
353 kvpairs.append(kvpair+' ')
354 projstring = ''.join(kvpairs)
355
356
357 projstring = projstring.replace('EPSG','epsg')
358 return _proj.Proj.__new__(self, projstring)
359
361
362 """
363 Calling a Proj class instance with the arguments lon, lat will
364 convert lon/lat (in degrees) to x/y native map projection
365 coordinates (in meters). If optional keyword 'inverse' is True
366 (default is False), the inverse transformation from x/y to
367 lon/lat is performed. If optional keyword 'radians' is True
368 (default is False) the units of lon/lat are radians instead of
369 degrees. If optional keyword 'errcheck' is True (default is
370 False) an exception is raised if the transformation is invalid.
371 If errcheck=False and the transformation is invalid, no
372 exception is raised and 1.e30 is returned.
373
374 Inputs should be doubles (they will be cast to doubles if they
375 are not, causing a slight performance hit).
376
377 Works with numpy and regular python array objects, python
378 sequences and scalars, but is fastest for array objects.
379 """
380 inverse = kw.get('inverse', False)
381 radians = kw.get('radians', False)
382 errcheck = kw.get('errcheck', False)
383
384
385
386
387
388
389
390
391 lon, lat = args
392
393 inx, xisfloat, xislist, xistuple = _copytobuffer(lon)
394 iny, yisfloat, yislist, yistuple = _copytobuffer(lat)
395
396 if inverse:
397 _proj.Proj._inv(self, inx, iny, radians=radians, errcheck=errcheck)
398 else:
399 _proj.Proj._fwd(self, inx, iny, radians=radians, errcheck=errcheck)
400
401 outx = _convertback(xisfloat,xislist,xistuple,inx)
402 outy = _convertback(yisfloat,yislist,xistuple,iny)
403 return outx, outy
404
406 """returns an equivalent Proj in the corresponding lon/lat
407 coordinates. (see pj_latlong_from_proj() in the Proj.4 C API)"""
408 return _proj.Proj.to_latlong(self)
409
411 """returns True if projection in geographic (lon/lat) coordinates"""
412 return _proj.Proj.is_latlong(self)
413
415 """returns True if projection in geocentric (x/y) coordinates"""
416 return _proj.Proj.is_geocent(self)
417
513
515 try:
516
517 return array('d',(float(x),)),True,False,False
518 except:
519 raise TypeError('input must be an array, list, tuple or scalar')
520
522 """
523 return a copy of x as an object that supports the python Buffer
524 API (python array if input is float, list or tuple, numpy array
525 if input is a numpy array). returns copyofx, isfloat, islist,
526 istuple (islist is True if input is a list, istuple is true if
527 input is a tuple, isfloat is true if input is a float).
528 """
529
530 isfloat = False; islist = False; istuple = False
531
532
533 if hasattr(x,'shape'):
534 if x.shape == ():
535 return _copytobuffer_return_scalar(x)
536 else:
537 try:
538
539
540
541 x.dtype.char
542
543
544
545 inx = x.copy(order="C").astype('d')
546
547 return inx,False,False,False
548 except:
549 try:
550
551
552 x.typecode()
553 inx = x.astype('d')
554
555 return inx,False,False,False
556 except:
557 raise TypeError('input must be an array, list, tuple or scalar')
558 else:
559
560 if hasattr(x, 'typecode'):
561
562 inx = array('d',x)
563
564
565 elif type(x) == list:
566 inx = array('d',x)
567 islist = True
568
569 elif type(x) == tuple:
570 inx = array('d',x)
571 istuple = True
572
573 else:
574 return _copytobuffer_return_scalar(x)
575 return inx,isfloat,islist,istuple
576
578
579 if isfloat:
580 return inx[0]
581 elif islist:
582 return inx.tolist()
583 elif istuple:
584 return tuple(inx)
585 else:
586 return inx
587
589
590 pjargs = []
591 for key,value in projparams.items():
592 pjargs.append('+'+key+"="+str(value)+' ')
593 return ''.join(pjargs)
594
595 -class Geod(_proj.Geod):
596 """
597 performs forward and inverse geodetic, or Great Circle,
598 computations. The forward computation (using the 'fwd' method)
599 involves determining latitude, longitude and back azimuth of a
600 computations. The forward computation (using the 'fwd' method)
601 involves determining latitude, longitude and back azimuth of a
602 terminus point given the latitude and longitude of an initial
603 point, plus azimuth and distance. The inverse computation (using
604 the 'inv' method) involves determining the forward and back
605 azimuths and distance given the latitudes and longitudes of an
606 initial and terminus point.
607 """
608 - def __new__(self, initstring=None, **kwargs):
609 """
610 initialize a Geod class instance.
611
612 Geodetic parameters for specifying the ellipsoid
613 can be given in a dictionary 'initparams', as keyword arguments,
614 or as as proj4 geod initialization string.
615 Following is a list of the ellipsoids that may be defined using the
616 'ellps' keyword (these are stored in the model variable pj_ellps)::
617
618 MERIT a=6378137.0 rf=298.257 MERIT 1983
619 SGS85 a=6378136.0 rf=298.257 Soviet Geodetic System 85
620 GRS80 a=6378137.0 rf=298.257222101 GRS 1980(IUGG, 1980)
621 IAU76 a=6378140.0 rf=298.257 IAU 1976
622 airy a=6377563.396 b=6356256.910 Airy 1830
623 APL4.9 a=6378137.0. rf=298.25 Appl. Physics. 1965
624 airy a=6377563.396 b=6356256.910 Airy 1830
625 APL4.9 a=6378137.0. rf=298.25 Appl. Physics. 1965
626 NWL9D a=6378145.0. rf=298.25 Naval Weapons Lab., 1965
627 mod_airy a=6377340.189 b=6356034.446 Modified Airy
628 andrae a=6377104.43 rf=300.0 Andrae 1876 (Den., Iclnd.)
629 aust_SA a=6378160.0 rf=298.25 Australian Natl & S. Amer. 1969
630 GRS67 a=6378160.0 rf=298.247167427 GRS 67(IUGG 1967)
631 bessel a=6377397.155 rf=299.1528128 Bessel 1841
632 bess_nam a=6377483.865 rf=299.1528128 Bessel 1841 (Namibia)
633 clrk66 a=6378206.4 b=6356583.8 Clarke 1866
634 clrk80 a=6378249.145 rf=293.4663 Clarke 1880 mod.
635 CPM a=6375738.7 rf=334.29 Comm. des Poids et Mesures 1799
636 delmbr a=6376428. rf=311.5 Delambre 1810 (Belgium)
637 engelis a=6378136.05 rf=298.2566 Engelis 1985
638 evrst30 a=6377276.345 rf=300.8017 Everest 1830
639 evrst48 a=6377304.063 rf=300.8017 Everest 1948
640 evrst56 a=6377301.243 rf=300.8017 Everest 1956
641 evrst69 a=6377295.664 rf=300.8017 Everest 1969
642 evrstSS a=6377298.556 rf=300.8017 Everest (Sabah & Sarawak)
643 fschr60 a=6378166. rf=298.3 Fischer (Mercury Datum) 1960
644 fschr60m a=6378155. rf=298.3 Modified Fischer 1960
645 fschr68 a=6378150. rf=298.3 Fischer 1968
646 helmert a=6378200. rf=298.3 Helmert 1906
647 hough a=6378270.0 rf=297. Hough
648 helmert a=6378200. rf=298.3 Helmert 1906
649 hough a=6378270.0 rf=297. Hough
650 intl a=6378388.0 rf=297. International 1909 (Hayford)
651 krass a=6378245.0 rf=298.3 Krassovsky, 1942
652 kaula a=6378163. rf=298.24 Kaula 1961
653 lerch a=6378139. rf=298.257 Lerch 1979
654 mprts a=6397300. rf=191. Maupertius 1738
655 new_intl a=6378157.5 b=6356772.2 New International 1967
656 plessis a=6376523. b=6355863. Plessis 1817 (France)
657 SEasia a=6378155.0 b=6356773.3205 Southeast Asia
658 walbeck a=6376896.0 b=6355834.8467 Walbeck
659 WGS60 a=6378165.0 rf=298.3 WGS 60
660 WGS66 a=6378145.0 rf=298.25 WGS 66
661 WGS72 a=6378135.0 rf=298.26 WGS 72
662 WGS84 a=6378137.0 rf=298.257223563 WGS 84
663 sphere a=6370997.0 b=6370997.0 Normal Sphere (r=6370997)
664
665 The parameters of the ellipsoid may also be set directly using
666 the 'a' (semi-major or equatorial axis radius) keyword, and
667 any one of the following keywords: 'b' (semi-minor,
668 or polar axis radius), 'e' (eccentricity), 'es' (eccentricity
669 squared), 'f' (flattening), or 'rf' (reciprocal flattening).
670
671 See the proj documentation (https://github.com/OSGeo/proj.4/wiki) for more
672 information about specifying ellipsoid parameters (specifically,
673 the chapter 'Specifying the Earth's figure' in the main Proj
674 users manual).
675
676 Example usage:
677
678 >>> from pyproj import Geod
679 >>> g = Geod(ellps='clrk66') # Use Clarke 1866 ellipsoid.
680 >>> # specify the lat/lons of some cities.
681 >>> boston_lat = 42.+(15./60.); boston_lon = -71.-(7./60.)
682 >>> portland_lat = 45.+(31./60.); portland_lon = -123.-(41./60.)
683 >>> newyork_lat = 40.+(47./60.); newyork_lon = -73.-(58./60.)
684 >>> london_lat = 51.+(32./60.); london_lon = -(5./60.)
685 >>> # compute forward and back azimuths, plus distance
686 >>> # between Boston and Portland.
687 >>> az12,az21,dist = g.inv(boston_lon,boston_lat,portland_lon,portland_lat)
688 >>> "%7.3f %6.3f %12.3f" % (az12,az21,dist)
689 '-66.531 75.654 4164192.708'
690 >>> # compute latitude, longitude and back azimuth of Portland,
691 >>> # given Boston lat/lon, forward azimuth and distance to Portland.
692 >>> endlon, endlat, backaz = g.fwd(boston_lon, boston_lat, az12, dist)
693 >>> "%6.3f %6.3f %13.3f" % (endlat,endlon,backaz)
694 '45.517 -123.683 75.654'
695 >>> # compute the azimuths, distances from New York to several
696 >>> # cities (pass a list)
697 >>> lons1 = 3*[newyork_lon]; lats1 = 3*[newyork_lat]
698 >>> lons2 = [boston_lon, portland_lon, london_lon]
699 >>> lats2 = [boston_lat, portland_lat, london_lat]
700 >>> az12,az21,dist = g.inv(lons1,lats1,lons2,lats2)
701 >>> for faz,baz,d in list(zip(az12,az21,dist)): "%7.3f %7.3f %9.3f" % (faz,baz,d)
702 ' 54.663 -123.448 288303.720'
703 '-65.463 79.342 4013037.318'
704 ' 51.254 -71.576 5579916.651'
705 >>> g2 = Geod('+ellps=clrk66') # use proj4 style initialization string
706 >>> az12,az21,dist = g2.inv(boston_lon,boston_lat,portland_lon,portland_lat)
707 >>> "%7.3f %6.3f %12.3f" % (az12,az21,dist)
708 '-66.531 75.654 4164192.708'
709 """
710
711
712 ellpsd = {}
713 if initstring is not None:
714 for kvpair in initstring.split():
715
716
717 if kvpair.find('=') == -1:
718 continue
719 k,v = kvpair.split('=')
720 k = k.lstrip('+')
721 if k in ['a','b','rf','f','es','e']:
722 v = float(v)
723 ellpsd[k] = v
724
725 kwargs = dict(list(kwargs.items()) + list(ellpsd.items()))
726 self.sphere = False
727 if 'ellps' in kwargs:
728
729 ellps_dict = pj_ellps[kwargs['ellps']]
730 a = ellps_dict['a']
731 if ellps_dict['description']=='Normal Sphere':
732 self.sphere = True
733 if 'b' in ellps_dict:
734 b = ellps_dict['b']
735 es = 1. - (b * b) / (a * a)
736 f = (a - b)/a
737 elif 'rf' in ellps_dict:
738 f = 1./ellps_dict['rf']
739 b = a*(1. - f)
740 es = 1. - (b * b) / (a * a)
741 else:
742
743
744
745
746
747
748 a = kwargs['a']
749 if 'b' in kwargs:
750 b = kwargs['b']
751 es = 1. - (b * b) / (a * a)
752 f = (a - b)/a
753 elif 'rf' in kwargs:
754 f = 1./kwargs['rf']
755 b = a*(1. - f)
756 es = 1. - (b * b) / (a * a)
757 elif 'f' in kwargs:
758 f = kwargs['f']
759 b = a*(1. - f)
760 es = 1. - (b/a)**2
761 elif 'es' in kwargs:
762 es = kwargs['es']
763 b = math.sqrt(a**2 - es*a**2)
764 f = (a - b)/a
765 elif 'e' in kwargs:
766 es = kwargs['e']**2
767 b = math.sqrt(a**2 - es*a**2)
768 f = (a - b)/a
769 else:
770 b = a
771 f = 0.
772 es = 0.
773
774
775 if math.fabs(f) < 1.e-8: self.sphere = True
776 self.a = a
777 self.b = b
778 self.f = f
779 self.es = es
780 return _proj.Geod.__new__(self, a, f)
781
782 - def fwd(self, lons, lats, az, dist, radians=False):
783 """
784 forward transformation - Returns longitudes, latitudes and back
785 azimuths of terminus points given longitudes (lons) and
786 latitudes (lats) of initial points, plus forward azimuths (az)
787 and distances (dist).
788 latitudes (lats) of initial points, plus forward azimuths (az)
789 and distances (dist).
790
791 Works with numpy and regular python array objects, python
792 sequences and scalars.
793
794 if radians=True, lons/lats and azimuths are radians instead of
795 degrees. Distances are in meters.
796 """
797
798 inx, xisfloat, xislist, xistuple = _copytobuffer(lons)
799 iny, yisfloat, yislist, yistuple = _copytobuffer(lats)
800 inz, zisfloat, zislist, zistuple = _copytobuffer(az)
801 ind, disfloat, dislist, distuple = _copytobuffer(dist)
802 _proj.Geod._fwd(self, inx, iny, inz, ind, radians=radians)
803
804 outx = _convertback(xisfloat,xislist,xistuple,inx)
805 outy = _convertback(yisfloat,yislist,xistuple,iny)
806 outz = _convertback(zisfloat,zislist,zistuple,inz)
807 return outx, outy, outz
808
809 - def inv(self,lons1,lats1,lons2,lats2,radians=False):
810 """
811 inverse transformation - Returns forward and back azimuths, plus
812 distances between initial points (specified by lons1, lats1) and
813 terminus points (specified by lons2, lats2).
814
815 Works with numpy and regular python array objects, python
816 sequences and scalars.
817
818 if radians=True, lons/lats and azimuths are radians instead of
819 degrees. Distances are in meters.
820 """
821
822 inx, xisfloat, xislist, xistuple = _copytobuffer(lons1)
823 iny, yisfloat, yislist, yistuple = _copytobuffer(lats1)
824 inz, zisfloat, zislist, zistuple = _copytobuffer(lons2)
825 ind, disfloat, dislist, distuple = _copytobuffer(lats2)
826 _proj.Geod._inv(self,inx,iny,inz,ind,radians=radians)
827
828 outx = _convertback(xisfloat,xislist,xistuple,inx)
829 outy = _convertback(yisfloat,yislist,xistuple,iny)
830 outz = _convertback(zisfloat,zislist,zistuple,inz)
831 return outx, outy, outz
832
833 - def npts(self, lon1, lat1, lon2, lat2, npts, radians=False):
834 """
835 Given a single initial point and terminus point (specified by
836 python floats lon1,lat1 and lon2,lat2), returns a list of
837 longitude/latitude pairs describing npts equally spaced
838 intermediate points along the geodesic between the initial and
839 terminus points.
840
841 if radians=True, lons/lats are radians instead of degrees.
842
843 Example usage:
844
845 >>> from pyproj import Geod
846 >>> g = Geod(ellps='clrk66') # Use Clarke 1866 ellipsoid.
847 >>> # specify the lat/lons of Boston and Portland.
848 >>> g = Geod(ellps='clrk66') # Use Clarke 1866 ellipsoid.
849 >>> # specify the lat/lons of Boston and Portland.
850 >>> boston_lat = 42.+(15./60.); boston_lon = -71.-(7./60.)
851 >>> portland_lat = 45.+(31./60.); portland_lon = -123.-(41./60.)
852 >>> # find ten equally spaced points between Boston and Portland.
853 >>> lonlats = g.npts(boston_lon,boston_lat,portland_lon,portland_lat,10)
854 >>> for lon,lat in lonlats: '%6.3f %7.3f' % (lat, lon)
855 '43.528 -75.414'
856 '44.637 -79.883'
857 '45.565 -84.512'
858 '46.299 -89.279'
859 '46.830 -94.156'
860 '47.149 -99.112'
861 '47.251 -104.106'
862 '47.136 -109.100'
863 '46.805 -114.051'
864 '46.262 -118.924'
865 >>> # test with radians=True (inputs/outputs in radians, not degrees)
866 >>> import math
867 >>> dg2rad = math.radians(1.)
868 >>> rad2dg = math.degrees(1.)
869 >>> lonlats = g.npts(dg2rad*boston_lon,dg2rad*boston_lat,dg2rad*portland_lon,dg2rad*portland_lat,10,radians=True)
870 >>> for lon,lat in lonlats: '%6.3f %7.3f' % (rad2dg*lat, rad2dg*lon)
871 '43.528 -75.414'
872 '44.637 -79.883'
873 '45.565 -84.512'
874 '46.299 -89.279'
875 '46.830 -94.156'
876 '47.149 -99.112'
877 '47.251 -104.106'
878 '47.136 -109.100'
879 '46.805 -114.051'
880 '46.262 -118.924'
881 """
882 lons, lats = _proj.Geod._npts(self, lon1, lat1, lon2, lat2, npts, radians=radians)
883 return list(zip(lons, lats))
884
886 """run the examples in the docstrings using the doctest module"""
887 import doctest, pyproj
888 doctest.testmod(pyproj,verbose=True)
889
890 if __name__ == "__main__": test()
891