GeographicLib 2.5
GeographicLib::AuxLatitude Class Reference

Conversions between auxiliary latitudes. More...

#include <GeographicLib/AuxLatitude.hpp>

Inheritance diagram for GeographicLib::AuxLatitude:
GeographicLib::DAuxLatitude

Public Types

enum  aux {
  GEOGRAPHIC , PARAMETRIC , GEOCENTRIC , RECTIFYING ,
  CONFORMAL , AUTHALIC , AUXNUMBER , PHI ,
  BETA , THETA , MU , CHI ,
  XI , COMMON , GEODETIC , REDUCED
}
 

Public Member Functions

 AuxLatitude (real a, real f)
 
AuxAngle Convert (int auxin, int auxout, const AuxAngle &zeta, bool exact=false) const
 
Math::real Convert (int auxin, int auxout, real zeta, bool exact=false) const
 
AuxAngle ToAuxiliary (int auxout, const AuxAngle &phi, real *diff=nullptr) const
 
AuxAngle FromAuxiliary (int auxin, const AuxAngle &zeta, int *niter=nullptr) const
 
Math::real RectifyingRadius (bool exact=false) const
 
Math::real AuthalicRadiusSquared (bool exact=false) const
 
Math::real EquatorialRadius () const
 
Math::real PolarSemiAxis () const
 
Math::real Flattening () const
 

Static Public Member Functions

static AuxLatitude axes (real a, real b)
 
static Math::real Clenshaw (bool sinp, real szeta, real czeta, const real c[], int K)
 
static const AuxLatitudeWGS84 ()
 

Static Public Attributes

static const int Lmax
 

Protected Member Functions

AuxAngle Parametric (const AuxAngle &phi, real *diff=nullptr) const
 
AuxAngle Geocentric (const AuxAngle &phi, real *diff=nullptr) const
 
AuxAngle Rectifying (const AuxAngle &phi, real *diff=nullptr) const
 
AuxAngle Conformal (const AuxAngle &phi, real *diff=nullptr) const
 
AuxAngle Authalic (const AuxAngle &phi, real *diff=nullptr) const
 

Detailed Description

Conversions between auxiliary latitudes.

This class is an implementation of the methods described in

The provides accurate conversions between geographic (phi, φ), parametric (beta, β), geocentric (theta, θ), rectifying (mu, μ), conformal (chi, χ), and authalic (xi, ξ) latitudes for an ellipsoid of revolution. A latitude is represented by the class AuxAngle in order to maintain precision close to the poles.

The class implements two methods for the conversion:

  • Direct evaluation of the defining equations, the exact method. These equations are formulated so as to preserve relative accuracy of the tangent of the latitude, ensuring high accuracy near the equator and the poles. Newton's method is used for those conversions that can't be expressed in closed form.
  • Expansions in powers of n, the third flattening, the series method. This delivers full accuracy for abs(f) ≤ 1/150. Here, f is the flattening of the ellipsoid.

The series method is the preferred method of conversion for any conversion involving μ, χ, or ξ, with abs(f) ≤ 1/150. The equations for the conversions between φ, β, and θ are sufficiently simple that the exact method should be used for such conversions and also for conversions with with abs(f) > 1/150.

Example of use:

// Example of using the GeographicLib::AuxLatitude class. See the paper
//
// - C. F. F. Karney,
// On auxiliary latitudes,
// Survey Review 56(395), 165--180 (2024).
// https://doi.org/10.1080/00396265.2023.2217604
// preprint: https://arxiv.org/abs/2212.05818
#include <iostream>
#include <iomanip>
#include <exception>
int main(int argc, const char* const argv[]) {
try {
typedef GeographicLib::AuxLatitude latitude;
typedef GeographicLib::AuxAngle angle;
if (argc != 3) {
std::cerr << "Usage: example-AuxLatitude <n> <base-lat>\n";
return 1;
}
double n = GeographicLib::Utility::fract<double>(std::string(argv[1]));
int auxin = GeographicLib::Utility::val<int>(std::string(argv[2]));
double a = 1+n, b = 1-n; // Equatorial radius and polar semi-axis
latitude aux(latitude::axes(a, b));
bool series = false; // Don't use series method
std::cout << std::setprecision(9) << std::fixed;
int m = 1;
for (int l = 0; l < 90*m; ++l) {
angle phi(angle::degrees((l+0.5)/m));
for (int auxout = 0; auxout < latitude::AUXNUMBER; ++auxout) {
angle eta = aux.Convert(auxin, auxout, phi, series);
std::cout << (auxout ? " " : "") << eta.degrees();
}
std::cout << "\n";
}
}
catch (const std::exception& e) {
std::cerr << "Caught exception: " << e.what() << "\n";
return 1;
}
}
Header for the GeographicLib::AuxLatitude class.
int main(int argc, const char *const argv[])
Header for GeographicLib::Utility class.
An accurate representation of angles.
Definition AuxAngle.hpp:47
Conversions between auxiliary latitudes.

Definition at line 75 of file AuxLatitude.hpp.

Member Enumeration Documentation

◆ aux

The floating-point type for real numbers. This just connects to the template parameters for the class. The different auxiliary latitudes.

Enumerator
GEOGRAPHIC 

Geographic latitude, phi, φ

PARAMETRIC 

Parametric latitude, beta, β

GEOCENTRIC 

Geocentric latitude, theta, θ

RECTIFYING 

Rectifying latitude, mu, μ

CONFORMAL 

Conformal latitude, chi, χ

AUTHALIC 

Authalic latitude, xi, ξ

AUXNUMBER 

The total number of auxiliary latitudes

PHI 

An alias for GEOGRAPHIC

BETA 

An alias for PARAMETRIC

THETA 

An alias for GEOCENTRIC

MU 

An alias for RECTIFYING

CHI 

An alias for CONFORMAL

XI 

An alias for AUTHALIC

COMMON 

An alias for GEOGRAPHIC

GEODETIC 

An alias for GEOGRAPHIC

REDUCED 

An alias for PARAMETRIC

Definition at line 86 of file AuxLatitude.hpp.

Constructor & Destructor Documentation

◆ AuxLatitude()

GeographicLib::AuxLatitude::AuxLatitude ( real  a,
real  f 
)

Constructor

Parameters
[in]aequatorial radius.
[in]fflattening of ellipsoid. Setting f = 0 gives a sphere. Negative f gives a prolate ellipsoid.
Exceptions
GeographicErrif a or (1 − f) a is not positive.
Note
the constructor does not precompute the coefficients for the Fourier series for the series conversions. These are computed and saved when first needed.

Definition at line 25 of file AuxLatitude.cpp.

References AUXNUMBER, and Lmax.

Member Function Documentation

◆ axes()

static AuxLatitude GeographicLib::AuxLatitude::axes ( real  a,
real  b 
)
inlinestatic

Construct and return an AuxLatitude object specified in terms of the semi-axes

Parameters
[in]aequatorial radius.
[in]bpolar semi-axis.
Exceptions
GeographicErrif a or b is not positive.

This allows a new AuxAngle to be initialized as an angle in radians with

static AuxLatitude axes(real a, real b)

Definition at line 195 of file AuxLatitude.hpp.

◆ Convert() [1/2]

AuxAngle GeographicLib::AuxLatitude::Convert ( int  auxin,
int  auxout,
const AuxAngle zeta,
bool  exact = false 
) const

Convert between any two auxiliary latitudes specified as AuxAngle.

Parameters
[in]auxinan AuxLatitude::aux indicating the type of auxiliary latitude zeta.
[in]auxoutan AuxLatitude::aux indicating the type of auxiliary latitude eta.
[in]zetathe input auxiliary latitude as an AuxAngle
[in]exactif true use the exact equations instead of the Taylor series [default false].
Returns
the output auxiliary latitude eta as an AuxAngle.

With exact = false, the Fourier coefficients for a specific auxin and auxout are computed and saved on the first call; the saved coefficients are used on subsequent calls. The series method is accurate for abs(f) ≤ 1/150; for other f, the exact method should be used.

Definition at line 297 of file AuxLatitude.cpp.

References Clenshaw(), FromAuxiliary(), Lmax, GeographicLib::AuxAngle::NaN(), GeographicLib::AuxAngle::normalized(), GeographicLib::AuxAngle::radians(), ToAuxiliary(), GeographicLib::AuxAngle::x(), and GeographicLib::AuxAngle::y().

Referenced by Convert(), GeographicLib::Rhumb::GenInverse(), GeographicLib::RhumbLine::GenPosition(), and main().

◆ Convert() [2/2]

Math::real GeographicLib::AuxLatitude::Convert ( int  auxin,
int  auxout,
real  zeta,
bool  exact = false 
) const

Convert between any two auxiliary latitudes specified in degrees.

Parameters
[in]auxinan AuxLatitude::aux indicating the type of auxiliary latitude zeta.
[in]auxoutan AuxLatitude::aux indicating the type of auxiliary latitude eta.
[in]zetathe input auxiliary latitude in degrees.
[in]exactif true use the exact equations instead of the Taylor series [default false].
Returns
the output auxiliary latitude eta in degrees.

With exact = false, the Fourier coefficients for a specific auxin and auxout are computed and saved on the first call; the saved coefficients are used on subsequent calls. The series method is accurate for abs(f) ≤ 1/150; for other f, the exact method should be used.

Definition at line 318 of file AuxLatitude.cpp.

References Convert(), GeographicLib::AuxAngle::degrees(), and GeographicLib::Math::td.

◆ ToAuxiliary()

AuxAngle GeographicLib::AuxLatitude::ToAuxiliary ( int  auxout,
const AuxAngle phi,
real *  diff = nullptr 
) const

Convert geographic latitude to an auxiliary latitude eta.

Parameters
[in]auxoutan AuxLatitude::aux indicating the auxiliary latitude returned.
[in]phithe geographic latitude.
[out]diffoptional pointer to the derivative d tan(eta) / d tan(phi).
Returns
the auxiliary latitude eta.

This uses the exact equations.

Definition at line 220 of file AuxLatitude.cpp.

References AUTHALIC, Authalic(), CONFORMAL, Conformal(), GEOCENTRIC, GEOGRAPHIC, GeographicLib::AuxAngle::NaN(), PARAMETRIC, Parametric(), RECTIFYING, and Rectifying().

Referenced by Convert(), and FromAuxiliary().

◆ FromAuxiliary()

AuxAngle GeographicLib::AuxLatitude::FromAuxiliary ( int  auxin,
const AuxAngle zeta,
int *  niter = nullptr 
) const

Convert an auxiliary latitude zeta to geographic latitude.

Parameters
[in]auxinan AuxLatitude::aux indicating the type of auxiliary latitude zeta.
[in]zetathe input auxiliary latitude.
[out]niteroptional pointer to the number of iterations.
Returns
the geographic latitude phi.

This uses the exact equations.

Definition at line 236 of file AuxLatitude.cpp.

References AUTHALIC, CONFORMAL, GeographicLib::AuxAngle::copyquadrant(), GEOCENTRIC, GEOGRAPHIC, GeographicLib::AuxAngle::NaN(), PARAMETRIC, RECTIFYING, GeographicLib::AuxAngle::tan(), ToAuxiliary(), GeographicLib::AuxAngle::x(), and GeographicLib::AuxAngle::y().

Referenced by Convert().

◆ RectifyingRadius()

Math::real GeographicLib::AuxLatitude::RectifyingRadius ( bool  exact = false) const

Return the rectifying radius.

Parameters
[in]exactif true use the exact expression instead of the Taylor series [default false].
Returns
the rectifying radius in the same units as a.

Definition at line 325 of file AuxLatitude.cpp.

References Lmax, GeographicLib::Math::pi(), GeographicLib::Math::polyval(), GeographicLib::EllipticFunction::RG(), and GeographicLib::Math::sq().

Referenced by GeographicLib::DAuxLatitude::DRectifying().

◆ AuthalicRadiusSquared()

Math::real GeographicLib::AuxLatitude::AuthalicRadiusSquared ( bool  exact = false) const

Return the authalic radius squared.

Parameters
[in]exactif true use the exact expression instead of the Taylor series [default false].
Returns
the authalic radius squared in the same units as a.

Definition at line 352 of file AuxLatitude.cpp.

References Lmax, GeographicLib::Math::polyval(), and GeographicLib::Math::sq().

Referenced by main().

◆ EquatorialRadius()

Math::real GeographicLib::AuxLatitude::EquatorialRadius ( ) const
inline
Returns
a the equatorial radius of the ellipsoid (meters).

Definition at line 284 of file AuxLatitude.hpp.

◆ PolarSemiAxis()

Math::real GeographicLib::AuxLatitude::PolarSemiAxis ( ) const
inline
Returns
b the polar semi-axis of the ellipsoid (meters).

Definition at line 288 of file AuxLatitude.hpp.

◆ Flattening()

Math::real GeographicLib::AuxLatitude::Flattening ( ) const
inline
Returns
f, the flattening of the ellipsoid.

Definition at line 292 of file AuxLatitude.hpp.

◆ Clenshaw()

static Math::real GeographicLib::AuxLatitude::Clenshaw ( bool  sinp,
real  szeta,
real  czeta,
const real  c[],
int  K 
)
static

Use Clenshaw to sum a Fouier series.

Parameters
[in]sinpif true sum the sine series, else sum the cosine series.
[in]szetasin(zeta).
[in]czetacos(zeta).
[in]cthe array of coefficients.
[in]Kthe number of coefficients.
Returns
the Clenshaw sum.

The result returned is \( \sum_0^{K-1} c_k \sin (2k+2)\zeta \), if sinp is true; with sinp false, replace sin by cos.

Referenced by Convert().

◆ WGS84()

const AuxLatitude & GeographicLib::AuxLatitude::WGS84 ( )
static

A global instantiation of Ellipsoid with the parameters for the WGS84 ellipsoid.

Definition at line 79 of file AuxLatitude.cpp.

References GeographicLib::Constants::WGS84_a(), and GeographicLib::Constants::WGS84_f().

◆ Parametric()

AuxAngle GeographicLib::AuxLatitude::Parametric ( const AuxAngle phi,
real *  diff = nullptr 
) const
protected

Convert geographic latitude to parametric latitude

Parameters
[in]phigeographic latitude.
[out]diffoptional pointer to the derivative d tan(beta) / d tan(phi).
Returns
beta, the parametric latitude

Definition at line 84 of file AuxLatitude.cpp.

References GeographicLib::AuxAngle::x(), and GeographicLib::AuxAngle::y().

Referenced by Authalic(), Conformal(), GeographicLib::DAuxLatitude::DRectifying(), Rectifying(), and ToAuxiliary().

◆ Geocentric()

AuxAngle GeographicLib::AuxLatitude::Geocentric ( const AuxAngle phi,
real *  diff = nullptr 
) const
protected

Convert geographic latitude to geocentric latitude

Parameters
[in]phigeographic latitude.
[out]diffoptional pointer to the derivative d tan(theta) / d tan(phi).
Returns
theta, the geocentric latitude.

Definition at line 89 of file AuxLatitude.cpp.

References GeographicLib::AuxAngle::x(), and GeographicLib::AuxAngle::y().

◆ Rectifying()

AuxAngle GeographicLib::AuxLatitude::Rectifying ( const AuxAngle phi,
real *  diff = nullptr 
) const
protected

Convert geographic latitude to rectifying latitude

Parameters
[in]phigeographic latitude.
[out]diffoptional pointer to the derivative d tan(mu) / d tan(phi).
Returns
mu, the rectifying latitude.

Definition at line 94 of file AuxLatitude.cpp.

References GeographicLib::AuxAngle::normalized(), Parametric(), GeographicLib::Math::pi(), GeographicLib::EllipticFunction::RD(), GeographicLib::EllipticFunction::RF(), GeographicLib::Math::sq(), GeographicLib::AuxAngle::tan(), GeographicLib::AuxAngle::x(), and GeographicLib::AuxAngle::y().

Referenced by GeographicLib::DAuxLatitude::DRectifying(), and ToAuxiliary().

◆ Conformal()

AuxAngle GeographicLib::AuxLatitude::Conformal ( const AuxAngle phi,
real *  diff = nullptr 
) const
protected

Convert geographic latitude to conformal latitude

Parameters
[in]phigeographic latitude.
[out]diffoptional pointer to the derivative d tan(chi) / d tan(phi).
Returns
chi, the conformal latitude.

Definition at line 140 of file AuxLatitude.cpp.

References GeographicLib::AuxAngle::normalized(), Parametric(), GeographicLib::AuxAngle::tan(), and GeographicLib::AuxAngle::x().

Referenced by ToAuxiliary().

◆ Authalic()

AuxAngle GeographicLib::AuxLatitude::Authalic ( const AuxAngle phi,
real *  diff = nullptr 
) const
protected

Convert geographic latitude to authalic latitude

Parameters
[in]phigeographic latitude.
[out]diffoptional pointer to the derivative d tan(xi) / d tan(phi).
Returns
xi, the authalic latitude.

Definition at line 198 of file AuxLatitude.cpp.

References GeographicLib::AuxAngle::normalized(), Parametric(), GeographicLib::Math::sq(), GeographicLib::AuxAngle::tan(), GeographicLib::AuxAngle::x(), and GeographicLib::AuxAngle::y().

Referenced by ToAuxiliary().

Member Data Documentation

◆ Lmax

const int GeographicLib::AuxLatitude::Lmax
static

The order of the series expansions. This is set at compile time to either 4, 6, or 8, by the preprocessor macro GEOGRAPHICLIB_AUXLATITUDE_ORDER.

Definition at line 316 of file AuxLatitude.hpp.

Referenced by AuthalicRadiusSquared(), AuxLatitude(), Convert(), GeographicLib::DAuxLatitude::DConvert(), and RectifyingRadius().


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