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.3804e+03, -6.9927e+03,  1.9941e+02],
         [ 1.9805e+00,  5.8818e-01,  7.1998e+00]],

        [[ 1.3693e+03, -6.9959e+03,  1.5917e+02],
         [ 1.9889e+00,  5.4514e-01,  7.2009e+00]],

        [[ 1.3919e+03, -6.9891e+03,  2.4137e+02],
         [ 1.9716e+00,  6.3304e-01,  7.1984e+00]],

        [[ 1.3802e+03, -6.9928e+03,  1.9855e+02],
         [ 1.9806e+00,  5.8726e-01,  7.1998e+00]],

        [[ 1.4227e+03, -6.9782e+03,  3.5449e+02],
         [ 1.9472e+00,  7.5389e-01,  7.1932e+00]],

        [[ 1.3498e+03, -7.0008e+03,  8.8896e+01],
         [ 2.0035e+00,  4.6996e-01,  7.2023e+00]],

        [[ 1.3533e+03, -7.0000e+03,  1.0125e+02],
         [ 2.0010e+00,  4.8318e-01,  7.2021e+00]],

        [[ 1.4273e+03, -6.9765e+03,  3.7138e+02],
         [ 1.9435e+00,  7.7193e-01,  7.1923e+00]],

        [[ 1.3452e+03, -7.0019e+03,  7.2326e+01],
         [ 2.0069e+00,  4.5223e-01,  7.2025e+00]],

        [[ 1.4423e+03, -6.9703e+03,  4.2694e+02],
         [ 1.9313e+00,  8.3124e-01,  7.1889e+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.1883e+02,  3.5291e+01,  4.3199e+02],
         [-9.1182e-02,  4.6190e-01, -1.3206e-02]],

        [[ 1.1934e+02,  3.2709e+01,  4.3206e+02],
         [-9.0451e-02,  4.6212e-01, -1.0541e-02]],

        [[ 1.1830e+02,  3.7982e+01,  4.3190e+02],
         [-9.1941e-02,  4.6165e-01, -1.5985e-02]],

        [[ 1.1884e+02,  3.5236e+01,  4.3199e+02],
         [-9.1166e-02,  4.6190e-01, -1.3150e-02]],

        [[ 1.1683e+02,  4.5234e+01,  4.3159e+02],
         [-9.3969e-02,  4.6090e-01, -2.3475e-02]],

        [[ 1.2021e+02,  2.8198e+01,  4.3214e+02],
         [-8.9167e-02,  4.6246e-01, -5.8878e-03]],

        [[ 1.2006e+02,  2.8991e+01,  4.3213e+02],
         [-8.9394e-02,  4.6240e-01, -6.7062e-03]],

        [[ 1.1661e+02,  4.6316e+01,  4.3154e+02],
         [-9.4270e-02,  4.6078e-01, -2.4593e-02]],

        [[ 1.2042e+02,  2.7134e+01,  4.3215e+02],
         [-8.8863e-02,  4.6253e-01, -4.7904e-03]],

        [[ 1.1588e+02,  4.9875e+01,  4.3133e+02],
         [-9.5256e-02,  4.6036e-01, -2.8271e-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.7180e+01, -3.6004e+01,  4.4363e+02],
         [-2.0204e-01, -4.2680e-01, -4.5892e-03]],

        [[ 2.4651e+01, -6.5092e+01,  4.4366e+02],
         [-4.1144e-01, -2.3154e-01, -1.0460e-02]],

        [[-1.0796e+02, -3.3918e+01,  4.4222e+02],
         [ 9.6697e-02, -4.8886e-01, -1.3957e-02]],

        [[-3.7333e+01, -1.0605e+02,  4.4122e+02],
         [ 4.5411e-01, -1.9882e-01, -6.7121e-03]],

        [[-7.8860e+01, -8.1715e+01,  4.4170e+02],
         [ 3.0359e-01, -3.9254e-01, -1.6069e-02]],

        [[ 4.6672e+01, -3.7112e+02, -2.3730e+02],
         [ 3.8606e-02,  2.4403e-01, -3.7421e-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([[[-2.7184e-05,  0.0000e+00,  0.0000e+00,  9.3720e+02,  1.9122e+03,
          -6.2828e+01,  1.9270e+03, -1.4037e+04,  7.0024e+03],
         [-1.4554e-04,  0.0000e+00,  0.0000e+00,  6.3800e+03,  4.3616e+02,
          -1.5214e+01,  4.2466e+02,  7.4512e+04,  1.3427e+03],
         [-1.8805e-04,  0.0000e+00,  0.0000e+00,  7.3696e+03,  6.8569e+03,
           1.4119e+01,  6.9006e+03,  1.9040e+02,  0.0000e+00],
         [-2.9565e-07,  0.0000e+00,  0.0000e+00,  9.2625e-01, -1.4080e+00,
          -7.0737e+00, -1.4165e+00,  1.0472e+01, -4.4273e-01],
         [-1.6639e-07,  0.0000e+00,  0.0000e+00,  4.3726e+00,  7.3645e+00,
          -1.3401e+00,  7.3870e+00,  3.3568e+00,  2.0088e+00],
         [-1.1234e-06,  0.0000e+00,  0.0000e+00,  6.0347e+00, -5.2318e-02,
           2.0517e+00, -6.7107e-02,  3.8330e+01,  0.0000e+00]],

        [[-2.2675e-05,  0.0000e+00,  0.0000e+00,  9.3589e+02,  1.9141e+03,
          -5.2993e+01,  1.9290e+03, -1.4051e+04,  7.0030e+03],
         [-1.2244e-04,  0.0000e+00,  0.0000e+00,  6.3739e+03,  4.2592e+02,
          -1.3351e+01,  4.1438e+02,  7.4508e+04,  1.3400e+03],
         [-1.5754e-04,  0.0000e+00,  0.0000e+00,  7.3612e+03,  6.8570e+03,
           1.1266e+01,  6.9007e+03,  1.3710e+02,  0.0000e+00],
         [-2.4923e-07,  0.0000e+00,  0.0000e+00,  9.3289e-01, -1.4050e+00,
          -7.0738e+00, -1.4136e+00,  1.0516e+01, -4.3201e-01],
         [-1.3913e-07,  0.0000e+00,  0.0000e+00,  4.3552e+00,  7.3651e+00,
          -1.3401e+00,  7.3877e+00,  3.1286e+00,  2.0108e+00],
         [-9.4630e-07,  0.0000e+00,  0.0000e+00,  6.0462e+00, -4.1794e-02,
           2.0517e+00, -5.6517e-02,  3.8333e+01,  0.0000e+00]],

        [[-6.6819e-05,  0.0000e+00,  0.0000e+00,  9.4742e+02,  1.8961e+03,
          -1.4278e+02,  1.9109e+03, -1.3921e+04,  6.9969e+03],
         [-3.3522e-04,  0.0000e+00,  0.0000e+00,  6.4302e+03,  5.1938e+02,
          -3.0361e+01,  5.0812e+02,  7.4561e+04,  1.3654e+03],
         [-4.4754e-04,  0.0000e+00,  0.0000e+00,  7.4373e+03,  6.8558e+03,
           3.7310e+01,  6.8994e+03,  6.2346e+02,  0.0000e+00],
         [-6.6920e-07,  0.0000e+00,  0.0000e+00,  8.7196e-01, -1.4317e+00,
          -7.0724e+00, -1.4403e+00,  1.0108e+01, -5.2984e-01],
         [-4.0160e-07,  0.0000e+00,  0.0000e+00,  4.5129e+00,  7.3588e+00,
          -1.3399e+00,  7.3809e+00,  5.2104e+00,  1.9919e+00],
         [-2.5591e-06,  0.0000e+00,  0.0000e+00,  5.9390e+00, -1.3787e-01,
           2.0514e+00, -1.5320e-01,  3.8284e+01,  0.0000e+00]],

        [[-1.2760e-04,  0.0000e+00,  0.0000e+00,  9.6007e+02,  1.8744e+03,
          -2.4905e+02,  1.8890e+03, -1.3772e+04,  6.9881e+03],
         [-5.9282e-04,  0.0000e+00,  0.0000e+00,  6.4995e+03,  6.2990e+02,
          -5.0495e+01,  6.1898e+02,  7.4657e+04,  1.3951e+03],
         [-8.2377e-04,  0.0000e+00,  0.0000e+00,  7.5256e+03,  6.8529e+03,
           6.8135e+01,  6.8962e+03,  1.1980e+03,  0.0000e+00],
         [-1.1551e-06,  0.0000e+00,  0.0000e+00,  7.9903e-01, -1.4630e+00,
          -7.0691e+00, -1.4717e+00,  9.6080e+00, -6.4554e-01],
         [-7.5131e-07,  0.0000e+00,  0.0000e+00,  4.6954e+00,  7.3496e+00,
          -1.3394e+00,  7.3712e+00,  7.6697e+00,  1.9691e+00],
         [-4.4549e-06,  0.0000e+00,  0.0000e+00,  5.8066e+00, -2.5158e-01,
           2.0505e+00, -2.6760e-01,  3.8163e+01,  0.0000e+00]],

        [[-3.7251e-05,  0.0000e+00,  0.0000e+00,  9.3999e+02,  1.9079e+03,
          -8.4182e+01,  1.9228e+03, -1.4005e+04,  7.0011e+03],
         [-1.9588e-04,  0.0000e+00,  0.0000e+00,  6.3933e+03,  4.5839e+02,
          -1.9259e+01,  4.4695e+02,  7.4523e+04,  1.3488e+03],
         [-2.5537e-04,  0.0000e+00,  0.0000e+00,  7.3878e+03,  6.8567e+03,
           2.0313e+01,  6.9004e+03,  3.0610e+02,  0.0000e+00],
         [-3.9608e-07,  0.0000e+00,  0.0000e+00,  9.1180e-01, -1.4143e+00,
          -7.0735e+00, -1.4229e+00,  1.0376e+01, -4.6600e-01],
         [-2.2685e-07,  0.0000e+00,  0.0000e+00,  4.4103e+00,  7.3631e+00,
          -1.3400e+00,  7.3855e+00,  3.8521e+00,  2.0043e+00],
         [-1.5076e-06,  0.0000e+00,  0.0000e+00,  6.0094e+00, -7.5167e-02,
           2.0516e+00, -9.0100e-02,  3.8321e+01,  0.0000e+00]],

        [[-1.4198e-04,  0.0000e+00,  0.0000e+00,  9.6266e+02,  1.8696e+03,
          -2.7208e+02,  1.8842e+03, -1.3741e+04,  6.9860e+03],
         [-6.4955e-04,  0.0000e+00,  0.0000e+00,  6.5148e+03,  6.5385e+02,
          -5.4859e+01,  6.4299e+02,  7.4683e+04,  1.4015e+03],
         [-9.1003e-04,  0.0000e+00,  0.0000e+00,  7.5445e+03,  6.8520e+03,
           7.4816e+01,  6.8953e+03,  1.3223e+03,  0.0000e+00],
         [-1.2588e-06,  0.0000e+00,  0.0000e+00,  7.8311e-01, -1.4697e+00,
          -7.0682e+00, -1.4784e+00,  9.4972e+00, -6.7060e-01],
         [-8.3262e-07,  0.0000e+00,  0.0000e+00,  4.7344e+00,  7.3474e+00,
          -1.3392e+00,  7.3688e+00,  8.2021e+00,  1.9641e+00],
         [-4.8638e-06,  0.0000e+00,  0.0000e+00,  5.7771e+00, -2.7623e-01,
           2.0502e+00, -2.9240e-01,  3.8128e+01,  0.0000e+00]],

        [[-2.3639e-04,  0.0000e+00,  0.0000e+00,  9.7706e+02,  1.8406e+03,
          -4.0934e+02,  1.8551e+03, -1.3563e+04,  6.9715e+03],
         [-9.9493e-04,  0.0000e+00,  0.0000e+00,  6.6091e+03,  7.9646e+02,
          -8.0867e+01,  7.8601e+02,  7.4873e+04,  1.4394e+03],
         [-1.4587e-03,  0.0000e+00,  0.0000e+00,  7.6551e+03,  6.8452e+03,
           1.1463e+02,  6.8882e+03,  2.0606e+03,  0.0000e+00],
         [-1.8648e-06,  0.0000e+00,  0.0000e+00,  6.8746e-01, -1.5094e+00,
          -7.0609e+00, -1.5183e+00,  8.8192e+00, -8.1986e-01],
         [-1.3573e-06,  0.0000e+00,  0.0000e+00,  4.9620e+00,  7.3322e+00,
          -1.3380e+00,  7.3530e+00,  1.1368e+01,  1.9337e+00],
         [-7.2855e-06,  0.0000e+00,  0.0000e+00,  5.5959e+00, -4.2310e-01,
           2.0482e+00, -4.4015e-01,  3.7854e+01,  0.0000e+00]],

        [[-2.1347e-05,  0.0000e+00,  0.0000e+00,  9.3551e+02,  1.9147e+03,
          -5.0063e+01,  1.9296e+03, -1.4056e+04,  7.0032e+03],
         [-1.1556e-04,  0.0000e+00,  0.0000e+00,  6.3721e+03,  4.2287e+02,
          -1.2796e+01,  4.1132e+02,  7.4506e+04,  1.3391e+03],
         [-1.4850e-04,  0.0000e+00,  0.0000e+00,  7.3586e+03,  6.8570e+03,
           1.0416e+01,  6.9007e+03,  1.2123e+02,  0.0000e+00],
         [-2.3538e-07,  0.0000e+00,  0.0000e+00,  9.3487e-01, -1.4041e+00,
          -7.0738e+00, -1.4127e+00,  1.0529e+01, -4.2882e-01],
         [-1.3108e-07,  0.0000e+00,  0.0000e+00,  4.3500e+00,  7.3653e+00,
          -1.3401e+00,  7.3879e+00,  3.0606e+00,  2.0114e+00],
         [-8.9351e-07,  0.0000e+00,  0.0000e+00,  6.0496e+00, -3.8659e-02,
           2.0517e+00, -5.3361e-02,  3.8334e+01,  0.0000e+00]],

        [[-7.5956e-05,  0.0000e+00,  0.0000e+00,  9.4953e+02,  1.8926e+03,
          -1.5984e+02,  1.9074e+03, -1.3896e+04,  6.9956e+03],
         [-3.7613e-04,  0.0000e+00,  0.0000e+00,  6.4412e+03,  5.3712e+02,
          -3.3593e+01,  5.2592e+02,  7.4574e+04,  1.3702e+03],
         [-5.0552e-04,  0.0000e+00,  0.0000e+00,  7.4516e+03,  6.8555e+03,
           4.2258e+01,  6.8990e+03,  7.1577e+02,  0.0000e+00],
         [-7.4801e-07,  0.0000e+00,  0.0000e+00,  8.6031e-01, -1.4367e+00,
          -7.0720e+00, -1.4454e+00,  1.0029e+01, -5.4842e-01],
         [-4.5490e-07,  0.0000e+00,  0.0000e+00,  4.5425e+00,  7.3574e+00,
          -1.3398e+00,  7.3795e+00,  5.6055e+00,  1.9883e+00],
         [-2.8643e-06,  0.0000e+00,  0.0000e+00,  5.9181e+00, -1.5612e-01,
           2.0513e+00, -1.7156e-01,  3.8269e+01,  0.0000e+00]],

        [[-8.5479e-05,  0.0000e+00,  0.0000e+00,  9.5163e+02,  1.8891e+03,
          -1.7716e+02,  1.9039e+03, -1.3872e+04,  6.9943e+03],
         [-4.1784e-04,  0.0000e+00,  0.0000e+00,  6.4523e+03,  5.5514e+02,
          -3.6874e+01,  5.4400e+02,  7.4588e+04,  1.3750e+03],
         [-5.6535e-04,  0.0000e+00,  0.0000e+00,  7.4661e+03,  6.8551e+03,
           4.7282e+01,  6.8985e+03,  8.0949e+02,  0.0000e+00],
         [-8.2773e-07,  0.0000e+00,  0.0000e+00,  8.4846e-01, -1.4419e+00,
          -7.0715e+00, -1.4505e+00,  9.9482e+00, -5.6728e-01],
         [-5.1013e-07,  0.0000e+00,  0.0000e+00,  4.5724e+00,  7.3560e+00,
          -1.3398e+00,  7.3780e+00,  6.0066e+00,  1.9846e+00],
         [-3.1740e-06,  0.0000e+00,  0.0000e+00,  5.8968e+00, -1.7466e-01,
           2.0511e+00, -1.9021e-01,  3.8252e+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.5642e-05,  0.0000e+00,  0.0000e+00, -1.3664e+03,  9.4478e+02,
           3.3895e+01,  9.2021e+02, -3.2021e+04, -6.4385e+03],
         [-2.2545e-04,  0.0000e+00,  0.0000e+00,  2.0294e+03, -4.9629e+02,
          -1.8026e+01, -5.4239e+02, -6.7916e+04,  3.0428e+03],
         [ 3.6253e-04,  0.0000e+00,  0.0000e+00, -1.3881e+04,  7.0426e+03,
          -3.5086e+00,  7.0257e+03,  4.5703e+01,  0.0000e+00],
         [-1.5704e-07,  0.0000e+00,  0.0000e+00,  3.0031e+00, -3.2061e+00,
           6.6833e+00, -3.1952e+00,  4.8505e+00,  5.7024e-01],
         [-1.3699e-07,  0.0000e+00,  0.0000e+00,  6.7659e+00, -6.7649e+00,
          -3.1531e+00, -6.7610e+00, -3.5072e+00,  9.6711e-01],
         [-5.5091e-07,  0.0000e+00,  0.0000e+00, -1.1228e+00, -9.4534e-02,
          -1.1154e+00, -4.0051e-02,  3.9013e+01,  0.0000e+00]],

        [[-8.0306e-06,  0.0000e+00,  0.0000e+00, -5.0992e+03,  5.2614e+02,
           2.2033e+00,  5.2422e+02, -6.5461e+04, -3.5128e+03],
         [-9.1220e-06,  0.0000e+00,  0.0000e+00, -1.1010e+03, -9.5769e+02,
          -7.6949e+00, -9.5969e+02, -3.7114e+04,  6.1947e+03],
         [ 2.5588e-05,  0.0000e+00,  0.0000e+00, -1.0039e+04,  7.0386e+03,
           1.1141e+00,  7.0437e+03, -1.1550e+02,  0.0000e+00],
         [-2.1186e-08,  0.0000e+00,  0.0000e+00,  5.0342e+00, -6.5178e+00,
           3.6462e+00, -6.5199e+00,  2.8621e+00,  1.0064e+00],
         [ 1.0290e-08,  0.0000e+00,  0.0000e+00,  1.9299e+00, -3.6955e+00,
          -6.4300e+00, -3.6972e+00, -5.3299e+00,  5.4971e-01],
         [-1.2481e-07,  0.0000e+00,  0.0000e+00,  5.1949e+00, -1.0590e-02,
          -1.1439e+00, -7.8453e-03,  3.9115e+01,  0.0000e+00]],

        [[-1.0210e-03,  0.0000e+00,  0.0000e+00,  6.8213e+02, -1.6197e+03,
           2.9051e+02, -1.6404e+03,  1.2950e+04, -6.7880e+03],
         [ 5.6088e-03,  0.0000e+00,  0.0000e+00, -6.9387e+03, -6.2658e+02,
           5.8482e+01, -6.2559e+02, -6.9854e+04, -1.3691e+03],
         [-1.8029e-04,  0.0000e+00,  0.0000e+00,  2.4669e+03,  6.7121e+03,
           7.1367e+01,  6.7913e+03,  1.3672e+03,  0.0000e+00],
         [ 5.3471e-06,  0.0000e+00,  0.0000e+00, -1.4288e+00,  1.5024e+00,
           7.2321e+00,  1.5134e+00, -8.0795e+00,  6.7804e-01],
         [ 1.7537e-06,  0.0000e+00,  0.0000e+00, -2.2066e+00, -7.4597e+00,
           1.3796e+00, -7.5028e+00, -8.4025e+00, -1.7767e+00],
         [-2.1914e-05,  0.0000e+00,  0.0000e+00,  7.1464e+00, -3.1671e-01,
           1.8681e+00, -3.2754e-01,  3.7514e+01,  0.0000e+00]],

        [[-4.1331e-03,  0.0000e+00,  0.0000e+00,  4.8139e+03, -2.9842e+02,
           1.6269e+02, -2.7500e+02,  6.5136e+04, -2.7085e+03],
         [-7.1212e-04,  0.0000e+00,  0.0000e+00,  6.8869e+02, -1.7380e+03,
           3.7592e+02, -1.7620e+03, -2.9324e+04, -6.3858e+03],
         [ 9.1307e-03,  0.0000e+00,  0.0000e+00, -1.0126e+04,  6.7235e+03,
           1.0007e+02,  6.7726e+03,  1.9440e+03,  0.0000e+00],
         [ 6.6237e-06,  0.0000e+00,  0.0000e+00, -5.3289e+00,  6.9847e+00,
           2.9636e+00,  7.0099e+00,  4.9130e+00,  1.9075e+00],
         [ 4.3232e-06,  0.0000e+00,  0.0000e+00,  8.6657e-01, -2.9542e+00,
           6.7120e+00, -2.9730e+00, -1.2479e+01, -2.9706e-01],
         [-2.6048e-05,  0.0000e+00,  0.0000e+00,  5.3598e+00, -4.7766e-01,
           1.8608e+00, -4.5003e-01,  3.7214e+01,  0.0000e+00]],

        [[-3.2274e-03,  0.0000e+00,  0.0000e+00,  4.9141e+03, -1.2074e+03,
           1.9085e+02, -1.2069e+03,  4.3002e+04, -5.4877e+03],
         [ 2.4299e-03,  0.0000e+00,  0.0000e+00, -3.5001e+03, -1.2403e+03,
           1.4899e+02, -1.2755e+03, -5.7011e+04, -4.2498e+03],
         [ 4.3177e-03,  0.0000e+00,  0.0000e+00, -7.0483e+03,  6.7279e+03,
           5.7591e+01,  6.8103e+03,  1.0943e+03,  0.0000e+00],
         [ 4.4628e-06,  0.0000e+00,  0.0000e+00, -3.4989e+00,  4.6611e+00,
           5.8553e+00,  4.6853e+00, -4.1652e+00,  1.3771e+00],
         [ 8.3629e-07,  0.0000e+00,  0.0000e+00,  1.8756e+00, -6.0077e+00,
           4.4551e+00, -6.0499e+00, -1.0302e+01, -1.3025e+00],
         [-1.5942e-05,  0.0000e+00,  0.0000e+00,  6.4941e+00, -2.9147e-01,
           1.8662e+00, -2.6673e-01,  3.7692e+01,  0.0000e+00]],

        [[-1.7827e-06,  0.0000e+00,  0.0000e+00, -6.0292e+02,  7.6122e+02,
           6.1426e+03,  7.5833e+02,  7.1336e+03,  3.8856e+03],
         [ 5.4816e-06,  0.0000e+00,  0.0000e+00, -4.1835e+03, -6.1837e+03,
           2.5800e+02, -6.1607e+03,  4.1623e+04, -6.4098e+02],
         [ 1.0834e-05,  0.0000e+00,  0.0000e+00,  5.9495e+03, -3.8233e+03,
           8.0417e+02, -3.8087e+03, -6.7945e+04,  0.0000e+00],
         [-5.1654e-08,  0.0000e+00,  0.0000e+00, -7.5270e-01,  6.4857e-01,
          -3.8556e+00,  6.4732e-01,  4.2900e+00,  6.2455e+00],
         [ 4.1567e-07,  0.0000e+00,  0.0000e+00,  6.3441e+00,  3.9251e+00,
          -1.6105e-01,  3.9181e+00, -3.3627e+01,  7.6815e-01],
         [ 2.6039e-07,  0.0000e+00,  0.0000e+00,  3.7111e+00, -6.2227e+00,
          -5.0294e-01, -6.2109e+00, -2.2023e+01,  0.0000e+00]]])