High level API¶
poliastro.twobody package¶
poliastro.twobody.angles module¶
Angles and anomalies.
-
poliastro.twobody.angles.D_to_nu(D)¶ True anomaly from parabolic eccentric anomaly.
Parameters: D (Quantity) – Eccentric anomaly. Returns: nu – True anomaly. Return type: Quantity Notes
Taken from Farnocchia, Davide, Davide Bracali Cioci, and Andrea Milani. “Robust resolution of Kepler’s equation in all eccentricity regimes.” Celestial Mechanics and Dynamical Astronomy 116, no. 1 (2013): 21-34.
-
poliastro.twobody.angles.nu_to_D(nu)¶ Parabolic eccentric anomaly from true anomaly.
Parameters: nu (Quantity) – True anomaly. Returns: D – Hyperbolic eccentric anomaly. Return type: Quantity Notes
Taken from Farnocchia, Davide, Davide Bracali Cioci, and Andrea Milani. “Robust resolution of Kepler’s equation in all eccentricity regimes.” Celestial Mechanics and Dynamical Astronomy 116, no. 1 (2013): 21-34.
-
poliastro.twobody.angles.nu_to_E(nu, ecc)¶ Eccentric anomaly from true anomaly.
New in version 0.4.0.
Parameters: Returns: E – Eccentric anomaly.
Return type:
-
poliastro.twobody.angles.nu_to_F(nu, ecc)¶ Hyperbolic eccentric anomaly from true anomaly.
Parameters: Returns: F – Hyperbolic eccentric anomaly.
Return type: Note
Taken from Curtis, H. (2013). Orbital mechanics for engineering students. 167
-
poliastro.twobody.angles.E_to_nu(E, ecc)¶ True anomaly from eccentric anomaly.
New in version 0.4.0.
Parameters: Returns: nu – True anomaly.
Return type:
-
poliastro.twobody.angles.F_to_nu(F, ecc)¶ True anomaly from hyperbolic eccentric anomaly.
Parameters: Returns: nu – True anomaly.
Return type:
-
poliastro.twobody.angles.M_to_E(M, ecc)¶ Eccentric anomaly from mean anomaly.
New in version 0.4.0.
Parameters: Returns: E – Eccentric anomaly.
Return type:
-
poliastro.twobody.angles.M_to_F(M, ecc)¶ Hyperbolic eccentric anomaly from mean anomaly.
Parameters: Returns: F – Hyperbolic eccentric anomaly.
Return type:
-
poliastro.twobody.angles.M_to_D(M, ecc)¶ Parabolic eccentric anomaly from mean anomaly.
Parameters: Returns: D – Parabolic eccentric anomaly.
Return type:
-
poliastro.twobody.angles.E_to_M(E, ecc)¶ Mean anomaly from eccentric anomaly.
New in version 0.4.0.
Parameters: Returns: M – Mean anomaly.
Return type:
-
poliastro.twobody.angles.F_to_M(F, ecc)¶ Mean anomaly from eccentric anomaly.
Parameters: Returns: M – Mean anomaly.
Return type:
-
poliastro.twobody.angles.D_to_M(D, ecc)¶ Mean anomaly from eccentric anomaly.
Parameters: Returns: M – Mean anomaly.
Return type:
-
poliastro.twobody.angles.M_to_nu(M, ecc, delta=0.01)¶ True anomaly from mean anomaly.
New in version 0.4.0.
Parameters: Returns: nu – True anomaly.
Return type: Examples
>>> M_to_nu(30.0 * u.deg, 0.06 * u.one) <Quantity 33.67328493 deg>
-
poliastro.twobody.angles.nu_to_M(nu, ecc, delta=0.01)¶ Mean anomaly from true anomaly.
New in version 0.4.0.
Parameters: Returns: M – Mean anomaly.
Return type:
poliastro.twobody.classical module¶
Functions to define orbits from classical orbital elements.
-
class
poliastro.twobody.classical.ClassicalState(attractor, p, ecc, inc, raan, argp, nu)¶ State defined by its classical orbital elements.
-
p¶ Semilatus rectum.
-
ecc¶ Eccentricity.
-
inc¶ Inclination.
-
raan¶ Right ascension of the ascending node.
-
argp¶ Argument of the perigee.
-
nu¶ True anomaly.
-
to_vectors()¶ Converts to position and velocity vector representation.
-
to_classical()¶ Converts to classical orbital elements representation.
-
to_equinoctial()¶ Converts to modified equinoctial elements representation.
-
poliastro.twobody.decorators module¶
Decorators.
-
poliastro.twobody.decorators.state_from_vector(func)¶ Changes signature to receive Orbit instead of state array.
Examples
>>> from poliastro.twobody.decorators import state_from_vector >>> @state_from_vector ... def func(_, ss): ... return ss.r, ss.v ... >>> func(0.0, [1, 2, 3, -1, -2, -3], 1.0) (<Quantity [1., 2., 3.] km>, <Quantity [-1., -2., -3.] km / s>)
Notes
Functions decorated with this will have poor performance.
poliastro.twobody.equinoctial module¶
Functions to define orbits from modified equinoctial orbital elements.
-
poliastro.twobody.equinoctial.mee2coe(p, f, g, h, k, L)¶ Converts from modified equinoctial orbital elements to classical orbital elements.
The definition of the modified equinoctial orbital elements is taken from [Walker, 1985].
Note
The conversion is always safe because arctan2 works also for 0, 0 arguments.
-
class
poliastro.twobody.equinoctial.ModifiedEquinoctialState(attractor, p, f, g, h, k, L)¶ -
p¶ Semimajor axis.
-
f¶ Second modified equinoctial element.
-
g¶ Third modified equinoctial element.
-
h¶ Fourth modified equinoctial element.
-
k¶ Fifth modified equinoctial element.
-
L¶ True longitude.
-
to_classical()¶ Converts to classical orbital elements representation.
-
poliastro.twobody.orbit module¶
-
exception
poliastro.twobody.orbit.TimeScaleWarning¶
-
class
poliastro.twobody.orbit.Orbit(state, epoch, plane)¶ 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.
-
state¶ Position and velocity or orbital elements.
-
epoch¶ Epoch of the orbit.
-
plane¶ Fundamental plane of the frame.
-
frame¶ Reference frame of the orbit.
New in version 0.11.0.
-
classmethod
from_vectors(attractor, r, v, epoch=<Time object: scale='tt' format='jyear_str' value=J2000.000>, plane=<Planes.EARTH_EQUATOR: 'Earth mean Equator and Equinox of epoch (J2000.0)'>)¶ Return Orbit from position and velocity vectors.
Parameters:
-
classmethod
from_classical(attractor, a, ecc, inc, raan, argp, nu, epoch=<Time object: scale='tt' format='jyear_str' value=J2000.000>, plane=<Planes.EARTH_EQUATOR: 'Earth mean Equator and Equinox of epoch (J2000.0)'>)¶ 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(attractor, p, f, g, h, k, L, epoch=<Time object: scale='tt' format='jyear_str' value=J2000.000>, plane=<Planes.EARTH_EQUATOR: 'Earth mean Equator and Equinox of epoch (J2000.0)'>)¶ 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_body_ephem(body, epoch=None)¶ Return osculating Orbit of a body at a given time.
-
classmethod
circular(attractor, alt, inc=<Quantity 0. deg>, raan=<Quantity 0. deg>, arglat=<Quantity 0. deg>, epoch=<Time object: scale='tt' format='jyear_str' value=J2000.000>, plane=<Planes.EARTH_EQUATOR: 'Earth mean Equator and Equinox of epoch (J2000.0)'>)¶ 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
parabolic(attractor, p, inc, raan, argp, nu, epoch=<Time object: scale='tt' format='jyear_str' value=J2000.000>, plane=<Planes.EARTH_EQUATOR: 'Earth mean Equator and Equinox of epoch (J2000.0)'>)¶ Return 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.
-
represent_as(representation)¶ 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. Examples
>>> from poliastro.examples import iss >>> from astropy.coordinates import CartesianRepresentation, SphericalRepresentation >>> iss.represent_as(CartesianRepresentation) <CartesianRepresentation (x, y, z) in km (859.07256, -4137.20368, 5295.56871) (has differentials w.r.t.: 's')> >>> iss.represent_as(CartesianRepresentation).xyz <Quantity [ 859.07256, -4137.20368, 5295.56871] km> >>> iss.represent_as(CartesianRepresentation).differentials['s'] <CartesianDifferential (d_x, d_y, d_z) in km / s (7.37289205, 2.08223573, 0.43999979)> >>> iss.represent_as(CartesianRepresentation).differentials['s'].d_xyz <Quantity [7.37289205, 2.08223573, 0.43999979] km / s> >>> iss.represent_as(SphericalRepresentation) <SphericalRepresentation (lon, lat, distance) in (rad, rad, km) (4.91712525, 0.89732339, 6774.76995296) (has differentials w.r.t.: 's')>
-
to_icrs()¶ Creates a new Orbit object with its coordinates transformed to ICRS.
Notice that, strictly speaking, the center of ICRS is the Solar System Barycenter and not the Sun, and therefore these orbits cannot be propagated in the context of the two body problem. Therefore, this function exists merely for practical purposes.
New in version 0.11.0.
-
propagate(value, method=<function mean_motion>, rtol=1e-10, **kwargs)¶ Propagates an orbit.
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.
Parameters: - value (Multiple options) – True anomaly values or time values. If given an angle, it will always propagate forward.
- rtol (float, optional) – Relative tolerance for the propagation algorithm, default to 1e-10.
- method (function, optional) – Method used for propagation
- **kwargs – parameters used in perturbation models
-
sample(values=None, method=<function mean_motion>)¶ Samples an orbit to some specified time values.
New in version 0.8.0.
Parameters: - values (Multiple options) – Number of interval points (default to 100), True anomaly values, Time values.
- method (function, optional) – Method used for propagation
Returns: positions – Array of x, y, z positions, with proper times as the frame attributes if supported.
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() # doctest: +ELLIPSIS <GCRS Coordinate ...> >>> iss.sample(10) # doctest: +ELLIPSIS <GCRS Coordinate ...> >>> iss.sample([0, 180] * u.deg) # doctest: +ELLIPSIS <GCRS Coordinate ...> >>> iss.sample([0, 10, 20] * u.minute) # doctest: +ELLIPSIS <GCRS Coordinate ...> >>> iss.sample([iss.epoch + iss.period / 2]) # doctest: +ELLIPSIS <GCRS Coordinate ...>
-
poliastro.twobody.propagation module¶
Propagation algorithms
-
poliastro.twobody.propagation.func_twobody(t0, u_, k, ad, ad_kwargs)¶ Differential equation for the initial value two body problem.
This function follows Cowell’s formulation.
Parameters:
-
poliastro.twobody.propagation.cowell(orbit, tof, rtol=1e-11, *, ad=None, **ad_kwargs)¶ Propagates orbit using Cowell’s formulation.
Parameters: - orbit (Orbit) – the Orbit object to propagate.
- ad (function(t0, u, k), optional) – Non Keplerian acceleration (km/s2), default to None.
- tof (Multiple options) – Time to propagate, float (s), Times to propagate, array of float (s).
- rtol (float, optional) – Maximum relative error permitted, default to 1e-10.
Raises: RuntimeError– If the algorithm didn’t converge.Note
This method uses a Dormand & Prince method of order 8(5,3) available in the
poliastro.integratorsmodule. If multiple tofs are provided, the method propagates to the maximum value and calculates the other values via dense output
-
poliastro.twobody.propagation.propagate(orbit, time_of_flight, *, method=<function mean_motion>, rtol=1e-10, **kwargs)¶ Propagate an orbit some time and return the result.
poliastro.twobody.rv module¶
Functions to define orbits from position and velocity vectors.
poliastro.twobody.thrust package¶
-
poliastro.twobody.thrust.change_a_inc.change_a_inc(k, a_0, a_f, inc_0, inc_f, f)¶ - Guidance law from the Edelbaum/Kéchichian theory, optimal transfer between circular inclined orbits
- (a_0, i_0) –> (a_f, i_f), ecc = 0.
Parameters: Notes
Edelbaum theory, reformulated by Kéchichian.
References
- Edelbaum, T. N. “Propulsion Requirements for Controllable Satellites”, 1961.
- Kéchichian, J. A. “Reformulation of Edelbaum’s Low-Thrust Transfer Problem Using Optimal Control Theory”, 1997.
Argument of perigee change, with formulas developed by Pollard.
References
- Pollard, J. E. “Simplified Approach for Assessment of Low-Thrust Elliptical Orbit Transfers”, 1997.
- Pollard, J. E. “Evaluation of Low-Thrust Orbital Maneuvers”, 1998.
-
poliastro.twobody.thrust.change_argp.change_argp(k, a, ecc, argp_0, argp_f, f)¶ Guidance law from the model. Thrust is aligned with an inertially fixed direction perpendicular to the semimajor axis of the orbit.
Parameters: f (float) – Magnitude of constant acceleration
Quasi optimal eccentricity-only change, with formulas developed by Pollard.
References
- Pollard, J. E. “Simplified Approach for Assessment of Low-Thrust Elliptical Orbit Transfers”, 1997.
-
poliastro.twobody.thrust.change_ecc_quasioptimal.change_ecc_quasioptimal(ss_0, ecc_f, f)¶ Guidance law from the model. Thrust is aligned with an inertially fixed direction perpendicular to the semimajor axis of the orbit.
Parameters:
Simultaneous eccentricity and inclination changes.
References
- Pollard, J. E. “Simplified Analysis of Low-Thrust Orbital Maneuvers”, 2000.
-
poliastro.twobody.thrust.change_inc_ecc.change_inc_ecc(ss_0, ecc_f, inc_f, f)¶ Guidance law from the model. Thrust is aligned with an inertially fixed direction perpendicular to the semimajor axis of the orbit.
Parameters:
poliastro.iod package¶
poliastro.iod.izzo module¶
Izzo’s algorithm for Lambert’s problem
-
poliastro.iod.izzo.lambert(k, r0, r, tof, M=0, numiter=35, rtol=1e-08)¶ Solves the Lambert problem using the Izzo algorithm.
New in version 0.5.0.
Parameters: - k (Quantity) – Gravitational constant of main attractor (km^3 / s^2).
- r0 (Quantity) – Initial position (km).
- r (Quantity) – Final position (km).
- tof (Quantity) – Time of flight (s).
- M (int, optional) – Number of full revolutions, default to 0.
- numiter (int, optional) – Maximum number of iterations, default to 35.
- rtol (float, optional) – Relative tolerance of the algorithm, default to 1e-8.
Yields: v0, v (tuple) – Pair of velocity solutions.
poliastro.iod.vallado module¶
Initial orbit determination.
-
poliastro.iod.vallado.lambert(k, r0, r, tof, short=True, numiter=35, rtol=1e-08)¶ Solves the Lambert problem.
New in version 0.3.0.
Parameters: - k (Quantity) – Gravitational constant of main attractor (km^3 / s^2).
- r0 (Quantity) – Initial position (km).
- r (Quantity) – Final position (km).
- tof (Quantity) – Time of flight (s).
- short (boolean, optional) – Find out the short path, default to True. If False, find long path.
- numiter (int, optional) – Maximum number of iterations, default to 35.
- rtol (float, optional) – Relative tolerance of the algorithm, default to 1e-8.
Raises: RuntimeError– If it was not possible to compute the orbit.Note
This uses the universal variable approach found in Battin, Mueller & White with the bisection iteration suggested by Vallado. Multiple revolutions not supported.
poliastro.neos package¶
Code related to NEOs.
Functions related to NEOs and different NASA APIs. All of them are coded as part of SOCIS 2017 proposal.
Notes
The orbits returned by the functions in this package are in the
HeliocentricEclipticJ2000 frame.
poliastro.neos.dastcom5 module¶
NEOs orbit from DASTCOM5 database.
-
poliastro.neos.dastcom5.asteroid_db()¶ Return complete DASTCOM5 asteroid database.
Returns: database – Database with custom dtype. Return type: numpy.ndarray
-
poliastro.neos.dastcom5.comet_db()¶ Return complete DASTCOM5 comet database.
Returns: database – Database with custom dtype. Return type: numpy.ndarray
-
poliastro.neos.dastcom5.orbit_from_name(name)¶ Return
Orbitgiven a name.Retrieve info from JPL DASTCOM5 database.
Parameters: name (str) – NEO name. Returns: orbit – NEO orbits. Return type: list (Orbit)
-
poliastro.neos.dastcom5.orbit_from_record(record)¶ Return
Orbitgiven a record.Retrieve info from JPL DASTCOM5 database.
Parameters: record (int) – Object record. Returns: orbit – NEO orbit. Return type: Orbit
-
poliastro.neos.dastcom5.record_from_name(name)¶ Search dastcom.idx and return logical records that match a given string.
Body name, SPK-ID, or alternative designations can be used.
Parameters: name (str) – Body name. Returns: records – DASTCOM5 database logical records matching str. Return type: list (int)
-
poliastro.neos.dastcom5.string_record_from_name(name)¶ Search dastcom.idx and return body full record.
Search DASTCOM5 index and return body records that match string, containing logical record, name, alternative designations, SPK-ID, etc.
Parameters: name (str) – Body name. Returns: lines – Body records Return type: list(str)
-
poliastro.neos.dastcom5.read_headers()¶ Read DASTCOM5 headers and return asteroid and comet headers.
Headers are two numpy arrays with custom dtype.
Returns: ast_header, com_header – DASTCOM5 headers. Return type: tuple (numpy.ndarray)
-
poliastro.neos.dastcom5.read_record(record)¶ Read DASTCOM5 record and return body data.
Body data consists of numpy array with custom dtype.
Parameters: record (int) – Body record. Returns: body_data – Body information. Return type: numpy.ndarray
-
poliastro.neos.dastcom5.download_dastcom5()¶ Downloads DASTCOM5 database.
Downloads and unzip DASTCOM5 file in default poliastro path (~/.poliastro).
-
poliastro.neos.dastcom5.entire_db()¶ Return complete DASTCOM5 database.
Merge asteroid and comet databases, only with fields related to orbital data, discarding the rest.
Returns: database – Database with custom dtype. Return type: numpy.ndarray
dastcom5 parameters¶
poliastro.neos.neows module¶
NEOs orbit from NEOWS and JPL SBDB
-
poliastro.neos.neows.orbit_from_spk_id(spk_id, api_key=None)¶ Return
Orbitgiven a SPK-ID.Retrieve info from NASA NeoWS API, and therefore it only works with NEAs (Near Earth Asteroids).
Parameters: Returns: orbit – NEA orbit.
Return type:
-
poliastro.neos.neows.spk_id_from_name(name)¶ Return SPK-ID number given a small-body name.
Retrieve and parse HTML from JPL Small Body Database to get SPK-ID.
Parameters: name (str) – Small-body object name. Wildcards “*” and/or “?” can be used. Returns: spk_id – SPK-ID number. Return type: str
poliastro.bodies module¶
Bodies of the Solar System.
Contains some predefined bodies of the Solar System:
- Sun (☉)
- Earth (♁)
- Moon (☾)
- Mercury (☿)
- Venus (♀)
- Mars (♂)
- Jupiter (♃)
- Saturn (♄)
- Uranus (⛢)
- Neptune (♆)
- Pluto (♇)
and a way to define new bodies (Body class).
Data references can be found in constants
-
class
poliastro.bodies.Body(parent, k, name, symbol=None, R=<Quantity 0. km>, **kwargs)¶ Class to represent a generic body.
poliastro.constants module¶
Astronomical and physics constants.
This module complements constants defined in astropy.constants, with gravitational paremeters and radii.
Note that GM_jupiter and GM_neptune are both referred to the whole planetary system gravitational parameter.
Unless otherwise specified, gravitational and mass parameters were obtained from:
- Luzum, Brian et al. “The IAU 2009 System of Astronomical Constants: The Report of the IAU Working Group on Numerical Standards for Fundamental Astronomy.” Celestial Mechanics and Dynamical Astronomy 110.4 (2011): 293–304. Crossref. Web. DOI: 10.1007/s10569-011-9352-4
radii were obtained from:
- Archinal, B. A. et al. “Report of the IAU Working Group on Cartographic Coordinates and Rotational Elements: 2009.” Celestial Mechanics and Dynamical Astronomy 109.2 (2010): 101–135. Crossref. Web. DOI: 10.1007/s10569-010-9320-4
J2 for the Sun was obtained from:
- https://hal.archives-ouvertes.fr/hal-00433235/document (New values of gravitational moments J2 and J4 deduced from helioseismology, Redouane Mecheri et al)
poliastro.coordinates module¶
Functions related to coordinate systems and transformations.
This module complements astropy.coordinates.
-
poliastro.coordinates.body_centered_to_icrs(r, v, source_body, epoch=<Time object: scale='tt' format='jyear_str' value=J2000.000>, rotate_meridian=False)¶ Converts position and velocity body-centered frame to ICRS.
Parameters: - r (Quantity) – Position vector in a body-centered reference frame.
- v (Quantity) – Velocity vector in a body-centered reference frame.
- source_body (Body) – Source body.
- epoch (Time, optional) – Epoch, default to J2000.
- rotate_meridian (bool, optional) – Whether to apply the rotation of the meridian too, default to False.
Returns: r, v – Position and velocity vectors in ICRS.
Return type:
-
poliastro.coordinates.icrs_to_body_centered(r, v, target_body, epoch=<Time object: scale='tt' format='jyear_str' value=J2000.000>, rotate_meridian=False)¶ Converts position and velocity in ICRS to body-centered frame.
Parameters: Returns: r, v – Position and velocity vectors in a body-centered reference frame.
Return type:
-
poliastro.coordinates.inertial_body_centered_to_pqw(r, v, source_body)¶ Converts position and velocity from inertial body-centered frame to perifocal frame.
Parameters: Returns: r_pqw, v_pqw – Position and velocity vectors in ICRS.
Return type:
-
poliastro.coordinates.transform(orbit, frame_orig, frame_dest)¶ Transforms Orbit from one frame to another.
Parameters: - orbit (Orbit) – Orbit to transform
- frame_orig (BaseCoordinateFrame) – Initial frame
- frame_dest (BaseCoordinateFrame) – Final frame
Returns: orbit – Orbit in the new frame
Return type: Orbit
poliastro.cli module¶
Command line functions.
poliastro.examples module¶
Example data.
-
poliastro.examples.iss= 6772 x 6790 km x 51.6 deg (GCRS) orbit around Earth (♁) at epoch 2013-03-18 12:00:00.000 (UTC)¶ ISS orbit example
Taken from Plyades (c) 2012 Helge Eichhorn (MIT License)
-
poliastro.examples.molniya= 6650 x 46550 km x 63.4 deg (GCRS) orbit around Earth (♁) at epoch J2000.000 (TT)¶ Molniya orbit example
-
poliastro.examples.soyuz_gto= 6628 x 42328 km x 6.0 deg (GCRS) orbit around Earth (♁) at epoch J2000.000 (TT)¶ Soyuz geostationary transfer orbit (GTO) example
Taken from Soyuz User’s Manual, issue 2 revision 0
-
poliastro.examples.churi= 1 x 6 AU x 7.0 deg (HCRS) orbit around Sun (☉) at epoch 2015-11-05 12:00:00.000 (UTC)¶ Comet 67P/Churyumov–Gerasimenko orbit example
poliastro.frames module¶
Coordinate frames definitions.
-
class
poliastro.frames.Planes¶ An enumeration.
-
class
poliastro.frames.HeliocentricEclipticJ2000(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶ Heliocentric ecliptic coordinates. These origin of the coordinates are the center of the sun, with the x axis pointing in the direction of the mean equinox of J2000 and the xy-plane in the plane of the ecliptic of J2000 (according to the IAU 1976/1980 obliquity model).
-
class
poliastro.frames.HCRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
class
poliastro.frames.MercuryICRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
class
poliastro.frames.VenusICRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
class
poliastro.frames.MarsICRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
class
poliastro.frames.JupiterICRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
class
poliastro.frames.SaturnICRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
class
poliastro.frames.UranusICRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
class
poliastro.frames.NeptuneICRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
class
poliastro.frames.PlutoICRS(*args, copy=True, representation_type=None, differential_type=None, **kwargs)¶
-
poliastro.frames.get_frame(attractor, plane, obstime=<Time object: scale='tt' format='jyear_str' value=J2000.000>)¶ Returns an appropriate reference frame from an attractor and a plane.
Available planes are Earth equator (parallel to GCRS) and Earth ecliptic. The fundamental direction of both is the equinox of epoch (J2000). An obstime is needed to properly locate the attractor.
Parameters:
poliastro.maneuver module¶
Orbital maneuvers.
-
class
poliastro.maneuver.Maneuver(*impulses)¶ Class to represent a Maneuver.
Each
Maneuverconsists on a list of impulses \(\Delta v_i\) (changes in velocity) each one applied at a certain instant \(t_i\). You can access them directly indexing theManeuverobject itself.>>> man = Maneuver((0 * u.s, [1, 0, 0] * u.km / u.s), ... (10 * u.s, [1, 0, 0] * u.km / u.s)) >>> man[0] (<Quantity 0. s>, <Quantity [1., 0., 0.] km / s>) >>> man.impulses[1] (<Quantity 10. s>, <Quantity [1., 0., 0.] km / s>)
-
__init__(*impulses)¶ Constructor.
Parameters: impulses (list) – List of pairs (delta_time, delta_velocity) Notes
TODO: Fix docstring, *args convention
-
classmethod
impulse(dv)¶ Single impulse at current time.
-
classmethod
hohmann(orbit_i, r_f)¶ Compute a Hohmann transfer between two circular orbits.
-
classmethod
bielliptic(orbit_i, r_b, r_f)¶ Compute a bielliptic transfer between two circular orbits.
-
get_total_time()¶ Returns total time of the maneuver.
-
get_total_cost()¶ Returns total cost of the maneuver.
-
poliastro.threebody package¶
poliastro.threebody.restricted module¶
Circular Restricted 3-Body Problem (CR3BP)
Includes the computation of Lagrange points
-
poliastro.threebody.restricted.lagrange_points(r12, m1, m2)¶ Computes the Lagrangian points of CR3BP.
Computes the Lagrangian points of CR3BP given the distance between two bodies and their masses. It uses the formulation found in Eq. (2.204) of Curtis, Howard. ‘Orbital mechanics for engineering students’. Elsevier, 3rd Edition.
Parameters: Returns: Distance of the Lagrangian points to the main body, projected on the axis main body - secondary body
Return type:
-
poliastro.threebody.restricted.lagrange_points_vec(m1, r1, m2, r2, n)¶ Computes the five Lagrange points in the CR3BP.
Returns the positions in the same frame of reference as r1 and r2 for the five Lagrangian points.
Parameters: Returns: Position of the Lagrange points: [L1, L2, L3, L4, L5] The positions are of type ~astropy.units.Quantity
Return type:
poliastro.threebody.flybys module¶
-
poliastro.threebody.flybys.compute_flyby(v_spacecraft, v_body, k, r_p, theta=<Quantity 0. deg>)¶ Computes outbound velocity after a flyby.
Parameters: - v_spacecraft (Quantity) – Velocity of the spacecraft, relative to the attractor of the body.
- v_body (Quantity) – Velocity of the body, relative to its attractor.
- k (Quantity) – Standard gravitational parameter of the body.
- r_p (Quantity) – Radius of periapsis, measured from the center of the body.
- theta (Quantity, optional) – Aim angle of the B vector, default to 0.
Returns: - v_spacecraft_out (~astropy.units.Quantity) – Outbound velocity of the spacecraft.
- delta (~astropy.units.Quantity) – Turn angle.
poliastro.patched_conics module¶
Patched Conics computations
Contains methods to compute interplanetary trajectories approximating the three body problem with Patched Conics.
-
poliastro.patched_conics.compute_soi(body, a=None)¶ Approximated radius of the Laplace Sphere of Influence (SOI) for a body.
Parameters: - body (~poliastro.bodies.Body) – Astronomical body which the SOI’s radius is computed for.
- a (float, optional) – Semimajor axis of the body’s orbit, default to None (will be computed from ephemerides).
Returns: Approximated radius of the Sphere of Influence (SOI) [m]
Return type:
poliastro.plotting module¶
Plotting utilities.
-
poliastro.plotting.plot(state, label=None, color=None, dark=False)¶ Plots an
Orbitin 2D.For more advanced tuning, use the
OrbitPlotterclass.
-
poliastro.plotting.plot3d(orbit, *, label=None, color=None, dark=False)¶ Plots an
Orbitin 3D.For more advanced tuning, use the
OrbitPlotter3Dclass.
-
class
poliastro.plotting.OrbitPlotter(ax=None, num_points=150, dark=False)¶ OrbitPlotter class.
This class holds the perifocal plane of the first
Orbitplotted in it usingplot(), so all following plots will be projected on that plane. Alternatively, you can callset_frame()to set the frame before plotting.-
__init__(ax=None, num_points=150, dark=False)¶ Constructor.
Parameters:
-
set_frame(p_vec, q_vec, w_vec)¶ Sets perifocal frame.
Raises: ValueError– If the vectors are not a set of mutually orthogonal unit vectors.
-
plot_trajectory(trajectory, *, label=None, color=None)¶ Plots a precomputed trajectory.
Parameters: trajectory (BaseRepresentation, BaseCoordinateFrame) – Trajectory to plot.
-
plot(orbit, label=None, color=None, method=<function mean_motion>)¶ Plots state and osculating orbit in their plane.
-
-
class
poliastro.plotting.OrbitPlotter3D(dark=False)¶ OrbitPlotter3D class.
-
__init__(dark=False)¶ Initialize self. See help(type(self)) for accurate signature.
-
poliastro.util module¶
Function helpers.
-
poliastro.util.circular_velocity(k, a)¶ Compute circular velocity for a given body (k) and semimajor axis (a).
-
poliastro.util.rotate(vector, angle, axis='z')¶ Rotates a vector around axis a right-handed positive angle.
This is just a convenience function around
astropy.coordinates.matrix_utilities.rotation_matrix().Parameters: Note
This performs a so-called active or alibi transformation: rotates the vector while the coordinate system remains unchanged. To do the opposite operation (passive or alias transformation) call the function as
rotate(vec, ax, -angle, unit)or use the convenience functiontransform(), see [1].References
[1] http://en.wikipedia.org/wiki/Rotation_matrix#Ambiguities
-
poliastro.util.transform(vector, angle, axis='z')¶ Rotates a coordinate system around axis a positive right-handed angle.
Note
This is a convenience function, equivalent to
rotate(vec, -angle, axis, unit). Refer to the documentation ofrotate()for further information.
-
poliastro.util.norm(vec)¶ Norm of a Quantity vector that respects units.
Parameters: vec (Quantity) – Vector with units.
-
poliastro.util.time_range(start, *, periods=50, spacing=None, end=None, format=None, scale=None)¶ Generates range of astronomical times.
New in version 0.8.0.
Parameters: - periods (int, optional) – Number of periods, default to 50.
- spacing (Time or Quantity, optional) – Spacing between periods, optional.
- end (Time or equivalent, optional) – End date.
Returns: Array of time values.
Return type: Time