public final class AnomalyUtils extends Object
Specific methods for elliptic or hyperbolic orbits can be used.
Specific methods are available depending if the eccentricity is expressed with one component (e) or two components
(ex / ey).
| Modifier and Type | Field and Description |
|---|---|
static int |
MAX_ITERATIONS_DEFAULT_VALUE
Default value for the maximum number of iterations for convergence algorithms.
|
| Modifier and Type | Method and Description |
|---|---|
static double |
eccentricToMeanElliptic(double e,
double eccentricAnomaly)
Convert eccentric anomaly → mean anomaly (elliptic orbit).
|
static double |
eccentricToMeanElliptic(double ex,
double ey,
double alphaE)
Convert eccentric latitude/longitude argument → mean latitude/longitude argument (elliptic orbit).
|
static double |
eccentricToMeanHyperbolic(double e,
double eccentricAnomaly)
Convert eccentric anomaly → mean anomaly (hyperbolic orbit).
|
static double |
eccentricToMeanHyperbolic(double ex,
double ey,
double alphaE)
Convert eccentric latitude/longitude argument → mean latitude/longitude argument (hyperbolic orbit).
|
static double |
eccentricToTrueElliptic(double e,
double eccentricAnomaly)
Convert eccentric anomaly → true anomaly (elliptic orbit).
|
static double |
eccentricToTrueElliptic(double ex,
double ey,
double alphaE)
Convert eccentric latitude/longitude argument → true latitude/longitude argument (elliptic orbit).
|
static double |
eccentricToTrueHyperbolic(double e,
double eccentricAnomaly)
Convert eccentric anomaly → true anomaly (hyperbolic orbit).
|
static double |
eccentricToTrueHyperbolic(double ex,
double ey,
double alphaE)
Convert eccentric latitude/longitude argument → true latitude/longitude argument (hyperbolic orbit).
|
static Pair<Double,PositionAngle> |
getInitializedValue(double meanVal,
double eccentricVal,
double trueVal)
Evaluate which of the three values is initialized ({@code !
|
static double[] |
initializeValues(double value,
PositionAngle type)
Return an array describing the given value according to the given type and
Double.NaN for the other
elements. |
static double |
meanToEccentricElliptic(double e,
double meanAnomaly)
Convert mean anomaly → eccentric anomaly (elliptic orbit).
|
static double |
meanToEccentricElliptic(double ex,
double ey,
double alphaM)
Convert mean latitude/longitude argument → eccentric latitude/longitude argument (elliptic orbit).
|
static double |
meanToEccentricHyperbolic(double e,
double meanAnomaly)
Convert mean anomaly → eccentric anomaly (hyperbolic orbit).
|
static double |
meanToEccentricHyperbolic(double ex,
double ey,
double alphaM)
Convert mean latitude/longitude argument → eccentric latitude/longitude argument (hyperbolic orbit).
|
static void |
setMaxIterations(int maxIterations)
Setter for the maximum number of iterations for convergence algorithms.
|
static double |
trueToEccentricElliptic(double e,
double trueAnomaly)
Convert true anomaly → eccentric anomaly (elliptic orbit).
|
static double |
trueToEccentricElliptic(double ex,
double ey,
double alphaV)
Convert true latitude/longitude argument → eccentric latitude/longitude argument (elliptic orbit).
|
static double |
trueToEccentricHyperbolic(double e,
double trueAnomaly)
Convert true anomaly → eccentric anomaly (hyperbolic orbit).
|
static double |
trueToEccentricHyperbolic(double ex,
double ey,
double alphaV)
Convert true latitude/longitude argument → eccentric latitude/longitude argument (hyperbolic orbit).
|
public static final int MAX_ITERATIONS_DEFAULT_VALUE
public static double meanToEccentricElliptic(double e,
double meanAnomaly)
The algorithm used here for solving Kepler equation has been published in: "Procedures for solving Kepler's Equation", A. W. Odell and R. H. Gooding, Celestial Mechanics 38 (1986) 307-334
e - eccentricity (e >= 0)meanAnomaly - mean anomaly (rad)public static double meanToEccentricElliptic(double ex,
double ey,
double alphaM)
ex - e cos(ω), first component of circular eccentricity vectorey - e sin(ω), second component of circular eccentricity vectoralphaM - = M + ω mean latitude/longitude argument (rad)ConvergenceException - if the algorithm fails to convergepublic static double meanToEccentricHyperbolic(double e,
double meanAnomaly)
The algorithm used here for solving hyperbolic Kepler equation is a naive initialization and classical Halley method for iterations.
e - eccentricity (e >= 0)meanAnomaly - mean anomaly (rad)ConvergenceException - if the algorithm fails to convergepublic static double meanToEccentricHyperbolic(double ex,
double ey,
double alphaM)
The algorithm used here for solving hyperbolic Kepler equation is a naive initialization and classical Halley method for iterations adapted to circular parameters.
ex - e cos(ω), first component of circular eccentricity vectorey - e sin(ω), second component of circular eccentricity vectoralphaM - mean latitude/longitude argument (rad)ConvergenceException - if the algorithm fails to convergepublic static double trueToEccentricElliptic(double e,
double trueAnomaly)
e - eccentricity (e >= 0)trueAnomaly - true anomaly (rad)public static double trueToEccentricElliptic(double ex,
double ey,
double alphaV)
ex - e cos(ω), first component of circular eccentricity vectorey - e sin(ω), second component of circular eccentricity vectoralphaV - true latitude/longitude argument (rad)public static double trueToEccentricHyperbolic(double e,
double trueAnomaly)
e - eccentricity (e >= 0)trueAnomaly - true anomaly (rad)public static double trueToEccentricHyperbolic(double ex,
double ey,
double alphaV)
ex - e cos(ω), first component of circular eccentricity vectorey - e sin(ω), second component of circular eccentricity vectoralphaV - true latitude/longitude argument (rad)public static double eccentricToMeanElliptic(double e,
double eccentricAnomaly)
e - eccentricity (e >= 0)eccentricAnomaly - eccentric anomaly (rad)public static double eccentricToMeanElliptic(double ex,
double ey,
double alphaE)
ex - e cos(ω), first component of circular eccentricity vectorey - e sin(ω), second component of circular eccentricity vectoralphaE - eccentric argument of latitude/longitude (rad)public static double eccentricToMeanHyperbolic(double e,
double eccentricAnomaly)
e - eccentricity (e >= 0)eccentricAnomaly - eccentric anomaly (rad)public static double eccentricToMeanHyperbolic(double ex,
double ey,
double alphaE)
ex - e cos(ω), first component of circular eccentricity vectorey - e sin(ω), second component of circular eccentricity vectoralphaE - eccentric argument of latitude/longitude (rad)public static double eccentricToTrueElliptic(double e,
double eccentricAnomaly)
e - eccentricity (e >= 0)eccentricAnomaly - eccentric anomaly (rad)public static double eccentricToTrueElliptic(double ex,
double ey,
double alphaE)
ex - e cos(ω), first component of circular eccentricity vectorey - e sin(ω), second component of circular eccentricity vectoralphaE - eccentric argument of latitude/longitude (rad)public static double eccentricToTrueHyperbolic(double e,
double eccentricAnomaly)
e - eccentricity (e >= 0)eccentricAnomaly - hyperbolic eccentric anomaly (rad)public static double eccentricToTrueHyperbolic(double ex,
double ey,
double alphaE)
ex - e cos(ω), first component of circular eccentricity vectorey - e sin(ω), second component of circular eccentricity vectoralphaE - eccentric argument of latitude/longitude (rad)public static double[] initializeValues(double value,
PositionAngle type)
Double.NaN for the other
elements.
Array index:
value - Valuetype - TypeIllegalArgumentException - if the provided value is Double.NaNpublic static Pair<Double,PositionAngle> getInitializedValue(double meanVal, double eccentricVal, double trueVal)
!= Double.NaN), then return it with its associated
type.
Note: the method purpose is to avoid unnecessary computations to use already initialized values.
In that matter computation methods should NEVER be called to provide the inputs.
To be more explicit through an example, don't call "getInitializedValue(getMeanAnomaly(),
getEccentricAnomaly(), getTrueAnomaly())", but instead "getInitializedValue(this.meanAnomaly,
this.eccentricAnomaly, this.trueAnomaly)"
meanVal - Mean valueeccentricVal - Eccentric valuetrueVal - True valuepublic static void setMaxIterations(int maxIterations)
maxIterations - the value to setCopyright © 2025 CNES. All rights reserved.