Partial Derivatives Computation via Autodiff#

In this notebook, we show how to use the autodiff feature of \(\partial\textrm{SGP4}\). Due to the fact that it is written in pytorch, it automatically supports automatic differentiation via torch.autograd.

In this notebook, we show how these partial derivatives can be constructed: for more advanced examples on how to use these gradients for practical applications, see the tutorials on state_transition_matrix_computation, covariance_propagation, graident_based_optimization, orbit_determination.

import dsgp4
import torch

We create a TLE object

#as always, first, we create a TLE object:
tle=[]
tle.append('0 COSMOS 2251 DEB')
tle.append('1 34454U 93036SX  22068.91971155  .00000319  00000-0  11812-3 0  9996')
tle.append('2 34454  74.0583 280.7094 0037596 327.9100  31.9764 14.35844873683320')
tle = dsgp4.tle.TLE(tle)
print(tle)
TLE(
0 COSMOS 2251 DEB
1 34454U 93036SX  22068.91971155  .00000319  00000-0  11812-3 0  9996
2 34454  74.0583 280.7094 0037596 327.9100  31.9764 14.35844873683320
)

Now, as shown in the tle_propagation tutorial, we can propagate the TLE. However, instead of using the standard API, we require torch.autograd to record the operations w.r.t. the time.

Partials with respect to time#

Let’s compute the partials of the \(\partial \textrm{SGP4}\) output w.r.t. the propagation times

Single TLEs#

Let’s first see the case of single TLEs, propagated at various future times:

#let's take a random tensor of 10 tsince elements, where we track the gradients:
tsince=torch.rand((10,),requires_grad=True)
#the state is then:
state_teme = dsgp4.propagate(tle,
                tsinces=tsince,
                initialized=False)
#now, we can see that the gradient is tracked:
print(state_teme)
tensor([[[ 1.4322e+03, -6.9745e+03,  3.8945e+02],
         [ 1.9396e+00,  7.9122e-01,  7.1912e+00]],

        [[ 1.4066e+03, -6.9842e+03,  2.9524e+02],
         [ 1.9600e+00,  6.9061e-01,  7.1962e+00]],

        [[ 1.3925e+03, -6.9889e+03,  2.4360e+02],
         [ 1.9711e+00,  6.3542e-01,  7.1983e+00]],

        [[ 1.4077e+03, -6.9838e+03,  2.9910e+02],
         [ 1.9592e+00,  6.9473e-01,  7.1960e+00]],

        [[ 1.4143e+03, -6.9814e+03,  3.2340e+02],
         [ 1.9540e+00,  7.2068e-01,  7.1948e+00]],

        [[ 1.3538e+03, -6.9999e+03,  1.0329e+02],
         [ 2.0005e+00,  4.8537e-01,  7.2020e+00]],

        [[ 1.3980e+03, -6.9871e+03,  2.6369e+02],
         [ 1.9668e+00,  6.5689e-01,  7.1975e+00]],

        [[ 1.4173e+03, -6.9803e+03,  3.3451e+02],
         [ 1.9515e+00,  7.3256e-01,  7.1943e+00]],

        [[ 1.4380e+03, -6.9721e+03,  4.1104e+02],
         [ 1.9348e+00,  8.1427e-01,  7.1899e+00]],

        [[ 1.4113e+03, -6.9825e+03,  3.1244e+02],
         [ 1.9563e+00,  7.0898e-01,  7.1954e+00]]],
       grad_fn=<TransposeBackward0>)

Let’s now retrieve the partial derivatives of the SGP4 output w.r.t. time.

Since the state is position and velocity (i.e., \([x,y,z,v_x,v_y,v_z]\)), these partials will be all the elements of type:

(1)#\[\begin{equation} \dfrac{d \pmb{x}}{d t}=[\dfrac{dx}{dt}, \dfrac{dy}{dt}, \dfrac{dz}{dt}, \dfrac{d^2x}{dt^2}, \dfrac{d^2y}{dt^2}, \dfrac{d^2z}{dt^2}]^T=[v_x, v_y, v_z, \dfrac{dv_x}{dt}, \dfrac{dv_y}{dt}, \dfrac{dv_z}{dt}]^T \end{equation}\]

Note

One thing to be careful about is that \(\partial\textrm{SGP4}\), mirroring the original \(\textrm{SGP4}\), takes the time in minutes, and returns the state in km and km/s. Hence, the derivatives will have dimensions coherent to these, and to return to SI, conversions have to be made.

partial_derivatives = torch.zeros_like(state_teme)
for i in [0,1]:
    for j in [0,1,2]:
        tsince.grad=None
        state_teme[:,i,j].backward(torch.ones_like(tsince),retain_graph=True)
        partial_derivatives[:,i,j] = tsince.grad

#let's print to screen the partials:
print(partial_derivatives)
tensor([[[ 1.1638e+02,  4.7473e+01,  4.3147e+02],
         [-9.4591e-02,  4.6065e-01, -2.5789e-02]],

        [[ 1.1760e+02,  4.1437e+01,  4.3177e+02],
         [-9.2910e-02,  4.6131e-01, -1.9552e-02]],

        [[ 1.1827e+02,  3.8126e+01,  4.3190e+02],
         [-9.1981e-02,  4.6164e-01, -1.6133e-02]],

        [[ 1.1755e+02,  4.1684e+01,  4.3176e+02],
         [-9.2979e-02,  4.6129e-01, -1.9807e-02]],

        [[ 1.1724e+02,  4.3241e+01,  4.3169e+02],
         [-9.3414e-02,  4.6112e-01, -2.1416e-02]],

        [[ 1.2003e+02,  2.9122e+01,  4.3212e+02],
         [-8.9431e-02,  4.6239e-01, -6.8413e-03]],

        [[ 1.1801e+02,  3.9414e+01,  4.3185e+02],
         [-9.2343e-02,  4.6151e-01, -1.7463e-02]],

        [[ 1.1709e+02,  4.3954e+01,  4.3166e+02],
         [-9.3613e-02,  4.6104e-01, -2.2152e-02]],

        [[ 1.1609e+02,  4.8856e+01,  4.3140e+02],
         [-9.4974e-02,  4.6048e-01, -2.7218e-02]],

        [[ 1.1738e+02,  4.2539e+01,  4.3172e+02],
         [-9.3218e-02,  4.6120e-01, -2.0690e-02]]])

Batch TLEs#

Let’s now see how it works for batch TLEs. The API is basically identical:

#we load 6 TLEs:
inp_file="""0 PSLV DEB
1 35350U 01049QJ  22068.76869562  .00000911  00000-0  24939-3 0  9998
2 35350  98.6033  64.7516 0074531  99.8340 261.1278 14.48029442457561
0 PSLV DEB *
1 35351U 01049QK  22066.70636923  .00002156  00000-0  63479-3 0  9999
2 35351  98.8179  29.5651 0005211  45.5944 314.5671 14.44732274457505
0 SL-18 DEB
1 35354U 93014BD  22068.76520028  .00021929  00000-0  20751-2 0  9995
2 35354  75.7302 100.7819 0059525 350.7978   9.2117 14.92216400847487
0 SL-18 DEB
1 35359U 93014BJ  22068.55187275  .00025514  00000-0  24908-2 0  9992
2 35359  75.7369 156.1582 0054843  50.5279 310.0745 14.91164684775759
0 SL-18 DEB
1 35360U 93014BK  22068.44021735  .00019061  00000-0  20292-2 0  9992
2 35360  75.7343 127.2487 0071107  32.5913 327.9635 14.86997880798827
0 METEOR 2-17 DEB
1 35364U 88005Y   22067.81503681  .00001147  00000-0  84240-3 0  9995
2 35364  82.5500  92.4124 0018834 303.2489 178.0638 13.94853833332534"""
lines=inp_file.splitlines()
#let's create the TLE objects
tles=[]
for i in range(0,len(lines),3):
    data=[]
    data.append(lines[i])
    data.append(lines[i+1])
    data.append(lines[i+2])
    tles.append(dsgp4.tle.TLE(data))
#we also create 9 random times, tracking the gradients:
tsinces=torch.rand((6,),requires_grad=True)
#let's initialize the TLEs:
_,tle_batch=dsgp4.initialize_tle(tles)
#let's propagate the batch:
state_teme = dsgp4.propagate_batch(tle_batch,
                tsinces=tsinces)

Now, let’s retrieve the partial of each TLE, at each propagated time, and store them into a Nx2x3 matrix:

#let's retrieve the partials w.r.t. time:
partial_derivatives = torch.zeros_like(state_teme)
for i in [0,1]:
    for j in [0,1,2]:
        tsinces.grad=None
        state_teme[:,i,j].backward(torch.ones_like(tsinces),retain_graph=True)
        partial_derivatives[:,i,j] = tsinces.grad

#let's print to screen the partials:
print(partial_derivatives)
tensor([[[ 5.8186e+01, -3.3876e+01,  4.4365e+02],
         [-2.0170e-01, -4.2694e-01, -2.1380e-03]],

        [[ 2.1731e+01, -6.6733e+01,  4.4357e+02],
         [-4.1163e-01, -2.3103e-01, -1.3947e-02]],

        [[-1.0886e+02, -2.9339e+01,  4.4232e+02],
         [ 9.5481e-02, -4.8922e-01, -8.9801e-03]],

        [[-2.4058e+01, -1.1180e+02,  4.4080e+02],
         [ 4.5540e-01, -1.9514e-01, -2.2070e-02]],

        [[-8.4714e+01, -7.4104e+01,  4.4192e+02],
         [ 3.0162e-01, -3.9422e-01, -5.8558e-03]],

        [[ 4.6044e+01, -3.7500e+02, -2.3124e+02],
         [ 3.9367e-02,  2.3788e-01, -3.7806e-01]]])

Partials with respect to TLE parameters#

Let’s now tackle the case in which we are interested in the partials of the \(\partial\textrm{SGP4}\) output w.r.t. the TLE parameters.

Single TLEs#

We first tackle the case of single TLE, propagated at multiple times:

In this case, we want the Jacobian of the output state, w.r.t. the following TLE parameters \(\textrm{TLE}=[n,e,i,\Omega,\omega,M,B^*,\dot{n},\ddot{n}]\), where:

  • \(n\) is the mean motion (also known as no_kozai in the original implementation) [rad/minute];

  • \(e\) is the eccentricity [-];

  • \(i\) is the inclination [rad];

  • \(\Omega\) is the right ascension of the ascending node [rad];

  • \(\omega\) is the argument of perigee [rad];

  • \(M\) is the mean anomaly [rad];

  • \(B^*\) is the Bstar parameter [1/earth radii]

  • \(\dot{n}\) mean motion first derivative [radians/\(\textrm{minute}^2\)]

  • \(\ddot{n}\) mean motion second derivative [radians/\(\textrm{minute}^2\)]

(2)#\[\begin{equation} \dfrac{\partial \pmb{x}}{\partial \textrm{TLE}}= \begin{bmatrix} \frac{\partial x}{\partial B^*} & \frac{\partial x}{\partial \dot{n}} & \frac{\partial x}{\partial \ddot{n}} & \frac{\partial x}{\partial e} & \frac{\partial x}{\partial \omega} & \frac{\partial x}{\partial i} & \frac{\partial x}{\partial M} & \frac{\partial x}{\partial n} & \frac{\partial x}{\partial \Omega} \\ \\ \frac{\partial y}{\partial B^*} & \frac{\partial y}{\partial \dot{n}} & \frac{\partial y}{\partial \ddot{n}} & \frac{\partial y}{\partial e} & \frac{\partial y}{\partial \omega} & \frac{\partial y}{\partial i} & \frac{\partial y}{\partial M} & \frac{\partial y}{\partial n} & \frac{\partial y}{\partial \Omega} \\ \\ \frac{\partial z}{\partial B^*} & \frac{\partial z}{\partial \dot{n}} & \frac{\partial z}{\partial \ddot{n}} & \frac{\partial z}{\partial e} & \frac{\partial z}{\partial \omega} & \frac{\partial z}{\partial i} & \frac{\partial z}{\partial M} & \frac{\partial z}{\partial n} & \frac{\partial z}{\partial \Omega} \\ \\ \frac{\partial v_x}{\partial B^*} & \frac{\partial v_x}{\partial \dot{n}} & \frac{\partial v_x}{\partial \ddot{n}} & \frac{\partial v_x}{\partial e} & \frac{\partial v_x}{\partial \omega} & \frac{\partial v_x}{\partial i} & \frac{\partial v_x}{\partial M} & \frac{\partial v_x}{\partial n} & \frac{\partial v_x}{\partial \Omega} \\ \\ \frac{\partial v_y}{\partial B^*} & \frac{\partial v_y}{\partial \dot{n}} & \frac{\partial v_y}{\partial \ddot{n}} & \frac{\partial v_y}{\partial e} & \frac{\partial v_y}{\partial \omega} & \frac{\partial v_y}{\partial i} & \frac{\partial v_y}{\partial M} & \frac{\partial v_y}{\partial n} & \frac{\partial v_y}{\partial \Omega} \\ \\ \frac{\partial v_z}{\partial B^*} & \frac{\partial v_z}{\partial \dot{n}} & \frac{\partial v_z}{\partial \ddot{n}} & \frac{\partial v_z}{\partial e} & \frac{\partial v_z}{\partial \omega} & \frac{\partial v_z}{\partial i} & \frac{\partial v_z}{\partial M} & \frac{\partial v_z}{\partial n} & \frac{\partial v_z}{\partial \Omega} \\ \\ \frac{\partial \dot{n}}{\partial B^*} & \frac{\partial \dot{n}}{\partial \dot{n}} & \frac{\partial \dot{n}}{\partial \ddot{n}} & \frac{\partial \dot{n}}{\partial e} & \frac{\partial \dot{n}}{\partial \omega} & \frac{\partial \dot{n}}{\partial i} & \frac{\partial \dot{n}}{\partial M} & \frac{\partial \dot{n}}{\partial n} & \frac{\partial \dot{n}}{\partial \Omega} \\ \\ \frac{\partial \ddot{n}}{\partial B^*} & \frac{\partial \ddot{n}}{\partial \dot{n}} & \frac{\partial \ddot{n}}{\partial \ddot{n}} & \frac{\partial \ddot{n}}{\partial e} & \frac{\partial \ddot{n}}{\partial \omega} & \frac{\partial \ddot{n}}{\partial i} & \frac{\partial \ddot{n}}{\partial M} & \frac{\partial \ddot{n}}{\partial n} & \frac{\partial \ddot{n}}{\partial \Omega} \end{bmatrix} \end{equation}\]
#as always, first, we create a TLE object:
tle=[]
tle.append('0 COSMOS 2251 DEB')
tle.append('1 34454U 93036SX  22068.91971155  .00000319  00000-0  11812-3 0  9996')
tle.append('2 34454  74.0583 280.7094 0037596 327.9100  31.9764 14.35844873683320')
tle = dsgp4.tle.TLE(tle)
print(tle)
tle_elements=dsgp4.initialize_tle(tle,with_grad=True)
TLE(
0 COSMOS 2251 DEB
1 34454U 93036SX  22068.91971155  .00000319  00000-0  11812-3 0  9996
2 34454  74.0583 280.7094 0037596 327.9100  31.9764 14.35844873683320
)
#let's select 10 random times:
tsince=torch.rand((10,))
#and let's propagate:
state_teme=dsgp4.propagate(tle,tsince)
#now we can build the partial derivatives matrix, of shape Nx6x9 (N is the number of tsince elements, 6 is the number of elements in the state vector, and 9 is the number of elements in the TLE):
partial_derivatives = torch.zeros((len(tsince),6,9))
for k in range(len(tsince)):
    for i in range(6):
        tle_elements.grad=None
        state_teme[k].flatten()[i].backward(retain_graph=True)
        partial_derivatives[k,i,:] = tle_elements.grad
#let's print them to screen:
print(partial_derivatives)
tensor([[[-5.8691e-05,  0.0000e+00,  0.0000e+00,  9.4548e+02,  1.8992e+03,
          -1.2722e+02,  1.9141e+03, -1.3943e+04,  6.9981e+03],
         [-2.9804e-04,  0.0000e+00,  0.0000e+00,  6.4203e+03,  5.0318e+02,
          -2.7413e+01,  4.9188e+02,  7.4550e+04,  1.3610e+03],
         [-3.9545e-04,  0.0000e+00,  0.0000e+00,  7.4242e+03,  6.8561e+03,
           3.2796e+01,  6.8997e+03,  5.3921e+02,  0.0000e+00],
         [-5.9703e-07,  0.0000e+00,  0.0000e+00,  8.8257e-01, -1.4271e+00,
          -7.0727e+00, -1.4357e+00,  1.0180e+01, -5.1289e-01],
         [-3.5393e-07,  0.0000e+00,  0.0000e+00,  4.4858e+00,  7.3600e+00,
          -1.3399e+00,  7.3822e+00,  4.8499e+00,  1.9952e+00],
         [-2.2802e-06,  0.0000e+00,  0.0000e+00,  5.9579e+00, -1.2122e-01,
           2.0514e+00, -1.3644e-01,  3.8296e+01,  0.0000e+00]],

        [[-7.2297e-05,  0.0000e+00,  0.0000e+00,  9.4870e+02,  1.8940e+03,
          -1.5306e+02,  1.9088e+03, -1.3906e+04,  6.9962e+03],
         [-3.5986e-04,  0.0000e+00,  0.0000e+00,  6.4368e+03,  5.3007e+02,
          -3.2309e+01,  5.1885e+02,  7.4568e+04,  1.3683e+03],
         [-4.8238e-04,  0.0000e+00,  0.0000e+00,  7.4459e+03,  6.8556e+03,
           4.0292e+01,  6.8991e+03,  6.7909e+02,  0.0000e+00],
         [-7.1674e-07,  0.0000e+00,  0.0000e+00,  8.6494e-01, -1.4347e+00,
          -7.0722e+00, -1.4434e+00,  1.0060e+01, -5.4104e-01],
         [-4.3359e-07,  0.0000e+00,  0.0000e+00,  4.5307e+00,  7.3580e+00,
          -1.3399e+00,  7.3801e+00,  5.4485e+00,  1.9897e+00],
         [-2.7431e-06,  0.0000e+00,  0.0000e+00,  5.9265e+00, -1.4887e-01,
           2.0513e+00, -1.6426e-01,  3.8275e+01,  0.0000e+00]],

        [[-9.4827e-05,  0.0000e+00,  0.0000e+00,  9.5362e+02,  1.8857e+03,
          -1.9375e+02,  1.9005e+03, -1.3849e+04,  6.9929e+03],
         [-4.5794e-04,  0.0000e+00,  0.0000e+00,  6.4631e+03,  5.7240e+02,
          -4.0017e+01,  5.6130e+02,  7.4602e+04,  1.3797e+03],
         [-6.2353e-04,  0.0000e+00,  0.0000e+00,  7.4799e+03,  6.8546e+03,
           5.2094e+01,  6.8981e+03,  8.9920e+02,  0.0000e+00],
         [-9.0377e-07,  0.0000e+00,  0.0000e+00,  8.3709e-01, -1.4467e+00,
          -7.0711e+00, -1.4554e+00,  9.8705e+00, -5.8534e-01],
         [-5.6408e-07,  0.0000e+00,  0.0000e+00,  4.6010e+00,  7.3546e+00,
          -1.3397e+00,  7.3765e+00,  6.3906e+00,  1.9810e+00],
         [-3.4701e-06,  0.0000e+00,  0.0000e+00,  5.8762e+00, -1.9241e-01,
           2.0510e+00, -2.0807e-01,  3.8234e+01,  0.0000e+00]],

        [[-9.8537e-05,  0.0000e+00,  0.0000e+00,  9.5439e+02,  1.8844e+03,
          -2.0022e+02,  1.8991e+03, -1.3840e+04,  6.9924e+03],
         [-4.7364e-04,  0.0000e+00,  0.0000e+00,  6.4673e+03,  5.7913e+02,
          -4.1244e+01,  5.6806e+02,  7.4608e+04,  1.3815e+03],
         [-6.4648e-04,  0.0000e+00,  0.0000e+00,  7.4853e+03,  6.8544e+03,
           5.3972e+01,  6.8979e+03,  9.3420e+02,  0.0000e+00],
         [-9.3337e-07,  0.0000e+00,  0.0000e+00,  8.3265e-01, -1.4486e+00,
          -7.0709e+00, -1.4573e+00,  9.8400e+00, -5.9239e-01],
         [-5.8541e-07,  0.0000e+00,  0.0000e+00,  4.6121e+00,  7.3541e+00,
          -1.3397e+00,  7.3759e+00,  6.5405e+00,  1.9796e+00],
         [-3.5856e-06,  0.0000e+00,  0.0000e+00,  5.8682e+00, -1.9933e-01,
           2.0510e+00, -2.1504e-01,  3.8227e+01,  0.0000e+00]],

        [[-2.4530e-04,  0.0000e+00,  0.0000e+00,  9.7823e+02,  1.8381e+03,
          -4.2129e+02,  1.8525e+03, -1.3548e+04,  6.9701e+03],
         [-1.0256e-03,  0.0000e+00,  0.0000e+00,  6.6175e+03,  8.0887e+02,
          -8.3132e+01,  7.9846e+02,  7.4893e+04,  1.4427e+03],
         [-1.5092e-03,  0.0000e+00,  0.0000e+00,  7.6646e+03,  6.8445e+03,
           1.1810e+02,  6.8874e+03,  2.1246e+03,  0.0000e+00],
         [-1.9166e-06,  0.0000e+00,  0.0000e+00,  6.7907e-01, -1.5129e+00,
          -7.0601e+00, -1.5218e+00,  8.7588e+00, -8.3285e-01],
         [-1.4062e-06,  0.0000e+00,  0.0000e+00,  4.9815e+00,  7.3308e+00,
          -1.3378e+00,  7.3514e+00,  1.1643e+01,  1.9310e+00],
         [-7.4951e-06,  0.0000e+00,  0.0000e+00,  5.5796e+00, -4.3589e-01,
           2.0480e+00, -4.5302e-01,  3.7825e+01,  0.0000e+00]],

        [[-7.6894e-05,  0.0000e+00,  0.0000e+00,  9.4974e+02,  1.8923e+03,
          -1.6156e+02,  1.9071e+03, -1.3894e+04,  6.9955e+03],
         [-3.8028e-04,  0.0000e+00,  0.0000e+00,  6.4423e+03,  5.3892e+02,
          -3.3920e+01,  5.2772e+02,  7.4575e+04,  1.3706e+03],
         [-5.1145e-04,  0.0000e+00,  0.0000e+00,  7.4531e+03,  6.8554e+03,
           4.2758e+01,  6.8989e+03,  7.2511e+02,  0.0000e+00],
         [-7.5597e-07,  0.0000e+00,  0.0000e+00,  8.5913e-01, -1.4373e+00,
          -7.0720e+00, -1.4459e+00,  1.0021e+01, -5.5030e-01],
         [-4.6035e-07,  0.0000e+00,  0.0000e+00,  4.5455e+00,  7.3573e+00,
          -1.3398e+00,  7.3794e+00,  5.6455e+00,  1.9879e+00],
         [-2.8952e-06,  0.0000e+00,  0.0000e+00,  5.9160e+00, -1.5797e-01,
           2.0512e+00, -1.7342e-01,  3.8267e+01,  0.0000e+00]],

        [[-1.2716e-04,  0.0000e+00,  0.0000e+00,  9.5999e+02,  1.8745e+03,
          -2.4833e+02,  1.8892e+03, -1.3773e+04,  6.9882e+03],
         [-5.9105e-04,  0.0000e+00,  0.0000e+00,  6.4990e+03,  6.2915e+02,
          -5.0358e+01,  6.1823e+02,  7.4657e+04,  1.3949e+03],
         [-8.2111e-04,  0.0000e+00,  0.0000e+00,  7.5251e+03,  6.8529e+03,
           6.7926e+01,  6.8962e+03,  1.1941e+03,  0.0000e+00],
         [-1.1519e-06,  0.0000e+00,  0.0000e+00,  7.9953e-01, -1.4627e+00,
          -7.0692e+00, -1.4715e+00,  9.6114e+00, -6.4475e-01],
         [-7.4880e-07,  0.0000e+00,  0.0000e+00,  4.6942e+00,  7.3497e+00,
          -1.3394e+00,  7.3713e+00,  7.6531e+00,  1.9692e+00],
         [-4.4421e-06,  0.0000e+00,  0.0000e+00,  5.8075e+00, -2.5081e-01,
           2.0505e+00, -2.6683e-01,  3.8164e+01,  0.0000e+00]],

        [[-7.3979e-05,  0.0000e+00,  0.0000e+00,  9.4908e+02,  1.8934e+03,
          -1.5619e+02,  1.9082e+03, -1.3901e+04,  6.9959e+03],
         [-3.6735e-04,  0.0000e+00,  0.0000e+00,  6.4388e+03,  5.3332e+02,
          -3.2900e+01,  5.2211e+02,  7.4571e+04,  1.3691e+03],
         [-4.9303e-04,  0.0000e+00,  0.0000e+00,  7.4485e+03,  6.8555e+03,
           4.1198e+01,  6.8991e+03,  6.9600e+02,  0.0000e+00],
         [-7.3116e-07,  0.0000e+00,  0.0000e+00,  8.6281e-01, -1.4357e+00,
          -7.0721e+00, -1.4443e+00,  1.0046e+01, -5.4444e-01],
         [-4.4339e-07,  0.0000e+00,  0.0000e+00,  4.5362e+00,  7.3577e+00,
          -1.3399e+00,  7.3798e+00,  5.5209e+00,  1.9891e+00],
         [-2.7990e-06,  0.0000e+00,  0.0000e+00,  5.9226e+00, -1.5221e-01,
           2.0513e+00, -1.6763e-01,  3.8272e+01,  0.0000e+00]],

        [[-3.4517e-05,  0.0000e+00,  0.0000e+00,  9.3924e+02,  1.9090e+03,
          -7.8461e+01,  1.9239e+03, -1.4014e+04,  7.0014e+03],
         [-1.8237e-04,  0.0000e+00,  0.0000e+00,  6.3897e+03,  4.5244e+02,
          -1.8175e+01,  4.4098e+02,  7.4520e+04,  1.3472e+03],
         [-2.3719e-04,  0.0000e+00,  0.0000e+00,  7.3829e+03,  6.8568e+03,
           1.8653e+01,  6.9004e+03,  2.7510e+02,  0.0000e+00],
         [-3.6922e-07,  0.0000e+00,  0.0000e+00,  9.1567e-01, -1.4126e+00,
          -7.0735e+00, -1.4212e+00,  1.0402e+01, -4.5977e-01],
         [-2.1048e-07,  0.0000e+00,  0.0000e+00,  4.4002e+00,  7.3635e+00,
          -1.3400e+00,  7.3859e+00,  3.7194e+00,  2.0055e+00],
         [-1.4047e-06,  0.0000e+00,  0.0000e+00,  6.0162e+00, -6.9045e-02,
           2.0516e+00, -8.3940e-02,  3.8324e+01,  0.0000e+00]],

        [[-3.6280e-05,  0.0000e+00,  0.0000e+00,  9.3972e+02,  1.9083e+03,
          -8.2157e+01,  1.9232e+03, -1.4008e+04,  7.0012e+03],
         [-1.9110e-04,  0.0000e+00,  0.0000e+00,  6.3920e+03,  4.5628e+02,
          -1.8876e+01,  4.4484e+02,  7.4522e+04,  1.3482e+03],
         [-2.4893e-04,  0.0000e+00,  0.0000e+00,  7.3860e+03,  6.8567e+03,
           1.9725e+01,  6.9004e+03,  2.9513e+02,  0.0000e+00],
         [-3.8658e-07,  0.0000e+00,  0.0000e+00,  9.1317e-01, -1.4137e+00,
          -7.0735e+00, -1.4223e+00,  1.0385e+01, -4.6380e-01],
         [-2.2104e-07,  0.0000e+00,  0.0000e+00,  4.4068e+00,  7.3632e+00,
          -1.3400e+00,  7.3857e+00,  3.8051e+00,  2.0047e+00],
         [-1.4712e-06,  0.0000e+00,  0.0000e+00,  6.0118e+00, -7.3000e-02,
           2.0516e+00, -8.7919e-02,  3.8322e+01,  0.0000e+00]]])

Batch TLEs:#

As for the time derivatives, the API stays practically identical:

#we load 6 TLEs:
inp_file="""0 PSLV DEB
1 35350U 01049QJ  22068.76869562  .00000911  00000-0  24939-3 0  9998
2 35350  98.6033  64.7516 0074531  99.8340 261.1278 14.48029442457561
0 PSLV DEB *
1 35351U 01049QK  22066.70636923  .00002156  00000-0  63479-3 0  9999
2 35351  98.8179  29.5651 0005211  45.5944 314.5671 14.44732274457505
0 SL-18 DEB
1 35354U 93014BD  22068.76520028  .00021929  00000-0  20751-2 0  9995
2 35354  75.7302 100.7819 0059525 350.7978   9.2117 14.92216400847487
0 SL-18 DEB
1 35359U 93014BJ  22068.55187275  .00025514  00000-0  24908-2 0  9992
2 35359  75.7369 156.1582 0054843  50.5279 310.0745 14.91164684775759
0 SL-18 DEB
1 35360U 93014BK  22068.44021735  .00019061  00000-0  20292-2 0  9992
2 35360  75.7343 127.2487 0071107  32.5913 327.9635 14.86997880798827
0 METEOR 2-17 DEB
1 35364U 88005Y   22067.81503681  .00001147  00000-0  84240-3 0  9995
2 35364  82.5500  92.4124 0018834 303.2489 178.0638 13.94853833332534"""
lines=inp_file.splitlines()
#let's create the TLE objects
tles=[]
for i in range(0,len(lines),3):
    data=[]
    data.append(lines[i])
    data.append(lines[i+1])
    data.append(lines[i+2])
    tles.append(dsgp4.tle.TLE(data))
#we also create 6 random times, tracking the gradients:
tsinces=torch.rand((6,))
#let's now initialize the TLEs, activating the gradient tracking for the TLE parameters:
tle_elements,tle_batch=dsgp4.initialize_tle(tles,with_grad=True)
#let's now propagate the batch of TLEs:
state_teme = dsgp4.propagate_batch(tle_batch,tsinces)

Finally, we can build the matrix that contains the partial of the SGP4 output w.r.t. the TLE parameters, for each TLE:

#now we can build the partial derivatives matrix, of shape Nx6x9 (N is the number of tsince elements, 6 is the number of elements in the state vector, and 9 is the number of elements in the TLE):
partial_derivatives = torch.zeros((len(tsinces),6,9))
for k in range(len(tsinces)):
    for i in range(6):
        tle_elements[k].grad=None
        state_teme[k].flatten()[i].backward(retain_graph=True)
        partial_derivatives[k,i,:] = tle_elements[k].grad
#let's print them to screen:
print(partial_derivatives)
tensor([[[-4.1427e-04,  0.0000e+00,  0.0000e+00, -1.2462e+03,  8.1817e+02,
           2.9610e+02,  7.9399e+02, -3.1886e+04, -6.4107e+03],
         [-1.9798e-03,  0.0000e+00,  0.0000e+00,  2.2957e+03, -7.6132e+02,
          -1.4174e+02, -8.0735e+02, -6.8169e+04,  3.0781e+03],
         [ 3.1036e-03,  0.0000e+00,  0.0000e+00, -1.3913e+04,  7.0329e+03,
          -4.7270e+01,  7.0182e+03,  1.5728e+03,  0.0000e+00],
         [-1.3947e-06,  0.0000e+00,  0.0000e+00,  3.1175e+00, -3.2448e+00,
           6.6760e+00, -3.2357e+00,  2.0305e+00,  8.4908e-01],
         [-1.2033e-06,  0.0000e+00,  0.0000e+00,  6.7932e+00, -6.7384e+00,
          -3.1499e+00, -6.7387e+00, -9.3853e+00,  8.3428e-01],
         [-4.9311e-06,  0.0000e+00,  0.0000e+00, -5.1465e-01, -4.0081e-01,
          -1.1143e+00, -3.4582e-01,  3.8725e+01,  0.0000e+00]],

        [[-2.7866e-06,  0.0000e+00,  0.0000e+00, -5.1025e+03,  5.3039e+02,
          -1.7839e-01,  5.2848e+02, -6.5463e+04, -3.5134e+03],
         [-3.1692e-06,  0.0000e+00,  0.0000e+00, -1.1022e+03, -9.5528e+02,
          -3.4949e+00, -9.5728e+02, -3.7110e+04,  6.1943e+03],
         [ 8.9004e-06,  0.0000e+00,  0.0000e+00, -1.0042e+04,  7.0386e+03,
           1.8614e+00,  7.0438e+03, -1.4105e+02,  0.0000e+00],
         [-7.3667e-09,  0.0000e+00,  0.0000e+00,  5.0399e+00, -6.5174e+00,
           3.6462e+00, -6.5195e+00,  2.9567e+00,  1.0039e+00],
         [ 3.5670e-09,  0.0000e+00,  0.0000e+00,  1.9344e+00, -3.6961e+00,
          -6.4300e+00, -3.6979e+00, -5.2762e+00,  5.5418e-01],
         [-4.3335e-08,  0.0000e+00,  0.0000e+00,  5.1877e+00, -5.4955e-03,
          -1.1439e+00, -2.7474e-03,  3.9115e+01,  0.0000e+00]],

        [[-1.1809e-03,  0.0000e+00,  0.0000e+00,  6.7181e+02, -1.6087e+03,
           3.4304e+02, -1.6293e+03,  1.2892e+04, -6.7828e+03],
         [ 6.6357e-03,  0.0000e+00,  0.0000e+00, -6.9552e+03, -6.8075e+02,
           6.8503e+01, -6.8008e+02, -6.9920e+04, -1.3820e+03],
         [-3.1695e-04,  0.0000e+00,  0.0000e+00,  2.5187e+03,  6.7096e+03,
           8.4937e+01,  6.7887e+03,  1.6393e+03,  0.0000e+00],
         [ 6.2789e-06,  0.0000e+00,  0.0000e+00, -1.3993e+00,  1.5165e+00,
           7.2294e+00,  1.5276e+00, -7.8240e+00,  7.3714e-01],
         [ 2.2516e-06,  0.0000e+00,  0.0000e+00, -2.3205e+00, -7.4542e+00,
           1.3791e+00, -7.4970e+00, -9.6093e+00, -1.7648e+00],
         [-2.5875e-05,  0.0000e+00,  0.0000e+00,  7.1164e+00, -3.7532e-01,
           1.8674e+00, -3.8683e-01,  3.7414e+01,  0.0000e+00]],

        [[-1.6710e-03,  0.0000e+00,  0.0000e+00,  5.0004e+03, -5.3535e+02,
           6.1957e+01, -5.1274e+02,  6.5059e+04, -2.7714e+03],
         [-3.1025e-04,  0.0000e+00,  0.0000e+00,  6.5516e+02, -1.6366e+03,
           1.4778e+02, -1.6599e+03, -2.8937e+04, -6.3714e+03],
         [ 3.7712e-03,  0.0000e+00,  0.0000e+00, -1.0301e+04,  6.7351e+03,
           3.6816e+01,  6.7832e+03,  6.7228e+02,  0.0000e+00],
         [ 2.8362e-06,  0.0000e+00,  0.0000e+00, -5.6490e+00,  6.9672e+00,
           2.9684e+00,  6.9901e+00, -3.8221e-01,  1.7969e+00],
         [ 1.5265e-06,  0.0000e+00,  0.0000e+00,  1.1164e+00, -3.0221e+00,
           6.7224e+00, -3.0405e+00, -1.0252e+01, -5.5473e-01],
         [-9.9864e-06,  0.0000e+00,  0.0000e+00,  4.9693e+00, -2.0506e-01,
           1.8636e+00, -1.7532e-01,  3.7622e+01,  0.0000e+00]],

        [[-3.4576e-03,  0.0000e+00,  0.0000e+00,  4.9057e+03, -1.1961e+03,
           2.0502e+02, -1.1956e+03,  4.2992e+04, -5.4843e+03],
         [ 2.6095e-03,  0.0000e+00,  0.0000e+00, -3.4956e+03, -1.2548e+03,
           1.5976e+02, -1.2901e+03, -5.7037e+04, -4.2529e+03],
         [ 4.6107e-03,  0.0000e+00,  0.0000e+00, -7.0325e+03,  6.7272e+03,
           6.2106e+01,  6.8097e+03,  1.1855e+03,  0.0000e+00],
         [ 4.7694e-06,  0.0000e+00,  0.0000e+00, -3.4812e+00,  4.6646e+00,
           5.8547e+00,  4.6889e+00, -3.9112e+00,  1.3929e+00],
         [ 9.3343e-07,  0.0000e+00,  0.0000e+00,  1.8445e+00, -6.0042e+00,
           4.4546e+00, -6.0465e+00, -1.0625e+01, -1.2903e+00],
         [-1.7132e-05,  0.0000e+00,  0.0000e+00,  6.5126e+00, -3.1092e-01,
           1.8661e+00, -2.8642e-01,  3.7666e+01,  0.0000e+00]],

        [[-2.4021e-05,  0.0000e+00,  0.0000e+00, -6.4060e+02,  7.9145e+02,
           5.9504e+03,  7.8850e+02,  7.3554e+03,  4.1808e+03],
         [ 1.3217e-04,  0.0000e+00,  0.0000e+00, -3.8880e+03, -5.9880e+03,
           2.4996e+02, -5.9654e+03,  4.0114e+04, -6.0336e+02],
         [ 1.3422e-04,  0.0000e+00,  0.0000e+00,  6.1419e+03, -4.1175e+03,
           7.7910e+02, -4.1022e+03, -6.9158e+04,  0.0000e+00],
         [-3.4035e-07,  0.0000e+00,  0.0000e+00, -8.1373e-01,  6.1055e-01,
          -4.1524e+00,  6.0944e-01,  4.9368e+00,  6.0475e+00],
         [ 2.6143e-06,  0.0000e+00,  0.0000e+00,  5.9487e+00,  4.2244e+00,
          -1.7375e-01,  4.2162e+00, -2.9134e+01,  7.9869e-01],
         [ 1.7554e-06,  0.0000e+00,  0.0000e+00,  4.3058e+00, -6.0282e+00,
          -5.4154e-01, -6.0170e+00, -2.8435e+01,  0.0000e+00]]])