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.4001e+03, -6.9865e+03,  2.7125e+02],
         [ 1.9652e+00,  6.6496e-01,  7.1972e+00]],

        [[ 1.3968e+03, -6.9875e+03,  2.5932e+02],
         [ 1.9677e+00,  6.5222e-01,  7.1977e+00]],

        [[ 1.3811e+03, -6.9925e+03,  2.0192e+02],
         [ 1.9799e+00,  5.9086e-01,  7.1997e+00]],

        [[ 1.3638e+03, -6.9974e+03,  1.3915e+02],
         [ 1.9931e+00,  5.2373e-01,  7.2014e+00]],

        [[ 1.4332e+03, -6.9741e+03,  3.9332e+02],
         [ 1.9387e+00,  7.9535e-01,  7.1910e+00]],

        [[ 1.3698e+03, -6.9957e+03,  1.6106e+02],
         [ 1.9885e+00,  5.4716e-01,  7.2009e+00]],

        [[ 1.4362e+03, -6.9728e+03,  4.0451e+02],
         [ 1.9363e+00,  8.0730e-01,  7.1903e+00]],

        [[ 1.4215e+03, -6.9787e+03,  3.5010e+02],
         [ 1.9482e+00,  7.4920e-01,  7.1935e+00]],

        [[ 1.3864e+03, -6.9909e+03,  2.2114e+02],
         [ 1.9759e+00,  6.1141e-01,  7.1991e+00]],

        [[ 1.4101e+03, -6.9830e+03,  3.0797e+02],
         [ 1.9573e+00,  7.0420e-01,  7.1956e+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.1791e+02,  3.9898e+01,  4.3183e+02],
         [-9.2479e-02,  4.6147e-01, -1.7963e-02]],

        [[ 1.1807e+02,  3.9133e+01,  4.3186e+02],
         [-9.2264e-02,  4.6154e-01, -1.7173e-02]],

        [[ 1.1880e+02,  3.5452e+01,  4.3199e+02],
         [-9.1227e-02,  4.6188e-01, -1.3373e-02]],

        [[ 1.1959e+02,  3.1424e+01,  4.3208e+02],
         [-9.0086e-02,  4.6222e-01, -9.2156e-03]],

        [[ 1.1632e+02,  4.7721e+01,  4.3146e+02],
         [-9.4660e-02,  4.6062e-01, -2.6045e-02]],

        [[ 1.1931e+02,  3.2830e+01,  4.3205e+02],
         [-9.0485e-02,  4.6211e-01, -1.0666e-02]],

        [[ 1.1618e+02,  4.8438e+01,  4.3142e+02],
         [-9.4858e-02,  4.6053e-01, -2.6786e-02]],

        [[ 1.1689e+02,  4.4953e+01,  4.3161e+02],
         [-9.3891e-02,  4.6093e-01, -2.3184e-02]],

        [[ 1.1855e+02,  3.6685e+01,  4.3195e+02],
         [-9.1575e-02,  4.6177e-01, -1.4645e-02]],

        [[ 1.1744e+02,  4.2252e+01,  4.3174e+02],
         [-9.3138e-02,  4.6123e-01, -2.0394e-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.2191e+01, -4.6490e+01,  4.4337e+02],
         [-2.0366e-01, -4.2596e-01, -1.6686e-02]],

        [[ 1.7526e+01, -6.9088e+01,  4.4341e+02],
         [-4.1187e-01, -2.3027e-01, -1.8967e-02]],

        [[-1.0935e+02, -2.6783e+01,  4.4236e+02],
         [ 9.4798e-02, -4.8940e-01, -6.2029e-03]],

        [[-4.1119e+01, -1.0438e+02,  4.4125e+02],
         [ 4.5366e-01, -1.9984e-01, -2.3245e-03]],

        [[-7.4109e+01, -8.7830e+01,  4.4139e+02],
         [ 3.0509e-01, -3.9105e-01, -2.4309e-02]],

        [[ 4.7924e+01, -3.6283e+02, -2.4956e+02],
         [ 3.7008e-02,  2.5648e-01, -3.6598e-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([[[-7.2095e-05,  0.0000e+00,  0.0000e+00,  9.4865e+02,  1.8941e+03,
          -1.5269e+02,  1.9089e+03, -1.3906e+04,  6.9962e+03],
         [-3.5895e-04,  0.0000e+00,  0.0000e+00,  6.4366e+03,  5.2968e+02,
          -3.2237e+01,  5.1846e+02,  7.4568e+04,  1.3681e+03],
         [-4.8110e-04,  0.0000e+00,  0.0000e+00,  7.4456e+03,  6.8556e+03,
           4.0183e+01,  6.8991e+03,  6.7706e+02,  0.0000e+00],
         [-7.1500e-07,  0.0000e+00,  0.0000e+00,  8.6520e-01, -1.4346e+00,
          -7.0722e+00, -1.4433e+00,  1.0062e+01, -5.4063e-01],
         [-4.3241e-07,  0.0000e+00,  0.0000e+00,  4.5301e+00,  7.3580e+00,
          -1.3399e+00,  7.3801e+00,  5.4398e+00,  1.9898e+00],
         [-2.7363e-06,  0.0000e+00,  0.0000e+00,  5.9269e+00, -1.4847e-01,
           2.0513e+00, -1.6386e-01,  3.8275e+01,  0.0000e+00]],

        [[-1.3433e-04,  0.0000e+00,  0.0000e+00,  9.6130e+02,  1.8721e+03,
          -2.5992e+02,  1.8868e+03, -1.3758e+04,  6.9871e+03],
         [-6.1955e-04,  0.0000e+00,  0.0000e+00,  6.5067e+03,  6.4120e+02,
          -5.2554e+01,  6.3031e+02,  7.4669e+04,  1.3982e+03],
         [-8.6426e-04,  0.0000e+00,  0.0000e+00,  7.5346e+03,  6.8525e+03,
           7.1287e+01,  6.8958e+03,  1.2566e+03,  0.0000e+00],
         [-1.2041e-06,  0.0000e+00,  0.0000e+00,  7.9153e-01, -1.4661e+00,
          -7.0687e+00, -1.4749e+00,  9.5558e+00, -6.5736e-01],
         [-7.8943e-07,  0.0000e+00,  0.0000e+00,  4.7138e+00,  7.3486e+00,
          -1.3393e+00,  7.3701e+00,  7.9209e+00,  1.9667e+00],
         [-4.6479e-06,  0.0000e+00,  0.0000e+00,  5.7927e+00, -2.6321e-01,
           2.0504e+00, -2.7930e-01,  3.8147e+01,  0.0000e+00]],

        [[-1.7689e-04,  0.0000e+00,  0.0000e+00,  9.6848e+02,  1.8584e+03,
          -3.2540e+02,  1.8730e+03, -1.3671e+04,  6.9807e+03],
         [-7.8219e-04,  0.0000e+00,  0.0000e+00,  6.5509e+03,  7.0926e+02,
          -6.4961e+01,  6.9857e+02,  7.4750e+04,  1.4163e+03],
         [-1.1161e-03,  0.0000e+00,  0.0000e+00,  7.5879e+03,  6.8497e+03,
           9.0283e+01,  6.8929e+03,  1.6096e+03,  0.0000e+00],
         [-1.4967e-06,  0.0000e+00,  0.0000e+00,  7.4611e-01, -1.4852e+00,
          -7.0657e+00, -1.4940e+00,  9.2374e+00, -7.2860e-01],
         [-1.0283e-06,  0.0000e+00,  0.0000e+00,  4.8238e+00,  7.3419e+00,
          -1.3388e+00,  7.3631e+00,  9.4332e+00,  1.9523e+00],
         [-5.8077e-06,  0.0000e+00,  0.0000e+00,  5.7079e+00, -3.3328e-01,
           2.0496e+00, -3.4980e-01,  3.8035e+01,  0.0000e+00]],

        [[-1.1353e-04,  0.0000e+00,  0.0000e+00,  9.5740e+02,  1.8792e+03,
          -2.2581e+02,  1.8938e+03, -1.3804e+04,  6.9902e+03],
         [-5.3591e-04,  0.0000e+00,  0.0000e+00,  6.4841e+03,  6.0574e+02,
          -4.6091e+01,  5.9474e+02,  7.4633e+04,  1.3886e+03],
         [-7.3844e-04,  0.0000e+00,  0.0000e+00,  7.5065e+03,  6.8537e+03,
           6.1393e+01,  6.8970e+03,  1.0725e+03,  0.0000e+00],
         [-1.0499e-06,  0.0000e+00,  0.0000e+00,  8.1506e-01, -1.4562e+00,
          -7.0700e+00, -1.4649e+00,  9.7189e+00, -6.2024e-01],
         [-6.7124e-07,  0.0000e+00,  0.0000e+00,  4.6559e+00,  7.3518e+00,
          -1.3395e+00,  7.3735e+00,  7.1323e+00,  1.9741e+00],
         [-4.0415e-06,  0.0000e+00,  0.0000e+00,  5.8361e+00, -2.2671e-01,
           2.0507e+00, -2.4258e-01,  3.8195e+01,  0.0000e+00]],

        [[-5.3907e-06,  0.0000e+00,  0.0000e+00,  9.3058e+02,  1.9219e+03,
          -1.3513e+01,  1.9369e+03, -1.4111e+04,  7.0053e+03],
         [-3.0150e-05,  0.0000e+00,  0.0000e+00,  6.3498e+03,  3.8481e+02,
          -5.8715e+00,  3.7315e+02,  7.4493e+04,  1.3287e+03],
         [-3.8132e-05,  0.0000e+00,  0.0000e+00,  7.3272e+03,  6.8571e+03,
          -1.8546e-01,  6.9009e+03, -7.6861e+01,  0.0000e+00],
         [-6.1884e-08,  0.0000e+00,  0.0000e+00,  9.5947e-01, -1.3932e+00,
          -7.0740e+00, -1.4017e+00,  1.0691e+01, -3.8897e-01],
         [-3.3418e-08,  0.0000e+00,  0.0000e+00,  4.2849e+00,  7.3675e+00,
          -1.3401e+00,  7.3903e+00,  2.2123e+00,  2.0190e+00],
         [-2.3423e-07,  0.0000e+00,  0.0000e+00,  6.0920e+00,  4.4947e-04,
           2.0517e+00, -1.4004e-02,  3.8340e+01,  0.0000e+00]],

        [[-2.2609e-04,  0.0000e+00,  0.0000e+00,  9.7568e+02,  1.8436e+03,
          -3.9536e+02,  1.8581e+03, -1.3581e+04,  6.9731e+03],
         [-9.5916e-04,  0.0000e+00,  0.0000e+00,  6.5993e+03,  7.8194e+02,
          -7.8218e+01,  7.7145e+02,  7.4851e+04,  1.4356e+03],
         [-1.4001e-03,  0.0000e+00,  0.0000e+00,  7.6440e+03,  6.8461e+03,
           1.1058e+02,  6.8890e+03,  1.9856e+03,  0.0000e+00],
         [-1.8040e-06,  0.0000e+00,  0.0000e+00,  6.9726e-01, -1.5054e+00,
          -7.0617e+00, -1.5143e+00,  8.8897e+00, -8.0467e-01],
         [-1.3007e-06,  0.0000e+00,  0.0000e+00,  4.9392e+00,  7.3339e+00,
          -1.3381e+00,  7.3547e+00,  1.1046e+01,  1.9368e+00],
         [-7.0401e-06,  0.0000e+00,  0.0000e+00,  5.6148e+00, -4.0814e-01,
           2.0485e+00, -4.2511e-01,  3.7887e+01,  0.0000e+00]],

        [[-1.5139e-05,  0.0000e+00,  0.0000e+00,  9.3364e+02,  1.9174e+03,
          -3.6145e+01,  1.9324e+03, -1.4077e+04,  7.0041e+03],
         [-8.2964e-05,  0.0000e+00,  0.0000e+00,  6.3636e+03,  4.0838e+02,
          -1.0159e+01,  3.9679e+02,  7.4501e+04,  1.3352e+03],
         [-1.0597e-04,  0.0000e+00,  0.0000e+00,  7.3467e+03,  6.8570e+03,
           6.3791e+00,  6.9008e+03,  4.5800e+01,  0.0000e+00],
         [-1.6948e-07,  0.0000e+00,  0.0000e+00,  9.4425e-01, -1.4000e+00,
          -7.0739e+00, -1.4085e+00,  1.0591e+01, -4.1365e-01],
         [-9.3292e-08,  0.0000e+00,  0.0000e+00,  4.3253e+00,  7.3662e+00,
          -1.3401e+00,  7.3888e+00,  2.7376e+00,  2.0143e+00],
         [-6.4264e-07,  0.0000e+00,  0.0000e+00,  6.0658e+00, -2.3767e-02,
           2.0517e+00, -3.8375e-02,  3.8337e+01,  0.0000e+00]],

        [[-5.1773e-05,  0.0000e+00,  0.0000e+00,  9.4377e+02,  1.9020e+03,
          -1.1367e+02,  1.9168e+03, -1.3962e+04,  6.9991e+03],
         [-2.6576e-04,  0.0000e+00,  0.0000e+00,  6.4118e+03,  4.8908e+02,
          -2.4845e+01,  4.7773e+02,  7.4541e+04,  1.3571e+03],
         [-3.5070e-04,  0.0000e+00,  0.0000e+00,  7.4128e+03,  6.8563e+03,
           2.8864e+01,  6.8999e+03,  4.6580e+02,  0.0000e+00],
         [-5.3395e-07,  0.0000e+00,  0.0000e+00,  8.9179e-01, -1.4231e+00,
          -7.0730e+00, -1.4317e+00,  1.0242e+01, -4.9813e-01],
         [-3.1315e-07,  0.0000e+00,  0.0000e+00,  4.4621e+00,  7.3610e+00,
          -1.3400e+00,  7.3833e+00,  4.5357e+00,  1.9981e+00],
         [-2.0371e-06,  0.0000e+00,  0.0000e+00,  5.9742e+00, -1.0671e-01,
           2.0515e+00, -1.2185e-01,  3.8305e+01,  0.0000e+00]],

        [[-1.1549e-04,  0.0000e+00,  0.0000e+00,  9.5778e+02,  1.8785e+03,
          -2.2908e+02,  1.8932e+03, -1.3800e+04,  6.9899e+03],
         [-5.4391e-04,  0.0000e+00,  0.0000e+00,  6.4863e+03,  6.0915e+02,
          -4.6712e+01,  5.9816e+02,  7.4636e+04,  1.3896e+03],
         [-7.5037e-04,  0.0000e+00,  0.0000e+00,  7.5092e+03,  6.8536e+03,
           6.2344e+01,  6.8969e+03,  1.0902e+03,  0.0000e+00],
         [-1.0648e-06,  0.0000e+00,  0.0000e+00,  8.1280e-01, -1.4571e+00,
          -7.0699e+00, -1.4658e+00,  9.7033e+00, -6.2381e-01],
         [-6.8241e-07,  0.0000e+00,  0.0000e+00,  4.6615e+00,  7.3515e+00,
          -1.3395e+00,  7.3732e+00,  7.2081e+00,  1.9734e+00],
         [-4.0998e-06,  0.0000e+00,  0.0000e+00,  5.8319e+00, -2.3022e-01,
           2.0507e+00, -2.4611e-01,  3.8191e+01,  0.0000e+00]],

        [[-1.9200e-04,  0.0000e+00,  0.0000e+00,  9.7080e+02,  1.8538e+03,
          -3.4747e+02,  1.8683e+03, -1.3642e+04,  6.9784e+03],
         [-8.3766e-04,  0.0000e+00,  0.0000e+00,  6.5660e+03,  7.3220e+02,
          -6.9144e+01,  7.2157e+02,  7.4780e+04,  1.4224e+03],
         [-1.2040e-03,  0.0000e+00,  0.0000e+00,  7.6057e+03,  6.8487e+03,
           9.6687e+01,  6.8917e+03,  1.7284e+03,  0.0000e+00],
         [-1.5942e-06,  0.0000e+00,  0.0000e+00,  7.3074e-01, -1.4916e+00,
          -7.0645e+00, -1.5004e+00,  9.1285e+00, -7.5261e-01],
         [-1.1124e-06,  0.0000e+00,  0.0000e+00,  4.8604e+00,  7.3394e+00,
          -1.3386e+00,  7.3605e+00,  9.9425e+00,  1.9475e+00],
         [-6.1974e-06,  0.0000e+00,  0.0000e+00,  5.6788e+00, -3.5690e-01,
           2.0492e+00, -3.7356e-01,  3.7991e+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([[[-3.8039e-04,  0.0000e+00,  0.0000e+00, -1.2569e+03,  8.2936e+02,
           2.7307e+02,  8.0515e+02, -3.1893e+04, -6.4136e+03],
         [-1.8234e-03,  0.0000e+00,  0.0000e+00,  2.2723e+03, -7.3807e+02,
          -1.3087e+02, -7.8411e+02, -6.8137e+04,  3.0752e+03],
         [ 2.8650e-03,  0.0000e+00,  0.0000e+00, -1.3911e+04,  7.0343e+03,
          -4.3427e+01,  7.0193e+03,  1.4391e+03,  0.0000e+00],
         [-1.2835e-06,  0.0000e+00,  0.0000e+00,  3.1082e+00, -3.2416e+00,
           6.6771e+00, -3.2324e+00,  2.2806e+00,  8.2462e-01],
         [-1.1088e-06,  0.0000e+00,  0.0000e+00,  6.7926e+00, -6.7412e+00,
          -3.1504e+00, -6.7411e+00, -8.8704e+00,  8.4602e-01],
         [-4.5337e-06,  0.0000e+00,  0.0000e+00, -5.6834e-01, -3.7390e-01,
          -1.1144e+00, -3.1894e-01,  3.8769e+01,  0.0000e+00]],

        [[-4.5400e-04,  0.0000e+00,  0.0000e+00, -4.8483e+03,  1.8465e+02,
           1.9277e+02,  1.8261e+02, -6.5510e+04, -3.4549e+03],
         [-4.6869e-04,  0.0000e+00,  0.0000e+00, -1.0099e+03, -1.1494e+03,
          -3.4377e+02, -1.1515e+03, -3.7504e+04,  6.2141e+03],
         [ 1.1838e-03,  0.0000e+00,  0.0000e+00, -9.7525e+03,  7.0274e+03,
          -5.8680e+01,  7.0327e+03,  1.9225e+03,  0.0000e+00],
         [-1.0088e-06,  0.0000e+00,  0.0000e+00,  4.5465e+00, -6.5393e+00,
           3.6403e+00, -6.5418e+00, -4.7363e+00,  1.2078e+00],
         [ 6.3273e-07,  0.0000e+00,  0.0000e+00,  1.5556e+00, -3.6351e+00,
          -6.4200e+00, -3.6369e+00, -9.5653e+00,  1.9102e-01],
         [-6.7419e-06,  0.0000e+00,  0.0000e+00,  5.7432e+00, -4.1827e-01,
          -1.1423e+00, -4.1584e-01,  3.8693e+01,  0.0000e+00]],

        [[-1.2301e-03,  0.0000e+00,  0.0000e+00,  6.6858e+02, -1.6052e+03,
           3.5969e+02, -1.6258e+03,  1.2874e+04, -6.7811e+03],
         [ 6.9620e-03,  0.0000e+00,  0.0000e+00, -6.9605e+03, -6.9791e+02,
           7.1678e+01, -6.9734e+02, -6.9942e+04, -1.3860e+03],
         [-3.6694e-04,  0.0000e+00,  0.0000e+00,  2.5351e+03,  6.7087e+03,
           8.9237e+01,  6.7877e+03,  1.7254e+03,  0.0000e+00],
         [ 6.5718e-06,  0.0000e+00,  0.0000e+00, -1.3899e+00,  1.5209e+00,
           7.2284e+00,  1.5321e+00, -7.7422e+00,  7.5586e-01],
         [ 2.4209e-06,  0.0000e+00,  0.0000e+00, -2.3565e+00, -7.4523e+00,
           1.3790e+00, -7.4951e+00, -9.9914e+00, -1.7609e+00],
         [-2.7129e-05,  0.0000e+00,  0.0000e+00,  7.1064e+00, -3.9389e-01,
           1.8672e+00, -4.0561e-01,  3.7379e+01,  0.0000e+00]],

        [[-1.0295e-03,  0.0000e+00,  0.0000e+00,  5.0485e+03, -5.9417e+02,
           3.6882e+01, -5.7176e+02,  6.5068e+04, -2.7865e+03],
         [-1.9392e-04,  0.0000e+00,  0.0000e+00,  6.4551e+02, -1.6110e+03,
           9.0992e+01, -1.6342e+03, -2.8853e+04, -6.3664e+03],
         [ 2.3333e-03,  0.0000e+00,  0.0000e+00, -1.0343e+04,  6.7365e+03,
           2.1074e+01,  6.7844e+03,  3.5430e+02,  0.0000e+00],
         [ 1.7681e-06,  0.0000e+00,  0.0000e+00, -5.7238e+00,  6.9613e+00,
           2.9690e+00,  6.9837e+00, -1.6954e+00,  1.7689e+00],
         [ 9.0765e-07,  0.0000e+00,  0.0000e+00,  1.1776e+00, -3.0383e+00,
           6.7236e+00, -3.0566e+00, -9.6838e+00, -6.1869e-01],
         [-6.0675e-06,  0.0000e+00,  0.0000e+00,  4.8679e+00, -1.3723e-01,
           1.8639e+00, -1.0699e-01,  3.7667e+01,  0.0000e+00]],

        [[-5.7140e-05,  0.0000e+00,  0.0000e+00,  5.0307e+03, -1.3567e+03,
           2.2320e+00, -1.3570e+03,  4.3190e+04, -5.5286e+03],
         [ 4.1725e-05,  0.0000e+00,  0.0000e+00, -3.5669e+03, -1.0461e+03,
           5.4761e+00, -1.0800e+03, -5.6749e+04, -4.2052e+03],
         [ 7.9508e-05,  0.0000e+00,  0.0000e+00, -7.2534e+03,  6.7332e+03,
          -2.5289e+00,  6.8147e+03, -1.2289e+02,  0.0000e+00],
         [ 8.1265e-08,  0.0000e+00,  0.0000e+00, -3.7254e+00,  4.6114e+00,
           5.8592e+00,  4.6338e+00, -7.5126e+00,  1.1658e+00],
         [ 6.8508e-09,  0.0000e+00,  0.0000e+00,  2.2837e+00, -6.0509e+00,
           4.4578e+00, -6.0921e+00, -5.9685e+00, -1.4647e+00],
         [-2.7037e-07,  0.0000e+00,  0.0000e+00,  6.2312e+00, -3.2569e-02,
           1.8673e+00, -4.6230e-03,  3.7857e+01,  0.0000e+00]],

        [[-7.7926e-06,  0.0000e+00,  0.0000e+00, -6.1704e+02,  7.7304e+02,
           6.0705e+03,  7.7013e+02,  7.2151e+03,  4.0001e+03],
         [ 3.4910e-05,  0.0000e+00,  0.0000e+00, -4.0678e+03, -6.1103e+03,
           2.5499e+02, -6.0874e+03,  4.1018e+04, -6.2671e+02],
         [ 4.5114e-05,  0.0000e+00,  0.0000e+00,  6.0200e+03, -3.9374e+03,
           7.9476e+02, -3.9225e+03, -6.8374e+04,  0.0000e+00],
         [-1.6094e-07,  0.0000e+00,  0.0000e+00, -7.7700e-01,  6.3415e-01,
          -3.9707e+00,  6.3295e-01,  4.5467e+00,  6.1712e+00],
         [ 1.2721e-06,  0.0000e+00,  0.0000e+00,  6.1991e+00,  4.0412e+00,
          -1.6597e-01,  4.0337e+00, -3.1967e+01,  7.8009e-01],
         [ 8.1866e-07,  0.0000e+00,  0.0000e+00,  3.9439e+00, -6.1498e+00,
          -5.1790e-01, -6.1382e+00, -2.4528e+01,  0.0000e+00]]])