public abstract class FieldOrbit<T extends org.hipparchus.RealFieldElement<T>> extends Object implements FieldPVCoordinatesProvider<T>, FieldTimeStamped<T>, FieldTimeShiftable<FieldOrbit<T>,T>, FieldTimeInterpolable<FieldOrbit<T>,T>
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 and Description |
---|---|
protected |
FieldOrbit(Frame frame,
FieldAbsoluteDate<T> date,
double mu)
Default constructor.
|
protected |
FieldOrbit(TimeStampedFieldPVCoordinates<T> FieldPVCoordinates,
Frame frame,
double mu)
Set the orbit from Cartesian parameters.
|
Modifier and Type | Method and Description |
---|---|
abstract void |
addKeplerContribution(PositionAngle type,
double gm,
T[] pDot)
Add the contribution of the Keplerian motion to parameters derivatives
|
protected abstract T[][] |
computeJacobianEccentricWrtCartesian()
Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
|
protected abstract T[][] |
computeJacobianMeanWrtCartesian()
Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
|
protected abstract T[][] |
computeJacobianTrueWrtCartesian()
Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
|
protected void |
fillHalfRow(T a,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v,
T[] row,
int j)
Fill a Jacobian half row with a single vector.
|
protected void |
fillHalfRow(T a1,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1,
T a2,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2,
T[] row,
int j)
Fill a Jacobian half row with a linear combination of vectors.
|
protected void |
fillHalfRow(T a1,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1,
T a2,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2,
T a3,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3,
T[] row,
int j)
Fill a Jacobian half row with a linear combination of vectors.
|
protected void |
fillHalfRow(T a1,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1,
T a2,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2,
T a3,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3,
T a4,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4,
T[] row,
int j)
Fill a Jacobian half row with a linear combination of vectors.
|
protected void |
fillHalfRow(T a1,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1,
T a2,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2,
T a3,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3,
T a4,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4,
T a5,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v5,
T[] row,
int j)
Fill a Jacobian half row with a linear combination of vectors.
|
protected void |
fillHalfRow(T a1,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1,
T a2,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2,
T a3,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3,
T a4,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4,
T a5,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v5,
T a6,
org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v6,
T[] row,
int j)
Fill a Jacobian half row with a linear combination of vectors.
|
abstract T |
getA()
Get the semi-major axis.
|
abstract T |
getADot()
Get the semi-major axis derivative.
|
FieldAbsoluteDate<T> |
getDate()
Get the date of orbital parameters.
|
abstract T |
getE()
Get the eccentricity.
|
abstract T |
getEDot()
Get the eccentricity derivative.
|
abstract T |
getEquinoctialEx()
Get the first component of the equinoctial eccentricity vector.
|
abstract T |
getEquinoctialExDot()
Get the first component of the equinoctial eccentricity vector.
|
abstract T |
getEquinoctialEy()
Get the second component of the equinoctial eccentricity vector.
|
abstract T |
getEquinoctialEyDot()
Get the second component of the equinoctial eccentricity vector.
|
Frame |
getFrame()
Get the frame in which the orbital parameters are defined.
|
abstract T |
getHx()
Get the first component of the inclination vector.
|
abstract T |
getHxDot()
Get the first component of the inclination vector derivative.
|
abstract T |
getHy()
Get the second component of the inclination vector.
|
abstract T |
getHyDot()
Get the second component of the inclination vector derivative.
|
abstract T |
getI()
Get the inclination.
|
abstract T |
getIDot()
Get the inclination derivative.
|
void |
getJacobianWrtCartesian(PositionAngle type,
T[][] jacobian)
Compute the Jacobian of the orbital parameters with respect to the Cartesian parameters.
|
void |
getJacobianWrtParameters(PositionAngle type,
T[][] jacobian)
Compute the Jacobian of the Cartesian parameters with respect to the orbital parameters.
|
T |
getKeplerianMeanMotion()
Get the Keplerian mean motion.
|
T |
getKeplerianPeriod()
Get the Keplerian period.
|
abstract T |
getLE()
Get the eccentric longitude argument.
|
abstract T |
getLEDot()
Get the eccentric longitude argument derivative.
|
abstract T |
getLM()
Get the mean longitude argument.
|
abstract T |
getLMDot()
Get the mean longitude argument derivative.
|
abstract T |
getLv()
Get the true longitude argument.
|
abstract T |
getLvDot()
Get the true longitude argument derivative.
|
double |
getMu() |
TimeStampedFieldPVCoordinates<T> |
getPVCoordinates()
Get the
TimeStampedPVCoordinates in definition frame. |
TimeStampedFieldPVCoordinates<T> |
getPVCoordinates(FieldAbsoluteDate<T> otherDate,
Frame otherFrame)
Get the
FieldPVCoordinates of the body in the selected frame. |
TimeStampedFieldPVCoordinates<T> |
getPVCoordinates(Frame outputFrame)
Get the
TimeStampedPVCoordinates in a specified frame. |
abstract OrbitType |
getType()
Get the orbit type.
|
abstract boolean |
hasDerivatives()
Check if orbit includes derivatives.
|
protected static <T extends org.hipparchus.RealFieldElement<T>> |
hasNonKeplerianAcceleration(FieldPVCoordinates<T> pva,
double mu)
Check if Cartesian coordinates include non-Keplerian acceleration.
|
protected abstract TimeStampedFieldPVCoordinates<T> |
initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters.
|
abstract FieldOrbit<T> |
shiftedBy(T dt)
Get a time-shifted orbit.
|
abstract Orbit |
toOrbit()
Transforms the FieldOrbit instance into an Orbit instance.
|
clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
shiftedBy
interpolate, interpolate
protected FieldOrbit(Frame frame, FieldAbsoluteDate<T> 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 FieldOrbit(TimeStampedFieldPVCoordinates<T> FieldPVCoordinates, Frame frame, double mu) throws IllegalArgumentException
The acceleration provided in FieldPVCoordinates
is accessible using
getPVCoordinates()
and getPVCoordinates(Frame)
. All other methods
use mu
and the position to compute the acceleration, including
shiftedBy(RealFieldElement)
and getPVCoordinates(FieldAbsoluteDate, Frame)
.
FieldPVCoordinates
- 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 <T extends org.hipparchus.RealFieldElement<T>> boolean hasNonKeplerianAcceleration(FieldPVCoordinates<T> pva, double mu)
T
- type of the field elementspva
- Cartesian coordinatesmu
- central attraction coefficientpublic abstract OrbitType getType()
public Frame getFrame()
public abstract Orbit toOrbit()
public abstract T getA()
Note that the semi-major axis is considered negative for hyperbolic orbits.
public abstract T getADot()
Note that the semi-major axis is considered negative for hyperbolic orbits.
If the orbit was created without derivatives, the value returned is null.
public abstract T getEquinoctialEx()
public abstract T getEquinoctialExDot()
If the orbit was created without derivatives, the value returned is null.
public abstract T getEquinoctialEy()
public abstract T getEquinoctialEyDot()
If the orbit was created without derivatives, the value returned is null.
public abstract T getHx()
public abstract T getHxDot()
If the orbit was created without derivatives, the value returned is null.
public abstract T getHy()
public abstract T getHyDot()
public abstract T getLE()
public abstract T getLEDot()
If the orbit was created without derivatives, the value returned is null.
public abstract T getLv()
public abstract T getLvDot()
If the orbit was created without derivatives, the value returned is null.
public abstract T getLM()
public abstract T getLMDot()
If the orbit was created without derivatives, the value returned is null.
public abstract T getE()
public abstract T getEDot()
If the orbit was created without derivatives, the value returned is null.
public abstract T getI()
If the orbit was created without derivatives, the value returned is null.
public abstract T getIDot()
public abstract boolean hasDerivatives()
getADot()
,
getEquinoctialExDot()
,
getEquinoctialEyDot()
,
getHxDot()
,
getHyDot()
,
getLEDot()
,
getLvDot()
,
getLMDot()
,
getEDot()
,
getIDot()
public double getMu()
public T getKeplerianPeriod()
The Keplerian period is computed directly from semi major axis and central acceleration constant.
public T getKeplerianMeanMotion()
The Keplerian mean motion is computed directly from semi major axis and central acceleration constant.
public FieldAbsoluteDate<T> getDate()
getDate
in interface FieldTimeStamped<T extends org.hipparchus.RealFieldElement<T>>
public TimeStampedFieldPVCoordinates<T> 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 TimeStampedFieldPVCoordinates<T> getPVCoordinates(FieldAbsoluteDate<T> otherDate, Frame otherFrame) throws OrekitException
FieldPVCoordinates
of the body in the selected frame.getPVCoordinates
in interface FieldPVCoordinatesProvider<T extends org.hipparchus.RealFieldElement<T>>
otherDate
- current dateotherFrame
- the frame where to define the positionOrekitException
- if position cannot be computed in given framepublic TimeStampedFieldPVCoordinates<T> getPVCoordinates()
TimeStampedPVCoordinates
in definition frame.getPVCoordinates(Frame)
protected abstract TimeStampedFieldPVCoordinates<T> initPVCoordinates()
public abstract FieldOrbit<T> shiftedBy(T dt)
The orbit can be slightly shifted to close dates. This shift is based on a simple Keplerian model. It is not intended as a replacement for proper orbit and attitude propagation but should be sufficient for small time shifts or coarse accuracy.
shiftedBy
in interface FieldTimeShiftable<FieldOrbit<T extends org.hipparchus.RealFieldElement<T>>,T extends org.hipparchus.RealFieldElement<T>>
dt
- time shift in secondspublic void getJacobianWrtCartesian(PositionAngle type, T[][] 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, T[][] 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 T[][] 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 T[][] 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 T[][] 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, T[] 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 void fillHalfRow(T a, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v, T[] 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 void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T[] 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 void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T[] 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 void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T[] 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 void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T a5, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v5, T[] 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 void fillHalfRow(T a1, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v1, T a2, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v2, T a3, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v3, T a4, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v4, T a5, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v5, T a6, org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> v6, T[] 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-2017 CS Systèmes d'information. All rights reserved.