Class Orbit
- java.lang.Object
-
- org.orekit.orbits.Orbit
-
- All Implemented Interfaces:
Serializable
,TimeShiftable<Orbit>
,TimeStamped
,PVCoordinatesProvider
- Direct Known Subclasses:
CartesianOrbit
,CircularOrbit
,EquinoctialOrbit
,KeplerianOrbit
public abstract class Orbit extends Object implements TimeStamped, TimeShiftable<Orbit>, Serializable, PVCoordinatesProvider
This class handles orbital parameters.For user convenience, both the Cartesian and the equinoctial elements are provided by this class, regardless of the canonical representation implemented in the derived class (which may be classical Keplerian elements for example).
The parameters are defined in a frame specified by the user. It is important to make sure this frame is consistent: it probably is inertial and centered on the central body. This information is used for example by some force models.
Instance of this class are guaranteed to be immutable.
- Author:
- Luc Maisonobe, Guylaine Prat, Fabien Maussion, Véronique Pommier-Maurussane
- See Also:
- Serialized Form
-
-
Constructor Summary
Constructors Modifier Constructor Description protected
Orbit(Frame frame, AbsoluteDate date, double mu)
Default constructor.protected
Orbit(TimeStampedPVCoordinates pvCoordinates, Frame frame, double mu)
Set the orbit from Cartesian parameters.
-
Method Summary
All Methods Static Methods Instance Methods Abstract Methods Concrete Methods Modifier and Type Method Description abstract void
addKeplerContribution(PositionAngleType type, double gm, double[] pDot)
Add the contribution of the Keplerian motion to parameters derivativesprotected abstract double[][]
computeJacobianEccentricWrtCartesian()
Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.protected abstract double[][]
computeJacobianMeanWrtCartesian()
Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.protected abstract double[][]
computeJacobianTrueWrtCartesian()
Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.protected static void
fillHalfRow(double a, Vector3D v, double[] row, int j)
Fill a Jacobian half row with a single vector.protected static void
fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.protected static void
fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.protected static void
fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double a4, Vector3D v4, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.protected static void
fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double a4, Vector3D v4, double a5, Vector3D v5, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.protected static void
fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double a4, Vector3D v4, double a5, Vector3D v5, double a6, Vector3D v6, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.abstract double
getA()
Get the semi-major axis.abstract double
getADot()
Get the semi-major axis derivative.AbsoluteDate
getDate()
Get the date of orbital parameters.abstract double
getE()
Get the eccentricity.abstract double
getEDot()
Get the eccentricity derivative.abstract double
getEquinoctialEx()
Get the first component of the equinoctial eccentricity vector.abstract double
getEquinoctialExDot()
Get the first component of the equinoctial eccentricity vector derivative.abstract double
getEquinoctialEy()
Get the second component of the equinoctial eccentricity vector.abstract double
getEquinoctialEyDot()
Get the second component of the equinoctial eccentricity vector derivative.Frame
getFrame()
Get the frame in which the orbital parameters are defined.abstract double
getHx()
Get the first component of the inclination vector.abstract double
getHxDot()
Get the first component of the inclination vector derivative.abstract double
getHy()
Get the second component of the inclination vector.abstract double
getHyDot()
Get the second component of the inclination vector derivative.abstract double
getI()
Get the inclination.abstract double
getIDot()
Get the inclination derivative.void
getJacobianWrtCartesian(PositionAngleType type, double[][] jacobian)
Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.void
getJacobianWrtParameters(PositionAngleType type, double[][] jacobian)
Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.double
getKeplerianMeanMotion()
Get the Keplerian mean motion.double
getKeplerianPeriod()
Get the Keplerian period.abstract double
getLE()
Get the eccentric longitude argument.abstract double
getLEDot()
Get the eccentric longitude argument derivative.abstract double
getLM()
Get the mean longitude argument.abstract double
getLMDot()
Get the mean longitude argument derivative.abstract double
getLv()
Get the true longitude argument.abstract double
getLvDot()
Get the true longitude argument derivative.double
getMeanAnomalyDotWrtA()
Get the derivative of the mean anomaly with respect to the semi major axis.double
getMu()
Get the central acceleration constant.Vector3D
getPosition()
Get the position in definition frame.Vector3D
getPosition(Frame outputFrame)
Get the position in a specified frame.Vector3D
getPosition(AbsoluteDate otherDate, Frame otherFrame)
Get the position of the body in the selected frame.TimeStampedPVCoordinates
getPVCoordinates()
Get theTimeStampedPVCoordinates
in definition frame.TimeStampedPVCoordinates
getPVCoordinates(Frame outputFrame)
Get theTimeStampedPVCoordinates
in a specified frame.TimeStampedPVCoordinates
getPVCoordinates(AbsoluteDate otherDate, Frame otherFrame)
Get thePVCoordinates
of the body in the selected frame.abstract OrbitType
getType()
Get the orbit type.boolean
hasDerivatives()
Check if orbit includes derivatives.protected static boolean
hasNonKeplerianAcceleration(PVCoordinates pva, double mu)
Check if Cartesian coordinates include non-Keplerian acceleration.protected abstract Vector3D
initPosition()
Compute the position coordinates from the canonical parameters.protected abstract TimeStampedPVCoordinates
initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters.boolean
isElliptical()
Returns true if and only if the orbit is elliptical i.e.abstract Orbit
shiftedBy(double dt)
Get a time-shifted orbit.-
Methods inherited from class java.lang.Object
clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
-
Methods inherited from interface org.orekit.time.TimeStamped
durationFrom
-
-
-
-
Constructor Detail
-
Orbit
protected Orbit(Frame frame, AbsoluteDate date, double mu) throws IllegalArgumentException
Default constructor. Build a new instance with arbitrary default elements.- Parameters:
frame
- the frame in which the parameters are defined (must be apseudo-inertial frame
)date
- date of the orbital parametersmu
- central attraction coefficient (m^3/s^2)- Throws:
IllegalArgumentException
- if frame is not apseudo-inertial frame
-
Orbit
protected Orbit(TimeStampedPVCoordinates pvCoordinates, Frame frame, double mu) throws IllegalArgumentException
Set the orbit from Cartesian parameters.The acceleration provided in
pvCoordinates
is accessible usinggetPVCoordinates()
andgetPVCoordinates(Frame)
. All other methods usemu
and the position to compute the acceleration, includingshiftedBy(double)
andgetPVCoordinates(AbsoluteDate, Frame)
.- Parameters:
pvCoordinates
- the position and velocity in the inertial frameframe
- the frame in which theTimeStampedPVCoordinates
are defined (must be apseudo-inertial frame
)mu
- central attraction coefficient (m^3/s^2)- Throws:
IllegalArgumentException
- if frame is not apseudo-inertial frame
-
-
Method Detail
-
hasNonKeplerianAcceleration
protected static boolean hasNonKeplerianAcceleration(PVCoordinates pva, double mu)
Check if Cartesian coordinates include non-Keplerian acceleration.- Parameters:
pva
- Cartesian coordinatesmu
- central attraction coefficient- Returns:
- true if Cartesian coordinates include non-Keplerian acceleration
-
isElliptical
public boolean isElliptical()
Returns true if and only if the orbit is elliptical i.e. has a non-negative semi-major axis.- Returns:
- true if getA() is strictly greater than 0
- Since:
- 12.0
-
getType
public abstract OrbitType getType()
Get the orbit type.- Returns:
- orbit type
-
getFrame
public Frame getFrame()
Get the frame in which the orbital parameters are defined.- Returns:
- frame in which the orbital parameters are defined
-
getA
public abstract double getA()
Get the semi-major axis.Note that the semi-major axis is considered negative for hyperbolic orbits.
- Returns:
- semi-major axis (m)
-
getADot
public abstract double getADot()
Get the semi-major axis derivative.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
.- Returns:
- semi-major axis derivative (m/s)
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getEquinoctialEx
public abstract double getEquinoctialEx()
Get the first component of the equinoctial eccentricity vector.- Returns:
- first component of the equinoctial eccentricity vector
-
getEquinoctialExDot
public abstract double getEquinoctialExDot()
Get the first component of the equinoctial eccentricity vector derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- first component of the equinoctial eccentricity vector derivative
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getEquinoctialEy
public abstract double getEquinoctialEy()
Get the second component of the equinoctial eccentricity vector.- Returns:
- second component of the equinoctial eccentricity vector
-
getEquinoctialEyDot
public abstract double getEquinoctialEyDot()
Get the second component of the equinoctial eccentricity vector derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- second component of the equinoctial eccentricity vector derivative
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getHx
public abstract double getHx()
Get the first component of the inclination vector.- Returns:
- first component of the inclination vector
-
getHxDot
public abstract double getHxDot()
Get the first component of the inclination vector derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- first component of the inclination vector derivative
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getHy
public abstract double getHy()
Get the second component of the inclination vector.- Returns:
- second component of the inclination vector
-
getHyDot
public abstract double getHyDot()
Get the second component of the inclination vector derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- second component of the inclination vector derivative
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getLE
public abstract double getLE()
Get the eccentric longitude argument.- Returns:
- E + ω + Ω eccentric longitude argument (rad)
-
getLEDot
public abstract double getLEDot()
Get the eccentric longitude argument derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- d(E + ω + Ω)/dt eccentric longitude argument derivative (rad/s)
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getLv
public abstract double getLv()
Get the true longitude argument.- Returns:
- v + ω + Ω true longitude argument (rad)
-
getLvDot
public abstract double getLvDot()
Get the true longitude argument derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- d(v + ω + Ω)/dt true longitude argument derivative (rad/s)
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getLM
public abstract double getLM()
Get the mean longitude argument.- Returns:
- M + ω + Ω mean longitude argument (rad)
-
getLMDot
public abstract double getLMDot()
Get the mean longitude argument derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- d(M + ω + Ω)/dt mean longitude argument derivative (rad/s)
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getE
public abstract double getE()
Get the eccentricity.- Returns:
- eccentricity
-
getEDot
public abstract double getEDot()
Get the eccentricity derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- eccentricity derivative
- Since:
- 9.0
- See Also:
hasDerivatives()
-
getI
public abstract double getI()
Get the inclination.- Returns:
- inclination (rad)
-
getIDot
public abstract double getIDot()
Get the inclination derivative.If the orbit was created without derivatives, the value returned is
Double.NaN
.- Returns:
- inclination derivative (rad/s)
- Since:
- 9.0
- See Also:
hasDerivatives()
-
hasDerivatives
public boolean hasDerivatives()
Check if orbit includes derivatives.- Returns:
- true if orbit includes derivatives
- Since:
- 9.0
- See Also:
getADot()
,getEquinoctialExDot()
,getEquinoctialEyDot()
,getHxDot()
,getHyDot()
,getLEDot()
,getLvDot()
,getLMDot()
,getEDot()
,getIDot()
-
getMu
public double getMu()
Get the central acceleration constant.- Returns:
- central acceleration constant
-
getKeplerianPeriod
public double getKeplerianPeriod()
Get the Keplerian period.The Keplerian period is computed directly from semi major axis and central acceleration constant.
- Returns:
- Keplerian period in seconds, or positive infinity for hyperbolic orbits
-
getKeplerianMeanMotion
public double getKeplerianMeanMotion()
Get the Keplerian mean motion.The Keplerian mean motion is computed directly from semi major axis and central acceleration constant.
- Returns:
- Keplerian mean motion in radians per second
-
getMeanAnomalyDotWrtA
public double getMeanAnomalyDotWrtA()
Get the derivative of the mean anomaly with respect to the semi major axis.- Returns:
- derivative of the mean anomaly with respect to the semi major axis
-
getDate
public AbsoluteDate getDate()
Get the date of orbital parameters.- Specified by:
getDate
in interfaceTimeStamped
- Returns:
- date of the orbital parameters
-
getPVCoordinates
public TimeStampedPVCoordinates getPVCoordinates(Frame outputFrame)
Get theTimeStampedPVCoordinates
in a specified frame.- Parameters:
outputFrame
- frame in which the position/velocity coordinates shall be computed- Returns:
- pvCoordinates in the specified output frame
- See Also:
getPVCoordinates()
-
getPVCoordinates
public TimeStampedPVCoordinates getPVCoordinates(AbsoluteDate otherDate, Frame otherFrame)
Get thePVCoordinates
of the body in the selected frame.- Specified by:
getPVCoordinates
in interfacePVCoordinatesProvider
- Parameters:
otherDate
- current dateotherFrame
- the frame where to define the position- Returns:
- time-stamped position/velocity of the body (m and m/s)
-
getPosition
public Vector3D getPosition(AbsoluteDate otherDate, Frame otherFrame)
Get the position of the body in the selected frame.- Specified by:
getPosition
in interfacePVCoordinatesProvider
- Parameters:
otherDate
- current dateotherFrame
- the frame where to define the position- Returns:
- position of the body (m and)
-
getPosition
public Vector3D getPosition(Frame outputFrame)
Get the position in a specified frame.- Parameters:
outputFrame
- frame in which the position coordinates shall be computed- Returns:
- position in the specified output frame
- Since:
- 12.0
- See Also:
getPosition()
-
getPosition
public Vector3D getPosition()
Get the position in definition frame.- Returns:
- position in the definition frame
- Since:
- 12.0
- See Also:
getPVCoordinates()
-
getPVCoordinates
public TimeStampedPVCoordinates getPVCoordinates()
Get theTimeStampedPVCoordinates
in definition frame.- Returns:
- pvCoordinates in the definition frame
- See Also:
getPVCoordinates(Frame)
-
initPosition
protected abstract Vector3D initPosition()
Compute the position coordinates from the canonical parameters.- Returns:
- computed position coordinates
- Since:
- 12.0
-
initPVCoordinates
protected abstract TimeStampedPVCoordinates initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters.- Returns:
- computed position/velocity coordinates
-
shiftedBy
public abstract Orbit shiftedBy(double dt)
Get a time-shifted orbit.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.
- Specified by:
shiftedBy
in interfaceTimeShiftable<Orbit>
- Parameters:
dt
- time shift in seconds- Returns:
- a new orbit, shifted with respect to the instance (which is immutable)
-
getJacobianWrtCartesian
public void getJacobianWrtCartesian(PositionAngleType type, double[][] jacobian)
Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.Element
jacobian[i][j]
is the derivative of parameter i of the orbit with respect to Cartesian coordinate j. This means each row corresponds to one orbital parameter whereas columns 0 to 5 correspond to the Cartesian coordinates x, y, z, xDot, yDot and zDot.- Parameters:
type
- type of the position angle to usejacobian
- placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix is larger than 6x6, only the 6x6 upper left corner will be modified
-
getJacobianWrtParameters
public void getJacobianWrtParameters(PositionAngleType type, double[][] jacobian)
Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.Element
jacobian[i][j]
is the derivative of Cartesian coordinate i of the orbit with respect to orbital parameter j. This means each row corresponds to one Cartesian coordinate x, y, z, xdot, ydot, zdot whereas columns 0 to 5 correspond to the orbital parameters.- Parameters:
type
- type of the position angle to usejacobian
- placeholder 6x6 (or larger) matrix to be filled with the Jacobian, if matrix is larger than 6x6, only the 6x6 upper left corner will be modified
-
computeJacobianMeanWrtCartesian
protected abstract double[][] computeJacobianMeanWrtCartesian()
Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.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.The array returned by this method will not be modified.
- Returns:
- 6x6 Jacobian matrix
- See Also:
computeJacobianEccentricWrtCartesian()
,computeJacobianTrueWrtCartesian()
-
computeJacobianEccentricWrtCartesian
protected abstract double[][] computeJacobianEccentricWrtCartesian()
Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.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.The array returned by this method will not be modified.
- Returns:
- 6x6 Jacobian matrix
- See Also:
computeJacobianMeanWrtCartesian()
,computeJacobianTrueWrtCartesian()
-
computeJacobianTrueWrtCartesian
protected abstract double[][] computeJacobianTrueWrtCartesian()
Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.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.The array returned by this method will not be modified.
- Returns:
- 6x6 Jacobian matrix
- See Also:
computeJacobianMeanWrtCartesian()
,computeJacobianEccentricWrtCartesian()
-
addKeplerContribution
public abstract void addKeplerContribution(PositionAngleType type, double gm, double[] pDot)
Add the contribution of the Keplerian motion to parameters derivativesThis method is used by integration-based propagators to evaluate the part of Keplerian motion to evolution of the orbital state.
- Parameters:
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)
-
fillHalfRow
protected static void fillHalfRow(double a, Vector3D v, double[] row, int j)
Fill a Jacobian half row with a single vector.- Parameters:
a
- coefficient of the vectorv
- vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected static void fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected static void fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectora3
- coefficient of the third vectorv3
- third vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected static void fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double a4, Vector3D v4, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectora3
- coefficient of the third vectorv3
- third vectora4
- coefficient of the fourth vectorv4
- fourth vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected static void fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double a4, Vector3D v4, double a5, Vector3D v5, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectora3
- coefficient of the third vectorv3
- third vectora4
- coefficient of the fourth vectorv4
- fourth vectora5
- coefficient of the fifth vectorv5
- fifth vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
fillHalfRow
protected static void fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double a4, Vector3D v4, double a5, Vector3D v5, double a6, Vector3D v6, double[] row, int j)
Fill a Jacobian half row with a linear combination of vectors.- Parameters:
a1
- coefficient of the first vectorv1
- first vectora2
- coefficient of the second vectorv2
- second vectora3
- coefficient of the third vectorv3
- third vectora4
- coefficient of the fourth vectorv4
- fourth vectora5
- coefficient of the fifth vectorv5
- fifth vectora6
- coefficient of the sixth vectorv6
- sixth vectorrow
- Jacobian matrix rowj
- index of the first element to set (row[j], row[j+1] and row[j+2] will all be set)
-
-