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.4325e+03, -6.9743e+03, 3.9084e+02],
[ 1.9393e+00, 7.9271e-01, 7.1911e+00]],
[[ 1.3737e+03, -6.9947e+03, 1.7493e+02],
[ 1.9856e+00, 5.6200e-01, 7.2005e+00]],
[[ 1.4198e+03, -6.9794e+03, 3.4383e+02],
[ 1.9495e+00, 7.4250e-01, 7.1938e+00]],
[[ 1.3989e+03, -6.9869e+03, 2.6669e+02],
[ 1.9662e+00, 6.6010e-01, 7.1974e+00]],
[[ 1.4198e+03, -6.9794e+03, 3.4384e+02],
[ 1.9495e+00, 7.4252e-01, 7.1938e+00]],
[[ 1.3602e+03, -6.9983e+03, 1.2617e+02],
[ 1.9958e+00, 5.0985e-01, 7.2016e+00]],
[[ 1.3902e+03, -6.9897e+03, 2.3492e+02],
[ 1.9729e+00, 6.2614e-01, 7.1986e+00]],
[[ 1.3555e+03, -6.9995e+03, 1.0936e+02],
[ 1.9993e+00, 4.9186e-01, 7.2019e+00]],
[[ 1.3654e+03, -6.9969e+03, 1.4486e+02],
[ 1.9919e+00, 5.2983e-01, 7.2013e+00]],
[[ 1.4277e+03, -6.9763e+03, 3.7290e+02],
[ 1.9432e+00, 7.7355e-01, 7.1922e+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.1636e+02, 4.7563e+01, 4.3147e+02],
[-9.4616e-02, 4.6064e-01, -2.5881e-02]],
[[ 1.1914e+02, 3.3720e+01, 4.3203e+02],
[-9.0738e-02, 4.6203e-01, -1.1585e-02]],
[[ 1.1697e+02, 4.4550e+01, 4.3163e+02],
[-9.3779e-02, 4.6098e-01, -2.2768e-02]],
[[ 1.1797e+02, 3.9606e+01, 4.3185e+02],
[-9.2397e-02, 4.6150e-01, -1.7661e-02]],
[[ 1.1697e+02, 4.4551e+01, 4.3163e+02],
[-9.3779e-02, 4.6098e-01, -2.2769e-02]],
[[ 1.1975e+02, 3.0591e+01, 4.3210e+02],
[-8.9849e-02, 4.6228e-01, -8.3564e-03]],
[[ 1.1838e+02, 3.7569e+01, 4.3192e+02],
[-9.1824e-02, 4.6169e-01, -1.5558e-02]],
[[ 1.1996e+02, 2.9512e+01, 4.3212e+02],
[-8.9542e-02, 4.6236e-01, -7.2429e-03]],
[[ 1.1952e+02, 3.1790e+01, 4.3208e+02],
[-9.0190e-02, 4.6219e-01, -9.5937e-03]],
[[ 1.1659e+02, 4.6413e+01, 4.3153e+02],
[-9.4297e-02, 4.6077e-01, -2.4693e-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.8500e+01, -3.3212e+01, 4.4365e+02],
[-2.0159e-01, -4.2698e-01, -1.3732e-03]],
[[ 1.8242e+01, -6.8688e+01, 4.4344e+02],
[-4.1183e-01, -2.3040e-01, -1.8113e-02]],
[[-1.0468e+02, -5.0116e+01, 4.4146e+02],
[ 1.0092e-01, -4.8717e-01, -3.1585e-02]],
[[-2.9611e+01, -1.0941e+02, 4.4103e+02],
[ 4.5492e-01, -1.9670e-01, -1.5651e-02]],
[[-7.6181e+01, -8.5170e+01, 4.4154e+02],
[ 3.0445e-01, -3.9171e-01, -2.0721e-02]],
[[ 4.7581e+01, -3.6518e+02, -2.4617e+02],
[ 3.7457e-02, 2.5304e-01, -3.6831e-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([[[-8.0801e-05, 0.0000e+00, 0.0000e+00, 9.5061e+02, 1.8908e+03,
-1.6871e+02, 1.9056e+03, -1.3884e+04, 6.9949e+03],
[-3.9746e-04, 0.0000e+00, 0.0000e+00, 6.4469e+03, 5.4635e+02,
-3.5273e+01, 5.3518e+02, 7.4581e+04, 1.3727e+03],
[-5.3604e-04, 0.0000e+00, 0.0000e+00, 7.4590e+03, 6.8553e+03,
4.4830e+01, 6.8988e+03, 7.6375e+02, 0.0000e+00],
[-7.8887e-07, 0.0000e+00, 0.0000e+00, 8.5425e-01, -1.4394e+00,
-7.0718e+00, -1.4480e+00, 9.9877e+00, -5.5808e-01],
[-4.8304e-07, 0.0000e+00, 0.0000e+00, 4.5578e+00, 7.3567e+00,
-1.3398e+00, 7.3787e+00, 5.8109e+00, 1.9864e+00],
[-3.0229e-06, 0.0000e+00, 0.0000e+00, 5.9072e+00, -1.6561e-01,
2.0512e+00, -1.8111e-01, 3.8260e+01, 0.0000e+00]],
[[-3.5563e-08, 0.0000e+00, 0.0000e+00, 9.2882e+02, 1.9245e+03,
-6.3346e-01, 1.9394e+03, -1.4130e+04, 7.0060e+03],
[-2.0128e-07, 0.0000e+00, 0.0000e+00, 6.3420e+03, 3.7140e+02,
-3.4316e+00, 3.5969e+02, 7.4489e+04, 1.3250e+03],
[-2.5311e-07, 0.0000e+00, 0.0000e+00, 7.3161e+03, 6.8571e+03,
-3.9212e+00, 6.9009e+03, -1.4667e+02, 0.0000e+00],
[-4.1424e-10, 0.0000e+00, 0.0000e+00, 9.6811e-01, -1.3894e+00,
-7.0740e+00, -1.3978e+00, 1.0748e+01, -3.7493e-01],
[-2.2124e-10, 0.0000e+00, 0.0000e+00, 4.2618e+00, 7.3682e+00,
-1.3401e+00, 7.3910e+00, 1.9132e+00, 2.0217e+00],
[-1.5663e-09, 0.0000e+00, 0.0000e+00, 6.1067e+00, 1.4230e-02,
2.0517e+00, -1.3554e-04, 3.8340e+01, 0.0000e+00]],
[[-6.2966e-05, 0.0000e+00, 0.0000e+00, 9.4651e+02, 1.8976e+03,
-1.3545e+02, 1.9124e+03, -1.3931e+04, 6.9975e+03],
[-3.1770e-04, 0.0000e+00, 0.0000e+00, 6.4256e+03, 5.1175e+02,
-2.8973e+01, 5.0048e+02, 7.4555e+04, 1.3633e+03],
[-4.2292e-04, 0.0000e+00, 0.0000e+00, 7.4311e+03, 6.8560e+03,
3.5184e+01, 6.8995e+03, 5.8378e+02, 0.0000e+00],
[-6.3525e-07, 0.0000e+00, 0.0000e+00, 8.7696e-01, -1.4295e+00,
-7.0726e+00, -1.4382e+00, 1.0142e+01, -5.2186e-01],
[-3.7904e-07, 0.0000e+00, 0.0000e+00, 4.5001e+00, 7.3594e+00,
-1.3399e+00, 7.3815e+00, 5.0406e+00, 1.9935e+00],
[-2.4278e-06, 0.0000e+00, 0.0000e+00, 5.9479e+00, -1.3003e-01,
2.0514e+00, -1.4531e-01, 3.8289e+01, 0.0000e+00]],
[[-4.2987e-05, 0.0000e+00, 0.0000e+00, 9.4151e+02, 1.9055e+03,
-9.6002e+01, 1.9204e+03, -1.3988e+04, 7.0003e+03],
[-2.2384e-04, 0.0000e+00, 0.0000e+00, 6.4007e+03, 4.7069e+02,
-2.1499e+01, 4.5929e+02, 7.4530e+04, 1.3522e+03],
[-2.9326e-04, 0.0000e+00, 0.0000e+00, 7.3978e+03, 6.8566e+03,
2.3741e+01, 6.9002e+03, 3.7013e+02, 0.0000e+00],
[-4.5146e-07, 0.0000e+00, 0.0000e+00, 9.0378e-01, -1.4178e+00,
-7.0733e+00, -1.4264e+00, 1.0322e+01, -4.7888e-01],
[-2.6105e-07, 0.0000e+00, 0.0000e+00, 4.4311e+00, 7.3623e+00,
-1.3400e+00, 7.3846e+00, 4.1262e+00, 2.0018e+00],
[-1.7200e-06, 0.0000e+00, 0.0000e+00, 5.9954e+00, -8.7814e-02,
2.0516e+00, -1.0283e-01, 3.8315e+01, 0.0000e+00]],
[[-1.3209e-05, 0.0000e+00, 0.0000e+00, 9.3305e+02, 1.9183e+03,
-3.1743e+01, 1.9333e+03, -1.4083e+04, 7.0043e+03],
[-7.2672e-05, 0.0000e+00, 0.0000e+00, 6.3609e+03, 4.0379e+02,
-9.3250e+00, 3.9219e+02, 7.4499e+04, 1.3339e+03],
[-9.2651e-05, 0.0000e+00, 0.0000e+00, 7.3429e+03, 6.8570e+03,
5.1022e+00, 6.9008e+03, 2.1941e+01, 0.0000e+00],
[-1.4860e-07, 0.0000e+00, 0.0000e+00, 9.4721e-01, -1.3987e+00,
-7.0739e+00, -1.4072e+00, 1.0611e+01, -4.0885e-01],
[-8.1493e-08, 0.0000e+00, 0.0000e+00, 4.3174e+00, 7.3664e+00,
-1.3401e+00, 7.3891e+00, 2.6355e+00, 2.0152e+00],
[-5.6325e-07, 0.0000e+00, 0.0000e+00, 6.0709e+00, -1.9057e-02,
2.0517e+00, -3.3635e-02, 3.8338e+01, 0.0000e+00]],
[[-1.8142e-04, 0.0000e+00, 0.0000e+00, 9.6918e+02, 1.8570e+03,
-3.3208e+02, 1.8716e+03, -1.3662e+04, 6.9800e+03],
[-7.9895e-04, 0.0000e+00, 0.0000e+00, 6.5555e+03, 7.1621e+02,
-6.6228e+01, 7.0553e+02, 7.4759e+04, 1.4182e+03],
[-1.1426e-03, 0.0000e+00, 0.0000e+00, 7.5933e+03, 6.8494e+03,
9.2221e+01, 6.8925e+03, 1.6456e+03, 0.0000e+00],
[-1.5263e-06, 0.0000e+00, 0.0000e+00, 7.4146e-01, -1.4871e+00,
-7.0654e+00, -1.4959e+00, 9.2046e+00, -7.3587e-01],
[-1.0536e-06, 0.0000e+00, 0.0000e+00, 4.8349e+00, 7.3411e+00,
-1.3388e+00, 7.3623e+00, 9.5874e+00, 1.9509e+00],
[-5.9258e-06, 0.0000e+00, 0.0000e+00, 5.6991e+00, -3.4043e-01,
2.0495e+00, -3.5699e-01, 3.8022e+01, 0.0000e+00]],
[[-2.6726e-05, 0.0000e+00, 0.0000e+00, 9.3706e+02, 1.9124e+03,
-6.1836e+01, 1.9272e+03, -1.4038e+04, 7.0025e+03],
[-1.4321e-04, 0.0000e+00, 0.0000e+00, 6.3794e+03, 4.3513e+02,
-1.5026e+01, 4.2362e+02, 7.4512e+04, 1.3425e+03],
[-1.8496e-04, 0.0000e+00, 0.0000e+00, 7.3687e+03, 6.8569e+03,
1.3831e+01, 6.9006e+03, 1.8502e+02, 0.0000e+00],
[-2.9097e-07, 0.0000e+00, 0.0000e+00, 9.2692e-01, -1.4077e+00,
-7.0737e+00, -1.4162e+00, 1.0477e+01, -4.4165e-01],
[-1.6363e-07, 0.0000e+00, 0.0000e+00, 4.3709e+00, 7.3646e+00,
-1.3401e+00, 7.3871e+00, 3.3338e+00, 2.0090e+00],
[-1.1056e-06, 0.0000e+00, 0.0000e+00, 6.0358e+00, -5.1256e-02,
2.0517e+00, -6.6038e-02, 3.8331e+01, 0.0000e+00]],
[[-5.8052e-05, 0.0000e+00, 0.0000e+00, 9.4533e+02, 1.8995e+03,
-1.2598e+02, 1.9143e+03, -1.3945e+04, 6.9982e+03],
[-2.9508e-04, 0.0000e+00, 0.0000e+00, 6.4196e+03, 5.0189e+02,
-2.7178e+01, 4.9059e+02, 7.4549e+04, 1.3606e+03],
[-3.9133e-04, 0.0000e+00, 0.0000e+00, 7.4232e+03, 6.8561e+03,
3.2437e+01, 6.8997e+03, 5.3250e+02, 0.0000e+00],
[-5.9127e-07, 0.0000e+00, 0.0000e+00, 8.8341e-01, -1.4267e+00,
-7.0728e+00, -1.4354e+00, 1.0185e+01, -5.1154e-01],
[-3.5017e-07, 0.0000e+00, 0.0000e+00, 4.4836e+00, 7.3601e+00,
-1.3400e+00, 7.3823e+00, 4.8211e+00, 1.9955e+00],
[-2.2580e-06, 0.0000e+00, 0.0000e+00, 5.9594e+00, -1.1989e-01,
2.0515e+00, -1.3511e-01, 3.8297e+01, 0.0000e+00]],
[[-9.7713e-06, 0.0000e+00, 0.0000e+00, 9.3198e+02, 1.9199e+03,
-2.3807e+01, 1.9348e+03, -1.4095e+04, 7.0048e+03],
[-5.4142e-05, 0.0000e+00, 0.0000e+00, 6.3561e+03, 3.9553e+02,
-7.8216e+00, 3.8390e+02, 7.4496e+04, 1.3316e+03],
[-6.8788e-05, 0.0000e+00, 0.0000e+00, 7.3361e+03, 6.8571e+03,
2.8003e+00, 6.9009e+03, -2.1069e+01, 0.0000e+00],
[-1.1089e-07, 0.0000e+00, 0.0000e+00, 9.5255e-01, -1.3963e+00,
-7.0740e+00, -1.4048e+00, 1.0646e+01, -4.0020e-01],
[-6.0409e-08, 0.0000e+00, 0.0000e+00, 4.3033e+00, 7.3669e+00,
-1.3401e+00, 7.3896e+00, 2.4513e+00, 2.0169e+00],
[-4.2006e-07, 0.0000e+00, 0.0000e+00, 6.0801e+00, -1.0565e-02,
2.0517e+00, -2.5089e-02, 3.8339e+01, 0.0000e+00]],
[[-1.7931e-04, 0.0000e+00, 0.0000e+00, 9.6886e+02, 1.8577e+03,
-3.2898e+02, 1.8722e+03, -1.3666e+04, 6.9803e+03],
[-7.9118e-04, 0.0000e+00, 0.0000e+00, 6.5534e+03, 7.1299e+02,
-6.5641e+01, 7.0230e+02, 7.4755e+04, 1.4173e+03],
[-1.1303e-03, 0.0000e+00, 0.0000e+00, 7.5908e+03, 6.8496e+03,
9.1323e+01, 6.8927e+03, 1.6289e+03, 0.0000e+00],
[-1.5126e-06, 0.0000e+00, 0.0000e+00, 7.4362e-01, -1.4862e+00,
-7.0655e+00, -1.4950e+00, 9.2198e+00, -7.3250e-01],
[-1.0418e-06, 0.0000e+00, 0.0000e+00, 4.8297e+00, 7.3415e+00,
-1.3388e+00, 7.3626e+00, 9.5159e+00, 1.9516e+00],
[-5.8711e-06, 0.0000e+00, 0.0000e+00, 5.7032e+00, -3.3712e-01,
2.0495e+00, -3.5366e-01, 3.8028e+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([[[-8.7316e-05, 0.0000e+00, 0.0000e+00, -1.3524e+03, 9.2991e+02,
6.4871e+01, 9.0539e+02, -3.2000e+04, -6.4358e+03],
[-4.2967e-04, 0.0000e+00, 0.0000e+00, 2.0608e+03, -5.2764e+02,
-3.2641e+01, -5.7373e+02, -6.7933e+04, 3.0472e+03],
[ 6.8895e-04, 0.0000e+00, 0.0000e+00, -1.3886e+04, 7.0421e+03,
-8.6782e+00, 7.0255e+03, 2.2650e+02, 0.0000e+00],
[-2.9977e-07, 0.0000e+00, 0.0000e+00, 3.0176e+00, -3.2110e+00,
6.6830e+00, -3.2003e+00, 4.5207e+00, 6.0321e-01],
[-2.6128e-07, 0.0000e+00, 0.0000e+00, 6.7715e+00, -6.7624e+00,
-3.1530e+00, -6.7590e+00, -4.2032e+00, 9.5151e-01],
[-1.0523e-06, 0.0000e+00, 0.0000e+00, -1.0514e+00, -1.3070e-01,
-1.1153e+00, -7.6140e-02, 3.9003e+01, 0.0000e+00]],
[[-3.2346e-04, 0.0000e+00, 0.0000e+00, -4.9164e+03, 2.8092e+02,
1.3914e+02, 2.7892e+02, -6.5456e+04, -3.4722e+03],
[-3.4292e-04, 0.0000e+00, 0.0000e+00, -1.0336e+03, -1.0957e+03,
-2.4919e+02, -1.0978e+03, -3.7372e+04, 6.2105e+03],
[ 8.9369e-04, 0.0000e+00, 0.0000e+00, -9.8360e+03, 7.0327e+03,
-4.1852e+01, 7.0379e+03, 1.3512e+03, 0.0000e+00],
[-7.5528e-07, 0.0000e+00, 0.0000e+00, 4.6897e+00, -6.5353e+00,
3.6431e+00, -6.5376e+00, -2.5933e+00, 1.1514e+00],
[ 4.4026e-07, 0.0000e+00, 0.0000e+00, 1.6631e+00, -3.6532e+00,
-6.4248e+00, -3.6550e+00, -8.3856e+00, 2.9214e-01],
[-4.8603e-06, 0.0000e+00, 0.0000e+00, 5.5956e+00, -3.0352e-01,
-1.1431e+00, -3.0100e-01, 3.8896e+01, 0.0000e+00]],
[[-2.1239e-04, 0.0000e+00, 0.0000e+00, 7.3083e+02, -1.6674e+03,
5.5532e+01, -1.6885e+03, 1.3230e+04, -6.8057e+03],
[ 1.0605e-03, 0.0000e+00, 0.0000e+00, -6.8753e+03, -3.8401e+02,
1.3660e+01, -3.8161e+02, -6.9669e+04, -1.3106e+03],
[ 4.0336e-05, 0.0000e+00, 0.0000e+00, 2.2327e+03, 6.7181e+03,
1.0668e+01, 6.7976e+03, 1.4416e+02, 0.0000e+00],
[ 1.0418e-06, 0.0000e+00, 0.0000e+00, -1.5560e+00, 1.4383e+00,
7.2389e+00, 1.4486e+00, -9.1725e+00, 4.1343e-01],
[ 2.0344e-07, 0.0000e+00, 0.0000e+00, -1.6907e+00, -7.4786e+00,
1.3807e+00, -7.5225e+00, -2.9885e+00, -1.8289e+00],
[-4.1679e-06, 0.0000e+00, 0.0000e+00, 7.2585e+00, -5.4554e-02,
1.8697e+00, -6.2305e-02, 3.7754e+01, 0.0000e+00]],
[[-2.7840e-03, 0.0000e+00, 0.0000e+00, 4.9165e+03, -4.3047e+02,
1.0660e+02, -4.0752e+02, 6.5071e+04, -2.7440e+03],
[-5.0181e-04, 0.0000e+00, 0.0000e+00, 6.7106e+02, -1.6818e+03,
2.4888e+02, -1.7054e+03, -2.9099e+04, -6.3789e+03],
[ 6.2292e-03, 0.0000e+00, 0.0000e+00, -1.0225e+04, 6.7311e+03,
6.4844e+01, 6.7796e+03, 1.2373e+03, 0.0000e+00],
[ 4.6158e-06, 0.0000e+00, 0.0000e+00, -5.5110e+00, 6.9761e+00,
2.9668e+00, 7.0001e+00, 1.9608e+00, 1.8462e+00],
[ 2.7043e-06, 0.0000e+00, 0.0000e+00, 1.0064e+00, -2.9925e+00,
6.7189e+00, -3.0111e+00, -1.1250e+01, -4.4068e-01],
[-1.7043e-05, 0.0000e+00, 0.0000e+00, 5.1458e+00, -3.2585e-01,
1.8627e+00, -2.9702e-01, 3.7486e+01, 0.0000e+00]],
[[-3.5487e-03, 0.0000e+00, 0.0000e+00, 4.9023e+03, -1.1916e+03,
2.1064e+02, -1.1911e+03, 4.2989e+04, -5.4830e+03],
[ 2.6808e-03, 0.0000e+00, 0.0000e+00, -3.4938e+03, -1.2606e+03,
1.6404e+02, -1.2959e+03, -5.7047e+04, -4.2541e+03],
[ 4.7261e-03, 0.0000e+00, 0.0000e+00, -7.0262e+03, 6.7269e+03,
6.3900e+01, 6.8094e+03, 1.2217e+03, 0.0000e+00],
[ 4.8903e-06, 0.0000e+00, 0.0000e+00, -3.4741e+00, 4.6660e+00,
5.8544e+00, 4.6904e+00, -3.8103e+00, 1.3992e+00],
[ 9.7338e-07, 0.0000e+00, 0.0000e+00, 1.8321e+00, -6.0027e+00,
4.4545e+00, -6.0451e+00, -1.0754e+01, -1.2854e+00],
[-1.7605e-05, 0.0000e+00, 0.0000e+00, 6.5198e+00, -3.1864e-01,
1.8660e+00, -2.9425e-01, 3.7656e+01, 0.0000e+00]],
[[-4.7826e-06, 0.0000e+00, 0.0000e+00, -6.1081e+02, 7.6790e+02,
6.1023e+03, 7.6500e+02, 7.1789e+03, 3.9502e+03],
[ 1.9071e-05, 0.0000e+00, 0.0000e+00, -4.1181e+03, -6.1426e+03,
2.5632e+02, -6.1197e+03, 4.1279e+04, -6.3297e+02],
[ 2.8167e-05, 0.0000e+00, 0.0000e+00, 5.9887e+03, -3.8877e+03,
7.9891e+02, -3.8729e+03, -6.8181e+04, 0.0000e+00],
[-1.1293e-07, 0.0000e+00, 0.0000e+00, -7.6651e-01, 6.4048e-01,
-3.9205e+00, 6.3926e-01, 4.4358e+00, 6.2039e+00],
[ 8.9963e-07, 0.0000e+00, 0.0000e+00, 6.2636e+00, 3.9906e+00,
-1.6382e-01, 3.9833e+00, -3.2703e+01, 7.7490e-01],
[ 5.7219e-07, 0.0000e+00, 0.0000e+00, 3.8428e+00, -6.1819e+00,
-5.1138e-01, -6.1702e+00, -2.3440e+01, 0.0000e+00]]])