public class EquinoctialOrbit extends Orbit
The parameters used internally are the equinoctial elements which can be related to Keplerian elements as follows:
a ex = e cos(ω + Ω) ey = e sin(ω + Ω) hx = tan(i/2) cos(Ω) hy = tan(i/2) sin(Ω) lv = v + ω + Ωwhere ω stands for the Perigee Argument and Ω stands for the Right Ascension of the Ascending Node.
The conversion equations from and to Keplerian elements given above hold only when both sides are unambiguously defined, i.e. when orbit is neither equatorial nor circular. When orbit is either equatorial or circular, the equinoctial parameters are still unambiguously defined whereas some Keplerian elements (more precisely ω and Ω) become ambiguous. For this reason, equinoctial parameters are the recommended way to represent orbits.
The instance EquinoctialOrbit
is guaranteed to be immutable.
Orbit
,
KeplerianOrbit
,
CircularOrbit
,
CartesianOrbit
,
Serialized FormConstructor and Description |
---|
EquinoctialOrbit(double a,
double ex,
double ey,
double hx,
double hy,
double l,
double aDot,
double exDot,
double eyDot,
double hxDot,
double hyDot,
double lDot,
PositionAngle type,
Frame frame,
AbsoluteDate date,
double mu)
Creates a new instance.
|
EquinoctialOrbit(double a,
double ex,
double ey,
double hx,
double hy,
double l,
PositionAngle type,
Frame frame,
AbsoluteDate date,
double mu)
Creates a new instance.
|
EquinoctialOrbit(Orbit op)
Constructor from any kind of orbital parameters.
|
EquinoctialOrbit(PVCoordinates pvCoordinates,
Frame frame,
AbsoluteDate date,
double mu)
Constructor from Cartesian parameters.
|
EquinoctialOrbit(TimeStampedPVCoordinates pvCoordinates,
Frame frame,
double mu)
Constructor from Cartesian parameters.
|
Modifier and Type | Method and Description |
---|---|
void |
addKeplerContribution(PositionAngle type,
double gm,
double[] pDot)
Add the contribution of the Keplerian motion to parameters derivatives
|
protected double[][] |
computeJacobianEccentricWrtCartesian()
Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
|
protected double[][] |
computeJacobianMeanWrtCartesian()
Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
|
protected double[][] |
computeJacobianTrueWrtCartesian()
Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
|
static double |
eccentricToMean(double lE,
double ex,
double ey)
Computes the mean longitude argument from the eccentric longitude argument.
|
static double |
eccentricToTrue(double lE,
double ex,
double ey)
Computes the true longitude argument from the eccentric longitude argument.
|
double |
getA()
Get the semi-major axis.
|
double |
getADot()
Get the semi-major axis derivative.
|
double |
getE()
Get the eccentricity.
|
double |
getEDot()
Get the eccentricity derivative.
|
double |
getEquinoctialEx()
Get the first component of the equinoctial eccentricity vector derivative.
|
double |
getEquinoctialExDot()
Get the first component of the equinoctial eccentricity vector.
|
double |
getEquinoctialEy()
Get the second component of the equinoctial eccentricity vector derivative.
|
double |
getEquinoctialEyDot()
Get the second component of the equinoctial eccentricity vector.
|
double |
getHx()
Get the first component of the inclination vector.
|
double |
getHxDot()
Get the first component of the inclination vector derivative.
|
double |
getHy()
Get the second component of the inclination vector.
|
double |
getHyDot()
Get the second component of the inclination vector derivative.
|
double |
getI()
Get the inclination.
|
double |
getIDot()
Get the inclination derivative.
|
double |
getL(PositionAngle type)
Get the longitude argument.
|
double |
getLDot(PositionAngle type)
Get the longitude argument derivative.
|
double |
getLE()
Get the eccentric longitude argument.
|
double |
getLEDot()
Get the eccentric longitude argument derivative.
|
double |
getLM()
Get the mean longitude argument.
|
double |
getLMDot()
Get the mean longitude argument derivative.
|
double |
getLv()
Get the true longitude argument.
|
double |
getLvDot()
Get the true longitude argument derivative.
|
OrbitType |
getType()
Get the orbit type.
|
protected TimeStampedPVCoordinates |
initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters.
|
EquinoctialOrbit |
interpolate(AbsoluteDate date,
Stream<Orbit> sample)
Get an interpolated instance.
|
static double |
meanToEccentric(double lM,
double ex,
double ey)
Computes the eccentric longitude argument from the mean longitude argument.
|
EquinoctialOrbit |
shiftedBy(double dt)
Get a time-shifted orbit.
|
String |
toString()
Returns a string representation of this equinoctial parameters object.
|
static double |
trueToEccentric(double lv,
double ex,
double ey)
Computes the eccentric longitude argument from the true longitude argument.
|
fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, getDate, getFrame, getJacobianWrtCartesian, getJacobianWrtParameters, getKeplerianMeanMotion, getKeplerianPeriod, getMu, getPVCoordinates, getPVCoordinates, getPVCoordinates, hasDerivatives, hasNonKeplerianAcceleration
clone, equals, finalize, getClass, hashCode, notify, notifyAll, wait, wait, wait
interpolate
public EquinoctialOrbit(double a, double ex, double ey, double hx, double hy, double l, PositionAngle type, Frame frame, AbsoluteDate date, double mu) throws IllegalArgumentException
a
- semi-major axis (m)ex
- e cos(ω + Ω), first component of eccentricity vectorey
- e sin(ω + Ω), second component of eccentricity vectorhx
- tan(i/2) cos(Ω), first component of inclination vectorhy
- tan(i/2) sin(Ω), second component of inclination vectorl
- (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)type
- type of longitude argumentframe
- the frame in which the parameters are defined
(must be a pseudo-inertial frame
)date
- date of the orbital parametersmu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if eccentricity is equal to 1 or larger or
if frame is not a pseudo-inertial frame
public EquinoctialOrbit(double a, double ex, double ey, double hx, double hy, double l, double aDot, double exDot, double eyDot, double hxDot, double hyDot, double lDot, PositionAngle type, Frame frame, AbsoluteDate date, double mu) throws IllegalArgumentException
a
- semi-major axis (m)ex
- e cos(ω + Ω), first component of eccentricity vectorey
- e sin(ω + Ω), second component of eccentricity vectorhx
- tan(i/2) cos(Ω), first component of inclination vectorhy
- tan(i/2) sin(Ω), second component of inclination vectorl
- (M or E or v) + ω + Ω, mean, eccentric or true longitude argument (rad)aDot
- semi-major axis derivative (m/s)exDot
- d(e cos(ω + Ω))/dt, first component of eccentricity vector derivativeeyDot
- d(e sin(ω + Ω))/dt, second component of eccentricity vector derivativehxDot
- d(tan(i/2) cos(Ω))/dt, first component of inclination vector derivativehyDot
- d(tan(i/2) sin(Ω))/dt, second component of inclination vector derivativelDot
- d(M or E or v) + ω + Ω)/dr, mean, eccentric or true longitude argument derivative (rad/s)type
- type of longitude argumentframe
- the frame in which the parameters are defined
(must be a pseudo-inertial frame
)date
- date of the orbital parametersmu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if eccentricity is equal to 1 or larger or
if frame is not a pseudo-inertial frame
public EquinoctialOrbit(TimeStampedPVCoordinates pvCoordinates, Frame frame, double mu) throws IllegalArgumentException
The acceleration provided in pvCoordinates
is accessible using
Orbit.getPVCoordinates()
and Orbit.getPVCoordinates(Frame)
. All other methods
use mu
and the position to compute the acceleration, including
shiftedBy(double)
and Orbit.getPVCoordinates(AbsoluteDate, Frame)
.
pvCoordinates
- the position, velocity and accelerationframe
- the frame in which are defined the PVCoordinates
(must be a pseudo-inertial frame
)mu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if eccentricity is equal to 1 or larger or
if frame is not a pseudo-inertial frame
public EquinoctialOrbit(PVCoordinates pvCoordinates, Frame frame, AbsoluteDate date, double mu) throws IllegalArgumentException
The acceleration provided in pvCoordinates
is accessible using
Orbit.getPVCoordinates()
and Orbit.getPVCoordinates(Frame)
. All other methods
use mu
and the position to compute the acceleration, including
shiftedBy(double)
and Orbit.getPVCoordinates(AbsoluteDate, Frame)
.
pvCoordinates
- the position end velocityframe
- the frame in which are defined the PVCoordinates
(must be a pseudo-inertial frame
)date
- date of the orbital parametersmu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if eccentricity is equal to 1 or larger or
if frame is not a pseudo-inertial frame
public EquinoctialOrbit(Orbit op)
op
- orbital parameters to copypublic OrbitType getType()
public double getA()
Note that the semi-major axis is considered negative for hyperbolic orbits.
public double getADot()
Note that the semi-major axis is considered negative for hyperbolic orbits.
If the orbit was created without derivatives, the value returned is Double.NaN
.
getADot
in class Orbit
Orbit.hasDerivatives()
public double getEquinoctialEx()
getEquinoctialEx
in class Orbit
public double getEquinoctialExDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getEquinoctialExDot
in class Orbit
Orbit.hasDerivatives()
public double getEquinoctialEy()
getEquinoctialEy
in class Orbit
public double getEquinoctialEyDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getEquinoctialEyDot
in class Orbit
Orbit.hasDerivatives()
public double getHx()
public double getHxDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getHxDot
in class Orbit
Orbit.hasDerivatives()
public double getHy()
public double getHyDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getHyDot
in class Orbit
Orbit.hasDerivatives()
public double getLv()
public double getLvDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getLvDot
in class Orbit
Orbit.hasDerivatives()
public double getLE()
public double getLEDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getLEDot
in class Orbit
Orbit.hasDerivatives()
public double getLM()
public double getLMDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getLMDot
in class Orbit
Orbit.hasDerivatives()
public double getL(PositionAngle type)
type
- type of the anglepublic double getLDot(PositionAngle type)
type
- type of the anglepublic static double eccentricToTrue(double lE, double ex, double ey)
lE
- = E + ω + Ω eccentric longitude argument (rad)ex
- first component of the eccentricity vectorey
- second component of the eccentricity vectorpublic static double trueToEccentric(double lv, double ex, double ey)
lv
- = v + ω + Ω true longitude argument (rad)ex
- first component of the eccentricity vectorey
- second component of the eccentricity vectorpublic static double meanToEccentric(double lM, double ex, double ey)
lM
- = M + ω + Ω mean longitude argument (rad)ex
- first component of the eccentricity vectorey
- second component of the eccentricity vectorpublic static double eccentricToMean(double lE, double ex, double ey)
lE
- = E + ω + Ω mean longitude argument (rad)ex
- first component of the eccentricity vectorey
- second component of the eccentricity vectorpublic double getE()
public double getEDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getEDot
in class Orbit
Orbit.hasDerivatives()
public double getI()
public double getIDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
getIDot
in class Orbit
Orbit.hasDerivatives()
protected TimeStampedPVCoordinates initPVCoordinates()
initPVCoordinates
in class Orbit
public EquinoctialOrbit shiftedBy(double dt)
The orbit can be slightly shifted to close dates. The shifting model is a Keplerian one if no derivatives are available in the orbit, or Keplerian plus quadratic effect of the non-Keplerian acceleration if derivatives are available. Shifting is not intended as a replacement for proper orbit propagation but should be sufficient for small time shifts or coarse accuracy.
shiftedBy
in interface TimeShiftable<Orbit>
shiftedBy
in class Orbit
dt
- time shift in secondspublic EquinoctialOrbit interpolate(AbsoluteDate date, Stream<Orbit> sample)
Note that the state of the current instance may not be used in the interpolation process, only its type and non interpolable fields are used (for example central attraction coefficient or frame when interpolating orbits). The interpolable fields taken into account are taken only from the states of the sample points. So if the state of the instance must be used, the instance should be included in the sample points.
Note that this method is designed for small samples only (say up to about 10-20 points) so it can be implemented using polynomial interpolation (typically Hermite interpolation). Using too much points may induce Runge's phenomenon and numerical problems (including NaN appearing).
The interpolated instance is created by polynomial Hermite interpolation on equinoctial elements, without derivatives (which means the interpolation falls back to Lagrange interpolation only).
As this implementation of interpolation is polynomial, it should be used only with small samples (about 10-20 points) in order to avoid Runge's phenomenon and numerical problems (including NaN appearing).
If orbit interpolation on large samples is needed, using the Ephemeris
class is a better way than using this
low-level interpolation. The Ephemeris class automatically handles selection of
a neighboring sub-sample with a predefined number of point from a large global sample
in a thread-safe way.
date
- interpolation datesample
- sample points on which interpolation should be doneprotected double[][] computeJacobianMeanWrtCartesian()
Element jacobian[i][j]
is the derivative of parameter i of the orbit with
respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
computeJacobianMeanWrtCartesian
in class Orbit
Orbit.computeJacobianEccentricWrtCartesian()
,
Orbit.computeJacobianTrueWrtCartesian()
protected double[][] computeJacobianEccentricWrtCartesian()
Element jacobian[i][j]
is the derivative of parameter i of the orbit with
respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
computeJacobianEccentricWrtCartesian
in class Orbit
Orbit.computeJacobianMeanWrtCartesian()
,
Orbit.computeJacobianTrueWrtCartesian()
protected double[][] computeJacobianTrueWrtCartesian()
Element jacobian[i][j]
is the derivative of parameter i of the orbit with
respect to Cartesian coordinate j. This means each row correspond to one orbital parameter
whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.
computeJacobianTrueWrtCartesian
in class Orbit
Orbit.computeJacobianMeanWrtCartesian()
,
Orbit.computeJacobianEccentricWrtCartesian()
public void addKeplerContribution(PositionAngle type, double gm, double[] pDot)
This method is used by integration-based propagators to evaluate the part of Keplerian motion to evolution of the orbital state.
addKeplerContribution
in class Orbit
type
- type of the position angle in the stategm
- attraction coefficient to usepDot
- array containing orbital state derivatives to update (the Keplerian
part must be added to the array components, as the array may already
contain some non-zero elements corresponding to non-Keplerian parts)Copyright © 2002-2021 CS GROUP. All rights reserved.