Serializable
, TimeInterpolable<Orbit>
, TimeShiftable<Orbit>
, TimeStamped
, PVCoordinatesProvider
CartesianOrbit
, CircularOrbit
, EquinoctialOrbit
, KeplerianOrbit
public abstract class Orbit extends Object implements TimeStamped, TimeShiftable<Orbit>, TimeInterpolable<Orbit>, Serializable, PVCoordinatesProvider
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.
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.
|
Modifier and Type | Method | Description |
---|---|---|
abstract void |
addKeplerContribution(PositionAngle type,
double gm,
double[] pDot) |
Add the contribution of the Keplerian motion to parameters derivatives
|
protected 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 derivative.
|
abstract double |
getEquinoctialExDot() |
Get the first component of the equinoctial eccentricity vector.
|
abstract double |
getEquinoctialEy() |
Get the second component of the equinoctial eccentricity vector derivative.
|
abstract double |
getEquinoctialEyDot() |
Get the second component of the equinoctial eccentricity vector.
|
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(PositionAngle type,
double[][] jacobian) |
Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
|
void |
getJacobianWrtParameters(PositionAngle 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 |
getMu() |
Get the central acceleration constant.
|
TimeStampedPVCoordinates |
getPVCoordinates() |
Get the
TimeStampedPVCoordinates in definition frame. |
TimeStampedPVCoordinates |
getPVCoordinates(Frame outputFrame) |
Get the
TimeStampedPVCoordinates in a specified frame. |
TimeStampedPVCoordinates |
getPVCoordinates(AbsoluteDate otherDate,
Frame otherFrame) |
Get the
PVCoordinates 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 TimeStampedPVCoordinates |
initPVCoordinates() |
Compute the position/velocity coordinates from the canonical parameters.
|
abstract Orbit |
shiftedBy(double dt) |
Get a time-shifted orbit.
|
clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
interpolate, interpolate
protected Orbit(Frame frame, AbsoluteDate date, double mu) throws IllegalArgumentException
frame
- the frame in which the parameters are defined
(must be a pseudo-inertial frame
)date
- date of the orbital parametersmu
- central attraction coefficient (m^3/s^2)IllegalArgumentException
- if frame is not a pseudo-inertial frame
protected Orbit(TimeStampedPVCoordinates pvCoordinates, Frame frame, double mu) throws IllegalArgumentException
The acceleration provided in pvCoordinates
is accessible using
getPVCoordinates()
and getPVCoordinates(Frame)
. All other methods
use mu
and the position to compute the acceleration, including
shiftedBy(double)
and getPVCoordinates(AbsoluteDate, Frame)
.
pvCoordinates
- the position and velocity in the inertial frameframe
- the frame in which the TimeStampedPVCoordinates
are defined
(must be a pseudo-inertial frame
)mu
- central attraction coefficient (m^3/s^2)IllegalArgumentException
- if frame is not a pseudo-inertial frame
protected static boolean hasNonKeplerianAcceleration(PVCoordinates pva, double mu)
pva
- Cartesian coordinatesmu
- central attraction coefficientpublic abstract OrbitType getType()
public Frame getFrame()
public abstract double getA()
Note that the semi-major axis is considered negative for hyperbolic orbits.
public abstract 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
.
hasDerivatives()
public abstract double getEquinoctialEx()
public abstract double getEquinoctialExDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public abstract double getEquinoctialEy()
public abstract double getEquinoctialEyDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public abstract double getHx()
public abstract double getHxDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public abstract double getHy()
public abstract double getHyDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public abstract double getLE()
public abstract double getLEDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public abstract double getLv()
public abstract double getLvDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public abstract double getLM()
public abstract double getLMDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public abstract double getE()
public abstract double getEDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public abstract double getI()
public abstract double getIDot()
If the orbit was created without derivatives, the value returned is Double.NaN
.
hasDerivatives()
public boolean hasDerivatives()
getADot()
,
getEquinoctialExDot()
,
getEquinoctialEyDot()
,
getHxDot()
,
getHyDot()
,
getLEDot()
,
getLvDot()
,
getLMDot()
,
getEDot()
,
getIDot()
public double getMu()
public double getKeplerianPeriod()
The Keplerian period is computed directly from semi major axis and central acceleration constant.
public double getKeplerianMeanMotion()
The Keplerian mean motion is computed directly from semi major axis and central acceleration constant.
public AbsoluteDate getDate()
getDate
in interface TimeStamped
public TimeStampedPVCoordinates getPVCoordinates(Frame outputFrame) throws OrekitException
TimeStampedPVCoordinates
in a specified frame.outputFrame
- frame in which the position/velocity coordinates shall be computedOrekitException
- if transformation between frames cannot be computedgetPVCoordinates()
public TimeStampedPVCoordinates getPVCoordinates(AbsoluteDate otherDate, Frame otherFrame) throws OrekitException
PVCoordinates
of the body in the selected frame.getPVCoordinates
in interface PVCoordinatesProvider
otherDate
- current dateotherFrame
- the frame where to define the positionOrekitException
- if position cannot be computed in given framepublic TimeStampedPVCoordinates getPVCoordinates()
TimeStampedPVCoordinates
in definition frame.getPVCoordinates(Frame)
protected abstract TimeStampedPVCoordinates initPVCoordinates()
public abstract Orbit 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>
dt
- time shift in secondspublic void getJacobianWrtCartesian(PositionAngle type, double[][] jacobian)
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.
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 modifiedpublic void getJacobianWrtParameters(PositionAngle type, double[][] jacobian)
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.
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 modifiedprotected abstract 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.
computeJacobianEccentricWrtCartesian()
,
computeJacobianTrueWrtCartesian()
protected abstract 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.
computeJacobianMeanWrtCartesian()
,
computeJacobianTrueWrtCartesian()
protected abstract 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.
computeJacobianMeanWrtCartesian()
,
computeJacobianEccentricWrtCartesian()
public abstract 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.
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)protected static void fillHalfRow(double a, Vector3D v, double[] row, int j)
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)protected static void fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double[] row, int j)
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)protected static void fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double[] row, int j)
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)protected static void fillHalfRow(double a1, Vector3D v1, double a2, Vector3D v2, double a3, Vector3D v3, double a4, Vector3D v4, double[] row, int j)
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)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)
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)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)
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)Copyright © 2002-2018 CS Systèmes d'information. All rights reserved.