public class CircularOrbit extends Orbit
The parameters used internally are the circular elements which can be related to keplerian elements as follows:
The conversion equations from and to keplerian elements given above hold only
when both sides are unambiguously defined, i.e. when orbit is neither equatorial
nor circular. When orbit is circular (but not equatorial), the circular
parameters are still unambiguously defined whereas some keplerian elements
(more precisely ω and Ω) become ambiguous. When orbit is equatorial,
neither the keplerian nor the circular parameters can be defined unambiguously.
equinoctial orbits
is the recommended way to represent
orbits.
The instance CircularOrbit
is guaranteed to be immutable.
Orbit
,
KeplerianOrbit
,
CartesianOrbit
,
EquinoctialOrbit
,
Serialized FormConstructor and Description |
---|
CircularOrbit(double a,
double ex,
double ey,
double i,
double raan,
double alpha,
PositionAngle type,
Frame frame,
AbsoluteDate date,
double mu)
Creates a new instance.
|
CircularOrbit(double a,
double ex,
double ey,
double i,
double raan,
double alpha,
PositionAngle type,
TimeStampedPVCoordinates pvCoordinates,
Frame frame,
double mu)
Creates a new instance.
|
CircularOrbit(Orbit op)
Constructor from any kind of orbital parameters.
|
CircularOrbit(PVCoordinates pvCoordinates,
Frame frame,
AbsoluteDate date,
double mu)
Constructor from cartesian parameters.
|
CircularOrbit(TimeStampedPVCoordinates pvCoordinates,
Frame frame,
double mu)
Constructor from cartesian parameters.
|
Modifier and Type | Method and Description |
---|---|
void |
addKeplerContribution(PositionAngle type,
double gm,
double[] pDot)
Add the contribution of the Keplerian motion to parameters derivatives
|
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 |
getAlpha(PositionAngle type)
Get the latitude argument.
|
double |
getAlphaE()
Get the eccentric latitude argument.
|
double |
getAlphaM()
Get the mean latitude argument.
|
double |
getAlphaV()
Get the true latitude argument.
|
double |
getCircularEx()
Get the first component of the circular eccentricity vector.
|
double |
getCircularEy()
Get the second component of the circular eccentricity vector.
|
double |
getE()
Get the eccentricity.
|
double |
getEquinoctialEx()
Get the first component of the equinoctial eccentricity vector.
|
double |
getEquinoctialEy()
Get the second component of the equinoctial eccentricity vector.
|
double |
getHx()
Get the first component of the inclination vector.
|
double |
getHy()
Get the second component of the inclination vector.
|
double |
getI()
Get the inclination.
|
double |
getLE()
Get the eccentric longitude argument.
|
double |
getLM()
Get the mean longitude argument.
|
double |
getLv()
Get the true longitude argument.
|
double |
getRightAscensionOfAscendingNode()
Get the right ascension of the ascending node.
|
OrbitType |
getType()
Get the orbit type.
|
protected TimeStampedPVCoordinates |
initPVCoordinates()
Compute the position/velocity coordinates from the canonical parameters.
|
CircularOrbit |
interpolate(AbsoluteDate date,
Collection<Orbit> sample)
Get an interpolated instance.
|
CircularOrbit |
shiftedBy(double dt)
Get a time-shifted orbit.
|
String |
toString()
Returns a string representation of this Orbit object.
|
fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, fillHalfRow, getDate, getFrame, getJacobianWrtCartesian, getJacobianWrtParameters, getKeplerianMeanMotion, getKeplerianPeriod, getMu, getPVCoordinates, getPVCoordinates, getPVCoordinates
public CircularOrbit(double a, double ex, double ey, double i, double raan, double alpha, PositionAngle type, Frame frame, AbsoluteDate date, double mu) throws IllegalArgumentException
a
- semi-major axis (m)ex
- e cos(ω), first component of circular eccentricity vectorey
- e sin(ω), second component of circular eccentricity vectori
- inclination (rad)raan
- right ascension of ascending node (Ω, rad)alpha
- an + ω, mean, eccentric or true latitude argument (rad)type
- type of latitude argumentframe
- the frame in which are defined the parameters
(must be a pseudo-inertial frame
)date
- date of the orbital parametersmu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if eccentricity is equal to 1 or larger or
if frame is not a pseudo-inertial frame
public CircularOrbit(double a, double ex, double ey, double i, double raan, double alpha, PositionAngle type, TimeStampedPVCoordinates pvCoordinates, Frame frame, double mu) throws IllegalArgumentException
a
- semi-major axis (m)ex
- e cos(ω), first component of circular eccentricity vectorey
- e sin(ω), second component of circular eccentricity vectori
- inclination (rad)raan
- right ascension of ascending node (Ω, rad)alpha
- an + ω, mean, eccentric or true latitude argument (rad)type
- type of latitude argumentpvCoordinates
- the PVCoordinates
in inertial frameframe
- the frame in which are defined the parameters
(must be a pseudo-inertial frame
)mu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if eccentricity is equal to 1 or larger or
if frame is not a pseudo-inertial frame
public CircularOrbit(TimeStampedPVCoordinates pvCoordinates, Frame frame, double mu) throws IllegalArgumentException
The acceleration provided in pvCoordinates
is accessible using
Orbit.getPVCoordinates()
and Orbit.getPVCoordinates(Frame)
. All other methods
use mu
and the position to compute the acceleration, including
shiftedBy(double)
and Orbit.getPVCoordinates(AbsoluteDate, Frame)
.
pvCoordinates
- the PVCoordinates
in inertial frameframe
- the frame in which are defined the PVCoordinates
(must be a pseudo-inertial frame
)mu
- central attraction coefficient (m³/s²)IllegalArgumentException
- if frame is not a pseudo-inertial frame
public CircularOrbit(PVCoordinates pvCoordinates, Frame frame, AbsoluteDate date, double mu) throws IllegalArgumentException
The acceleration provided in pvCoordinates
is accessible using
Orbit.getPVCoordinates()
and Orbit.getPVCoordinates(Frame)
. All other methods
use mu
and the position to compute the acceleration, including
shiftedBy(double)
and Orbit.getPVCoordinates(AbsoluteDate, Frame)
.
pvCoordinates
- the PVCoordinates
in inertial frameframe
- the frame in which are defined the PVCoordinates
(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 CircularOrbit(Orbit op)
op
- orbital parameters to copypublic OrbitType getType()
public double getA()
Note that the semi-major axis is considered negative for hyperbolic orbits.
public double getEquinoctialEx()
getEquinoctialEx
in class Orbit
public double getEquinoctialEy()
getEquinoctialEy
in class Orbit
public double getCircularEx()
public double getCircularEy()
public double getHx()
public double getHy()
public double getAlphaV()
public double getAlpha(PositionAngle type)
type
- type of the anglepublic double getAlphaE()
public double getAlphaM()
public double getE()
public double getI()
public double getRightAscensionOfAscendingNode()
public double getLv()
public double getLE()
public double getLM()
protected TimeStampedPVCoordinates initPVCoordinates()
initPVCoordinates
in class Orbit
public CircularOrbit shiftedBy(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.
shiftedBy
in interface TimeShiftable<Orbit>
shiftedBy
in class Orbit
dt
- time shift in secondspublic CircularOrbit 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.
Note that this method is designed for small samples only (say up to about 10-20 points) so it can be implemented using polynomial interpolation (typically Hermite interpolation). Using too much points may induce Runge's phenomenon and numerical problems (including NaN appearing).
The interpolated instance is created by polynomial Hermite interpolation on circular elements, without derivatives (which means the interpolation falls back to Lagrange interpolation only).
As this implementation of interpolation is polynomial, it should be used only with small samples (about 10-20 points) in order to avoid Runge's phenomenon and numerical problems (including NaN appearing).
If orbit interpolation on large samples is needed, using the Ephemeris
class is a better way than using this
low-level interpolation. The Ephemeris class automatically handles selection of
a neighboring sub-sample with a predefined number of point from a large global sample
in a thread-safe way.
date
- interpolation datesample
- sample points on which interpolation should be doneprotected 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 Orbit
Orbit.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 Orbit
Orbit.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 Orbit
Orbit.computeJacobianMeanWrtCartesian()
,
Orbit.computeJacobianEccentricWrtCartesian()
public 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.
addKeplerContribution
in class Orbit
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)Copyright © 2002-2015 CS Systèmes d'information. All rights reserved.