Geod

pyproj.Geod

class pyproj.Geod(initstring=None, **kwargs)[source]

Bases: pyproj._geod.Geod

performs forward and inverse geodetic, or Great Circle, computations. The forward computation (using the ‘fwd’ method) involves determining latitude, longitude and back azimuth of a computations. The forward computation (using the ‘fwd’ method) involves determining latitude, longitude and back azimuth of a terminus point given the latitude and longitude of an initial point, plus azimuth and distance. The inverse computation (using the ‘inv’ method) involves determining the forward and back azimuths and distance given the latitudes and longitudes of an initial and terminus point.

initstring

The string form of the user input used to create the Geod.

Type

str

sphere

If True, it is a sphere.

Type

bool

a

The ellipsoid equatorial radius, or semi-major axis.

Type

float

b

The ellipsoid polar radius, or semi-minor axis.

Type

float

es

The ‘eccentricity’ of the ellipse, squared (1-b2/a2).

Type

float

f

The ellipsoid ‘flattening’ parameter ( (a-b)/a ).

Type

float

__init__(initstring=None, **kwargs)[source]

initialize a Geod class instance.

Geodetic parameters for specifying the ellipsoid can be given in a dictionary ‘initparams’, as keyword arguments, or as as proj geod initialization string.

You can get a dictionary of ellipsoids using get_ellps_map() or with the variable pyproj.pj_ellps.

The parameters of the ellipsoid may also be set directly using the ‘a’ (semi-major or equatorial axis radius) keyword, and any one of the following keywords: ‘b’ (semi-minor, or polar axis radius), ‘e’ (eccentricity), ‘es’ (eccentricity squared), ‘f’ (flattening), or ‘rf’ (reciprocal flattening).

See the proj documentation (https://proj.org) for more information about specifying ellipsoid parameters.

Example usage:

>>> from pyproj import Geod
>>> g = Geod(ellps='clrk66') # Use Clarke 1866 ellipsoid.
>>> # specify the lat/lons of some cities.
>>> boston_lat = 42.+(15./60.); boston_lon = -71.-(7./60.)
>>> portland_lat = 45.+(31./60.); portland_lon = -123.-(41./60.)
>>> newyork_lat = 40.+(47./60.); newyork_lon = -73.-(58./60.)
>>> london_lat = 51.+(32./60.); london_lon = -(5./60.)
>>> # compute forward and back azimuths, plus distance
>>> # between Boston and Portland.
>>> az12,az21,dist = g.inv(boston_lon,boston_lat,portland_lon,portland_lat)
>>> "%7.3f %6.3f %12.3f" % (az12,az21,dist)
'-66.531 75.654  4164192.708'
>>> # compute latitude, longitude and back azimuth of Portland,
>>> # given Boston lat/lon, forward azimuth and distance to Portland.
>>> endlon, endlat, backaz = g.fwd(boston_lon, boston_lat, az12, dist)
>>> "%6.3f  %6.3f %13.3f" % (endlat,endlon,backaz)
'45.517  -123.683        75.654'
>>> # compute the azimuths, distances from New York to several
>>> # cities (pass a list)
>>> lons1 = 3*[newyork_lon]; lats1 = 3*[newyork_lat]
>>> lons2 = [boston_lon, portland_lon, london_lon]
>>> lats2 = [boston_lat, portland_lat, london_lat]
>>> az12,az21,dist = g.inv(lons1,lats1,lons2,lats2)
>>> for faz, baz, d in list(zip(az12,az21,dist)):
...     "%7.3f %7.3f %9.3f" % (faz, baz, d)
' 54.663 -123.448 288303.720'
'-65.463  79.342 4013037.318'
' 51.254 -71.576 5579916.651'
>>> g2 = Geod('+ellps=clrk66') # use proj4 style initialization string
>>> az12,az21,dist = g2.inv(boston_lon,boston_lat,portland_lon,portland_lat)
>>> "%7.3f %6.3f %12.3f" % (az12,az21,dist)
'-66.531 75.654  4164192.708'
fwd(lons, lats, az, dist, radians=False)[source]

forward transformation - Returns longitudes, latitudes and back azimuths of terminus points given longitudes (lons) and latitudes (lats) of initial points, plus forward azimuths (az) and distances (dist). latitudes (lats) of initial points, plus forward azimuths (az) and distances (dist).

Works with numpy and regular python array objects, python sequences and scalars.

if radians=True, lons/lats and azimuths are radians instead of degrees. Distances are in meters.

geometry_area_perimeter(geometry, radians=False)[source]

New in version 2.3.0.

A simple interface for computing the area (meters^2) and perimeter (meters) of a geodesic polygon as a shapely geometry.

Warning

Only simple polygons (which are not self-intersecting) are allowed.

Note

lats should be in the range [-90 deg, 90 deg].

There’s no need to “close” the polygon by repeating the first vertex. The area returned is signed with counter-clockwise traversal being treated as positive.

If it is a Polygon, it will return the area and exterior perimeter. It will subtract the area of the interior holes. If it is a MultiPolygon or MultiLine, it will return the sum of the areas and perimeters of all geometries.

Example usage:

>>> from pyproj import Geod
>>> from shapely.geometry import LineString, Point, Polygon
>>> geod = Geod(ellps="WGS84")
>>> poly_area, poly_perimeter = geod.geometry_area_perimeter(
...     Polygon(
...         LineString([
...             Point(1, 1), Point(1, 10), Point(10, 10), Point(10, 1)
...         ]),
...         holes=[LineString([Point(1, 2), Point(3, 4), Point(5, 2)])],
...     )
... )
>>> "{:.3f} {:.3f}".format(poly_area, poly_perimeter)
'-944373881400.339 3979008.036'
Parameters
  • geometry (shapely.geometry.BaseGeometry) – The geometry to calculate the area and perimeter from.

  • radians (bool, optional) –

    If True, the input data is assumed to be in radians.

    Returns

  • -------

  • float) ((float,) – The geodesic area (meters^2) and permimeter (meters) of the polygon.

geometry_length(geometry, radians=False)[source]

New in version 2.3.0.

Returns the geodesic length (meters) of the shapely geometry.

If it is a Polygon, it will return the sum of the lengths along the perimeter. If it is a MultiPolygon or MultiLine, it will return the sum of the lengths.

Example usage:

>>> from pyproj import Geod
>>> from shapely.geometry import Point, LineString
>>> line_string = LineString([Point(1, 2), Point(3, 4)])
>>> geod = Geod(ellps="WGS84")
>>> "{:.3f}".format(geod.geometry_length(line_string))
'313588.397'
Parameters
  • geometry (shapely.geometry.BaseGeometry) – The geometry to calculate the length from.

  • radians (bool, optional) – If True, the input data is assumed to be in radians.

Returns

float

Return type

The total geodesic length of the geometry (meters)

inv(lons1, lats1, lons2, lats2, radians=False)[source]

inverse transformation - Returns forward and back azimuths, plus distances between initial points (specified by lons1, lats1) and terminus points (specified by lons2, lats2).

Works with numpy and regular python array objects, python sequences and scalars.

if radians=True, lons/lats and azimuths are radians instead of degrees. Distances are in meters.

line_length(lons, lats, radians=False)[source]

New in version 2.3.0.

Calculate the total distance between points along a line.

>>> from pyproj import Geod
>>> geod = Geod('+a=6378137 +f=0.0033528106647475126')
>>> lats = [-72.9, -71.9, -74.9, -74.3, -77.5, -77.4, -71.7, -65.9, -65.7,
...         -66.6, -66.9, -69.8, -70.0, -71.0, -77.3, -77.9, -74.7]
>>> lons = [-74, -102, -102, -131, -163, 163, 172, 140, 113,
...         88, 59, 25, -4, -14, -33, -46, -61]
>>> total_length = geod.line_length(lons, lats)
>>> "{:.3f}".format(total_length)
'14259605.611'
Parameters
  • lons (array, numpy.ndarray, list, tuple, or scalar) – The longitude points along a line.

  • lats (array, numpy.ndarray, list, tuple, or scalar) – The latitude points along a line.

  • radians (bool, optional) – If True, the input data is assumed to be in radians.

Returns

float

Return type

The total length of the line.

line_lengths(lons, lats, radians=False)[source]

New in version 2.3.0.

Calculate the distances between points along a line.

>>> from pyproj import Geod
>>> geod = Geod(ellps="WGS84")
>>> lats = [-72.9, -71.9, -74.9]
>>> lons = [-74, -102, -102]
>>> for line_length in geod.line_lengths(lons, lats):
...     "{:.3f}".format(line_length)
'943065.744'
'334805.010'
Parameters
  • lons (array, numpy.ndarray, list, tuple, or scalar) – The longitude points along a line.

  • lats (array, numpy.ndarray, list, tuple, or scalar) – The latitude points along a line.

  • radians (bool, optional) – If True, the input data is assumed to be in radians.

Returns

The total length of the line.

Return type

array, numpy.ndarray, list, tuple, or scalar

npts(lon1, lat1, lon2, lat2, npts, radians=False)[source]

Given a single initial point and terminus point (specified by python floats lon1,lat1 and lon2,lat2), returns a list of longitude/latitude pairs describing npts equally spaced intermediate points along the geodesic between the initial and terminus points.

if radians=True, lons/lats are radians instead of degrees.

Example usage:

>>> from pyproj import Geod
>>> g = Geod(ellps='clrk66') # Use Clarke 1866 ellipsoid.
>>> # specify the lat/lons of Boston and Portland.
>>> boston_lat = 42.+(15./60.); boston_lon = -71.-(7./60.)
>>> portland_lat = 45.+(31./60.); portland_lon = -123.-(41./60.)
>>> # find ten equally spaced points between Boston and Portland.
>>> lonlats = g.npts(boston_lon,boston_lat,portland_lon,portland_lat,10)
>>> for lon,lat in lonlats: '%6.3f  %7.3f' % (lat, lon)
'43.528  -75.414'
'44.637  -79.883'
'45.565  -84.512'
'46.299  -89.279'
'46.830  -94.156'
'47.149  -99.112'
'47.251  -104.106'
'47.136  -109.100'
'46.805  -114.051'
'46.262  -118.924'
>>> # test with radians=True (inputs/outputs in radians, not degrees)
>>> import math
>>> dg2rad = math.radians(1.)
>>> rad2dg = math.degrees(1.)
>>> lonlats = g.npts(
...    dg2rad*boston_lon,
...    dg2rad*boston_lat,
...    dg2rad*portland_lon,
...    dg2rad*portland_lat,
...    10,
...    radians=True
... )
>>> for lon,lat in lonlats: '%6.3f  %7.3f' % (rad2dg*lat, rad2dg*lon)
'43.528  -75.414'
'44.637  -79.883'
'45.565  -84.512'
'46.299  -89.279'
'46.830  -94.156'
'47.149  -99.112'
'47.251  -104.106'
'47.136  -109.100'
'46.805  -114.051'
'46.262  -118.924'
polygon_area_perimeter(lons, lats, radians=False)[source]

New in version 2.3.0.

A simple interface for computing the area (meters^2) and perimeter (meters) of a geodesic polygon.

Warning

Only simple polygons (which are not self-intersecting) are allowed.

Note

lats should be in the range [-90 deg, 90 deg].

There’s no need to “close” the polygon by repeating the first vertex. The area returned is signed with counter-clockwise traversal being treated as positive.

Example usage:

>>> from pyproj import Geod
>>> geod = Geod('+a=6378137 +f=0.0033528106647475126')
>>> lats = [-72.9, -71.9, -74.9, -74.3, -77.5, -77.4, -71.7, -65.9, -65.7,
...         -66.6, -66.9, -69.8, -70.0, -71.0, -77.3, -77.9, -74.7]
>>> lons = [-74, -102, -102, -131, -163, 163, 172, 140, 113,
...         88, 59, 25, -4, -14, -33, -46, -61]
>>> poly_area, poly_perimeter = geod.polygon_area_perimeter(lons, lats)
>>> "{:.1f} {:.1f}".format(poly_area, poly_perimeter)
'13376856682207.4 14710425.4'
Parameters
  • lons (array, numpy.ndarray, list, tuple, or scalar) – An array of longitude values.

  • lats (array, numpy.ndarray, list, tuple, or scalar) – An array of latitude values.

  • radians (bool, optional) – If True, the input data is assumed to be in radians.

Returns

The geodesic area (meters^2) and permimeter (meters) of the polygon.

Return type

(float, float)