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 CoordinateTransformerPtr = std::unique_ptr<CoordinateTransformer>¶
Shared pointer of CoordinateTransformer
-
using CoordinateTransformerNNPtr = util::nn<CoordinateTransformerPtr>¶
Non-null shared pointer of CoordinateTransformer
-
using TransformationPtr = std::shared_ptr<Transformation>¶
Shared pointer of Transformation
-
using TransformationNNPtr = util::nn<TransformationPtr>¶
Non-null shared pointer of Transformation
-
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
-
typedef std::shared_ptr<PointMotionOperation> PointMotionOperationPtr¶
Shared pointer of PointMotionOperation
-
typedef util::nn<PointMotionOperationPtr> PointMotionOperationNNPtr¶
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 std::list<crs::GeodeticCRSNNPtr> findCandidateGeodCRSForDatum(const io::AuthorityFactoryPtr &authFactory, const crs::GeodeticCRS *crs, const datum::GeodeticReferenceFrameNNPtr &datum)¶
-
static bool isRegularVerticalGridMethod(int methodEPSGCode, bool &reverseOffsetSign)¶
-
static double negate(double val)¶
-
static CoordinateOperationPtr createApproximateInverseIfPossible(const Transformation *op)¶
-
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.
-
std::string shortName¶
-
class CoordinateTransformer¶
- #include <coordinateoperation.hpp>
Coordinate transformer.
Performs coordinate transformation of coordinate tuplies.
- Since
9.3
Public Functions
-
PJ_COORD transform(PJ_COORD coord)¶
Transforms a coordinate tuple.
PJ_COORD is a union of many structures. In the context of this method, it is prudent to only use the v[] array, with the understanding that the expected input values should be passed in the order and the unit of the successive axis of the input CRS. Similarly the values returned in the v[] array of the output PJ_COORD are in the order and the unit of the successive axis of the output CRS. For coordinate operations involving a time-dependent operation, coord.v[3] is the decimal year of the coordinate epoch of the input (or HUGE_VAL to indicate none)
If an error occurs, HUGE_VAL is returned in the .v[0] member of the output coordinate tuple.
Example how to transform coordinates from EPSG:4326 (WGS 84 latitude/longitude) to EPSG:32631 (WGS 84 / UTM zone 31N).
auto authFactory = AuthorityFactory::create(DatabaseContext::create(), std::string()); auto coord_op_ctxt = CoordinateOperationContext::create( authFactory, nullptr, 0.0); auto authFactoryEPSG = AuthorityFactory::create(DatabaseContext::create(), "EPSG"); auto list = CoordinateOperationFactory::create()->createOperations( authFactoryEPSG->createCoordinateReferenceSystem("4326"), authFactoryEPSG->createCoordinateReferenceSystem("32631"), coord_op_ctxt); ASSERT_TRUE(!list.empty()); PJ_CONTEXT* ctx = proj_context_create(); auto transformer = list[0]->coordinateTransformer(ctx); PJ_COORD c; c.v[0] = 49; // latitude in degree c.v[1] = 2; // longitude in degree c.v[2] = 0; c.v[3] = HUGE_VAL; c = transformer->transform(c); EXPECT_NEAR(c.v[0], 426857.98771728, 1e-8); // easting in metre EXPECT_NEAR(c.v[1], 5427937.52346492, 1e-8); // northing in metre proj_context_destroy(ctx);
-
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.
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.
- Returns:
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.
- Returns:
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.
- Returns:
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.
- Returns:
target CRS, or null.
-
const crs::CRSPtr &interpolationCRS() const¶
Return the interpolation CRS of this coordinate operation.
- Returns:
interpolation CRS, or null.
-
const util::optional<common::DataEpoch> &sourceCoordinateEpoch() const¶
Return the source epoch of coordinates.
- Returns:
source epoch of coordinates, or empty.
-
const util::optional<common::DataEpoch> &targetCoordinateEpoch() const¶
Return the target epoch of coordinates.
- Returns:
target epoch of coordinates, or empty.
-
CoordinateTransformerNNPtr coordinateTransformer(PJ_CONTEXT *ctx) const¶
Return a coordinate transformer for this operation.
The returned coordinate transformer is tied to the provided context, and should only be called by the thread "owning" the passed context. It should not be used after the context has been destroyed.
- Since
9.3
- Parameters:
ctx -- Execution context to which the transformer will be tied to. If null, the default context will be used (only sfe for single-threaded applications).
- Throws:
UnsupportedOperationException -- if the transformer cannot be instantiated.
- Returns:
a new CoordinateTransformer instance.
-
virtual CoordinateOperationNNPtr inverse() const = 0¶
Return the inverse of the coordinate operation.
- Throws:
-
virtual std::set<GridDescription> gridsNeeded(const io::DatabaseContextPtr &databaseContext, bool considerKnownGridsAsAvailable) const = 0¶
Return grids needed by an operation.
-
bool isPROJInstantiable(const io::DatabaseContextPtr &databaseContext, bool considerKnownGridsAsAvailable) 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
-
static 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 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 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.
- Returns:
code, or 0 if not found
Public Static Functions
-
static OperationParameterNNPtr create(const util::PropertyMap &properties)¶
Instantiate a OperationParameter.
- Parameters:
properties -- See General properties. At minimum the name should be defined.
- Returns:
a new OperationParameter.
-
static const char *getNameForEPSGCode(int epsg_code) noexcept¶
Return the name of a parameter designed by its EPSG code.
- Returns:
name, or nullptr if not found
-
int getEPSGCode()¶
-
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
-
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
Public Functions
-
const common::Measure &value()¶
Returns the value as a Measure (assumes type() == Type::MEASURE)
- Returns:
the value as a Measure.
-
const std::string &stringValue()¶
Returns the value as a string (assumes type() == Type::STRING)
- Returns:
the value as a string.
-
const std::string &valueFile()¶
Returns the value as a filename (assumes type() == Type::FILENAME)
- Returns:
the value as a filename.
-
int integerValue()¶
Returns the value as a integer (assumes type() == Type::INTEGER)
- Returns:
the value as a integer.
-
bool booleanValue()¶
Returns the value as a boolean (assumes type() == Type::BOOLEAN)
- Returns:
the value as a boolean.
Public Static Functions
-
static ParameterValueNNPtr create(const common::Measure &measureIn)¶
Instantiate a ParameterValue from a Measure (i.e. a value associated with a unit)
- Returns:
a new ParameterValue.
-
static ParameterValueNNPtr create(const char *stringValueIn)¶
Instantiate a ParameterValue from a string value.
- Returns:
a new ParameterValue.
-
static ParameterValueNNPtr create(const std::string &stringValueIn)¶
Instantiate a ParameterValue from a string value.
- Returns:
a new ParameterValue.
-
static ParameterValueNNPtr create(int integerValueIn)¶
Instantiate a ParameterValue from a integer value.
- Returns:
a new ParameterValue.
-
static ParameterValueNNPtr create(bool booleanValueIn)¶
Instantiate a ParameterValue from a boolean value.
- Returns:
a new ParameterValue.
-
static ParameterValueNNPtr createFilename(const std::string &stringValueIn)¶
Instantiate a ParameterValue from a filename.
- Returns:
a new ParameterValue.
-
const common::Measure &value()¶
-
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 ¶meter()¶
Return the parameter (definition)
- Returns:
the parameter (definition).
-
const ParameterValueNNPtr ¶meterValue()¶
Return the parameter value.
- Returns:
the parameter value.
Public Static Functions
-
static OperationParameterValueNNPtr create(const OperationParameterNNPtr ¶meterIn, const ParameterValueNNPtr &valueIn)¶
Instantiate a OperationParameterValue.
- Parameters:
parameterIn -- Parameter (definition).
valueIn -- Parameter value.
- Returns:
a new OperationParameterValue.
-
const OperationParameterNNPtr ¶meter()¶
-
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.
- Returns:
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.
- Returns:
the formula citation, or empty.
-
const std::vector<GeneralOperationParameterNNPtr> ¶meters()¶
Return the parameters of this operation method.
- Returns:
the parameters.
-
int getEPSGCode()¶
Return the EPSG code, either directly, or through the name.
- Returns:
code, or 0 if not found
Public Static Functions
-
static OperationMethodNNPtr create(const util::PropertyMap &properties, const std::vector<GeneralOperationParameterNNPtr> ¶meters)¶
Instantiate a operation method from a vector of GeneralOperationParameter.
- Parameters:
properties -- See General properties. At minimum the name should be defined.
parameters -- Vector of GeneralOperationParameterNNPtr.
- Returns:
a new OperationMethod.
-
static OperationMethodNNPtr create(const util::PropertyMap &properties, const std::vector<OperationParameterNNPtr> ¶meters)¶
Instantiate a operation method from a vector of OperationParameter.
- Parameters:
properties -- See General properties. At minimum the name should be defined.
parameters -- Vector of OperationParameterNNPtr.
- Returns:
a new OperationMethod.
-
const util::optional<std::string> &formula()¶
-
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 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> ¶meterValues()¶
Return the parameter values.
- Returns:
the parameter values.
-
const OperationMethodNNPtr &method()¶
Return the operation method associated to the operation.
- Returns:
the operation method.
-
const ParameterValuePtr ¶meterValue(const std::string ¶mName, int epsg_code = 0) const noexcept¶
Return the parameter value corresponding to a parameter name or EPSG code.
- Parameters:
paramName -- the parameter name (or empty, in which case epsg_code should be non zero)
epsg_code -- the parameter EPSG code (possibly zero)
- Returns:
the value, or nullptr if not found.
-
const ParameterValuePtr ¶meterValue(int epsg_code) const noexcept¶
Return the parameter value corresponding to a EPSG code.
- Parameters:
epsg_code -- the parameter EPSG code
- Returns:
the value, or nullptr if not found.
-
const common::Measure ¶meterValueMeasure(const std::string ¶mName, int epsg_code = 0) const noexcept¶
Return the parameter value, as a measure, corresponding to a parameter name or EPSG code.
- Parameters:
paramName -- the parameter name (or empty, in which case epsg_code should be non zero)
epsg_code -- the parameter EPSG code (possibly zero)
- Returns:
the measure, or the empty Measure() object if not found.
-
const common::Measure ¶meterValueMeasure(int epsg_code) const noexcept¶
Return the parameter value, as a measure, corresponding to a EPSG code.
- Parameters:
epsg_code -- the parameter EPSG code
- Returns:
the measure, or the empty Measure() object if not found.
-
virtual std::set<GridDescription> gridsNeeded(const io::DatabaseContextPtr &databaseContext, bool considerKnownGridsAsAvailable) const override¶
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.
-
TransformationNNPtr substitutePROJAlternativeGridNames(io::DatabaseContextNNPtr databaseContext) const¶
Return an equivalent transformation to the current one, but using PROJ alternative grid names.
Public Static Functions
-
static 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.
- 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).
- Returns:
the new instance
-
const std::vector<GeneralParameterValueNNPtr> ¶meterValues()¶
-
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
Projection parameters¶
Co-latitude of cone axis¶
The rotation applied to spherical coordinates for the oblique projection, measured on the conformal sphere in the plane of the meridian of origin.
EPSG:1036
Latitude of natural origin/Center Latitude¶
The latitude of the point from which the values of both the geographical coordinates on the ellipsoid and the grid coordinates on the projection are deemed to increment or decrement for computational purposes. Alternatively it may be considered as the latitude of the point which in the absence of application of false coordinates has grid coordinates of (0,0).
EPSG:8801
Longitude of natural origin/Central Meridian¶
The longitude of the point from which the values of both the geographical coordinates on the ellipsoid and the grid coordinates on the projection are deemed to increment or decrement for computational purposes. Alternatively it may be considered as the longitude of the point which in the absence of application of false coordinates has grid coordinates of (0,0). Sometimes known as "central meridian (CM)".
EPSG:8802
Scale Factor¶
The factor by which the map grid is reduced or enlarged during the projection process, defined by its value at the natural origin.
EPSG:8805
False Easting¶
Since the natural origin may be at or near the centre of the projection and under normal coordinate circumstances would thus give rise to negative coordinates over parts of the mapped area, this origin is usually given false coordinates which are large enough to avoid this inconvenience. The False Easting, FE, is the value assigned to the abscissa (east or west) axis of the projection grid at the natural origin.
EPSG:8806
False Northing¶
Since the natural origin may be at or near the centre of the projection and under normal coordinate circumstances would thus give rise to negative coordinates over parts of the mapped area, this origin is usually given false coordinates which are large enough to avoid this inconvenience. The False Northing, FN, is the value assigned to the ordinate (north or south) axis of the projection grid at the natural origin.
EPSG:8807
Latitude of projection centre¶
For an oblique projection, this is the latitude of the point at which the azimuth of the central line is defined.
EPSG:8811
Longitude of projection centre¶
For an oblique projection, this is the longitude of the point at which the azimuth of the central line is defined.
EPSG:8812
Azimuth of initial line¶
The azimuthal direction (north zero, east of north being positive) of the great circle which is the centre line of an oblique projection. The azimuth is given at the projection centre.
EPSG:8813
Angle from Rectified to Skew Grid¶
The angle at the natural origin of an oblique projection through which the natural coordinate reference system is rotated to make the projection north axis parallel with true north.
EPSG:8814
Scale factor on initial line¶
The factor by which the map grid is reduced or enlarged during the projection process, defined by its value at the projection center.
EPSG:8815
Easting at projection centre¶
The easting value assigned to the projection centre.
EPSG:8816
Northing at projection centre¶
The northing value assigned to the projection centre.
EPSG:8817
Latitude of pseudo standard¶
parallel
Latitude of the parallel on which the conic or cylindrical projection is based. This latitude is not geographic, but is defined on the conformal sphere AFTER its rotation to obtain the oblique aspect of the projection.
EPSG:8818
Scale factor on pseudo¶
standard parallel
The factor by which the map grid is reduced or enlarged during the projection process, defined by its value at the pseudo-standard parallel. EPSG:8819
Latitude of false origin¶
The latitude of the point which is not the natural origin and at which grid coordinate values false easting and false northing are defined.
EPSG:8821
Longitude of false origin¶
The longitude of the point which is not the natural origin and at which grid coordinate values false easting and false northing are defined.
EPSG:8822
Latitude of 1st standard parallel¶
For a conic projection with two standard parallels, this is the latitude of one of the parallels of intersection of the cone with the ellipsoid. It is normally but not necessarily that nearest to the pole. Scale is true along this parallel.
EPSG:8823
Latitude of 2nd standard parallel¶
For a conic projection with two standard parallels, this is the latitude of one of the parallels at which the cone intersects with the ellipsoid. It is normally but not necessarily that nearest to the equator. Scale is true along this parallel.
EPSG:8824
Easting of false origin¶
The easting value assigned to the false origin.
EPSG:8826
Northing of false origin¶
The northing value assigned to the false origin.
EPSG:8827
Latitude of standard parallel¶
For polar aspect azimuthal projections, the parallel on which the scale factor is defined to be unity.
EPSG:8832
Longitude of origin¶
For polar aspect azimuthal projections, the meridian along which the northing axis increments and also across which parallels of latitude increment towards the north pole.
EPSG:8833
Public Functions
-
virtual CoordinateOperationNNPtr inverse() const override¶
Return the inverse of the coordinate operation.
- Throws:
-
bool isUTM(int &zone, bool &north) const¶
Return whether a conversion is a Universal Transverse Mercator conversion.
- Parameters:
zone -- [out] UTM zone number between 1 and 60.
north -- [out] true for UTM northern hemisphere, false for UTM southern hemisphere.
- Returns:
true if it is a UTM conversion.
-
ConversionNNPtr identify() const¶
Return a Conversion object where some parameters are better identified.
- Returns:
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
- Parameters:
targetEPSGCode -- EPSG code of the target method.
- Returns:
new conversion, or nullptr
Public Static Functions
-
static ConversionNNPtr create(const util::PropertyMap &properties, const OperationMethodNNPtr &methodIn, const std::vector<GeneralParameterValueNNPtr> &values)¶
Instantiate a Conversion from a vector of GeneralParameterValue.
- Parameters:
properties -- See General properties. At minimum the name should be defined.
methodIn -- the operation method.
values -- the values.
- Throws:
- Returns:
a new Conversion.
-
static ConversionNNPtr create(const util::PropertyMap &propertiesConversion, const util::PropertyMap &propertiesOperationMethod, const std::vector<OperationParameterNNPtr> ¶meters, const std::vector<ParameterValueNNPtr> &values)¶
Instantiate a Conversion and its OperationMethod.
- 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()
- Throws:
- Returns:
a new Conversion.
-
static 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.
- 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.
- Returns:
a new Conversion.
-
static ConversionNNPtr createTransverseMercator(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createGaussSchreiberTransverseMercator(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createTransverseMercatorSouthOriented(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createTunisiaMappingGrid(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- Deprecated:
. Use createTunisiaMiningGrid() instead
Note
There is currently no implementation of the method formulas in PROJ.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createTunisiaMiningGrid(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Tunisia Mining Grid projection method.
This method is defined as EPSG:9816.
- Since
9.2
Note
There is currently no implementation of the method formulas in PROJ.
- 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
- Returns:
a new Conversion.
-
static 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
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createLambertConicConformal_1SP(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createLambertConicConformal_1SP_VariantB(const util::PropertyMap &properties, const common::Angle &latitudeNatOrigin, const common::Scale &scale, const common::Angle &latitudeFalseOrigin, const common::Angle &longitudeFalseOrigin, const common::Length &eastingFalseOrigin, const common::Length &northingFalseOrigin)¶
Instantiate a conversion based on the Lambert Conic Conformal 1SP Variant B projection method.
This method is defined as EPSG:1102.
- Since
9.2.1
- 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
scale -- See Scale Factor
latitudeFalseOrigin -- See Latitude of false origin
longitudeFalseOrigin -- See Longitude of false origin
eastingFalseOrigin -- See Easting of false origin
northingFalseOrigin -- See Northing of false origin
- Returns:
a new Conversion.
-
static 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
- 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
- Returns:
a new Conversion.
-
static 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.
- 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.
- Returns:
a new Conversion.
-
static 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.
Note
the order of arguments is conformant with the corresponding EPSG mode and different than OGRSpatialReference::setLCCB() of GDAL <= 2.3
Warning
The formulas used currently in PROJ are, incorrectly, the ones of the regular LCC_2SP method.
- 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
- Returns:
a new Conversion.
-
static 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 Azimuthal Equidistant projection method.
This method is defined as EPSG:1125.
- 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
- Returns:
a new Conversion.
-
static 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 method.
This method is defined as EPSG:9831.
- 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
- Returns:
a new Conversion.
-
static 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.
- 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
- Returns:
a new Conversion.
-
static 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 GeographicCRS.
- 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
- Returns:
a new Conversion.
-
static 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createCassiniSoldner(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createEquidistantConic(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 Equidistant Conic projection method.
This method is defined as EPSG:1119.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createEckertI(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Eckert I projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createEckertII(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Eckert II projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createEckertIII(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Eckert III projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createEckertIV(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Eckert IV projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createEckertV(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Eckert V projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createEckertVI(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Eckert VI projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static 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.
- 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
- Returns:
a new Conversion.
-
static 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createGall(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Gall (Stereographic) projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createGoodeHomolosine(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Goode Homolosine projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createInterruptedGoodeHomolosine(const util::PropertyMap &properties, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createGeostationarySatelliteSweepX(const util::PropertyMap &properties, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createGeostationarySatelliteSweepY(const util::PropertyMap &properties, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createGnomonic(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Gnomonic projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static 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.
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
latitudeProjectionCentre -- See Latitude 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
- Returns:
a new Conversion.
-
static 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.
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
latitudeProjectionCentre -- See Latitude 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
- Returns:
a new Conversion.
-
static 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.
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
latitudeProjectionCentre -- See Latitude 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
- Returns:
a new Conversion.
-
static 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.
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
latitudeProjectionCentre -- See Latitude 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createInternationalMapWorldPolyconic(const util::PropertyMap &properties, const common::Angle ¢erLong, 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
- 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
- Returns:
a new Conversion.
-
static 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.
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
latitudeProjectionCentre -- See Latitude 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
- Returns:
a new Conversion.
-
static 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.
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
latitudeProjectionCentre -- See Latitude 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
- Returns:
a new Conversion.
-
static 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createMillerCylindrical(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Miller Cylindrical projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createMercatorVariantA(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, const common::Scale &scale, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Mercator (variant A) projection method.
This is the A 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createMercatorVariantB(const util::PropertyMap &properties, const common::Angle &latitudeFirstParallel, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Mercator (variant B) projection method.
This is the B 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createPopularVisualisationPseudoMercator(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createMercatorSpherical(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Mercator projection method, using its spherical formulation.
When used with an ellipsoid, the radius used is the radius of the conformal sphere at centerLat.
This method is defined as EPSG:1026.
- Since
9.3
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createMollweide(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Mollweide projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createNewZealandMappingGrid(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createObliqueStereographic(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createOrthographic(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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
Before PROJ 7.2, only the spherical formulation was implemented.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createAmericanPolyconic(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createPolarStereographicVariantA(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createRobinson(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Robinson projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createSinusoidal(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Sinusoidal projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createStereographic(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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().- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createVanDerGrinten(const util::PropertyMap &properties, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createWagnerI(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Wagner I projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createWagnerII(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Wagner II projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createWagnerIII(const util::PropertyMap &properties, const common::Angle &latitudeTrueScale, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Wagner III projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createWagnerIV(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Wagner IV projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createWagnerV(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Wagner V projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createWagnerVI(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Wagner VI projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createWagnerVII(const util::PropertyMap &properties, const common::Angle ¢erLong, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Wagner VII projection method.
There is no equivalent in EPSG.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createQuadrilateralizedSphericalCube(const util::PropertyMap &properties, const common::Angle ¢erLat, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static 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.
- 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.
- Returns:
a new Conversion.
-
static ConversionNNPtr createEqualEarth(const util::PropertyMap &properties, const common::Angle ¢erLong, 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.
- 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
- Returns:
a new Conversion.
-
static ConversionNNPtr createVerticalPerspective(const util::PropertyMap &properties, const common::Angle &topoOriginLat, const common::Angle &topoOriginLong, const common::Length &topoOriginHeight, const common::Length &viewPointHeight, const common::Length &falseEasting, const common::Length &falseNorthing)¶
Instantiate a conversion based on the Vertical Perspective projection method.
This method is defined as EPSG:9838.
The PROJ implementation of the EPSG Vertical Perspective has the current limitations with respect to the method described in EPSG:
it is a 2D-only method, ignoring the ellipsoidal height of the point to project.
it has only a spherical development.
the height of the topocentric origin is ignored, and thus assumed to be 0.
For completeness, PROJ adds the falseEasting and falseNorthing parameter, which are not described in EPSG. They should usually be set to 0.
- Since
6.3
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
topoOriginLat -- Latitude of topocentric origin
topoOriginLong -- Longitude of topocentric origin
topoOriginHeight -- Ellipsoidal height of topocentric origin. Ignored by PROJ (that is assumed to be 0)
viewPointHeight -- Viewpoint height with respect to the topocentric/mapping plane. In the case where topoOriginHeight = 0, this is the height above the ellipsoid surface at topoOriginLat, topoOriginLong.
falseEasting -- See False Easting . (not in EPSG)
falseNorthing -- See False Northing . (not in EPSG)
- Returns:
a new Conversion.
-
static ConversionNNPtr createPoleRotationGRIBConvention(const util::PropertyMap &properties, const common::Angle &southPoleLatInUnrotatedCRS, const common::Angle &southPoleLongInUnrotatedCRS, const common::Angle &axisRotation)¶
Instantiate a conversion based on the Pole Rotation method, using the conventions of the GRIB 1 and GRIB 2 data formats.
Those are mentioned in the Note 2 of https://www.nco.ncep.noaa.gov/pmb/docs/grib2/grib2_doc/grib2_temp3-1.shtml
Several conventions for the pole rotation method exists. The parameters provided in this method are remapped to the PROJ ob_tran operation with:
Another implementation of that convention is also in the netcdf-java library: https://github.com/Unidata/netcdf-java/blob/3ce72c0cd167609ed8c69152bb4a004d1daa9273/cdm/core/src/main/java/ucar/unidata/geoloc/projection/RotatedLatLon.java
The PROJ implementation of this method assumes a spherical ellipsoid.
- Since
7.0
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
southPoleLatInUnrotatedCRS -- Latitude of the point from the unrotated CRS, expressed in the unrotated CRS, that will become the south pole of the rotated CRS.
southPoleLongInUnrotatedCRS -- Longitude of the point from the unrotated CRS, expressed in the unrotated CRS, that will become the south pole of the rotated CRS.
axisRotation -- The angle of rotation about the new polar axis (measured clockwise when looking from the southern to the northern pole) of the coordinate system, assuming the new axis to have been obtained by first rotating the sphere through southPoleLongInUnrotatedCRS degrees about the geographic polar axis and then rotating through (90 + southPoleLatInUnrotatedCRS) degrees so that the southern pole moved along the (previously rotated) Greenwich meridian.
- Returns:
a new Conversion.
-
static ConversionNNPtr createPoleRotationNetCDFCFConvention(const util::PropertyMap &properties, const common::Angle &gridNorthPoleLatitude, const common::Angle &gridNorthPoleLongitude, const common::Angle &northPoleGridLongitude)¶
Instantiate a conversion based on the Pole Rotation method, using the conventions of the netCDF CF convention for the netCDF format.
Those are mentioned in the Note 2 of https://cfconventions.org/Data/cf-conventions/cf-conventions-1.8/cf-conventions.html#_rotated_pole
Several conventions for the pole rotation method exists. The parameters provided in this method are remapped to the PROJ ob_tran operation with:
Another implementation of that convention is also in the netcdf-java library: https://github.com/Unidata/netcdf-java/blob/3ce72c0cd167609ed8c69152bb4a004d1daa9273/cdm/core/src/main/java/ucar/unidata/geoloc/projection/RotatedPole.java
The PROJ implementation of this method assumes a spherical ellipsoid.
- Since
8.2
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
gridNorthPoleLatitude -- True latitude of the north pole of the rotated grid
gridNorthPoleLongitude -- True longitude of the north pole of the rotated grid.
northPoleGridLongitude -- Longitude of the true north pole in the rotated grid.
- Returns:
a new Conversion.
-
static 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 [DEPRECATED].
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
factor -- Conversion factor
- Returns:
a new Conversion.
-
static ConversionNNPtr createChangeVerticalUnit(const util::PropertyMap &properties)¶
Instantiate a conversion based on the Change of Vertical Unit method (without explicit conversion factor)
This method is defined as EPSG:1104.
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
- Returns:
a new Conversion.
-
static ConversionNNPtr createHeightDepthReversal(const util::PropertyMap &properties)¶
Instantiate a conversion based on the Height Depth Reversal method.
This method is defined as EPSG:1068.
- Since
6.3
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
- Returns:
a new Conversion.
-
static 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 for 2D or EPSG:9844 for Geographic3D horizontal.
- Parameters:
is3D -- Whether this should apply on 3D geographicCRS
- Returns:
a new Conversion.
-
static ConversionNNPtr createGeographicGeocentric(const util::PropertyMap &properties)¶
Instantiate a conversion based on the Geographic/Geocentric method.
This method is defined as EPSG:9602.
- Parameters:
properties -- See General properties of the conversion. If the name is not provided, it is automatically set.
- Returns:
a new Conversion.
-
static ConversionNNPtr createGeographic2DOffsets(const util::PropertyMap &properties, const common::Angle &offsetLat, const common::Angle &offsetLong)¶
Instantiate a conversion with method Geographic 2D offsets.
This method is defined as EPSG:9619.
- Parameters:
properties -- See General properties of the conversion. At minimum the name should be defined.
offsetLat -- Latitude offset to add.
offsetLong -- Longitude offset to add.
- Returns:
new conversion.
-
static ConversionNNPtr createGeographic3DOffsets(const util::PropertyMap &properties, const common::Angle &offsetLat, const common::Angle &offsetLong, const common::Length &offsetHeight)¶
Instantiate a conversion with method Geographic 3D offsets.
This method is defined as EPSG:9660.
- Parameters:
properties -- See General properties of the Conversion. At minimum the name should be defined.
offsetLat -- Latitude offset to add.
offsetLong -- Longitude offset to add.
offsetHeight -- Height offset to add.
- Returns:
new Conversion.
-
static ConversionNNPtr createGeographic2DWithHeightOffsets(const util::PropertyMap &properties, const common::Angle &offsetLat, const common::Angle &offsetLong, const common::Length &offsetHeight)¶
Instantiate a conversion with method Geographic 2D with height offsets.
This method is defined as EPSG:9618.
- Parameters:
properties -- See General properties of the Conversion. At minimum the name should be defined.
offsetLat -- Latitude offset to add.
offsetLong -- Longitude offset to add.
offsetHeight -- Geoid undulation to add.
- Returns:
new Conversion.
-
static ConversionNNPtr createVerticalOffset(const util::PropertyMap &properties, const common::Length &offsetHeight)¶
Instantiate a conversion with method Vertical Offset.
This method is defined as EPSG:9616.
- Parameters:
properties -- See General properties of the Conversion. At minimum the name should be defined.
offsetHeight -- Geoid undulation to add.
- Returns:
new Conversion.
-
virtual CoordinateOperationNNPtr inverse() const override¶
-
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.
- Returns:
the source CRS.
-
const crs::CRSNNPtr &targetCRS()¶
Return the target crs::CRS of the transformation.
- Returns:
the target CRS.
-
virtual CoordinateOperationNNPtr inverse() const override¶
Return the inverse of the coordinate operation.
- Throws:
Public Static Functions
-
static 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.
- 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).
- Throws:
- Returns:
new Transformation.
-
static 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> ¶meters, const std::vector<ParameterValueNNPtr> &values, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)¶
Instantiate a transformation and its OperationMethod.
- 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).
- Throws:
- Returns:
new Transformation.
-
static 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.
- 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).
- Returns:
new Transformation.
-
static 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.
- 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).
- Returns:
new Transformation.
-
static 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.
- 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).
- Returns:
new Transformation.
-
static 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.
- 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).
- Returns:
new Transformation.
-
static 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.
- 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).
- Throws:
- Returns:
new Transformation.
-
static 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
- 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()
- Throws:
- Returns:
new Transformation.
-
static 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.
- 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).
- Returns:
new Transformation.
-
static 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 also
createAbridgedMolodensky() for a related method.
- 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).
- Throws:
- Returns:
new Transformation.
-
static 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 also
createdMolodensky() for a related method.
- 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).
- Throws:
- Returns:
new Transformation.
-
static 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.
- 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).
- Returns:
new Transformation.
-
static 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.
- 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).
- Returns:
new Transformation.
-
static 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.
- Parameters:
properties -- See General properties of the Transformation. At minimum the name should be defined.
sourceCRSIn -- Source CRS.
targetCRSIn -- Target CRS.
offset -- Longitude offset to add.
- Returns:
new Transformation.
-
static TransformationNNPtr createGeographic2DOffsets(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Angle &offsetLat, const common::Angle &offsetLong, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)¶
Instantiate a transformation with method Geographic 2D offsets.
This method is defined as EPSG:9619.
- 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.
offsetLong -- Longitude offset to add.
accuracies -- Vector of positional accuracy (might be empty).
- Returns:
new Transformation.
-
static TransformationNNPtr createGeographic3DOffsets(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Angle &offsetLat, const common::Angle &offsetLong, 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.
- 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.
offsetLong -- Longitude offset to add.
offsetHeight -- Height offset to add.
accuracies -- Vector of positional accuracy (might be empty).
- Returns:
new Transformation.
-
static TransformationNNPtr createGeographic2DWithHeightOffsets(const util::PropertyMap &properties, const crs::CRSNNPtr &sourceCRSIn, const crs::CRSNNPtr &targetCRSIn, const common::Angle &offsetLat, const common::Angle &offsetLong, 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.
- 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.
offsetLong -- Longitude offset to add.
offsetHeight -- Geoid undulation to add.
accuracies -- Vector of positional accuracy (might be empty).
- Returns:
new Transformation.
-
static 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.
- 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).
- Returns:
new Transformation.
-
static 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 [DEPRECATED].
- 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).
- Returns:
a new Transformation.
-
const crs::CRSNNPtr &sourceCRS()¶
-
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
Public Functions
-
const crs::CRSNNPtr &sourceCRS()¶
Return the source crs::CRS of the operation.
- Returns:
the source CRS.
-
virtual CoordinateOperationNNPtr inverse() const override¶
Return the inverse of the coordinate operation.
- Throws:
-
PointMotionOperationNNPtr substitutePROJAlternativeGridNames(io::DatabaseContextNNPtr databaseContext) const¶
Return an equivalent transformation to the current one, but using PROJ alternative grid names.
Public Static Functions
-
static PointMotionOperationNNPtr create(const util::PropertyMap &properties, const crs::CRSNNPtr &crsIn, const OperationMethodNNPtr &methodIn, const std::vector<GeneralParameterValueNNPtr> &values, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)¶
Instantiate a point motion operation from a vector of GeneralParameterValue.
- Parameters:
properties -- See General properties. At minimum the name should be defined.
crsIn -- Source and target CRS.
methodIn -- Operation method.
values -- Vector of GeneralOperationParameterNNPtr.
accuracies -- Vector of positional accuracy (might be empty).
- Throws:
- Returns:
new PointMotionOperation.
-
static PointMotionOperationNNPtr create(const util::PropertyMap &propertiesOperation, const crs::CRSNNPtr &crsIn, const util::PropertyMap &propertiesOperationMethod, const std::vector<OperationParameterNNPtr> ¶meters, const std::vector<ParameterValueNNPtr> &values, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)¶
Instantiate a point motion operation and its OperationMethod.
- Parameters:
propertiesOperation -- The General properties of the PointMotionOperation. At minimum the name should be defined.
crsIn -- Source and target CRS.
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).
- Throws:
- Returns:
new PointMotionOperation.
-
const crs::CRSNNPtr &sourceCRS()¶
-
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.
- Returns:
the operation steps.
-
virtual CoordinateOperationNNPtr inverse() const override¶
Return the inverse of the coordinate operation.
- Throws:
-
virtual std::set<GridDescription> gridsNeeded(const io::DatabaseContextPtr &databaseContext, bool considerKnownGridsAsAvailable) const override¶
Return grids needed by an operation.
Public Static Functions
-
static ConcatenatedOperationNNPtr create(const util::PropertyMap &properties, const std::vector<CoordinateOperationNNPtr> &operationsIn, const std::vector<metadata::PositionalAccuracyNNPtr> &accuracies)¶
Instantiate a ConcatenatedOperation.
- 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).
- Throws:
- Returns:
new Transformation.
-
static 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-emptiness of the intersection of the extents of the operations
- Throws:
-
const std::vector<CoordinateOperationNNPtr> &operations() const¶
-
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 class 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:
-
enumerator NONE¶
Ignore CRS extent
-
enumerator BOTH¶
Test coordinate operation extent against both CRS extent.
-
enumerator INTERSECTION¶
Test coordinate operation extent against the intersection of both CRS extent.
-
enumerator SMALLEST¶
Test coordinate operation against the smallest of both CRS extent.
-
enumerator NONE¶
-
enum class SpatialCriterion¶
Spatial criterion to restrict candidate operations.
Values:
-
enumerator STRICT_CONTAINMENT¶
The area of validity of transforms should strictly contain the are of interest.
-
enumerator PARTIAL_INTERSECTION¶
The area of validity of transforms should at least intersect the area of interest.
-
enumerator STRICT_CONTAINMENT¶
-
enum class GridAvailabilityUse¶
Describe how grid availability is used.
Values:
-
enumerator USE_FOR_SORTING¶
Grid availability is only used for sorting results. Operations where some grids are missing will be sorted last.
-
enumerator DISCARD_OPERATION_IF_MISSING_GRID¶
Completely discard an operation if a required grid is missing.
-
enumerator IGNORE_GRID_AVAILABILITY¶
Ignore grid availability at all. Results will be presented as if all grids were available.
-
enumerator KNOWN_AVAILABLE¶
Results will be presented as if grids known to PROJ (that is registered in the grid_alternatives table of its database) were available. Used typically when networking is enabled.
-
enumerator USE_FOR_SORTING¶
Public Functions
-
const io::AuthorityFactoryPtr &getAuthorityFactory() const¶
Return the authority factory, 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 setAllowBallparkTransformations(bool allow)¶
Set whether ballpark transformations are allowed.
-
bool getAllowBallparkTransformations() const¶
Return whether ballpark transformations are allowed.
-
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.
-
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.
-
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.
-
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 strategy, all potential C candidates will be used if there is no direct transformation.
-
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.
-
void setSourceCoordinateEpoch(const util::optional<common::DataEpoch> &epoch)¶
Set the source coordinate epoch.
-
const util::optional<common::DataEpoch> &getSourceCoordinateEpoch() const¶
Return the source coordinate epoch.
-
void setTargetCoordinateEpoch(const util::optional<common::DataEpoch> &epoch)¶
Set the target coordinate epoch.
-
const util::optional<common::DataEpoch> &getTargetCoordinateEpoch() const¶
Return the target coordinate epoch.
-
CoordinateOperationContextNNPtr clone() const¶
Clone a coordinate operation context.
- Since
9.2
- Returns:
a new context.
Public Static Functions
-
static 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 operations will be searched only in that authority namespace.
- 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.
- Returns:
a new context.
-
enum class SourceTargetCRSExtentUse¶
-
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.
- Parameters:
sourceCRS -- source CRS.
targetCRS -- source CRS.
- Returns:
a CoordinateOperation or nullptr.
-
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.
When one of the source or target CRS has a vertical component but not the other one, the one that has no vertical component is automatically promoted to a 3D version, where its vertical axis is the ellipsoidal height in metres, using the ellipsoid of the base geodetic CRS.
- Parameters:
sourceCRS -- source CRS.
targetCRS -- target CRS.
context -- Search context.
- Returns:
a list
-
std::vector<CoordinateOperationNNPtr> createOperations(const coordinates::CoordinateMetadataNNPtr &sourceCoordinateMetadata, const crs::CRSNNPtr &targetCRS, const CoordinateOperationContextNNPtr &context) const¶
Find a list of CoordinateOperation from a source coordinate metadata to targetCRS.
- Since
9.2
- Parameters:
sourceCoordinateMetadata -- source CoordinateMetadata.
targetCRS -- target CRS.
context -- Search context.
- Returns:
a list
-
std::vector<CoordinateOperationNNPtr> createOperations(const crs::CRSNNPtr &sourceCRS, const coordinates::CoordinateMetadataNNPtr &targetCoordinateMetadata, const CoordinateOperationContextNNPtr &context) const¶
Find a list of CoordinateOperation from a source CRS to a target coordinate metadata.
- Since
9.2
- Parameters:
sourceCRS -- source CRS.
targetCoordinateMetadata -- target CoordinateMetadata.
context -- Search context.
- Returns:
a list
-
std::vector<CoordinateOperationNNPtr> createOperations(const coordinates::CoordinateMetadataNNPtr &sourceCoordinateMetadata, const coordinates::CoordinateMetadataNNPtr &targetCoordinateMetadata, const CoordinateOperationContextNNPtr &context) const¶
Find a list of CoordinateOperation from a source coordinate metadata to a target coordinate metadata.
Both source_crs and target_crs can be a CoordinateMetadata with an associated coordinate epoch, to perform changes of coordinate epochs. Note however than this is in practice limited to use of velocity grids inside the same dynamic CRS.
- Since
9.4
- Parameters:
sourceCoordinateMetadata -- source CoordinateMetadata.
targetCoordinateMetadata -- target CoordinateMetadata.
context -- Search context.
- Returns:
a list
Public Static Functions
-
static CoordinateOperationFactoryNNPtr create()¶
Instantiate a CoordinateOperationFactory.
-
CoordinateOperationPtr createOperation(const crs::CRSNNPtr &sourceCRS, const crs::CRSNNPtr &targetCRS) const¶
-
typedef std::shared_ptr<CoordinateOperation> CoordinateOperationPtr¶