-
Notifications
You must be signed in to change notification settings - Fork 73
/
orbital.py
1334 lines (1105 loc) · 46.8 KB
/
orbital.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
#!/usr/bin/env python
# -*- coding: utf-8 -*-
# Copyright (c) 2011, 2012, 2013.
# Author(s):
# Esben S. Nielsen <esn@dmi.dk>
# Adam Dybbroe <adam.dybbroe@smhi.se>
# This program is free software: you can redistribute it and/or modify
# it under the terms of the GNU General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
# This program is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU General Public License for more details.
# You should have received a copy of the GNU General Public License
# along with this program. If not, see <http://www.gnu.org/licenses/>.
"""Module for computing the orbital parameters of satellites.
"""
from datetime import datetime, timedelta
import numpy as np
from pyorbital import tlefile
from pyorbital import astronomy
ECC_EPS = 1.0e-6 # Too low for computing further drops.
ECC_LIMIT_LOW = -1.0e-3
ECC_LIMIT_HIGH = 1.0 - ECC_EPS # Too close to 1
ECC_ALL = 1.0e-4
EPS_COS = 1.5e-12
EPS_SIN = 1.0e-12
EPS_COSIO = 1.5e-12
NR_EPS = 1.0e-12
CK2 = 5.413080e-4
CK4 = 0.62098875e-6
E6A = 1.0e-6
QOMS2T = 1.88027916e-9
S = 1.01222928
S0 = 78.0
XJ3 = -0.253881e-5
#XKE = 0.743669161e-1
XKE = 7.43669161331734132e-2
XKMPER = 6378.135
XMNPDA = 1440.0
#MFACTOR = 7.292115E-5
AE = 1.0
SECDAY = 8.6400E4
F = 1 / 298.257223563 # Earth flattening WGS-84
A = 6378.137 # WGS84 Equatorial radius
SGDP4_ZERO_ECC = 0
SGDP4_NEAR_SIMP = 1
SGDP4_NEAR_NORM = 2
SGDP4_DEEP_NORM = 3
SGDP4_DEEP_RESN = 4
SGDP4_DEEP_SYNC = 5
KS = AE * (1.0 + S0 / XKMPER)
A3OVK2 = (-XJ3 / CK2) * AE ** 3
ZNS = 1.19459e-5
C1SS = 2.9864797e-6
ZES = 0.01675
ZNL = 1.5835218e-4
C1L = 4.7968065e-7
ZEL = 0.0549
ZCOSIS = 0.91744867
ZSINIS = 0.39785416
ZCOSGS = 0.1945905
ZSINGS = -0.98088458
ROOT22 = 1.7891679e-6
ROOT32 = 3.7393792e-7
ROOT44 = 7.3636953e-9
ROOT52 = 1.1428639e-7
ROOT54 = 2.1765803e-9
G22 = 5.7686396
G32 = 0.95240898
G44 = 1.8014998
G52 = 1.0508330
G54 = 4.4108898
THDT = 4.37526908801129966e-3
STEP = 720.0
class OrbitalError(Exception):
pass
class Orbital(object):
"""Class for orbital computations.
The *satellite* parameter is the name of the satellite to work on and is
used to retreive the right TLE data for internet or from *tle_file* in case
it is provided.
"""
def __init__(self, satellite, tle_file=None, line1=None, line2=None):
satellite = satellite.upper()
self.satellite_name = satellite
self.tle = tlefile.read(satellite, tle_file=tle_file,
line1=line1, line2=line2)
self.orbit_elements = OrbitElements(self.tle)
self._sgdp4 = _SGDP4(self.orbit_elements)
pos_epoch, vel_epoch = self.get_position(self.tle.epoch,
normalize=False)
"""
if np.abs(pos_epoch[2]) > 1 or not vel_epoch[2] > 0:
# Epoch not at ascending node
self.orbit_elements.an_time = self.get_last_an_time(self.tle.epoch)
else:
# Epoch at ascending node (z < 1 km) and positive v_z
self.orbit_elements.an_time = self.tle.epoch
self.orbit_elements.an_period = self.orbit_elements.an_time - \
self.get_last_an_time(self.orbit_elements.an_time)"""
def __str__(self):
return self.satellite_name + " " + str(self.tle)
def get_last_an_time(self, utc_time):
"""Calculate time of last ascending node relative to the
specified time
"""
# Propagate backwards to ascending node
dt = timedelta(minutes=10)
t_old = utc_time
t_new = t_old - dt
pos0, vel0 = self.get_position(t_old, normalize=False)
pos1, vel1 = self.get_position(t_new, normalize=False)
while not (pos0[2] > 0 and pos1[2] < 0):
pos0, vel0 = pos1, vel1
t_old = t_new
t_new = t_old - dt
pos1, vel1 = self.get_position(t_new, normalize=False)
# Return if z within 1 km of an
if np.abs(pos0[2]) < 1:
return t_old
elif np.abs(pos1[2]) < 1:
return t_new
# Bisect to z within 1 km
while np.abs(pos1[2]) > 1:
pos0, vel0 = pos1, vel1
dt = (t_old - t_new) / 2
t_mid = t_old - dt
pos1, vel1 = self.get_position(t_mid, normalize=False)
if pos1[2] > 0:
t_old = t_mid
else:
t_new = t_mid
return t_mid
def get_position(self, utc_time, normalize=True):
"""Get the cartesian position and velocity from the satellite.
"""
kep = self._sgdp4.propagate(utc_time)
pos, vel = kep2xyz(kep)
if normalize:
pos /= XKMPER
vel /= XKMPER * XMNPDA / SECDAY
return pos, vel
def get_lonlatalt(self, utc_time):
"""Calculate sublon, sublat and altitude of satellite.
http://celestrak.com/columns/v02n03/
"""
(pos_x, pos_y, pos_z), (vel_x, vel_y, vel_z) = self.get_position(utc_time, normalize=True)
lon = ((np.arctan2(pos_y * XKMPER, pos_x * XKMPER) - astronomy.gmst(utc_time))
% (2 * np.pi))
lon = np.where(lon > np.pi, lon - np.pi * 2, lon)
lon = np.where(lon <= -np.pi, lon + np.pi *2, lon)
r = np.sqrt(pos_x ** 2 + pos_y ** 2)
lat = np.arctan2(pos_z, r)
e2 = F * (2 - F)
while True:
lat2 = lat
c = 1/(np.sqrt(1 - e2 * (np.sin(lat2) ** 2)))
lat = np.arctan2(pos_z + c * e2 *np.sin(lat2), r)
if np.all(abs(lat - lat2) < 1e-10):
break
alt = r / np.cos(lat)- c;
alt *= A
return np.rad2deg(lon), np.rad2deg(lat), alt
def find_aos(self, time, lon, lat):
pass
def find_aol(self, time, lon, lat):
pass
def get_observer_look(self, time, lon, lat, alt):
"""Calculate observers look angle to a satellite.
http://celestrak.com/columns/v02n02/
time: Observation time (datetime object)
lon: Longitude of observer position on ground
lat: Latitude of observer position on ground
alt: Altitude above sea-level (geoid) of observer position on ground
Return: (Azimuth, Elevation)
"""
(pos_x, pos_y, pos_z), (vel_x, vel_y, vel_z) = self.get_position(time, normalize=False)
(opos_x, opos_y, opos_z), (ovel_x, ovel_y, ovel_z) = \
astronomy.observer_position(time, lon, lat, alt)
lon = np.deg2rad(lon)
lat = np.deg2rad(lat)
theta = (astronomy.gmst(time) + lon) % (2 * np.pi)
rx = pos_x - opos_x
ry = pos_y - opos_y
rz = pos_z - opos_z
rvx = vel_x - ovel_x
rvy = vel_y - ovel_y
rvz = vel_z - ovel_z
sin_lat = np.sin(lat)
cos_lat = np.cos(lat)
sin_theta = np.sin(theta)
cos_theta = np.cos(theta)
top_s = sin_lat * cos_theta * rx + sin_lat * sin_theta * ry - cos_lat * rz
top_e = -sin_theta * rx + cos_theta * ry
top_z = cos_lat * cos_theta * rx + cos_lat * sin_theta * ry + sin_lat * rz
az = np.arctan(-top_e / top_s)
if top_s > 0:
az = az + np.pi
if az < 0:
az = az + 2 * np.pi
rg = np.sqrt(rx * rx + ry * ry + rz * rz)
el = np.arcsin(top_z / rg)
w = (rx * rvx + ry * rvy + rz * rvz) / rg
return np.rad2deg(az), np.rad2deg(el)
def get_orbit_number(self, utc_time, tbus_style=False):
"""Calculate orbit number at specified time.
Optionally use TBUS-style orbit numbering (TLE orbit number + 1)
"""
dt = astronomy._days(utc_time - self.orbit_elements.an_time)
orbit_period = astronomy._days(self.orbit_elements.an_period)
orbit = int(self.tle.orbit + dt / orbit_period +
self.tle.mean_motion_derivative * dt**2 +
self.tle.mean_motion_sec_derivative * dt**3)
if tbus_style:
orbit += 1
return orbit
def get_zenith_overpass(self, utc_time, obslon, obslat):
"""Get the time when the satellite is highest on the horizon relative
to the observer position on ground given by *obslon*,*obslat* closest
in time (in the future) to the time given by *utc_time*. Using the
method of gradient ascent with variable step parameter (decreasing in
size as we apprach 90 degrees zenith angle.
"""
# First check if the elevation is above zero. If not continue until it
# is, using larger steps when the ange is far from zero and shorter
# ones when the elevation gets closer to zero:
el_start = self.get_observer_look(utc_time, obslon, obslat, 0.0)[1]
if el_start < 0:
start_time = utc_time
elev = el_start
idx = 0
NIDX = 100
while elev < 0 and idx < NIDX:
var_scale = np.exp(np.square(np.sin(elev * np.pi/180.)))
t_step = timedelta(seconds = (300 * var_scale))
start_time = start_time + t_step
elev = self.get_observer_look(start_time, obslon, obslat, 0.0)[1]
idx = idx + 1
#print idx, start_time, var_scale, elev
utc_time = start_time - t_step
precision = timedelta(seconds=0.01)
sec_step = 0.5
t_step = timedelta(seconds=sec_step/2.0)
# Local derivative:
def fprime(timex):
el0 = self.get_observer_look(timex - t_step,
obslon, obslat, 0.0)[1]
el1 = self.get_observer_look(timex + t_step,
obslon, obslat, 0.0)[1]
return el0, (el1 - el0) / sec_step
tx0 = utc_time - timedelta(seconds=1.0)
tx1 = utc_time
idx = 0
NIDX = 100
eps = 1000. # Step size scaling
while abs(tx1 - tx0) > precision and idx < NIDX:
tx0 = tx1
fpr = fprime(tx0)
#var_scale = abs(90.0-fpr[0])
var_scale = np.abs(np.cos(fpr[0] * np.pi/180.))
tx1 = tx0 + timedelta(seconds = (eps * var_scale * fpr[1]))
#print idx, tx0, tx1, fpr
idx = idx + 1
if abs(tx1 - tx0) <= precision and idx < NIDX:
return tx1, idx
else:
return None
def get_risetime(self, utc_time, obslon, obslat, **kwargs):
"""Get the risetime of the satellite closest in time in the future to
the time *utc_time* for a reception station at *obslon*, *obslat*
position on the ground
"""
retv = self.get_zenith_overpass(utc_time, obslon, obslat)
if retv:
zenith_time = retv[0]
else:
zenith_time = utc_time
one_minute = timedelta(seconds = 60)
return self._get_time_at_horizon(zenith_time - one_minute,
obslon, obslat, **kwargs)
def get_falltime(self, utc_time, obslon, obslat, **kwargs):
"""Get the falltime of the satellite closest in time in the future to
the time *utc_time* for a reception station at *obslon*, *obslat*
position on the ground
"""
retv = self.get_zenith_overpass(utc_time, obslon, obslat)
if retv:
zenith_time = retv[0]
else:
zenith_time = utc_time
one_minute = timedelta(seconds = 60)
return self._get_time_at_horizon(zenith_time + one_minute,
obslon, obslat, **kwargs)
def _get_time_at_horizon(self, utc_time, obslon, obslat, **kwargs):
"""Get the time closest in time to *utc_time* when the
satellite is at the horizon relative to the position of an observer on
ground (altitude = 0)
"""
if "precision" in kwargs:
precision = kwargs['precision']
else:
precision = timedelta(seconds=0.001)
if "max_iterations" in kwargs:
nmax_iter = kwargs["max_iterations"]
else:
nmax_iter = 100
sec_step = 0.5
t_step = timedelta(seconds=sec_step/2.0)
# Local derivative:
def fprime(timex):
el0 = self.get_observer_look(timex - t_step,
obslon, obslat, 0.0)[1]
el1 = self.get_observer_look(timex + t_step,
obslon, obslat, 0.0)[1]
return el0, (abs(el1) - abs(el0)) / sec_step
tx0 = utc_time - timedelta(seconds=1.0)
tx1 = utc_time
idx = 0
#eps = 500.
eps = 100.
while abs(tx1 - tx0) > precision and idx < nmax_iter:
tx0 = tx1
fpr = fprime(tx0)
# When the elevation is high the scale is high, and when
# the elevation is low the scale is low
#var_scale = np.abs(np.sin(fpr[0] * np.pi/180.))
#var_scale = np.sqrt(var_scale)
var_scale = np.abs(fpr[0])
tx1 = tx0 - timedelta(seconds = (eps * var_scale * fpr[1]))
idx = idx + 1
#print idx, tx0, tx1, var_scale, fpr
if abs(tx1 - utc_time) < precision and idx < 2:
tx1 = tx1 + timedelta(seconds=1.0)
if abs(tx1 - tx0) <= precision and idx < nmax_iter:
return tx1
else:
return None
class OrbitElements(object):
"""Class holding the orbital elements.
"""
def __init__(self, tle):
self.epoch = tle.epoch
self.excentricity = tle.excentricity
self.inclination = np.deg2rad(tle.inclination)
self.right_ascension = np.deg2rad(tle.right_ascension)
self.arg_perigee = np.deg2rad(tle.arg_perigee)
self.mean_anomaly = np.deg2rad(tle.mean_anomaly)
self.mean_motion = tle.mean_motion * (np.pi * 2 / XMNPDA)
self.mean_motion_derivative = tle.mean_motion_derivative * np.pi * 2 / XMNPDA ** 2
self.mean_motion_sec_derivative = tle.mean_motion_sec_derivative * np.pi * 2 / XMNPDA ** 3
self.bstar = tle.bstar * AE
n_0 = self.mean_motion
k_e = XKE
k_2 = CK2
i_0 = self.inclination
e_0 = self.excentricity
a_1 = (k_e / n_0) ** (2.0/3)
delta_1 = ((3/2.0) * (k_2 / a_1**2) * ((3 * np.cos(i_0)**2 - 1) /
(1 - e_0**2)**(2.0/3)))
a_0 = a_1 * (1 - delta_1/3 - delta_1**2 - (134.0/81) * delta_1**3)
delta_0 = ((3/2.0) * (k_2 / a_0**2) * ((3 * np.cos(i_0)**2 - 1) /
(1 - e_0**2)**(2.0/3)))
# original mean motion
n_0pp = n_0 / (1 + delta_0)
self.original_mean_motion = n_0pp
# semi major axis
a_0pp = a_0 / (1 - delta_0)
self.semi_major_axis = a_0pp
self.period = np.pi * 2 / n_0pp
self.perigee = (a_0pp * (1 - e_0) / AE - AE) * XKMPER
self.right_ascension_lon = (self.right_ascension
- astronomy.gmst(self.epoch))
if self.right_ascension_lon > np.pi:
self.right_ascension_lon -= 2 * np.pi
class _SGDP4(object):
"""Class for the SGDP4 computations.
"""
def __init__(self, orbit_elements):
self.mode = None
perigee = orbit_elements.perigee
self.eo = orbit_elements.excentricity
print 'eo', self.eo
self.xincl = orbit_elements.inclination
self.xno = orbit_elements.original_mean_motion
k_2 = CK2
k_4 = CK4
k_e = XKE
self.bstar = orbit_elements.bstar
self.omegao = orbit_elements.arg_perigee
self.xmo = orbit_elements.mean_anomaly
self.xnodeo = orbit_elements.right_ascension
self.t_0 = orbit_elements.epoch
self.xn_0 = orbit_elements.mean_motion
A30 = -XJ3 * AE**3
if not(0 <= self.eo < ECC_LIMIT_HIGH):
raise OrbitalError('Eccentricity out of range: %e' % self.eo)
elif not((0.0035 * 2 * np.pi / XMNPDA) < self.xn_0 < (18 * 2 * np.pi / XMNPDA)):
raise OrbitalError('Mean motion out of range: %e' % self.xn_0)
elif not(0 < self.xincl < np.pi):
raise OrbitalError('Inclination out of range: %e' % self.xincl)
self.cosIO = np.cos(self.xincl)
self.sinIO = np.sin(self.xincl)
theta2 = self.cosIO**2
theta4 = theta2 ** 2
self.x3thm1 = 3.0 * theta2 - 1.0
self.x1mth2 = 1.0 - theta2
self.x7thm1 = 7.0 * theta2 - 1.0
a1 = (XKE / self.xn_0) ** (2. / 3)
betao2 = 1.0 - self.eo**2
betao = np.sqrt(betao2)
temp0 = 1.5 * CK2 * self.x3thm1 / (betao * betao2)
del1 = temp0 / (a1**2)
a0 = a1 * (1.0 - del1 * (1.0 / 3.0 + del1 * (1.0 + del1 * 134.0 / 81.0)))
del0 = temp0 / (a0**2)
self.xnodp = self.xn_0 / (1.0 + del0)
self.aodp = (a0 / (1.0 - del0))
self.perigee = (self.aodp * (1.0 - self.eo) - AE) * XKMPER
self.apogee = (self.aodp * (1.0 + self.eo) - AE) * XKMPER
self.period = (2 * np.pi * 1440.0 / XMNPDA) / self.xnodp
if self.eo == 0:
self.mode = self.SGDP4_ZERO_ECC
return
if self.period >= 225:
# Deep-Space model
self.mode = SGDP4_DEEP_NORM
elif self.perigee < 220:
# Near-space, simplified equations
self.mode = SGDP4_NEAR_SIMP
else:
# Near-space, normal equations
self.mode = SGDP4_NEAR_NORM
if self.perigee < 156:
s4 = self.perigee - 78
if s4 < 20:
s4 = 20
qoms24 = ((120 - s4) * (AE / XKMPER))**4
s4 = (s4 / XKMPER + AE)
else:
s4 = KS
qoms24 = QOMS2T
pinvsq = 1.0 / (self.aodp**2 * betao2**2)
tsi = 1.0 / (self.aodp - s4)
self.eta = self.aodp * self.eo * tsi
etasq = self.eta**2
eeta = self.eo * self.eta
psisq = np.abs(1.0 - etasq)
coef = qoms24 * tsi**4
coef_1 = coef / psisq**3.5
self.c2 = (coef_1 * self.xnodp * (self.aodp *
(1.0 + 1.5 * etasq + eeta * (4.0 + etasq)) +
(0.75 * CK2) * tsi / psisq * self.x3thm1 *
(8.0 + 3.0 * etasq * (8.0 + etasq))))
self.c1 = self.bstar * self.c2
self.c4 = (2.0 * self.xnodp * coef_1 * self.aodp * betao2 * (self.eta *
(2.0 + 0.5 * etasq) + self.eo * (0.5 + 2.0 *
etasq) - (2.0 * CK2) * tsi / (self.aodp * psisq) * (-3.0 *
self.x3thm1 * (1.0 - 2.0 * eeta + etasq *
(1.5 - 0.5 * eeta)) + 0.75 * self.x1mth2 * (2.0 *
etasq - eeta * (1.0 + etasq)) * np.cos(2.0 * self.omegao))))
self.c5, self.c3, self.omgcof = 0.0, 0.0, 0.0
if self.mode == SGDP4_NEAR_NORM:
self.c5 = (2.0 * coef_1 * self.aodp * betao2 *
(1.0 + 2.75 * (etasq + eeta) + eeta * etasq))
if self.eo > ECC_ALL:
self.c3 = coef * tsi * A3OVK2 * self.xnodp * AE * self.sinIO / self.eo
self.omgcof = self.bstar * self.c3 * np.cos(self.omegao)
temp1 = 3.0 * CK2 * pinvsq * self.xnodp
temp2 = temp1 * CK2 * pinvsq
temp3 = 1.25 * CK4 * pinvsq**2 * self.xnodp
self.xmdot = (self.xnodp + (0.5 * temp1 * betao * self.x3thm1 + 0.0625 *
temp2 * betao * (13.0 - 78.0 * theta2 +
137.0 * theta4)))
x1m5th = 1.0 - 5.0 * theta2
self.omgdot = (-0.5 * temp1 * x1m5th + 0.0625 * temp2 *
(7.0 - 114.0 * theta2 + 395.0 * theta4) +
temp3 * (3.0 - 36.0 * theta2 + 49.0 * theta4))
xhdot1 = -temp1 * self.cosIO
self.xnodot = (xhdot1 + (0.5 * temp2 * (4.0 - 19.0 * theta2) +
2.0 * temp3 * (3.0 - 7.0 * theta2)) * self.cosIO)
if self.eo > ECC_ALL:
self.xmcof = (-(2. / 3) * AE) * coef * self.bstar / eeta
else:
self.xmcof = 0.0
self.xnodcf = 3.5 * betao2 * xhdot1 * self.c1
self.t2cof = 1.5 * self.c1
# Check for possible divide-by-zero for X/(1+cos(xincl)) when calculating xlcof */
temp0 = 1.0 + self.cosIO
if np.abs(temp0) < EPS_COS:
temp0 = np.sign(temp0) * EPS_COS
self.xlcof = 0.125 * A3OVK2 * self.sinIO * (3.0 + 5.0 * self.cosIO) / temp0
self.aycof = 0.25 * A3OVK2 * self.sinIO
self.cosXMO = np.cos(self.xmo)
self.sinXMO = np.sin(self.xmo)
self.delmo = (1.0 + self.eta * self.cosXMO)**3
self._deep_space = None
if self.mode == SGDP4_NEAR_NORM:
c1sq = self.c1**2
self.d2 = 4.0 * self.aodp * tsi * c1sq
temp0 = self.d2 * tsi * self.c1 / 3.0
self.d3 = (17.0 * self.aodp + s4) * temp0
self.d4 = 0.5 * temp0 * self.aodp * tsi * (221.0 * self.aodp + 31.0 * s4) * self.c1
self.t3cof = self.d2 + 2.0 * c1sq
self.t4cof = 0.25 * (3.0 * self.d3 + self.c1 * (12.0 * self.d2 + 10.0 * c1sq))
self.t5cof = (0.2 * (3.0 * self.d4 + 12.0 * self.c1 * self.d3 + 6.0 * self.d2**2 +
15.0 * c1sq * (2.0 * self.d2 + c1sq)))
elif self.mode == SGDP4_DEEP_NORM:
self._deep_space = _DeepSpace(self.t_0, self.omegao, self.xnodeo,
self.xmo, self.eo, self.xincl,
self.aodp, self.xmdot,
self.omgdot, self.xnodot,
self.xnodp)
#raise NotImplementedError('Deep space calculations not supported')
def propagate(self, utc_time):
kep = {}
ts = astronomy._days(utc_time - self.t_0) * XMNPDA
print 'ts', ts
em = self.eo
xinc = self.xincl
xmp = self.xmo + self.xmdot * ts
xnode = self.xnodeo + ts * (self.xnodot + ts * self.xnodcf)
omega = self.omegao + self.omgdot * ts
if self.mode == SGDP4_ZERO_ECC:
raise NotImplementedError('Mode SGDP4_ZERO_ECC not implemented')
elif self.mode == SGDP4_NEAR_SIMP:
raise NotImplementedError('Mode "Near-space, simplified equations"'
' not implemented')
elif self.mode == SGDP4_NEAR_NORM:
delm = self.xmcof * ((1.0 + self.eta * np.cos(xmp))**3 - self.delmo)
temp0 = ts * self.omgcof + delm
xmp += temp0
omega -= temp0
tempa = 1.0 - (ts * (self.c1 + ts * (self.d2 + ts * (self.d3 + ts * self.d4))))
tempe = self.bstar * (self.c4 * ts + self.c5 * (np.sin(xmp) - self.sinXMO))
templ = ts * ts * (self.t2cof + ts * (self.t3cof + ts * (self.t4cof + ts * self.t5cof)))
a = self.aodp * tempa**2
e = em - tempe
xl = xmp + omega + xnode + self.xnodp * templ
x3thm1 = self.x3thm1
x1mth2 = self.x1mth2
x7thm1 = self.x7thm1
xlcof = self.xlcof
aycof = self.aycof
sinIO = self.sinIO
cosIO = self.cosIO
else:
#raise NotImplementedError('Deep space calculations not supported')
tempa = 1.0 - ts * self.c1
tempe = self.bstar * ts * self.c4
templ = ts * ts * self.t2cof
xn = self.xnodp
xmp, omega, xnode, em, xinc, xn = \
self._deep_space.dpsec(xmp, omega, xnode, em, xinc, xn, ts)
print 'xmp', xmp, 'omega', omega, 'xnode', xnode, 'em', em, 'xinc', xinc, 'xn', xn
a = (XKE / xn) ** (2. / 3.) * tempa ** 2
e = em - tempe
xmam = xmp + self.xnodp * templ
print 'xnodp', self.xnodp, 'templ', templ
print 'a', a, 'e', e, 'xmam', xmam
e, xinc, omega, xnode, xmam = \
self._deep_space.dpper(e, xinc, omega, xnode, xmam, ts)
print 'dpper'
print 'e', e, 'xinc', xinc, 'omega', omega, 'xnode', xnode, 'xmam', xmam
if xinc < 0:
xinc = -xinc
xnode += np.pi
omega -= np.pi
xl = xmam + omega + xnode
# Re-compute the perturbed values
sinIO = np.sin(xinc)
cosIO = np.cos(xinc)
theta2 = cosIO * cosIO
x3thm1 = 3.0 * theta2 - 1.0
x1mth2 = 1.0 - theta2
x7thm1 = 7.0 * theta2 - 1.0
# Check for possible divide-by-zero for X/(1+cosIO) when calculating xlcof
temp0 = 1.0 + cosIO
temp0 = np.where(np.abs(temp0) < EPS_COSIO, np.sign(temp0) * EPS_COSIO, temp0)
xlcof = 0.125 * A3OVK2 * sinIO * (3.0 + 5.0 * cosIO) / temp0
aycof = 0.25 * A3OVK2 * sinIO
if np.any(a < 1):
raise Exception('Satellite crashed at time %s', utc_time)
elif np.any(e < ECC_LIMIT_LOW):
raise ValueError('Satellite modified eccentricity to low: %e < %e'
% (e, ECC_LIMIT_LOW))
e = np.where(e < ECC_EPS, ECC_EPS, e)
e = np.where(e > ECC_LIMIT_HIGH, ECC_LIMIT_HIGH, e)
beta2 = 1.0 - e ** 2
# Long period periodics
sinOMG = np.sin(omega)
cosOMG = np.cos(omega)
temp0 = 1.0 / (a * beta2)
axn = e * cosOMG
ayn = e * sinOMG + temp0 * aycof
print 'axn', axn, 'ayn', ayn
xlt = xl + temp0 * xlcof * axn
print 'xlt', xlt, 'xnode', xnode
elsq = axn ** 2 + ayn ** 2
print 'elsq', elsq
if np.any(elsq >= 1):
raise Exception('e**2 >= 1 at %s', utc_time)
kep['ecc'] = np.sqrt(elsq)
print 'xlt - xnode', xlt - xnode
epw = np.fmod(xlt - xnode, 2 * np.pi)
print 'pre epw', epw
# needs a copy in case of an array
capu = np.array(epw)
maxnr = kep['ecc']
for i in range(10):
sinEPW = np.sin(epw)
cosEPW = np.cos(epw)
ecosE = axn * cosEPW + ayn * sinEPW
esinE = axn * sinEPW - ayn * cosEPW
f = capu - epw + esinE
if np.all(np.abs(f) < NR_EPS):
break
df = 1.0 - ecosE
# 1st order Newton-Raphson correction.
nr = f / df
# 2nd order Newton-Raphson correction.
nr = np.where(np.logical_and(i == 0, np.abs(nr) > 1.25 * maxnr),
np.sign(nr) * maxnr,
f / (df + 0.5*esinE*nr))
epw += nr
print 'post epw', epw
print 'sinEPW', sinEPW
print 'cosEPW', cosEPW
print 'ecosE', ecosE
print 'axn * cosEPW', axn * cosEPW
print 'ayn * sinEPW', ayn * sinEPW
# Short period preliminary quantities
temp0 = 1.0 - elsq
betal = np.sqrt(temp0)
pl = a * temp0
r = a * (1.0 - ecosE)
print 'r', r
invR = 1.0 / r
temp2 = a * invR
print 'temp2', temp2
temp3 = 1.0 / (1.0 + betal)
print 'temp3', temp3
cosu = temp2 * (cosEPW - axn + ayn * esinE * temp3)
sinu = temp2 * (sinEPW - ayn - axn * esinE * temp3)
print 'cosu', cosu, 'sinu', sinu
u = np.arctan2(sinu, cosu)
sin2u = 2.0 * sinu * cosu
cos2u = 2.0 * cosu**2 - 1.0
temp0 = 1.0 / pl
temp1 = CK2 * temp0
temp2 = temp1 * temp0
# Update for short term periodics to position terms.
rk = r * (1.0 - 1.5 * temp2 * betal * x3thm1) + 0.5 * temp1 * x1mth2 * cos2u
uk = u - 0.25 * temp2 * x7thm1 * sin2u
print 'u', u, 'temp2', temp2, 'x7thm1', x7thm1, 'sin2u', sin2u
xnodek = xnode + 1.5 * temp2 * cosIO * sin2u
xinck = xinc + 1.5 * temp2 * cosIO * sinIO * cos2u
if np.any(rk < 1):
raise Exception('Satellite crashed at time %s', utc_time)
temp0 = np.sqrt(a)
temp2 = XKE / (a * temp0)
rdotk = ((XKE * temp0 * esinE * invR -temp2 * temp1 * x1mth2 * sin2u) *
(XKMPER / AE * XMNPDA / 86400.0))
rfdotk = ((XKE * np.sqrt(pl) * invR + temp2 * temp1 *
(x1mth2 * cos2u + 1.5 * x3thm1)) *
(XKMPER / AE * XMNPDA / 86400.0))
kep['radius'] = rk * XKMPER / AE
kep['theta'] = uk
kep['eqinc'] = xinck
kep['ascn'] = xnodek
kep['argp'] = omega
kep['smjaxs'] = a * XKMPER / AE
kep['rdotk'] = rdotk
kep['rfdotk'] = rfdotk
print 'kep:', kep
return kep
class _DeepSpace(object):
def __init__(self, epoch, omegao, xnodeo, xmo, orb_eo, orb_xincl, aodp,
xlldot, omgdot, xnodot, xnodp):
def mod2pi(a):
b = np.fmod(a, 2 * np.pi)
if b < 0:
return b + 2 * np.pi
else:
return b
eo = eq = orb_eo
xincl = orb_xincl
# Decide on direct or Lyddane Lunar-Solar perturbations.
self.ilsd = False
if xincl >= 0.2:
self.ilsd = True
# Drop some terms below 3 deg inclination.
self.ishq = False
if xincl >= 0.052359877:
self.ishq = True
sinomo = np.sin(omegao)
cosomo = np.cos(omegao)
sinq = np.sin(xnodeo)
cosq = np.cos(xnodeo)
siniq = np.sin(xincl)
cosiq = np.cos(xincl)
if np.abs(siniq) <= EPS_SIN:
siniq = np.sign(siniq) * EPS_SIN
cosiq2 = cosiq ** 2
siniq2 = siniq ** 2
print cosiq2, siniq2
ao = aodp
self.omgdt = omgdot
eqsq = eo ** 2
bsq = 1 - eqsq
rteqsq = np.sqrt(bsq)
print 'rteqsq', rteqsq
self.thgr = astronomy.gmst(epoch)
print 'thgr', self.thgr
xnq = xnodp
aqnv = 1. / ao
xmao = xmo
xpidot = self.omgdt + xnodot
print 'xpidot', xpidot
self.omegaq = omegao
# Initialize lunar terms.
# Note: the Dundee reference code d50 is off by 1 day
jd1900 = astronomy.jdays1900(epoch)
jd1900 += 1 #TODO: compability?
print 'epoch', epoch
print 'd1900', jd1900
xnodce = 4.523602 - jd1900 * 9.2422029e-4
temp0 = np.fmod(xnodce, 2 * np.pi)
stem = np.sin(temp0)
ctem = np.cos(temp0)
zcosil = 0.91375164 - ctem * 0.03568096
zsinil = np.sqrt(1.0 - zcosil * zcosil)
zsinhl = stem * 0.089683511 / zsinil
zcoshl = np.sqrt(1.0 - zsinhl * zsinhl)
print 'zcoshl', zcoshl
c = jd1900 * 0.2299715 + 4.7199672
gam = jd1900 * 0.001944368 + 5.8351514
print 'c', c
print 'gam', gam
self.zmol = mod2pi(c - gam)
print 'zmol', self.zmol
zx = stem * 0.39785416 / zsinil
zy = zcoshl * ctem + zsinhl * 0.91744867 * stem
zx = np.arctan2(zx, zy)
zx = np.fmod(gam + zx - xnodce, 2 * np.pi)
print 'zx', zx
zsingl = np.sin(zx)
zcosgl = np.cos(zx)
self.zmos = mod2pi(jd1900 * 0.017201977 + 6.2565837)
print 'zsingl', zsingl
print 'zcosgl', zcosgl
print 'zmos', self.zmos
# Do solar terms
zcosg = ZCOSGS
zsing = ZSINGS
zcosi = 0.91744867
zsini = 0.39785416
zcosh = cosq
zsinh = sinq
cc = C1SS
zn = ZNS
ze = ZES
zmo = self.zmos
xnoi = 1. / xnq
for ls in range(2):
a1 = zcosg * zcosh + zsing * zcosi * zsinh
a3 = -zsing * zcosh + zcosg * zcosi * zsinh
a7 = -zcosg * zsinh + zsing * zcosi * zcosh
a8 = zsing * zsini
a9 = zsing * zsinh + zcosg * zcosi * zcosh
a10 = zcosg * zsini
a2 = cosiq * a7 + siniq * a8
a4 = cosiq * a9 + siniq * a10
a5 = -siniq * a7 + cosiq * a8
a6 = -siniq * a9 + cosiq * a10
x1 = a1 * cosomo + a2 * sinomo
x2 = a3 * cosomo + a4 * sinomo
x3 = -a1 * sinomo + a2 * cosomo
x4 = -a3 * sinomo + a4 * cosomo
x5 = a5 * sinomo
x6 = a6 * sinomo
x7 = a5 * cosomo
x8 = a6 * cosomo
z31 = x1 * 12.0 * x1 - x3 * 3.0 * x3
z32 = x1 * 24.0 * x2 - x3 * 6.0 * x4
z33 = x2 * 12.0 * x2 - x4 * 3.0 * x4
z1 = (a1 * a1 + a2 * a2) * 3.0 + z31 * eqsq
z2 = (a1 * a3 + a2 * a4) * 6.0 + z32 * eqsq
z3 = (a3 * a3 + a4 * a4) * 3.0 + z33 * eqsq
z11 = a1 * -6.0 * a5 + eqsq * (x1 * -24.0 * x7 - x3 * 6.0 * x5)
z12 = ((a1 * a6 + a3 * a5) * -6.0 + eqsq * ((x2 * x7
+ x1 * x8) * -24.0 - (x3 * x6 + x4 * x5) * 6.0))
z13 = a3 * -6.0 * a6 + eqsq * (x2 * -24.0 * x8 - x4 * 6.0 * x6)
z21 = a2 * 6.0 * a5 + eqsq * (x1 * 24.0 * x5 - x3 * 6.0 * x7)
z22 = ((a4 * a5 + a2 * a6) * 6.0 + eqsq * ((x2 * x5 + x1 * x6)
* 24.0 - (x4 * x7 + x3 * x8) * 6.0))
z23 = a4 * 6.0 * a6 + eqsq * (x2 * 24.0 * x6 - x4 * 6.0 * x8)
z1 = z1 + z1 + bsq * z31
z2 = z2 + z2 + bsq * z32
z3 = z3 + z3 + bsq * z33
s3 = cc * xnoi
s2 = s3 * -0.5 / rteqsq
s4 = s3 * rteqsq
s1 = eq * -15.0 * s4
s5 = x1 * x3 + x2 * x4
s6 = x2 * x3 + x1 * x4
s7 = x2 * x4 - x1 * x3
se = s1 * zn * s5
si = s2 * zn * (z11 + z13)
sl = -zn * s3 * (z1 + z3 - 14.0 - eqsq * 6.0)
sgh = s4 * zn * (z31 + z33 - 6.0)
print 'sgh', sgh
shdq = 0
if self.ishq:
sh = -zn * s2 * (z21 + z23)
shdq = sh / siniq
self.ee2 = s1 * 2.0 * s6
self.e3 = s1 * 2.0 * s7
self.xi2 = s2 * 2.0 * z12
self.xi3 = s2 * 2.0 * (z13 - z11)