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.3306e+03, -7.0050e+03, 2.0089e+01],
[ 2.0176e+00, 3.9631e-01, 7.2028e+00]],
[[ 1.3913e+03, -6.9894e+03, 2.3896e+02],
[ 1.9721e+00, 6.3046e-01, 7.1985e+00]],
[[ 1.3671e+03, -6.9965e+03, 1.5114e+02],
[ 1.9906e+00, 5.3656e-01, 7.2011e+00]],
[[ 1.3910e+03, -6.9894e+03, 2.3807e+02],
[ 1.9723e+00, 6.2951e-01, 7.1985e+00]],
[[ 1.3843e+03, -6.9915e+03, 2.1368e+02],
[ 1.9774e+00, 6.0343e-01, 7.1994e+00]],
[[ 1.3945e+03, -6.9883e+03, 2.5096e+02],
[ 1.9695e+00, 6.4328e-01, 7.1980e+00]],
[[ 1.4288e+03, -6.9759e+03, 3.7678e+02],
[ 1.9423e+00, 7.7770e-01, 7.1920e+00]],
[[ 1.3292e+03, -7.0052e+03, 1.5159e+01],
[ 2.0186e+00, 3.9103e-01, 7.2028e+00]],
[[ 1.3559e+03, -6.9994e+03, 1.1088e+02],
[ 1.9990e+00, 4.9348e-01, 7.2019e+00]],
[[ 1.3639e+03, -6.9973e+03, 1.3969e+02],
[ 1.9930e+00, 5.2431e-01, 7.2014e+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:
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.2106e+02, 2.3778e+01, 4.3217e+02],
[-8.7901e-02, 4.6275e-01, -1.3308e-03]],
[[ 1.1833e+02, 3.7828e+01, 4.3191e+02],
[-9.1897e-02, 4.6167e-01, -1.5825e-02]],
[[ 1.1944e+02, 3.2194e+01, 4.3207e+02],
[-9.0305e-02, 4.6216e-01, -1.0010e-02]],
[[ 1.1834e+02, 3.7771e+01, 4.3191e+02],
[-9.1881e-02, 4.6167e-01, -1.5766e-02]],
[[ 1.1865e+02, 3.6206e+01, 4.3196e+02],
[-9.1440e-02, 4.6182e-01, -1.4151e-02]],
[[ 1.1817e+02, 3.8597e+01, 4.3188e+02],
[-9.2113e-02, 4.6159e-01, -1.6619e-02]],
[[ 1.1654e+02, 4.6662e+01, 4.3152e+02],
[-9.4366e-02, 4.6074e-01, -2.4950e-02]],
[[ 1.2112e+02, 2.3462e+01, 4.3217e+02],
[-8.7810e-02, 4.6277e-01, -1.0043e-03]],
[[ 1.1994e+02, 2.9609e+01, 4.3212e+02],
[-8.9570e-02, 4.6236e-01, -7.3435e-03]],
[[ 1.1958e+02, 3.1459e+01, 4.3208e+02],
[-9.0096e-02, 4.6222e-01, -9.2519e-03]]])
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.3076e+01, -4.4639e+01, 4.4344e+02],
[-2.0339e-01, -4.2613e-01, -1.4548e-02]],
[[ 9.9186e+00, -7.3327e+01, 4.4297e+02],
[-4.1218e-01, -2.2883e-01, -2.8037e-02]],
[[-1.0729e+02, -3.7270e+01, 4.4211e+02],
[ 9.7582e-02, -4.8857e-01, -1.7603e-02]],
[[-2.5585e+01, -1.1114e+02, 4.4087e+02],
[ 4.5527e-01, -1.9557e-01, -2.0306e-02]],
[[-7.7650e+01, -8.3277e+01, 4.4164e+02],
[ 3.0398e-01, -3.9217e-01, -1.8171e-02]],
[[ 4.7077e+01, -3.6852e+02, -2.4123e+02],
[ 3.8101e-02, 2.4803e-01, -3.7163e-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_kozaiin 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\)]
#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.0287e-04, 0.0000e+00, 0.0000e+00, 9.7241e+02, 1.8505e+03,
-3.6303e+02, 1.8650e+03, -1.3622e+04, 6.9767e+03],
[-8.7694e-04, 0.0000e+00, 0.0000e+00, 6.5768e+03, 7.4835e+02,
-7.2091e+01, 7.3777e+02, 7.4802e+04, 1.4267e+03],
[-1.2669e-03, 0.0000e+00, 0.0000e+00, 7.6182e+03, 6.8478e+03,
1.0120e+02, 6.8909e+03, 1.8120e+03, 0.0000e+00],
[-1.6626e-06, 0.0000e+00, 0.0000e+00, 7.1988e-01, -1.4961e+00,
-7.0637e+00, -1.5049e+00, 9.0514e+00, -7.6951e-01],
[-1.1726e-06, 0.0000e+00, 0.0000e+00, 4.8861e+00, 7.3377e+00,
-1.3385e+00, 7.3587e+00, 1.0301e+01, 1.9440e+00],
[-6.4714e-06, 0.0000e+00, 0.0000e+00, 5.6581e+00, -3.7354e-01,
2.0490e+00, -3.9030e-01, 3.7959e+01, 0.0000e+00]],
[[-2.2946e-04, 0.0000e+00, 0.0000e+00, 9.7613e+02, 1.8426e+03,
-3.9996e+02, 1.8571e+03, -1.3575e+04, 6.9726e+03],
[-9.7091e-04, 0.0000e+00, 0.0000e+00, 6.6025e+03, 7.8671e+02,
-7.9089e+01, 7.7624e+02, 7.4858e+04, 1.4368e+03],
[-1.4193e-03, 0.0000e+00, 0.0000e+00, 7.6477e+03, 6.8458e+03,
1.1191e+02, 6.8888e+03, 2.0103e+03, 0.0000e+00],
[-1.8240e-06, 0.0000e+00, 0.0000e+00, 6.9404e-01, -1.5067e+00,
-7.0615e+00, -1.5156e+00, 8.8665e+00, -8.0966e-01],
[-1.3193e-06, 0.0000e+00, 0.0000e+00, 4.9467e+00, 7.3334e+00,
-1.3381e+00, 7.3542e+00, 1.1152e+01, 1.9358e+00],
[-7.1209e-06, 0.0000e+00, 0.0000e+00, 5.6086e+00, -4.1306e-01,
2.0484e+00, -4.3005e-01, 3.7876e+01, 0.0000e+00]],
[[-1.0691e-04, 0.0000e+00, 0.0000e+00, 9.5609e+02, 1.8815e+03,
-2.1462e+02, 1.8962e+03, -1.3820e+04, 6.9912e+03],
[-5.0864e-04, 0.0000e+00, 0.0000e+00, 6.4767e+03, 5.9411e+02,
-4.3973e+01, 5.8308e+02, 7.4622e+04, 1.3855e+03],
[-6.9799e-04, 0.0000e+00, 0.0000e+00, 7.4972e+03, 6.8540e+03,
5.8149e+01, 6.8974e+03, 1.0120e+03, 0.0000e+00],
[-9.9904e-07, 0.0000e+00, 0.0000e+00, 8.2275e-01, -1.4529e+00,
-7.0704e+00, -1.4616e+00, 9.7720e+00, -6.0807e-01],
[-6.3343e-07, 0.0000e+00, 0.0000e+00, 4.6368e+00, 7.3528e+00,
-1.3396e+00, 7.3746e+00, 6.8737e+00, 1.9765e+00],
[-3.8423e-06, 0.0000e+00, 0.0000e+00, 5.8501e+00, -2.1474e-01,
2.0508e+00, -2.3054e-01, 3.8209e+01, 0.0000e+00]],
[[-8.2863e-05, 0.0000e+00, 0.0000e+00, 9.5106e+02, 1.8901e+03,
-1.7245e+02, 1.9048e+03, -1.3878e+04, 6.9946e+03],
[-4.0647e-04, 0.0000e+00, 0.0000e+00, 6.4493e+03, 5.5024e+02,
-3.5981e+01, 5.3908e+02, 7.4584e+04, 1.3737e+03],
[-5.4898e-04, 0.0000e+00, 0.0000e+00, 7.4622e+03, 6.8552e+03,
4.5915e+01, 6.8987e+03, 7.8399e+02, 0.0000e+00],
[-8.0607e-07, 0.0000e+00, 0.0000e+00, 8.5169e-01, -1.4405e+00,
-7.0717e+00, -1.4491e+00, 9.9703e+00, -5.6215e-01],
[-4.9499e-07, 0.0000e+00, 0.0000e+00, 4.5643e+00, 7.3564e+00,
-1.3398e+00, 7.3784e+00, 5.8975e+00, 1.9856e+00],
[-3.0897e-06, 0.0000e+00, 0.0000e+00, 5.9026e+00, -1.6961e-01,
2.0512e+00, -1.8513e-01, 3.8257e+01, 0.0000e+00]],
[[-8.2963e-05, 0.0000e+00, 0.0000e+00, 9.5109e+02, 1.8900e+03,
-1.7263e+02, 1.9048e+03, -1.3878e+04, 6.9946e+03],
[-4.0691e-04, 0.0000e+00, 0.0000e+00, 6.4494e+03, 5.5043e+02,
-3.6016e+01, 5.3927e+02, 7.4584e+04, 1.3738e+03],
[-5.4960e-04, 0.0000e+00, 0.0000e+00, 7.4623e+03, 6.8552e+03,
4.5967e+01, 6.8987e+03, 7.8496e+02, 0.0000e+00],
[-8.0690e-07, 0.0000e+00, 0.0000e+00, 8.5156e-01, -1.4405e+00,
-7.0717e+00, -1.4492e+00, 9.9694e+00, -5.6235e-01],
[-4.9557e-07, 0.0000e+00, 0.0000e+00, 4.5646e+00, 7.3564e+00,
-1.3398e+00, 7.3784e+00, 5.9017e+00, 1.9855e+00],
[-3.0930e-06, 0.0000e+00, 0.0000e+00, 5.9024e+00, -1.6981e-01,
2.0512e+00, -1.8533e-01, 3.8256e+01, 0.0000e+00]],
[[-7.3898e-05, 0.0000e+00, 0.0000e+00, 9.4906e+02, 1.8934e+03,
-1.5603e+02, 1.9082e+03, -1.3902e+04, 6.9959e+03],
[-3.6699e-04, 0.0000e+00, 0.0000e+00, 6.4387e+03, 5.3317e+02,
-3.2872e+01, 5.2195e+02, 7.4571e+04, 1.3691e+03],
[-4.9251e-04, 0.0000e+00, 0.0000e+00, 7.4484e+03, 6.8555e+03,
4.1154e+01, 6.8991e+03, 6.9518e+02, 0.0000e+00],
[-7.3046e-07, 0.0000e+00, 0.0000e+00, 8.6291e-01, -1.4356e+00,
-7.0721e+00, -1.4443e+00, 1.0047e+01, -5.4428e-01],
[-4.4292e-07, 0.0000e+00, 0.0000e+00, 4.5359e+00, 7.3578e+00,
-1.3399e+00, 7.3798e+00, 5.5174e+00, 1.9891e+00],
[-2.7963e-06, 0.0000e+00, 0.0000e+00, 5.9228e+00, -1.5205e-01,
2.0513e+00, -1.6746e-01, 3.8272e+01, 0.0000e+00]],
[[-1.0982e-05, 0.0000e+00, 0.0000e+00, 9.3236e+02, 1.9193e+03,
-2.6615e+01, 1.9343e+03, -1.4091e+04, 7.0046e+03],
[-6.0695e-05, 0.0000e+00, 0.0000e+00, 6.3578e+03, 3.9845e+02,
-8.3535e+00, 3.8683e+02, 7.4497e+04, 1.3324e+03],
[-7.7208e-05, 0.0000e+00, 0.0000e+00, 7.3385e+03, 6.8571e+03,
3.6148e+00, 6.9009e+03, -5.8517e+00, 0.0000e+00],
[-1.2424e-07, 0.0000e+00, 0.0000e+00, 9.5066e-01, -1.3971e+00,
-7.0740e+00, -1.4057e+00, 1.0633e+01, -4.0326e-01],
[-6.7841e-08, 0.0000e+00, 0.0000e+00, 4.3083e+00, 7.3667e+00,
-1.3401e+00, 7.3894e+00, 2.5164e+00, 2.0163e+00],
[-4.7073e-07, 0.0000e+00, 0.0000e+00, 6.0769e+00, -1.3569e-02,
2.0517e+00, -2.8113e-02, 3.8339e+01, 0.0000e+00]],
[[-2.1373e-04, 0.0000e+00, 0.0000e+00, 9.7396e+02, 1.8473e+03,
-3.7828e+02, 1.8618e+03, -1.3602e+04, 6.9750e+03],
[-9.1564e-04, 0.0000e+00, 0.0000e+00, 6.5874e+03, 7.6420e+02,
-7.4982e+01, 7.5366e+02, 7.4825e+04, 1.4309e+03],
[-1.3293e-03, 0.0000e+00, 0.0000e+00, 7.6304e+03, 6.8470e+03,
1.0562e+02, 6.8900e+03, 1.8939e+03, 0.0000e+00],
[-1.7295e-06, 0.0000e+00, 0.0000e+00, 7.0922e-01, -1.5005e+00,
-7.0628e+00, -1.5093e+00, 8.9753e+00, -7.8610e-01],
[-1.2326e-06, 0.0000e+00, 0.0000e+00, 4.9112e+00, 7.3359e+00,
-1.3383e+00, 7.3568e+00, 1.0653e+01, 1.9406e+00],
[-6.7399e-06, 0.0000e+00, 0.0000e+00, 5.6377e+00, -3.8987e-01,
2.0487e+00, -4.0672e-01, 3.7926e+01, 0.0000e+00]],
[[-2.5914e-05, 0.0000e+00, 0.0000e+00, 9.3683e+02, 1.9127e+03,
-6.0075e+01, 1.9276e+03, -1.4041e+04, 7.0026e+03],
[-1.3907e-04, 0.0000e+00, 0.0000e+00, 6.3783e+03, 4.3329e+02,
-1.4692e+01, 4.2178e+02, 7.4511e+04, 1.3420e+03],
[-1.7948e-04, 0.0000e+00, 0.0000e+00, 7.3672e+03, 6.8569e+03,
1.3320e+01, 6.9006e+03, 1.7548e+02, 0.0000e+00],
[-2.8266e-07, 0.0000e+00, 0.0000e+00, 9.2811e-01, -1.4071e+00,
-7.0737e+00, -1.4157e+00, 1.0484e+01, -4.3973e-01],
[-1.5872e-07, 0.0000e+00, 0.0000e+00, 4.3677e+00, 7.3647e+00,
-1.3401e+00, 7.3872e+00, 3.2929e+00, 2.0093e+00],
[-1.0738e-06, 0.0000e+00, 0.0000e+00, 6.0379e+00, -4.9372e-02,
2.0517e+00, -6.4142e-02, 3.8331e+01, 0.0000e+00]],
[[-6.8646e-05, 0.0000e+00, 0.0000e+00, 9.4785e+02, 1.8954e+03,
-1.4623e+02, 1.9102e+03, -1.3916e+04, 6.9967e+03],
[-3.4348e-04, 0.0000e+00, 0.0000e+00, 6.4324e+03, 5.2296e+02,
-3.1014e+01, 5.1172e+02, 7.4563e+04, 1.3663e+03],
[-4.5919e-04, 0.0000e+00, 0.0000e+00, 7.4402e+03, 6.8557e+03,
3.8310e+01, 6.8993e+03, 6.4211e+02, 0.0000e+00],
[-6.8516e-07, 0.0000e+00, 0.0000e+00, 8.6961e-01, -1.4327e+00,
-7.0723e+00, -1.4414e+00, 1.0092e+01, -5.3360e-01],
[-4.1229e-07, 0.0000e+00, 0.0000e+00, 4.5189e+00, 7.3585e+00,
-1.3399e+00, 7.3807e+00, 5.2903e+00, 1.9912e+00],
[-2.6208e-06, 0.0000e+00, 0.0000e+00, 5.9348e+00, -1.4156e-01,
2.0513e+00, -1.5691e-01, 3.8281e+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([[[-2.5379e-04, 0.0000e+00, 0.0000e+00, -1.2976e+03, 8.7194e+02,
1.8516e+02, 8.4761e+02, -3.1930e+04, -6.4238e+03],
[-1.2303e-03, 0.0000e+00, 0.0000e+00, 2.1829e+03, -6.4927e+02,
-8.9391e+01, -6.9532e+02, -6.8033e+04, 3.0638e+03],
[ 1.9502e-03, 0.0000e+00, 0.0000e+00, -1.3902e+04, 7.0385e+03,
-2.8753e+01, 7.0228e+03, 9.2787e+02, 0.0000e+00],
[-8.6312e-07, 0.0000e+00, 0.0000e+00, 3.0716e+00, -3.2291e+00,
6.6805e+00, -3.2192e+00, 3.2313e+00, 7.3117e-01],
[-7.4893e-07, 0.0000e+00, 0.0000e+00, 6.7873e+00, -6.7510e+00,
-3.1519e+00, -6.7496e+00, -6.9021e+00, 8.9071e-01],
[-3.0395e-06, 0.0000e+00, 0.0000e+00, -7.7286e-01, -2.7119e-01,
-1.1150e+00, -2.1636e-01, 3.8904e+01, 0.0000e+00]],
[[-1.0017e-04, 0.0000e+00, 0.0000e+00, -5.0428e+03, 4.5239e+02,
4.3437e+01, 4.5045e+02, -6.5438e+04, -3.5011e+03],
[-1.1142e-04, 0.0000e+00, 0.0000e+00, -1.0796e+03, -9.9941e+02,
-8.0411e+01, -1.0014e+03, -3.7179e+04, 6.2005e+03],
[ 3.0596e-04, 0.0000e+00, 0.0000e+00, -9.9795e+03, 7.0380e+03,
-1.1823e+01, 7.0432e+03, 3.2678e+02, 0.0000e+00],
[-2.5488e-07, 0.0000e+00, 0.0000e+00, 4.9337e+00, -6.5241e+00,
3.6459e+00, -6.5263e+00, 1.2223e+00, 1.0502e+00],
[ 1.3071e-07, 0.0000e+00, 0.0000e+00, 1.8508e+00, -3.6834e+00,
-6.4295e+00, -3.6851e+00, -6.2563e+00, 4.7225e-01],
[-1.5402e-06, 0.0000e+00, 0.0000e+00, 5.3191e+00, -9.8789e-02,
-1.1439e+00, -9.6110e-02, 3.9094e+01, 0.0000e+00]],
[[-9.7121e-04, 0.0000e+00, 0.0000e+00, 6.8529e+02, -1.6230e+03,
2.7465e+02, -1.6437e+03, 1.2968e+04, -6.7895e+03],
[ 5.2997e-03, 0.0000e+00, 0.0000e+00, -6.9339e+03, -6.1022e+02,
5.5457e+01, -6.0914e+02, -6.9836e+04, -1.3652e+03],
[-1.4533e-04, 0.0000e+00, 0.0000e+00, 2.4512e+03, 6.7128e+03,
6.7270e+01, 6.7920e+03, 1.2849e+03, 0.0000e+00],
[ 5.0635e-06, 0.0000e+00, 0.0000e+00, -1.4376e+00, 1.4981e+00,
7.2329e+00, 1.5091e+00, -8.1558e+00, 6.6020e-01],
[ 1.6142e-06, 0.0000e+00, 0.0000e+00, -2.1720e+00, -7.4613e+00,
1.3797e+00, -7.5044e+00, -8.0379e+00, -1.7803e+00],
[-2.0718e-05, 0.0000e+00, 0.0000e+00, 7.1552e+00, -2.9902e-01,
1.8683e+00, -3.0964e-01, 3.7540e+01, 0.0000e+00]],
[[-2.3451e-04, 0.0000e+00, 0.0000e+00, 5.1076e+03, -6.6545e+02,
6.4626e+00, -6.4326e+02, 6.5094e+04, -2.8044e+03],
[-4.4865e-05, 0.0000e+00, 0.0000e+00, 6.3311e+02, -1.5797e+03,
2.2104e+01, -1.6028e+03, -2.8758e+04, -6.3597e+03],
[ 5.3396e-04, 0.0000e+00, 0.0000e+00, -1.0392e+04, 6.7375e+03,
1.9762e+00, 6.7850e+03, -3.1759e+01, 0.0000e+00],
[ 4.0810e-07, 0.0000e+00, 0.0000e+00, -5.8119e+00, 6.9534e+00,
2.9694e+00, 6.9751e+00, -3.2854e+00, 1.7349e+00],
[ 1.9784e-07, 0.0000e+00, 0.0000e+00, 1.2512e+00, -3.0576e+00,
6.7242e+00, -3.0758e+00, -8.9879e+00, -6.9618e-01],
[-1.3585e-06, 0.0000e+00, 0.0000e+00, 4.7426e+00, -5.4954e-02,
1.8640e+00, -2.4120e-02, 3.7691e+01, 0.0000e+00]],
[[-4.7561e-03, 0.0000e+00, 0.0000e+00, 4.8581e+03, -1.1313e+03,
2.8618e+02, -1.1304e+03, 4.2948e+04, -5.4644e+03],
[ 3.6408e-03, 0.0000e+00, 0.0000e+00, -3.4713e+03, -1.3379e+03,
2.2152e+02, -1.3738e+03, -5.7197e+04, -4.2703e+03],
[ 6.2207e-03, 0.0000e+00, 0.0000e+00, -6.9414e+03, 6.7221e+03,
8.7978e+01, 6.8049e+03, 1.7067e+03, 0.0000e+00],
[ 6.4639e-06, 0.0000e+00, 0.0000e+00, -3.3775e+00, 4.6840e+00,
5.8505e+00, 4.7091e+00, -2.4494e+00, 1.4834e+00],
[ 1.5845e-06, 0.0000e+00, 0.0000e+00, 1.6653e+00, -5.9830e+00,
4.4516e+00, -6.0257e+00, -1.2472e+01, -1.2198e+00],
[-2.3979e-05, 0.0000e+00, 0.0000e+00, 6.6149e+00, -4.2237e-01,
1.8648e+00, -3.9929e-01, 3.7484e+01, 0.0000e+00]],
[[-3.6326e-06, 0.0000e+00, 0.0000e+00, -6.0807e+02, 7.6560e+02,
6.1163e+03, 7.6270e+02, 7.1630e+03, 3.9279e+03],
[ 1.3495e-05, 0.0000e+00, 0.0000e+00, -4.1406e+03, -6.1569e+03,
2.5690e+02, -6.1339e+03, 4.1397e+04, -6.3575e+02],
[ 2.1597e-05, 0.0000e+00, 0.0000e+00, 5.9750e+03, -3.8654e+03,
8.0074e+02, -3.8507e+03, -6.8098e+04, 0.0000e+00],
[-9.1641e-08, 0.0000e+00, 0.0000e+00, -7.6176e-01, 6.4329e-01,
-3.8981e+00, 6.4206e-01, 4.3857e+00, 6.2184e+00],
[ 7.3260e-07, 0.0000e+00, 0.0000e+00, 6.2917e+00, 3.9680e+00,
-1.6286e-01, 3.9608e+00, -3.3025e+01, 7.7257e-01],
[ 4.6351e-07, 0.0000e+00, 0.0000e+00, 3.7974e+00, -6.1961e+00,
-5.0846e-01, -6.1843e+00, -2.2952e+01, 0.0000e+00]]])