T
- type of the field elementspublic class FieldCartesianOrbit<T extends org.hipparchus.CalculusFieldElement<T>> extends FieldOrbit<T>
The parameters used internally are the Cartesian coordinates:
PVCoordinates
.
Note that the implementation of this class delegates all non-Cartesian related
computations (getA()
, getEquinoctialEx()
, ...) to an underlying
instance of the EquinoctialOrbit
class. This implies that using this class
only for analytical computations which are always based on non-Cartesian
parameters is perfectly possible but somewhat sub-optimal.
The instance CartesianOrbit
is guaranteed to be immutable.
Orbit
,
KeplerianOrbit
,
CircularOrbit
,
EquinoctialOrbit
Constructor and Description |
---|
FieldCartesianOrbit(org.hipparchus.Field<T> field,
CartesianOrbit op)
Constructor from Field and CartesianOrbit.
|
FieldCartesianOrbit(org.hipparchus.Field<T> field,
Orbit op)
Constructor from Field and Orbit.
|
FieldCartesianOrbit(FieldOrbit<T> op)
Constructor from any kind of orbital parameters.
|
FieldCartesianOrbit(FieldPVCoordinates<T> pvaCoordinates,
Frame frame,
FieldAbsoluteDate<T> date,
T mu)
Constructor from Cartesian parameters.
|
FieldCartesianOrbit(TimeStampedFieldPVCoordinates<T> pvaCoordinates,
Frame frame,
T mu)
Constructor from Cartesian parameters.
|
Modifier and Type | Method and Description |
---|---|
void |
addKeplerContribution(PositionAngleType type,
T gm,
T[] pDot)
Add the contribution of the Keplerian motion to parameters derivatives
|
protected T[][] |
computeJacobianEccentricWrtCartesian()
Compute the Jacobian of the orbital parameters with eccentric angle with respect to the Cartesian parameters.
|
protected T[][] |
computeJacobianMeanWrtCartesian()
Compute the Jacobian of the orbital parameters with mean angle with respect to the Cartesian parameters.
|
protected T[][] |
computeJacobianTrueWrtCartesian()
Compute the Jacobian of the orbital parameters with true angle with respect to the Cartesian parameters.
|
T |
getA()
Get the semi-major axis.
|
T |
getADot()
Get the semi-major axis derivative.
|
T |
getE()
Get the eccentricity.
|
T |
getEDot()
Get the eccentricity derivative.
|
T |
getEquinoctialEx()
Get the first component of the equinoctial eccentricity vector.
|
T |
getEquinoctialExDot()
Get the first component of the equinoctial eccentricity vector.
|
T |
getEquinoctialEy()
Get the second component of the equinoctial eccentricity vector.
|
T |
getEquinoctialEyDot()
Get the second component of the equinoctial eccentricity vector.
|
T |
getHx()
Get the first component of the inclination vector.
|
T |
getHxDot()
Get the first component of the inclination vector derivative.
|
T |
getHy()
Get the second component of the inclination vector.
|
T |
getHyDot()
Get the second component of the inclination vector derivative.
|
T |
getI()
Get the inclination.
|
T |
getIDot()
Get the inclination derivative.
|
T |
getLE()
Get the eccentric longitude argument.
|
T |
getLEDot()
Get the eccentric longitude argument derivative.
|
T |
getLM()
Get the mean longitude argument.
|
T |
getLMDot()
Get the mean longitude argument derivative.
|
T |
getLv()
Get the true longitude argument.
|
T |
getLvDot()
Get the true longitude argument derivative.
|
OrbitType |
getType()
Get the orbit type.
|
boolean |
hasDerivatives()
Check if orbit includes derivatives.
|
protected org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> |
initPosition()
Compute the position coordinates from the canonical parameters.
|
protected TimeStampedFieldPVCoordinates<T> |
initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters.
|
FieldCartesianOrbit<T> |
shiftedBy(double dt)
Get a time-shifted instance.
|
FieldCartesianOrbit<T> |
shiftedBy(T dt)
Get a time-shifted orbit.
|
CartesianOrbit |
toOrbit()
Transforms the FieldOrbit instance into an Orbit instance.
|
String |
toString()
Returns a string representation of this Orbit object.
|
fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, getDate, getField, getFrame, getJacobianWrtCartesian, getJacobianWrtParameters, getKeplerianMeanMotion, getKeplerianPeriod, getMeanAnomalyDotWrtA, getMu, getOne, getPosition, getPosition, getPVCoordinates, getPVCoordinates, getPVCoordinates, getZero, hasNonKeplerianAcceleration, isElliptical
clone, equals, finalize, getClass, hashCode, notify, notifyAll, wait, wait, wait
getPosition
durationFrom
public FieldCartesianOrbit(TimeStampedFieldPVCoordinates<T> pvaCoordinates, Frame frame, T mu) throws IllegalArgumentException
The acceleration provided in pvCoordinates
is accessible using
FieldOrbit.getPVCoordinates()
and FieldOrbit.getPVCoordinates(Frame)
. All other methods
use mu
and the position to compute the acceleration, including
shiftedBy(CalculusFieldElement)
and FieldOrbit.getPVCoordinates(FieldAbsoluteDate, Frame)
.
pvaCoordinates
- the position, velocity and acceleration of the satellite.frame
- the frame in which the PVCoordinates
are defined
(must be a pseudo-inertial frame
)mu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if frame is not a pseudo-inertial frame
public FieldCartesianOrbit(FieldPVCoordinates<T> pvaCoordinates, Frame frame, FieldAbsoluteDate<T> date, T mu) throws IllegalArgumentException
The acceleration provided in pvCoordinates
is accessible using
FieldOrbit.getPVCoordinates()
and FieldOrbit.getPVCoordinates(Frame)
. All other methods
use mu
and the position to compute the acceleration, including
shiftedBy(CalculusFieldElement)
and FieldOrbit.getPVCoordinates(FieldAbsoluteDate, Frame)
.
pvaCoordinates
- the position and velocity of the satellite.frame
- the frame in which the PVCoordinates
are defined
(must be a pseudo-inertial frame
)date
- date of the orbital parametersmu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if frame is not a pseudo-inertial frame
public FieldCartesianOrbit(FieldOrbit<T> op)
op
- orbital parameters to copypublic FieldCartesianOrbit(org.hipparchus.Field<T> field, CartesianOrbit op)
Build a FieldCartesianOrbit from non-Field CartesianOrbit.
field
- CalculusField to base object onop
- non-field orbit with only "constant" termspublic OrbitType getType()
getType
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getA()
Note that the semi-major axis is considered negative for hyperbolic orbits.
getA
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public 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.
getADot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getE()
getE
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getEDot()
If the orbit was created without derivatives, the value returned is null.
getEDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getI()
If the orbit was created without derivatives, the value returned is null.
getI
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getIDot()
getIDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getEquinoctialEx()
getEquinoctialEx
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getEquinoctialExDot()
If the orbit was created without derivatives, the value returned is null.
getEquinoctialExDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getEquinoctialEy()
getEquinoctialEy
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getEquinoctialEyDot()
If the orbit was created without derivatives, the value returned is null.
getEquinoctialEyDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getHx()
getHx
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getHxDot()
If the orbit was created without derivatives, the value returned is null.
getHxDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getHy()
getHy
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getHyDot()
getHyDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getLv()
getLv
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getLvDot()
If the orbit was created without derivatives, the value returned is null.
getLvDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getLE()
getLE
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getLEDot()
If the orbit was created without derivatives, the value returned is null.
getLEDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getLM()
getLM
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public T getLMDot()
If the orbit was created without derivatives, the value returned is null.
getLMDot
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public boolean hasDerivatives()
hasDerivatives
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
FieldOrbit.getADot()
,
FieldOrbit.getEquinoctialExDot()
,
FieldOrbit.getEquinoctialEyDot()
,
FieldOrbit.getHxDot()
,
FieldOrbit.getHyDot()
,
FieldOrbit.getLEDot()
,
FieldOrbit.getLvDot()
,
FieldOrbit.getLMDot()
,
FieldOrbit.getEDot()
,
FieldOrbit.getIDot()
protected org.hipparchus.geometry.euclidean.threed.FieldVector3D<T> initPosition()
initPosition
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
protected TimeStampedFieldPVCoordinates<T> initPVCoordinates()
initPVCoordinates
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
public FieldCartesianOrbit<T> shiftedBy(double dt)
dt
- time shift in secondspublic FieldCartesianOrbit<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.CalculusFieldElement<T>>,T extends org.hipparchus.CalculusFieldElement<T>>
shiftedBy
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
dt
- time shift in secondsprotected T[][] computeJacobianMeanWrtCartesian()
FieldOrbit
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 FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
FieldOrbit.computeJacobianEccentricWrtCartesian()
,
FieldOrbit.computeJacobianTrueWrtCartesian()
protected T[][] computeJacobianEccentricWrtCartesian()
FieldOrbit
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 FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
FieldOrbit.computeJacobianMeanWrtCartesian()
,
FieldOrbit.computeJacobianTrueWrtCartesian()
protected T[][] computeJacobianTrueWrtCartesian()
FieldOrbit
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 FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
FieldOrbit.computeJacobianMeanWrtCartesian()
,
FieldOrbit.computeJacobianEccentricWrtCartesian()
public void addKeplerContribution(PositionAngleType type, T gm, T[] 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 FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
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)public String toString()
public CartesianOrbit toOrbit()
FieldOrbit
toOrbit
in class FieldOrbit<T extends org.hipparchus.CalculusFieldElement<T>>
Copyright © 2002-2023 CS GROUP. All rights reserved.