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.4432e+03, -6.9699e+03, 4.3030e+02],
[ 1.9306e+00, 8.3482e-01, 7.1887e+00]],
[[ 1.3291e+03, -7.0052e+03, 1.4813e+01],
[ 2.0187e+00, 3.9066e-01, 7.2028e+00]],
[[ 1.3573e+03, -6.9990e+03, 1.1562e+02],
[ 1.9980e+00, 4.9855e-01, 7.2018e+00]],
[[ 1.4416e+03, -6.9705e+03, 4.2442e+02],
[ 1.9319e+00, 8.2855e-01, 7.1890e+00]],
[[ 1.4430e+03, -6.9700e+03, 4.2957e+02],
[ 1.9307e+00, 8.3404e-01, 7.1887e+00]],
[[ 1.4384e+03, -6.9719e+03, 4.1254e+02],
[ 1.9345e+00, 8.1587e-01, 7.1898e+00]],
[[ 1.3962e+03, -6.9878e+03, 2.5682e+02],
[ 1.9683e+00, 6.4954e-01, 7.1978e+00]],
[[ 1.3908e+03, -6.9895e+03, 2.3729e+02],
[ 1.9724e+00, 6.2867e-01, 7.1985e+00]],
[[ 1.3337e+03, -7.0044e+03, 3.0909e+01],
[ 2.0154e+00, 4.0789e-01, 7.2028e+00]],
[[ 1.3367e+03, -7.0037e+03, 4.1839e+01],
[ 2.0132e+00, 4.1959e-01, 7.2027e+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.1584e+02, 5.0089e+01, 4.3132e+02],
[-9.5315e-02, 4.6033e-01, -2.8493e-02]],
[[ 1.2112e+02, 2.3439e+01, 4.3217e+02],
[-8.7804e-02, 4.6277e-01, -9.8138e-04]],
[[ 1.1988e+02, 2.9913e+01, 4.3211e+02],
[-8.9656e-02, 4.6233e-01, -7.6573e-03]],
[[ 1.1591e+02, 4.9713e+01, 4.3134e+02],
[-9.5211e-02, 4.6038e-01, -2.8104e-02]],
[[ 1.1585e+02, 5.0043e+01, 4.3132e+02],
[-9.5302e-02, 4.6034e-01, -2.8445e-02]],
[[ 1.1607e+02, 4.8953e+01, 4.3139e+02],
[-9.5001e-02, 4.6047e-01, -2.7318e-02]],
[[ 1.1810e+02, 3.8973e+01, 4.3187e+02],
[-9.2219e-02, 4.6156e-01, -1.7007e-02]],
[[ 1.1835e+02, 3.7721e+01, 4.3191e+02],
[-9.1867e-02, 4.6168e-01, -1.5715e-02]],
[[ 1.2093e+02, 2.4474e+01, 4.3217e+02],
[-8.8101e-02, 4.6270e-01, -2.0474e-03]],
[[ 1.2079e+02, 2.5176e+01, 4.3216e+02],
[-8.8302e-02, 4.6266e-01, -2.7713e-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.2179e+01, -4.6516e+01, 4.4337e+02],
[-2.0367e-01, -4.2595e-01, -1.6715e-02]],
[[ 1.0497e+01, -7.3005e+01, 4.4301e+02],
[-4.1216e-01, -2.2895e-01, -2.7348e-02]],
[[-1.0487e+02, -4.9166e+01, 4.4152e+02],
[ 1.0068e-01, -4.8729e-01, -3.0551e-02]],
[[-3.3037e+01, -1.0792e+02, 4.4113e+02],
[ 4.5458e-01, -1.9765e-01, -1.1688e-02]],
[[-8.6371e+01, -7.1936e+01, 4.4194e+02],
[ 3.0103e-01, -3.9466e-01, -2.9538e-03]],
[[ 4.7268e+01, -3.6727e+02, -2.4310e+02],
[ 3.7859e-02, 2.4993e-01, -3.7039e-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.2839e-04, 0.0000e+00, 0.0000e+00, 9.7599e+02, 1.8429e+03,
-3.9849e+02, 1.8574e+03, -1.3577e+04, 6.9727e+03],
[-9.6716e-04, 0.0000e+00, 0.0000e+00, 6.6015e+03, 7.8519e+02,
-7.8812e+01, 7.7472e+02, 7.4856e+04, 1.4364e+03],
[-1.4131e-03, 0.0000e+00, 0.0000e+00, 7.6465e+03, 6.8459e+03,
1.1149e+02, 6.8888e+03, 2.0024e+03, 0.0000e+00],
[-1.8177e-06, 0.0000e+00, 0.0000e+00, 6.9507e-01, -1.5063e+00,
-7.0616e+00, -1.5152e+00, 8.8739e+00, -8.0807e-01],
[-1.3133e-06, 0.0000e+00, 0.0000e+00, 4.9443e+00, 7.3335e+00,
-1.3381e+00, 7.3543e+00, 1.1118e+01, 1.9361e+00],
[-7.0951e-06, 0.0000e+00, 0.0000e+00, 5.6105e+00, -4.1149e-01,
2.0484e+00, -4.2848e-01, 3.7880e+01, 0.0000e+00]],
[[-5.9600e-06, 0.0000e+00, 0.0000e+00, 9.3076e+02, 1.9217e+03,
-1.4863e+01, 1.9366e+03, -1.4109e+04, 7.0053e+03],
[-3.3294e-05, 0.0000e+00, 0.0000e+00, 6.3506e+03, 3.8622e+02,
-6.1272e+00, 3.7456e+02, 7.4493e+04, 1.3291e+03],
[-4.2133e-05, 0.0000e+00, 0.0000e+00, 7.3284e+03, 6.8571e+03,
2.0607e-01, 6.9009e+03, -6.9545e+01, 0.0000e+00],
[-6.8317e-08, 0.0000e+00, 0.0000e+00, 9.5856e-01, -1.3936e+00,
-7.0740e+00, -1.4021e+00, 1.0685e+01, -3.9045e-01],
[-3.6935e-08, 0.0000e+00, 0.0000e+00, 4.2873e+00, 7.3674e+00,
-1.3401e+00, 7.3902e+00, 2.2436e+00, 2.0187e+00],
[-2.5861e-07, 0.0000e+00, 0.0000e+00, 6.0904e+00, -9.9486e-04,
2.0517e+00, -1.5458e-02, 3.8340e+01, 0.0000e+00]],
[[-2.4774e-04, 0.0000e+00, 0.0000e+00, 9.7854e+02, 1.8374e+03,
-4.2453e+02, 1.8518e+03, -1.3544e+04, 6.9697e+03],
[-1.0340e-03, 0.0000e+00, 0.0000e+00, 6.6198e+03, 8.1223e+02,
-8.3746e+01, 8.0183e+02, 7.4898e+04, 1.4436e+03],
[-1.5230e-03, 0.0000e+00, 0.0000e+00, 7.6672e+03, 6.8443e+03,
1.1904e+02, 6.8872e+03, 2.1420e+03, 0.0000e+00],
[-1.9306e-06, 0.0000e+00, 0.0000e+00, 6.7680e-01, -1.5138e+00,
-7.0599e+00, -1.5227e+00, 8.7423e+00, -8.3637e-01],
[-1.4195e-06, 0.0000e+00, 0.0000e+00, 4.9868e+00, 7.3304e+00,
-1.3378e+00, 7.3510e+00, 1.1718e+01, 1.9303e+00],
[-7.5518e-06, 0.0000e+00, 0.0000e+00, 5.5752e+00, -4.3936e-01,
2.0479e+00, -4.5650e-01, 3.7817e+01, 0.0000e+00]],
[[-1.4731e-04, 0.0000e+00, 0.0000e+00, 9.6359e+02, 1.8678e+03,
-2.8045e+02, 1.8825e+03, -1.3730e+04, 6.9852e+03],
[-6.7025e-04, 0.0000e+00, 0.0000e+00, 6.5205e+03, 6.6255e+02,
-5.6444e+01, 6.5172e+02, 7.4693e+04, 1.4039e+03],
[-9.4178e-04, 0.0000e+00, 0.0000e+00, 7.5514e+03, 6.8517e+03,
7.7244e+01, 6.8949e+03, 1.3674e+03, 0.0000e+00],
[-1.2964e-06, 0.0000e+00, 0.0000e+00, 7.7732e-01, -1.4721e+00,
-7.0678e+00, -1.4809e+00, 9.4567e+00, -6.7970e-01],
[-8.6264e-07, 0.0000e+00, 0.0000e+00, 4.7485e+00, 7.3465e+00,
-1.3392e+00, 7.3680e+00, 8.3954e+00, 1.9622e+00],
[-5.0122e-06, 0.0000e+00, 0.0000e+00, 5.7664e+00, -2.8518e-01,
2.0501e+00, -3.0141e-01, 3.8114e+01, 0.0000e+00]],
[[-1.4681e-04, 0.0000e+00, 0.0000e+00, 9.6351e+02, 1.8680e+03,
-2.7966e+02, 1.8826e+03, -1.3731e+04, 6.9852e+03],
[-6.6830e-04, 0.0000e+00, 0.0000e+00, 6.5199e+03, 6.6173e+02,
-5.6295e+01, 6.5090e+02, 7.4692e+04, 1.4036e+03],
[-9.3878e-04, 0.0000e+00, 0.0000e+00, 7.5507e+03, 6.8517e+03,
7.7015e+01, 6.8950e+03, 1.3632e+03, 0.0000e+00],
[-1.2928e-06, 0.0000e+00, 0.0000e+00, 7.7786e-01, -1.4719e+00,
-7.0679e+00, -1.4807e+00, 9.4605e+00, -6.7885e-01],
[-8.5981e-07, 0.0000e+00, 0.0000e+00, 4.7472e+00, 7.3466e+00,
-1.3392e+00, 7.3681e+00, 8.3772e+00, 1.9624e+00],
[-4.9983e-06, 0.0000e+00, 0.0000e+00, 5.7674e+00, -2.8434e-01,
2.0501e+00, -3.0056e-01, 3.8116e+01, 0.0000e+00]],
[[-1.0787e-04, 0.0000e+00, 0.0000e+00, 9.5628e+02, 1.8811e+03,
-2.1626e+02, 1.8958e+03, -1.3817e+04, 6.9910e+03],
[-5.1263e-04, 0.0000e+00, 0.0000e+00, 6.4778e+03, 5.9581e+02,
-4.4283e+01, 5.8479e+02, 7.4624e+04, 1.3860e+03],
[-7.0388e-04, 0.0000e+00, 0.0000e+00, 7.4986e+03, 6.8540e+03,
5.8624e+01, 6.8974e+03, 1.0209e+03, 0.0000e+00],
[-1.0065e-06, 0.0000e+00, 0.0000e+00, 8.2162e-01, -1.4534e+00,
-7.0703e+00, -1.4621e+00, 9.7642e+00, -6.0985e-01],
[-6.3893e-07, 0.0000e+00, 0.0000e+00, 4.6396e+00, 7.3526e+00,
-1.3396e+00, 7.3744e+00, 6.9115e+00, 1.9762e+00],
[-3.8715e-06, 0.0000e+00, 0.0000e+00, 5.8481e+00, -2.1650e-01,
2.0508e+00, -2.3231e-01, 3.8207e+01, 0.0000e+00]],
[[-2.3238e-04, 0.0000e+00, 0.0000e+00, 9.7653e+02, 1.8418e+03,
-4.0393e+02, 1.8563e+03, -1.3570e+04, 6.9721e+03],
[-9.8107e-04, 0.0000e+00, 0.0000e+00, 6.6053e+03, 7.9084e+02,
-7.9842e+01, 7.8038e+02, 7.4865e+04, 1.4379e+03],
[-1.4359e-03, 0.0000e+00, 0.0000e+00, 7.6508e+03, 6.8456e+03,
1.1306e+02, 6.8885e+03, 2.0316e+03, 0.0000e+00],
[-1.8413e-06, 0.0000e+00, 0.0000e+00, 6.9126e-01, -1.5079e+00,
-7.0612e+00, -1.5168e+00, 8.8465e+00, -8.1398e-01],
[-1.3353e-06, 0.0000e+00, 0.0000e+00, 4.9532e+00, 7.3329e+00,
-1.3380e+00, 7.3537e+00, 1.1244e+01, 1.9349e+00],
[-7.1906e-06, 0.0000e+00, 0.0000e+00, 5.6032e+00, -4.1731e-01,
2.0483e+00, -4.3433e-01, 3.7867e+01, 0.0000e+00]],
[[-2.3819e-04, 0.0000e+00, 0.0000e+00, 9.7730e+02, 1.8401e+03,
-4.1177e+02, 1.8546e+03, -1.3560e+04, 6.9712e+03],
[-1.0012e-03, 0.0000e+00, 0.0000e+00, 6.6108e+03, 7.9898e+02,
-8.1328e+01, 7.8854e+02, 7.4877e+04, 1.4401e+03],
[-1.4689e-03, 0.0000e+00, 0.0000e+00, 7.6571e+03, 6.8451e+03,
1.1534e+02, 6.8880e+03, 2.0736e+03, 0.0000e+00],
[-1.8753e-06, 0.0000e+00, 0.0000e+00, 6.8576e-01, -1.5101e+00,
-7.0607e+00, -1.5190e+00, 8.8070e+00, -8.2250e-01],
[-1.3672e-06, 0.0000e+00, 0.0000e+00, 4.9660e+00, 7.3319e+00,
-1.3379e+00, 7.3527e+00, 1.1424e+01, 1.9331e+00],
[-7.3281e-06, 0.0000e+00, 0.0000e+00, 5.5926e+00, -4.2570e-01,
2.0482e+00, -4.4277e-01, 3.7848e+01, 0.0000e+00]],
[[-2.4610e-04, 0.0000e+00, 0.0000e+00, 9.7833e+02, 1.8378e+03,
-4.2235e+02, 1.8523e+03, -1.3547e+04, 6.9700e+03],
[-1.0284e-03, 0.0000e+00, 0.0000e+00, 6.6183e+03, 8.0997e+02,
-8.3333e+01, 7.9956e+02, 7.4895e+04, 1.4430e+03],
[-1.5137e-03, 0.0000e+00, 0.0000e+00, 7.6654e+03, 6.8444e+03,
1.1841e+02, 6.8874e+03, 2.1303e+03, 0.0000e+00],
[-1.9212e-06, 0.0000e+00, 0.0000e+00, 6.7833e-01, -1.5132e+00,
-7.0600e+00, -1.5221e+00, 8.7534e+00, -8.3400e-01],
[-1.4105e-06, 0.0000e+00, 0.0000e+00, 4.9832e+00, 7.3306e+00,
-1.3378e+00, 7.3513e+00, 1.1668e+01, 1.9307e+00],
[-7.5137e-06, 0.0000e+00, 0.0000e+00, 5.5782e+00, -4.3703e-01,
2.0480e+00, -4.5416e-01, 3.7822e+01, 0.0000e+00]],
[[-2.1251e-04, 0.0000e+00, 0.0000e+00, 9.7379e+02, 1.8476e+03,
-3.7658e+02, 1.8621e+03, -1.3605e+04, 6.9752e+03],
[-9.1130e-04, 0.0000e+00, 0.0000e+00, 6.5862e+03, 7.6243e+02,
-7.4659e+01, 7.5189e+02, 7.4822e+04, 1.4304e+03],
[-1.3223e-03, 0.0000e+00, 0.0000e+00, 7.6290e+03, 6.8471e+03,
1.0513e+02, 6.8901e+03, 1.8848e+03, 0.0000e+00],
[-1.7220e-06, 0.0000e+00, 0.0000e+00, 7.1041e-01, -1.5000e+00,
-7.0629e+00, -1.5089e+00, 8.9838e+00, -7.8425e-01],
[-1.2259e-06, 0.0000e+00, 0.0000e+00, 4.9084e+00, 7.3361e+00,
-1.3383e+00, 7.3570e+00, 1.0613e+01, 1.9410e+00],
[-6.7099e-06, 0.0000e+00, 0.0000e+00, 5.6400e+00, -3.8804e-01,
2.0488e+00, -4.0489e-01, 3.7930e+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([[[-5.2684e-04, 0.0000e+00, 0.0000e+00, -1.2110e+03, 7.8163e+02,
3.7115e+02, 7.5754e+02, -3.1868e+04, -6.4007e+03],
[-2.4933e-03, 0.0000e+00, 0.0000e+00, 2.3721e+03, -8.3704e+02,
-1.7715e+02, -8.8309e+02, -6.8284e+04, 3.0873e+03],
[ 3.8779e-03, 0.0000e+00, 0.0000e+00, -1.3918e+04, 7.0279e+03,
-5.9797e+01, 7.0138e+03, 2.0073e+03, 0.0000e+00],
[-1.7603e-06, 0.0000e+00, 0.0000e+00, 3.1464e+00, -3.2548e+00,
6.6718e+00, -3.2463e+00, 1.2121e+00, 9.2878e-01],
[-1.5116e-06, 0.0000e+00, 0.0000e+00, 6.7924e+00, -6.7287e+00,
-3.1480e+00, -6.7302e+00, -1.1061e+01, 7.9593e-01],
[-6.2441e-06, 0.0000e+00, 0.0000e+00, -3.3939e-01, -4.8852e-01,
-1.1136e+00, -4.3345e-01, 3.8556e+01, 0.0000e+00]],
[[-2.3368e-04, 0.0000e+00, 0.0000e+00, -4.9656e+03, 3.4874e+02,
1.0131e+02, 3.4677e+02, -6.5437e+04, -3.4840e+03],
[-2.5246e-04, 0.0000e+00, 0.0000e+00, -1.0512e+03, -1.0578e+03,
-1.8248e+02, -1.0598e+03, -3.7289e+04, 6.2071e+03],
[ 6.7202e-04, 0.0000e+00, 0.0000e+00, -9.8936e+03, 7.0354e+03,
-2.9983e+01, 7.0407e+03, 9.4683e+02, 0.0000e+00],
[-5.6469e-07, 0.0000e+00, 0.0000e+00, 4.7879e+00, -6.5314e+00,
3.6445e+00, -6.5338e+00, -1.0837e+00, 1.1115e+00],
[ 3.1279e-07, 0.0000e+00, 0.0000e+00, 1.7380e+00, -3.6655e+00,
-6.4273e+00, -3.6673e+00, -7.5476e+00, 3.6338e-01],
[-3.5422e-06, 0.0000e+00, 0.0000e+00, 5.4883e+00, -2.2260e-01,
-1.1435e+00, -2.2001e-01, 3.8999e+01, 0.0000e+00]],
[[-4.6817e-04, 0.0000e+00, 0.0000e+00, 7.1597e+02, -1.6535e+03,
1.2518e+02, -1.6744e+03, 1.3143e+04, -6.8014e+03],
[ 2.4026e-03, 0.0000e+00, 0.0000e+00, -6.8924e+03, -4.5595e+02,
2.6945e+01, -4.5397e+02, -6.9706e+04, -1.3281e+03],
[ 4.1332e-05, 0.0000e+00, 0.0000e+00, 2.3024e+03, 6.7172e+03,
2.8660e+01, 6.7966e+03, 5.0730e+02, 0.0000e+00],
[ 2.3406e-06, 0.0000e+00, 0.0000e+00, -1.5191e+00, 1.4575e+00,
7.2379e+00, 1.4680e+00, -8.8571e+00, 4.9191e-01],
[ 5.4750e-07, 0.0000e+00, 0.0000e+00, -1.8445e+00, -7.4740e+00,
1.3805e+00, -7.5177e+00, -4.5956e+00, -1.8137e+00],
[-9.4309e-06, 0.0000e+00, 0.0000e+00, 7.2292e+00, -1.3226e-01,
1.8695e+00, -1.4092e-01, 3.7718e+01, 0.0000e+00]],
[[-3.6075e-03, 0.0000e+00, 0.0000e+00, 4.8539e+03, -3.5052e+02,
1.4058e+02, -3.2728e+02, 6.5104e+04, -2.7226e+03],
[-6.3343e-04, 0.0000e+00, 0.0000e+00, 6.8206e+02, -1.7159e+03,
3.2583e+02, -1.7398e+03, -2.9232e+04, -6.3834e+03],
[ 8.0119e-03, 0.0000e+00, 0.0000e+00, -1.0166e+04, 6.7268e+03,
8.6180e+01, 6.7757e+03, 1.6659e+03, 0.0000e+00],
[ 5.8631e-06, 0.0000e+00, 0.0000e+00, -5.4018e+00, 6.9817e+00,
2.9650e+00, 7.0064e+00, 3.7482e+00, 1.8834e+00],
[ 3.6664e-06, 0.0000e+00, 0.0000e+00, 9.2192e-01, -2.9695e+00,
6.7150e+00, -2.9882e+00, -1.1998e+01, -3.5372e-01],
[-2.2480e-05, 0.0000e+00, 0.0000e+00, 5.2765e+00, -4.1781e-01,
1.8616e+00, -3.8970e-01, 3.7335e+01, 0.0000e+00]],
[[-7.3754e-04, 0.0000e+00, 0.0000e+00, 5.0057e+03, -1.3256e+03,
4.1705e+01, -1.3257e+03, 4.3142e+04, -5.5206e+03],
[ 5.4188e-04, 0.0000e+00, 0.0000e+00, -3.5519e+03, -1.0868e+03,
3.5508e+01, -1.1210e+03, -5.6793e+04, -4.2150e+03],
[ 1.0184e-03, 0.0000e+00, 0.0000e+00, -7.2112e+03, 6.7327e+03,
1.0052e+01, 6.8145e+03, 1.3215e+02, 0.0000e+00],
[ 1.0434e-06, 0.0000e+00, 0.0000e+00, -3.6796e+00, 4.6223e+00,
5.8590e+00, 4.6451e+00, -6.8175e+00, 1.2101e+00],
[ 1.0931e-07, 0.0000e+00, 0.0000e+00, 2.1992e+00, -6.0425e+00,
4.4577e+00, -6.0839e+00, -6.8796e+00, -1.4309e+00],
[-3.5223e-06, 0.0000e+00, 0.0000e+00, 6.2888e+00, -8.6740e-02,
1.8673e+00, -5.9454e-02, 3.7849e+01, 0.0000e+00]],
[[-1.3036e-05, 0.0000e+00, 0.0000e+00, -6.2593e+02, 7.8016e+02,
6.0252e+03, 7.7724e+02, 7.2674e+03, 4.0696e+03],
[ 6.4906e-05, 0.0000e+00, 0.0000e+00, -3.9981e+03, -6.0642e+03,
2.5309e+02, -6.0413e+03, 4.0663e+04, -6.1784e+02],
[ 7.4163e-05, 0.0000e+00, 0.0000e+00, 6.0654e+03, -4.0067e+03,
7.8885e+02, -3.9917e+03, -6.8661e+04, 0.0000e+00],
[-2.2903e-07, 0.0000e+00, 0.0000e+00, -7.9139e-01, 6.2519e-01,
-4.0406e+00, 6.2402e-01, 4.6992e+00, 6.1245e+00],
[ 1.7905e-06, 0.0000e+00, 0.0000e+00, 6.1059e+00, 4.1117e+00,
-1.6896e-01, 4.1039e+00, -3.0908e+01, 7.8729e-01],
[ 1.1713e-06, 0.0000e+00, 0.0000e+00, 4.0841e+00, -6.1039e+00,
-5.2700e-01, -6.0925e+00, -2.6040e+01, 0.0000e+00]]])