This repository has been archived by the owner on Oct 14, 2023. It is now read-only.
-
-
Notifications
You must be signed in to change notification settings - Fork 282
/
fixed.py
374 lines (306 loc) · 11.1 KB
/
fixed.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
import numpy as np
from astropy import units as u
from astropy.coordinates import (
HCRS,
ITRS,
BaseRADecFrame,
FunctionTransform,
TimeAttribute,
frame_transform_graph,
)
from astropy.coordinates.builtin_frames.utils import DEFAULT_OBSTIME
from astropy.coordinates.matrix_utilities import rotation_matrix
from poliastro.bodies import (
Jupiter,
Mars,
Mercury,
Moon,
Neptune,
Saturn,
Sun,
Uranus,
Venus,
)
from poliastro.constants import J2000
from .equatorial import (
JupiterICRS,
MarsICRS,
MercuryICRS,
MoonICRS,
NeptuneICRS,
SaturnICRS,
UranusICRS,
VenusICRS,
)
__all__ = [
"SunFixed",
"MercuryFixed",
"VenusFixed",
"ITRS",
"MarsFixed",
"JupiterFixed",
"SaturnFixed",
"UranusFixed",
"NeptuneFixed",
]
class _PlanetaryFixed(BaseRADecFrame):
obstime = TimeAttribute(default=DEFAULT_OBSTIME)
def __new__(cls, *args, **kwargs):
frame_transform_graph.transform(FunctionTransform, cls, cls.equatorial)(
cls.to_equatorial
)
frame_transform_graph.transform(FunctionTransform, cls.equatorial, cls)(
cls.from_equatorial
)
return super().__new__(cls)
@staticmethod
def to_equatorial(fixed_coo, equatorial_frame):
# TODO replace w/ something smart (Sun/Earth special cased)
if fixed_coo.body == Sun:
assert type(equatorial_frame) == HCRS
else:
assert fixed_coo.body == equatorial_frame.body
r = fixed_coo.cartesian
ra, dec, W = fixed_coo.rot_elements_at_epoch(equatorial_frame.obstime)
r = r.transform(rotation_matrix(-W, "z"))
r_trans1 = r.transform(rotation_matrix(-(90 * u.deg - dec), "x"))
data = r_trans1.transform(rotation_matrix(-(90 * u.deg + ra), "z"))
return equatorial_frame.realize_frame(data)
@staticmethod
def from_equatorial(equatorial_coo, fixed_frame):
# TODO replace w/ something smart (Sun/Earth special cased)
if fixed_frame.body == Sun:
assert type(equatorial_coo) == HCRS
else:
assert equatorial_coo.body == fixed_frame.body
r = equatorial_coo.cartesian
ra, dec, W = fixed_frame.rot_elements_at_epoch(fixed_frame.obstime)
r_trans2 = r.transform(rotation_matrix(90 * u.deg + ra, "z"))
r_f = r_trans2.transform(rotation_matrix(90 * u.deg - dec, "x"))
r_f = r_f.transform(rotation_matrix(W, "z"))
return fixed_frame.realize_frame(r_f)
@classmethod
def rot_elements_at_epoch(cls, epoch):
"""Provides rotational elements at epoch.
Provides north pole of body and angle to prime meridian.
Parameters
----------
epoch : ~astropy.time.Time, optional
Epoch, default to J2000.
Returns
-------
ra, dec, W: tuple (~astropy.units.Quantity)
Right ascension and declination of north pole, and angle of the prime meridian.
"""
T = (epoch.tdb - J2000).to(u.day).value / 36525
d = (epoch.tdb - J2000).to(u.day).value
return cls._rot_elements_at_epoch(T, d)
@staticmethod
def _rot_elements_at_epoch(T, d):
raise NotImplementedError
class SunFixed(_PlanetaryFixed):
body = Sun
equatorial = HCRS
@staticmethod
def _rot_elements_at_epoch(T, d):
ra = 286.13 * u.deg
dec = 63.87 * u.deg
W = (84.176 + 14.1844000 * d) * u.deg
return ra, dec, W
class MercuryFixed(_PlanetaryFixed):
body = Mercury
equatorial = MercuryICRS
@staticmethod
def _rot_elements_at_epoch(T, d):
M1 = (174.7910857 + 4.092335 * d) * u.deg
M2 = (349.5821714 + 8.184670 * d) * u.deg
M3 = (164.3732571 + 12.277005 * d) * u.deg
M4 = (339.1643429 + 16.369340 * d) * u.deg
M5 = (153.9554286 + 20.461675 * d) * u.deg
ra = (281.0103 - 0.0328 * T) * u.deg
dec = (61.45 - 0.005 * T) * u.deg
W = (329.5988 + 6.1385108 * d) * u.deg + (
0.01067257 * np.sin(M1.to("rad").value)
- 0.00112309 * np.sin(M2.to("rad").value)
- 0.00011040 * np.sin(M3.to("rad").value)
- 0.00002539 * np.sin(M4.to("rad").value)
- 0.00000571 * np.sin(M5.to("rad").value)
) * u.deg
return ra, dec, W
class VenusFixed(_PlanetaryFixed):
body = Venus
equatorial = VenusICRS
@staticmethod
def _rot_elements_at_epoch(T, d):
ra = 272.76 * u.deg
dec = 67.16 * u.deg
W = (160.20 - 1.4813688 * d) * u.deg
return ra, dec, W
class MarsFixed(_PlanetaryFixed):
body = Mars
equatorial = MarsICRS
@staticmethod
def _rot_elements_at_epoch(T, d):
M1 = (198.991226 + 19139.4819985 * T) * u.deg
M2 = (226.292679 + 38280.8511281 * T) * u.deg
M3 = (249.663391 + 57420.7251593 * T) * u.deg
M4 = (266.183510 + 76560.6367950 * T) * u.deg
M5 = (79.398797 + 0.5042615 * T) * u.deg
ra = (
317.269202
- 0.10927547 * T
+ 0.000068 * np.sin(M1.to("rad").value)
+ 0.000238 * np.sin(M2.to("rad").value)
+ 0.000052 * np.sin(M3.to("rad").value)
+ 0.000009 * np.sin(M4.to("rad").value)
+ 0.419057 * np.sin(M5.to("rad").value)
) * u.deg
K1 = (122.433576 + 19139.9407476 * T) * u.deg
K2 = (43.058401 + 38280.8753272 * T) * u.deg
K3 = (57.663379 + 57420.7517205 * T) * u.deg
K4 = (79.476401 + 76560.6495004 * T) * u.deg
K5 = (166.325722 + 0.5042615 * T) * u.deg
dec = (
54.432516
- 0.05827105 * T
+ 0.000051 * np.cos(K1.to("rad").value)
+ 0.000141 * np.cos(K2.to("rad").value)
+ 0.000031 * np.cos(K3.to("rad").value)
+ 0.000005 * np.cos(K4.to("rad").value)
+ 1.591274 * np.cos(K5.to("rad").value)
) * u.deg
J1 = (129.071773 + 19140.0328244 * T) * u.deg
J2 = (36.352167 + 38281.0473591 * T) * u.deg
J3 = (56.668646 + 57420.9295360 * T) * u.deg
J4 = (67.364003 + 76560.2552215 * T) * u.deg
J5 = (104.792680 + 95700.4387578 * T) * u.deg
J6 = (95.391654 + 0.5042615 * T) * u.deg
W = (
176.049863
+ 350.891982443297 * d
+ 0.000145 * np.sin(J1.to("rad").value)
+ 0.000157 * np.sin(J2.to("rad").value)
+ 0.000040 * np.sin(J3.to("rad").value)
+ 0.000001 * np.sin(J4.to("rad").value)
+ 0.000001 * np.sin(J5.to("rad").value)
+ 0.584542 * np.sin(J6.to("rad").value)
) * u.deg
return ra, dec, W
class JupiterFixed(_PlanetaryFixed):
body = Jupiter
equatorial = JupiterICRS
@staticmethod
def _rot_elements_at_epoch(T, d):
Ja = (99.360714 + 4850.4046 * T) * u.deg
Jb = (175.895369 + 1191.9605 * T) * u.deg
Jc = (300.323162 + 262.5475 * T) * u.deg
Jd = (114.012305 + 6070.2476 * T) * u.deg
Je = (49.511251 + 64.3000 * T) * u.deg
ra = (
268.056595
- 0.006499 * T
+ 0.000117 * np.sin(Ja.to("rad").value)
+ 0.000938 * np.sin(Jb.to("rad").value)
+ 0.001432 * np.sin(Jc.to("rad").value)
+ 0.000030 * np.sin(Jd.to("rad").value)
+ 0.002150 * np.sin(Je.to("rad").value)
) * u.deg
dec = (
64.495303
+ 0.002413 * T
+ 0.000050 * np.cos(Ja.to("rad").value)
+ 0.000404 * np.cos(Jb.to("rad").value)
+ 0.000617 * np.cos(Jc.to("rad").value)
- 0.000013 * np.cos(Jd.to("rad").value)
+ 0.000926 * np.cos(Je.to("rad").value)
) * u.deg
W = (284.95 + 870.536 * d) * u.deg
return ra, dec, W
class SaturnFixed(_PlanetaryFixed):
body = Saturn
equatorial = SaturnICRS
@staticmethod
def _rot_elements_at_epoch(T, d):
ra = (40.589 - 0.036 * T) * u.deg
dec = (83.537 - 0.004 * T) * u.deg
W = (38.90 + 810.7939024 * d) * u.deg
return ra, dec, W
class UranusFixed(_PlanetaryFixed):
body = Uranus
equatorial = UranusICRS
@staticmethod
def _rot_elements_at_epoch(T, d):
ra = 257.311 * u.deg
dec = -15.175 * u.deg
W = (203.81 - 501.1600928 * d) * u.deg
return ra, dec, W
class NeptuneFixed(_PlanetaryFixed):
body = Neptune
equatorial = NeptuneICRS
@staticmethod
def _rot_elements_at_epoch(T, d):
N = (357.85 + 52.316 * T) * u.deg
ra = (299.36 + 0.70 * np.sin(N.to("rad").value)) * u.deg
dec = (43.46 - 0.51 * np.cos(N.to("rad").value)) * u.deg
W = (249.978 + 541.1397757 * d - 0.48 * np.sin(N.to("rad").value)) * u.deg
return ra, dec, W
class MoonFixed(_PlanetaryFixed):
body = Moon
equatorial = MoonICRS
@staticmethod
def _rot_elements_at_epoch(T, d):
E1 = (125.045 - 0.0529921 * d) * u.deg
E2 = (250.089 - 0.1059842 * d) * u.deg
E3 = (260.008 + 13.0120009 * d) * u.deg
E4 = (176.625 + 13.3407154 * d) * u.deg
E5 = (357.529 + 0.9856003 * d) * u.deg
E6 = (311.589 + 26.4057084 * d) * u.deg
E7 = (134.963 + 13.0649930 * d) * u.deg
E8 = (276.617 + 0.3287146 * d) * u.deg
E9 = (34.226 + 1.7484877 * d) * u.deg
E10 = (15.134 - 0.1589763 * d) * u.deg
E11 = (119.743 + 0.0036096 * d) * u.deg
E12 = (239.961 + 0.1643573 * d) * u.deg
E13 = (25.053 + 12.9590088 * d) * u.deg
ra = (
269.9949
+ 0.0031 * T
- 3.8787 * np.sin(E1.to("rad").value)
- 0.1204 * np.sin(E2.to("rad").value)
+ 0.0700 * np.sin(E3.to("rad").value)
- 0.0172 * np.sin(E4.to("rad").value)
+ 0.0072 * np.sin(E6.to("rad").value)
- 0.0052 * np.sin(E10.to("rad").value)
+ 0.0043 * np.sin(E13.to("rad").value)
) * u.deg
dec = (
66.5392
+ 0.0130 * T
+ 1.5419 * np.cos(E1.to("rad").value)
+ 0.0239 * np.cos(E2.to("rad").value)
- 0.0278 * np.cos(E3.to("rad").value)
+ 0.0068 * np.cos(E4.to("rad").value)
- 0.0029 * np.cos(E6.to("rad").value)
+ 0.0009 * np.cos(E7.to("rad").value)
+ 0.0008 * np.cos(E10.to("rad").value)
- 0.0009 * np.cos(E13.to("rad").value)
) * u.deg
W = (
38.321
+ 13.17635815 * d
- 1.4e-12 * d ** 2
+ 3.5610 * np.sin(E1.to("rad").value)
+ 0.1208 * np.sin(E2.to("rad").value)
- 0.0642 * np.sin(E3.to("rad").value)
+ 0.0158 * np.sin(E4.to("rad").value)
+ 0.0252 * np.sin(E5.to("rad").value)
- 0.0066 * np.sin(E6.to("rad").value)
- 0.0047 * np.sin(E7.to("rad").value)
- 0.0046 * np.sin(E8.to("rad").value)
+ 0.0028 * np.sin(E9.to("rad").value)
+ 0.0052 * np.sin(E10.to("rad").value)
+ 0.0040 * np.sin(E11.to("rad").value)
+ 0.0019 * np.sin(E12.to("rad").value)
- 0.0044 * np.sin(E13.to("rad").value)
) * u.deg
return ra, dec, W