OSSIM - Open Source Software Image Map  Version 1.9.0 (20180803)
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ossimVanDerGrintenProjection Class Reference

#include <ossimVanDerGrintenProjection.h>

Inheritance diagram for ossimVanDerGrintenProjection:
ossimMapProjection ossimProjection ossimObject ossimErrorStatusInterface ossimReferenced

Public Member Functions

 ossimVanDerGrintenProjection (const ossimEllipsoid &ellipsoid=ossimEllipsoid(), const ossimGpt &origin=ossimGpt())
 
 ossimVanDerGrintenProjection (const ossimEllipsoid &ellipsoid, const ossimGpt &origin, const double falseEasting, const double falseNorthing)
 
 ~ossimVanDerGrintenProjection ()
 
virtual ossimObjectdup () const
 
virtual ossimGpt inverse (const ossimDpt &eastingNorthing) const
 Will take a point in meters and convert it to ground. More...
 
virtual ossimDpt forward (const ossimGpt &latLon) const
 All map projections will convert the world coordinate to an easting northing (Meters). More...
 
virtual void update ()
 
void setFalseEasting (double falseEasting)
 
void setFalseNorthing (double falseNorthing)
 
void setFalseEastingNorthing (double falseEasting, double falseNorthing)
 
void setDefaults ()
 
void setCentralMeridian (double centralMeridian)
 
double getFalseEasting () const
 
double getFalseNorthing () const
 
virtual bool saveState (ossimKeywordlist &kwl, const char *prefix=0) const
 
virtual bool loadState (const ossimKeywordlist &kwl, const char *prefix=0)
 
virtual bool operator== (const ossimProjection &projection) const
 Returns TRUE if principal parameters are within epsilon tolerance. More...
 
- Public Member Functions inherited from ossimMapProjection
 ossimMapProjection (const ossimEllipsoid &ellipsoid=ossimEllipsoid(), const ossimGpt &origin=ossimGpt())
 
 ossimMapProjection (const ossimMapProjection &src)
 
virtual ossimGpt origin () const
 
virtual ossimDpt worldToLineSample (const ossimGpt &worldPoint) const
 
virtual void worldToLineSample (const ossimGpt &worldPoint, ossimDpt &lineSample) const
 
virtual ossimGpt lineSampleToWorld (const ossimDpt &projectedPoint) const
 
virtual void lineSampleToWorld (const ossimDpt &projectedPoint, ossimGpt &gpt) const
 
virtual void lineSampleHeightToWorld (const ossimDpt &lineSampPt, const double &heightAboveEllipsoid, ossimGpt &worldPt) const
 This is the pure virtual that projects the image point to the given elevation above ellipsoid, thereby bypassing reference to a DEM. More...
 
virtual void lineSampleToEastingNorthing (const ossimDpt &liineSample, ossimDpt &eastingNorthing) const
 
virtual void eastingNorthingToLineSample (const ossimDpt &eastingNorthing, ossimDpt &lineSample) const
 
virtual void eastingNorthingToWorld (const ossimDpt &eastingNorthing, ossimGpt &worldPt) const
 
virtual double getStandardParallel1 () const
 Derived classes should implement as needed. More...
 
virtual double getStandardParallel2 () const
 Derived classes should implement as needed. More...
 
virtual void setPcsCode (ossim_uint32 pcsCode)
 
virtual ossim_uint32 getPcsCode () const
 Returns the EPSG PCS code or 32767 if the projection is a custom (non-EPSG) projection. More...
 
virtual ossimString getProjectionName () const
 Returns the projection name. More...
 
virtual double getA () const
 ACCESS METHODS: More...
 
virtual double getB () const
 
virtual double getF () const
 
virtual ossimDpt getMetersPerPixel () const
 
virtual const ossimDptgetDecimalDegreesPerPixel () const
 Returns decimal degrees per pixel as an ossimDpt with "x" representing longitude and "y" representing latitude. More...
 
virtual const ossimDptgetUlEastingNorthing () const
 
virtual const ossimGptgetUlGpt () const
 
virtual const ossimDatumgetDatum () const
 
const ossimEllipsoidgetEllipsoid () const
 
const ossimGptgetOrigin () const
 
virtual bool isGeographic () const
 
virtual void applyScale (const ossimDpt &scale, bool recenterTiePoint)
 Applies scale to theDeltaLonPerPixel, theDeltaLatPerPixel and theMetersPerPixel data members (eg: theDeltaLonPerPixel *= scale.x). More...
 
virtual void setEllipsoid (const ossimEllipsoid &ellipsoid)
 SET METHODS: More...
 
virtual void setAB (double a, double b)
 
virtual void setDatum (const ossimDatum *datum)
 Sets theDatum to datum. More...
 
virtual void setOrigin (const ossimGpt &origin)
 Sets theOrigin to origin. More...
 
virtual void setMetersPerPixel (const ossimDpt &gsd)
 
virtual void setDecimalDegreesPerPixel (const ossimDpt &gsd)
 
virtual void setUlTiePoints (const ossimGpt &gpt)
 
virtual void setUlTiePoints (const ossimDpt &eastingNorthing)
 
virtual void setUlEastingNorthing (const ossimDpt &ulEastingNorthing)
 
virtual void setUlGpt (const ossimGpt &ulGpt)
 
virtual void assign (const ossimProjection &aProjection)
 
virtual std::ostream & print (std::ostream &out) const
 Prints data members to stream. More...
 
virtual void computeDegreesPerPixel ()
 Computes the approximate resolution in degrees/pixel. More...
 
virtual void computeMetersPerPixel ()
 This will go from the ground point and give you an approximate meters per pixel. More...
 
void setMatrix (double rotation, const ossimDpt &scale, const ossimDpt &translation)
 
void setMatrixScale (const ossimDpt &scale)
 
void setMatrixRotation (double rotation)
 
void setMatrixTranslation (const ossimDpt &translation)
 
void snapTiePointTo (ossim_float64 multiple, ossimUnitType unitType)
 Utility method to snap the tie point to some multiple. More...
 
void snapTiePointToOrigin ()
 
void setElevationLookupFlag (bool flag)
 
bool getElevationLookupFlag () const
 
ossimUnitType getModelTransformUnitType () const
 
void setModelTransformUnitType (ossimUnitType unit)
 
bool hasModelTransform () const
 
virtual bool isAffectedByElevation () const
 Implementation of pure virtual ossimProjection::isAffectedByElevation method. More...
 
void setProjectionUnits (ossimUnitType units)
 
ossimUnitType getProjectionUnits () const
 OSSIM considers all map projection coordinates (including false eastings and northings) to be in meters. More...
 
virtual bool isEqualTo (const ossimObject &obj, ossimCompareType compareType=OSSIM_COMPARE_FULL) const
 
- Public Member Functions inherited from ossimProjection
 ossimProjection ()
 
virtual ~ossimProjection ()
 
virtual void getRoundTripError (const ossimDpt &imagePoint, ossimDpt &errorResult) const
 
virtual void getRoundTripError (const ossimGpt &groundPoint, ossimDpt &errorResult) const
 
virtual void getGroundClipPoints (ossimGeoPolygon &gpts) const
 
virtual bool operator!= (const ossimProjection &projection) const
 
- Public Member Functions inherited from ossimObject
 ossimObject ()
 
virtual ~ossimObject ()
 
virtual ossimString getShortName () const
 
virtual ossimString getLongName () const
 
virtual ossimString getDescription () const
 
virtual ossimString getClassName () const
 
virtual RTTItypeid getType () const
 
virtual bool canCastTo (ossimObject *obj) const
 
virtual bool canCastTo (const RTTItypeid &id) const
 
virtual bool canCastTo (const ossimString &parentClassName) const
 
virtual void accept (ossimVisitor &visitor)
 
- Public Member Functions inherited from ossimReferenced
 ossimReferenced ()
 
 ossimReferenced (const ossimReferenced &)
 
ossimReferencedoperator= (const ossimReferenced &)
 
void ref () const
 increment the reference count by one, indicating that this object has another pointer which is referencing it. More...
 
void unref () const
 decrement the reference count by one, indicating that a pointer to this object is referencing it. More...
 
void unref_nodelete () const
 decrement the reference count by one, indicating that a pointer to this object is referencing it. More...
 
int referenceCount () const
 
- Public Member Functions inherited from ossimErrorStatusInterface
 ossimErrorStatusInterface ()
 
virtual ~ossimErrorStatusInterface ()
 
virtual ossimErrorCode getErrorStatus () const
 
virtual ossimString getErrorStatusString () const
 
virtual void setErrorStatus (ossimErrorCode error_status) const
 
virtual void setErrorStatus () const
 
virtual void clearErrorStatus () const
 
bool hasError () const
 

Protected Member Functions

long Set_Van_der_Grinten_Parameters (double a, double f, double Central_Meridian, double False_Easting, double False_Northing)
 
void Get_Van_der_Grinten_Parameters (double *a, double *f, double *Central_Meridian, double *False_Easting, double *False_Northing) const
 
long Convert_Geodetic_To_Van_der_Grinten (double Latitude, double Longitude, double *Easting, double *Northing) const
 
long Convert_Van_der_Grinten_To_Geodetic (double Easting, double Northing, double *Latitude, double *Longitude) const
 
- Protected Member Functions inherited from ossimMapProjection
virtual ~ossimMapProjection ()
 
void updateFromTransform ()
 
- Protected Member Functions inherited from ossimReferenced
virtual ~ossimReferenced ()
 

Protected Attributes

double Grin_a
 
double Grin_f
 
double es2
 
double es4
 
double es6
 
double Ra
 
double PI_Ra
 
double Grin_Origin_Long
 
double Grin_False_Easting
 
double Grin_False_Northing
 
- Protected Attributes inherited from ossimMapProjection
ossimEllipsoid theEllipsoid
 This method verifies that the projection parameters match the current pcs code. More...
 
ossimGpt theOrigin
 
const ossimDatumtheDatum
 This is only set if we want to have built in datum shifting. More...
 
ossimDpt theMetersPerPixel
 Holds the number of meters per pixel. More...
 
ossimDpt theDegreesPerPixel
 Hold the decimal degrees per pixel. More...
 
ossimGpt theUlGpt
 Hold tie point in decimal degrees. More...
 
ossimDpt theUlEastingNorthing
 Hold tie point as easting northing. More...
 
ossimDpt theFalseEastingNorthing
 Hold the false easting northing. More...
 
ossim_uint32 thePcsCode
 Projection Coordinate System(PCS) code. More...
 
bool theElevationLookupFlag
 
ossimMatrix4x4 theModelTransform
 
ossimMatrix4x4 theInverseModelTransform
 
ossimUnitType theModelTransformUnitType
 
ossimUnitType theProjectionUnits
 Linear units of the projection as indicated in the projection's specification: More...
 
- Protected Attributes inherited from ossimErrorStatusInterface
ossimErrorCode theErrorStatus
 

Detailed Description

Definition at line 17 of file ossimVanDerGrintenProjection.h.

Constructor & Destructor Documentation

◆ ossimVanDerGrintenProjection() [1/2]

ossimVanDerGrintenProjection::ossimVanDerGrintenProjection ( const ossimEllipsoid ellipsoid = ossimEllipsoid(),
const ossimGpt origin = ossimGpt() 
)

Definition at line 54 of file ossimVanDerGrintenProjection.cpp.

References setDefaults(), and update().

56  :ossimMapProjection(ellipsoid, origin)
57 {
58  setDefaults();
59  update();
60 }
virtual ossimGpt origin() const
ossimMapProjection(const ossimEllipsoid &ellipsoid=ossimEllipsoid(), const ossimGpt &origin=ossimGpt())

◆ ossimVanDerGrintenProjection() [2/2]

ossimVanDerGrintenProjection::ossimVanDerGrintenProjection ( const ossimEllipsoid ellipsoid,
const ossimGpt origin,
const double  falseEasting,
const double  falseNorthing 
)

Definition at line 62 of file ossimVanDerGrintenProjection.cpp.

References Grin_False_Easting, Grin_False_Northing, and update().

66  :ossimMapProjection(ellipsoid, origin)
67 {
68  Grin_False_Easting = falseEasting;
69  Grin_False_Northing = falseNorthing;
70 
71  update();
72 }
virtual ossimGpt origin() const
ossimMapProjection(const ossimEllipsoid &ellipsoid=ossimEllipsoid(), const ossimGpt &origin=ossimGpt())

◆ ~ossimVanDerGrintenProjection()

ossimVanDerGrintenProjection::~ossimVanDerGrintenProjection ( )
inline

Definition at line 26 of file ossimVanDerGrintenProjection.h.

26 {}

Member Function Documentation

◆ Convert_Geodetic_To_Van_der_Grinten()

long ossimVanDerGrintenProjection::Convert_Geodetic_To_Van_der_Grinten ( double  Latitude,
double  Longitude,
double *  Easting,
double *  Northing 
) const
protected

The function Convert_Geodetic_To_Van_der_Grinten converts geodetic (latitude and longitude) coordinates to Van Der Grinten projection easting, and northing coordinates, according to the current ellipsoid and Van Der Grinten projection parameters. If any errors occur, the error code(s) are returned by the function, otherwise GRIN_NO_ERROR is returned.

Latitude : Latitude (phi) in radians (input) Longitude : Longitude (lambda) in radians (input) Easting : Easting (X) in meters (output) Northing : Northing (Y) in meters (output)

Definition at line 278 of file ossimVanDerGrintenProjection.cpp.

References FLOAT_EQ, Grin_False_Easting, Grin_False_Northing, GRIN_NO_ERROR, Grin_Origin_Long, M_PI, MAX_LAT, PI_Ra, Ra, and TWO_OVER_PI.

Referenced by forward().

283 { /* BEGIN Convert_Geodetic_To_Van_der_Grinten */
284 /*
285  * The function Convert_Geodetic_To_Van_der_Grinten converts geodetic (latitude and
286  * longitude) coordinates to Van Der Grinten projection (easting and northing)
287  * coordinates, according to the current ellipsoid and Van Der Grinten projection
288  * parameters. If any errors occur, the error code(s) are returned by the
289  * function, otherwise GRIN_NO_ERROR is returned.
290  *
291  * Latitude : Latitude (phi) in radians (input)
292  * Longitude : Longitude (lambda) in radians (input)
293  * Easting : Easting (X) in meters (output)
294  * Northing : Northing (Y) in meters (output)
295  */
296 
297  double dlam; /* Longitude - Central Meridan */
298  double aa, aasqr;
299  double gg;
300  double pp, ppsqr;
301  double gg_MINUS_ppsqr, ppsqr_PLUS_aasqr;
302  double in_theta;
303  double theta;
304  double sin_theta, cos_theta;
305  double qq;
306  long Error_Code = GRIN_NO_ERROR;
307 
308 // if ((Latitude < -PI_OVER_2) || (Latitude > PI_OVER_2))
309 // { /* Latitude out of range */
310 // Error_Code |= GRIN_LAT_ERROR;
311 // }
312 // if ((Longitude < -M_PI) || (Longitude > TWO_PI))
313 // { /* Longitude out of range */
314 // Error_Code|= GRIN_LON_ERROR;
315 // }
316 
317  if (!Error_Code)
318  { /* no errors */
319 
320  dlam = Longitude - Grin_Origin_Long;
321 // if (dlam > M_PI)
322 // {
323 // dlam -= TWO_PI;
324 // }
325 // if (dlam < -M_PI)
326 // {
327 // dlam += TWO_PI;
328 // }
329 
330  if (Latitude == 0.0)
331  {
332  *Easting = Ra * dlam + Grin_False_Easting;
333  *Northing = 0.0;
334  }
335  else if (dlam == 0.0 || FLOAT_EQ(Latitude,MAX_LAT,.00001) || FLOAT_EQ(Latitude,-MAX_LAT,.00001))
336  {
337  in_theta = fabs(TWO_OVER_PI * Latitude);
338 
339  if (in_theta > 1.0)
340  in_theta = 1.0;
341  else if (in_theta < -1.0)
342  in_theta = -1.0;
343 
344  theta = asin(in_theta);
345  *Easting = 0.0;
346  *Northing = PI_Ra * tan(theta / 2) + Grin_False_Northing;
347  if (Latitude < 0.0)
348  *Northing *= -1.0;
349  }
350  else
351  {
352  aa = 0.5 * fabs(M_PI / dlam - dlam / M_PI);
353  in_theta = fabs(TWO_OVER_PI * Latitude);
354 
355  if (in_theta > 1.0)
356  in_theta = 1.0;
357  else if (in_theta < -1.0)
358  in_theta = -1.0;
359 
360  theta = asin(in_theta);
361  sin_theta = sin(theta);
362  cos_theta = cos(theta);
363  gg = cos_theta / (sin_theta + cos_theta - 1);
364  pp = gg * (2 / sin_theta - 1);
365  aasqr = aa * aa;
366  ppsqr = pp * pp;
367  gg_MINUS_ppsqr = gg - ppsqr;
368  ppsqr_PLUS_aasqr = ppsqr + aasqr;
369  qq = aasqr + gg;
370  *Easting = PI_Ra * (aa * (gg_MINUS_ppsqr) +
371  sqrt(aasqr * (gg_MINUS_ppsqr) * (gg_MINUS_ppsqr) -
372  (ppsqr_PLUS_aasqr) * (gg * gg - ppsqr))) /
373  (ppsqr_PLUS_aasqr) + Grin_False_Easting;
374  if (dlam < 0.0)
375  *Easting *= -1.0;
376  *Northing = PI_Ra * (pp * qq - aa * sqrt ((aasqr + 1) * (ppsqr_PLUS_aasqr) - qq * qq)) /
377  (ppsqr_PLUS_aasqr) + Grin_False_Northing;
378  if (Latitude < 0.0)
379  *Northing *= -1.0;
380  }
381  }
382  return (Error_Code);
383 
384 } /* END OF Convert_Geodetic_To_Van_der_Grinten */
#define M_PI
const double TWO_OVER_PI
#define GRIN_NO_ERROR
#define FLOAT_EQ(x, v, epsilon)

◆ Convert_Van_der_Grinten_To_Geodetic()

long ossimVanDerGrintenProjection::Convert_Van_der_Grinten_To_Geodetic ( double  Easting,
double  Northing,
double *  Latitude,
double *  Longitude 
) const
protected

The function Convert_Van_der_Grinten_To_Geodetic converts Van Der Grinten projection easting and northing coordinates to geodetic (latitude and longitude) coordinates, according to the current ellipsoid and Van Der Grinten projection coordinates. If any errors occur, the error code(s) are returned by the function, otherwise GRIN_NO_ERROR is returned.

Easting : Easting (X) in meters (input) Northing : Northing (Y) in meters (input) Latitude : Latitude (phi) in radians (output) Longitude : Longitude (lambda) in radians (output)

Definition at line 387 of file ossimVanDerGrintenProjection.cpp.

References Grin_False_Easting, Grin_False_Northing, GRIN_NO_ERROR, Grin_Origin_Long, M_PI, MAX_LAT, ONE_THIRD, PI_OVER_3, and PI_Ra.

Referenced by inverse().

391 { /* BEGIN Convert_Van_der_Grinten_To_Geodetic */
392 /*
393  * The function Convert_Van_der_Grinten_To_Geodetic converts Grinten projection
394  * (easting and northing) coordinates to geodetic (latitude and longitude)
395  * coordinates, according to the current ellipsoid and Grinten projection
396  * coordinates. If any errors occur, the error code(s) are returned by the
397  * function, otherwise GRIN_NO_ERROR is returned.
398  *
399  * Easting : Easting (X) in meters (input)
400  * Northing : Northing (Y) in meters (input)
401  * Latitude : Latitude (phi) in radians (output)
402  * Longitude : Longitude (lambda) in radians (output)
403  */
404 
405  double dx, dy;
406  double xx, xxsqr;
407  double yy, yysqr, two_yysqr;
408  double xxsqr_PLUS_yysqr;
409  double c1;
410  double c2;
411  double c3, c3sqr;
412  double c2_OVER_3c3;
413  double dd;
414  double a1;
415  double m1;
416  double i;
417  double theta1;
418 // double temp;
419 // const double epsilon = 1.0e-2;
420 
421  long Error_Code = GRIN_NO_ERROR;
422 
423 // if ((Easting > (Grin_False_Easting + PI_Ra + epsilon)) ||
424 // (Easting < (Grin_False_Easting - PI_Ra - epsilon)))
425 // { /* Easting out of range */
426 // Error_Code |= GRIN_EASTING_ERROR;
427 // }
428 // if ((Northing > (Grin_False_Northing + PI_Ra + epsilon)) ||
429 // (Northing < (Grin_False_Northing - PI_Ra - epsilon)))
430 // { /* Northing out of range */
431 // Error_Code |= GRIN_NORTHING_ERROR;
432 // }
433 // if (!Error_Code)
434 // {
435 // temp = sqrt(Easting * Easting + Northing * Northing);
436 
437 // if ((temp > (Grin_False_Easting + PI_Ra + epsilon)) ||
438 // (temp > (Grin_False_Northing + PI_Ra + epsilon)) ||
439 // (temp < (Grin_False_Easting - PI_Ra - epsilon)) ||
440 // (temp < (Grin_False_Northing - PI_Ra - epsilon)))
441 // { /* Point is outside of projection area */
442 // Error_Code |= GRIN_RADIUS_ERROR;
443 // }
444 // }
445 
446  if (!Error_Code)
447  {
448  dy = Northing - Grin_False_Northing;
449  dx = Easting - Grin_False_Easting;
450  xx = dx / PI_Ra;
451  yy = dy / PI_Ra;
452  xxsqr = xx * xx;
453  yysqr = yy * yy;
454  xxsqr_PLUS_yysqr = xxsqr + yysqr;
455  two_yysqr = 2 * yysqr;
456 
457  if (Northing == 0.0)
458  *Latitude = 0.0;
459 
460  else
461  {
462  c1 = - fabs(yy) * (1 + xxsqr_PLUS_yysqr);
463  c2 = c1 - two_yysqr + xxsqr;
464  c3 = - 2 * c1 + 1 + two_yysqr + (xxsqr_PLUS_yysqr) * (xxsqr_PLUS_yysqr);
465  c2_OVER_3c3 = c2 / (3.0 * c3);
466  c3sqr = c3 * c3;
467  dd = yysqr / c3 + ((2 * c2 * c2 * c2) / (c3sqr * c3) - (9 * c1 * c2) / (c3sqr)) / 27;
468  a1 = (c1 - c2 * c2_OVER_3c3) /c3;
469  m1 = 2 * sqrt(-ONE_THIRD * a1);
470  i = 3 * dd/ (a1 * m1);
471  if ((i > 1.0)||(i < -1.0))
472  *Latitude = MAX_LAT;
473  else
474  {
475  theta1 = ONE_THIRD * acos(3 * dd / (a1 * m1));
476  *Latitude = M_PI * (-m1 * cos(theta1 + PI_OVER_3) - c2_OVER_3c3);
477  }
478  }
479  if (Northing < 0.0)
480  *Latitude *= -1.0;
481 
482  if (xx == 0.0)
483  *Longitude = Grin_Origin_Long;
484  else
485  {
486  *Longitude = M_PI * (xxsqr_PLUS_yysqr - 1 +
487  sqrt(1 + (2 * xxsqr - two_yysqr) + (xxsqr_PLUS_yysqr) * (xxsqr_PLUS_yysqr))) /
488  (2 * xx) + Grin_Origin_Long;
489  }
490 // if (*Latitude > PI_OVER_2) /* force distorted values to 90, -90 degrees */
491 // *Latitude = PI_OVER_2;
492 // else if (*Latitude < -PI_OVER_2)
493 // *Latitude = -PI_OVER_2;
494 
495 // if (*Longitude > M_PI)
496 // *Longitude -= TWO_PI;
497 // if (*Longitude < -M_PI)
498 // *Longitude += TWO_PI;
499 
500 // if (*Longitude > M_PI) /* force distorted values to 180, -180 degrees */
501 // *Longitude = M_PI;
502 // else if (*Longitude < -M_PI)
503 // *Longitude = -M_PI;
504 
505  }
506  return (Error_Code);
507 
508 } /* END OF Convert_Van_der_Grinten_To_Geodetic */
#define M_PI
const double ONE_THIRD
#define GRIN_NO_ERROR
const double PI_OVER_3

◆ dup()

virtual ossimObject* ossimVanDerGrintenProjection::dup ( ) const
inlinevirtual

Implements ossimProjection.

Definition at line 27 of file ossimVanDerGrintenProjection.h.

28  {
29  return new ossimVanDerGrintenProjection(*this);
30  }
ossimVanDerGrintenProjection(const ossimEllipsoid &ellipsoid=ossimEllipsoid(), const ossimGpt &origin=ossimGpt())

◆ forward()

ossimDpt ossimVanDerGrintenProjection::forward ( const ossimGpt worldPoint) const
virtual

All map projections will convert the world coordinate to an easting northing (Meters).

Implements ossimMapProjection.

Definition at line 136 of file ossimVanDerGrintenProjection.cpp.

References ossimGpt::changeDatum(), ossimDatum::code(), Convert_Geodetic_To_Van_der_Grinten(), ossimGpt::datum(), ossimGpt::latr(), ossimGpt::lonr(), and ossimMapProjection::theDatum.

137 {
138  double easting = 0.0;
139  double northing = 0.0;
140  ossimGpt gpt = latLon;
141 
142  if (theDatum)
143  {
144  if (theDatum->code() != latLon.datum()->code())
145  {
146  gpt.changeDatum(theDatum); // Shift to our datum.
147  }
148  }
149 
151  gpt.lonr(),
152  &easting,
153  &northing);
154  return ossimDpt(easting, northing);
155 }
virtual const ossimString & code() const
Definition: ossimDatum.h:57
void changeDatum(const ossimDatum *datum)
This will actually perform a shift.
Definition: ossimGpt.cpp:316
double lonr() const
Returns the longitude in radian measure.
Definition: ossimGpt.h:76
long Convert_Geodetic_To_Van_der_Grinten(double Latitude, double Longitude, double *Easting, double *Northing) const
double latr() const
latr().
Definition: ossimGpt.h:66
const ossimDatum * theDatum
This is only set if we want to have built in datum shifting.

◆ Get_Van_der_Grinten_Parameters()

void ossimVanDerGrintenProjection::Get_Van_der_Grinten_Parameters ( double *  a,
double *  f,
double *  Central_Meridian,
double *  False_Easting,
double *  False_Northing 
) const
protected

The function Get_Van_der_Grinten_Parameters returns the current ellipsoid parameters, and Van Der Grinten projection parameters.

a : Semi-major axis of ellipsoid, in meters (output) f : Flattening of ellipsoid (output) Central_Meridian : Longitude in radians at the center of (output) the projection False_Easting : A coordinate value in meters assigned to the central meridian of the projection. (output) False_Northing : A coordinate value in meters assigned to the origin latitude of the projection (output)

Definition at line 247 of file ossimVanDerGrintenProjection.cpp.

References Grin_a, Grin_f, Grin_False_Easting, Grin_False_Northing, and Grin_Origin_Long.

253 { /* BEGIN Get_Van_der_Grinten_Parameters */
254 /*
255  * The function Get_Van_der_Grinten_Parameters returns the current ellipsoid
256  * parameters, and Van Der Grinten projection parameters.
257  *
258  * a : Semi-major axis of ellipsoid, in meters (output)
259  * f : Flattening of ellipsoid (output)
260  * Central_Meridian : Longitude in radians at the center of (output)
261  * the projection
262  * False_Easting : A coordinate value in meters assigned to the
263  * central meridian of the projection. (output)
264  * False_Northing : A coordinate value in meters assigned to the
265  * origin latitude of the projection (output)
266  */
267 
268  *a = Grin_a;
269  *f = Grin_f;
270  *Central_Meridian = Grin_Origin_Long;
271  *False_Easting = Grin_False_Easting;
272  *False_Northing = Grin_False_Northing;
273 
274  return;
275 } /* END OF Get_Van_der_Grinten_Parameters */

◆ getFalseEasting()

double ossimVanDerGrintenProjection::getFalseEasting ( ) const
inlinevirtual
Returns
The false easting.

Reimplemented from ossimMapProjection.

Definition at line 60 of file ossimVanDerGrintenProjection.h.

◆ getFalseNorthing()

double ossimVanDerGrintenProjection::getFalseNorthing ( ) const
inlinevirtual
Returns
The false northing.

Reimplemented from ossimMapProjection.

Definition at line 61 of file ossimVanDerGrintenProjection.h.

◆ inverse()

ossimGpt ossimVanDerGrintenProjection::inverse ( const ossimDpt projectedPoint) const
virtual

Will take a point in meters and convert it to ground.

Implements ossimMapProjection.

Definition at line 123 of file ossimVanDerGrintenProjection.cpp.

References Convert_Van_der_Grinten_To_Geodetic(), DEG_PER_RAD, ossimMapProjection::theDatum, ossimDpt::x, and ossimDpt::y.

124 {
125  double lat = 0.0;
126  double lon = 0.0;
127 
128  Convert_Van_der_Grinten_To_Geodetic(eastingNorthing.x,
129  eastingNorthing.y,
130  &lat,
131  &lon);
132 
133  return ossimGpt(lat*DEG_PER_RAD, lon*DEG_PER_RAD, 0.0, theDatum);
134 }
#define DEG_PER_RAD
long Convert_Van_der_Grinten_To_Geodetic(double Easting, double Northing, double *Latitude, double *Longitude) const
const ossimDatum * theDatum
This is only set if we want to have built in datum shifting.

◆ loadState()

bool ossimVanDerGrintenProjection::loadState ( const ossimKeywordlist kwl,
const char *  prefix = 0 
)
virtual

Method to the load (recreate) the state of an object from a keyword list. Return true if ok or false on error.

Reimplemented from ossimMapProjection.

Definition at line 162 of file ossimVanDerGrintenProjection.cpp.

References ossimKeywordlist::find(), Grin_False_Easting, Grin_False_Northing, ossimMapProjection::loadState(), setDefaults(), STATIC_TYPE_NAME, ossimMapProjection::theFalseEastingNorthing, ossimKeywordNames::TYPE_KW, update(), ossimDpt::x, and ossimDpt::y.

164 {
165  bool flag = ossimMapProjection::loadState(kwl, prefix);
166  const char* type = kwl.find(prefix, ossimKeywordNames::TYPE_KW);
167 
168  setDefaults();
169 
171  {
174  }
175 
176  update();
177 
178  return flag;
179 
180 }
virtual bool loadState(const ossimKeywordlist &kwl, const char *prefix=0)
Method to the load (recreate) the state of an object from a keyword list.
const char * find(const char *key) const
double y
Definition: ossimDpt.h:165
static const char * TYPE_KW
#define STATIC_TYPE_NAME(T)
Definition: ossimRtti.h:325
double x
Definition: ossimDpt.h:164
ossimDpt theFalseEastingNorthing
Hold the false easting northing.

◆ operator==()

bool ossimVanDerGrintenProjection::operator== ( const ossimProjection projection) const
virtual

Returns TRUE if principal parameters are within epsilon tolerance.

Reimplemented from ossimMapProjection.

Definition at line 513 of file ossimVanDerGrintenProjection.cpp.

References ossim::almostEqual(), and Grin_Origin_Long.

514 {
515  if (!ossimMapProjection::operator==(proj))
516  return false;
517 
519  dynamic_cast<const ossimVanDerGrintenProjection*>(&proj);
520  if (!p) return false;
521 
523 
524  return true;
525 }
bool almostEqual(T x, T y, T tolerance=FLT_EPSILON)
Definition: ossimCommon.h:53

◆ saveState()

bool ossimVanDerGrintenProjection::saveState ( ossimKeywordlist kwl,
const char *  prefix = 0 
) const
virtual

Method to save the state of an object to a keyword list. Return true if ok or false on error.

Reimplemented from ossimMapProjection.

Definition at line 157 of file ossimVanDerGrintenProjection.cpp.

References ossimMapProjection::saveState().

158 {
159  return ossimMapProjection::saveState(kwl, prefix);
160 }
virtual bool saveState(ossimKeywordlist &kwl, const char *prefix=0) const
Method to save the state of an object to a keyword list.

◆ Set_Van_der_Grinten_Parameters()

long ossimVanDerGrintenProjection::Set_Van_der_Grinten_Parameters ( double  a,
double  f,
double  Central_Meridian,
double  False_Easting,
double  False_Northing 
)
protected

The function Set_Van_der_Grinten_Parameters receives the ellipsoid parameters and Van Der Grinten projcetion parameters as inputs, and sets the corresponding state variables. If any errors occur, the error code(s) are returned by the function, otherwise GRIN_NO_ERROR is returned.

a : Semi-major axis of ellipsoid, in meters (input) f : Flattening of ellipsoid (input) Central_Meridian : Longitude in radians at the center of (input) the projection False_Easting : A coordinate value in meters assigned to the central meridian of the projection. (input) False_Northing : A coordinate value in meters assigned to the origin latitude of the projection (input)

Definition at line 188 of file ossimVanDerGrintenProjection.cpp.

References es2, es4, es6, Grin_a, Grin_f, Grin_False_Easting, Grin_False_Northing, GRIN_NO_ERROR, Grin_Origin_Long, M_PI, PI_Ra, and Ra.

Referenced by update().

194 { /* BEGIN Set_Van_der_Grinten_Parameters */
195 /*
196  * The function Set_Van_der_Grinten_Parameters receives the ellipsoid parameters and
197  * projection parameters as inputs, and sets the corresponding state
198  * variables. If any errors occur, the error code(s) are returned by the function,
199  * otherwise Grin_NO_ERROR is returned.
200  *
201  * a : Semi-major axis of ellipsoid, in meters (input)
202  * f : Flattening of ellipsoid (input)
203  * Central_Meridian : Longitude in radians at the center of (input)
204  * the projection
205  * False_Easting : A coordinate value in meters assigned to the
206  * central meridian of the projection. (input)
207  * False_Northing : A coordinate value in meters assigned to the
208  * origin latitude of the projection (input)
209  */
210 
211 // double inv_f = 1 / f;
212  long Error_Code = GRIN_NO_ERROR;
213 
214 // if (a <= 0.0)
215 // { /* Semi-major axis must be greater than zero */
216 // Error_Code |= GRIN_A_ERROR;
217 // }
218 // if ((inv_f < 250) || (inv_f > 350))
219 // { /* Inverse flattening must be between 250 and 350 */
220 // Error_Code |= GRIN_INV_F_ERROR;
221 // }
222 // if ((Central_Meridian < -M_PI) || (Central_Meridian > TWO_PI))
223 // { /* origin longitude out of range */
224 // Error_Code |= GRIN_CENT_MER_ERROR;
225 // }
226  if (!Error_Code)
227  { /* no errors */
228  Grin_a = a;
229  Grin_f = f;
230  es2 = 2 * Grin_f - Grin_f * Grin_f;
231  es4 = es2 * es2;
232  es6 = es4 * es2;
233  /* spherical radius */
234  Ra = Grin_a * (1.0 - es2 / 6.0 - 17.0 * es4 / 360.0 - 67.0 * es6 /3024.0);
235  PI_Ra = M_PI * Ra;
236 // if (Central_Meridian > M_PI)
237 // Central_Meridian -= TWO_PI;
238  Grin_Origin_Long = Central_Meridian;
239  Grin_False_Easting = False_Easting;
240  Grin_False_Northing = False_Northing;
241 
242  } /* END OF if(!Error_Code) */
243  return (Error_Code);
244 } /* END OF Set_Van_der_Grinten_Parameters */
#define M_PI
#define GRIN_NO_ERROR

◆ setCentralMeridian()

void ossimVanDerGrintenProjection::setCentralMeridian ( double  centralMeridian)

Definition at line 108 of file ossimVanDerGrintenProjection.cpp.

References Grin_Origin_Long, and update().

109 {
110  Grin_Origin_Long = centralMeridian;
111  update();
112 }

◆ setDefaults()

void ossimVanDerGrintenProjection::setDefaults ( )

◆ setFalseEasting()

void ossimVanDerGrintenProjection::setFalseEasting ( double  falseEasting)

SetFalseEasting. The value is in meters. Update is then called so we can pre-compute paramters

Definition at line 88 of file ossimVanDerGrintenProjection.cpp.

References Grin_False_Easting, and update().

89 {
90  Grin_False_Easting = falseEasting;
91 
92  update();
93 }

◆ setFalseEastingNorthing()

void ossimVanDerGrintenProjection::setFalseEastingNorthing ( double  falseEasting,
double  falseNorthing 
)

Sets both false easting and northing values. The values are expected to be in meters. Update is then called so we can pre-compute paramters

Definition at line 114 of file ossimVanDerGrintenProjection.cpp.

References Grin_False_Easting, Grin_False_Northing, and update().

116 {
117  Grin_False_Easting = falseEasting;
118  Grin_False_Northing = falseNorthing;
119 
120  update();
121 }

◆ setFalseNorthing()

void ossimVanDerGrintenProjection::setFalseNorthing ( double  falseNorthing)

SetFalseNorthing. The value is in meters. Update is then called so we can pre-compute paramters

Definition at line 95 of file ossimVanDerGrintenProjection.cpp.

References Grin_False_Northing, and update().

96 {
97  Grin_False_Northing = falseNorthing;
98 
99  update();
100 }

◆ update()

void ossimVanDerGrintenProjection::update ( )
virtual

Reimplemented from ossimMapProjection.

Definition at line 74 of file ossimVanDerGrintenProjection.cpp.

References ossimEllipsoid::getA(), ossimEllipsoid::getFlattening(), Grin_False_Easting, Grin_False_Northing, ossimGpt::lonr(), Set_Van_der_Grinten_Parameters(), ossimMapProjection::theEllipsoid, ossimMapProjection::theFalseEastingNorthing, ossimMapProjection::theOrigin, ossimMapProjection::update(), ossimDpt::x, and ossimDpt::y.

Referenced by loadState(), ossimVanDerGrintenProjection(), setCentralMeridian(), setFalseEasting(), setFalseEastingNorthing(), and setFalseNorthing().

75 {
78  theOrigin.lonr(),
81 
84 
86 }
double y
Definition: ossimDpt.h:165
const double & getA() const
double lonr() const
Returns the longitude in radian measure.
Definition: ossimGpt.h:76
double x
Definition: ossimDpt.h:164
const double & getFlattening() const
ossimEllipsoid theEllipsoid
This method verifies that the projection parameters match the current pcs code.
ossimDpt theFalseEastingNorthing
Hold the false easting northing.
long Set_Van_der_Grinten_Parameters(double a, double f, double Central_Meridian, double False_Easting, double False_Northing)

Member Data Documentation

◆ es2

double ossimVanDerGrintenProjection::es2
mutableprotected

Definition at line 87 of file ossimVanDerGrintenProjection.h.

Referenced by Set_Van_der_Grinten_Parameters().

◆ es4

double ossimVanDerGrintenProjection::es4
mutableprotected

Definition at line 88 of file ossimVanDerGrintenProjection.h.

Referenced by Set_Van_der_Grinten_Parameters().

◆ es6

double ossimVanDerGrintenProjection::es6
mutableprotected

Definition at line 89 of file ossimVanDerGrintenProjection.h.

Referenced by Set_Van_der_Grinten_Parameters().

◆ Grin_a

double ossimVanDerGrintenProjection::Grin_a
mutableprotected

◆ Grin_f

double ossimVanDerGrintenProjection::Grin_f
mutableprotected

◆ Grin_False_Easting

double ossimVanDerGrintenProjection::Grin_False_Easting
mutableprotected

◆ Grin_False_Northing

double ossimVanDerGrintenProjection::Grin_False_Northing
mutableprotected

◆ Grin_Origin_Long

double ossimVanDerGrintenProjection::Grin_Origin_Long
mutableprotected

◆ PI_Ra

double ossimVanDerGrintenProjection::PI_Ra
mutableprotected

◆ Ra

double ossimVanDerGrintenProjection::Ra
mutableprotected

The documentation for this class was generated from the following files: