GDAL
|
This class represents an OpenGIS Spatial Reference System, and contains methods for converting between this object organization and well known text (WKT) format. More...
#include <ogr_spatialref.h>
Public Member Functions | |
OGRSpatialReference (const char *=nullptr) | |
Constructor. | |
OGRSpatialReference (const OGRSpatialReference &) | |
Copy constructor. | |
OGRSpatialReference (OGRSpatialReference &&) | |
Move constructor. | |
virtual | ~OGRSpatialReference () |
OGRSpatialReference destructor. | |
OGRSpatialReference & | operator= (const OGRSpatialReference &) |
Assignment operator. | |
OGRSpatialReference & | operator= (OGRSpatialReference &&) |
Move assignment operator. | |
OGRSpatialReference & | AssignAndSetThreadSafe (const OGRSpatialReference &) |
Assignment method, with thread-safety. | |
int | Reference () |
Increments the reference count by one. | |
int | Dereference () |
Decrements the reference count by one. | |
int | GetReferenceCount () const |
Fetch current reference count. | |
void | Release () |
Decrements the reference count by one, and destroy if zero. | |
const char * | GetName () const |
Return the CRS name. | |
OGRSpatialReference * | Clone () const |
Make a duplicate of this OGRSpatialReference. | |
OGRSpatialReference * | CloneGeogCS () const |
Make a duplicate of the GEOGCS node of this OGRSpatialReference object. | |
void | dumpReadable () |
Dump pretty wkt to stdout, mostly for debugging. | |
OGRErr | exportToWkt (char **) const |
Convert this SRS into WKT 1 format. | |
OGRErr | exportToWkt (char **ppszWKT, const char *const *papszOptions) const |
Convert this SRS into a WKT string. | |
std::string | exportToWkt (const char *const *papszOptions=nullptr) const |
Convert this SRS into a WKT string. | |
OGRErr | exportToPrettyWkt (char **, int=FALSE) const |
Convert this SRS into a nicely formatted WKT 1 string for display to a person. | |
OGRErr | exportToPROJJSON (char **, const char *const *papszOptions) const |
Convert this SRS into a PROJJSON string. | |
OGRErr | exportToProj4 (char **) const |
Export coordinate system in PROJ.4 legacy format. | |
OGRErr | exportToPCI (char **, char **, double **) const |
Export coordinate system in PCI projection definition. | |
OGRErr | exportToUSGS (long *, long *, double **, long *) const |
Export coordinate system in USGS GCTP projection definition. | |
OGRErr | exportToXML (char **, const char *=nullptr) const |
Export coordinate system in XML format. | |
OGRErr | exportToPanorama (long *, long *, long *, long *, double *) const |
Export coordinate system in "Panorama" GIS projection definition. | |
OGRErr | exportVertCSToPanorama (int *) const |
Export vertical coordinate system to "Panorama" GIS projection definition. | |
OGRErr | exportToERM (char *pszProj, char *pszDatum, char *pszUnits) |
Convert coordinate system to ERMapper format. | |
OGRErr | exportToMICoordSys (char **) const |
Export coordinate system in Mapinfo style CoordSys format. | |
OGRErr | exportToCF1 (char **ppszGridMappingName, char ***ppapszKeyValues, char **ppszUnits, CSLConstList papszOptions) const |
Export a CRS to netCDF CF-1 definitions. | |
OGRErr | importFromWkt (char **) |
Import from WKT string. | |
OGRErr | importFromWkt (const char **) |
Import from WKT string. | |
OGRErr | importFromWkt (const char *) |
Import from WKT string. | |
OGRErr | importFromProj4 (const char *) |
Import PROJ coordinate string. | |
OGRErr | importFromEPSG (int) |
Initialize SRS based on EPSG geographic, projected or vertical CRS code. | |
OGRErr | importFromEPSGA (int) |
Initialize SRS based on EPSG geographic, projected or vertical CRS code. | |
OGRErr | importFromESRI (char **) |
Import coordinate system from ESRI .prj format(s). | |
OGRErr | importFromPCI (const char *, const char *=nullptr, const double *=nullptr) |
Import coordinate system from PCI projection definition. | |
OGRErr | importFromUSGS (long iProjSys, long iZone, double *padfPrjParams, long iDatum, int nUSGSAngleFormat=USGS_ANGLE_PACKEDDMS) |
Import coordinate system from USGS projection definition. | |
OGRErr | importFromPanorama (long, long, long, double *, bool bNorth=true) |
Import coordinate system from "Panorama" GIS projection definition. | |
OGRErr | importVertCSFromPanorama (int) |
Import vertical coordinate system from "Panorama" GIS projection definition. | |
OGRErr | importFromOzi (const char *const *papszLines) |
Import coordinate system from OziExplorer projection definition. | |
OGRErr | importFromWMSAUTO (const char *pszAutoDef) |
Initialize from WMSAUTO string. | |
OGRErr | importFromXML (const char *) |
Import coordinate system from XML format (GML only currently). | |
OGRErr | importFromDict (const char *pszDict, const char *pszCode) |
Read SRS from WKT dictionary. | |
OGRErr | importFromURN (const char *) |
Initialize from OGC URN. | |
OGRErr | importFromCRSURL (const char *) |
Initialize from OGC URL. | |
OGRErr | importFromERM (const char *pszProj, const char *pszDatum, const char *pszUnits) |
Create OGR WKT from ERMapper projection definitions. | |
OGRErr | importFromUrl (const char *) |
Set spatial reference from a URL. | |
OGRErr | importFromMICoordSys (const char *) |
Import Mapinfo style CoordSys definition. | |
OGRErr | importFromCF1 (CSLConstList papszKeyValues, const char *pszUnits) |
Import a CRS from netCDF CF-1 definitions. | |
OGRErr | morphToESRI () |
Convert in place to ESRI WKT format. | |
OGRErr | morphFromESRI () |
Convert in place from ESRI WKT format. | |
OGRSpatialReference * | convertToOtherProjection (const char *pszTargetProjection, const char *const *papszOptions=nullptr) const |
Convert to another equivalent projection. | |
OGRErr | Validate () const |
Validate CRS imported with importFromWkt() or with modified with direct node manipulations. | |
OGRErr | StripVertical () |
Convert a compound cs into a horizontal CS. | |
bool | StripTOWGS84IfKnownDatumAndAllowed () |
Remove TOWGS84 information if the CRS has a known horizontal datum and this is allowed by the user. | |
bool | StripTOWGS84IfKnownDatum () |
Remove TOWGS84 information if the CRS has a known horizontal datum. | |
int | EPSGTreatsAsLatLong () const |
This method returns TRUE if this geographic coordinate system should be treated as having lat/long coordinate ordering. | |
int | EPSGTreatsAsNorthingEasting () const |
This method returns TRUE if this projected coordinate system should be treated as having northing/easting coordinate ordering. | |
int | GetAxesCount () const |
Return the number of axis of the coordinate system of the CRS. | |
const char * | GetAxis (const char *pszTargetKey, int iAxis, OGRAxisOrientation *peOrientation, double *pdfConvFactor=nullptr) const |
Fetch the orientation of one axis. | |
OGRErr | SetAxes (const char *pszTargetKey, const char *pszXAxisName, OGRAxisOrientation eXAxisOrientation, const char *pszYAxisName, OGRAxisOrientation eYAxisOrientation) |
Set the axes for a coordinate system. | |
OSRAxisMappingStrategy | GetAxisMappingStrategy () const |
Return the data axis to CRS axis mapping strategy. | |
void | SetAxisMappingStrategy (OSRAxisMappingStrategy) |
Set the data axis to CRS axis mapping strategy. | |
const std::vector< int > & | GetDataAxisToSRSAxisMapping () const |
Return the data axis to SRS axis mapping. | |
OGRErr | SetDataAxisToSRSAxisMapping (const std::vector< int > &mapping) |
Set a custom data axis to CRS axis mapping. | |
OGR_SRSNode * | GetRoot () |
Return root node. | |
const OGR_SRSNode * | GetRoot () const |
Return root node. | |
void | SetRoot (OGR_SRSNode *) |
Set the root SRS node. | |
OGR_SRSNode * | GetAttrNode (const char *) |
Find named node in tree. | |
const OGR_SRSNode * | GetAttrNode (const char *) const |
Find named node in tree. | |
const char * | GetAttrValue (const char *, int=0) const |
Fetch indicated attribute of named node. | |
OGRErr | SetNode (const char *, const char *) |
Set attribute value in spatial reference. | |
OGRErr | SetNode (const char *, double) |
Set attribute value in spatial reference. | |
OGRErr | SetLinearUnitsAndUpdateParameters (const char *pszName, double dfInMeters, const char *pszUnitAuthority=nullptr, const char *pszUnitCode=nullptr) |
Set the linear units for the projection. | |
OGRErr | SetLinearUnits (const char *pszName, double dfInMeters) |
Set the linear units for the projection. | |
OGRErr | SetTargetLinearUnits (const char *pszTargetKey, const char *pszName, double dfInMeters, const char *pszUnitAuthority=nullptr, const char *pszUnitCode=nullptr) |
Set the linear units for the projection. | |
double | GetLinearUnits (char **) const |
Fetch linear projection units. | |
double | GetLinearUnits (const char **=nullptr) const |
Fetch linear projection units. | |
double | GetTargetLinearUnits (const char *pszTargetKey, char **ppszRetName) const |
Fetch linear units for target. | |
double | GetTargetLinearUnits (const char *pszTargetKey, const char **ppszRetName=nullptr) const |
Fetch linear units for target. | |
OGRErr | SetAngularUnits (const char *pszName, double dfInRadians) |
Set the angular units for the geographic coordinate system. | |
double | GetAngularUnits (char **) const |
Fetch angular geographic coordinate system units. | |
double | GetAngularUnits (const char **=nullptr) const |
Fetch angular geographic coordinate system units. | |
double | GetPrimeMeridian (char **) const |
Fetch prime meridian info. | |
double | GetPrimeMeridian (const char **=nullptr) const |
Fetch prime meridian info. | |
bool | IsEmpty () const |
Return if the SRS is not set. | |
int | IsGeographic () const |
Check if geographic coordinate system. | |
int | IsDerivedGeographic () const |
Check if the CRS is a derived geographic coordinate system. | |
int | IsProjected () const |
Check if projected coordinate system. | |
int | IsDerivedProjected () const |
Check if the CRS is a derived projected coordinate system. | |
int | IsGeocentric () const |
Check if geocentric coordinate system. | |
bool | IsDynamic () const |
Check if a CRS is a dynamic CRS. | |
bool | HasPointMotionOperation () const |
Check if a CRS has at least an associated point motion operation. | |
int | IsLocal () const |
Check if local coordinate system. | |
int | IsVertical () const |
Check if vertical coordinate system. | |
int | IsCompound () const |
Check if coordinate system is compound. | |
int | IsSameGeogCS (const OGRSpatialReference *) const |
Do the GeogCS'es match? | |
int | IsSameGeogCS (const OGRSpatialReference *, const char *const *papszOptions) const |
Do the GeogCS'es match? | |
int | IsSameVertCS (const OGRSpatialReference *) const |
Do the VertCS'es match? | |
int | IsSame (const OGRSpatialReference *) const |
Do these two spatial references describe the same system ? | |
int | IsSame (const OGRSpatialReference *, const char *const *papszOptions) const |
Do these two spatial references describe the same system ? | |
void | Clear () |
Wipe current definition. | |
OGRErr | SetLocalCS (const char *) |
Set the user visible LOCAL_CS name. | |
OGRErr | SetProjCS (const char *) |
Set the user visible PROJCS name. | |
OGRErr | SetProjection (const char *) |
Set a projection name. | |
OGRErr | SetGeocCS (const char *pszGeocName) |
Set the user visible GEOCCS name. | |
OGRErr | SetGeogCS (const char *pszGeogName, const char *pszDatumName, const char *pszEllipsoidName, double dfSemiMajor, double dfInvFlattening, const char *pszPMName=nullptr, double dfPMOffset=0.0, const char *pszUnits=nullptr, double dfConvertToRadians=0.0) |
Set geographic coordinate system. | |
OGRErr | SetWellKnownGeogCS (const char *) |
Set a GeogCS based on well known name. | |
OGRErr | CopyGeogCSFrom (const OGRSpatialReference *poSrcSRS) |
Copy GEOGCS from another OGRSpatialReference. | |
OGRErr | SetVertCS (const char *pszVertCSName, const char *pszVertDatumName, int nVertDatumClass=2005) |
Set the user visible VERT_CS name. | |
OGRErr | SetCompoundCS (const char *pszName, const OGRSpatialReference *poHorizSRS, const OGRSpatialReference *poVertSRS) |
Setup a compound coordinate system. | |
void | SetCoordinateEpoch (double dfCoordinateEpoch) |
Set the coordinate epoch, as decimal year. | |
double | GetCoordinateEpoch () const |
Return the coordinate epoch, as decimal year. | |
OGRErr | PromoteTo3D (const char *pszName) |
"Promotes" a 2D CRS to a 3D CRS one. | |
OGRErr | DemoteTo2D (const char *pszName) |
"Demote" a 3D CRS to a 2D CRS one. | |
OGRErr | SetFromUserInput (const char *) |
Set spatial reference from various text formats. | |
OGRErr | SetFromUserInput (const char *, CSLConstList papszOptions) |
Set spatial reference from various text formats. | |
OGRErr | SetTOWGS84 (double, double, double, double=0.0, double=0.0, double=0.0, double=0.0) |
Set the Bursa-Wolf conversion to WGS84. | |
OGRErr | GetTOWGS84 (double *padfCoef, int nCoeff=7) const |
Fetch TOWGS84 parameters, if available. | |
OGRErr | AddGuessedTOWGS84 () |
Try to add a a 3-parameter or 7-parameter Helmert transformation to WGS84. | |
double | GetSemiMajor (OGRErr *=nullptr) const |
Get spheroid semi major axis (in metres starting with GDAL 3.0) | |
double | GetSemiMinor (OGRErr *=nullptr) const |
Get spheroid semi minor axis. | |
double | GetInvFlattening (OGRErr *=nullptr) const |
Get spheroid inverse flattening. | |
double | GetEccentricity () const |
Get spheroid eccentricity. | |
double | GetSquaredEccentricity () const |
Get spheroid squared eccentricity. | |
OGRErr | SetAuthority (const char *pszTargetKey, const char *pszAuthority, int nCode) |
Set the authority for a node. | |
OGRErr | AutoIdentifyEPSG () |
Set EPSG authority info if possible. | |
OGRSpatialReferenceH * | FindMatches (char **papszOptions, int *pnEntries, int **ppanMatchConfidence) const |
Try to identify a match between the passed SRS and a related SRS in a catalog. | |
OGRSpatialReference * | FindBestMatch (int nMinimumMatchConfidence=90, const char *pszPreferredAuthority="EPSG", CSLConstList papszOptions=nullptr) const |
Try to identify the best match between the passed SRS and a related SRS in a catalog. | |
int | GetEPSGGeogCS () const |
Try to establish what the EPSG code for this coordinate systems GEOGCS might be. | |
const char * | GetAuthorityCode (const char *pszTargetKey) const |
Get the authority code for a node. | |
const char * | GetAuthorityName (const char *pszTargetKey) const |
Get the authority name for a node. | |
char * | GetOGCURN () const |
Get a OGC URN string describing the CRS, when possible. | |
bool | GetAreaOfUse (double *pdfWestLongitudeDeg, double *pdfSouthLatitudeDeg, double *pdfEastLongitudeDeg, double *pdfNorthLatitudeDeg, const char **ppszAreaName) const |
Return the area of use of the CRS. | |
const char * | GetExtension (const char *pszTargetKey, const char *pszName, const char *pszDefault=nullptr) const |
Fetch extension value. | |
OGRErr | SetExtension (const char *pszTargetKey, const char *pszName, const char *pszValue) |
Set extension value. | |
int | FindProjParm (const char *pszParameter, const OGR_SRSNode *poPROJCS=nullptr) const |
Return the child index of the named projection parameter on its parent PROJCS node. | |
OGRErr | SetProjParm (const char *, double) |
Set a projection parameter value. | |
double | GetProjParm (const char *, double=0.0, OGRErr *=nullptr) const |
Fetch a projection parameter value. | |
OGRErr | SetNormProjParm (const char *, double) |
Set a projection parameter with a normalized value. | |
double | GetNormProjParm (const char *, double=0.0, OGRErr *=nullptr) const |
Fetch a normalized projection parameter value. | |
OGRErr | SetACEA (double dfStdP1, double dfStdP2, double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Albers Conic Equal Area. | |
OGRErr | SetAE (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Azimuthal Equidistant. | |
OGRErr | SetBonne (double dfStdP1, double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Bonne. | |
OGRErr | SetCEA (double dfStdP1, double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Cylindrical Equal Area. | |
OGRErr | SetCS (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Cassini-Soldner. | |
OGRErr | SetEC (double dfStdP1, double dfStdP2, double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Equidistant Conic. | |
OGRErr | SetEckert (int nVariation, double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Eckert I. | |
OGRErr | SetEckertIV (double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Eckert IV. | |
OGRErr | SetEckertVI (double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Eckert VI. | |
OGRErr | SetEquirectangular (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Equirectangular. | |
OGRErr | SetEquirectangular2 (double dfCenterLat, double dfCenterLong, double dfPseudoStdParallel1, double dfFalseEasting, double dfFalseNorthing) |
Equirectangular generalized form : | |
OGRErr | SetGEOS (double dfCentralMeridian, double dfSatelliteHeight, double dfFalseEasting, double dfFalseNorthing) |
Geostationary Satellite. | |
OGRErr | SetGH (double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Goode Homolosine. | |
OGRErr | SetIGH () |
Interrupted Goode Homolosine. | |
OGRErr | SetGS (double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Gall Stereographic. | |
OGRErr | SetGaussSchreiberTMercator (double dfCenterLat, double dfCenterLong, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Gauss Schreiber Transverse Mercator. | |
OGRErr | SetGnomonic (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Gnomonic. | |
OGRErr | SetHOM (double dfCenterLat, double dfCenterLong, double dfAzimuth, double dfRectToSkew, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Hotine Oblique Mercator. | |
OGRErr | SetHOM2PNO (double dfCenterLat, double dfLat1, double dfLong1, double dfLat2, double dfLong2, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Hotine Oblique Mercator 2 points. | |
OGRErr | SetHOMAC (double dfCenterLat, double dfCenterLong, double dfAzimuth, double dfRectToSkew, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Hotine Oblique Mercator Azimuth Center / Variant B. | |
OGRErr | SetLOM (double dfCenterLat, double dfCenterLong, double dfAzimuth, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Laborde Oblique Mercator. | |
OGRErr | SetIWMPolyconic (double dfLat1, double dfLat2, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
International Map of the World Polyconic. | |
OGRErr | SetKrovak (double dfCenterLat, double dfCenterLong, double dfAzimuth, double dfPseudoStdParallelLat, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Krovak Oblique Conic Conformal. | |
OGRErr | SetLAEA (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Lambert Azimuthal Equal-Area. | |
OGRErr | SetLCC (double dfStdP1, double dfStdP2, double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Lambert Conformal Conic. | |
OGRErr | SetLCC1SP (double dfCenterLat, double dfCenterLong, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Lambert Conformal Conic 1SP. | |
OGRErr | SetLCCB (double dfStdP1, double dfStdP2, double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Lambert Conformal Conic (Belgium) | |
OGRErr | SetMC (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Miller Cylindrical. | |
OGRErr | SetMercator (double dfCenterLat, double dfCenterLong, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Mercator 1SP. | |
OGRErr | SetMercator2SP (double dfStdP1, double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Mercator 2SP. | |
OGRErr | SetMollweide (double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Mollweide. | |
OGRErr | SetNZMG (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
New Zealand Map Grid. | |
OGRErr | SetOS (double dfOriginLat, double dfCMeridian, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Oblique Stereographic. | |
OGRErr | SetOrthographic (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Orthographic. | |
OGRErr | SetPolyconic (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Polyconic. | |
OGRErr | SetPS (double dfCenterLat, double dfCenterLong, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Polar Stereographic. | |
OGRErr | SetRobinson (double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Robinson. | |
OGRErr | SetSinusoidal (double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Sinusoidal. | |
OGRErr | SetStereographic (double dfCenterLat, double dfCenterLong, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Stereographic. | |
OGRErr | SetSOC (double dfLatitudeOfOrigin, double dfCentralMeridian, double dfFalseEasting, double dfFalseNorthing) |
Swiss Oblique Cylindrical. | |
OGRErr | SetTM (double dfCenterLat, double dfCenterLong, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Transverse Mercator. | |
OGRErr | SetTMVariant (const char *pszVariantName, double dfCenterLat, double dfCenterLong, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Transverse Mercator variants. | |
OGRErr | SetTMG (double dfCenterLat, double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
Tunesia Mining Grid | |
OGRErr | SetTMSO (double dfCenterLat, double dfCenterLong, double dfScale, double dfFalseEasting, double dfFalseNorthing) |
Transverse Mercator (South Oriented) | |
OGRErr | SetTPED (double dfLat1, double dfLong1, double dfLat2, double dfLong2, double dfFalseEasting, double dfFalseNorthing) |
Two Point Equidistant. | |
OGRErr | SetVDG (double dfCenterLong, double dfFalseEasting, double dfFalseNorthing) |
VanDerGrinten. | |
OGRErr | SetUTM (int nZone, int bNorth=TRUE) |
Universal Transverse Mercator. | |
int | GetUTMZone (int *pbNorth=nullptr) const |
Get utm zone information. | |
OGRErr | SetWagner (int nVariation, double dfCenterLat, double dfFalseEasting, double dfFalseNorthing) |
Wagner I – VII. | |
OGRErr | SetQSC (double dfCenterLat, double dfCenterLong) |
Quadrilateralized Spherical Cube. | |
OGRErr | SetSCH (double dfPegLat, double dfPegLong, double dfPegHeading, double dfPegHgt) |
Spherical, Cross-track, Height. | |
OGRErr | SetVerticalPerspective (double dfTopoOriginLat, double dfTopoOriginLon, double dfTopoOriginHeight, double dfViewPointHeight, double dfFalseEasting, double dfFalseNorthing) |
Vertical Perspective / Near-sided Perspective. | |
OGRErr | SetDerivedGeogCRSWithPoleRotationGRIBConvention (const char *pszCRSName, double dfSouthPoleLat, double dfSouthPoleLon, double dfAxisRotation) |
Pole rotation (GRIB convention) | |
OGRErr | SetDerivedGeogCRSWithPoleRotationNetCDFCFConvention (const char *pszCRSName, double dfGridNorthPoleLat, double dfGridNorthPoleLon, double dfNorthPoleGridLon) |
Pole rotation (netCDF CF convention) | |
OGRErr | SetStatePlane (int nZone, int bNAD83=TRUE, const char *pszOverrideUnitName=nullptr, double dfOverrideUnit=0.0) |
State Plane. | |
OGRErr | ImportFromESRIStatePlaneWKT (int nCode, const char *pszDatumName, const char *pszUnitsName, int nPCSCode, const char *pszCRSName=nullptr) |
ImportFromESRIStatePlaneWKT. | |
OGRErr | ImportFromESRIWisconsinWKT (const char *pszPrjName, double dfCentralMeridian, double dfLatOfOrigin, const char *pszUnitsName, const char *pszCRSName=nullptr) |
ImportFromESRIWisconsinWKT. | |
Static Public Member Functions | |
static void | DestroySpatialReference (OGRSpatialReference *poSRS) |
OGRSpatialReference destructor. | |
static CSLConstList | SET_FROM_USER_INPUT_LIMITATIONS_get () |
Return OGRSpatialReference::SET_FROM_USER_INPUT_LIMITATIONS. | |
static int | IsAngularParameter (const char *) |
Is the passed projection parameter an angular one? | |
static int | IsLongitudeParameter (const char *) |
Is the passed projection parameter an angular longitude (relative to a prime meridian)? | |
static int | IsLinearParameter (const char *) |
Is the passed projection parameter an linear one measured in meters or some similar linear measure. | |
static OGRSpatialReference * | GetWGS84SRS () |
Returns an instance of a SRS object with WGS84 WKT. | |
static OGRSpatialReferenceH | ToHandle (OGRSpatialReference *poSRS) |
Convert a OGRSpatialReference* to a OGRSpatialReferenceH. | |
static OGRSpatialReference * | FromHandle (OGRSpatialReferenceH hSRS) |
Convert a OGRSpatialReferenceH to a OGRSpatialReference*. | |
Static Public Attributes | |
static const char *const | SET_FROM_USER_INPUT_LIMITATIONS [] |
Limitations for OGRSpatialReference::SetFromUserInput(). | |
This class represents an OpenGIS Spatial Reference System, and contains methods for converting between this object organization and well known text (WKT) format.
This object is reference counted as one instance of the object is normally shared between many OGRGeometry objects.
Normally application code can fetch needed parameter values for this SRS using GetAttrValue(), but in special cases the underlying parse tree (or OGR_SRSNode objects) can be accessed more directly.
See the tutorial for more information on how to use this class.
Consult also the OGC WKT Coordinate System Issues page for implementation details of WKT in OGR.
|
explicit |
Constructor.
This constructor takes an optional string argument which if passed should be a WKT representation of an SRS. Passing this is equivalent to not passing it, and then calling importFromWkt() with the WKT string.
Note that newly created objects are given a reference count of one.
Starting with GDAL 3.0, coordinates associated with a OGRSpatialReference object are assumed to be in the order of the axis of the CRS definition (which for example means latitude first, longitude second for geographic CRS belonging to the EPSG authority). It is possible to define a data axis to CRS axis mapping strategy with the SetAxisMappingStrategy() method.
Starting with GDAL 3.5, the OSR_DEFAULT_AXIS_MAPPING_STRATEGY configuration option can be set to "TRADITIONAL_GIS_ORDER" / "AUTHORITY_COMPLIANT" (the later being the default value when the option is not set) to control the value of the data axis to CRS axis mapping strategy when a OSRSpatialReference object is created. Calling SetAxisMappingStrategy() will override this default value.
The C function OSRNewSpatialReference() does the same thing as this constructor.
pszWKT | well known text definition to which the object should be initialized, or NULL (the default). |
OGRSpatialReference::OGRSpatialReference | ( | const OGRSpatialReference & | oOther | ) |
OGRSpatialReference::OGRSpatialReference | ( | OGRSpatialReference && | oOther | ) |
Move constructor.
oOther | other spatial reference |
|
virtual |
OGRSpatialReference destructor.
The C function OSRDestroySpatialReference() does the same thing as this method. Preferred C++ method : OGRSpatialReference::DestroySpatialReference()
OGRErr OGRSpatialReference::AddGuessedTOWGS84 | ( | ) |
Try to add a a 3-parameter or 7-parameter Helmert transformation to WGS84.
This method try to attach a 3-parameter or 7-parameter Helmert transformation to WGS84 when there is one and only one such method available for the CRS. Note: this is more restrictive to how GDAL < 3 worked.
This method is the same as the C function OSRAddGuessedTOWGS84().
OGRSpatialReference & OGRSpatialReference::AssignAndSetThreadSafe | ( | const OGRSpatialReference & | oSource | ) |
Assignment method, with thread-safety.
Same as an assignment operator, but asking also that the *this instance becomes thread-safe.
oSource | SRS to assign to *this |
OGRErr OGRSpatialReference::AutoIdentifyEPSG | ( | ) |
Set EPSG authority info if possible.
This method inspects a WKT definition, and adds EPSG authority nodes where an aspect of the coordinate system can be easily and safely corresponded with an EPSG identifier. In practice, this method will evolve over time. In theory it can add authority nodes for any object (i.e. spheroid, datum, GEOGCS, units, and PROJCS) that could have an authority node. Mostly this is useful to inserting appropriate PROJCS codes for common formulations (like UTM n WGS84).
If it success the OGRSpatialReference is updated in place, and the method return OGRERR_NONE. If the method fails to identify the general coordinate system OGRERR_UNSUPPORTED_SRS is returned but no error message is posted via CPLError().
This method is the same as the C function OSRAutoIdentifyEPSG().
Since GDAL 2.3, the FindMatches() method can also be used for improved matching by researching the EPSG catalog.
void OGRSpatialReference::Clear | ( | ) |
Wipe current definition.
Returns OGRSpatialReference to a state with no definition, as it exists when first created. It does not affect reference counts.
OGRSpatialReference * OGRSpatialReference::Clone | ( | ) | const |
Make a duplicate of this OGRSpatialReference.
This method is the same as the C function OSRClone().
OGRSpatialReference * OGRSpatialReference::CloneGeogCS | ( | ) | const |
Make a duplicate of the GEOGCS node of this OGRSpatialReference object.
OGRSpatialReference * OGRSpatialReference::convertToOtherProjection | ( | const char * | pszTargetProjection, |
const char *const * | papszOptions = nullptr |
||
) | const |
Convert to another equivalent projection.
Currently implemented:
pszTargetProjection | target projection. |
papszOptions | lists of options. None supported currently. |
OGRErr OGRSpatialReference::CopyGeogCSFrom | ( | const OGRSpatialReference * | poSrcSRS | ) |
Copy GEOGCS from another OGRSpatialReference.
The GEOGCS information is copied into this OGRSpatialReference from another. If this object has a PROJCS root already, the GEOGCS is installed within it, otherwise it is installed as the root.
poSrcSRS | the spatial reference to copy the GEOGCS information from. |
OGRErr OGRSpatialReference::DemoteTo2D | ( | const char * | pszName | ) |
"Demote" a 3D CRS to a 2D CRS one.
pszName | New name for the CRS. If set to NULL, the previous name will be used. |
int OGRSpatialReference::Dereference | ( | ) |
Decrements the reference count by one.
The method does the same thing as the C function OSRDereference().
|
static |
OGRSpatialReference destructor.
This static method will destroy a OGRSpatialReference. It is equivalent to calling delete on the object, but it ensures that the deallocation is properly executed within the OGR libraries heap on platforms where this can matter (win32).
This function is the same as OSRDestroySpatialReference()
poSRS | the object to delete |
int OGRSpatialReference::EPSGTreatsAsLatLong | ( | ) | const |
This method returns TRUE if this geographic coordinate system should be treated as having lat/long coordinate ordering.
Currently this returns TRUE for all geographic coordinate systems with axes set defining it as lat, long (prior to GDAL 3.10, it also checked that the CRS had belonged to EPSG authority, but this check has now been removed).
FALSE will be returned for all coordinate systems that are not geographic, or whose axes ordering is not latitude, longitude.
This method is the same as the C function OSREPSGTreatsAsLatLong().
int OGRSpatialReference::EPSGTreatsAsNorthingEasting | ( | ) | const |
This method returns TRUE if this projected coordinate system should be treated as having northing/easting coordinate ordering.
Currently this returns TRUE for all projected coordinate systems with axes set defining it as northing, easting (prior to GDAL 3.10, it also checked that the CRS had belonged to EPSG authority, but this check has now been removed).
FALSE will be returned for all coordinate systems that are not projected, or whose axes ordering is not northing, easting.
This method is the same as the C function EPSGTreatsAsNorthingEasting().
OGRErr OGRSpatialReference::exportToCF1 | ( | char ** | ppszGridMappingName, |
char *** | ppapszKeyValues, | ||
char ** | ppszUnits, | ||
CSLConstList | papszOptions | ||
) | const |
Export a CRS to netCDF CF-1 definitions.
http://cfconventions.org/cf-conventions/cf-conventions.html#appendix-grid-mappings
This function is the equivalent of the C function OSRExportToCF1().
[out] | ppszGridMappingName | Pointer to the suggested name for the grid mapping variable. ppszGridMappingName may be nullptr. *ppszGridMappingName should be freed with CPLFree(). |
[out] | ppapszKeyValues | Pointer to a null-terminated list of key/value pairs, to write into the grid mapping variable. ppapszKeyValues may be nullptr. *ppapszKeyValues should be freed with CSLDestroy() Values may be of type string, double or a list of 2 double values (comma separated). |
[out] | ppszUnits | Pointer to the value of the "units" attribute of the X/Y arrays. ppszGridMappingName may be nullptr. *ppszUnits should be freed with CPLFree(). |
[in] | papszOptions | Options. Currently none supported |
OGRErr OGRSpatialReference::exportToERM | ( | char * | pszProj, |
char * | pszDatum, | ||
char * | pszUnits | ||
) |
Convert coordinate system to ERMapper format.
pszProj | 32 character buffer to receive projection name. |
pszDatum | 32 character buffer to receive datum name. |
pszUnits | 32 character buffer to receive units name. |
OGRErr OGRSpatialReference::exportToMICoordSys | ( | char ** | ppszResult | ) | const |
Export coordinate system in Mapinfo style CoordSys format.
Note that the returned WKT string should be freed with CPLFree() when no longer needed. It is the responsibility of the caller.
This method is the same as the C function OSRExportToMICoordSys().
ppszResult | pointer to which dynamically allocated Mapinfo CoordSys definition will be assigned. |
OGRErr OGRSpatialReference::exportToPanorama | ( | long * | piProjSys, |
long * | piDatum, | ||
long * | piEllips, | ||
long * | piZone, | ||
double * | padfPrjParams | ||
) | const |
Export coordinate system in "Panorama" GIS projection definition.
This method is the equivalent of the C function OSRExportToPanorama().
piProjSys | Pointer to variable, where the projection system code will be returned. |
piDatum | Pointer to variable, where the coordinate system code will be returned. |
piEllips | Pointer to variable, where the spheroid code will be returned. |
piZone | Pointer to variable, where the zone for UTM projection system will be returned. |
padfPrjParams | an existing 7 double buffer into which the projection parameters will be placed. See importFromPanorama() for the list of parameters. |
OGRErr OGRSpatialReference::exportToPCI | ( | char ** | ppszProj, |
char ** | ppszUnits, | ||
double ** | ppadfPrjParams | ||
) | const |
Export coordinate system in PCI projection definition.
Converts the loaded coordinate reference system into PCI projection definition to the extent possible. The strings returned in ppszProj, ppszUnits and ppadfPrjParams array should be deallocated by the caller with CPLFree() when no longer needed.
LOCAL_CS coordinate systems are not translatable. An empty string will be returned along with OGRERR_NONE.
This method is the equivalent of the C function OSRExportToPCI().
ppszProj | pointer to which dynamically allocated PCI projection definition will be assigned. |
ppszUnits | pointer to which dynamically allocated units definition will be assigned. |
ppadfPrjParams | pointer to which dynamically allocated array of 17 projection parameters will be assigned. See importFromPCI() for the list of parameters. |
OGRErr OGRSpatialReference::exportToPrettyWkt | ( | char ** | ppszResult, |
int | bSimplify = FALSE |
||
) | const |
Convert this SRS into a nicely formatted WKT 1 string for display to a person.
Consult also the OGC WKT Coordinate System Issues page for implementation details of WKT 1 in OGR.
Note that the returned WKT string should be freed with CPLFree() when no longer needed. It is the responsibility of the caller.
The WKT version can be overridden by using the OSR_WKT_FORMAT configuration option. Valid values are the one of the FORMAT option of exportToWkt( char ** ppszResult, const char* const* papszOptions ) const
This method is the same as the C function OSRExportToPrettyWkt().
ppszResult | the resulting string is returned in this pointer. |
bSimplify | TRUE if the AXIS, AUTHORITY and EXTENSION nodes should be stripped off. |
OGRErr OGRSpatialReference::exportToProj4 | ( | char ** | ppszProj4 | ) | const |
Export coordinate system in PROJ.4 legacy format.
Converts the loaded coordinate reference system into PROJ format to the extent possible. The string returned in ppszProj4 should be deallocated by the caller with CPLFree() when no longer needed.
LOCAL_CS coordinate systems are not translatable. An empty string will be returned along with OGRERR_NONE.
Special processing for Transverse Mercator: Starting with GDAL 3.0, if the OSR_USE_APPROX_TMERC configuration option is set to YES, the PROJ definition built from the SRS will use the +approx flag for the tmerc and utm projection methods, rather than the more accurate method.
Starting with GDAL 3.0.3, this method will try to add a +towgs84 parameter, if there's none attached yet to the SRS and if the SRS has a EPSG code. See the AddGuessedTOWGS84() method for how this +towgs84 parameter may be added. This automatic addition may be disabled by setting the OSR_ADD_TOWGS84_ON_EXPORT_TO_PROJ4 configuration option to NO.
This method is the equivalent of the C function OSRExportToProj4().
ppszProj4 | pointer to which dynamically allocated PROJ definition will be assigned. |
OGRErr OGRSpatialReference::exportToPROJJSON | ( | char ** | ppszResult, |
const char *const * | papszOptions | ||
) | const |
Convert this SRS into a PROJJSON string.
Note that the returned JSON string should be freed with CPLFree() when no longer needed. It is the responsibility of the caller.
ppszResult | the resulting string is returned in this pointer. |
papszOptions | NULL terminated list of options, or NULL. Currently supported options are
|
OGRErr OGRSpatialReference::exportToUSGS | ( | long * | piProjSys, |
long * | piZone, | ||
double ** | ppadfPrjParams, | ||
long * | piDatum | ||
) | const |
Export coordinate system in USGS GCTP projection definition.
This method is the equivalent of the C function OSRExportToUSGS().
piProjSys | Pointer to variable, where the projection system code will be returned. |
piZone | Pointer to variable, where the zone for UTM and State Plane projection systems will be returned. |
ppadfPrjParams | Pointer to which dynamically allocated array of 15 projection parameters will be assigned. See importFromUSGS() for the list of parameters. Caller responsible to free this array. |
piDatum | Pointer to variable, where the datum code will be returned. |
OGRErr OGRSpatialReference::exportToWkt | ( | char ** | ppszResult | ) | const |
Convert this SRS into WKT 1 format.
Consult also the OGC WKT Coordinate System Issues page for implementation details of WKT 1 in OGR.
Note that the returned WKT string should be freed with CPLFree() when no longer needed. It is the responsibility of the caller.
The WKT version can be overridden by using the OSR_WKT_FORMAT configuration option. Valid values are the one of the FORMAT option of exportToWkt( char ** ppszResult, const char* const* papszOptions ) const
This method is the same as the C function OSRExportToWkt().
ppszResult | the resulting string is returned in this pointer. |
OGRErr OGRSpatialReference::exportToWkt | ( | char ** | ppszResult, |
const char *const * | papszOptions | ||
) | const |
Convert this SRS into a WKT string.
Note that the returned WKT string should be freed with CPLFree() when no longer needed. It is the responsibility of the caller.
Consult also the OGC WKT Coordinate System Issues page for implementation details of WKT 1 in OGR.
ppszResult | the resulting string is returned in this pointer. |
papszOptions | NULL terminated list of options, or NULL. Currently supported options are
|
Starting with GDAL 3.0.3, if the OSR_ADD_TOWGS84_ON_EXPORT_TO_WKT1 configuration option is set to YES, when exporting to WKT1_GDAL, this method will try to add a TOWGS84[] node, if there's none attached yet to the SRS and if the SRS has a EPSG code. See the AddGuessedTOWGS84() method for how this TOWGS84[] node may be added.
std::string OGRSpatialReference::exportToWkt | ( | const char *const * | papszOptions = nullptr | ) | const |
Convert this SRS into a WKT string.
Consult also the OGC WKT Coordinate System Issues page for implementation details of WKT 1 in OGR.
papszOptions | NULL terminated list of options, or NULL. Currently supported options are
|
If the OSR_ADD_TOWGS84_ON_EXPORT_TO_WKT1 configuration option is set to YES, when exporting to WKT1_GDAL, this method will try to add a TOWGS84[] node, if there's none attached yet to the SRS and if the SRS has a EPSG code. See the AddGuessedTOWGS84() method for how this TOWGS84[] node may be added.
OGRErr OGRSpatialReference::exportToXML | ( | char ** | ppszRawXML, |
const char * | pszDialect = nullptr |
||
) | const |
Export coordinate system in XML format.
Converts the loaded coordinate reference system into XML format to the extent possible. The string returned in ppszRawXML should be deallocated by the caller with CPLFree() when no longer needed.
LOCAL_CS coordinate systems are not translatable. An empty string will be returned along with OGRERR_NONE.
This method is the equivalent of the C function OSRExportToXML().
ppszRawXML | pointer to which dynamically allocated XML definition will be assigned. |
pszDialect | currently ignored. The dialect used is GML based. |
OGRSpatialReference * OGRSpatialReference::FindBestMatch | ( | int | nMinimumMatchConfidence = 90 , |
const char * | pszPreferredAuthority = "EPSG" , |
||
CSLConstList | papszOptions = nullptr |
||
) | const |
Try to identify the best match between the passed SRS and a related SRS in a catalog.
This is a wrapper over OGRSpatialReference::FindMatches() that takes care of filtering its output. Only matches whose confidence is greater or equal to nMinimumMatchConfidence will be considered. If there is a single match, it is returned. If there are several matches, only return the one under the pszPreferredAuthority, if there is a single one under that authority.
nMinimumMatchConfidence | Minimum match confidence (value between 0 and 100). If set to 0, 90 is used. |
pszPreferredAuthority | Preferred CRS authority. If set to nullptr, "EPSG" is used. |
papszOptions | NULL terminated list of options or NULL. No option is defined at time of writing. |
OGRSpatialReferenceH * OGRSpatialReference::FindMatches | ( | char ** | papszOptions, |
int * | pnEntries, | ||
int ** | ppanMatchConfidence | ||
) | const |
Try to identify a match between the passed SRS and a related SRS in a catalog.
Matching may be partial, or may fail. Returned entries will be sorted by decreasing match confidence (first entry has the highest match confidence).
The exact way matching is done may change in future versions. Starting with GDAL 3.0, it relies on PROJ' proj_identify() function.
This method is the same as OSRFindMatches().
papszOptions | NULL terminated list of options or NULL |
pnEntries | Output parameter. Number of values in the returned array. |
ppanMatchConfidence | Output parameter (or NULL). *ppanMatchConfidence will be allocated to an array of *pnEntries whose values between 0 and 100 indicate the confidence in the match. 100 is the highest confidence level. The array must be freed with CPLFree(). |
int OGRSpatialReference::FindProjParm | ( | const char * | pszParameter, |
const OGR_SRSNode * | poPROJCS = nullptr |
||
) | const |
Return the child index of the named projection parameter on its parent PROJCS node.
pszParameter | projection parameter to look for |
poPROJCS | projection CS node to look in. If NULL is passed, the PROJCS node of the SpatialReference object will be searched. |
|
inlinestatic |
Convert a OGRSpatialReferenceH to a OGRSpatialReference*.
double OGRSpatialReference::GetAngularUnits | ( | char ** | ppszName | ) | const |
Fetch angular geographic coordinate system units.
If no units are available, a value of "degree" and SRS_UA_DEGREE_CONV will be assumed. This method only checks directly under the GEOGCS node for units.
This method does the same thing as the C function OSRGetAngularUnits().
ppszName | a pointer to be updated with the pointer to the units name. The returned value remains internal to the OGRSpatialReference and should not be freed, or modified. It may be invalidated on the next OGRSpatialReference call. |
double OGRSpatialReference::GetAngularUnits | ( | const char ** | ppszName = nullptr | ) | const |
Fetch angular geographic coordinate system units.
If no units are available, a value of "degree" and SRS_UA_DEGREE_CONV will be assumed. This method only checks directly under the GEOGCS node for units.
This method does the same thing as the C function OSRGetAngularUnits().
ppszName | a pointer to be updated with the pointer to the units name. The returned value remains internal to the OGRSpatialReference and should not be freed, or modified. It may be invalidated on the next OGRSpatialReference call. |
bool OGRSpatialReference::GetAreaOfUse | ( | double * | pdfWestLongitudeDeg, |
double * | pdfSouthLatitudeDeg, | ||
double * | pdfEastLongitudeDeg, | ||
double * | pdfNorthLatitudeDeg, | ||
const char ** | ppszAreaName | ||
) | const |
Return the area of use of the CRS.
This method is the same as the OSRGetAreaOfUse() function.
pdfWestLongitudeDeg | Pointer to a double to receive the western-most longitude, expressed in degree. Might be NULL. If the returned value is -1000, the bounding box is unknown. |
pdfSouthLatitudeDeg | Pointer to a double to receive the southern-most latitude, expressed in degree. Might be NULL. If the returned value is -1000, the bounding box is unknown. |
pdfEastLongitudeDeg | Pointer to a double to receive the eastern-most longitude, expressed in degree. Might be NULL. If the returned value is -1000, the bounding box is unknown. |
pdfNorthLatitudeDeg | Pointer to a double to receive the northern-most latitude, expressed in degree. Might be NULL. If the returned value is -1000, the bounding box is unknown. |
ppszAreaName | Pointer to a string to receive the name of the area of use. Might be NULL. Note that *ppszAreaName is short-lived and might be invalidated by further calls. |
OGR_SRSNode * OGRSpatialReference::GetAttrNode | ( | const char * | pszNodePath | ) |
Find named node in tree.
This method does a pre-order traversal of the node tree searching for a node with this exact value (case insensitive), and returns it. Leaf nodes are not considered, under the assumption that they are just attribute value nodes.
If a node appears more than once in the tree (such as UNIT for instance), the first encountered will be returned. Use GetNode() on a subtree to be more specific.
pszNodePath | the name of the node to search for. May contain multiple components such as "GEOGCS|UNIT". |
const OGR_SRSNode * OGRSpatialReference::GetAttrNode | ( | const char * | pszNodePath | ) | const |
Find named node in tree.
This method does a pre-order traversal of the node tree searching for a node with this exact value (case insensitive), and returns it. Leaf nodes are not considered, under the assumption that they are just attribute value nodes.
If a node appears more than once in the tree (such as UNIT for instance), the first encountered will be returned. Use GetNode() on a subtree to be more specific.
pszNodePath | the name of the node to search for. May contain multiple components such as "GEOGCS|UNIT". |
const char * OGRSpatialReference::GetAttrValue | ( | const char * | pszNodeName, |
int | iAttr = 0 |
||
) | const |
Fetch indicated attribute of named node.
This method uses GetAttrNode() to find the named node, and then extracts the value of the indicated child. Thus a call to GetAttrValue("UNIT",1) would return the second child of the UNIT node, which is normally the length of the linear unit in meters.
This method does the same thing as the C function OSRGetAttrValue().
pszNodeName | the tree node to look for (case insensitive). |
iAttr | the child of the node to fetch (zero based). |
const char * OGRSpatialReference::GetAuthorityCode | ( | const char * | pszTargetKey | ) | const |
Get the authority code for a node.
This method is used to query an AUTHORITY[] node from within the WKT tree, and fetch the code value.
While in theory values may be non-numeric, for the EPSG authority all code values should be integral.
This method is the same as the C function OSRGetAuthorityCode().
pszTargetKey | the partial or complete path to the node to get an authority from. i.e. "PROJCS", "GEOGCS", "GEOGCS|UNIT" or NULL to search for an authority node on the root element. |
const char * OGRSpatialReference::GetAuthorityName | ( | const char * | pszTargetKey | ) | const |
Get the authority name for a node.
This method is used to query an AUTHORITY[] node from within the WKT tree, and fetch the authority name value.
The most common authority is "EPSG".
This method is the same as the C function OSRGetAuthorityName().
pszTargetKey | the partial or complete path to the node to get an authority from. i.e. "PROJCS", "GEOGCS", "GEOGCS|UNIT" or NULL to search for an authority node on the root element. |
int OGRSpatialReference::GetAxesCount | ( | ) | const |
Return the number of axis of the coordinate system of the CRS.
const char * OGRSpatialReference::GetAxis | ( | const char * | pszTargetKey, |
int | iAxis, | ||
OGRAxisOrientation * | peOrientation, | ||
double * | pdfConvUnit = nullptr |
||
) | const |
Fetch the orientation of one axis.
Fetches the request axis (iAxis - zero based) from the indicated portion of the coordinate system (pszTargetKey) which should be either "GEOGCS" or "PROJCS".
No CPLError is issued on routine failures (such as not finding the AXIS).
This method is equivalent to the C function OSRGetAxis().
pszTargetKey | the coordinate system part to query ("PROJCS" or "GEOGCS"). |
iAxis | the axis to query (0 for first, 1 for second, 2 for third). |
peOrientation | location into which to place the fetch orientation, may be NULL. |
pdfConvUnit | (GDAL >= 3.4) Location into which to place axis conversion factor. May be NULL. Only set if pszTargetKey == NULL |
OSRAxisMappingStrategy OGRSpatialReference::GetAxisMappingStrategy | ( | ) | const |
Return the data axis to CRS axis mapping strategy.
double OGRSpatialReference::GetCoordinateEpoch | ( | ) | const |
Return the coordinate epoch, as decimal year.
In a dynamic CRS, coordinates of a point on the surface of the Earth may change with time. To be unambiguous the coordinates must always be qualified with the epoch at which they are valid. The coordinate epoch is not necessarily the epoch at which the observation was collected.
Pedantically the coordinate epoch of an observation belongs to the observation, and not to the CRS, however it is often more practical to bind it to the CRS. The coordinate epoch should be specified for dynamic CRS (see IsDynamic())
This method is the same as the OSRGetCoordinateEpoch() function.
const std::vector< int > & OGRSpatialReference::GetDataAxisToSRSAxisMapping | ( | ) | const |
Return the data axis to SRS axis mapping.
The number of elements of the vector will be the number of axis of the CRS. Values start at 1.
If m = GetDataAxisToSRSAxisMapping(), then m[0] is the data axis number for the first axis of the CRS.
double OGRSpatialReference::GetEccentricity | ( | ) | const |
Get spheroid eccentricity.
int OGRSpatialReference::GetEPSGGeogCS | ( | ) | const |
Try to establish what the EPSG code for this coordinate systems GEOGCS might be.
Returns -1 if no reasonable guess can be made.
const char * OGRSpatialReference::GetExtension | ( | const char * | pszTargetKey, |
const char * | pszName, | ||
const char * | pszDefault = nullptr |
||
) | const |
Fetch extension value.
Fetch the value of the named EXTENSION item for the identified target node.
pszTargetKey | the name or path to the parent node of the EXTENSION. |
pszName | the name of the extension being fetched. |
pszDefault | the value to return if the extension is not found. |
double OGRSpatialReference::GetInvFlattening | ( | OGRErr * | pnErr = nullptr | ) | const |
Get spheroid inverse flattening.
This method does the same thing as the C function OSRGetInvFlattening().
pnErr | if non-NULL set to OGRERR_FAILURE if no inverse flattening can be found. |
double OGRSpatialReference::GetLinearUnits | ( | char ** | ppszName | ) | const |
Fetch linear projection units.
If no units are available, a value of "Meters" and 1.0 will be assumed. This method only checks directly under the PROJCS, GEOCCS, GEOGCS or LOCAL_CS node for units. When called on a Geographic 3D CRS the vertical axis units will be returned.
This method does the same thing as the C function OSRGetLinearUnits()
ppszName | a pointer to be updated with the pointer to the units name. The returned value remains internal to the OGRSpatialReference and should not be freed, or modified. It may be invalidated on the next OGRSpatialReference call. |
double OGRSpatialReference::GetLinearUnits | ( | const char ** | ppszName = nullptr | ) | const |
Fetch linear projection units.
If no units are available, a value of "Meters" and 1.0 will be assumed. This method only checks directly under the PROJCS, GEOCCS or LOCAL_CS node for units.
This method does the same thing as the C function OSRGetLinearUnits()
ppszName | a pointer to be updated with the pointer to the units name. The returned value remains internal to the OGRSpatialReference and should not be freed, or modified. It may be invalidated on the next OGRSpatialReference call. |
const char * OGRSpatialReference::GetName | ( | ) | const |
Return the CRS name.
The returned value is only short lived and should not be used after other calls to methods on this object.
double OGRSpatialReference::GetNormProjParm | ( | const char * | pszName, |
double | dfDefaultValue = 0.0 , |
||
OGRErr * | pnErr = nullptr |
||
) | const |
Fetch a normalized projection parameter value.
This method is the same as GetProjParm() except that the value of the parameter is "normalized" into degrees or meters depending on whether it is linear or angular.
This method is the same as the C function OSRGetNormProjParm().
pszName | the name of the parameter to fetch, from the set of SRS_PP codes in ogr_srs_api.h. |
dfDefaultValue | the value to return if this parameter doesn't exist. |
pnErr | place to put error code on failure. Ignored if NULL. |
char * OGRSpatialReference::GetOGCURN | ( | ) | const |
Get a OGC URN string describing the CRS, when possible.
This method assumes that the CRS has a top-level identifier, or is a compound CRS whose horizontal and vertical parts have a top-level identifier.
double OGRSpatialReference::GetPrimeMeridian | ( | char ** | ppszName | ) | const |
Fetch prime meridian info.
Returns the offset of the prime meridian from greenwich in degrees, and the prime meridian name (if requested). If no PRIMEM value exists in the coordinate system definition a value of "Greenwich" and an offset of 0.0 is assumed.
If the prime meridian name is returned, the pointer is to an internal copy of the name. It should not be freed, altered or depended on after the next OGR call.
This method is the same as the C function OSRGetPrimeMeridian().
ppszName | return location for prime meridian name. If NULL, name is not returned. |
double OGRSpatialReference::GetPrimeMeridian | ( | const char ** | ppszName = nullptr | ) | const |
Fetch prime meridian info.
Returns the offset of the prime meridian from greenwich in degrees, and the prime meridian name (if requested). If no PRIMEM value exists in the coordinate system definition a value of "Greenwich" and an offset of 0.0 is assumed.
If the prime meridian name is returned, the pointer is to an internal copy of the name. It should not be freed, altered or depended on after the next OGR call.
This method is the same as the C function OSRGetPrimeMeridian().
ppszName | return location for prime meridian name. If NULL, name is not returned. |
double OGRSpatialReference::GetProjParm | ( | const char * | pszName, |
double | dfDefaultValue = 0.0 , |
||
OGRErr * | pnErr = nullptr |
||
) | const |
Fetch a projection parameter value.
NOTE: This code should be modified to translate non degree angles into degrees based on the GEOGCS unit. This has not yet been done.
This method is the same as the C function OSRGetProjParm().
pszName | the name of the parameter to fetch, from the set of SRS_PP codes in ogr_srs_api.h. |
dfDefaultValue | the value to return if this parameter doesn't exist. |
pnErr | place to put error code on failure. Ignored if NULL. |
int OGRSpatialReference::GetReferenceCount | ( | ) | const |
Fetch current reference count.
double OGRSpatialReference::GetSemiMajor | ( | OGRErr * | pnErr = nullptr | ) | const |
Get spheroid semi major axis (in metres starting with GDAL 3.0)
This method does the same thing as the C function OSRGetSemiMajor().
pnErr | if non-NULL set to OGRERR_FAILURE if semi major axis can be found. |
double OGRSpatialReference::GetSemiMinor | ( | OGRErr * | pnErr = nullptr | ) | const |
Get spheroid semi minor axis.
This method does the same thing as the C function OSRGetSemiMinor().
pnErr | if non-NULL set to OGRERR_FAILURE if semi minor axis can be found. |
double OGRSpatialReference::GetSquaredEccentricity | ( | ) | const |
Get spheroid squared eccentricity.
double OGRSpatialReference::GetTargetLinearUnits | ( | const char * | pszTargetKey, |
char ** | ppszName | ||
) | const |
Fetch linear units for target.
If no units are available, a value of "Meters" and 1.0 will be assumed.
This method does the same thing as the C function OSRGetTargetLinearUnits()
pszTargetKey | the key to look on. i.e. "PROJCS" or "VERT_CS". Might be NULL, in which case PROJCS will be implied (and if not found, LOCAL_CS, GEOCCS and VERT_CS are looked up) |
ppszName | a pointer to be updated with the pointer to the units name. The returned value remains internal to the OGRSpatialReference and should not be freed, or modified. It may be invalidated on the next OGRSpatialReference call. ppszName can be set to NULL. |
double OGRSpatialReference::GetTargetLinearUnits | ( | const char * | pszTargetKey, |
const char ** | ppszName = nullptr |
||
) | const |
Fetch linear units for target.
If no units are available, a value of "Meters" and 1.0 will be assumed.
This method does the same thing as the C function OSRGetTargetLinearUnits()
pszTargetKey | the key to look on. i.e. "PROJCS" or "VERT_CS". Might be NULL, in which case PROJCS will be implied (and if not found, LOCAL_CS, GEOCCS, GEOGCS and VERT_CS are looked up) |
ppszName | a pointer to be updated with the pointer to the units name. The returned value remains internal to the OGRSpatialReference and should not be freed, or modified. It may be invalidated on the next OGRSpatialReference call. ppszName can be set to NULL. |
OGRErr OGRSpatialReference::GetTOWGS84 | ( | double * | padfCoeff, |
int | nCoeffCount = 7 |
||
) | const |
Fetch TOWGS84 parameters, if available.
The parameters have the same meaning as EPSG transformation 9606 (Position Vector 7-param. transformation).
padfCoeff | array into which up to 7 coefficients are placed. |
nCoeffCount | size of padfCoeff - defaults to 7. |
int OGRSpatialReference::GetUTMZone | ( | int * | pbNorth = nullptr | ) | const |
Get utm zone information.
This is the same as the C function OSRGetUTMZone().
In SWIG bindings (Python, Java, etc) the GetUTMZone() method returns a zone which is negative in the southern hemisphere instead of having the pbNorth flag used in the C and C++ interface.
pbNorth | pointer to in to set to TRUE if northern hemisphere, or FALSE if southern. |
|
static |
Returns an instance of a SRS object with WGS84 WKT.
Note: the instance will have SetAxisMappingStrategy(OAMS_TRADITIONAL_GIS_ORDER)
The reference counter of the returned object is not increased by this operation.
bool OGRSpatialReference::HasPointMotionOperation | ( | ) | const |
Check if a CRS has at least an associated point motion operation.
Some CRS are not formally declared as dynamic, but may behave as such in practice due to the prsence of point motion operation, to perform coordinate epoch changes within the CRS. Typically NAD83(CSRS)v7
OGRErr OGRSpatialReference::importFromCF1 | ( | CSLConstList | papszKeyValues, |
const char * | pszUnits | ||
) |
Import a CRS from netCDF CF-1 definitions.
http://cfconventions.org/cf-conventions/cf-conventions.html#appendix-grid-mappings
This function is the equivalent of the C function OSRImportFromCF1().
papszKeyValues | Key/value pairs from the grid mapping variable. Multi-valued parameters (typically "standard_parallel") should be comma separated. |
pszUnits | Value of the "units" attribute of the X/Y arrays. May be nullptr |
OGRErr OGRSpatialReference::importFromCRSURL | ( | const char * | pszURL | ) |
Initialize from OGC URL.
Initializes this spatial reference from a coordinate system defined by an OGC URL prefixed with "http://opengis.net/def/crs" per best practice paper 11-135. Currently EPSG and OGC authority values are supported, including OGC auto codes, but not including CRS1 or CRS88 (NAVD88).
This method is also supported through SetFromUserInput() which can normally be used for URLs.
pszURL | the URL string. |
OGRErr OGRSpatialReference::importFromDict | ( | const char * | pszDictFile, |
const char * | pszCode | ||
) |
Read SRS from WKT dictionary.
This method will attempt to find the indicated coordinate system identity in the indicated dictionary file. If found, the WKT representation is imported and used to initialize this OGRSpatialReference.
More complete information on the format of the dictionary files can be found in the epsg.wkt file in the GDAL data tree. The dictionary files are searched for in the "GDAL" domain using CPLFindFile(). Normally this results in searching /usr/local/share/gdal or somewhere similar.
This method is the same as the C function OSRImportFromDict().
pszDictFile | the name of the dictionary file to load. |
pszCode | the code to lookup in the dictionary. |
OGRErr OGRSpatialReference::importFromEPSG | ( | int | nCode | ) |
Initialize SRS based on EPSG geographic, projected or vertical CRS code.
This method will initialize the spatial reference based on the passed in EPSG CRS code found in the PROJ database.
This method is the same as the C function OSRImportFromEPSG().
Before GDAL 3.0.3, this method would try to attach a 3-parameter or 7-parameter Helmert transformation to WGS84 when there is one and only one such method available for the CRS. This behavior might not always be desirable, so starting with GDAL 3.0.3, this is no longer done unless the OSR_ADD_TOWGS84_ON_IMPORT_FROM_EPSG configuration option is set to YES.
nCode | a GCS or PCS code from the horizontal coordinate system table. |
OGRErr OGRSpatialReference::importFromEPSGA | ( | int | nCode | ) |
Initialize SRS based on EPSG geographic, projected or vertical CRS code.
This method will initialize the spatial reference based on the passed in EPSG CRS code found in the PROJ database.
Since GDAL 3.0, this method is identical to importFromEPSG().
Before GDAL 3.0.3, this method would try to attach a 3-parameter or 7-parameter Helmert transformation to WGS84 when there is one and only one such method available for the CRS. This behavior might not always be desirable, so starting with GDAL 3.0.3, this is no longer done unless the OSR_ADD_TOWGS84_ON_IMPORT_FROM_EPSG configuration option is set to YES. The AddGuessedTOWGS84() method can also be used for that purpose.
The method will also by default substitute a deprecated EPSG code by its non-deprecated replacement. If this behavior is not desired, the OSR_USE_NON_DEPRECATED configuration option can be set to NO.
This method is the same as the C function OSRImportFromEPSGA().
nCode | a CRS code. |
OGRErr OGRSpatialReference::importFromERM | ( | const char * | pszProj, |
const char * | pszDatum, | ||
const char * | pszUnits | ||
) |
Create OGR WKT from ERMapper projection definitions.
Generates an OGRSpatialReference definition from an ERMapper datum and projection name. Based on the ecw_cs.wkt dictionary file from gdal/data.
pszProj | the projection name, such as "NUTM11" or "GEOGRAPHIC". |
pszDatum | the datum name, such as "NAD83". |
pszUnits | the linear units "FEET" or "METERS". |
OGRErr OGRSpatialReference::importFromESRI | ( | char ** | papszPrj | ) |
Import coordinate system from ESRI .prj format(s).
This function will read the text loaded from an ESRI .prj file, and translate it into an OGRSpatialReference definition. This should support many (but by no means all) old style (Arc/Info 7.x) .prj files, as well as the newer pseudo-OGC WKT .prj files. Note that new style .prj files are in OGC WKT format, but require some manipulation to correct datum names, and units on some projection parameters. This is addressed within importFromESRI() by an automatic call to morphFromESRI().
Currently only GEOGRAPHIC, UTM, STATEPLANE, GREATBRITIAN_GRID, ALBERS, EQUIDISTANT_CONIC, TRANSVERSE (mercator), POLAR, LAMBERT (Conic Conformal), LAMBERT_AZIMUTHAL, MERCATOR and POLYCONIC projections are supported from old style files.
At this time there is no equivalent exportToESRI() method. Writing old style .prj files is not supported by OGRSpatialReference. However the morphToESRI() and exportToWkt() methods can be used to generate output suitable to write to new style (Arc 8) .prj files.
This function is the equivalent of the C function OSRImportFromESRI().
papszPrj | NULL terminated list of strings containing the definition. |
OGRErr OGRSpatialReference::importFromMICoordSys | ( | const char * | pszCoordSys | ) |
Import Mapinfo style CoordSys definition.
The OGRSpatialReference is initialized from the passed Mapinfo style CoordSys definition string.
This method is the equivalent of the C function OSRImportFromMICoordSys().
pszCoordSys | Mapinfo style CoordSys definition string. |
OGRErr OGRSpatialReference::importFromOzi | ( | const char *const * | papszLines | ) |
Import coordinate system from OziExplorer projection definition.
This method will import projection definition in style, used by OziExplorer software.
papszLines | Map file lines. This is an array of strings containing the whole OziExplorer .MAP file. The array is terminated by a NULL pointer. |
OGRErr OGRSpatialReference::importFromPanorama | ( | long | iProjSys, |
long | iDatum, | ||
long | iEllips, | ||
double * | padfPrjParams, | ||
bool | bNorth = true |
||
) |
Import coordinate system from "Panorama" GIS projection definition.
This method will import projection definition in style, used by "Panorama" GIS.
This function is the equivalent of the C function OSRImportFromPanorama().
iProjSys | Input projection system code, used in GIS "Panorama". |
Supported Projections are:
iDatum | Input coordinate system. |
Supported Datums are:
iEllips | Input spheroid. |
Supported Spheroids are:
padfPrjParams | Array of 8 coordinate system parameters: |
bNorth | If northern hemisphere true, else false. Defaults to true. |
Particular projection uses different parameters, unused ones may be set to zero. If NULL supplied instead of array pointer default values will be used (i.e., zeroes).
OGRErr OGRSpatialReference::importFromPCI | ( | const char * | pszProj, |
const char * | pszUnits = nullptr , |
||
const double * | padfPrjParams = nullptr |
||
) |
Import coordinate system from PCI projection definition.
PCI software uses 16-character string to specify coordinate system and datum/ellipsoid. You should supply at least this string to the importFromPCI() function.
This function is the equivalent of the C function OSRImportFromPCI().
pszProj | NULL terminated string containing the definition. Looks like "pppppppppppp Ennn" or "pppppppppppp Dnnn", where "pppppppppppp" is a projection code, "Ennn" is an ellipsoid code, "Dnnn" a datum code. |
pszUnits | Grid units code ("DEGREE" or "METRE"). If NULL "METRE" will be used. |
padfPrjParams | Array of 17 coordinate system parameters: |
[0] Spheroid semi major axis [1] Spheroid semi minor axis [2] Reference Longitude [3] Reference Latitude [4] First Standard Parallel [5] Second Standard Parallel [6] False Easting [7] False Northing [8] Scale Factor [9] Height above sphere surface [10] Longitude of 1st point on center line [11] Latitude of 1st point on center line [12] Longitude of 2nd point on center line [13] Latitude of 2nd point on center line [14] Azimuth east of north for center line [15] Landsat satellite number [16] Landsat path number
Particular projection uses different parameters, unused ones may be set to zero. If NULL is supplied instead of an array pointer, default values will be used (i.e., zeroes).
OGRErr OGRSpatialReference::importFromProj4 | ( | const char * | pszProj4 | ) |
Import PROJ coordinate string.
The OGRSpatialReference is initialized from the passed PROJs style coordinate system string.
Example: pszProj4 = "+proj=utm +zone=11 +datum=WGS84"
It is also possible to import "+init=epsg:n" style definitions. Those are a legacy syntax that should be avoided in the future. In particular they will result in CRS objects whose axis order might not correspond to the official EPSG axis order.
This method is the equivalent of the C function OSRImportFromProj4().
pszProj4 | the PROJ style string. |
OGRErr OGRSpatialReference::importFromUrl | ( | const char * | pszUrl | ) |
Set spatial reference from a URL.
This method will download the spatial reference at a given URL and feed it into SetFromUserInput for you.
This method does the same thing as the OSRImportFromUrl() function.
pszUrl | text definition to try to deduce SRS from. |
OGRErr OGRSpatialReference::importFromURN | ( | const char * | pszURN | ) |
Initialize from OGC URN.
Initializes this spatial reference from a coordinate system defined by an OGC URN prefixed with "urn:ogc:def:crs:" per recommendation paper 06-023r1. Currently EPSG and OGC authority values are supported, including OGC auto codes, but not including CRS1 or CRS88 (NAVD88).
This method is also support through SetFromUserInput() which can normally be used for URNs.
pszURN | the urn string. |
OGRErr OGRSpatialReference::importFromUSGS | ( | long | iProjSys, |
long | iZone, | ||
double * | padfPrjParams, | ||
long | iDatum, | ||
int | nUSGSAngleFormat = USGS_ANGLE_PACKEDDMS |
||
) |
Import coordinate system from USGS projection definition.
This method will import projection definition in style, used by USGS GCTP software. GCTP operates on angles in packed DMS format (see CPLDecToPackedDMS() function for details), so all angle values (latitudes, longitudes, azimuths, etc.) specified in the padfPrjParams array should be in the packed DMS format, unless bAnglesInPackedDMSFormat is set to FALSE.
This function is the equivalent of the C function OSRImportFromUSGS(). Note that the bAnglesInPackedDMSFormat parameter is only present in the C++ method. The C function assumes bAnglesInPackedFormat = TRUE.
iProjSys | Input projection system code, used in GCTP. |
iZone | Input zone for UTM and State Plane projection systems. For Southern Hemisphere UTM use a negative zone code. iZone ignored for all other projections. |
padfPrjParams | Array of 15 coordinate system parameters. These parameters differs for different projections. |
Projection Transformation Package Projection Parameters: ---------------------------------------------------------------------------- | Array Element Code & Projection Id |--------------------------------------------------- | 0 | 1 | 2 | 3 | 4 | 5 |6 | 7 ---------------------------------------------------------------------------- 0 Geographic | | | | | | | | 1 U T M |Lon/Z |Lat/Z | | | | | | 2 State Plane | | | | | | | | 3 Albers Equal Area |SMajor|SMinor|STDPR1|STDPR2|CentMer|OriginLat|FE|FN 4 Lambert Conformal C |SMajor|SMinor|STDPR1|STDPR2|CentMer|OriginLat|FE|FN 5 Mercator |SMajor|SMinor| | |CentMer|TrueScale|FE|FN 6 Polar Stereographic |SMajor|SMinor| | |LongPol|TrueScale|FE|FN 7 Polyconic |SMajor|SMinor| | |CentMer|OriginLat|FE|FN 8 Equid. Conic A |SMajor|SMinor|STDPAR| |CentMer|OriginLat|FE|FN Equid. Conic B |SMajor|SMinor|STDPR1|STDPR2|CentMer|OriginLat|FE|FN 9 Transverse Mercator |SMajor|SMinor|Factor| |CentMer|OriginLat|FE|FN 10 Stereographic |Sphere| | | |CentLon|CenterLat|FE|FN 11 Lambert Azimuthal |Sphere| | | |CentLon|CenterLat|FE|FN 12 Azimuthal |Sphere| | | |CentLon|CenterLat|FE|FN 13 Gnomonic |Sphere| | | |CentLon|CenterLat|FE|FN 14 Orthographic |Sphere| | | |CentLon|CenterLat|FE|FN 15 Gen. Vert. Near Per |Sphere| |Height| |CentLon|CenterLat|FE|FN 16 Sinusoidal |Sphere| | | |CentMer| |FE|FN 17 Equirectangular |Sphere| | | |CentMer|TrueScale|FE|FN 18 Miller Cylindrical |Sphere| | | |CentMer| |FE|FN 19 Van der Grinten |Sphere| | | |CentMer|OriginLat|FE|FN 20 Hotin Oblique Merc A |SMajor|SMinor|Factor| | |OriginLat|FE|FN Hotin Oblique Merc B |SMajor|SMinor|Factor|AziAng|AzmthPt|OriginLat|FE|FN 21 Robinson |Sphere| | | |CentMer| |FE|FN 22 Space Oblique Merc A |SMajor|SMinor| |IncAng|AscLong| |FE|FN Space Oblique Merc B |SMajor|SMinor|Satnum|Path | | |FE|FN 23 Alaska Conformal |SMajor|SMinor| | | | |FE|FN 24 Interrupted Goode |Sphere| | | | | | | 25 Mollweide |Sphere| | | |CentMer| |FE|FN 26 Interrupt Mollweide |Sphere| | | | | | | 27 Hammer |Sphere| | | |CentMer| |FE|FN 28 Wagner IV |Sphere| | | |CentMer| |FE|FN 29 Wagner VII |Sphere| | | |CentMer| |FE|FN 30 Oblated Equal Area |Sphere| |Shapem|Shapen|CentLon|CenterLat|FE|FN ---------------------------------------------------------------------------- ---------------------------------------------------- | Array Element | Code & Projection Id |--------------------------- | 8 | 9 | 10 | 11 | 12 | ---------------------------------------------------- 0 Geographic | | | | | | 1 U T M | | | | | | 2 State Plane | | | | | | 3 Albers Equal Area | | | | | | 4 Lambert Conformal C | | | | | | 5 Mercator | | | | | | 6 Polar Stereographic | | | | | | 7 Polyconic | | | | | | 8 Equid. Conic A |zero | | | | | Equid. Conic B |one | | | | | 9 Transverse Mercator | | | | | | 10 Stereographic | | | | | | 11 Lambert Azimuthal | | | | | | 12 Azimuthal | | | | | | 13 Gnomonic | | | | | | 14 Orthographic | | | | | | 15 Gen. Vert. Near Per | | | | | | 16 Sinusoidal | | | | | | 17 Equirectangular | | | | | | 18 Miller Cylindrical | | | | | | 19 Van der Grinten | | | | | | 20 Hotin Oblique Merc A |Long1|Lat1|Long2|Lat2|zero| Hotin Oblique Merc B | | | | |one | 21 Robinson | | | | | | 22 Space Oblique Merc A |PSRev|LRat|PFlag| |zero| Space Oblique Merc B | | | | |one | 23 Alaska Conformal | | | | | | 24 Interrupted Goode | | | | | | 25 Mollweide | | | | | | 26 Interrupt Mollweide | | | | | | 27 Hammer | | | | | | 28 Wagner IV | | | | | | 29 Wagner VII | | | | | | 30 Oblated Equal Area |Angle| | | | | ---------------------------------------------------- where Lon/Z Longitude of any point in the UTM zone or zero. If zero, a zone code must be specified. Lat/Z Latitude of any point in the UTM zone or zero. If zero, a zone code must be specified. SMajor Semi-major axis of ellipsoid. If zero, Clarke 1866 in meters is assumed. SMinor Eccentricity squared of the ellipsoid if less than zero, if zero, a spherical form is assumed, or if greater than zero, the semi-minor axis of ellipsoid. Sphere Radius of reference sphere. If zero, 6370997 meters is used. STDPAR Latitude of the standard parallel STDPR1 Latitude of the first standard parallel STDPR2 Latitude of the second standard parallel CentMer Longitude of the central meridian OriginLat Latitude of the projection origin FE False easting in the same units as the semi-major axis FN False northing in the same units as the semi-major axis TrueScale Latitude of true scale LongPol Longitude down below pole of map Factor Scale factor at central meridian (Transverse Mercator) or center of projection (Hotine Oblique Mercator) CentLon Longitude of center of projection CenterLat Latitude of center of projection Height Height of perspective point Long1 Longitude of first point on center line (Hotine Oblique Mercator, format A) Long2 Longitude of second point on center line (Hotine Oblique Mercator, format A) Lat1 Latitude of first point on center line (Hotine Oblique Mercator, format A) Lat2 Latitude of second point on center line (Hotine Oblique Mercator, format A) AziAng Azimuth angle east of north of center line (Hotine Oblique Mercator, format B) AzmthPt Longitude of point on central meridian where azimuth occurs (Hotine Oblique Mercator, format B) IncAng Inclination of orbit at ascending node, counter-clockwise from equator (SOM, format A) AscLong Longitude of ascending orbit at equator (SOM, format A) PSRev Period of satellite revolution in minutes (SOM, format A) LRat Landsat ratio to compensate for confusion at northern end of orbit (SOM, format A -- use 0.5201613) PFlag End of path flag for Landsat: 0 = start of path, 1 = end of path (SOM, format A) Satnum Landsat Satellite Number (SOM, format B) Path Landsat Path Number (Use WRS-1 for Landsat 1, 2 and 3 and WRS-2 for Landsat 4, 5 and 6.) (SOM, format B) Shapem Oblated Equal Area oval shape parameter m Shapen Oblated Equal Area oval shape parameter n Angle Oblated Equal Area oval rotation angle Array elements 13 and 14 are set to zero. All array elements with blank fields are set to zero too.
iDatum | Input spheroid. |
If the datum code is negative, the first two values in the parameter array (param) are used to define the values as follows:
If padfPrjParams[0] is a non-zero value and padfPrjParams[1] is greater than one, the semimajor axis is set to padfPrjParams[0] and the semiminor axis is set to padfPrjParams[1].
If padfPrjParams[0] is nonzero and padfPrjParams[1] is greater than zero but less than or equal to one, the semimajor axis is set to padfPrjParams[0] and the semiminor axis is computed from the eccentricity squared value padfPrjParams[1]:
semiminor = sqrt(1.0 - ES)semimajor
where
ES = eccentricity squared
If padfPrjParams[0] is nonzero and padfPrjParams[1] is equal to zero, the semimajor axis and semiminor axis are set to padfPrjParams[0].
If padfPrjParams[0] equals zero and padfPrjParams[1] is greater than zero, the default Clarke 1866 is used to assign values to the semimajor axis and semiminor axis.
If padfPrjParams[0] and padfPrjParams[1] equals zero, the semimajor axis is set to 6370997.0 and the semiminor axis is set to zero.
If a datum code is zero or greater, the semimajor and semiminor axis are defined by the datum code as found in the following table:
Supported Datums are:
nUSGSAngleFormat | one of USGS_ANGLE_DECIMALDEGREES, USGS_ANGLE_PACKEDDMS, or USGS_ANGLE_RADIANS (default is USGS_ANGLE_PACKEDDMS). |
OGRErr OGRSpatialReference::importFromWkt | ( | char ** | ppszInput | ) |
Import from WKT string.
This method will wipe the existing SRS definition, and reassign it based on the contents of the passed WKT string. Only as much of the input string as needed to construct this SRS is consumed from the input string, and the input string pointer is then updated to point to the remaining (unused) input.
Consult also the OGC WKT Coordinate System Issues page for implementation details of WKT in OGR.
This method is the same as the C function OSRImportFromWkt().
ppszInput | Pointer to pointer to input. The pointer is updated to point to remaining unused input text. |
OGRErr OGRSpatialReference::importFromWkt | ( | const char * | pszInput | ) |
Import from WKT string.
This method will wipe the existing SRS definition, and reassign it based on the contents of the passed WKT string. Only as much of the input string as needed to construct this SRS is consumed from the input string, and the input string pointer is then updated to point to the remaining (unused) input.
Consult also the OGC WKT Coordinate System Issues page for implementation details of WKT in OGR.
pszInput | Input WKT |
OGRErr OGRSpatialReference::importFromWkt | ( | const char ** | ppszInput | ) |
Import from WKT string.
This method will wipe the existing SRS definition, and reassign it based on the contents of the passed WKT string. Only as much of the input string as needed to construct this SRS is consumed from the input string, and the input string pointer is then updated to point to the remaining (unused) input.
Starting with PROJ 9.2, if invoked on a COORDINATEMETADATA[] construct, the CRS contained in it will be used to fill the OGRSpatialReference object, and the coordinate epoch potentially present used as the coordinate epoch property of the OGRSpatialReference object.
Consult also the OGC WKT Coordinate System Issues page for implementation details of WKT in OGR.
This method is the same as the C function OSRImportFromWkt().
ppszInput | Pointer to pointer to input. The pointer is updated to point to remaining unused input text. |
OGRErr OGRSpatialReference::importFromWMSAUTO | ( | const char * | pszDefinition | ) |
Initialize from WMSAUTO string.
Note that the WMS 1.3 specification does not include the units code, while apparently earlier specs do. We try to guess around this.
pszDefinition | the WMSAUTO string |
OGRErr OGRSpatialReference::importFromXML | ( | const char * | pszXML | ) |
Import coordinate system from XML format (GML only currently).
This method is the same as the C function OSRImportFromXML()
pszXML | XML string to import |
OGRErr OGRSpatialReference::importVertCSFromPanorama | ( | int | iVCS | ) |
Import vertical coordinate system from "Panorama" GIS projection definition.
iVCS | Input vertical coordinate system ID. |
Supported VCS are:
|
static |
Is the passed projection parameter an angular one?
int OGRSpatialReference::IsCompound | ( | ) | const |
Check if coordinate system is compound.
This method is the same as the C function OSRIsCompound().
int OGRSpatialReference::IsDerivedGeographic | ( | ) | const |
Check if the CRS is a derived geographic coordinate system.
(for example a rotated long/lat grid)
This method is the same as the C function OSRIsDerivedGeographic().
int OGRSpatialReference::IsDerivedProjected | ( | ) | const |
Check if the CRS is a derived projected coordinate system.
This method is the same as the C function OSRIsDerivedGeographic().
bool OGRSpatialReference::IsDynamic | ( | ) | const |
Check if a CRS is a dynamic CRS.
A dynamic CRS relies on a dynamic datum, that is a datum that is not plate-fixed.
This method is the same as the C function OSRIsDynamic().
int OGRSpatialReference::IsGeocentric | ( | ) | const |
Check if geocentric coordinate system.
This method is the same as the C function OSRIsGeocentric().
int OGRSpatialReference::IsGeographic | ( | ) | const |
Check if geographic coordinate system.
This method is the same as the C function OSRIsGeographic().
|
static |
Is the passed projection parameter an linear one measured in meters or some similar linear measure.
int OGRSpatialReference::IsLocal | ( | ) | const |
Check if local coordinate system.
This method is the same as the C function OSRIsLocal().
|
static |
Is the passed projection parameter an angular longitude (relative to a prime meridian)?
int OGRSpatialReference::IsProjected | ( | ) | const |
Check if projected coordinate system.
This method is the same as the C function OSRIsProjected().
int OGRSpatialReference::IsSame | ( | const OGRSpatialReference * | poOtherSRS | ) | const |
Do these two spatial references describe the same system ?
poOtherSRS | the SRS being compared to. |
int OGRSpatialReference::IsSame | ( | const OGRSpatialReference * | poOtherSRS, |
const char *const * | papszOptions | ||
) | const |
Do these two spatial references describe the same system ?
This also takes into account the data axis to CRS axis mapping by default
poOtherSRS | the SRS being compared to. |
papszOptions | options. NULL or NULL terminated list of options. Currently supported options are:
|
int OGRSpatialReference::IsSameGeogCS | ( | const OGRSpatialReference * | poOther | ) | const |
Do the GeogCS'es match?
This method is the same as the C function OSRIsSameGeogCS().
poOther | the SRS being compared against. |
int OGRSpatialReference::IsSameGeogCS | ( | const OGRSpatialReference * | poOther, |
const char *const * | papszOptions | ||
) | const |
Do the GeogCS'es match?
This method is the same as the C function OSRIsSameGeogCS().
poOther | the SRS being compared against. |
papszOptions | options. ignored |
int OGRSpatialReference::IsSameVertCS | ( | const OGRSpatialReference * | poOther | ) | const |
Do the VertCS'es match?
This method is the same as the C function OSRIsSameVertCS().
poOther | the SRS being compared against. |
int OGRSpatialReference::IsVertical | ( | ) | const |
Check if vertical coordinate system.
This method is the same as the C function OSRIsVertical().
OGRErr OGRSpatialReference::morphFromESRI | ( | ) |
Convert in place from ESRI WKT format.
The value notes of this coordinate system are modified in various manners to adhere more closely to the WKT standard. This mostly involves translating a variety of ESRI names for projections, arguments and datums to "standard" names, as defined by Adam Gawne-Cain's reference translation of EPSG to WKT for the CT specification.
This does the same as the C function OSRMorphFromESRI().
OGRErr OGRSpatialReference::morphToESRI | ( | ) |
Convert in place to ESRI WKT format.
The value nodes of this coordinate system are modified in various manners more closely map onto the ESRI concept of WKT format. This includes renaming a variety of projections and arguments, and stripping out nodes note recognised by ESRI (like AUTHORITY and AXIS).
This does the same as the C function OSRMorphToESRI().
OGRSpatialReference & OGRSpatialReference::operator= | ( | const OGRSpatialReference & | oSource | ) |
Assignment operator.
oSource | SRS to assign to *this |
OGRSpatialReference & OGRSpatialReference::operator= | ( | OGRSpatialReference && | oSource | ) |
Move assignment operator.
oSource | SRS to assign to *this |
OGRErr OGRSpatialReference::PromoteTo3D | ( | const char * | pszName | ) |
"Promotes" a 2D CRS to a 3D CRS one.
The new axis will be ellipsoidal height, oriented upwards, and with metre units.
pszName | New name for the CRS. If set to NULL, the previous name will be used. |
int OGRSpatialReference::Reference | ( | ) |
Increments the reference count by one.
The reference count is used keep track of the number of OGRGeometry objects referencing this SRS.
The method does the same thing as the C function OSRReference().
void OGRSpatialReference::Release | ( | ) |
Decrements the reference count by one, and destroy if zero.
The method does the same thing as the C function OSRRelease().
OGRErr OGRSpatialReference::SetAngularUnits | ( | const char * | pszUnitsName, |
double | dfInRadians | ||
) |
Set the angular units for the geographic coordinate system.
This method creates a UNIT subnode with the specified values as a child of the GEOGCS node.
This method does the same as the C function OSRSetAngularUnits().
pszUnitsName | the units name to be used. Some preferred units names can be found in ogr_srs_api.h such as SRS_UA_DEGREE. |
dfInRadians | the value to multiple by an angle in the indicated units to transform to radians. Some standard conversion factors can be found in ogr_srs_api.h. |
OGRErr OGRSpatialReference::SetAuthority | ( | const char * | pszTargetKey, |
const char * | pszAuthority, | ||
int | nCode | ||
) |
Set the authority for a node.
This method is the same as the C function OSRSetAuthority().
pszTargetKey | the partial or complete path to the node to set an authority on. i.e. "PROJCS", "GEOGCS" or "GEOGCS|UNIT". |
pszAuthority | authority name, such as "EPSG". |
nCode | code for value with this authority. |
OGRErr OGRSpatialReference::SetAxes | ( | const char * | pszTargetKey, |
const char * | pszXAxisName, | ||
OGRAxisOrientation | eXAxisOrientation, | ||
const char * | pszYAxisName, | ||
OGRAxisOrientation | eYAxisOrientation | ||
) |
Set the axes for a coordinate system.
Set the names, and orientations of the axes for either a projected (PROJCS) or geographic (GEOGCS) coordinate system.
This method is equivalent to the C function OSRSetAxes().
pszTargetKey | either "PROJCS" or "GEOGCS", must already exist in SRS. |
pszXAxisName | name of first axis, normally "Long" or "Easting". |
eXAxisOrientation | normally OAO_East. |
pszYAxisName | name of second axis, normally "Lat" or "Northing". |
eYAxisOrientation | normally OAO_North. |
void OGRSpatialReference::SetAxisMappingStrategy | ( | OSRAxisMappingStrategy | strategy | ) |
Set the data axis to CRS axis mapping strategy.
Starting with GDAL 3.5, the OSR_DEFAULT_AXIS_MAPPING_STRATEGY configuration option can be set to "TRADITIONAL_GIS_ORDER" / "AUTHORITY_COMPLIANT" (the later being the default value when the option is not set) to control the value of the data axis to CRS axis mapping strategy when a OSRSpatialReference object is created. Calling SetAxisMappingStrategy() will override this default value.
See OGRSpatialReference::GetAxisMappingStrategy()
OGRErr OGRSpatialReference::SetCompoundCS | ( | const char * | pszName, |
const OGRSpatialReference * | poHorizSRS, | ||
const OGRSpatialReference * | poVertSRS | ||
) |
Setup a compound coordinate system.
This method is the same as the C function OSRSetCompoundCS().
This method is replace the current SRS with a COMPD_CS coordinate system consisting of the passed in horizontal and vertical coordinate systems.
pszName | the name of the compound coordinate system. |
poHorizSRS | the horizontal SRS (PROJCS or GEOGCS). |
poVertSRS | the vertical SRS (VERT_CS). |
void OGRSpatialReference::SetCoordinateEpoch | ( | double | dfCoordinateEpoch | ) |
Set the coordinate epoch, as decimal year.
In a dynamic CRS, coordinates of a point on the surface of the Earth may change with time. To be unambiguous the coordinates must always be qualified with the epoch at which they are valid. The coordinate epoch is not necessarily the epoch at which the observation was collected.
Pedantically the coordinate epoch of an observation belongs to the observation, and not to the CRS, however it is often more practical to bind it to the CRS. The coordinate epoch should be specified for dynamic CRS (see IsDynamic())
This method is the same as the OSRSetCoordinateEpoch() function.
dfCoordinateEpoch | Coordinate epoch as decimal year (e.g. 2021.3) |
OGRErr OGRSpatialReference::SetDataAxisToSRSAxisMapping | ( | const std::vector< int > & | mapping | ) |
Set a custom data axis to CRS axis mapping.
The number of elements of the mapping vector should be the number of axis of the CRS (as returned by GetAxesCount()) (although this method does not check that, beyond checking there are at least 2 elements, so that this method and setting the CRS can be done in any order). This is taken into account by OGRCoordinateTransformation to transform the order of coordinates to the order expected by the CRS before transformation, and back to the data order after transformation.
The mapping[i] value (one based) represents the data axis number for the i(th) axis of the CRS. A negative value can also be used to ask for a sign reversal during coordinate transformation (to deal with northing vs southing, easting vs westing, heights vs depths).
When used with OGRCoordinateTransformation,
mapping=[2,1] typically expresses the inversion of axis between the data axis and the CRS axis for a 2D CRS.
Automatically implies SetAxisMappingStrategy(OAMS_CUSTOM)
This is the same as the C function OSRSetDataAxisToSRSAxisMapping().
mapping | The new data axis to CRS axis mapping. |
OGRErr OGRSpatialReference::SetExtension | ( | const char * | pszTargetKey, |
const char * | pszName, | ||
const char * | pszValue | ||
) |
Set extension value.
Set the value of the named EXTENSION item for the identified target node.
pszTargetKey | the name or path to the parent node of the EXTENSION. |
pszName | the name of the extension being fetched. |
pszValue | the value to set |
OGRErr OGRSpatialReference::SetFromUserInput | ( | const char * | pszDefinition | ) |
Set spatial reference from various text formats.
This method will examine the provided input, and try to deduce the format, and then use it to initialize the spatial reference system. It may take the following forms:
It is expected that this method will be extended in the future to support XML and perhaps a simplified "minilanguage" for indicating common UTM and State Plane definitions.
This method is intended to be flexible, but by its nature it is imprecise as it must guess information about the format intended. When possible applications should call the specific method appropriate if the input is known to be in a particular format.
This method does the same thing as the OSRSetFromUserInput() function.
pszDefinition | text definition to try to deduce SRS from. |
OGRErr OGRSpatialReference::SetFromUserInput | ( | const char * | pszDefinition, |
CSLConstList | papszOptions | ||
) |
Set spatial reference from various text formats.
This method will examine the provided input, and try to deduce the format, and then use it to initialize the spatial reference system. It may take the following forms:
It is expected that this method will be extended in the future to support XML and perhaps a simplified "minilanguage" for indicating common UTM and State Plane definitions.
This method is intended to be flexible, but by its nature it is imprecise as it must guess information about the format intended. When possible applications should call the specific method appropriate if the input is known to be in a particular format.
This method does the same thing as the OSRSetFromUserInput() and OSRSetFromUserInputEx() functions.
pszDefinition | text definition to try to deduce SRS from. |
papszOptions | NULL terminated list of options, or NULL. |
OGRErr OGRSpatialReference::SetGeocCS | ( | const char * | pszName | ) |
Set the user visible GEOCCS name.
This method is the same as the C function OSRSetGeocCS().
This method will ensure a GEOCCS node is created as the root, and set the provided name on it. If used on a GEOGCS coordinate system, the DATUM and PRIMEM nodes from the GEOGCS will be transferred over to the GEOGCS.
pszName | the user visible name to assign. Not used as a key. |
OGRErr OGRSpatialReference::SetGeogCS | ( | const char * | pszGeogName, |
const char * | pszDatumName, | ||
const char * | pszSpheroidName, | ||
double | dfSemiMajor, | ||
double | dfInvFlattening, | ||
const char * | pszPMName = nullptr , |
||
double | dfPMOffset = 0.0 , |
||
const char * | pszAngularUnits = nullptr , |
||
double | dfConvertToRadians = 0.0 |
||
) |
Set geographic coordinate system.
This method is used to set the datum, ellipsoid, prime meridian and angular units for a geographic coordinate system. It can be used on its own to establish a geographic spatial reference, or applied to a projected coordinate system to establish the underlying geographic coordinate system.
This method does the same as the C function OSRSetGeogCS().
pszGeogName | user visible name for the geographic coordinate system (not to serve as a key). |
pszDatumName | key name for this datum. The OpenGIS specification lists some known values, and otherwise EPSG datum names with a standard transformation are considered legal keys. |
pszSpheroidName | user visible spheroid name (not to serve as a key) |
dfSemiMajor | the semi major axis of the spheroid. |
dfInvFlattening | the inverse flattening for the spheroid. This can be computed from the semi minor axis as 1/f = 1.0 / (1.0 - semiminor/semimajor). |
pszPMName | the name of the prime meridian (not to serve as a key) If this is NULL a default value of "Greenwich" will be used. |
dfPMOffset | the longitude of Greenwich relative to this prime meridian. Always in Degrees |
pszAngularUnits | the angular units name (see ogr_srs_api.h for some standard names). If NULL a value of "degrees" will be assumed. |
dfConvertToRadians | value to multiply angular units by to transform them to radians. A value of SRS_UA_DEGREE_CONV will be used if pszAngularUnits is NULL. |
OGRErr OGRSpatialReference::SetHOM | ( | double | dfCenterLat, |
double | dfCenterLong, | ||
double | dfAzimuth, | ||
double | dfRectToSkew, | ||
double | dfScale, | ||
double | dfFalseEasting, | ||
double | dfFalseNorthing | ||
) |
Hotine Oblique Mercator.
Set a Hotine Oblique Mercator projection using azimuth angle.
This projection corresponds to EPSG projection method 9812, also sometimes known as hotine oblique mercator (variant A)..
This method does the same thing as the C function OSRSetHOM().
dfCenterLat | Latitude of the projection origin. |
dfCenterLong | Longitude of the projection origin. |
dfAzimuth | Azimuth, measured clockwise from North, of the projection centerline. |
dfRectToSkew | Angle from Rectified to Skew Grid |
dfScale | Scale factor applies to the projection origin. |
dfFalseEasting | False easting. |
dfFalseNorthing | False northing. |
OGRErr OGRSpatialReference::SetHOM2PNO | ( | double | dfCenterLat, |
double | dfLat1, | ||
double | dfLong1, | ||
double | dfLat2, | ||
double | dfLong2, | ||
double | dfScale, | ||
double | dfFalseEasting, | ||
double | dfFalseNorthing | ||
) |
Hotine Oblique Mercator 2 points.
Set a Hotine Oblique Mercator projection using two points on projection centerline.
This method does the same thing as the C function OSRSetHOM2PNO().
dfCenterLat | Latitude of the projection origin. |
dfLat1 | Latitude of the first point on center line. |
dfLong1 | Longitude of the first point on center line. |
dfLat2 | Latitude of the second point on center line. |
dfLong2 | Longitude of the second point on center line. |
dfScale | Scale factor applies to the projection origin. |
dfFalseEasting | False easting. |
dfFalseNorthing | False northing. |
OGRErr OGRSpatialReference::SetHOMAC | ( | double | dfCenterLat, |
double | dfCenterLong, | ||
double | dfAzimuth, | ||
double | dfRectToSkew, | ||
double | dfScale, | ||
double | dfFalseEasting, | ||
double | dfFalseNorthing | ||
) |
Hotine Oblique Mercator Azimuth Center / Variant B.
Set an Hotine Oblique Mercator Azimuth Center projection using azimuth angle.
This projection corresponds to EPSG projection method 9815, also sometimes known as hotine oblique mercator (variant B).
This method does the same thing as the C function OSRSetHOMAC().
dfCenterLat | Latitude of the projection origin. |
dfCenterLong | Longitude of the projection origin. |
dfAzimuth | Azimuth, measured clockwise from North, of the projection centerline. |
dfRectToSkew | Angle from Rectified to Skew Grid |
dfScale | Scale factor applies to the projection origin. |
dfFalseEasting | False easting. |
dfFalseNorthing | False northing. |
OGRErr OGRSpatialReference::SetKrovak | ( | double | dfCenterLat, |
double | dfCenterLong, | ||
double | dfAzimuth, | ||
double | dfPseudoStdParallel1, | ||
double | dfScale, | ||
double | dfFalseEasting, | ||
double | dfFalseNorthing | ||
) |
Krovak Oblique Conic Conformal.
Krovak east-north projection.
Note that dfAzimuth and dfPseudoStdParallel1 are ignored when exporting to PROJ and should be respectively set to 30.28813972222222 and 78.5
OGRErr OGRSpatialReference::SetLinearUnits | ( | const char * | pszUnitsName, |
double | dfInMeters | ||
) |
Set the linear units for the projection.
This method creates a UNIT subnode with the specified values as a child of the PROJCS, GEOCCS, GEOGCS or LOCAL_CS node. When called on a Geographic 3D CRS the vertical axis units will be set.
This method does the same as the C function OSRSetLinearUnits().
pszUnitsName | the units name to be used. Some preferred units names can be found in ogr_srs_api.h such as SRS_UL_METER, SRS_UL_FOOT and SRS_UL_US_FOOT. |
dfInMeters | the value to multiple by a length in the indicated units to transform to meters. Some standard conversion factors can be found in ogr_srs_api.h. |
OGRErr OGRSpatialReference::SetLinearUnitsAndUpdateParameters | ( | const char * | pszName, |
double | dfInMeters, | ||
const char * | pszUnitAuthority = nullptr , |
||
const char * | pszUnitCode = nullptr |
||
) |
Set the linear units for the projection.
This method creates a UNIT subnode with the specified values as a child of the PROJCS or LOCAL_CS node. It works the same as the SetLinearUnits() method, but it also updates all existing linear projection parameter values from the old units to the new units.
pszName | the units name to be used. Some preferred units names can be found in ogr_srs_api.h such as SRS_UL_METER, SRS_UL_FOOT and SRS_UL_US_FOOT. |
dfInMeters | the value to multiple by a length in the indicated units to transform to meters. Some standard conversion factors can be found in ogr_srs_api.h. |
pszUnitAuthority | Unit authority name. Or nullptr |
pszUnitCode | Unit code. Or nullptr |
OGRErr OGRSpatialReference::SetLocalCS | ( | const char * | pszName | ) |
Set the user visible LOCAL_CS name.
This method is the same as the C function OSRSetLocalCS().
This method will ensure a LOCAL_CS node is created as the root, and set the provided name on it. It must be used before SetLinearUnits().
pszName | the user visible name to assign. Not used as a key. |
OGRErr OGRSpatialReference::SetLOM | ( | double | dfCenterLat, |
double | dfCenterLong, | ||
double | dfAzimuth, | ||
double | dfScale, | ||
double | dfFalseEasting, | ||
double | dfFalseNorthing | ||
) |
Laborde Oblique Mercator.
Set a Laborde Oblique Mercator projection.
dfCenterLat | Latitude of the projection origin. |
dfCenterLong | Longitude of the projection origin. |
dfAzimuth | Azimuth, measured clockwise from North, of the projection centerline. |
dfScale | Scale factor on the initiali line |
dfFalseEasting | False easting. |
dfFalseNorthing | False northing. |
OGRErr OGRSpatialReference::SetNode | ( | const char * | pszNodePath, |
const char * | pszNewNodeValue | ||
) |
Set attribute value in spatial reference.
Missing intermediate nodes in the path will be created if not already in existence. If the attribute has no children one will be created and assigned the value otherwise the zeroth child will be assigned the value.
This method does the same as the C function OSRSetAttrValue().
pszNodePath | full path to attribute to be set. For instance "PROJCS|GEOGCS|UNIT". |
pszNewNodeValue | value to be assigned to node, such as "meter". This may be NULL if you just want to force creation of the intermediate path. |
OGRErr OGRSpatialReference::SetNode | ( | const char * | pszNodePath, |
double | dfValue | ||
) |
Set attribute value in spatial reference.
Missing intermediate nodes in the path will be created if not already in existence. If the attribute has no children one will be created and assigned the value otherwise the zeroth child will be assigned the value.
This method does the same as the C function OSRSetAttrValue().
pszNodePath | full path to attribute to be set. For instance "PROJCS|GEOGCS|UNIT". |
dfValue | value to be assigned to node. |
OGRErr OGRSpatialReference::SetNormProjParm | ( | const char * | pszName, |
double | dfValue | ||
) |
Set a projection parameter with a normalized value.
This method is the same as SetProjParm() except that the value of the parameter passed in is assumed to be in "normalized" form (decimal degrees for angular values, meters for linear values. The values are converted in a form suitable for the GEOGCS and linear units in effect.
This method is the same as the C function OSRSetNormProjParm().
pszName | the parameter name, which should be selected from the macros in ogr_srs_api.h, such as SRS_PP_CENTRAL_MERIDIAN. |
dfValue | value to assign. |
OGRErr OGRSpatialReference::SetProjCS | ( | const char * | pszName | ) |
Set the user visible PROJCS name.
This method is the same as the C function OSRSetProjCS().
This method will ensure a PROJCS node is created as the root, and set the provided name on it. If used on a GEOGCS coordinate system, the GEOGCS node will be demoted to be a child of the new PROJCS root.
pszName | the user visible name to assign. Not used as a key. |
OGRErr OGRSpatialReference::SetProjection | ( | const char * | pszProjection | ) |
Set a projection name.
This method is the same as the C function OSRSetProjection().
pszProjection | the projection name, which should be selected from the macros in ogr_srs_api.h, such as SRS_PT_TRANSVERSE_MERCATOR. |
OGRErr OGRSpatialReference::SetProjParm | ( | const char * | pszParamName, |
double | dfValue | ||
) |
Set a projection parameter value.
Adds a new PARAMETER under the PROJCS with the indicated name and value.
This method is the same as the C function OSRSetProjParm().
Please check https://gdal.org/proj_list pages for legal parameter names for specific projections.
pszParamName | the parameter name, which should be selected from the macros in ogr_srs_api.h, such as SRS_PP_CENTRAL_MERIDIAN. |
dfValue | value to assign. |
OGRErr OGRSpatialReference::SetPS | ( | double | dfCenterLat, |
double | dfCenterLong, | ||
double | dfScale, | ||
double | dfFalseEasting, | ||
double | dfFalseNorthing | ||
) |
Polar Stereographic.
Sets a Polar Stereographic projection.
Two variants are possible:
void OGRSpatialReference::SetRoot | ( | OGR_SRSNode * | poNewRoot | ) |
Set the root SRS node.
If the object has an existing tree of OGR_SRSNodes, they are destroyed as part of assigning the new root. Ownership of the passed OGR_SRSNode is is assumed by the OGRSpatialReference.
poNewRoot | object to assign as root. |
OGRErr OGRSpatialReference::SetStatePlane | ( | int | nZone, |
int | bNAD83 = TRUE , |
||
const char * | pszOverrideUnitName = nullptr , |
||
double | dfOverrideUnit = 0.0 |
||
) |
State Plane.
Set State Plane projection definition.
This will attempt to generate a complete definition of a state plane zone based on generating the entire SRS from the EPSG tables. If the EPSG tables are unavailable, it will produce a stubbed LOCAL_CS definition and return OGRERR_FAILURE.
This method is the same as the C function OSRSetStatePlaneWithUnits().
nZone | State plane zone number, in the USGS numbering scheme (as distinct from the Arc/Info and Erdas numbering scheme. |
bNAD83 | TRUE if the NAD83 zone definition should be used or FALSE if the NAD27 zone definition should be used. |
pszOverrideUnitName | Linear unit name to apply overriding the legal definition for this zone. |
dfOverrideUnit | Linear unit conversion factor to apply overriding the legal definition for this zone. |
OGRErr OGRSpatialReference::SetTargetLinearUnits | ( | const char * | pszTargetKey, |
const char * | pszUnitsName, | ||
double | dfInMeters, | ||
const char * | pszUnitAuthority = nullptr , |
||
const char * | pszUnitCode = nullptr |
||
) |
Set the linear units for the projection.
This method creates a UNIT subnode with the specified values as a child of the target node.
This method does the same as the C function OSRSetTargetLinearUnits().
pszTargetKey | the keyword to set the linear units for. i.e. "PROJCS" or "VERT_CS" |
pszUnitsName | the units name to be used. Some preferred units names can be found in ogr_srs_api.h such as SRS_UL_METER, SRS_UL_FOOT and SRS_UL_US_FOOT. |
dfInMeters | the value to multiple by a length in the indicated units to transform to meters. Some standard conversion factors can be found in ogr_srs_api.h. |
pszUnitAuthority | Unit authority name. Or nullptr |
pszUnitCode | Unit code. Or nullptr |
OGRErr OGRSpatialReference::SetTOWGS84 | ( | double | dfDX, |
double | dfDY, | ||
double | dfDZ, | ||
double | dfEX = 0.0 , |
||
double | dfEY = 0.0 , |
||
double | dfEZ = 0.0 , |
||
double | dfPPM = 0.0 |
||
) |
Set the Bursa-Wolf conversion to WGS84.
This will create the TOWGS84 node as a child of the DATUM. It will fail if there is no existing DATUM node. It will replace an existing TOWGS84 node if there is one.
The parameters have the same meaning as EPSG transformation 9606 (Position Vector 7-param. transformation).
This method is the same as the C function OSRSetTOWGS84().
dfDX | X child in meters. |
dfDY | Y child in meters. |
dfDZ | Z child in meters. |
dfEX | X rotation in arc seconds (optional, defaults to zero). |
dfEY | Y rotation in arc seconds (optional, defaults to zero). |
dfEZ | Z rotation in arc seconds (optional, defaults to zero). |
dfPPM | scaling factor (parts per million). |
OGRErr OGRSpatialReference::SetUTM | ( | int | nZone, |
int | bNorth = TRUE |
||
) |
Universal Transverse Mercator.
Set UTM projection definition.
This will generate a projection definition with the full set of transverse mercator projection parameters for the given UTM zone. If no PROJCS[] description is set yet, one will be set to look like "UTM Zone %d, {Northern, Southern} Hemisphere".
This method is the same as the C function OSRSetUTM().
nZone | UTM zone. |
bNorth | TRUE for northern hemisphere, or FALSE for southern hemisphere. |
OGRErr OGRSpatialReference::SetVertCS | ( | const char * | pszVertCSName, |
const char * | pszVertDatumName, | ||
int | nVertDatumType = 2005 |
||
) |
Set the user visible VERT_CS name.
This method is the same as the C function OSRSetVertCS().
This method will ensure a VERT_CS node is created if needed. If the existing coordinate system is GEOGCS or PROJCS rooted, then it will be turned into a COMPD_CS.
pszVertCSName | the user visible name of the vertical coordinate system. Not used as a key. |
pszVertDatumName | the user visible name of the vertical datum. It is helpful if this matches the EPSG name. |
nVertDatumType | the OGC vertical datum type. Ignored |
OGRErr OGRSpatialReference::SetWellKnownGeogCS | ( | const char * | pszName | ) |
Set a GeogCS based on well known name.
This may be called on an empty OGRSpatialReference to make a geographic coordinate system, or on something with an existing PROJCS node to set the underlying geographic coordinate system of a projected coordinate system.
The following well known text values are currently supported, Except for "EPSG:n", the others are without dependency on EPSG data files:
pszName | name of well known geographic coordinate system. |
bool OGRSpatialReference::StripTOWGS84IfKnownDatum | ( | ) |
Remove TOWGS84 information if the CRS has a known horizontal datum.
bool OGRSpatialReference::StripTOWGS84IfKnownDatumAndAllowed | ( | ) |
Remove TOWGS84 information if the CRS has a known horizontal datum and this is allowed by the user.
The default behavior is to remove TOWGS84 information if the CRS has a known horizontal datum. This can be disabled by setting the OSR_STRIP_TOWGS84 configuration option to NO.
OGRErr OGRSpatialReference::StripVertical | ( | ) |
Convert a compound cs into a horizontal CS.
If this SRS is of type COMPD_CS[] then the vertical CS and the root COMPD_CS nodes are stripped resulting and only the horizontal coordinate system portion remains (normally PROJCS, GEOGCS or LOCAL_CS).
If this is not a compound coordinate system then nothing is changed.
This method is the same as the C function OSRStripVertical().
|
inlinestatic |
Convert a OGRSpatialReference* to a OGRSpatialReferenceH.
OGRErr OGRSpatialReference::Validate | ( | ) | const |
Validate CRS imported with importFromWkt() or with modified with direct node manipulations.
Otherwise the CRS should be always valid.
This method attempts to verify that the spatial reference system is well formed, and consists of known tokens. The validation is not comprehensive.
This method is the same as the C function OSRValidate().
|
static |
Limitations for OGRSpatialReference::SetFromUserInput().
Currently ALLOW_NETWORK_ACCESS=NO and ALLOW_FILE_ACCESS=NO.