|
||||||||||
| PREV CLASS NEXT CLASS | FRAMES NO FRAMES | |||||||||
| SUMMARY: NESTED | FIELD | CONSTR | METHOD | DETAIL: FIELD | CONSTR | METHOD | |||||||||
java.lang.Objectorg.orekit.orbits.Orbit
org.orekit.orbits.EquatorialOrbit
public final class EquatorialOrbit
This class handles non circular equatorial orbital parameters.
The parameters used internally are the equatorial elements (see EquatorialParameters
for more information.
The instance EquatorialOrbit is guaranteed to be immutable.
Orbit,
CircularOrbit,
CartesianOrbit,
EquinoctialOrbit,
Serialized Form| Constructor Summary | |
|---|---|
EquatorialOrbit(double a,
double e,
double pomega,
double ix,
double iy,
double anomaly,
PositionAngle type,
Frame frame,
AbsoluteDate date,
double mu)
Creates a new instance. |
|
EquatorialOrbit(IOrbitalParameters parameters,
Frame frame,
AbsoluteDate date)
Creates a new instance. |
|
EquatorialOrbit(Orbit op)
Constructor from any kind of orbital parameters. |
|
EquatorialOrbit(PVCoordinates pvCoordinates,
Frame frame,
AbsoluteDate date,
double mu)
Constructor from cartesian parameters. |
|
| Method Summary | |
|---|---|
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. |
double |
getA()
Get the semi-major axis. |
double |
getAnomaly(PositionAngle type)
Get the anomaly. |
double |
getE()
Get the eccentricity. |
double |
getEccentricAnomaly()
Get the eccentric anomaly. |
EquatorialParameters |
getEquatorialParameters()
Getter for underlying equatorial parameters. |
double |
getEquinoctialEx()
Get the first component of the eccentricity vector. |
double |
getEquinoctialEy()
Get the second component of the eccentricity vector. |
double |
getHx()
Get hx = ix / (2 * cos(i/2)), where ix is the first component of the inclination vector. |
double |
getHy()
Get hy = iy / (2 * cos(i/2)), where iy is the second component of the inclination vector. |
double |
getI()
Get the inclination angle. |
double |
getIx()
Get the first component of the inclination vector. |
double |
getIy()
Get the second component of the inclination vector. |
double |
getLE()
Get the eccentric latitude argument. |
double |
getLM()
Get the mean latitude argument. |
double |
getLv()
Get the true latitude argument. |
double |
getMeanAnomaly()
Get the mean anomaly. |
IOrbitalParameters |
getParameters()
Get underlying orbital parameters. |
double |
getPomega()
Get the longitude of the periapsis (ω + Ω). |
double |
getTrueAnomaly()
Get the true anomaly. |
OrbitType |
getType()
Get the orbit type. |
protected PVCoordinates |
initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters. |
EquatorialOrbit |
interpolate(AbsoluteDate date,
Collection<Orbit> sample)
Get an interpolated instance. |
protected void |
orbitAddKeplerContribution(PositionAngle type,
double gm,
double[] pDot)
Add the contribution of the Keplerian motion to parameters derivatives |
protected EquatorialOrbit |
orbitShiftedBy(double dt)
Get a time-shifted orbit. |
String |
toString()
Returns a string representation of this non circular equatorial orbital parameters object. |
| Methods inherited from class java.lang.Object |
|---|
clone, equals, finalize, getClass, hashCode, notify, notifyAll, wait, wait, wait |
| Constructor Detail |
|---|
public EquatorialOrbit(IOrbitalParameters parameters,
Frame frame,
AbsoluteDate date)
parameters - orbital parametersframe - the frame in which the parameters are defined
(must be a pseudo-inertial frame)date - date of the orbital parameters
public EquatorialOrbit(double a,
double e,
double pomega,
double ix,
double iy,
double anomaly,
PositionAngle type,
Frame frame,
AbsoluteDate date,
double mu)
throws IllegalArgumentException
a - semi-major axis (m)e - eccentricitypomega - ω + Ω (rad)ix - 2 sin(i/2) cos(Ω), first component of inclination vectoriy - 2 sin(i/2) sin(Ω), second component of inclination vectoranomaly - (M or E or v) = anomaly mean, eccentric or true anomaly (rad)type - type of anomalyframe - the frame in which the parameters are defineddate - date of the orbital parametersmu - central attraction coefficient (m3/s2)
IllegalArgumentException - if orbit is hyperbolic
IllegalArgumentException - if orbit mismatch with conic type
IllegalArgumentException - if inclination vector is not valid, meaning ix^2 + iy^2 > 4
public EquatorialOrbit(PVCoordinates pvCoordinates,
Frame frame,
AbsoluteDate date,
double mu)
throws IllegalArgumentException
pvCoordinates - the PVCoordinates of the satelliteframe - the frame in which are defined the PVCoordinatesdate - date of the orbital parametersmu - central attraction coefficient (m3/s2)
IllegalArgumentException - if orbit is hyperbolicpublic EquatorialOrbit(Orbit op)
op - orbital parameters to copy| Method Detail |
|---|
public IOrbitalParameters getParameters()
getParameters in class Orbitpublic EquatorialParameters getEquatorialParameters()
public OrbitType getType()
getType in class Orbitpublic double getA()
getA in class Orbitpublic double getE()
getE in class Orbitpublic double getI()
getI in class Orbitpublic double getAnomaly(PositionAngle type)
type - type of the angle
public double getPomega()
public double getTrueAnomaly()
public double getEccentricAnomaly()
public double getMeanAnomaly()
public double getEquinoctialEx()
getEquinoctialEx in class Orbitpublic double getEquinoctialEy()
getEquinoctialEy in class Orbitpublic double getIx()
public double getIy()
public double getHx()
getHx in class Orbitpublic double getHy()
getHy in class Orbitpublic double getLv()
getLv in class Orbitpublic double getLE()
getLE in class Orbitpublic double getLM()
getLM in class Orbitprotected PVCoordinates initPVCoordinates()
initPVCoordinates in class Orbitprotected EquatorialOrbit orbitShiftedBy(double 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.
orbitShiftedBy in class Orbitdt - time shift in seconds
public EquatorialOrbit interpolate(AbsoluteDate date,
Collection<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.
The interpolated instance is created by polynomial Hermite interpolation on Keplerian elements, without derivatives (which means the interpolation falls back to Lagrange interpolation only).
date - interpolation datesample - sample points on which interpolation should be done
protected 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 OrbitOrbit.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 OrbitOrbit.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 OrbitOrbit.computeJacobianMeanWrtCartesian(),
Orbit.computeJacobianEccentricWrtCartesian()
protected void orbitAddKeplerContribution(PositionAngle type,
double gm,
double[] pDot)
This method is used by numerical propagators to evaluate the part of Keplerrian motion to evolution of the orbital state.
orbitAddKeplerContribution in class Orbittype - 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()
toString in class Object
|
||||||||||
| PREV CLASS NEXT CLASS | FRAMES NO FRAMES | |||||||||
| SUMMARY: NESTED | FIELD | CONSTR | METHOD | DETAIL: FIELD | CONSTR | METHOD | |||||||||