-
Notifications
You must be signed in to change notification settings - Fork 14
/
Copy pathtess_stars2px.py
2030 lines (1919 loc) · 73.5 KB
/
tess_stars2px.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
"""
tess_stars2px.py - High precision TESS pointing tool.
Convert target coordinates given in Right Ascension and Declination to
TESS detector pixel coordinates for the prime mission TESS observing
sectors (Year 1 & 2), Extendend mission Years 3-5.
Can also query MAST to obtain detector
pixel coordinates for a star by TIC ID or common star name (must be online for this option).
USAGE to display command line arguments:
python tess_stars2px.py -h
AUTHORS: Original programming in C and focal plane geometry solutions
by Alan Levine (MIT)
This python translation by Christopher J. Burke (MIT)
Testing and focal plane geometry refinements by Michael Fausnaugh &
Roland Vanderspek (MIT)
Testing by Thomas Barclay (NASA Goddard) &
Jessica Roberts (Univ. of Colorado)
Sesame queries by Brett Morris (UW)
Proxy Support added by Dishendra Mishra
Deprecation warnings correctsion by Ethan Kruse
Updates by Tyler Pritchard, Christina Hedges
VERSION: 0.9.1
WHAT'S NEW:
-Year 8 pointings for Sectors 97-107 now available
-Year 7 pointings for Sectors 84-96 now available
-Year 6 pointings for Sectors 70-83 now available
-Deprecation Warning Corrections
-Year 5 pointings (Sectors 56-69) now available
-Added Sector Pointing override file input
Supports mission planning as well as enabling the user
to speed up the calculation by only searching a subset of all sectors
-Bug correction for aberration. Only impacts if you were using
aberration flag (-a) WITHOUT the single sector flag (-s). In other words,
does not affect users that did not use --aberrate or -aberrate with -s
NOTES:
-Pointing table is for TESS Year 1 - 5 (Sectors 1-69)
-Pointing table is unofficial, and the pointings may change.
-See https://tess.mit.edu/observations/ for latest TESS pointing table
-Pointing prediction algorithm is similar to internally at MIT for
target management. However, hard coded focal plane geometry is not
up to date and may contain inaccurate results.
-Testing shows pointing with this tool should be accurate to better than
a pixel, but without including aberration effects, ones algorithm
adopted for centroiding highly assymmetric point-spread function
at edge of
camera, and by-eye source location, a 2 pixel accuracy estimate is
warranted. There is an approximate aberration option now available
-The output pixel coordinates assume the ds9 convention with
1,1 being the middle of the lower left corner pixel.
-No corrections for velocity aberration are calculated by default.
Potentially more accurate
results can be obtained if the target RA and Declination coordinates
have aberration effects applied. An aberration approximation is available
by using the -a flag. The aberration approximation assumes Earth motion without
taking into account the TESS spacecraft motion around Earth.
-For proposals to the TESS science office or directors discretionary time,
please consult the TESS prediction webtool available at
https://heasarc.gsfc.nasa.gov/cgi-bin/tess/webtess/wtv.py
for official identification of 'observable' targets. However,
if your proposal depends on a single or few targets, then this tool is
helpful to further refine the likelihood of the target being available
on the detectors.
-The calibrated FFI fits file release at MAST and calibrated by
NASA Ames SPOC will have WCS information available to
supplant this code. The WCS generation is independent of the
focal plane geometry model employed in this code and will give
different results at the pixel level. However, the WCS information
is not available until the FFI files are released, whereas
this code can predict positions in advance of data release.
-Hard coded focal plane geometry parameters from rfpg5_c1kb.txt
TODOS:
-Time dependent Focal plane geometry
DEPENDENCIES:
python 3+
astropy
numpy
SPECIAL THANKS TO:
Includes code from the python MAST query examples
https://mast.stsci.edu/api/v0/pyex.html
IMPLEMENTATION DETAILS:
In summary, the code begins with a space craft bore site pointing in RA,
Dec, and roll angle. A series of Euler angle translation matrices
are calculated based upon the space craft bore site. Next the target
coordinates in RA and Dec are translated to the space craft bore site
frame. Next, the target coordinates are translated to each of the four
TESS camera frames. Once target coordinates are translated to the
camera frame the radial position of the target relative to the camera
center is checked to see if it is potentially in the camera field of view.
If so, the focal plane position is calculated using a radial polynomial
model with a constant term and terms the even powers (2nd ,4th , and 8th).
Rotations are applied to convert the on sky positions to the detector
readout directions.
MIT License
Copyright (c) 2018 Christopher J Burke
Permission is hereby granted, free of charge, to any person obtaining a copy
of this software and associated documentation files (the "Software"), to deal
in the Software without restriction, including without limitation the rights
to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
copies of the Software, and to permit persons to whom the Software is
furnished to do so, subject to the following conditions:
The above copyright notice and this permission notice shall be included in all
copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
SOFTWARE.
"""
import numpy as np
import os
import argparse
from astropy.coordinates import SkyCoord
import astropy.units as u
from astropy.time import Time
import sys
import datetime
import json
try: # Python 3.x
from urllib.parse import quote as urlencode
from urllib.request import urlretrieve
from urllib.parse import urlparse
except ImportError: # Python 2.x
from urllib import pathname2url as urlencode
from urllib import urlretrieve
from urlparse import urlparse
try: # Python 3.x
import http.client as httplib
except ImportError: # Python 2.x
import httplib
import scipy.optimize as opt
import base64
max_sector = 107
class Levine_FPG:
"""Al Levine Focal Plane Geometry Methods
Translated from starspx6.c
INPUT:
sc_ra_dec_roll = numpy array of the SpaceCraft boresite (sc Z-axis)
ra, dec, and roll [deg]
The roll angle is in RA, Dec space clockwise relative to the celestial
pole. roll angle = 0 [deg] implies space craft X-axis points N celestial (increasing dec)
roll angle = 90 [deg] implies sc X-axis points towards increasing/decreasing (?) RA
*** In practice there is a separate fpg file for each of the four cameras ***
rmat1[3,3] = is the rotation matrix from ra&dec to spacecraft boresite coords
rmat4[NCAM,3,3] - is the rotation matrix from ra&dec to NCAM coords
"""
parm_dict_list = [{}, {}, {}, {}]
NCAM = 4 # Number of Cameras
NCCD = 4 # Number of CCDs per Camera
def __init__(self, sc_ra_dec_roll=None, fpg_file_list=None):
self.eulcam = np.zeros((self.NCAM, 3), dtype=np.double)
self.optcon = np.zeros((self.NCAM, 6), dtype=np.double)
self.ccdxy0 = np.zeros((self.NCAM, self.NCCD, 2), dtype=np.double)
self.pixsz = np.zeros((self.NCAM, self.NCCD, 2), dtype=np.double)
self.ccdang = np.zeros((self.NCAM, self.NCCD), dtype=np.double)
self.ccdtilt = np.zeros((self.NCAM, self.NCCD, 2), dtype=np.double)
self.asymang = np.zeros((self.NCAM,), dtype=np.double)
self.asymfac = np.zeros((self.NCAM,), dtype=np.double)
self.rmat1 = np.zeros((3, 3), dtype=np.double)
self.rmat4 = np.zeros((self.NCAM, 3, 3), dtype=np.double)
self.havePointing = False
# Read in the fpg parameter files
self.read_all_levine_fpg_files(fpg_file_list)
# Generate rotation matrices if ra dec and roll values given
if not sc_ra_dec_roll is None:
# go from sky to spacecraft
self.sky_to_sc_mat(sc_ra_dec_roll)
# Go from spacecraft to each camera's coords
for icam in range(self.NCAM):
cureul = self.eulcam[icam, :]
rmat2 = self.sc_to_cam_mat(cureul)
self.rmat4[icam] = np.matmul(rmat2, self.rmat1)
self.havePointing = True
def read_all_levine_fpg_files(self, fpg_file_list=None):
default_fpg_file_list = [
"fpg_pars.txt-",
"fpg_pars.txt-",
"fpg_pars.txt-",
"fpg_pars.txt-",
]
# For each camera read in the separate fpg parameter file
for icam in range(self.NCAM):
if fpg_file_list == None:
fpg_file = default_fpg_file_list[icam]
else:
fpg_file = fpg_file_list[icam]
self.read_levine_fpg_file(icam, fpg_file)
# We now have parameters for all 4 cameras in the parm_dict_list
# parse the dictionary values into the working numpy arrays
for icam in range(self.NCAM):
pd = self.parm_dict_list[icam]
self.eulcam[icam][0] = pd["ang1_cam1"]
self.eulcam[icam][1] = pd["ang2_cam1"]
self.eulcam[icam][2] = pd["ang3_cam1"]
self.optcon[icam][0] = pd["fl_cam1"]
self.optcon[icam][1] = pd["opt_coef1_cam1"]
self.optcon[icam][2] = pd["opt_coef2_cam1"]
self.optcon[icam][3] = pd["opt_coef3_cam1"]
self.optcon[icam][4] = pd["opt_coef4_cam1"]
self.optcon[icam][5] = pd["opt_coef5_cam1"]
self.asymang[icam] = pd["asymang_cam1"]
self.asymfac[icam] = pd["asymfac_cam1"]
for iccd in range(self.NCCD):
self.ccdxy0[icam][iccd][0] = pd["x0_ccd{0:1d}_cam1".format(iccd + 1)]
self.ccdxy0[icam][iccd][1] = pd["y0_ccd{0:1d}_cam1".format(iccd + 1)]
self.pixsz[icam][iccd][0] = pd["pix_x_ccd{0:1d}_cam1".format(iccd + 1)]
self.pixsz[icam][iccd][1] = pd["pix_y_ccd{0:1d}_cam1".format(iccd + 1)]
self.ccdang[icam][iccd] = pd["ang_ccd{0:1d}_cam1".format(iccd + 1)]
self.ccdtilt[icam][iccd][0] = pd[
"tilt_x_ccd{0:1d}_cam1".format(iccd + 1)
]
self.ccdtilt[icam][iccd][1] = pd[
"tilt_y_ccd{0:1d}_cam1".format(iccd + 1)
]
def read_levine_fpg_file(self, icam, fpg_file):
gotParm = False
parm_dict = {}
if os.path.isfile(fpg_file):
try:
fpin = open(fpg_file, "r")
# Read in parameters
dtypeseq = ["U20", "i4", "f16"]
dataBlock = np.genfromtxt(fpin, dtype=dtypeseq)
parm_keys = dataBlock["f0"]
parm_fitted_flags = dataBlock["f1"]
parm_values = dataBlock["f2"]
# Now build dictionary of the parameters
for i in range(len(parm_keys)):
parm_dict[parm_keys[i]] = parm_values[i]
self.parm_dict_list[icam] = parm_dict
gotParm = True
print("Successful Focal Plane Geometry Read From {0}".format(fpg_file))
except:
print(
"Could not open {0}! Using Hard-coded Focal Plane Geometry from Levine_FPG read_levine_fpg_file()".format(
fpg_file
)
)
# If anything goes wrong with reading in parameters revert to hard coded version
# or file was never given and default_fpg_file does not exist
if not gotParm:
# print('Using Hard-coded Focal Plane Geometry from Levine_FPG read_levine_fpg_file')
# *** For now this hard code is just a filler need to actually fill in values for all cameras separately
# to prepare parameters for dictionary
# awk -v q="'" -v d=":" '{print q $1 q d $3 ",\"}' rfpg5_c1kb.txt
if icam == 0:
parm_dict = {
"ang1_cam1": 0.101588,
"ang2_cam1": -36.022035,
"ang3_cam1": 90.048315,
"fl_cam1": 145.948116,
"opt_coef1_cam1": 1.00000140,
"opt_coef2_cam1": 0.24779006,
"opt_coef3_cam1": -0.22681254,
"opt_coef4_cam1": 10.78243356,
"opt_coef5_cam1": -34.97817276,
"asymang_cam1": 0.00000000,
"asymfac_cam1": 1.00000000,
"x0_ccd1_cam1": 31.573417,
"y0_ccd1_cam1": 31.551637,
"pix_x_ccd1_cam1": 0.015000,
"pix_y_ccd1_cam1": 0.015000,
"ang_ccd1_cam1": 179.980833,
"tilt_x_ccd1_cam1": 0.000000,
"tilt_y_ccd1_cam1": 0.000000,
"x0_ccd2_cam1": -0.906060,
"y0_ccd2_cam1": 31.536148,
"pix_x_ccd2_cam1": 0.015000,
"pix_y_ccd2_cam1": 0.015000,
"ang_ccd2_cam1": 180.000000,
"tilt_x_ccd2_cam1": 0.000000,
"tilt_y_ccd2_cam1": 0.000000,
"x0_ccd3_cam1": -31.652818,
"y0_ccd3_cam1": -31.438350,
"pix_x_ccd3_cam1": 0.015000,
"pix_y_ccd3_cam1": 0.015000,
"ang_ccd3_cam1": -0.024851,
"tilt_x_ccd3_cam1": 0.000000,
"tilt_y_ccd3_cam1": 0.000000,
"x0_ccd4_cam1": 0.833161,
"y0_ccd4_cam1": -31.458180,
"pix_x_ccd4_cam1": 0.015000,
"pix_y_ccd4_cam1": 0.015000,
"ang_ccd4_cam1": 0.001488,
"tilt_x_ccd4_cam1": 0.000000,
"tilt_y_ccd4_cam1": 0.000000,
}
if icam == 1:
parm_dict = {
"ang1_cam1": -0.179412,
"ang2_cam1": -12.017260,
"ang3_cam1": 90.046500,
"fl_cam1": 145.989933,
"opt_coef1_cam1": 1.00000140,
"opt_coef2_cam1": 0.24069345,
"opt_coef3_cam1": 0.15391120,
"opt_coef4_cam1": 4.05433503,
"opt_coef5_cam1": 3.43136895,
"asymang_cam1": 0.00000000,
"asymfac_cam1": 1.00000000,
"x0_ccd1_cam1": 31.653635,
"y0_ccd1_cam1": 31.470291,
"pix_x_ccd1_cam1": 0.015000,
"pix_y_ccd1_cam1": 0.015000,
"ang_ccd1_cam1": 180.010890,
"tilt_x_ccd1_cam1": 0.000000,
"tilt_y_ccd1_cam1": 0.000000,
"x0_ccd2_cam1": -0.827405,
"y0_ccd2_cam1": 31.491388,
"pix_x_ccd2_cam1": 0.015000,
"pix_y_ccd2_cam1": 0.015000,
"ang_ccd2_cam1": 180.000000,
"tilt_x_ccd2_cam1": 0.000000,
"tilt_y_ccd2_cam1": 0.000000,
"x0_ccd3_cam1": -31.543794,
"y0_ccd3_cam1": -31.550699,
"pix_x_ccd3_cam1": 0.015000,
"pix_y_ccd3_cam1": 0.015000,
"ang_ccd3_cam1": -0.006624,
"tilt_x_ccd3_cam1": 0.000000,
"tilt_y_ccd3_cam1": 0.000000,
"x0_ccd4_cam1": 0.922834,
"y0_ccd4_cam1": -31.557268,
"pix_x_ccd4_cam1": 0.015000,
"pix_y_ccd4_cam1": 0.015000,
"ang_ccd4_cam1": -0.015464,
"tilt_x_ccd4_cam1": 0.000000,
"tilt_y_ccd4_cam1": 0.000000,
}
if icam == 2:
parm_dict = {
"ang1_cam1": 0.066596,
"ang2_cam1": 12.007750,
"ang3_cam1": -89.889085,
"fl_cam1": 146.006602,
"opt_coef1_cam1": 1.00000140,
"opt_coef2_cam1": 0.23452229,
"opt_coef3_cam1": 0.33552009,
"opt_coef4_cam1": 1.92009863,
"opt_coef5_cam1": 12.48880182,
"asymang_cam1": 0.00000000,
"asymfac_cam1": 1.00000000,
"x0_ccd1_cam1": 31.615486,
"y0_ccd1_cam1": 31.413644,
"pix_x_ccd1_cam1": 0.015000,
"pix_y_ccd1_cam1": 0.015000,
"ang_ccd1_cam1": 179.993948,
"tilt_x_ccd1_cam1": 0.000000,
"tilt_y_ccd1_cam1": 0.000000,
"x0_ccd2_cam1": -0.832993,
"y0_ccd2_cam1": 31.426621,
"pix_x_ccd2_cam1": 0.015000,
"pix_y_ccd2_cam1": 0.015000,
"ang_ccd2_cam1": 180.000000,
"tilt_x_ccd2_cam1": 0.000000,
"tilt_y_ccd2_cam1": 0.000000,
"x0_ccd3_cam1": -31.548296,
"y0_ccd3_cam1": -31.606976,
"pix_x_ccd3_cam1": 0.015000,
"pix_y_ccd3_cam1": 0.015000,
"ang_ccd3_cam1": 0.000298,
"tilt_x_ccd3_cam1": 0.000000,
"tilt_y_ccd3_cam1": 0.000000,
"x0_ccd4_cam1": 0.896018,
"y0_ccd4_cam1": -31.569542,
"pix_x_ccd4_cam1": 0.015000,
"pix_y_ccd4_cam1": 0.015000,
"ang_ccd4_cam1": -0.006464,
"tilt_x_ccd4_cam1": 0.000000,
"tilt_y_ccd4_cam1": 0.000000,
}
if icam == 3:
parm_dict = {
"ang1_cam1": 0.030756,
"ang2_cam1": 35.978116,
"ang3_cam1": -89.976802,
"fl_cam1": 146.039793,
"opt_coef1_cam1": 1.00000140,
"opt_coef2_cam1": 0.23920416,
"opt_coef3_cam1": 0.13349450,
"opt_coef4_cam1": 4.77768896,
"opt_coef5_cam1": -1.75114744,
"asymang_cam1": 0.00000000,
"asymfac_cam1": 1.00000000,
"x0_ccd1_cam1": 31.575820,
"y0_ccd1_cam1": 31.316510,
"pix_x_ccd1_cam1": 0.015000,
"pix_y_ccd1_cam1": 0.015000,
"ang_ccd1_cam1": 179.968217,
"tilt_x_ccd1_cam1": 0.000000,
"tilt_y_ccd1_cam1": 0.000000,
"x0_ccd2_cam1": -0.890877,
"y0_ccd2_cam1": 31.363511,
"pix_x_ccd2_cam1": 0.015000,
"pix_y_ccd2_cam1": 0.015000,
"ang_ccd2_cam1": 180.000000,
"tilt_x_ccd2_cam1": 0.000000,
"tilt_y_ccd2_cam1": 0.000000,
"x0_ccd3_cam1": -31.630470,
"y0_ccd3_cam1": -31.716942,
"pix_x_ccd3_cam1": 0.015000,
"pix_y_ccd3_cam1": 0.015000,
"ang_ccd3_cam1": -0.024359,
"tilt_x_ccd3_cam1": 0.000000,
"tilt_y_ccd3_cam1": 0.000000,
"x0_ccd4_cam1": 0.824159,
"y0_ccd4_cam1": -31.728751,
"pix_x_ccd4_cam1": 0.015000,
"pix_y_ccd4_cam1": 0.015000,
"ang_ccd4_cam1": -0.024280,
"tilt_x_ccd4_cam1": 0.000000,
"tilt_y_ccd4_cam1": 0.000000,
}
self.parm_dict_list[icam] = parm_dict
def sky_to_sc_mat(self, sc_ra_dec_roll):
"""Calculate the rotation matrix that will convert a vector in ra&dec
into the spacecraft boresite frame
"""
deg2rad = np.pi / 180.0
# Define the 3 euler angles of rotation
xeul = np.zeros((3,), dtype=np.double)
xeul[0] = deg2rad * sc_ra_dec_roll[0]
xeul[1] = np.pi / 2.0 - deg2rad * sc_ra_dec_roll[1]
xeul[2] = deg2rad * sc_ra_dec_roll[2] + np.pi
# Generate the rotation matrix from the 3 euler angles
self.rmat1 = self.eulerm323(xeul)
def sc_to_cam_mat(self, eul):
"""Calculate the rotation matrix that will convert a vector in spacecraft
into the a camera's coords
"""
deg2rad = np.pi / 180.0
# Generate the rotation matrix from the 3 euler angles
xeul = deg2rad * eul
return self.eulerm323(xeul)
def eulerm323(self, eul):
mat1 = self.rotm1(2, eul[0])
mat2 = self.rotm1(1, eul[1])
mata = np.matmul(mat2, mat1)
mat1 = self.rotm1(2, eul[2])
rmat = np.matmul(mat1, mata)
return rmat
def rotm1(self, ax, ang):
mat = np.zeros((3, 3), dtype=np.double)
n1 = ax
n2 = np.mod((n1 + 1), 3)
n3 = np.mod((n2 + 1), 3)
sinang = np.sin(ang)
cosang = np.cos(ang)
mat[n1][n1] = 1.0
mat[n2][n2] = cosang
mat[n3][n3] = cosang
mat[n2][n3] = sinang
mat[n3][n2] = -sinang
return mat
def sphereToCart(self, ras, decs):
"""Convert 3d spherical coordinates to cartesian"""
deg2rad = np.pi / 180.0
rarads = deg2rad * ras
decrads = deg2rad * decs
sinras = np.sin(rarads)
cosras = np.cos(rarads)
sindecs = np.sin(decrads)
cosdecs = np.cos(decrads)
vec0s = cosras * cosdecs
vec1s = sinras * cosdecs
vec2s = sindecs
return vec0s, vec1s, vec2s
def cartToSphere(self, vec):
ra = 0.0
dec = 0.0
norm = np.sqrt(np.sum(vec * vec))
if norm > 0.0:
dec = np.arcsin(vec[2] / norm)
if (not vec[0] == 0.0) or (not vec[1] == 0.0):
ra = np.arctan2(vec[1], vec[0])
ra = np.mod(ra, 2.0 * np.pi)
return ra, dec
def star_in_fov(self, lng, lat):
deg2rad = np.pi / 180.0
inView = False
if lat > 70.0:
vec0, vec1, vec2 = self.sphereToCart(lng, lat)
vec = np.array([vec0, vec1, vec2], dtype=np.double)
norm = np.sqrt(np.sum(vec * vec))
if norm > 0.0:
vec = vec / norm
xlen = np.abs(np.arctan(vec[0] / vec[2]))
ylen = np.abs(np.arctan(vec[1] / vec[2]))
if (xlen <= (12.5 * deg2rad)) and (ylen <= (12.5 * deg2rad)):
inView = True
return inView
def optics_fp(self, icam, lng_deg, lat_deg):
deg2rad = np.pi / 180.0
thetar = np.pi / 2.0 - (lat_deg * deg2rad)
tanth = np.tan(thetar)
cphi = np.cos(deg2rad * lng_deg)
sphi = np.sin(deg2rad * lng_deg)
rfp0 = self.optcon[icam][0] * tanth
noptcon = len(self.optcon[icam])
ii = np.arange(1, noptcon)
rfp = np.sum(self.optcon[icam][1:] * np.power(tanth, 2.0 * (ii - 1)))
xytmp = np.zeros((2,), dtype=np.double)
xytmp[0] = -cphi * rfp0 * rfp
xytmp[1] = -sphi * rfp0 * rfp
return self.make_az_asym(icam, xytmp)
def fp_optics(self, icam, xyfp):
deg2rad = np.pi / 180.0
xy = self.rev_az_asym(icam, xyfp)
rfp_times_rfp0 = np.sqrt(xy[0] * xy[0] + xy[1] * xy[1])
phirad = np.arctan2(-xy[1], -xy[0])
phideg = phirad / deg2rad
if phideg < 0.0:
phideg += 360.0
thetarad = self.tanth_of_r(icam, rfp_times_rfp0)
thetadeg = thetarad / deg2rad
lng_deg = phideg
lat_deg = 90.0 - thetadeg
return lng_deg, lat_deg
def r_of_tanth(self, icam, z):
tanth = np.tan(z)
rfp0 = self.optcon[icam][0] * tanth
noptcon = len(self.optcon[icam])
ii = np.arange(1, noptcon)
rfp = np.sum(self.optcon[icam][1:] * np.power(tanth, 2.0 * (ii - 1)))
return rfp0 * rfp
def tanth_of_r(self, icam, rprp0):
if np.abs(rprp0) > 1.0e-10:
c0 = self.optcon[icam][0]
zi = np.arctan(np.sqrt(rprp0) / c0)
def minFunc(z, icam, rp):
rtmp = self.r_of_tanth(icam, z)
return (rtmp - rprp0) * (rtmp - rprp0)
optResult = opt.minimize(
minFunc,
[zi],
args=(icam, rprp0),
method="Nelder-Mead",
tol=1.0e-10,
options={"maxiter": 500},
)
# print(optResult)
return optResult.x[0]
else:
return 0.0
def make_az_asym(self, icam, xy):
xyp = self.xyrotate(self.asymang[icam], xy)
xypa = np.zeros_like(xyp)
xypa[0] = self.asymfac[icam] * xyp[0]
xypa[1] = xyp[1]
xyout = self.xyrotate(-self.asymang[icam], xypa)
return xyout
def rev_az_asym(self, icam, xyin):
xyp = self.xyrotate(self.asymang[icam], xyin)
xypa = np.zeros_like(xyp)
xypa[0] = xyp[0] / self.asymfac[icam]
xypa[1] = xyp[1]
xyout = self.xyrotate(-self.asymang[icam], xypa)
return xyout
def xyrotate(self, angle_deg, xin):
deg2rad = np.pi / 180.0
ca = np.cos(deg2rad * angle_deg)
sa = np.sin(deg2rad * angle_deg)
xyout = np.zeros_like(xin)
xyout[0] = ca * xin[0] + sa * xin[1]
xyout[1] = -sa * xin[0] + ca * xin[1]
return xyout
def mm_to_pix(self, icam, xy):
"""Convert focal plane to pixel location also need to add in the
auxillary pixels added into FFIs
"""
CCDWD_T = 2048
CCDHT_T = 2058
ROWA = 44
ROWB = 44
COLDK_T = 20
xya = np.copy(xy)
xyb = np.zeros_like(xya)
ccdpx = np.zeros_like(xya)
fitpx = np.zeros_like(xya)
if xya[0] >= 0.0:
if xya[1] >= 0.0:
iccd = 0
xyb[0] = xya[0] - self.ccdxy0[icam][iccd][0]
xyb[1] = xya[1] - self.ccdxy0[icam][iccd][1]
xyccd = self.xyrotate(self.ccdang[icam][iccd], xyb)
ccdpx[0] = (xyccd[0] / self.pixsz[icam][iccd][0]) - 0.5
ccdpx[1] = (xyccd[1] / self.pixsz[icam][iccd][1]) - 0.5
fitpx[0] = (CCDWD_T - ccdpx[0]) + CCDWD_T + 2 * ROWA + ROWB - 1.0
fitpx[1] = (CCDHT_T - ccdpx[1]) + CCDHT_T + 2 * COLDK_T - 1.0
else:
iccd = 3
xyb[0] = xya[0] - self.ccdxy0[icam][iccd][0]
xyb[1] = xya[1] - self.ccdxy0[icam][iccd][1]
xyccd = self.xyrotate(self.ccdang[icam][iccd], xyb)
ccdpx[0] = (xyccd[0] / self.pixsz[icam][iccd][0]) - 0.5
ccdpx[1] = (xyccd[1] / self.pixsz[icam][iccd][1]) - 0.5
fitpx[0] = ccdpx[0] + CCDWD_T + 2 * ROWA + ROWB
fitpx[1] = ccdpx[1]
else:
if xya[1] >= 0.0:
iccd = 1
xyb[0] = xya[0] - self.ccdxy0[icam][iccd][0]
xyb[1] = xya[1] - self.ccdxy0[icam][iccd][1]
xyccd = self.xyrotate(self.ccdang[icam][iccd], xyb)
ccdpx[0] = (xyccd[0] / self.pixsz[icam][iccd][0]) - 0.5
ccdpx[1] = (xyccd[1] / self.pixsz[icam][iccd][1]) - 0.5
fitpx[0] = (CCDWD_T - ccdpx[0]) + ROWA - 1.0
fitpx[1] = (CCDHT_T - ccdpx[1]) + CCDHT_T + 2 * COLDK_T - 1.0
else:
iccd = 2
xyb[0] = xya[0] - self.ccdxy0[icam][iccd][0]
xyb[1] = xya[1] - self.ccdxy0[icam][iccd][1]
xyccd = self.xyrotate(self.ccdang[icam][iccd], xyb)
ccdpx[0] = (xyccd[0] / self.pixsz[icam][iccd][0]) - 0.5
ccdpx[1] = (xyccd[1] / self.pixsz[icam][iccd][1]) - 0.5
fitpx[0] = ccdpx[0] + ROWA
fitpx[1] = ccdpx[1]
return iccd, ccdpx, fitpx
def mm_to_pix_single_ccd(self, icam, xy, iccd):
"""Convert focal plane to pixel location also need to add in the
auxillary pixels added into FFIs
"""
CCDWD_T = 2048
CCDHT_T = 2058
ROWA = 44
ROWB = 44
COLDK_T = 20
xya = np.copy(xy)
xyb = np.zeros_like(xya)
ccdpx = np.zeros_like(xya)
fitpx = np.zeros_like(xya)
if iccd == 0:
xyb[0] = xya[0] - self.ccdxy0[icam][iccd][0]
xyb[1] = xya[1] - self.ccdxy0[icam][iccd][1]
xyccd = self.xyrotate(self.ccdang[icam][iccd], xyb)
ccdpx[0] = (xyccd[0] / self.pixsz[icam][iccd][0]) - 0.5
ccdpx[1] = (xyccd[1] / self.pixsz[icam][iccd][1]) - 0.5
fitpx[0] = (CCDWD_T - ccdpx[0]) + CCDWD_T + 2 * ROWA + ROWB - 1.0
fitpx[1] = (CCDHT_T - ccdpx[1]) + CCDHT_T + 2 * COLDK_T - 1.0
if iccd == 3:
xyb[0] = xya[0] - self.ccdxy0[icam][iccd][0]
xyb[1] = xya[1] - self.ccdxy0[icam][iccd][1]
xyccd = self.xyrotate(self.ccdang[icam][iccd], xyb)
ccdpx[0] = (xyccd[0] / self.pixsz[icam][iccd][0]) - 0.5
ccdpx[1] = (xyccd[1] / self.pixsz[icam][iccd][1]) - 0.5
fitpx[0] = ccdpx[0] + CCDWD_T + 2 * ROWA + ROWB
fitpx[1] = ccdpx[1]
if iccd == 1:
xyb[0] = xya[0] - self.ccdxy0[icam][iccd][0]
xyb[1] = xya[1] - self.ccdxy0[icam][iccd][1]
xyccd = self.xyrotate(self.ccdang[icam][iccd], xyb)
ccdpx[0] = (xyccd[0] / self.pixsz[icam][iccd][0]) - 0.5
ccdpx[1] = (xyccd[1] / self.pixsz[icam][iccd][1]) - 0.5
fitpx[0] = (CCDWD_T - ccdpx[0]) + ROWA - 1.0
fitpx[1] = (CCDHT_T - ccdpx[1]) + CCDHT_T + 2 * COLDK_T - 1.0
if iccd == 2:
xyb[0] = xya[0] - self.ccdxy0[icam][iccd][0]
xyb[1] = xya[1] - self.ccdxy0[icam][iccd][1]
xyccd = self.xyrotate(self.ccdang[icam][iccd], xyb)
ccdpx[0] = (xyccd[0] / self.pixsz[icam][iccd][0]) - 0.5
ccdpx[1] = (xyccd[1] / self.pixsz[icam][iccd][1]) - 0.5
fitpx[0] = ccdpx[0] + ROWA
fitpx[1] = ccdpx[1]
return ccdpx, fitpx
def pix_to_mm_single_ccd(self, icam, ccdpx, iccd):
"""convert pixel to mm focal plane position"""
xyccd = np.zeros_like(ccdpx)
xyccd[0] = (ccdpx[0] + 0.5) * self.pixsz[icam][iccd][0]
xyccd[1] = (ccdpx[1] + 0.5) * self.pixsz[icam][iccd][1]
xyb = self.xyrotate(-self.ccdang[icam][iccd], xyccd)
xya = np.zeros_like(xyb)
xya[0] = xyb[0] + self.ccdxy0[icam][iccd][0]
xya[1] = xyb[1] + self.ccdxy0[icam][iccd][1]
return xya
def radec2pix(self, ras, decs):
"""After the rotation matrices are defined to the actual
ra and dec to pixel coords mapping
"""
nStar = len(ras)
inCamera = np.array([], dtype=int)
ccdNum = np.array([], dtype=int)
fitsxpos = np.array([], dtype=np.double)
fitsypos = np.array([], dtype=np.double)
ccdxpos = np.array([], dtype=np.double)
ccdypos = np.array([], dtype=np.double)
deg2rad = np.pi / 180.0
if self.havePointing == True:
# Convert ra and dec spherical coords to cartesian
vec0s, vec1s, vec2s = self.sphereToCart(ras, decs)
for i in range(nStar):
curVec = np.array([vec0s[i], vec1s[i], vec2s[i]], dtype=np.double)
# Find the new vector in all cameras
for j in range(self.NCAM):
# Do the rotation from ra dec coords to camera coords
camVec = np.matmul(self.rmat4[j], curVec)
# Get the longitude and latitude of camera coords position
lng, lat = self.cartToSphere(camVec)
lng = lng / deg2rad
lat = lat / deg2rad
if self.star_in_fov(lng, lat):
# Get the xy focal plane position in mm
xyfp = self.optics_fp(j, lng, lat)
# Convert mm to pixels
iccd, ccdpx, fitpx = self.mm_to_pix(j, xyfp)
inCamera = np.append(
inCamera, j + 1
) # Als code is base 0 convert to base 1
ccdNum = np.append(ccdNum, iccd + 1) # ""
fitsxpos = np.append(fitsxpos, fitpx[0])
fitsypos = np.append(fitsypos, fitpx[1])
ccdxpos = np.append(ccdxpos, ccdpx[0])
ccdypos = np.append(ccdypos, ccdpx[1])
else:
print("Spacecraft Pointing Not specified!")
return inCamera, ccdNum, fitsxpos, fitsypos, ccdxpos, ccdypos
def radec2pix_nocheck_single(self, ras, decs, cam, iccd):
"""
ra and dec to pixel coords mapping
With no checks and assuming a single target and detector
Supports minimizing for reverse mode
"""
deg2rad = np.pi / 180.0
# Convert ra and dec spherical coords to cartesian
vec0s, vec1s, vec2s = self.sphereToCart(ras, decs)
curVec = np.array([vec0s, vec1s, vec2s], dtype=np.double)
j = cam
# Do the rotation from ra dec coords to camera coords
camVec = np.matmul(self.rmat4[j], curVec)
# Get the longitude and latitude of camera coords position
lng, lat = self.cartToSphere(camVec)
lng = lng / deg2rad
lat = lat / deg2rad
# Get the xy focal plane position in mm
xyfp = self.optics_fp(j, lng, lat)
# Convert mm to pixels
ccdpx, fitpx = self.mm_to_pix_single_ccd(j, xyfp, iccd)
ccdNum = iccd + 1
fitsxpos = fitpx[0]
fitsypos = fitpx[1]
ccdxpos = ccdpx[0]
ccdypos = ccdpx[1]
return ccdNum, fitsxpos, fitsypos, ccdxpos, ccdypos, lat
def pix2radec_nocheck_single(self, cam, iccd, ccdpx):
"""
Reverse the transform going from pixel coords to Ra & Dec
"""
deg2rad = np.pi / 180.0
# Convert pixels to mm
xyfp = self.pix_to_mm_single_ccd(cam, ccdpx, iccd)
lng_deg, lat_deg = self.fp_optics(cam, xyfp)
vcam0, vcam1, vcam2 = self.sphereToCart(lng_deg, lat_deg)
vcam = np.array([vcam0, vcam1, vcam2], dtype=np.double)
curVec = np.matmul(np.transpose(self.rmat4[cam]), vcam)
ra, dec = self.cartToSphere(curVec)
return ra / deg2rad, dec / deg2rad
class TESS_Spacecraft_Pointing_Data:
# Hard coded spacecraft pointings by Sector
# When adding sectors the arg2 needs to end +1 from sector
# due to the np.arange function ending at arg2-1
sectors = np.arange(1, max_sector + 1, dtype=int)
# Arrays are broken up into the following sectors:
# Line 1: Sectors 1-5 Start Year 1
# Line 2: Secotrs 6-9
# Line 3: Sectors 10-13 End Year 1
# Line 4: Sectors 14-17 Start Year 2
# Line 5: Sectors 18-22
# Line 6: Sectors 23-26 End Year 2
# Line 7: Sectors 27-30 Star Year 3
# Line 8: Sectors 31-34
# Line 9: Sectors 35-38
# Line 10: Sectors 39 End Year 3
# Line 11: Sectors 40-43 Star Year 4
# Line 12: S 44-47
# Line 13: S 48-51
# Line 14: S 52-55 END YEAR 4
# Line 15: S 56-59 START YEAR 5
# Line 16: S 60-63
# Line 17: S 64-67
# Line 18: S 68-69 END Year 5
# Line 19: S 70-72 START Year 6
# Line 20: S 73-76
# Line 21: S 77-80
# Line 22: S 81-83 END Year 6
# Line 23: S 84-87 START Year 7
# Line 24: S 88-90
# Line 25: S 91-93
# Line 26: S 94-96 END Year 7
# Line 27: S 97-99 START Year 8
# Line 28: S 100-102
# Line 29: S 103-105
# Line 30: S 106-107 END Year 8
### NOTE IF you add Sectors be sure to update the allowed range
### for sectors in argparse arguments!!!
ras = np.array(
[
352.6844,
16.5571,
36.3138,
55.0070,
73.5382,
92.0096,
110.2559,
128.1156,
145.9071,
165.0475,
189.1247,
229.5885,
298.6671,
276.7169,
280.3985,
282.4427,
351.2381,
16.1103,
60.2026,
129.3867,
171.7951,
197.1008,
217.2879,
261.4516,
265.6098,
270.1381,
326.8525,
357.2944,
18.9190,
38.3564,
57.6357,
77.1891,
96.5996,
115.2951,
133.2035,
150.9497,
170.2540,
195.7176,
242.1981,
273.0766,
277.6209,
13.0140,
49.5260,
89.6066,
130.2960,
157.6997,
143.3807,
179.4254,
202.6424,
221.8575,
239.4257,
266.3618,
270.8126,
290.1210,
307.8655,
324.2778,
344.2275,
9.3118,
52.9755,
125.6742,
118.0446,
135.2412,
153.0613,
173.2653,
201.6239,
259.1702,
326.7691,
359.2829,
20.0449,
24.0414,
77.3449,
133.7631,
80.6709,
261.2194,
254.9290,
253.5335,
255.8590,
260.4232,
266.0595,
275.6978,
292.5709,
309.1299,
325.8933,
343.8717,
5.8776,
42.6741,
97.9629,
116.8315,
135.8409,
155.0705,
227.0409,
311.9633,
260.5165,
323.4700,
355.7011,
16.7112,
41.6611,
80.1703,
157.2084,
180.5408,
209.0351,
249.2003,
297.5655,
335.4735,
2.5814,
25.0163,
45.8675,
],
dtype=float,
)
decs = np.array(
[
-64.8531,
-54.0160,
-44.2590,
-36.6420,
-31.9349,
-30.5839,
-32.6344,
-37.7370,
-45.3044,
-54.8165,
-65.5369,
-75.1256,
-76.3281,
62.4756,
64.0671,
66.1422,
57.8456,
67.9575,
76.2343,
75.2520,
65.1924,
53.7434,
43.8074,
63.1181,
61.9383,
61.5637,
-72.4265,
-63.0056,
-52.8296,
-43.3178,
-35.7835,
-31.3957,