operation namespace

namespace operation

Coordinate operations (relationship between any two coordinate reference systems).

osgeo.proj.operation namespace

This covers Conversion, Transformation, PointMotionOperation or ConcatenatedOperation.

Typedefs

typedef std::shared_ptr<CoordinateOperation> CoordinateOperationPtr

Shared pointer of CoordinateOperation

typedef util::nn<CoordinateOperationPtr> CoordinateOperationNNPtr

Non-null shared pointer of CoordinateOperation

using GeneralOperationParameterPtr = std::shared_ptr<GeneralOperationParameter>

Shared pointer of GeneralOperationParameter

using GeneralOperationParameterNNPtr = util::nn<GeneralOperationParameterPtr>

Non-null shared pointer of GeneralOperationParameter

using OperationParameterPtr = std::shared_ptr<OperationParameter>

Shared pointer of OperationParameter

using OperationParameterNNPtr = util::nn<OperationParameterPtr>

Non-null shared pointer of OperationParameter

using GeneralParameterValuePtr = std::shared_ptr<GeneralParameterValue>

Shared pointer of GeneralParameterValue

using GeneralParameterValueNNPtr = util::nn<GeneralParameterValuePtr>

Non-null shared pointer of GeneralParameterValue

using ParameterValuePtr = std::shared_ptr<ParameterValue>

Shared pointer of ParameterValue

using ParameterValueNNPtr = util::nn<ParameterValuePtr>

Non-null shared pointer of ParameterValue

using OperationParameterValuePtr = std::shared_ptr<OperationParameterValue>

Shared pointer of OperationParameterValue

using OperationParameterValueNNPtr = util::nn<OperationParameterValuePtr>

Non-null shared pointer of OperationParameterValue

using OperationMethodPtr = std::shared_ptr<OperationMethod>

Shared pointer of OperationMethod

using OperationMethodNNPtr = util::nn<OperationMethodPtr>

Non-null shared pointer of OperationMethod

using SingleOperationPtr = std::shared_ptr<SingleOperation>

Shared pointer of SingleOperation

using SingleOperationNNPtr = util::nn<SingleOperationPtr>

Non-null shared pointer of SingleOperation

typedef std::shared_ptr<Conversion> ConversionPtr

Shared pointer of Conversion

typedef util::nn<ConversionPtr> ConversionNNPtr

Non-null shared pointer of Conversion

using TransformationPtr = std::shared_ptr<Transformation>

Shared pointer of Transformation

using TransformationNNPtr = util::nn<TransformationPtr>

Non-null shared pointer of Transformation

using PointMotionOperationPtr = std::shared_ptr<PointMotionOperation>

Shared pointer of PointMotionOperation

using PointMotionOperationNNPtr = util::nn<PointMotionOperationPtr>

Non-null shared pointer of PointMotionOperation

using ConcatenatedOperationPtr = std::shared_ptr<ConcatenatedOperation>

Shared pointer of ConcatenatedOperation

using ConcatenatedOperationNNPtr = util::nn<ConcatenatedOperationPtr>

Non-null shared pointer of ConcatenatedOperation

using CoordinateOperationContextPtr = std::unique_ptr<CoordinateOperationContext>

Unique pointer of CoordinateOperationContext

using CoordinateOperationContextNNPtr = util::nn<CoordinateOperationContextPtr>

Non-null unique pointer of CoordinateOperationContext

using CoordinateOperationFactoryPtr = std::unique_ptr<CoordinateOperationFactory>

Unique pointer of CoordinateOperationFactory

using CoordinateOperationFactoryNNPtr = util::nn<CoordinateOperationFactoryPtr>

Non-null unique pointer of CoordinateOperationFactory

Functions

static const char *getCRSQualifierStr(const crs::CRSPtr &crs)
static std::string buildOpName(const char *opType, const crs::CRSPtr &source, const crs::CRSPtr &target)
static util::PropertyMap createPropertiesForInverse(const CoordinateOperation *op, bool derivedFrom, bool approximateInversion)
static bool isTimeDependent(const std::string &methodName)
static double negate(double val)
static CoordinateOperationPtr createApproximateInverseIfPossible(const Transformation *op)
static void exportSourceCRSAndTargetCRSToWKT(const CoordinateOperation *co, io::WKTFormatter *formatter)
static crs::CRSNNPtr getResolvedCRS(const crs::CRSNNPtr &crs, const CoordinateOperationContextNNPtr &context)
class ConcatenatedOperation : public osgeo::proj::operation::CoordinateOperation
#include <coordinateoperation.hpp>

An ordered sequence of two or more single coordinate operations (SingleOperation).

The sequence of coordinate operations is constrained by the requirement that the source coordinate reference system of step n+1 shall be the same as the target coordinate reference system of step n.

Remark

Implements ConcatenatedOperation from ISO_19111_2019

Public Functions

const std::vector<CoordinateOperationNNPtr> &operations() const

Return the operation steps of the concatenated operation.

Return

the operation steps.

CoordinateOperationNNPtr inverse() const

Return the inverse of the coordinate operation.

Exceptions

std::set<GridDescription> gridsNeeded(const io::DatabaseContextPtr &databaseContext) const

Return grids needed by an operation.

Public Static Functions

ConcatenatedOperationNNPtr create(const util::PropertyMap &properties, const std::vector<CoordinateOperationNNPtr> &operationsIn, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a ConcatenatedOperation.

Return

new Transformation.

Parameters
  • properties: See general_properties. At minimum the name should be defined.

  • operationsIn: Vector of the CoordinateOperation steps.

  • accuracies: Vector of positional accuracy (might be empty).

Exceptions

CoordinateOperationNNPtr createComputeMetadata(const std::vector<CoordinateOperationNNPtr> &operationsIn, bool checkExtent)

Instantiate a ConcatenatedOperation, or return a single coordinate operation.

This computes its accuracy from the sum of its member operations, its extent

Parameters
  • operationsIn: Vector of the CoordinateOperation steps.

  • checkExtent: Whether we should check the non-emptyness of the intersection of the extents of the operations

Exceptions

class Conversion : public osgeo::proj::operation::SingleOperation
#include <coordinateoperation.hpp>

A mathematical operation on coordinates in which the parameter values are defined rather than empirically derived.

Application of the coordinate conversion introduces no error into output coordinates. The best-known example of a coordinate conversion is a map projection. For coordinate conversions the output coordinates are referenced to the same datum as are the input coordinates.

Coordinate conversions forming a component of a derived CRS have a source crs::CRS and a target crs::CRS that are NOT specified through the source and target associations, but through associations from crs::DerivedCRS to crs::SingleCRS.

Remark

Implements Conversion from ISO_19111_2019

Public Functions

CoordinateOperationNNPtr inverse() const

Return the inverse of the coordinate operation.

Exceptions

bool isUTM(int &zone, bool &north) const

Return whether a conversion is a Universal Transverse Mercator conversion.

Return

true if it is a UTM conversion.

Parameters
  • [out] zone: UTM zone number between 1 and 60.

  • [out] north: true for UTM northern hemisphere, false for UTM southern hemisphere.

ConversionNNPtr identify() const

Return a Conversion object where some parameters are better identified.

Return

a new Conversion.

ConversionPtr convertToOtherMethod(int targetEPSGCode) const

Return an equivalent projection.

Currently implemented:

  • EPSG_CODE_METHOD_MERCATOR_VARIANT_A (1SP) to EPSG_CODE_METHOD_MERCATOR_VARIANT_B (2SP)

  • EPSG_CODE_METHOD_MERCATOR_VARIANT_B (2SP) to EPSG_CODE_METHOD_MERCATOR_VARIANT_A (1SP)

  • EPSG_CODE_METHOD_LAMBERT_CONIC_CONFORMAL_1SP to EPSG_CODE_METHOD_LAMBERT_CONIC_CONFORMAL_2SP

  • EPSG_CODE_METHOD_LAMBERT_CONIC_CONFORMAL_2SP to EPSG_CODE_METHOD_LAMBERT_CONIC_CONFORMAL_1SP

Return

new conversion, or nullptr

Parameters
  • targetEPSGCode: EPSG code of the target method.

Public Static Functions

ConversionNNPtr create(const util::PropertyMap &properties, const OperationMethodNNPtr &methodIn, const std::vector<GeneralParameterValueNNPtr> &values)

Instantiate a Conversion from a vector of GeneralParameterValue.

Return

a new Conversion.

Parameters
  • properties: See general_properties. At minimum the name should be defined.

  • methodIn: the operation method.

  • values: the values.

Exceptions

ConversionNNPtr create(const util::PropertyMap &propertiesConversion, const util::PropertyMap &propertiesOperationMethod, const std::vector<OperationParameterNNPtr> &parameters, const std::vector<ParameterValueNNPtr> &values)

Instantiate a Conversion and its OperationMethod.

Return

a new Conversion.

Parameters
  • propertiesConversion: See general_properties of the conversion. At minimum the name should be defined.

  • propertiesOperationMethod: See general_properties of the operation method. At minimum the name should be defined.

  • parameters: the operation parameters.

  • values: the operation values. Constraint: values.size() == parameters.size()

Exceptions

ConversionNNPtr createUTM(const util::PropertyMap &properties, int zone, bool north)

Instantiate a Universal Transverse Mercator conversion.

UTM is a family of conversions, of EPSG codes from 16001 to 16060 for the northern hemisphere, and 17001 to 17060 for the southern hemisphere, based on the Transverse Mercator projection method.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • zone: UTM zone number between 1 and 60.

  • north: true for UTM northern hemisphere, false for UTM southern hemisphere.

ConversionNNPtr createTransverseMercator(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Transverse Mercator projection method.

This method is defined as EPSG:9807

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • scale: See Scale Factor

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createGaussSchreiberTransverseMercator(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Gauss Schreiber Transverse Mercator projection method.

This method is also known as Gauss-Laborde Reunion.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • scale: See Scale Factor

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createTransverseMercatorSouthOriented(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Transverse Mercator South Orientated projection method.

This method is defined as EPSG:9808

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • scale: See Scale Factor

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createTwoPointEquidistant(const util::PropertyMap &properties, const common::Angle &latitudeFirstPoint, const common::Angle &longitudeFirstPoint, const common::Angle &latitudeSecondPoint, const common::Angle &longitudeSeconPoint, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Two Point Equidistant projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFirstPoint: Latitude of first point.

  • longitudeFirstPoint: Longitude of first point.

  • latitudeSecondPoint: Latitude of second point.

  • longitudeSeconPoint: Longitude of second point.

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createTunisiaMappingGrid(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Tunisia Mapping Grid projection method.

This method is defined as EPSG:9816

Note

There is currently no implementation of the method formulas in PROJ.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createAlbersEqualArea(const util::PropertyMap &properties, const common::Angle &latitudeFalseOrigin, const common::Angle &longitudeFalseOrigin, const common::Angle &latitudeFirstParallel, const common::Angle &latitudeSecondParallel, const common::Length &eastingFalseOrigin, const common::Length &northingFalseOrigin)

Instantiate a conversion based on the Albers Conic Equal Area projection method.

This method is defined as EPSG:9822

Note

the order of arguments is conformant with the corresponding EPSG mode and different than OGRSpatialReference::setACEA() of GDAL <= 2.3

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFalseOrigin: See Latitude of false origin

  • longitudeFalseOrigin: See Longitude of false origin

  • latitudeFirstParallel: See Latitude of 1st standard parallel

  • latitudeSecondParallel: See Latitude of 2nd standard parallel

  • eastingFalseOrigin: See Easting of false origin

  • northingFalseOrigin: See Northing of false origin

ConversionNNPtr createLambertConicConformal_1SP(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Lambert Conic Conformal 1SP projection method.

This method is defined as EPSG:9801

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • scale: See Scale Factor

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createLambertConicConformal_2SP(const util::PropertyMap &properties, const common::Angle &latitudeFalseOrigin, const common::Angle &longitudeFalseOrigin, const common::Angle &latitudeFirstParallel, const common::Angle &latitudeSecondParallel, const common::Length &eastingFalseOrigin, const common::Length &northingFalseOrigin)

Instantiate a conversion based on the Lambert Conic Conformal (2SP) projection method.

This method is defined as EPSG:9802

Note

the order of arguments is conformant with the corresponding EPSG mode and different than OGRSpatialReference::setLCC() of GDAL <= 2.3

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFalseOrigin: See Latitude of false origin

  • longitudeFalseOrigin: See Longitude of false origin

  • latitudeFirstParallel: See Latitude of 1st standard parallel

  • latitudeSecondParallel: See Latitude of 2nd standard parallel

  • eastingFalseOrigin: See Easting of false origin

  • northingFalseOrigin: See Northing of false origin

ConversionNNPtr createLambertConicConformal_2SP_Michigan(const util::PropertyMap &properties, const common::Angle &latitudeFalseOrigin, const common::Angle &longitudeFalseOrigin, const common::Angle &latitudeFirstParallel, const common::Angle &latitudeSecondParallel, const common::Length &eastingFalseOrigin, const common::Length &northingFalseOrigin, const common::Scale &ellipsoidScalingFactor)

Instantiate a conversion based on the Lambert Conic Conformal (2SP Michigan) projection method.

This method is defined as EPSG:1051

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFalseOrigin: See Latitude of false origin

  • longitudeFalseOrigin: See Longitude of false origin

  • latitudeFirstParallel: See Latitude of 1st standard parallel

  • latitudeSecondParallel: See Latitude of 2nd standard parallel

  • eastingFalseOrigin: See Easting of false origin

  • northingFalseOrigin: See Northing of false origin

  • ellipsoidScalingFactor: Ellipsoid scaling factor.

ConversionNNPtr createLambertConicConformal_2SP_Belgium(const util::PropertyMap &properties, const common::Angle &latitudeFalseOrigin, const common::Angle &longitudeFalseOrigin, const common::Angle &latitudeFirstParallel, const common::Angle &latitudeSecondParallel, const common::Length &eastingFalseOrigin, const common::Length &northingFalseOrigin)

Instantiate a conversion based on the Lambert Conic Conformal (2SP Belgium) projection method.

This method is defined as EPSG:9803

Warning

The formulas used currently in PROJ are, incorrectly, the ones of the regular LCC_2SP method.

Note

the order of arguments is conformant with the corresponding EPSG mode and different than OGRSpatialReference::setLCCB() of GDAL <= 2.3

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFalseOrigin: See Latitude of false origin

  • longitudeFalseOrigin: See Longitude of false origin

  • latitudeFirstParallel: See Latitude of 1st standard parallel

  • latitudeSecondParallel: See Latitude of 2nd standard parallel

  • eastingFalseOrigin: See Easting of false origin

  • northingFalseOrigin: See Northing of false origin

ConversionNNPtr createAzimuthalEquidistant(const util::PropertyMap &properties, const common::Angle &latitudeNatOrigin, const common::Angle &longitudeNatOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Modified Azimuthal Equidistant projection method.

This method is defined as EPSG:9832

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeNatOrigin: See Latitude of natural origin/Center Latitude

  • longitudeNatOrigin: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createGuamProjection(const util::PropertyMap &properties, const common::Angle &latitudeNatOrigin, const common::Angle &longitudeNatOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Guam Projection projection method.

This method is defined as EPSG:9831

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeNatOrigin: See Latitude of natural origin/Center Latitude

  • longitudeNatOrigin: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createBonne(const util::PropertyMap &properties, const common::Angle &latitudeNatOrigin, const common::Angle &longitudeNatOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Bonne projection method.

This method is defined as EPSG:9827

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeNatOrigin: See Latitude of natural origin/Center Latitude . PROJ calls its the standard parallel 1.

  • longitudeNatOrigin: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createLambertCylindricalEqualAreaSpherical(const util::PropertyMap &properties, const common::Angle &latitudeFirstParallel, const common::Angle &longitudeNatOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Lambert Cylindrical Equal Area (Spherical) projection method.

This method is defined as EPSG:9834

Warning

The PROJ cea computation code would select the ellipsoidal form if a non-spherical ellipsoid is used for the base GeographicalCRS.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFirstParallel: See Latitude of 1st standard parallel.

  • longitudeNatOrigin: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createLambertCylindricalEqualArea(const util::PropertyMap &properties, const common::Angle &latitudeFirstParallel, const common::Angle &longitudeNatOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Lambert Cylindrical Equal Area (ellipsoidal form) projection method.

This method is defined as EPSG:9835

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFirstParallel: See Latitude of 1st standard parallel.

  • longitudeNatOrigin: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createCassiniSoldner(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Cassini-Soldner projection method.

This method is defined as EPSG:9806

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEquidistantConic(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Angle &latitudeFirstParallel, const common::Angle &latitudeSecondParallel, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Equidistant Conic projection method.

There is no equivalent in EPSG.

Note

Although not found in EPSG, the order of arguments is conformant with the “spirit” of EPSG and different than OGRSpatialReference::setEC() of GDAL <= 2.3 *

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • latitudeFirstParallel: See Latitude of 1st standard parallel

  • latitudeSecondParallel: See Latitude of 2nd standard parallel

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEckertI(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Eckert I projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEckertII(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Eckert II projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEckertIII(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Eckert III projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEckertIV(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Eckert IV projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEckertV(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Eckert V projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEckertVI(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Eckert VI projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEquidistantCylindrical(const util::PropertyMap &properties, const common::Angle &latitudeFirstParallel, const common::Angle &longitudeNatOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Equidistant Cylindrical projection method.

This is also known as the Equirectangular method, and in the particular case where the latitude of first parallel is 0.

This method is defined as EPSG:1028

Note

This is the equivalent OGRSpatialReference::SetEquirectangular2( 0.0, latitudeFirstParallel, falseEasting, falseNorthing ) of GDAL <= 2.3, where the lat_0 / center_latitude parameter is forced to 0.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFirstParallel: See Latitude of 1st standard parallel.

  • longitudeNatOrigin: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createEquidistantCylindricalSpherical(const util::PropertyMap &properties, const common::Angle &latitudeFirstParallel, const common::Angle &longitudeNatOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Equidistant Cylindrical (Spherical) projection method.

This is also known as the Equirectangular method, and in the particular case where the latitude of first parallel is 0.

This method is defined as EPSG:1029

Note

This is the equivalent OGRSpatialReference::SetEquirectangular2( 0.0, latitudeFirstParallel, falseEasting, falseNorthing ) of GDAL <= 2.3, where the lat_0 / center_latitude parameter is forced to 0.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFirstParallel: See Latitude of 1st standard parallel.

  • longitudeNatOrigin: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createGall(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Gall (Stereographic) projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createGoodeHomolosine(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Goode Homolosine projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createInterruptedGoodeHomolosine(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Interrupted Goode Homolosine projection method.

There is no equivalent in EPSG.

Note

OGRSpatialReference::SetIGH() of GDAL <= 2.3 assumes the 3 projection parameters to be zero and this is the nominal case.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createGeostationarySatelliteSweepX(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &height, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Geostationary Satellite View projection method, with the sweep angle axis of the viewing instrument being x.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • height: Height of the view point above the Earth.

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createGeostationarySatelliteSweepY(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &height, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Geostationary Satellite View projection method, with the sweep angle axis of the viewing instrument being y.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • height: Height of the view point above the Earth.

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createGnomonic(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Gnomonic projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createHotineObliqueMercatorVariantA(const util::PropertyMap &properties, const common::Angle &latitudeProjectionCentre, const common::Angle &longitudeProjectionCentre, const common::Angle &azimuthInitialLine, const common::Angle &angleFromRectifiedToSkrewGrid, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Hotine Oblique Mercator (Variant A) projection method.

This is the variant with the no_uoff parameter, which corresponds to GDAL >=2.3 Hotine_Oblique_Mercator projection. In this variant, the false grid coordinates are defined at the intersection of the initial line and the aposphere (the equator on one of the intermediate surfaces inherent in the method), that is at the natural origin of the coordinate system).

This method is defined as EPSG:9812

Note

In the case where azimuthInitialLine = angleFromRectifiedToSkrewGrid = 90deg, this maps to the Swiss Oblique Mercator formulas.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeProjectionCentre: See Latitute of projection centre

  • longitudeProjectionCentre: See Longitude of projection centre

  • azimuthInitialLine: See Azimuth of initial line

  • angleFromRectifiedToSkrewGrid: See Angle from Rectified to Skew Grid

  • scale: See Scale factor on initial line

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createHotineObliqueMercatorVariantB(const util::PropertyMap &properties, const common::Angle &latitudeProjectionCentre, const common::Angle &longitudeProjectionCentre, const common::Angle &azimuthInitialLine, const common::Angle &angleFromRectifiedToSkrewGrid, const common::Scale &scale, const common::Length &eastingProjectionCentre, const common::Length &northingProjectionCentre)

Instantiate a conversion based on the Hotine Oblique Mercator (Variant B) projection method.

This is the variant without the no_uoff parameter, which corresponds to GDAL >=2.3 Hotine_Oblique_Mercator_Azimuth_Center projection. In this variant, the false grid coordinates are defined at the projection centre.

This method is defined as EPSG:9815

Note

In the case where azimuthInitialLine = angleFromRectifiedToSkrewGrid = 90deg, this maps to the Swiss Oblique Mercator formulas.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeProjectionCentre: See Latitute of projection centre

  • longitudeProjectionCentre: See Longitude of projection centre

  • azimuthInitialLine: See Azimuth of initial line

  • angleFromRectifiedToSkrewGrid: See Angle from Rectified to Skew Grid

  • scale: See Scale factor on initial line

  • eastingProjectionCentre: See Easting at projection centre

  • northingProjectionCentre: See Northing at projection centre

ConversionNNPtr createHotineObliqueMercatorTwoPointNaturalOrigin(const util::PropertyMap &properties, const common::Angle &latitudeProjectionCentre, const common::Angle &latitudePoint1, const common::Angle &longitudePoint1, const common::Angle &latitudePoint2, const common::Angle &longitudePoint2, const common::Scale &scale, const common::Length &eastingProjectionCentre, const common::Length &northingProjectionCentre)

Instantiate a conversion based on the Hotine Oblique Mercator Two Point Natural Origin projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeProjectionCentre: See Latitute of projection centre

  • latitudePoint1: Latitude of point 1.

  • longitudePoint1: Latitude of point 1.

  • latitudePoint2: Latitude of point 2.

  • longitudePoint2: Longitude of point 2.

  • scale: See Scale factor on initial line

  • eastingProjectionCentre: See Easting at projection centre

  • northingProjectionCentre: See Northing at projection centre

ConversionNNPtr createLabordeObliqueMercator(const util::PropertyMap &properties, const common::Angle &latitudeProjectionCentre, const common::Angle &longitudeProjectionCentre, const common::Angle &azimuthInitialLine, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Laborde Oblique Mercator projection method.

This method is defined as EPSG:9813

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeProjectionCentre: See Latitute of projection centre

  • longitudeProjectionCentre: See Longitude of projection centre

  • azimuthInitialLine: See Azimuth of initial line

  • scale: See Scale factor on initial line

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createInternationalMapWorldPolyconic(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Angle &latitudeFirstParallel, const common::Angle &latitudeSecondParallel, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the International Map of the World Polyconic projection method.

There is no equivalent in EPSG.

Note

the order of arguments is conformant with the corresponding EPSG mode and different than OGRSpatialReference::SetIWMPolyconic() of GDAL <= 2.3

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • latitudeFirstParallel: See Latitude of 1st standard parallel

  • latitudeSecondParallel: See Latitude of 2nd standard parallel

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createKrovakNorthOriented(const util::PropertyMap &properties, const common::Angle &latitudeProjectionCentre, const common::Angle &longitudeOfOrigin, const common::Angle &colatitudeConeAxis, const common::Angle &latitudePseudoStandardParallel, const common::Scale &scaleFactorPseudoStandardParallel, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Krovak (north oriented) projection method.

This method is defined as EPSG:1041

The coordinates are returned in the “GIS friendly” order: easting, northing. This method is similar to createKrovak(), except that the later returns projected values as southing, westing, where southing(Krovak) = -northing(Krovak_North) and westing(Krovak) = -easting(Krovak_North).

Note

The current PROJ implementation of Krovak hard-codes colatitudeConeAxis = 30deg17’17.30311” and latitudePseudoStandardParallel = 78deg30’N, which are the values used for the ProjectedCRS S-JTSK (Ferro) / Krovak East North (EPSG:5221). It also hard-codes the parameters of the Bessel ellipsoid typically used for Krovak.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeProjectionCentre: See Latitute of projection centre

  • longitudeOfOrigin: See Longitude of origin

  • colatitudeConeAxis: See Co-latitude of cone axis

  • latitudePseudoStandardParallel: See Latitude of pseudo standard

  • scaleFactorPseudoStandardParallel: See Scale factor on pseudo

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createKrovak(const util::PropertyMap &properties, const common::Angle &latitudeProjectionCentre, const common::Angle &longitudeOfOrigin, const common::Angle &colatitudeConeAxis, const common::Angle &latitudePseudoStandardParallel, const common::Scale &scaleFactorPseudoStandardParallel, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Krovak projection method.

This method is defined as EPSG:9819

The coordinates are returned in the historical order: southing, westing This method is similar to createKrovakNorthOriented(), except that the later returns projected values as easting, northing, where easting(Krovak_North) = -westing(Krovak) and northing(Krovak_North) = -southing(Krovak).

Note

The current PROJ implementation of Krovak hard-codes colatitudeConeAxis = 30deg17’17.30311” and latitudePseudoStandardParallel = 78deg30’N, which are the values used for the ProjectedCRS S-JTSK (Ferro) / Krovak East North (EPSG:5221). It also hard-codes the parameters of the Bessel ellipsoid typically used for Krovak.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeProjectionCentre: See Latitute of projection centre

  • longitudeOfOrigin: See Longitude of origin

  • colatitudeConeAxis: See Co-latitude of cone axis

  • latitudePseudoStandardParallel: See Latitude of pseudo standard

  • scaleFactorPseudoStandardParallel: See Scale factor on pseudo

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createLambertAzimuthalEqualArea(const util::PropertyMap &properties, const common::Angle &latitudeNatOrigin, const common::Angle &longitudeNatOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Lambert Azimuthal Equal Area projection method.

This method is defined as EPSG:9820

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeNatOrigin: See Latitude of natural origin/Center Latitude

  • longitudeNatOrigin: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createMillerCylindrical(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Miller Cylindrical projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createMercatorVariantA(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Mercator projection method.

This is the variant, also known as Mercator (1SP), defined with the scale factor. Note that latitude of natural origin (centerLat) is a parameter, but unused in the transformation formulas.

This method is defined as EPSG:9804

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude . Should be 0.

  • centerLong: See Longitude of natural origin/Central Meridian

  • scale: See Scale Factor

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createMercatorVariantB(const util::PropertyMap &properties, const common::Angle &latitudeFirstParallel, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Mercator projection method.

This is the variant, also known as Mercator (2SP), defined with the latitude of the first standard parallel (the second standard parallel is implicitly the opposite value). The latitude of natural origin is fixed to zero.

This method is defined as EPSG:9805

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeFirstParallel: See Latitude of 1st standard parallel

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createPopularVisualisationPseudoMercator(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Popular Visualisation Pseudo Mercator projection method.

Also known as WebMercator. Mostly/only used for Projected CRS EPSG:3857 (WGS 84 / Pseudo-Mercator)

This method is defined as EPSG:1024

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude . Usually 0

  • centerLong: See Longitude of natural origin/Central Meridian . Usually 0

  • falseEasting: See False Easting . Usually 0

  • falseNorthing: See False Northing . Usually 0

ConversionNNPtr createMollweide(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Mollweide projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createNewZealandMappingGrid(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the New Zealand Map Grid projection method.

This method is defined as EPSG:9811

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createObliqueStereographic(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Oblique Stereographic (Alternative) projection method.

This method is defined as EPSG:9809

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • scale: See Scale Factor

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createOrthographic(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Orthographic projection method.

This method is defined as EPSG:9840

Note

At the time of writing, PROJ only implements the spherical formulation

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createAmericanPolyconic(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the American Polyconic projection method.

This method is defined as EPSG:9818

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createPolarStereographicVariantA(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Polar Stereographic (Variant A) projection method.

This method is defined as EPSG:9810

This is the variant of polar stereographic defined with a scale factor.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude . Should be 90 deg ou -90 deg.

  • centerLong: See Longitude of natural origin/Central Meridian

  • scale: See Scale Factor

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createPolarStereographicVariantB(const util::PropertyMap &properties, const common::Angle &latitudeStandardParallel, const common::Angle &longitudeOfOrigin, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Polar Stereographic (Variant B) projection method.

This method is defined as EPSG:9829

This is the variant of polar stereographic defined with a latitude of standard parallel.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeStandardParallel: See Latitude of standard parallel

  • longitudeOfOrigin: See Longitude of origin

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createRobinson(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Robinson projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createSinusoidal(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Sinusoidal projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createStereographic(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Stereographic projection method.

There is no equivalent in EPSG. This method implements the original “Oblique

Stereographic” method described in “Snyder’s Map Projections - A Working

manual”, which is different from the “Oblique Stereographic (alternative”) method implemented in

createObliqueStereographic().

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • scale: See Scale Factor

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createVanDerGrinten(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Van der Grinten projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createWagnerI(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Wagner I projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createWagnerII(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Wagner II projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createWagnerIII(const util::PropertyMap &properties, const common::Angle &latitudeTrueScale, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Wagner III projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • latitudeTrueScale: Latitude of true scale.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createWagnerIV(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Wagner IV projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createWagnerV(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Wagner V projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createWagnerVI(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Wagner VI projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createWagnerVII(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Wagner VII projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createQuadrilateralizedSphericalCube(const util::PropertyMap &properties, const common::Angle &centerLat, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Quadrilateralized Spherical Cube projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLat: See Latitude of natural origin/Center Latitude

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createSphericalCrossTrackHeight(const util::PropertyMap &properties, const common::Angle &pegPointLat, const common::Angle &pegPointLong, const common::Angle &pegPointHeading, const common::Length &pegPointHeight)

Instantiate a conversion based on the Spherical Cross-Track Height projection method.

There is no equivalent in EPSG.

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • pegPointLat: Peg point latitude.

  • pegPointLong: Peg point longitude.

  • pegPointHeading: Peg point heading.

  • pegPointHeight: Peg point height.

ConversionNNPtr createEqualEarth(const util::PropertyMap &properties, const common::Angle &centerLong, const common::Length &falseEasting, const common::Length &falseNorthing)

Instantiate a conversion based on the Equal Earth projection method.

This method is defined as EPSG:1078

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • centerLong: See Longitude of natural origin/Central Meridian

  • falseEasting: See False Easting

  • falseNorthing: See False Northing

ConversionNNPtr createChangeVerticalUnit(const util::PropertyMap &properties, const common::Scale &factor)

Instantiate a conversion based on the Change of Vertical Unit method.

This method is defined as EPSG:1069

Return

a new Conversion.

Parameters

ConversionNNPtr createAxisOrderReversal(bool is3D)

Instantiate a conversion based on the Axis order reversal method.

This swaps the longitude, latitude axis.

This method is defined as EPSG:9843, or for 3D as EPSG:9844

Return

a new Conversion.

Parameters
  • is3D: Whether this should apply on 3D geographicCRS

ConversionNNPtr createGeographicGeocentric(const util::PropertyMap &properties)

Instantiate a conversion based on the Geographic/Geocentric method.

This method is defined as EPSG:9602,

Return

a new Conversion.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

class CoordinateOperation : public osgeo::proj::common::ObjectUsage, public osgeo::proj::io::IPROJStringExportable, public osgeo::proj::io::IJSONExportable
#include <coordinateoperation.hpp>

Abstract class for a mathematical operation on coordinates.

A mathematical operation:

  • on coordinates that transforms or converts them from one coordinate reference system to another coordinate reference system

  • or that describes the change of coordinate values within one coordinate reference system due to the motion of the point between one coordinate epoch and another coordinate epoch.

Many but not all coordinate operations (from CRS A to CRS B) also uniquely define the inverse coordinate operation (from CRS B to CRS A). In some cases, the coordinate operation method algorithm for the inverse coordinate operation is the same as for the forward algorithm, but the signs of some coordinate operation parameter values have to be reversed. In other cases, different algorithms are required for the forward and inverse coordinate operations, but the same coordinate operation parameter values are used. If (some) entirely different parameter values are needed, a different coordinate operation shall be defined.

Remark

Implements CoordinateOperation from ISO_19111_2019

Subclassed by osgeo::proj::operation::ConcatenatedOperation, osgeo::proj::operation::SingleOperation

Public Functions

const util::optional<std::string> &operationVersion() const

Return the version of the coordinate transformation (i.e. instantiation due to the stochastic nature of the parameters).

Mandatory when describing a coordinate transformation or point motion operation, and should not be supplied for a coordinate conversion.

Return

version or empty.

const std::vector<metadata::PositionalAccuracyNNPtr> &coordinateOperationAccuracies() const

Return estimate(s) of the impact of this coordinate operation on point accuracy.

Gives position error estimates for target coordinates of this coordinate operation, assuming no errors in source coordinates.

Return

estimate(s) or empty vector.

const crs::CRSPtr sourceCRS() const

Return the source CRS of this coordinate operation.

This should not be null, expect for of a derivingConversion of a DerivedCRS when the owning DerivedCRS has been destroyed.

Return

source CRS, or null.

const crs::CRSPtr targetCRS() const

Return the target CRS of this coordinate operation.

This should not be null, expect for of a derivingConversion of a DerivedCRS when the owning DerivedCRS has been destroyed.

Return

target CRS, or null.

const crs::CRSPtr &interpolationCRS() const

Return the interpolation CRS of this coordinate operation.

Return

interpolation CRS, or null.

const util::optional<common::DataEpoch> &sourceCoordinateEpoch() const

Return the source epoch of coordinates.

Return

source epoch of coordinates, or empty.

const util::optional<common::DataEpoch> &targetCoordinateEpoch() const

Return the target epoch of coordinates.

Return

target epoch of coordinates, or empty.

virtual CoordinateOperationNNPtr inverse() const = 0

Return the inverse of the coordinate operation.

Exceptions

virtual std::set<GridDescription> gridsNeeded(const io::DatabaseContextPtr &databaseContext) const = 0

Return grids needed by an operation.

bool isPROJInstantiable(const io::DatabaseContextPtr &databaseContext) const

Return whether a coordinate operation can be instantiated as a PROJ pipeline, checking in particular that referenced grids are available.

bool hasBallparkTransformation() const

Return whether a coordinate operation has a “ballpark” transformation, that is a very approximate one, due to lack of more accurate transformations.

Typically a null geographic offset between two horizontal datum, or a null vertical offset (or limited to unit changes) between two vertical datum. Errors of several tens to one hundred meters might be expected, compared to more accurate transformations.

CoordinateOperationNNPtr normalizeForVisualization() const

Return a variation of the current coordinate operation whose axis order is the one expected for visualization purposes.

Public Static Attributes

const std::string OPERATION_VERSION_KEY

Key to set the operation version of a operation::CoordinateOperation.

The value is to be provided as a string.

class CoordinateOperationContext
#include <coordinateoperation.hpp>

Context in which a coordinate operation is to be used.

Remark

Implements [CoordinateOperationFactory https://sis.apache.org/apidocs/org/apache/sis/referencing/operation/CoordinateOperationContext.html] from Apache SIS

Public Types

enum SourceTargetCRSExtentUse

Specify how source and target CRS extent should be used to restrict candidate operations (only taken into account if no explicit area of interest is specified.

Values:

NONE

Ignore CRS extent

BOTH

Test coordinate operation extent against both CRS extent.

INTERSECTION

Test coordinate operation extent against the intersection of both CRS extent.

SMALLEST

Test coordinate operation against the smallest of both CRS extent.

enum SpatialCriterion

Spatial criterion to restrict candidate operations.

Values:

STRICT_CONTAINMENT

The area of validity of transforms should strictly contain the are of interest.

PARTIAL_INTERSECTION

The area of validity of transforms should at least intersect the area of interest.

enum GridAvailabilityUse

Describe how grid availability is used.

Values:

USE_FOR_SORTING

Grid availability is only used for sorting results. Operations where some grids are missing will be sorted last.

DISCARD_OPERATION_IF_MISSING_GRID

Completely discard an operation if a required grid is missing.

IGNORE_GRID_AVAILABILITY

Ignore grid availability at all. Results will be presented as if all grids were available.

enum IntermediateCRSUse

Describe if and how intermediate CRS should be used

Values:

ALWAYS

Always search for intermediate CRS.

IF_NO_DIRECT_TRANSFORMATION

Only attempt looking for intermediate CRS if there is no direct transformation available.

NEVER

Public Functions

const io::AuthorityFactoryPtr &getAuthorityFactory() const

Return the authority factory, or null.

const metadata::ExtentPtr &getAreaOfInterest() const

Return the desired area of interest, or null.

void setAreaOfInterest(const metadata::ExtentPtr &extent)

Set the desired area of interest, or null.

double getDesiredAccuracy() const

Return the desired accuracy (in metre), or 0.

void setDesiredAccuracy(double accuracy)

Set the desired accuracy (in metre), or 0.

void setSourceAndTargetCRSExtentUse(SourceTargetCRSExtentUse use)

Set how source and target CRS extent should be used when considering if a transformation can be used (only takes effect if no area of interest is explicitly defined).

The default is CoordinateOperationContext::SourceTargetCRSExtentUse::SMALLEST.

CoordinateOperationContext::SourceTargetCRSExtentUse getSourceAndTargetCRSExtentUse() const

Return how source and target CRS extent should be used when considering if a transformation can be used (only takes effect if no area of interest is explicitly defined).

The default is CoordinateOperationContext::SourceTargetCRSExtentUse::SMALLEST.

void setSpatialCriterion(SpatialCriterion criterion)

Set the spatial criterion to use when comparing the area of validity of coordinate operations with the area of interest / area of validity of source and target CRS.

The default is STRICT_CONTAINMENT.

CoordinateOperationContext::SpatialCriterion getSpatialCriterion() const

Return the spatial criterion to use when comparing the area of validity of coordinate operations with the area of interest / area of validity of source and target CRS.

The default is STRICT_CONTAINMENT.

void setUsePROJAlternativeGridNames(bool usePROJNames)

Set whether PROJ alternative grid names should be substituted to the official authority names.

This only has effect is an authority factory with a non-null database context has been attached to this context.

If set to false, it is still possible to obtain later the substitution by using io::PROJStringFormatter::create() with a non-null database context.

The default is true.

bool getUsePROJAlternativeGridNames() const

Return whether PROJ alternative grid names should be substituted to the official authority names.

The default is true.

void setDiscardSuperseded(bool discard)

Set whether transformations that are superseded (but not deprecated) should be discarded.

The default is true.

bool getDiscardSuperseded() const

Return whether transformations that are superseded (but not deprecated) should be discarded.

The default is true.

void setGridAvailabilityUse(GridAvailabilityUse use)

Set how grid availability is used.

The default is USE_FOR_SORTING.

CoordinateOperationContext::GridAvailabilityUse getGridAvailabilityUse() const

Return how grid availability is used.

The default is USE_FOR_SORTING.

void setAllowUseIntermediateCRS(IntermediateCRSUse use)

Set whether an intermediate pivot CRS can be used for researching coordinate operations between a source and target CRS.

Concretely if in the database there is an operation from A to C (or C to A), and another one from C to B (or B to C), but no direct operation between A and B, setting this parameter to ALWAYS/IF_NO_DIRECT_TRANSFORMATION, allow chaining both operations.

The current implementation is limited to researching one intermediate step.

By default, with the IF_NO_DIRECT_TRANSFORMATION stratgey, all potential C candidates will be used if there is no direct tranformation.

CoordinateOperationContext::IntermediateCRSUse getAllowUseIntermediateCRS() const

Return whether an intermediate pivot CRS can be used for researching coordinate operations between a source and target CRS.

Concretely if in the database there is an operation from A to C (or C to A), and another one from C to B (or B to C), but no direct operation between A and B, setting this parameter to ALWAYS/IF_NO_DIRECT_TRANSFORMATION, allow chaining both operations.

The default is IF_NO_DIRECT_TRANSFORMATION.

void setIntermediateCRS(const std::vector<std::pair<std::string, std::string>> &intermediateCRSAuthCodes)

Restrict the potential pivot CRSs that can be used when trying to build a coordinate operation between two CRS that have no direct operation.

Parameters
  • intermediateCRSAuthCodes: a vector of (auth_name, code) that can be used as potential pivot RS

const std::vector<std::pair<std::string, std::string>> &getIntermediateCRS() const

Return the potential pivot CRSs that can be used when trying to build a coordinate operation between two CRS that have no direct operation.

Public Static Functions

CoordinateOperationContextNNPtr create(const io::AuthorityFactoryPtr &authorityFactory, const metadata::ExtentPtr &extent, double accuracy)

Creates a context for a coordinate operation.

If a non null authorityFactory is provided, the resulting context should not be used simultaneously by more than one thread.

If authorityFactory->getAuthority() is the empty string, then coordinate operations from any authority will be searched, with the restrictions set in the authority_to_authority_preference database table. If authorityFactory->getAuthority() is set to “any”, then coordinate operations from any authority will be searched If authorityFactory->getAuthority() is a non-empty string different of “any”, then coordinate operatiosn will be searched only in that authority namespace.

Return

a new context.

Parameters
  • authorityFactory: Authority factory, or null if no database lookup is allowed. Use io::authorityFactory::create(context, std::string()) to allow all authorities to be used.

  • extent: Area of interest, or null if none is known.

  • accuracy: Maximum allowed accuracy in metre, as specified in or 0 to get best accuracy.

class CoordinateOperationFactory
#include <coordinateoperation.hpp>

Creates coordinate operations. This factory is capable to find coordinate transformations or conversions between two coordinate reference systems.

Remark

Implements (partially) CoordinateOperationFactory from GeoAPI

Public Functions

CoordinateOperationPtr createOperation(const crs::CRSNNPtr &sourceCRS, const crs::CRSNNPtr &targetCRS) const

Find a CoordinateOperation from sourceCRS to targetCRS.

This is a helper of createOperations(), using a coordinate operation context with no authority factory (so no catalog searching is done), no desired accuracy and no area of interest. This returns the first operation of the result set of createOperations(), or null if none found.

Return

a CoordinateOperation or nullptr.

Parameters
  • sourceCRS: source CRS.

  • targetCRS: source CRS.

std::vector<CoordinateOperationNNPtr> createOperations(const crs::CRSNNPtr &sourceCRS, const crs::CRSNNPtr &targetCRS, const CoordinateOperationContextNNPtr &context) const

Find a list of CoordinateOperation from sourceCRS to targetCRS.

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.

Return

a list

Parameters
  • sourceCRS: source CRS.

  • targetCRS: source CRS.

  • context: Search context.

Public Static Functions

CoordinateOperationFactoryNNPtr create()

Instantiate a CoordinateOperationFactory.

class GeneralOperationParameter : public osgeo::proj::common::IdentifiedObject
#include <coordinateoperation.hpp>

Abstract class modelling a parameter value (OperationParameter) or group of parameters.

Remark

Implements GeneralOperationParameter from ISO_19111_2019

Subclassed by osgeo::proj::operation::OperationParameter

class GeneralParameterValue : public osgeo::proj::util::BaseObject, public osgeo::proj::io::IWKTExportable, public osgeo::proj::io::IJSONExportable, public osgeo::proj::util::IComparable
#include <coordinateoperation.hpp>

Abstract class modelling a parameter value (OperationParameterValue) or group of parameter values.

Remark

Implements GeneralParameterValue from ISO_19111_2019

Subclassed by osgeo::proj::operation::OperationParameterValue

struct GridDescription
#include <coordinateoperation.hpp>

Grid description.

Public Members

std::string shortName

Grid short filename

std::string fullName

Grid full path name (if found)

std::string packageName

Package name (or empty)

std::string url

Grid URL (if packageName is empty), or package URL (or empty)

bool directDownload

Whether url can be fetched directly.

bool openLicense

Whether the grid is released with an open license.

bool available

Whether GRID is available.

class InvalidOperation : public osgeo::proj::util::Exception
#include <coordinateoperation.hpp>

Exception that can be thrown when an invalid operation is attempted to be constructed.

class OperationMethod : public osgeo::proj::common::IdentifiedObject, public osgeo::proj::io::IJSONExportable
#include <coordinateoperation.hpp>

The method (algorithm or procedure) used to perform the coordinate operation.

For a projection method, this contains the name of the projection method and the name of the projection parameters.

Remark

Implements OperationMethod from ISO_19111_2019

Public Functions

const util::optional<std::string> &formula()

Return the formula(s) or procedure used by this coordinate operation method.

This may be a reference to a publication (in which case use formulaCitation()).

Note that the operation method may not be analytic, in which case this attribute references or contains the procedure, not an analytic formula.

Return

the formula, or empty.

const util::optional<metadata::Citation> &formulaCitation()

Return a reference to a publication giving the formula(s) or procedure used by the coordinate operation method.

Return

the formula citation, or empty.

const std::vector<GeneralOperationParameterNNPtr> &parameters()

Return the parameters of this operation method.

Return

the parameters.

int getEPSGCode()

Return the EPSG code, either directly, or through the name.

Return

code, or 0 if not found

Public Static Functions

OperationMethodNNPtr create(const util::PropertyMap &properties, const std::vector<GeneralOperationParameterNNPtr> &parameters)

Instantiate a operation method from a vector of GeneralOperationParameter.

Return

a new OperationMethod.

Parameters
  • properties: See general_properties. At minimum the name should be defined.

  • parameters: Vector of GeneralOperationParameterNNPtr.

OperationMethodNNPtr create(const util::PropertyMap &properties, const std::vector<OperationParameterNNPtr> &parameters)

Instantiate a operation method from a vector of OperationParameter.

Return

a new OperationMethod.

Parameters
  • properties: See general_properties. At minimum the name should be defined.

  • parameters: Vector of OperationParameterNNPtr.

class OperationParameter : public osgeo::proj::operation::GeneralOperationParameter
#include <coordinateoperation.hpp>

The definition of a parameter used by a coordinate operation method.

Most parameter values are numeric, but other types of parameter values are possible.

Remark

Implements OperationParameter from ISO_19111_2019

Public Functions

int getEPSGCode()

Return the EPSG code, either directly, or through the name.

Return

code, or 0 if not found

Public Static Functions

OperationParameterNNPtr create(const util::PropertyMap &properties)

Instantiate a OperationParameter.

Return

a new OperationParameter.

Parameters

const char *getNameForEPSGCode(int epsg_code)

Return the name of a parameter designed by its EPSG code.

Return

name, or nullptr if not found

class OperationParameterValue : public osgeo::proj::operation::GeneralParameterValue
#include <coordinateoperation.hpp>

A parameter value, ordered sequence of values, or reference to a file of parameter values.

This combines a OperationParameter with the corresponding ParameterValue.

Remark

Implements OperationParameterValue from ISO_19111_2019

Public Functions

const OperationParameterNNPtr &parameter()

Return the parameter (definition)

Return

the parameter (definition).

const ParameterValueNNPtr &parameterValue()

Return the parameter value.

Return

the parameter value.

Public Static Functions

OperationParameterValueNNPtr create(const OperationParameterNNPtr &parameterIn, const ParameterValueNNPtr &valueIn)

Instantiate a OperationParameterValue.

Return

a new OperationParameterValue.

Parameters
  • parameterIn: Parameter (definition).

  • valueIn: Parameter value.

class ParameterValue : public osgeo::proj::util::BaseObject, public osgeo::proj::io::IWKTExportable, public osgeo::proj::util::IComparable
#include <coordinateoperation.hpp>

The value of the coordinate operation parameter.

Most parameter values are numeric, but other types of parameter values are possible.

Remark

Implements ParameterValue from ISO_19111_2019

Public Types

enum Type

Type of the value.

Values:

MEASURE

Measure (i.e. value with a unit)

STRING

String

INTEGER

Integer

BOOLEAN

Boolean

FILENAME

Filename

Public Functions

const ParameterValue::Type &type()

Returns the type of a parameter value.

Return

the type.

const common::Measure &value()

Returns the value as a Measure (assumes type() == Type::MEASURE)

Return

the value as a Measure.

const std::string &stringValue()

Returns the value as a string (assumes type() == Type::STRING)

Return

the value as a string.

const std::string &valueFile()

Returns the value as a filename (assumes type() == Type::FILENAME)

Return

the value as a filename.

int integerValue()

Returns the value as a integer (assumes type() == Type::INTEGER)

Return

the value as a integer.

bool booleanValue()

Returns the value as a boolean (assumes type() == Type::BOOLEAN)

Return

the value as a boolean.

Public Static Functions

ParameterValueNNPtr create(const common::Measure &measureIn)

Instantiate a ParameterValue from a Measure (i.e. a value associated with a unit)

Return

a new ParameterValue.

ParameterValueNNPtr create(const char *stringValueIn)

Instantiate a ParameterValue from a string value.

Return

a new ParameterValue.

ParameterValueNNPtr create(const std::string &stringValueIn)

Instantiate a ParameterValue from a string value.

Return

a new ParameterValue.

ParameterValueNNPtr create(int integerValueIn)

Instantiate a ParameterValue from a integer value.

Return

a new ParameterValue.

ParameterValueNNPtr create(bool booleanValueIn)

Instantiate a ParameterValue from a boolean value.

Return

a new ParameterValue.

ParameterValueNNPtr createFilename(const std::string &stringValueIn)

Instantiate a ParameterValue from a filename.

Return

a new ParameterValue.

class PointMotionOperation : public osgeo::proj::operation::SingleOperation
#include <coordinateoperation.hpp>

A mathematical operation that describes the change of coordinate values within one coordinate reference system due to the motion of the point between one coordinate epoch and another coordinate epoch.

The motion is due to tectonic plate movement or deformation.

Remark

Implements PointMotionOperation from ISO_19111_2019

class SingleOperation : public virtual osgeo::proj::operation::CoordinateOperation
#include <coordinateoperation.hpp>

A single (not concatenated) coordinate operation (CoordinateOperation)

Remark

Implements SingleOperation from ISO_19111_2019

Subclassed by osgeo::proj::operation::Conversion, osgeo::proj::operation::PointMotionOperation, osgeo::proj::operation::Transformation

Public Functions

const std::vector<GeneralParameterValueNNPtr> &parameterValues()

Return the parameter values.

Return

the parameter values.

const OperationMethodNNPtr &method()

Return the operation method associated to the operation.

Return

the operation method.

const ParameterValuePtr &parameterValue(const std::string &paramName, int epsg_code = 0) const

Return the parameter value corresponding to a parameter name or EPSG code.

Return

the value, or nullptr if not found.

Parameters
  • paramName: the parameter name (or empty, in which case epsg_code should be non zero)

  • epsg_code: the parameter EPSG code (possibly zero)

const ParameterValuePtr &parameterValue(int epsg_code) const

Return the parameter value corresponding to a EPSG code.

Return

the value, or nullptr if not found.

Parameters
  • epsg_code: the parameter EPSG code

const common::Measure &parameterValueMeasure(const std::string &paramName, int epsg_code = 0) const

Return the parameter value, as a measure, corresponding to a parameter name or EPSG code.

Return

the measure, or the empty Measure() object if not found.

Parameters
  • paramName: the parameter name (or empty, in which case epsg_code should be non zero)

  • epsg_code: the parameter EPSG code (possibly zero)

const common::Measure &parameterValueMeasure(int epsg_code) const

Return the parameter value, as a measure, corresponding to a EPSG code.

Return

the measure, or the empty Measure() object if not found.

Parameters
  • epsg_code: the parameter EPSG code

std::set<GridDescription> gridsNeeded(const io::DatabaseContextPtr &databaseContext) const

Return grids needed by an operation.

std::list<std::string> validateParameters() const

Validate the parameters used by a coordinate operation.

Return whether the method is known or not, or a list of missing or extra parameters for the operations recognized by this implementation.

Public Static Functions

SingleOperationNNPtr createPROJBased(const util::PropertyMap &properties, const std::string &PROJString, const crs::CRSPtr &sourceCRS, const crs::CRSPtr &targetCRS, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies = std::vector<metadata::PositionalAccuracyNNPtr>())

Instantiate a PROJ-based single operation.

Note

The operation might internally be a pipeline chaining several operations. The use of the SingleOperation modeling here is mostly to be able to get the PROJ string as a parameter.

Return

the new instance

Parameters
  • properties: Properties

  • PROJString: the PROJ string.

  • sourceCRS: source CRS (might be null).

  • targetCRS: target CRS (might be null).

  • accuracies: Vector of positional accuracy (might be empty).

class Transformation : public osgeo::proj::operation::SingleOperation
#include <coordinateoperation.hpp>

A mathematical operation on coordinates in which parameters are empirically derived from data containing the coordinates of a series of points in both coordinate reference systems.

This computational process is usually “over-determined”, allowing derivation of error (or accuracy) estimates for the coordinate transformation. Also, the stochastic nature of the parameters may result in multiple (different) versions of the same coordinate transformations between the same source and target CRSs. Any single coordinate operation in which the input and output coordinates are referenced to different datums (reference frames) will be a coordinate transformation.

Remark

Implements Transformation from ISO_19111_2019

Public Functions

const crs::CRSNNPtr &sourceCRS()

Return the source crs::CRS of the transformation.

Return

the source CRS.

const crs::CRSNNPtr &targetCRS()

Return the target crs::CRS of the transformation.

Return

the target CRS.

CoordinateOperationNNPtr inverse() const

Return the inverse of the coordinate operation.

Exceptions

TransformationNNPtr substitutePROJAlternativeGridNames(io::DatabaseContextNNPtr databaseContext) const

Return an equivalent transformation to the current one, but using PROJ alternative grid names.

Public Static Functions

TransformationNNPtr create(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const crs::CRSPtr &interpolationCRSIn, const OperationMethodNNPtr &methodIn, const std::vector<GeneralParameterValueNNPtr> &values, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation from a vector of GeneralParameterValue.

Return

new Transformation.

Parameters
  • properties: See general_properties. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • interpolationCRSIn: Interpolation CRS (might be null)

  • methodIn: Operation method.

  • values: Vector of GeneralOperationParameterNNPtr.

  • accuracies: Vector of positional accuracy (might be empty).

Exceptions

TransformationNNPtr create(const util::PropertyMap &propertiesTransformation, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const crs::CRSPtr &interpolationCRSIn, const util::PropertyMap &propertiesOperationMethod, const std::vector<OperationParameterNNPtr> &parameters, const std::vector<ParameterValueNNPtr> &values, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation ands its OperationMethod.

Return

new Transformation.

Parameters
  • propertiesTransformation: The general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • interpolationCRSIn: Interpolation CRS (might be null)

  • propertiesOperationMethod: The general_properties of the OperationMethod. At minimum the name should be defined.

  • parameters: Vector of parameters of the operation method.

  • values: Vector of ParameterValueNNPtr. Constraint: values.size() == parameters.size()

  • accuracies: Vector of positional accuracy (might be empty).

Exceptions

TransformationNNPtr createGeocentricTranslations(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, double translationXMetre, double translationYMetre, double translationZMetre, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with Geocentric Translations method.

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • translationXMetre: Value of the Translation_X parameter (in metre).

  • translationYMetre: Value of the Translation_Y parameter (in metre).

  • translationZMetre: Value of the Translation_Z parameter (in metre).

  • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createPositionVector(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, double translationXMetre, double translationYMetre, double translationZMetre, double rotationXArcSecond, double rotationYArcSecond, double rotationZArcSecond, double scaleDifferencePPM, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with Position vector transformation method.

This is similar to createCoordinateFrameRotation(), except that the sign of the rotation terms is inverted.

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • translationXMetre: Value of the Translation_X parameter (in metre).

  • translationYMetre: Value of the Translation_Y parameter (in metre).

  • translationZMetre: Value of the Translation_Z parameter (in metre).

  • rotationXArcSecond: Value of the Rotation_X parameter (in arc-second).

  • rotationYArcSecond: Value of the Rotation_Y parameter (in arc-second).

  • rotationZArcSecond: Value of the Rotation_Z parameter (in arc-second).

  • scaleDifferencePPM: Value of the Scale_Difference parameter (in parts-per-million).

  • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createCoordinateFrameRotation(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, double translationXMetre, double translationYMetre, double translationZMetre, double rotationXArcSecond, double rotationYArcSecond, double rotationZArcSecond, double scaleDifferencePPM, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with Coordinate Frame Rotation method.

This is similar to createPositionVector(), except that the sign of the rotation terms is inverted.

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • translationXMetre: Value of the Translation_X parameter (in metre).

  • translationYMetre: Value of the Translation_Y parameter (in metre).

  • translationZMetre: Value of the Translation_Z parameter (in metre).

  • rotationXArcSecond: Value of the Rotation_X parameter (in arc-second).

  • rotationYArcSecond: Value of the Rotation_Y parameter (in arc-second).

  • rotationZArcSecond: Value of the Rotation_Z parameter (in arc-second).

  • scaleDifferencePPM: Value of the Scale_Difference parameter (in parts-per-million).

  • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createTimeDependentPositionVector(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, double translationXMetre, double translationYMetre, double translationZMetre, double rotationXArcSecond, double rotationYArcSecond, double rotationZArcSecond, double scaleDifferencePPM, double rateTranslationX, double rateTranslationY, double rateTranslationZ, double rateRotationX, double rateRotationY, double rateRotationZ, double rateScaleDifference, double referenceEpochYear, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with Time Dependent position vector transformation method.

This is similar to createTimeDependentCoordinateFrameRotation(), except that the sign of the rotation terms is inverted.

This method is defined as EPSG:1053

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • translationXMetre: Value of the Translation_X parameter (in metre).

  • translationYMetre: Value of the Translation_Y parameter (in metre).

  • translationZMetre: Value of the Translation_Z parameter (in metre).

  • rotationXArcSecond: Value of the Rotation_X parameter (in arc-second).

  • rotationYArcSecond: Value of the Rotation_Y parameter (in arc-second).

  • rotationZArcSecond: Value of the Rotation_Z parameter (in arc-second).

  • scaleDifferencePPM: Value of the Scale_Difference parameter (in parts-per-million).

  • rateTranslationX: Value of the rate of change of X-axis translation (in metre/year)

  • rateTranslationY: Value of the rate of change of Y-axis translation (in metre/year)

  • rateTranslationZ: Value of the rate of change of Z-axis translation (in metre/year)

  • rateRotationX: Value of the rate of change of X-axis rotation (in arc-second/year)

  • rateRotationY: Value of the rate of change of Y-axis rotation (in arc-second/year)

  • rateRotationZ: Value of the rate of change of Z-axis rotation (in arc-second/year)

  • rateScaleDifference: Value of the rate of change of scale difference (in PPM/year)

  • referenceEpochYear: Parameter reference epoch (in decimal year)

  • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createTimeDependentCoordinateFrameRotation(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, double translationXMetre, double translationYMetre, double translationZMetre, double rotationXArcSecond, double rotationYArcSecond, double rotationZArcSecond, double scaleDifferencePPM, double rateTranslationX, double rateTranslationY, double rateTranslationZ, double rateRotationX, double rateRotationY, double rateRotationZ, double rateScaleDifference, double referenceEpochYear, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with Time Dependent Position coordinate frame rotation transformation method.

This is similar to createTimeDependentPositionVector(), except that the sign of the rotation terms is inverted.

This method is defined as EPSG:1056

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • translationXMetre: Value of the Translation_X parameter (in metre).

  • translationYMetre: Value of the Translation_Y parameter (in metre).

  • translationZMetre: Value of the Translation_Z parameter (in metre).

  • rotationXArcSecond: Value of the Rotation_X parameter (in arc-second).

  • rotationYArcSecond: Value of the Rotation_Y parameter (in arc-second).

  • rotationZArcSecond: Value of the Rotation_Z parameter (in arc-second).

  • scaleDifferencePPM: Value of the Scale_Difference parameter (in parts-per-million).

  • rateTranslationX: Value of the rate of change of X-axis translation (in metre/year)

  • rateTranslationY: Value of the rate of change of Y-axis translation (in metre/year)

  • rateTranslationZ: Value of the rate of change of Z-axis translation (in metre/year)

  • rateRotationX: Value of the rate of change of X-axis rotation (in arc-second/year)

  • rateRotationY: Value of the rate of change of Y-axis rotation (in arc-second/year)

  • rateRotationZ: Value of the rate of change of Z-axis rotation (in arc-second/year)

  • rateScaleDifference: Value of the rate of change of scale difference (in PPM/year)

  • referenceEpochYear: Parameter reference epoch (in decimal year)

  • accuracies: Vector of positional accuracy (might be empty).

Exceptions

TransformationNNPtr createTOWGS84(const crs::CRSNNPtr &sourceCRSIn, const std::vector<double> &TOWGS84Parameters)

Instantiate a transformation from TOWGS84 parameters.

This is a helper of createPositionVector() with the source CRS being the GeographicCRS of sourceCRSIn, and the target CRS being EPSG:4326

Return

new Transformation.

Parameters
  • sourceCRSIn: Source CRS.

  • TOWGS84Parameters: The vector of 3 double values (Translation_X,_Y,_Z) or 7 double values (Translation_X,_Y,_Z, Rotation_X,_Y,_Z, Scale_Difference) passed to createPositionVector()

Exceptions

TransformationNNPtr createNTv2(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const std::string &filename, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with NTv2 method.

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • filename: NTv2 filename.

  • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createMolodensky(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, double translationXMetre, double translationYMetre, double translationZMetre, double semiMajorAxisDifferenceMetre, double flattingDifference, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with Molodensky method.

This method is defined as

EPSG:9604
See

createAbridgedMolodensky() for a related method.

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • translationXMetre: Value of the Translation_X parameter (in metre).

  • translationYMetre: Value of the Translation_Y parameter (in metre).

  • translationZMetre: Value of the Translation_Z parameter (in metre).

  • semiMajorAxisDifferenceMetre: The difference between the semi-major axis values of the ellipsoids used in the target and source CRS (in metre).

  • flattingDifference: The difference between the flattening values of the ellipsoids used in the target and source CRS.

  • accuracies: Vector of positional accuracy (might be empty).

Exceptions

TransformationNNPtr createAbridgedMolodensky(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, double translationXMetre, double translationYMetre, double translationZMetre, double semiMajorAxisDifferenceMetre, double flattingDifference, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with Abridged Molodensky method.

This method is defined as

EPSG:9605
See

createdMolodensky() for a related method.

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • translationXMetre: Value of the Translation_X parameter (in metre).

  • translationYMetre: Value of the Translation_Y parameter (in metre).

  • translationZMetre: Value of the Translation_Z parameter (in metre).

  • semiMajorAxisDifferenceMetre: The difference between the semi-major axis values of the ellipsoids used in the target and source CRS (in metre).

  • flattingDifference: The difference between the flattening values of the ellipsoids used in the target and source CRS.

  • accuracies: Vector of positional accuracy (might be empty).

Exceptions

TransformationNNPtr createGravityRelatedHeightToGeographic3D(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const crs::CRSPtr &interpolationCRSIn, const std::string &filename, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation from GravityRelatedHeight to Geographic3D.

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • interpolationCRSIn: Interpolation CRS. (might be null)

  • filename: GRID filename.

  • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createVERTCON(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const std::string &filename, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with method VERTCON.

Return

new Transformation.

Parameters
  • properties: See general_properties of the Transformation. At minimum the name should be defined.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • filename: GRID filename.

  • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createLongitudeRotation(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Angle &offset)

Instantiate a transformation with method Longitude rotation.

This method is defined as EPSG:9601

TransformationNNPtr createGeographic2DOffsets(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Angle &offsetLat, const common::Angle &offsetLon, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with method Geographic 2D offsets.

This method is defined as EPSG:9619

  • Return

    new Transformation.

    Parameters
    • properties: See general_properties of the Transformation. At minimum the name should be defined.

    • sourceCRSIn: Source CRS.

    • targetCRSIn: Target CRS.

    • offsetLat: Latitude offset to add.

    • offsetLon: Longitude offset to add.

    • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createGeographic3DOffsets(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Angle &offsetLat, const common::Angle &offsetLon, const common::Length &offsetHeight, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with method Geographic 3D offsets.

This method is defined as EPSG:9660

  • Return

    new Transformation.

    Parameters
    • properties: See general_properties of the Transformation. At minimum the name should be defined.

    • sourceCRSIn: Source CRS.

    • targetCRSIn: Target CRS.

    • offsetLat: Latitude offset to add.

    • offsetLon: Longitude offset to add.

    • offsetHeight: Height offset to add.

    • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createGeographic2DWithHeightOffsets(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Angle &offsetLat, const common::Angle &offsetLon, const common::Length &offsetHeight, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with method Geographic 2D with height offsets.

This method is defined as EPSG:9618

  • Return

    new Transformation.

    Parameters
    • properties: See general_properties of the Transformation. At minimum the name should be defined.

    • sourceCRSIn: Source CRS.

    • targetCRSIn: Target CRS.

    • offsetLat: Latitude offset to add.

    • offsetLon: Longitude offset to add.

    • offsetHeight: Geoid undulation to add.

    • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createVerticalOffset(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Length &offsetHeight, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation with method Vertical Offset.

This method is defined as EPSG:9616

  • Return

    new Transformation.

    Parameters
    • properties: See general_properties of the Transformation. At minimum the name should be defined.

    • sourceCRSIn: Source CRS.

    • targetCRSIn: Target CRS.

    • offsetHeight: Geoid undulation to add.

    • accuracies: Vector of positional accuracy (might be empty).

TransformationNNPtr createChangeVerticalUnit(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Scale &factor, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)

Instantiate a transformation based on the Change of Vertical Unit method.

This method is defined as EPSG:1069

Return

a new Transformation.

Parameters
  • properties: See general_properties of the conversion. If the name is not provided, it is automatically set.

  • sourceCRSIn: Source CRS.

  • targetCRSIn: Target CRS.

  • factor: Conversion factor

  • accuracies: Vector of positional accuracy (might be empty).