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.3621e+03, -6.9978e+03, 1.3319e+02],
[ 1.9943e+00, 5.1735e-01, 7.2015e+00]],
[[ 1.3398e+03, -7.0031e+03, 5.2983e+01],
[ 2.0109e+00, 4.3152e-01, 7.2026e+00]],
[[ 1.4317e+03, -6.9747e+03, 3.8766e+02],
[ 1.9400e+00, 7.8930e-01, 7.1913e+00]],
[[ 1.3921e+03, -6.9891e+03, 2.4202e+02],
[ 1.9714e+00, 6.3373e-01, 7.1984e+00]],
[[ 1.3997e+03, -6.9866e+03, 2.6990e+02],
[ 1.9655e+00, 6.6353e-01, 7.1973e+00]],
[[ 1.3398e+03, -7.0031e+03, 5.2718e+01],
[ 2.0110e+00, 4.3124e-01, 7.2026e+00]],
[[ 1.4259e+03, -6.9770e+03, 3.6618e+02],
[ 1.9447e+00, 7.6638e-01, 7.1926e+00]],
[[ 1.4197e+03, -6.9794e+03, 3.4333e+02],
[ 1.9496e+00, 7.4197e-01, 7.1938e+00]],
[[ 1.3385e+03, -7.0034e+03, 4.8278e+01],
[ 2.0119e+00, 4.2648e-01, 7.2027e+00]],
[[ 1.4024e+03, -6.9856e+03, 2.7986e+02],
[ 1.9633e+00, 6.7417e-01, 7.1969e+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.1966e+02, 3.1041e+01, 4.3209e+02],
[-8.9977e-02, 4.6225e-01, -8.8211e-03]],
[[ 1.2066e+02, 2.5891e+01, 4.3216e+02],
[-8.8508e-02, 4.6262e-01, -3.5094e-03]],
[[ 1.1640e+02, 4.7359e+01, 4.3148e+02],
[-9.4559e-02, 4.6066e-01, -2.5670e-02]],
[[ 1.1829e+02, 3.8024e+01, 4.3190e+02],
[-9.1952e-02, 4.6165e-01, -1.6028e-02]],
[[ 1.1793e+02, 3.9812e+01, 4.3184e+02],
[-9.2455e-02, 4.6148e-01, -1.7874e-02]],
[[ 1.2066e+02, 2.5874e+01, 4.3216e+02],
[-8.8503e-02, 4.6262e-01, -3.4918e-03]],
[[ 1.1668e+02, 4.5983e+01, 4.3156e+02],
[-9.4178e-02, 4.6082e-01, -2.4249e-02]],
[[ 1.1698e+02, 4.4519e+01, 4.3163e+02],
[-9.3770e-02, 4.6098e-01, -2.2736e-02]],
[[ 1.2071e+02, 2.5589e+01, 4.3216e+02],
[-8.8421e-02, 4.6263e-01, -3.1978e-03]],
[[ 1.1780e+02, 4.0451e+01, 4.3181e+02],
[-9.2634e-02, 4.6141e-01, -1.8533e-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.6271e+01, -3.7924e+01, 4.4360e+02],
[-2.0235e-01, -4.2667e-01, -6.8009e-03]],
[[ 1.0156e+01, -7.3195e+01, 4.4299e+02],
[-4.1217e-01, -2.2888e-01, -2.7755e-02]],
[[-1.0880e+02, -2.9612e+01, 4.4232e+02],
[ 9.5554e-02, -4.8920e-01, -9.2767e-03]],
[[-3.4161e+01, -1.0743e+02, 4.4116e+02],
[ 4.5446e-01, -1.9796e-01, -1.0386e-02]],
[[-7.4121e+01, -8.7815e+01, 4.4139e+02],
[ 3.0508e-01, -3.9105e-01, -2.4288e-02]],
[[ 4.5850e+01, -3.7617e+02, -2.2938e+02],
[ 3.9597e-02, 2.3599e-01, -3.7921e-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\)]
#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([[[-1.9316e-04, 0.0000e+00, 0.0000e+00, 9.7097e+02, 1.8534e+03,
-3.4915e+02, 1.8680e+03, -1.3640e+04, 6.9782e+03],
[-8.4188e-04, 0.0000e+00, 0.0000e+00, 6.5672e+03, 7.3394e+02,
-6.9461e+01, 7.2331e+02, 7.4782e+04, 1.4229e+03],
[-1.2108e-03, 0.0000e+00, 0.0000e+00, 7.6071e+03, 6.8486e+03,
9.7172e+01, 6.8917e+03, 1.7374e+03, 0.0000e+00],
[-1.6016e-06, 0.0000e+00, 0.0000e+00, 7.2957e-01, -1.4921e+00,
-7.0644e+00, -1.5009e+00, 9.1203e+00, -7.5443e-01],
[-1.1188e-06, 0.0000e+00, 0.0000e+00, 4.8632e+00, 7.3393e+00,
-1.3386e+00, 7.3603e+00, 9.9811e+00, 1.9471e+00],
[-6.2269e-06, 0.0000e+00, 0.0000e+00, 5.6765e+00, -3.5869e-01,
2.0492e+00, -3.7536e-01, 3.7988e+01, 0.0000e+00]],
[[-1.2926e-05, 0.0000e+00, 0.0000e+00, 9.3297e+02, 1.9184e+03,
-3.1094e+01, 1.9334e+03, -1.4084e+04, 7.0043e+03],
[-7.1155e-05, 0.0000e+00, 0.0000e+00, 6.3605e+03, 4.0312e+02,
-9.2020e+00, 3.9151e+02, 7.4499e+04, 1.3337e+03],
[-9.0691e-05, 0.0000e+00, 0.0000e+00, 7.3424e+03, 6.8570e+03,
4.9139e+00, 6.9008e+03, 1.8422e+01, 0.0000e+00],
[-1.4551e-07, 0.0000e+00, 0.0000e+00, 9.4765e-01, -1.3985e+00,
-7.0739e+00, -1.4070e+00, 1.0614e+01, -4.0814e-01],
[-7.9759e-08, 0.0000e+00, 0.0000e+00, 4.3163e+00, 7.3665e+00,
-1.3401e+00, 7.3892e+00, 2.6204e+00, 2.0154e+00],
[-5.5154e-07, 0.0000e+00, 0.0000e+00, 6.0717e+00, -1.8362e-02,
2.0517e+00, -3.2935e-02, 3.8338e+01, 0.0000e+00]],
[[-6.9267e-05, 0.0000e+00, 0.0000e+00, 9.4800e+02, 1.8952e+03,
-1.4740e+02, 1.9100e+03, -1.3914e+04, 6.9966e+03],
[-3.4627e-04, 0.0000e+00, 0.0000e+00, 6.4332e+03, 5.2418e+02,
-3.1235e+01, 5.1294e+02, 7.4564e+04, 1.3667e+03],
[-4.6314e-04, 0.0000e+00, 0.0000e+00, 7.4412e+03, 6.8557e+03,
3.8649e+01, 6.8993e+03, 6.4843e+02, 0.0000e+00],
[-6.9056e-07, 0.0000e+00, 0.0000e+00, 8.6881e-01, -1.4331e+00,
-7.0723e+00, -1.4417e+00, 1.0087e+01, -5.3487e-01],
[-4.1591e-07, 0.0000e+00, 0.0000e+00, 4.5209e+00, 7.3584e+00,
-1.3399e+00, 7.3806e+00, 5.3173e+00, 1.9909e+00],
[-2.6417e-06, 0.0000e+00, 0.0000e+00, 5.9334e+00, -1.4281e-01,
2.0513e+00, -1.5816e-01, 3.8280e+01, 0.0000e+00]],
[[-1.4103e-04, 0.0000e+00, 0.0000e+00, 9.6250e+02, 1.8699e+03,
-2.7059e+02, 1.8845e+03, -1.3743e+04, 6.9861e+03],
[-6.4586e-04, 0.0000e+00, 0.0000e+00, 6.5138e+03, 6.5229e+02,
-5.4575e+01, 6.4143e+02, 7.4681e+04, 1.4011e+03],
[-9.0438e-04, 0.0000e+00, 0.0000e+00, 7.5433e+03, 6.8521e+03,
7.4382e+01, 6.8954e+03, 1.3142e+03, 0.0000e+00],
[-1.2521e-06, 0.0000e+00, 0.0000e+00, 7.8415e-01, -1.4693e+00,
-7.0683e+00, -1.4780e+00, 9.5044e+00, -6.6897e-01],
[-8.2728e-07, 0.0000e+00, 0.0000e+00, 4.7319e+00, 7.3475e+00,
-1.3392e+00, 7.3690e+00, 8.1675e+00, 1.9644e+00],
[-4.8373e-06, 0.0000e+00, 0.0000e+00, 5.7791e+00, -2.7463e-01,
2.0503e+00, -2.9079e-01, 3.8130e+01, 0.0000e+00]],
[[-2.1867e-04, 0.0000e+00, 0.0000e+00, 9.7466e+02, 1.8458e+03,
-3.8515e+02, 1.8603e+03, -1.3594e+04, 6.9742e+03],
[-9.3311e-04, 0.0000e+00, 0.0000e+00, 6.5921e+03, 7.7133e+02,
-7.6283e+01, 7.6082e+02, 7.4835e+04, 1.4328e+03],
[-1.3577e-03, 0.0000e+00, 0.0000e+00, 7.6359e+03, 6.8466e+03,
1.0762e+02, 6.8896e+03, 1.9308e+03, 0.0000e+00],
[-1.7595e-06, 0.0000e+00, 0.0000e+00, 7.0441e-01, -1.5025e+00,
-7.0624e+00, -1.5113e+00, 8.9409e+00, -7.9357e-01],
[-1.2599e-06, 0.0000e+00, 0.0000e+00, 4.9225e+00, 7.3351e+00,
-1.3382e+00, 7.3560e+00, 1.0811e+01, 1.9391e+00],
[-6.8607e-06, 0.0000e+00, 0.0000e+00, 5.6285e+00, -3.9722e-01,
2.0486e+00, -4.1411e-01, 3.7911e+01, 0.0000e+00]],
[[-2.0382e-04, 0.0000e+00, 0.0000e+00, 9.7255e+02, 1.8502e+03,
-3.6437e+02, 1.8647e+03, -1.3620e+04, 6.9765e+03],
[-8.8034e-04, 0.0000e+00, 0.0000e+00, 6.5777e+03, 7.4975e+02,
-7.2346e+01, 7.3917e+02, 7.4804e+04, 1.4271e+03],
[-1.2724e-03, 0.0000e+00, 0.0000e+00, 7.6193e+03, 6.8478e+03,
1.0159e+02, 6.8908e+03, 1.8192e+03, 0.0000e+00],
[-1.6685e-06, 0.0000e+00, 0.0000e+00, 7.1894e-01, -1.4965e+00,
-7.0636e+00, -1.5053e+00, 9.0447e+00, -7.7098e-01],
[-1.1779e-06, 0.0000e+00, 0.0000e+00, 4.8883e+00, 7.3375e+00,
-1.3384e+00, 7.3585e+00, 1.0332e+01, 1.9437e+00],
[-6.4951e-06, 0.0000e+00, 0.0000e+00, 5.6563e+00, -3.7498e-01,
2.0490e+00, -3.9175e-01, 3.7956e+01, 0.0000e+00]],
[[-2.2435e-04, 0.0000e+00, 0.0000e+00, 9.7544e+02, 1.8441e+03,
-3.9298e+02, 1.8586e+03, -1.3584e+04, 6.9734e+03],
[-9.5306e-04, 0.0000e+00, 0.0000e+00, 6.5976e+03, 7.7946e+02,
-7.7766e+01, 7.6897e+02, 7.4847e+04, 1.4349e+03],
[-1.3901e-03, 0.0000e+00, 0.0000e+00, 7.6421e+03, 6.8462e+03,
1.0989e+02, 6.8892e+03, 1.9728e+03, 0.0000e+00],
[-1.7936e-06, 0.0000e+00, 0.0000e+00, 6.9893e-01, -1.5047e+00,
-7.0619e+00, -1.5136e+00, 8.9016e+00, -8.0207e-01],
[-1.2912e-06, 0.0000e+00, 0.0000e+00, 4.9353e+00, 7.3342e+00,
-1.3382e+00, 7.3550e+00, 1.0991e+01, 1.9373e+00],
[-6.9982e-06, 0.0000e+00, 0.0000e+00, 5.6180e+00, -4.0559e-01,
2.0485e+00, -4.2254e-01, 3.7893e+01, 0.0000e+00]],
[[-8.3525e-05, 0.0000e+00, 0.0000e+00, 9.5121e+02, 1.8898e+03,
-1.7364e+02, 1.9046e+03, -1.3877e+04, 6.9946e+03],
[-4.0935e-04, 0.0000e+00, 0.0000e+00, 6.4501e+03, 5.5148e+02,
-3.6208e+01, 5.4033e+02, 7.4585e+04, 1.3740e+03],
[-5.5312e-04, 0.0000e+00, 0.0000e+00, 7.4632e+03, 6.8551e+03,
4.6262e+01, 6.8986e+03, 7.9045e+02, 0.0000e+00],
[-8.1156e-07, 0.0000e+00, 0.0000e+00, 8.5087e-01, -1.4408e+00,
-7.0716e+00, -1.4495e+00, 9.9647e+00, -5.6345e-01],
[-4.9882e-07, 0.0000e+00, 0.0000e+00, 4.5664e+00, 7.3563e+00,
-1.3398e+00, 7.3783e+00, 5.9252e+00, 1.9853e+00],
[-3.1111e-06, 0.0000e+00, 0.0000e+00, 5.9012e+00, -1.7089e-01,
2.0512e+00, -1.8642e-01, 3.8255e+01, 0.0000e+00]],
[[-2.0803e-04, 0.0000e+00, 0.0000e+00, 9.7315e+02, 1.8489e+03,
-3.7030e+02, 1.8635e+03, -1.3613e+04, 6.9759e+03],
[-8.9537e-04, 0.0000e+00, 0.0000e+00, 6.5818e+03, 7.5591e+02,
-7.3470e+01, 7.4535e+02, 7.4813e+04, 1.4287e+03],
[-1.2966e-03, 0.0000e+00, 0.0000e+00, 7.6240e+03, 6.8475e+03,
1.0331e+02, 6.8905e+03, 1.8511e+03, 0.0000e+00],
[-1.6946e-06, 0.0000e+00, 0.0000e+00, 7.1480e-01, -1.4982e+00,
-7.0632e+00, -1.5070e+00, 9.0151e+00, -7.7743e-01],
[-1.2011e-06, 0.0000e+00, 0.0000e+00, 4.8981e+00, 7.3369e+00,
-1.3384e+00, 7.3578e+00, 1.0469e+01, 1.9424e+00],
[-6.5995e-06, 0.0000e+00, 0.0000e+00, 5.6484e+00, -3.8133e-01,
2.0489e+00, -3.9813e-01, 3.7943e+01, 0.0000e+00]],
[[-3.0418e-05, 0.0000e+00, 0.0000e+00, 9.3811e+02, 1.9108e+03,
-6.9776e+01, 1.9257e+03, -1.4027e+04, 7.0020e+03],
[-1.6190e-04, 0.0000e+00, 0.0000e+00, 6.3843e+03, 4.4339e+02,
-1.6530e+01, 4.3191e+02, 7.4516e+04, 1.3447e+03],
[-2.0980e-04, 0.0000e+00, 0.0000e+00, 7.3755e+03, 6.8568e+03,
1.6134e+01, 6.9005e+03, 2.2805e+02, 0.0000e+00],
[-3.2838e-07, 0.0000e+00, 0.0000e+00, 9.2155e-01, -1.4100e+00,
-7.0736e+00, -1.4186e+00, 1.0441e+01, -4.5030e-01],
[-1.8587e-07, 0.0000e+00, 0.0000e+00, 4.3849e+00, 7.3640e+00,
-1.3401e+00, 7.3865e+00, 3.5180e+00, 2.0073e+00],
[-1.2485e-06, 0.0000e+00, 0.0000e+00, 6.0265e+00, -5.9752e-02,
2.0517e+00, -7.4588e-02, 3.8328e+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.1568e-04, 0.0000e+00, 0.0000e+00, -1.2145e+03, 7.8521e+02,
3.6380e+02, 7.6112e+02, -3.1869e+04, -6.4017e+03],
[-2.4428e-03, 0.0000e+00, 0.0000e+00, 2.3646e+03, -8.2963e+02,
-1.7368e+02, -8.7568e+02, -6.8272e+04, 3.0864e+03],
[ 3.8023e-03, 0.0000e+00, 0.0000e+00, -1.3918e+04, 7.0285e+03,
-5.8571e+01, 7.0143e+03, 1.9649e+03, 0.0000e+00],
[-1.7243e-06, 0.0000e+00, 0.0000e+00, 3.1436e+00, -3.2538e+00,
6.6723e+00, -3.2453e+00, 1.2924e+00, 9.2098e-01],
[-1.4815e-06, 0.0000e+00, 0.0000e+00, 6.7927e+00, -6.7297e+00,
-3.1482e+00, -6.7310e+00, -1.0897e+01, 7.9969e-01],
[-6.1144e-06, 0.0000e+00, 0.0000e+00, -3.5656e-01, -4.7994e-01,
-1.1137e+00, -4.2487e-01, 3.8574e+01, 0.0000e+00]],
[[-1.4052e-04, 0.0000e+00, 0.0000e+00, -5.0189e+03, 4.2069e+02,
6.1150e+01, 4.1874e+02, -6.5434e+04, -3.4960e+03],
[-1.5491e-04, 0.0000e+00, 0.0000e+00, -1.0707e+03, -1.0173e+03,
-1.1165e+02, -1.0193e+03, -3.7211e+04, 6.2027e+03],
[ 4.2141e-04, 0.0000e+00, 0.0000e+00, -9.9535e+03, 7.0374e+03,
-1.7381e+01, 7.0426e+03, 5.1667e+02, 0.0000e+00],
[-3.5199e-07, 0.0000e+00, 0.0000e+00, 4.8897e+00, -6.5265e+00,
3.6455e+00, -6.5288e+00, 5.1712e-01, 1.0690e+00],
[ 1.8480e-07, 0.0000e+00, 0.0000e+00, 1.8165e+00, -3.6780e+00,
-6.4290e+00, -3.6798e+00, -6.6527e+00, 4.3895e-01],
[-2.1510e-06, 0.0000e+00, 0.0000e+00, 5.3715e+00, -1.3668e-01,
-1.1438e+00, -1.3403e-01, 3.9073e+01, 0.0000e+00]],
[[-6.6793e-04, 0.0000e+00, 0.0000e+00, 7.0404e+02, -1.6419e+03,
1.8236e+02, -1.6628e+03, 1.3075e+04, -6.7972e+03],
[ 3.5078e-03, 0.0000e+00, 0.0000e+00, -6.9074e+03, -5.1498e+02,
3.7852e+01, -5.1335e+02, -6.9747e+04, -1.3424e+03],
[ 3.9749e-07, 0.0000e+00, 0.0000e+00, 2.3595e+03, 6.7159e+03,
4.3430e+01, 6.7952e+03, 8.0511e+02, 0.0000e+00],
[ 3.3928e-06, 0.0000e+00, 0.0000e+00, -1.4883e+00, 1.4731e+00,
7.2364e+00, 1.4838e+00, -8.5927e+00, 5.5631e-01],
[ 9.0267e-07, 0.0000e+00, 0.0000e+00, -1.9703e+00, -7.4696e+00,
1.3803e+00, -7.5131e+00, -5.9136e+00, -1.8010e+00],
[-1.3751e-05, 0.0000e+00, 0.0000e+00, 7.2026e+00, -1.9605e-01,
1.8691e+00, -2.0547e-01, 3.7666e+01, 0.0000e+00]],
[[-3.2387e-03, 0.0000e+00, 0.0000e+00, 4.8820e+03, -3.8657e+02,
1.2526e+02, -3.6346e+02, 6.5087e+04, -2.7323e+03],
[-5.7567e-04, 0.0000e+00, 0.0000e+00, 6.7722e+02, -1.7006e+03,
2.9114e+02, -1.7243e+03, -2.9171e+04, -6.3815e+03],
[ 7.2178e-03, 0.0000e+00, 0.0000e+00, -1.0193e+04, 6.7289e+03,
7.6563e+01, 6.7776e+03, 1.4729e+03, 0.0000e+00],
[ 5.3125e-06, 0.0000e+00, 0.0000e+00, -5.4514e+00, 6.9793e+00,
2.9658e+00, 7.0037e+00, 2.9421e+00, 1.8667e+00],
[ 3.2257e-06, 0.0000e+00, 0.0000e+00, 9.6010e-01, -2.9799e+00,
6.7169e+00, -2.9986e+00, -1.1662e+01, -3.9294e-01],
[-2.0023e-05, 0.0000e+00, 0.0000e+00, 5.2179e+00, -3.7635e-01,
1.8621e+00, -3.4792e-01, 3.7408e+01, 0.0000e+00]],
[[-5.5013e-03, 0.0000e+00, 0.0000e+00, 4.8308e+03, -1.0932e+03,
3.3372e+02, -1.0921e+03, 4.2932e+04, -5.4521e+03],
[ 4.2476e-03, 0.0000e+00, 0.0000e+00, -3.4582e+03, -1.3865e+03,
2.5769e+02, -1.4228e+03, -5.7303e+04, -4.2801e+03],
[ 7.1096e-03, 0.0000e+00, 0.0000e+00, -6.8873e+03, 6.7184e+03,
1.0313e+02, 6.8014e+03, 2.0108e+03, 0.0000e+00],
[ 7.4066e-06, 0.0000e+00, 0.0000e+00, -3.3154e+00, 4.6949e+00,
5.8475e+00, 4.7204e+00, -1.5880e+00, 1.5363e+00],
[ 2.0409e-06, 0.0000e+00, 0.0000e+00, 1.5596e+00, -5.9699e+00,
4.4493e+00, -6.0128e+00, -1.3548e+01, -1.1784e+00],
[-2.8011e-05, 0.0000e+00, 0.0000e+00, 6.6721e+00, -4.8765e-01,
1.8639e+00, -4.6543e-01, 3.7348e+01, 0.0000e+00]],
[[-2.2730e-06, 0.0000e+00, 0.0000e+00, -6.0440e+02, 7.6248e+02,
6.1351e+03, 7.5959e+02, 7.1420e+03, 3.8978e+03],
[ 7.4537e-06, 0.0000e+00, 0.0000e+00, -4.1711e+03, -6.1760e+03,
2.5769e+02, -6.1530e+03, 4.1557e+04, -6.3947e+02],
[ 1.3718e-05, 0.0000e+00, 0.0000e+00, 5.9568e+03, -3.8355e+03,
8.0319e+02, -3.8208e+03, -6.7989e+04, 0.0000e+00],
[-6.3153e-08, 0.0000e+00, 0.0000e+00, -7.5532e-01, 6.4705e-01,
-3.8679e+00, 6.4580e-01, 4.3177e+00, 6.2377e+00],
[ 5.0724e-07, 0.0000e+00, 0.0000e+00, 6.3291e+00, 3.9375e+00,
-1.6157e-01, 3.9304e+00, -3.3454e+01, 7.6943e-01],
[ 3.1867e-07, 0.0000e+00, 0.0000e+00, 3.7360e+00, -6.2151e+00,
-5.0453e-01, -6.2032e+00, -2.2292e+01, 0.0000e+00]]])