poliastro.twobody.orbit.scalar
¶
Module Contents¶
Classes¶
Position and velocity of a body with respect to an attractor |
Attributes¶
- poliastro.twobody.orbit.scalar.ORBIT_FORMAT = {r_p:.0f} x {r_a:.0f} x {inc:.1f} ({frame}) orbit around {body} at epoch {epoch} ({scale})¶
- poliastro.twobody.orbit.scalar.ORBIT_NO_FRAME_FORMAT = {r_p:.0f} x {r_a:.0f} x {inc:.1f} orbit around {body} at epoch {epoch} ({scale})¶
- class poliastro.twobody.orbit.scalar.Orbit(state, epoch)¶
Bases:
poliastro.twobody.orbit.creation.OrbitCreationMixin
Position and velocity of a body with respect to an attractor at a given time (epoch).
Regardless of how the Orbit is created, the implicit reference system is an inertial one. For the specific case of the Solar System, this can be assumed to be the International Celestial Reference System or ICRS.
- property attractor(self)¶
Main attractor.
- property epoch(self)¶
Epoch of the orbit.
- property plane(self)¶
Fundamental plane of the frame.
- r(self)¶
Position vector.
- v(self)¶
Velocity vector.
- a(self)¶
Semimajor axis.
- p(self)¶
Semilatus rectum.
- r_p(self)¶
Radius of pericenter.
- r_a(self)¶
Radius of apocenter.
- ecc(self)¶
Eccentricity.
- inc(self)¶
Inclination.
- raan(self)¶
Right ascension of the ascending node.
- argp(self)¶
Argument of the perigee.
- property nu(self)¶
True anomaly.
- f(self)¶
Second modified equinoctial element.
- g(self)¶
Third modified equinoctial element.
- h(self)¶
Fourth modified equinoctial element.
- k(self)¶
Fifth modified equinoctial element.
- L(self)¶
True longitude.
- period(self)¶
Period of the orbit.
- n(self)¶
Mean motion.
- energy(self)¶
Specific energy.
- e_vec(self)¶
Eccentricity vector.
- h_vec(self)¶
Specific angular momentum vector.
- h_mag(self)¶
Specific angular momentum.
- arglat(self)¶
Argument of latitude.
- t_p(self)¶
Elapsed time since latest perifocal passage.
- get_frame(self)¶
Get equivalent reference frame of the orbit.
New in version 0.14.0.
- change_attractor(self, new_attractor, force=False)¶
Changes orbit attractor.
Only changes from attractor to parent or the other way around are allowed.
- Parameters
new_attractor (poliastro.bodies.Body) – Desired new attractor.
force (bool) – If True, changes attractor even if physically has no-sense.
- Returns
ss – Orbit with new attractor
- Return type
poliastro.twobody.orbit.Orbit
- change_plane(self, plane)¶
Changes fundamental plane.
- Parameters
plane (Planes) – Fundamental plane of the frame.
- represent_as(self, representation, differential_class=None)¶
Converts the orbit to a specific representation.
New in version 0.11.0.
- Parameters
representation (BaseRepresentation) – Representation object to use. It must be a class, not an instance.
differential_class (BaseDifferential, optional) – Class in which the differential should be represented, default to None.
Examples
>>> from poliastro.examples import iss >>> from astropy.coordinates import SphericalRepresentation >>> iss.represent_as(CartesianRepresentation) <CartesianRepresentation (x, y, z) in km (859.07256, -4137.20368, 5295.56871)> >>> iss.represent_as(CartesianRepresentation).xyz <Quantity [ 859.07256, -4137.20368, 5295.56871] km> >>> iss.represent_as(CartesianRepresentation, CartesianDifferential).differentials['s'] <CartesianDifferential (d_x, d_y, d_z) in km / s (7.37289205, 2.08223573, 0.43999979)> >>> iss.represent_as(CartesianRepresentation, CartesianDifferential).differentials['s'].d_xyz <Quantity [7.37289205, 2.08223573, 0.43999979] km / s> >>> iss.represent_as(SphericalRepresentation, CartesianDifferential) <SphericalRepresentation (lon, lat, distance) in (rad, rad, km) (4.91712525, 0.89732339, 6774.76995296) (has differentials w.r.t.: 's')>
- rv(self)¶
Position and velocity vectors.
- classical(self)¶
Classical orbital elements.
- pqw(self)¶
Perifocal frame (PQW) vectors.
- __str__(self)¶
Return str(self).
- __repr__(self)¶
Return repr(self).
- propagate(self, value, method=FarnocchiaPropagator())¶
Propagates an orbit a specified time.
If value is true anomaly, propagate orbit to this anomaly and return the result. Otherwise, if time is provided, propagate this Orbit some time and return the result.
- time_to_anomaly(self, value)¶
Returns time required to be in a specific true anomaly.
- propagate_to_anomaly(self, value)¶
Propagates an orbit to a specific true anomaly.
- to_ephem(self, strategy=TrueAnomalyBounds())¶
Samples Orbit to return an ephemerides.
New in version 0.17.0.
- sample(self, values=100, *, min_anomaly=None, max_anomaly=None)¶
Samples an orbit to some specified time values.
New in version 0.8.0.
- Parameters
values (int) – Number of interval points (default to 100).
min_anomaly (Quantity, optional) – Anomaly limits to sample the orbit. For elliptic orbits the default will be \(E \in \left[0, 2 \pi \right]\), and for hyperbolic orbits it will be \(\nu \in \left[-\nu_c, \nu_c \right]\), where \(\nu_c\) is either the current true anomaly or a value that corresponds to \(r = 3p\).
max_anomaly (Quantity, optional) – Anomaly limits to sample the orbit. For elliptic orbits the default will be \(E \in \left[0, 2 \pi \right]\), and for hyperbolic orbits it will be \(\nu \in \left[-\nu_c, \nu_c \right]\), where \(\nu_c\) is either the current true anomaly or a value that corresponds to \(r = 3p\).
- Returns
positions – Array of x, y, z positions.
- Return type
Notes
When specifying a number of points, the initial and final position is present twice inside the result (first and last row). This is more useful for plotting.
Examples
>>> from astropy import units as u >>> from poliastro.examples import iss >>> iss.sample() <CartesianRepresentation (x, y, z) in km ... >>> iss.sample(10) <CartesianRepresentation (x, y, z) in km ...
- apply_maneuver(self, maneuver, intermediate=False)¶
Returns resulting Orbit after applying maneuver to self.
Optionally return intermediate states (default to False).
- plot(self, label=None, use_3d=False, interactive=False)¶
Plots the orbit.
- Parameters
label (str, optional) – Label for the orbit, defaults to empty.
use_3d (bool, optional) – Produce a 3D plot, default to False.
interactive (bool, optional) – Produce an interactive (rather than static) image of the orbit, default to False. This option requires Plotly properly installed and configured for your environment.
- elevation(self, lat, theta, h)¶
Elevation
- Parameters
lat (astropy.units.Quantity) – Latitude of the observation point on the attractor.
theta (astropy.units.Quantity) – Local sideral time of the observation point on the attractor.
h (astropy.units.Quantity) – Height of the station above the attractor.
- Returns
elevation – Elevation of the orbit with respect to a location on attractor, in units of radian.
- Return type
Notes
Local sideral time needs to be precomputed. If Earth is the attractor, it can be computed using poliastro.earth.util.get_local_sidereal_time.
- classmethod from_vectors(cls, attractor, r, v, epoch=J2000, plane=Planes.EARTH_EQUATOR)¶
Return Orbit from position and velocity vectors.
- classmethod from_coords(cls, attractor, coord, plane=Planes.EARTH_EQUATOR)¶
Creates an Orbit from an attractor and astropy SkyCoord or BaseCoordinateFrame instance.
This method accepts position and velocity in any reference frame unlike Orbit.from_vector which can accept inputs in only inertial reference frame centred at attractor. Also note that the frame information is lost after creation of the orbit and only the inertial reference frame at body centre will be used for all purposes.
- Parameters
attractor (Body) – Main attractor
coord (SkyCoord or BaseCoordinateFrame) – Position and velocity vectors in any reference frame. Note that coord must have a representation and its differential with respect to time.
plane (Planes, optional) – Final orbit plane, default to Earth Equator.
- classmethod from_classical(cls, attractor, a, ecc, inc, raan, argp, nu, epoch=J2000, plane=Planes.EARTH_EQUATOR)¶
Return Orbit from classical orbital elements.
- Parameters
attractor (Body) – Main attractor.
a (Quantity) – Semi-major axis.
ecc (Quantity) – Eccentricity.
inc (Quantity) – Inclination
raan (Quantity) – Right ascension of the ascending node.
argp (Quantity) – Argument of the pericenter.
nu (Quantity) – True anomaly.
epoch (Time, optional) – Epoch, default to J2000.
plane (Planes) – Fundamental plane of the frame.
- classmethod from_equinoctial(cls, attractor, p, f, g, h, k, L, epoch=J2000, plane=Planes.EARTH_EQUATOR)¶
Return Orbit from modified equinoctial elements.
- Parameters
attractor (Body) – Main attractor.
p (Quantity) – Semilatus rectum.
f (Quantity) – Second modified equinoctial element.
g (Quantity) – Third modified equinoctial element.
h (Quantity) – Fourth modified equinoctial element.
k (Quantity) – Fifth modified equinoctial element.
L (Quantity) – True longitude.
epoch (Time, optional) – Epoch, default to J2000.
plane (Planes) – Fundamental plane of the frame.
- classmethod from_ephem(cls, attractor, ephem, epoch)¶
Create osculating orbit from ephemerides at a given epoch.
This will assume that the Ephem coordinates are expressed with respect the given body.
- classmethod from_sbdb(cls, name, **kwargs)¶
Return osculating Orbit by using SBDB from Astroquery.
- Parameters
name (str) – Name of the body to make the request.
**kwargs – Extra kwargs for astroquery.
- Returns
ss – Orbit corresponding to body_name
- Return type
poliastro.twobody.orbit.Orbit
Examples
>>> from poliastro.twobody.orbit import Orbit >>> apophis_orbit = Orbit.from_sbdb('apophis')
- classmethod circular(cls, attractor, alt, inc=0 * u.deg, raan=0 * u.deg, arglat=0 * u.deg, epoch=J2000, plane=Planes.EARTH_EQUATOR)¶
Return circular Orbit.
- Parameters
attractor (Body) – Main attractor.
alt (Quantity) – Altitude over surface.
inc (Quantity, optional) – Inclination, default to 0 deg (equatorial orbit).
raan (Quantity, optional) – Right ascension of the ascending node, default to 0 deg.
arglat (Quantity, optional) – Argument of latitude, default to 0 deg.
epoch (Time, optional) – Epoch, default to J2000.
plane (Planes) – Fundamental plane of the frame.
- classmethod stationary(cls, attractor)¶
Return the stationary orbit for the given attractor and its rotational speed.
- classmethod synchronous(cls, attractor, period_mul=1 * u.one, ecc=0 * u.one, inc=0 * u.deg, argp=0 * u.deg, arglat=0 * u.deg, raan=0 * u.deg, epoch=J2000, plane=Planes.EARTH_EQUATOR)¶
Returns an orbit where the orbital period equals the rotation rate of the orbited body. The synchronous altitude for any central body can directly be obtained from Kepler’s Third Law by setting the orbit period Psync, equal to the rotation period of the central body relative to the fixed stars D*. In order to obtain this, it’s important to match orbital period with sidereal rotation period.
- Parameters
attractor (Body) – Main attractor.
period_mul (Quantity) – Multiplier, default to 1 to indicate that the period of the body is equal to the sidereal rotational period of the body being orbited, 0.5 a period equal to half the average rotational period of the body being orbited, indicates a semi-synchronous orbit.
ecc (Quantity) – Eccentricity,default to 0 as a stationary orbit.
inc (Quantity) – Inclination,default to 0 deg.
raan (Quantity) – Right ascension of the ascending node,default to 0 deg.
argp (Quantity) – Argument of the pericenter,default to 0 deg.
arglat (Quantity, optional) – Argument of latitude, default to 0 deg.
epoch (Time, optional) – Epoch, default to J2000.
plane (Planes) – Fundamental plane of the frame.
- Returns
New orbit.
- Return type
- Raises
ValueError – If the pericenter is smaller than the attractor’s radius.
Notes
Thus:
\[ \begin{align}\begin{aligned}\begin{split}P_{s y n c}=D^{*} \\\end{split}\\\begin{split}a_{s y n c}=\left(\mu / 4 \pi^{2}\right)^{1 / 3}\left(D^{*}\right)^{2 / 3}\\\end{split}\\\begin{split}H_{s y n c}=a_{s y n c} - R_{p l a n e t}\\\end{split}\end{aligned}\end{align} \]
- classmethod heliosynchronous(cls, attractor, a=None, ecc=None, inc=None, raan=0 * u.deg, argp=0 * u.deg, nu=0 * u.deg, epoch=J2000, plane=Planes.EARTH_EQUATOR)¶
Creates a Sun-Synchronous orbit.
These orbits make use of the J2 perturbation to precess in order to be always towards Sun. At least two parameters of the set {a, ecc, inc} are needed in order to solve for these kind of orbits.
Relationships among them are given by:
\[\begin{split}\begin{align} a &= \left (\frac{-3R_{\bigoplus}J_{2}\sqrt{\mu}\cos(i)}{2\dot{\Omega}(1-e^2)^2} \right ) ^ {\frac{2}{7}}\\ e &= \sqrt{1 - \sqrt{\frac{-3R_{\bigoplus}J_{2}\sqrt{\mu}cos(i)}{2a^{\frac{7}{2}}\dot{\Omega}}}}\\ i &= \arccos{\left ( \frac{-2a^{\frac{7}{2}}\dot{\Omega}(1-e^2)^2}{3R_{\bigoplus}J_{2}\sqrt{\mu}} \right )}\\ \end{align}\end{split}\]- Parameters
attractor (SolarSystemPlanet) – Attractor.
a (Quantity) – Semi-major axis.
ecc (Quantity) – Eccentricity.
inc (Quantity) – Inclination.
raan (Quantity) – Right ascension of the ascending node.
argp (Quantity) – Argument of the pericenter.
nu (Quantity) – True anomaly.
epoch (Time, optional) – Epoch, default to J2000.
plane (Planes) – Fundamental plane of the frame.
- classmethod parabolic(cls, attractor, p, inc, raan, argp, nu, epoch=J2000, plane=Planes.EARTH_EQUATOR)¶
Return a parabolic Orbit.
- Parameters
attractor (Body) – Main attractor.
p (Quantity) – Semilatus rectum or parameter.
inc (Quantity, optional) – Inclination.
raan (Quantity) – Right ascension of the ascending node.
argp (Quantity) – Argument of the pericenter.
nu (Quantity) – True anomaly.
epoch (Time, optional) – Epoch, default to J2000.
plane (Planes) – Fundamental plane of the frame.
- classmethod frozen(cls, attractor, alt, inc=None, argp=None, raan=0 * u.deg, arglat=0 * u.deg, ecc=None, epoch=J2000, plane=Planes.EARTH_EQUATOR)¶
Return a frozen Orbit. If any of the given arguments results in an impossibility, some values will be overwritten
To achieve frozen orbit these two equations have to be set to zero.
\[\dfrac {d\overline {e}}{dt}=\dfrac {-3\overline {n}J_{3}R^{3}_{E}\sin \left( \overline {i}\right) }{2a^{3}\left( 1-\overline {e}^{2}\right) ^{2}}\left( 1-\dfrac {5}{4}\sin ^{2}\overline {i}\right) \cos \overline {w}\]\[\dfrac {d\overline {\omega }}{dt}=\dfrac {3\overline {n}J_{2}R^{2}_{E}}{a^{2}\left( 1-\overline {e}^{2}\right) ^{2}}\left( 1-\dfrac {5}{4}\sin ^{2}\overline {i}\right) \left[ 1+\dfrac {J_{3}R_{E}}{2J_{2}\overline {a}\left( 1-\overline {e}^{2}\right) }\left( \dfrac {\sin ^{2}\overline {i}-\overline {e}\cos ^{2}\overline {i}}{\sin \overline {i}}\right) \dfrac {\sin \overline {w}}{\overline {e}}\right]\]The first approach would be to nullify following term to zero:
\[( 1-\dfrac {5}{4}\sin ^{2})\]For which one obtains the so-called critical inclinations: i = 63.4349 or 116.5651 degrees. To escape the inclination requirement, the argument of periapsis can be set to w = 90 or 270 degrees to nullify the second equation. Then, one should nullify the right-hand side of the first equation, which yields an expression that correlates the inclination of the object and the eccentricity of the orbit:
\[\overline {e}=-\dfrac {J_{3}R_{E}}{2J_{2}\overline {a}\left( 1-\overline {e}^{2}\right) }\left( \dfrac {\sin ^{2}\overline {i}-\overline {e}\cos ^{2} \overline {i}}{\sin \overline {i}}\right)\]Assuming that e is negligible compared to J2, it can be shown that:
\[\overline {e}\approx -\dfrac {J_{3}R_{E}}{2J_{2}\overline {a}}\sin \overline {i}\]The implementation is divided in the following cases:
When the user gives a negative altitude, the method will raise a ValueError
When the attractor has not defined J2 or J3, the method will raise an AttributeError
When the attractor has J2/J3 outside of range 1 to 10 , the method will raise an NotImplementedError. Special case for Venus.See “Extension of the critical inclination” by Xiaodong Liu, Hexi Baoyin, and Xingrui Ma
If argp is not given or the given argp is a critical value:
if eccentricity is none and inclination is none, the inclination is set with a critical value and the eccentricity is obtained from the last formula mentioned
if only eccentricity is none, we calculate this value with the last formula mentioned
if only inclination is none ,we calculate this value with the formula for eccentricity with critical argp.
If inc is not given or the given inc is critical:
if the argp and the eccentricity is given we keep these values to create the orbit
if the eccentricity is given we keep this value, if not, default to the eccentricity of the Moon’s orbit around the Earth
if it’s not possible to create an orbit with the the argp and the inclination given, both of them are set to the critical values and the eccentricity is calculate with the last formula
- Parameters
attractor (Body) – Main attractor.
alt (Quantity) – Altitude over surface.
inc (Quantity, optional) – Inclination, default to critical value.
argp (Quantity, optional) – Argument of the pericenter, default to critical value.
raan (Quantity, optional) – Right ascension of the ascending node, default to 0 deg.
arglat (Quantity, optional) – Argument of latitude, default to 0 deg.
ecc (Quantity) – Eccentricity, default to the eccentricity of the Moon’s orbit around the Earth
epoch (Time, optional) – Epoch, default to J2000.
plane (Planes) – Fundamental plane of the frame.