Transformer#

The pyproj.Transformer has the capabilities of performing 2D, 3D, and 4D (time) transformations. It can do anything that the PROJ command line programs proj, cs2cs, and cct can do. This means that it allows translation between any pair of definable coordinate systems, including support for datum transformation.

Warning

The axis order may be swapped if the source and destination CRS’s are defined as having the first coordinate component point in a northerly direction (See PROJ FAQ on axis order). You can check the axis order with the pyproj.crs.CRS class. If you prefer to keep your axis order as always x,y, you can use the always_xy option when creating the pyproj.transformer.Transformer.

pyproj.Transformer#

class pyproj.transformer.Transformer(transformer_maker: TransformerMaker | None = None)[source]#

The Transformer class is for facilitating re-using transforms without needing to re-create them. The goal is to make repeated transforms faster.

Additionally, it provides multiple methods for initialization.

New in version 2.1.0.

property accuracy: float#

Expected accuracy of the transformation. -1 if unknown.

Type:

float

property area_of_use: AreaOfUse#

New in version 2.3.0.

Returns:

The area of use object with associated attributes.

Return type:

AreaOfUse

property definition: str#

Definition of the projection.

Type:

str

property description: str#

Description of the projection.

Type:

str

static from_crs(crs_from: Any, crs_to: Any, always_xy: bool = False, area_of_interest: AreaOfInterest | None = None, authority: str | None = None, accuracy: float | None = None, allow_ballpark: bool | None = None, force_over: bool = False, only_best: bool | None = None) Transformer[source]#

Make a Transformer from a pyproj.crs.CRS or input used to create one.

See:

New in version 2.2.0: always_xy

New in version 2.3.0: area_of_interest

New in version 3.1.0: authority, accuracy, allow_ballpark

New in version 3.4.0: force_over

New in version 3.5.0: only_best

Parameters:
  • crs_from (pyproj.crs.CRS or input used to create one) – Projection of input data.

  • crs_to (pyproj.crs.CRS or input used to create one) – Projection of output data.

  • always_xy (bool, default=False) – If true, the transform method will accept as input and return as output coordinates using the traditional GIS order, that is longitude, latitude for geographic CRS and easting, northing for most projected CRS.

  • area_of_interest (AreaOfInterest, optional) – The area of interest to help select the transformation.

  • authority (str, optional) – When not specified, coordinate operations from any authority will be searched, with the restrictions set in the authority_to_authority_preference database table related to the authority of the source/target CRS themselves. If authority is set to “any”, then coordinate operations from any authority will be searched. If authority is a non-empty string different from “any”, then coordinate operations will be searched only in that authority namespace (e.g. EPSG).

  • accuracy (float, optional) – The minimum desired accuracy (in metres) of the candidate coordinate operations.

  • allow_ballpark (bool, optional) – Set to False to disallow the use of Ballpark transformation in the candidate coordinate operations. Default is to allow.

  • force_over (bool, default=False) – If True, it will to force the +over flag on the transformation. Requires PROJ 9+.

  • only_best (bool, optional) – Can be set to True to cause PROJ to error out if the best transformation known to PROJ and usable by PROJ if all grids known and usable by PROJ were accessible, cannot be used. Best transformation should be understood as the transformation returned by proj_get_suggested_operation() if all known grids were accessible (either locally or through network). Note that the default value for this option can be also set with the PROJ_ONLY_BEST_DEFAULT environment variable, or with the only_best_default setting of proj.ini. The only_best kwarg overrides the default value if set. Requires PROJ 9.2+.

Return type:

Transformer

static from_pipeline(proj_pipeline: str) Transformer[source]#

Make a Transformer from a PROJ pipeline string.

The pipeline operator

See:

New in version 3.1.0: AUTH:CODE string support (e.g. EPSG:1671)

Allowed input:
  • a PROJ string

  • a WKT string

  • a PROJJSON string

  • an object code (e.g. “EPSG:1671” “urn:ogc:def:coordinateOperation:EPSG::1671”)

  • an object name. e.g “ITRF2014 to ETRF2014 (1)”. In that case as uniqueness is not guaranteed, heuristics are applied to determine the appropriate best match.

  • a OGC URN combining references for concatenated operations (e.g. “urn:ogc:def:coordinateOperation,coordinateOperation:EPSG::3895, coordinateOperation:EPSG::1618”)

Parameters:

proj_pipeline (str) – Projection pipeline string.

Return type:

Transformer

static from_proj(proj_from: Any, proj_to: Any, always_xy: bool = False, area_of_interest: AreaOfInterest | None = None) Transformer[source]#

Make a Transformer from a pyproj.Proj or input used to create one.

Deprecated since version 3.4.1: from_crs() is preferred.

New in version 2.2.0: always_xy

New in version 2.3.0: area_of_interest

Parameters:
  • proj_from (pyproj.Proj or input used to create one) – Projection of input data.

  • proj_to (pyproj.Proj or input used to create one) – Projection of output data.

  • always_xy (bool, default=False) – If true, the transform method will accept as input and return as output coordinates using the traditional GIS order, that is longitude, latitude for geographic CRS and easting, northing for most projected CRS.

  • area_of_interest (AreaOfInterest, optional) – The area of interest to help select the transformation.

Return type:

Transformer

get_last_used_operation() Transformer[source]#

New in version 3.4.0.

Note

Requires PROJ 9.1+

See: proj_trans_get_last_used_operation()

Returns:

The operation used in the transform call.

Return type:

Transformer

property has_inverse: bool#

True if an inverse mapping exists.

Type:

bool

is_exact_same(other: Any) bool[source]#

Check if the Transformer objects are the exact same. If it is not a Transformer, then it returns False.

Parameters:

other (Any) –

Return type:

bool

property is_network_enabled: bool#

New in version 3.0.0.

Returns:

If the network is enabled.

Return type:

bool

itransform(points: Any, switch: bool = False, time_3rd: bool = False, radians: bool = False, errcheck: bool = False, direction: TransformDirection | str = TransformDirection.FORWARD) Iterator[Iterable][source]#

Iterator/generator version of the function pyproj.Transformer.transform.

See: proj_trans_generic()

New in version 2.1.1: errcheck

New in version 2.2.0: direction

Parameters:
  • points (list) – List of point tuples.

  • switch (bool, default=False) – If True x, y or lon,lat coordinates of points are switched to y, x or lat, lon. Default is False.

  • time_3rd (bool, default=False) – If the input coordinates are 3 dimensional and the 3rd dimension is time.

  • radians (bool, default=False) – If True, will expect input data to be in radians and will return radians if the projection is geographic. Otherwise, it uses degrees. Ignored for pipeline transformations with pyproj 2, but will work in pyproj 3.

  • errcheck (bool, default=False) – If True, an exception is raised if the errors are found in the process. If False, inf is returned for errors.

  • direction (pyproj.enums.TransformDirection, optional) – The direction of the transform. Default is pyproj.enums.TransformDirection.FORWARD.

Example

>>> from pyproj import Transformer
>>> transformer = Transformer.from_crs(4326, 2100)
>>> points = [(22.95, 40.63), (22.81, 40.53), (23.51, 40.86)]
>>> for pt in transformer.itransform(points): '{:.3f} {:.3f}'.format(*pt)
'2221638.801 2637034.372'
'2212924.125 2619851.898'
'2238294.779 2703763.736'
>>> pipeline_str = (
...     "+proj=pipeline +step +proj=longlat +ellps=WGS84 "
...     "+step +proj=unitconvert +xy_in=rad +xy_out=deg"
... )
>>> pipe_trans = Transformer.from_pipeline(pipeline_str)
>>> for pt in pipe_trans.itransform([(2.1, 0.001)]):
...     '{:.3f} {:.3f}'.format(*pt)
'2.100 0.001'
>>> transproj = Transformer.from_crs(
...     {"proj":'geocent', "ellps":'WGS84', "datum":'WGS84'},
...     "EPSG:4326",
...     always_xy=True,
... )
>>> for pt in transproj.itransform(
...     [(-2704026.010, -4253051.810, 3895878.820)],
...     radians=True,
... ):
...     '{:.3f} {:.3f} {:.3f}'.format(*pt)
'-2.137 0.661 -20.531'
>>> transprojr = Transformer.from_crs(
...     "EPSG:4326",
...     {"proj":'geocent', "ellps":'WGS84', "datum":'WGS84'},
...     always_xy=True,
... )
>>> for pt in transprojr.itransform(
...     [(-2.137, 0.661, -20.531)],
...     radians=True
... ):
...     '{:.3f} {:.3f} {:.3f}'.format(*pt)
'-2704214.394 -4254414.478 3894270.731'
>>> transproj_eq = Transformer.from_crs(
...     'EPSG:4326',
...     '+proj=longlat +datum=WGS84 +no_defs +type=crs',
...     always_xy=True,
... )
>>> for pt in transproj_eq.itransform([(-2.137, 0.661)]):
...     '{:.3f} {:.3f}'.format(*pt)
'-2.137 0.661'
property name: str#

Name of the projection.

Type:

str

property operations: tuple[CoordinateOperation] | None#

New in version 2.4.0.

Returns:

The operations in a concatenated operation.

Return type:

tuple[CoordinateOperation]

property remarks: str#

New in version 2.4.0.

Returns:

Remarks about object.

Return type:

str

property scope: str#

New in version 2.4.0.

Returns:

Scope of object.

Return type:

str

property source_crs: CRS | None#

New in version 3.3.0.

Returns:

The source CRS of a CoordinateOperation.

Return type:

CRS | None

property target_crs: CRS | None#

New in version 3.3.0.

Returns:

The target CRS of a CoordinateOperation.

Return type:

CRS | None

to_json(pretty: bool = False, indentation: int = 2) str[source]#

Convert the projection to a JSON string.

New in version 2.4.0.

Parameters:
  • pretty (bool, default=False) – If True, it will set the output to be a multiline string.

  • indentation (int, default=2) – If pretty is True, it will set the width of the indentation.

Returns:

The JSON string.

Return type:

str

to_json_dict() dict[source]#

Convert the projection to a JSON dictionary.

New in version 2.4.0.

Returns:

The JSON dictionary.

Return type:

dict

to_proj4(version: ProjVersion | str = ProjVersion.PROJ_5, pretty: bool = False) str[source]#

Convert the projection to a PROJ string.

New in version 3.1.0.

Parameters:
Returns:

The PROJ string.

Return type:

str

to_wkt(version: WktVersion | str = WktVersion.WKT2_2019, pretty: bool = False) str[source]#

Convert the projection to a WKT string.

Version options:
  • WKT2_2015

  • WKT2_2015_SIMPLIFIED

  • WKT2_2019

  • WKT2_2019_SIMPLIFIED

  • WKT1_GDAL

  • WKT1_ESRI

Parameters:
Returns:

The WKT string.

Return type:

str

transform(xx: Any, yy: Any, radians: bool = False, errcheck: bool = False, direction: TransformDirection | str = TransformDirection.FORWARD, inplace: bool = False) tuple[Any, Any][source]#
transform(xx: Any, yy: Any, zz: Any, radians: bool = False, errcheck: bool = False, direction: TransformDirection | str = TransformDirection.FORWARD, inplace: bool = False) tuple[Any, Any, Any]
transform(xx: Any, yy: Any, zz: Any, tt: Any, radians: bool = False, errcheck: bool = False, direction: TransformDirection | str = TransformDirection.FORWARD, inplace: bool = False) tuple[Any, Any, Any, Any]

Transform points between two coordinate systems.

See: proj_trans_generic()

New in version 2.1.1: errcheck

New in version 2.2.0: direction

New in version 3.2.0: inplace

Accepted numeric scalar or array:

Parameters:
  • xx (scalar or array) – Input x coordinate(s).

  • yy (scalar or array) – Input y coordinate(s).

  • zz (scalar or array, optional) – Input z coordinate(s).

  • tt (scalar or array, optional) – Input time coordinate(s).

  • radians (bool, default=False) – If True, will expect input data to be in radians and will return radians if the projection is geographic. Otherwise, it uses degrees. Ignored for pipeline transformations with pyproj 2, but will work in pyproj 3.

  • errcheck (bool, default=False) – If True, an exception is raised if the errors are found in the process. If False, inf is returned for errors.

  • direction (pyproj.enums.TransformDirection, optional) – The direction of the transform. Default is pyproj.enums.TransformDirection.FORWARD.

  • inplace (bool, default=False) – If True, will attempt to write the results to the input array instead of returning a new array. This will fail if the input is not an array in C order with the double data type.

Example

>>> from pyproj import Transformer
>>> transformer = Transformer.from_crs("EPSG:4326", "EPSG:3857")
>>> x3, y3 = transformer.transform(33, 98)
>>> f"{x3:.3f}  {y3:.3f}"
'10909310.098  3895303.963'
>>> pipeline_str = (
...     "+proj=pipeline +step +proj=longlat +ellps=WGS84 "
...     "+step +proj=unitconvert +xy_in=rad +xy_out=deg"
... )
>>> pipe_trans = Transformer.from_pipeline(pipeline_str)
>>> xt, yt = pipe_trans.transform(2.1, 0.001)
>>> f"{xt:.3f}  {yt:.3f}"
'2.100  0.001'
>>> transproj = Transformer.from_crs(
...     {"proj":'geocent', "ellps":'WGS84', "datum":'WGS84'},
...     "EPSG:4326",
...     always_xy=True,
... )
>>> xpj, ypj, zpj = transproj.transform(
...     -2704026.010,
...     -4253051.810,
...     3895878.820,
...     radians=True,
... )
>>> f"{xpj:.3f} {ypj:.3f} {zpj:.3f}"
'-2.137 0.661 -20.531'
>>> transprojr = Transformer.from_crs(
...     "EPSG:4326",
...     {"proj":'geocent', "ellps":'WGS84', "datum":'WGS84'},
...     always_xy=True,
... )
>>> xpjr, ypjr, zpjr = transprojr.transform(xpj, ypj, zpj, radians=True)
>>> f"{xpjr:.3f} {ypjr:.3f} {zpjr:.3f}"
'-2704026.010 -4253051.810 3895878.820'
>>> transformer = Transformer.from_crs("EPSG:4326", 4326)
>>> xeq, yeq = transformer.transform(33, 98)
>>> f"{xeq:.0f}  {yeq:.0f}"
'33  98'
transform_bounds(left: float, bottom: float, right: float, top: float, densify_pts: int = 21, radians: bool = False, errcheck: bool = False, direction: TransformDirection | str = TransformDirection.FORWARD) tuple[float, float, float, float][source]#

New in version 3.1.0.

See: proj_trans_bounds()

Transform boundary densifying the edges to account for nonlinear transformations along these edges and extracting the outermost bounds.

If the destination CRS is geographic and right < left then the bounds crossed the antimeridian. In this scenario there are two polygons, one on each side of the antimeridian. The first polygon should be constructed with (left, bottom, 180, top) and the second with (-180, bottom, top, right).

To construct the bounding polygons with shapely:

def bounding_polygon(left, bottom, right, top):
    if right < left:
        return shapely.geometry.MultiPolygon(
            [
                shapely.geometry.box(left, bottom, 180, top),
                shapely.geometry.box(-180, bottom, right, top),
            ]
        )
    return shapely.geometry.box(left, bottom, right, top)
Parameters:
  • left (float) – Minimum bounding coordinate of the first axis in source CRS (or the target CRS if using the reverse direction).

  • bottom (float) – Minimum bounding coordinate of the second axis in source CRS. (or the target CRS if using the reverse direction).

  • right (float) – Maximum bounding coordinate of the first axis in source CRS. (or the target CRS if using the reverse direction).

  • top (float) – Maximum bounding coordinate of the second axis in source CRS. (or the target CRS if using the reverse direction).

  • densify_points (uint, default=21) – Number of points to add to each edge to account for nonlinear edges produced by the transform process. Large numbers will produce worse performance.

  • radians (bool, default=False) – If True, will expect input data to be in radians and will return radians if the projection is geographic. Otherwise, it uses degrees.

  • errcheck (bool, default=False) – If True, an exception is raised if the errors are found in the process. If False, inf is returned for errors.

  • direction (pyproj.enums.TransformDirection, optional) – The direction of the transform. Default is pyproj.enums.TransformDirection.FORWARD.

Returns:

left, bottom, right, top – Outermost coordinates in target coordinate reference system.

Return type:

float

pyproj.transformer.TransformerGroup#

class pyproj.transformer.TransformerGroup(crs_from: Any, crs_to: Any, always_xy: bool = False, area_of_interest: AreaOfInterest | None = None, authority: str | None = None, accuracy: float | None = None, allow_ballpark: bool = True, allow_superseded: bool = False)[source]#

The TransformerGroup is a set of possible transformers from one CRS to another.

New in version 2.3.0.

Warning

CoordinateOperation and Transformer objects returned are not thread-safe.

From PROJ docs:

The operations are sorted with the most relevant ones first: by
descending area (intersection of the transformation area with the
area of interest, or intersection of the transformation with the
area of use of the CRS), and by increasing accuracy. Operations
with unknown accuracy are sorted last, whatever their area.
__init__(crs_from: Any, crs_to: Any, always_xy: bool = False, area_of_interest: AreaOfInterest | None = None, authority: str | None = None, accuracy: float | None = None, allow_ballpark: bool = True, allow_superseded: bool = False) None[source]#

Get all possible transformations from a pyproj.crs.CRS or input used to create one.

New in version 3.4.0: authority, accuracy, allow_ballpark

New in version 3.6.0: allow_superseded

Parameters:
  • crs_from (pyproj.crs.CRS or input used to create one) – Projection of input data.

  • crs_to (pyproj.crs.CRS or input used to create one) – Projection of output data.

  • always_xy (bool, default=False) – If true, the transform method will accept as input and return as output coordinates using the traditional GIS order, that is longitude, latitude for geographic CRS and easting, northing for most projected CRS.

  • area_of_interest (AreaOfInterest, optional) – The area of interest to help order the transformations based on the best operation for the area.

  • authority (str, optional) – When not specified, coordinate operations from any authority will be searched, with the restrictions set in the authority_to_authority_preference database table related to the authority of the source/target CRS themselves. If authority is set to “any”, then coordinate operations from any authority will be searched. If authority is a non-empty string different from “any”, then coordinate operations will be searched only in that authority namespace (e.g. EPSG).

  • accuracy (float, optional) – The minimum desired accuracy (in metres) of the candidate coordinate operations.

  • allow_ballpark (bool, default=True) – Set to False to disallow the use of Ballpark transformation in the candidate coordinate operations. Default is to allow.

  • allow_superseded (bool, default=False) – Set to True to allow the use of superseded (but not deprecated) transformations in the candidate coordinate operations. Default is to disallow.

property best_available: bool#

If True, the best possible transformer is available.

Type:

bool

download_grids(directory: str | Path | None = None, open_license: bool = True, verbose: bool = False) None[source]#

New in version 3.0.0.

Download missing grids that can be downloaded automatically.

Warning

There are cases where the URL to download the grid is missing. In those cases, you can enable enable Debugging Internal PROJ and perform a transformation. The logs will show the grids PROJ searches for.

Parameters:
  • directory (str or Path, optional) – The directory to download the grids to. Defaults to pyproj.datadir.get_user_data_dir()

  • open_license (bool, default=True) – If True, will only download grids with an open license.

  • verbose (bool, default=False) – If True, will print information about grids downloaded.

property transformers: list[Transformer]#

list[Transformer]: List of available Transformer associated with the transformation.

property unavailable_operations: list[CoordinateOperation]#

list[pyproj.crs.CoordinateOperation]: List of pyproj.crs.CoordinateOperation that are not available due to missing grids.

pyproj.transform#

pyproj.transformer.transform(p1: Any, p2: Any, x: Any, y: Any, z: Any | None = None, tt: Any | None = None, radians: bool = False, errcheck: bool = False, always_xy: bool = False)[source]#

New in version 2.2.0: always_xy

Deprecated since version 2.6.1: This function is deprecated. See: Upgrading to pyproj 2 from pyproj 1

x2, y2, z2 = transform(p1, p2, x1, y1, z1)

Transform points between two coordinate systems defined by the Proj instances p1 and p2.

The points x1,y1,z1 in the coordinate system defined by p1 are transformed to x2,y2,z2 in the coordinate system defined by p2.

z1 is optional, if it is not set it is assumed to be zero (and only x2 and y2 are returned). If the optional keyword ‘radians’ is True (default is False), then all input and output coordinates will be in radians instead of the default of degrees for geographic input/output projections. If the optional keyword ‘errcheck’ is set to True an exception is raised if the transformation is invalid. By default errcheck=False and inf is returned for an invalid transformation (and no exception is raised). If always_xy is toggled, the transform method will accept as input and return as output coordinates using the traditional GIS order, that is longitude, latitude for geographic CRS and easting, northing for most projected CRS.

In addition to converting between cartographic and geographic projection coordinates, this function can take care of datum shifts (which cannot be done using the __call__ method of the Proj instances). It also allows for one of the coordinate systems to be geographic (proj = ‘latlong’).

x,y and z can be numpy or regular python arrays, python lists/tuples or scalars. Arrays are fastest. For projections in geocentric coordinates, values of x and y are given in meters. z is always meters.

pyproj.itransform#

pyproj.transformer.itransform(p1: Any, p2: Any, points: Iterable[Iterable], switch: bool = False, time_3rd: bool = False, radians: bool = False, errcheck: bool = False, always_xy: bool = False)[source]#

New in version 2.2.0: always_xy

Deprecated since version 2.6.1: This function is deprecated. See: Upgrading to pyproj 2 from pyproj 1

points2 = itransform(p1, p2, points1) Iterator/generator version of the function pyproj.transform. Transform points between two coordinate systems defined by the Proj instances p1 and p2. This function can be used as an alternative to pyproj.transform when there is a need to transform a big number of coordinates lazily, for example when reading and processing from a file. Points1 is an iterable/generator of coordinates x1,y1(,z1) or lon1,lat1(,z1) in the coordinate system defined by p1. Points2 is an iterator that returns tuples of x2,y2(,z2) or lon2,lat2(,z2) coordinates in the coordinate system defined by p2. z are provided optionally.

Points1 can be:
  • a tuple/list of tuples/lists i.e. for 2d points: [(xi,yi),(xj,yj),….(xn,yn)]

  • a Nx3 or Nx2 2d numpy array where N is the point number

  • a generator of coordinates (xi,yi) for 2d points or (xi,yi,zi) for 3d

If optional keyword ‘switch’ is True (default is False) then x, y or lon,lat coordinates of points are switched to y, x or lat, lon. If the optional keyword ‘radians’ is True (default is False), then all input and output coordinates will be in radians instead of the default of degrees for geographic input/output projections. If the optional keyword ‘errcheck’ is set to True an exception is raised if the transformation is invalid. By default errcheck=False and inf is returned for an invalid transformation (and no exception is raised). If always_xy is toggled, the transform method will accept as input and return as output coordinates using the traditional GIS order, that is longitude, latitude for geographic CRS and easting, northing for most projected CRS.

Example usage:

>>> from pyproj import Proj, itransform
>>> # projection 1: WGS84
>>> # (defined by epsg code 4326)
>>> p1 = Proj('epsg:4326', preserve_units=False)
>>> # projection 2: GGRS87 / Greek Grid
>>> p2 = Proj('epsg:2100', preserve_units=False)
>>> # Three points with coordinates lon, lat in p1
>>> points = [(22.95, 40.63), (22.81, 40.53), (23.51, 40.86)]
>>> # transform this point to projection 2 coordinates.
>>> for pt in itransform(p1,p2,points, always_xy=True): '%6.3f %7.3f' % pt
'411050.470 4497928.574'
'399060.236 4486978.710'
'458553.243 4523045.485'
>>> for pt in itransform(4326, 4326, [(30, 60)]):
...     '{:.0f} {:.0f}'.format(*pt)
'30 60'