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 (float) – Eccentric anomaly (rad).
  • ecc (float) – Eccentricity.
Returns:

nu – True anomaly (rad).

Return type:

float

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 (float) – True anomaly (rad).
  • ecc (float) – Eccentricity (>1).
Returns:

D – Hyperbolic eccentric anomaly.

Return type:

float

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:
  • nu (float) – True anomaly (rad).
  • ecc (float) – Eccentricity.
Returns:

E – Eccentric anomaly.

Return type:

float

poliastro.twobody.angles.nu_to_F(nu, ecc)

Hyperbolic eccentric anomaly from true anomaly.

Parameters:
  • nu (float) – True anomaly (rad).
  • ecc (float) – Eccentricity (>1).
Returns:

F – Hyperbolic eccentric anomaly.

Return type:

float

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:
  • E (float) – Eccentric anomaly (rad).
  • ecc (float) – Eccentricity.
Returns:

nu – True anomaly (rad).

Return type:

float

poliastro.twobody.angles.F_to_nu(F, ecc)

True anomaly from hyperbolic eccentric anomaly.

Parameters:
  • F (float) – Hyperbolic eccentric anomaly (rad).
  • ecc (float) – Eccentricity (>1).
Returns:

nu – True anomaly (rad).

Return type:

float

poliastro.twobody.angles.M_to_E(M, ecc)

Eccentric anomaly from mean anomaly.

New in version 0.4.0.

Parameters:
  • M (float) – Mean anomaly (rad).
  • ecc (float) – Eccentricity.
Returns:

E – Eccentric anomaly.

Return type:

float

poliastro.twobody.angles.M_to_F(M, ecc)

Hyperbolic eccentric anomaly from mean anomaly.

Parameters:
  • M (float) – Mean anomaly (rad).
  • ecc (float) – Eccentricity (>1).
Returns:

F – Hyperbolic eccentric anomaly.

Return type:

float

poliastro.twobody.angles.M_to_D(M, ecc)

Parabolic eccentric anomaly from mean anomaly.

Parameters:
  • M (float) – Mean anomaly (rad).
  • ecc (float) – Eccentricity (>1).
Returns:

D – Parabolic eccentric anomaly.

Return type:

float

poliastro.twobody.angles.E_to_M(E, ecc)

Mean anomaly from eccentric anomaly.

New in version 0.4.0.

Parameters:
  • E (float) – Eccentric anomaly (rad).
  • ecc (float) – Eccentricity.
Returns:

M – Mean anomaly (rad).

Return type:

float

poliastro.twobody.angles.F_to_M(F, ecc)

Mean anomaly from eccentric anomaly.

Parameters:
  • F (float) – Hyperbolic eccentric anomaly (rad).
  • ecc (float) – Eccentricity (>1).
Returns:

M – Mean anomaly (rad).

Return type:

float

poliastro.twobody.angles.D_to_M(D, ecc)

Mean anomaly from eccentric anomaly.

Parameters:
  • D (float) – Parabolic eccentric anomaly (rad).
  • ecc (float) – Eccentricity.
Returns:

M – Mean anomaly (rad).

Return type:

float

poliastro.twobody.angles.M_to_nu(M, ecc, delta=0.01)

True anomaly from mean anomaly.

New in version 0.4.0.

Parameters:
  • M (float) – Mean anomaly (rad).
  • ecc (float) – Eccentricity.
  • delta (float (optional)) – threshold of near-parabolic regime definition (from Davide Farnocchia et al)
Returns:

nu – True anomaly (rad).

Return type:

float

Examples

>>> nu = M_to_nu(np.radians(30.0), 0.06)
>>> np.rad2deg(nu)
33.673284930211658
poliastro.twobody.angles.nu_to_M(nu, ecc, delta=0.01)

Mean anomaly from true anomaly.

New in version 0.4.0.

Parameters:
  • nu (float) – True anomaly (rad).
  • ecc (float) – Eccentricity.
Returns:

M – Mean anomaly (rad).

Return type:

float

poliastro.twobody.angles.fp_angle(nu, ecc)

Flight path angle.

New in version 0.4.0.

Parameters:
  • nu (float) – True anomaly (rad).
  • ecc (float) – Eccentricity.

Note

Algorithm taken from Vallado 2007, pp. 113.

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)

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.

classmethod from_vectors(attractor, r, v, epoch=<Time object: scale='tdb' format='jyear_str' value=J2000.000>)

Return Orbit from position and velocity vectors.

Parameters:
  • attractor (Body) – Main attractor.
  • r (Quantity) – Position vector wrt attractor center.
  • v (Quantity) – Velocity vector.
  • epoch (Time, optional) – Epoch, default to J2000.
classmethod from_classical(attractor, a, ecc, inc, raan, argp, nu, epoch=<Time object: scale='tdb' format='jyear_str' value=J2000.000>)

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.
classmethod from_equinoctial(attractor, p, f, g, h, k, L, epoch=<Time object: scale='tdb' format='jyear_str' value=J2000.000>)

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.
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='tdb' format='jyear_str' value=J2000.000>)

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.
classmethod parabolic(attractor, p, inc, raan, argp, nu, epoch=<Time object: scale='tdb' format='jyear_str' value=J2000.000>)

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.
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, Time values.
  • 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:

A tuple containing Time and Position vector in each given value.

Return type:

(Time, CartesianRepresentation)

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()
>>> iss.sample(10)
>>> iss.sample([0, 180] * u.deg)
>>> iss.sample([0, 10, 20] * u.minute)
>>> iss.sample([iss.epoch + iss.period / 2])
apply_maneuver(maneuver, intermediate=False)

Returns resulting Orbit after applying maneuver to self.

Optionally return intermediate states (default to False).

Parameters:
  • maneuver (Maneuver) – Maneuver to apply.
  • intermediate (bool, optional) – Return intermediate states, default to False.

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:
  • t0 (float) – Time.
  • u (ndarray) – Six component state vector [x, y, z, vx, vy, vz] (km, km/s).
  • k (float) – Standard gravitational parameter.
  • ad (function(t0, u, k)) – Non Keplerian acceleration (km/s2).
  • ad_kwargs (optional) – perturbation parameters passed to ad
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.integrators module. 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.

class poliastro.twobody.rv.RVState(attractor, r, v)

State defined by its position and velocity vectors.

r

Position vector.

v

Velocity vector.

to_vectors()

Converts to position and velocity vector representation.

to_classical()

Converts to classical orbital elements representation.

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 Orbit given 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 Orbit given 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='DEMO_KEY')

Return Orbit given a SPK-ID.

Retrieve info from NASA NeoWS API, and therefore it only works with NEAs (Near Earth Asteroids).

Parameters:
  • spk_id (str) – SPK-ID number, which is given to each body by JPL.
  • api_key (str) – NASA OPEN APIs key (default: DEMO_KEY)
Returns:

orbit – NEA orbit.

Return type:

Orbit

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.neos.neows.orbit_from_name(name, api_key='DEMO_KEY')

Return Orbit given a name.

Retrieve info from NASA NeoWS API, and therefore it only works with NEAs (Near Earth Asteroids).

Parameters:
  • name (str) – NEA name.
  • api_key (str) – NASA OPEN APIs key (default: DEMO_KEY)
Returns:

orbit – NEA orbit.

Return type:

Orbit

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.

__init__(parent, k, name, symbol=None, R=<Quantity 0. km>, **kwargs)

Constructor.

Parameters:
  • parent (Body) – Central body.
  • k (Quantity) – Standard gravitational parameter.
  • name (str) – Name of the body.
  • symbol (str, optional) – Symbol for the body.
  • R (Quantity, optional) – Radius of the 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:

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='tdb' 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:

tuple (Quantity)

poliastro.coordinates.icrs_to_body_centered(r, v, target_body, epoch=<Time object: scale='tdb' format='jyear_str' value=J2000.000>, rotate_meridian=False)

Converts position and velocity in ICRS to body-centered frame.

Parameters:
  • r (Quantity) – Position vector in ICRS.
  • v (Quantity) – Velocity vector in ICRS.
  • target_body (Body) – Target 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 a body-centered reference frame.

Return type:

tuple (Quantity)

poliastro.coordinates.inertial_body_centered_to_pqw(r, v, source_body)

Converts position and velocity from inertial body-centered frame to perifocal frame.

Parameters:
  • r (Quantity) – Position vector in a inertial body-centered reference frame.
  • v (Quantity) – Velocity vector in a inertial body-centered reference frame.
  • source_body (Body) – Source body.
Returns:

r_pqw, v_pqw – Position and velocity vectors in ICRS.

Return type:

tuple (Quantity)

poliastro.coordinates.transform(orbit, frame_orig, frame_dest)

Transforms Orbit from one frame to another.

Parameters:
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 orbit around Earth (♁)

ISS orbit example

Taken from Plyades (c) 2012 Helge Eichhorn (MIT License)

poliastro.examples.molniya = 6650 x 46550 km x 63.4 deg orbit around Earth (♁)

Molniya orbit example

poliastro.examples.soyuz_gto = 6628 x 42328 km x 6.0 deg orbit around Earth (♁)

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 orbit around Sun (☉)

Comet 67P/Churyumov–Gerasimenko orbit example

poliastro.frames module

Coordinate frames definitions.

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).

poliastro.maneuver module

Orbital maneuvers.

class poliastro.maneuver.Maneuver(*impulses)

Class to represent a Maneuver.

Each Maneuver consists 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 the Maneuver object 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.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:

astropy.units.quantity.Quantity

poliastro.plotting module

Plotting utilities.

poliastro.plotting.plot(state, label=None, color=None)

Plots an Orbit in 2D.

For more advanced tuning, use the OrbitPlotter class.

poliastro.plotting.plot3d(orbit, *, label=None, color=None)

Plots an Orbit in 3D.

For more advanced tuning, use the OrbitPlotter3D class.

class poliastro.plotting.OrbitPlotter(ax=None, num_points=150, dark=False)

OrbitPlotter class.

This class holds the perifocal plane of the first Orbit plotted in it using plot(), so all following plots will be projected on that plane. Alternatively, you can call set_frame() to set the frame before plotting.

__init__(ax=None, num_points=150, dark=False)

Constructor.

Parameters:
  • ax (Axes) – Axes in which to plot. If not given, new ones will be created.
  • num_points (int, optional) – Number of points to use in plots, default to 150.
  • dark (bool, optional) – If set as True, plots the orbit in Dark mode.
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 (CartesianRepresentation) – Trajectory to plot.
set_attractor(attractor)

Sets plotting attractor.

Parameters:attractor (Body) – Central body.
plot(orbit, label=None, color=None, method=<function mean_motion>)

Plots state and osculating orbit in their plane.

class poliastro.plotting.OrbitPlotter3D

OrbitPlotter3D class.

__init__()

Initialize self. See help(type(self)) for accurate signature.

class poliastro.plotting.OrbitPlotter2D

OrbitPlotter2D class.

New in version 0.9.0.

__init__()

Initialize self. See help(type(self)) for accurate signature.

poliastro.plotting.plot_solar_system(outer=True, epoch=None)

Plots the whole solar system in one single call.

New in version 0.9.0.

Parameters:
  • outer (bool, optional) – Whether to print the outer Solar System, default to True.
  • epoch (Time, optional) – Epoch value of the plot, default to J2000.

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:
  • vector (Quantity) – Dimension 3 vector.
  • angle (Quantity) – Angle of rotation.
  • axis (str, optional) – Either ‘x’, ‘y’ or ‘z’.

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 function transform(), 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 of rotate() 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