Cowell’s formulation

For cases where we only study the gravitational forces, solving the Kepler’s equation is enough to propagate the orbit forward in time. However, when we want to take perturbations that deviate from Keplerian forces into account, we need a more complex method to solve our initial value problem: one of them is Cowell’s formulation.

In this formulation we write the two body differential equation separating the Keplerian and the perturbation accelerations:

\[\ddot{\mathbb{r}} = -\frac{\mu}{|\mathbb{r}|^3} \mathbb{r} + \mathbb{a}_d\]

For an in-depth exploration of this topic, still to be integrated in poliastro, check out this Master thesis.

An earlier version of this notebook allowed for more flexibility and interactivity, but was considerably more complex. Future versions of poliastro and plotly might bring back part of that functionality, depending on user feedback. You can still download the older version here.

First example

Let’s setup a very simple example with constant acceleration to visualize the effects on the orbit:

[1]:
from astropy import time
from astropy import units as u

import numpy as np

from poliastro.bodies import Earth
from poliastro.core.propagation import func_twobody
from poliastro.examples import iss
from poliastro.plotting import OrbitPlotter
from poliastro.plotting.orbit.backends import Plotly3D
from poliastro.twobody import Orbit
from poliastro.twobody.propagation import CowellPropagator
from poliastro.twobody.sampling import EpochsArray
from poliastro.util import norm

To provide an acceleration depending on an extra parameter, we can use closures like this one:

[2]:
accel = 2e-5
[3]:
def constant_accel_factory(accel):
    def constant_accel(t0, u, k):
        v = u[3:]
        norm_v = (v[0] ** 2 + v[1] ** 2 + v[2] ** 2) ** 0.5
        return accel * v / norm_v

    return constant_accel
[4]:
def f(t0, state, k):
    du_kep = func_twobody(t0, state, k)
    ax, ay, az = constant_accel_factory(accel)(t0, state, k)
    du_ad = np.array([0, 0, 0, ax, ay, az])

    return du_kep + du_ad
[5]:
times = np.linspace(0, 10 * iss.period, 500)
times
[5]:
$[0,~111.36212,~222.72424,~\dots,~55346.973,~55458.335,~55569.697] \; \mathrm{s}$
[6]:
ephem = iss.to_ephem(
    EpochsArray(iss.epoch + times, method=CowellPropagator(rtol=1e-11, f=f)),
)

And we plot the results:

[7]:
frame = OrbitPlotter(backend=Plotly3D())

frame.set_attractor(Earth)
frame.plot_ephem(ephem, label="ISS")
[7]:
(Scatter3d({
     'line': {'color': 'rgb(31, 119, 180)', 'dash': 'solid', 'width': 5},
     'mode': 'lines',
     'name': 'ISS',
     'showlegend': True,
     'x': array([ 859.07256   , 1671.25960505, 2457.15602056, ..., 1397.77703636,
                 2093.29692966, 2776.01690913]),
     'y': array([-4137.20368   , -3873.06052242, -3547.36485445, ..., -5585.63104917,
                 -5360.76635309, -5102.43872781]),
     'z': array([5295.56871   , 5302.40908447, 5225.0884093 , ..., 7244.2012053 ,
                 7251.65801693, 7213.95475356])
 }),
 None)

Error checking

[8]:
def state_to_vector(ss):
    r, v = ss.rv()
    x, y, z = r.to_value(u.km)
    vx, vy, vz = v.to_value(u.km / u.s)
    return np.array([x, y, z, vx, vy, vz])
[9]:
k = Earth.k.to(u.km**3 / u.s**2).value
[10]:
rtol = 1e-13
full_periods = 2
[11]:
u0 = state_to_vector(iss)
tf = (2 * full_periods + 1) * iss.period / 2

u0, tf
[11]:
(array([ 8.59072560e+02, -4.13720368e+03,  5.29556871e+03,  7.37289205e+00,
         2.08223573e+00,  4.39999794e-01]),
 <Quantity 13892.42425291 s>)
[12]:
iss_f_kep = iss.propagate(tf)
[13]:
iss_f_num = iss.propagate(tf, method=CowellPropagator(rtol=rtol))
[14]:
iss_f_num.r, iss_f_kep.r
[14]:
(<Quantity [ -835.92108005,  4151.60692532, -5303.60427969] km>,
 <Quantity [ -835.92108005,  4151.60692532, -5303.60427969] km>)
[15]:
assert np.allclose(iss_f_num.r, iss_f_kep.r, rtol=rtol, atol=1e-08 * u.km)
assert np.allclose(
    iss_f_num.v, iss_f_kep.v, rtol=rtol, atol=1e-08 * u.km / u.s
)
[16]:
assert np.allclose(iss_f_num.a, iss_f_kep.a, rtol=rtol, atol=1e-08 * u.km)
assert np.allclose(iss_f_num.ecc, iss_f_kep.ecc, rtol=rtol)
assert np.allclose(iss_f_num.inc, iss_f_kep.inc, rtol=rtol, atol=1e-08 * u.rad)
assert np.allclose(
    iss_f_num.raan, iss_f_kep.raan, rtol=rtol, atol=1e-08 * u.rad
)
assert np.allclose(
    iss_f_num.argp, iss_f_kep.argp, rtol=rtol, atol=1e-08 * u.rad
)
assert np.allclose(iss_f_num.nu, iss_f_kep.nu, rtol=rtol, atol=1e-08 * u.rad)

Numerical validation

According to [Edelbaum, 1961], a coplanar, semimajor axis change with tangent thrust is defined by:

\[\frac{\operatorname{d}\!a}{a_0} = 2 \frac{F}{m V_0}\operatorname{d}\!t, \qquad \frac{\Delta{V}}{V_0} = \frac{1}{2} \frac{\Delta{a}}{a_0}\]

So let’s create a new circular orbit and perform the necessary checks, assuming constant mass and thrust (i.e. constant acceleration):

[17]:
orb = Orbit.circular(Earth, 500 << u.km)
tof = 20 * orb.period

ad = constant_accel_factory(1e-7)


def f(t0, state, k):
    du_kep = func_twobody(t0, state, k)
    ax, ay, az = ad(t0, state, k)
    du_ad = np.array([0, 0, 0, ax, ay, az])

    return du_kep + du_ad


orb_final = orb.propagate(tof, method=CowellPropagator(f=f))
[18]:
da_a0 = (orb_final.a - orb.a) / orb.a
da_a0
[18]:
$2.989621 \times 10^{-6} \; \mathrm{\frac{km}{m}}$
[19]:
dv_v0 = abs(norm(orb_final.v) - norm(orb.v)) / norm(orb.v)
2 * dv_v0
[19]:
$0.0029960538 \; \mathrm{}$
[20]:
np.allclose(da_a0, 2 * dv_v0, rtol=1e-2)
[20]:
True

This means we successfully validated the model against an extremely simple orbit transfer with an approximate analytical solution. Notice that the final eccentricity, as originally noticed by Edelbaum, is nonzero:

[21]:
orb_final.ecc
[21]:
$6.6621426 \times 10^{-6} \; \mathrm{}$

References

  • [Edelbaum, 1961] “Propulsion requirements for controllable satellites”