-
Notifications
You must be signed in to change notification settings - Fork 45
/
exotic.py
2202 lines (1842 loc) · 100 KB
/
exotic.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
# ########################################################################### #
# Copyright (c) 2019-2020, California Institute of Technology.
# All rights reserved. Based on Government Sponsored Research under
# contracts NNN12AA01C, NAS7-1407 and/or NAS7-03001.
#
# Redistribution and use in source and binary forms, with or without
# modification, are permitted provided that the following conditions
# are met:
# 1. Redistributions of source code must retain the above copyright
# notice, this list of conditions and the following disclaimer.
# 2. Redistributions in binary form must reproduce the above copyright
# notice, this list of conditions and the following disclaimer in
# the documentation and/or other materials provided with the
# distribution.
# 3. Neither the name of the California Institute of
# Technology (Caltech), its operating division the Jet Propulsion
# Laboratory (JPL), the National Aeronautics and Space
# Administration (NASA), nor the names of its contributors may be
# used to endorse or promote products derived from this software
# without specific prior written permission.
#
# THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
# "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
# LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
# A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE CALIFORNIA
# INSTITUTE OF TECHNOLOGY BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
# SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
# TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
# PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
# LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
# NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
# SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
#
# ########################################################################### #
# EXOplanet Transit Interpretation Code (EXOTIC)
# # NOTE: See companion file version.py for version info.
# ########################################################################### #
# -- IMPORTS START ------------------------------------------------------------
# ########## IMPORTS -- PRELOAD ANIMATION START ##########
try: # animation
from animate import *
except ImportError:
from .animate import *
if __name__ == "__main__":
print("Importing Python Packages - please wait.")
animate_toggle(True)
# ########## IMPORTS -- PRELOAD ANIMATION END ##########
# preload to limit import warnings
import warnings
from astropy.utils.exceptions import AstropyDeprecationWarning
warnings.simplefilter('ignore', category=AstropyDeprecationWarning)
# standard imports
import argparse
# Image alignment import
import astroalign as aa
aa.PIXEL_TOL = 1
# aa.NUM_NEAREST_NEIGHBORS=10
# astropy imports
import astropy.units as u
from astropy.coordinates import SkyCoord, EarthLocation, AltAz
from astropy.io import fits
import astropy.time
from astropy.visualization import astropy_mpl_style
from astropy.wcs import WCS, FITSFixedWarning
from astroquery.simbad import Simbad
from astroquery.gaia import Gaia
# UTC to BJD converter import
from barycorrpy import utc_tdb
# julian conversion imports
import dateutil.parser as dup
from pathlib import Path
import pyvo as vo
import logging
from logging.handlers import TimedRotatingFileHandler
from matplotlib.animation import FuncAnimation
# Pyplot imports
import matplotlib.pyplot as plt
# from numba import njit
import numpy as np
# photometry
from photutils import CircularAperture
# scipy imports
from scipy.optimize import least_squares
from scipy.stats import mode
from scipy.signal import savgol_filter
from scipy.ndimage import binary_erosion
from skimage.util import view_as_windows
from skimage.transform import SimilarityTransform
# error handling for scraper
from tenacity import retry, stop_after_delay
# ########## EXOTIC imports ##########
try: # light curve numerics
from .api.elca import lc_fitter, binner, transit, get_phase
except ImportError: # package import
from api.elca import lc_fitter, binner, transit, get_phase
try: # output files
from inputs import Inputs
except ImportError: # package import
from .inputs import Inputs
try: # plate solution
from .api.plate_solution import PlateSolution
except ImportError: # package import
from api.plate_solution import PlateSolution
try: # nonlinear limb darkening numerics
from .api.gael_ld import createldgrid
except ImportError: # package import
from api.gael_ld import createldgrid
try: # filters
from .api.filters import fwhm
except ImportError: # package import
from api.filters import fwhm
try: # nea
from .api.nea import NASAExoplanetArchive
except ImportError: # package import
from api.nea import NASAExoplanetArchive
try: # output files
from output_files import OutputFiles
except ImportError: # package import
from .output_files import OutputFiles
try: # plots
from plots import plot_fov, plot_centroids, plot_obs_stats, plot_final_lightcurve, plot_flux
except ImportError: # package import
from .plots import plot_fov, plot_centroids, plot_obs_stats, plot_final_lightcurve, plot_flux
try: # tools
from utils import round_to_2, user_input
except ImportError: # package import
from .utils import round_to_2, user_input
try: # simple version
from .version import __version__
except ImportError: # package import
from version import __version__
animate_toggle() # CLOSE PRELOAD ANIMATION
# -- IMPORTS END --------------------------------------------------------------
# SETTINGS
plt.style.use(astropy_mpl_style)
# To increase memory allocation for EXOTIC; allows for more fits files
# resource.setrlimit(resource.RLIMIT_STACK, (resource.RLIM_INFINITY, resource.RLIM_INFINITY))
# ################### END PROPERTIES/SETTINGS ############################### #
# logging -- https://docs.python.org/3/library/logging.html
log = logging.getLogger(__name__)
def log_info(string, warn=False, error=False):
if error:
print(f"\033[91m {string}\033[00m")
elif warn:
print(f"\033[33m {string}\033[00m")
else:
print(string)
log.debug(string)
return True
def sigma_clip(ogdata, sigma=3, dt=21, po=2):
nanmask = np.isnan(ogdata)
if po < dt <= len(ogdata):
mdata = savgol_filter(ogdata[~nanmask], window_length=dt, polyorder=po)
# mdata = median_filter(ogdata[~nanmask], dt)
res = ogdata[~nanmask] - mdata
std = np.nanmedian([np.nanstd(np.random.choice(res, 25)) for i in range(100)])
# std = np.nanstd(res) # biased from large outliers
sigmask = np.abs(res) > sigma * std
nanmask[~nanmask] = sigmask
return nanmask
# Get Julian time, don't need to divide by 2 since assume mid-EXPOSURE
# Find separate funciton in code that does julian conversion to BJD_TDB
# Method that gets and returns the julian time of the observation
def getJulianTime(header):
exptime_offset = 0
exp = header.get('EXPTIME') # checking for variation in .fits header format
if exp:
exp = exp
else:
exp = header.get('EXPOSURE')
# Grab the BJD first
if 'BJD_TDB' in header:
julianTime = float(header['BJD_TDB'])
# If the time is from the beginning of the observation, then need to calculate mid-exposure time
if "start" in header.comments['BJD_TDB']:
exptime_offset = exp / 2. / 60. / 60. / 24. # assume exptime is in seconds for now
elif "BJD_TBD" in header: # looks like TheSky misspells BJD_TDB as BJD_TBD -> hardcoding this common misspelling
julianTime = float(header['BJD_TBD'])
# If the time is from the beginning of the observation, then need to calculate mid-exposure time
if "start" in header.comments['BJD_TBD']:
exptime_offset = exp / 2. / 60. / 60. / 24. # assume exptime is in seconds for now
elif 'BJD' in header:
julianTime = float(header['BJD'])
# If the time is from the beginning of the observation, then need to calculate mid-exposure time
if "start" in header.comments['BJD']:
exptime_offset = exp / 2. / 60. / 60. / 24. # assume exptime is in seconds for now
# then the DATE-OBS
elif "UT-OBS" in header:
gDateTime = header['UT-OBS'] # gets the gregorian date and time from the fits file header
dt = dup.parse(gDateTime)
atime = astropy.time.Time(dt)
julianTime = atime.jd
# If the time is from the beginning of the observation, then need to calculate mid-exposure time
if "start" in header.comments['UT-OBS']:
exptime_offset = exp / 2. / 60. / 60. / 24. # assume exptime is in seconds for now
# Then Julian Date
elif 'JULIAN' in header:
julianTime = float(header['JULIAN'])
# If the time is from the beginning of the observation, then need to calculate mid-exposure time
if "start" in header.comments['JULIAN']:
exptime_offset = exp / 2. / 60. / 60. / 24. # assume exptime is in seconds for now
# Then MJD-OBS last, as in the MicroObservatory headers, it has less precision
elif ("MJD-OBS" in header) and ("epoch" not in header.comments['MJD-OBS']):
julianTime = float(header["MJD-OBS"]) + 2400000.5
# If the time is from the beginning of the observation, then need to calculate mid-exposure time
if "start" in header.comments['MJD-OBS']:
exptime_offset = exp / 2. / 60. / 60. / 24. # assume exptime is in seconds for now
else:
if 'T' in header['DATE-OBS']:
gDateTime = header['DATE-OBS'] # gets the gregorian date and time from the fits file header
else:
gDateTime = f"{header['DATE-OBS']}T{header['TIME-OBS']}"
dt = dup.parse(gDateTime)
atime = astropy.time.Time(dt)
julianTime = atime.jd
# If the time is from the beginning of the observation, then need to calculate mid-exposure time
if "start" in header.comments['DATE-OBS']:
exptime_offset = exp / 2. / 60. / 60. / 24. # assume exptime is in seconds for now
# If the mid-exposure time is given in the fits header, then no offset is needed to calculate the mid-exposure time
return julianTime + exptime_offset
# Method that gets and returns the airmass from the fits file (Really the Altitude)
def getAirMass(image_header, ra, dec, lati, longit, elevation):
# Grab airmass from image header; if not listed, calculate it from TELALT; if that isn't listed, then calculate it the hard way
if 'AIRMASS' in image_header:
am = float(image_header['AIRMASS'])
elif 'TELALT' in image_header:
alt = float(image_header['TELALT']) # gets the airmass from the fits file header in (sec(z)) (Secant of the zenith angle)
cosam = np.cos((np.pi / 180) * (90.0 - alt))
am = 1 / cosam
else:
# pointing = SkyCoord(str(astropy.coordinates.Angle(raStr+" hours").deg)+" "+str(astropy.coordinates.Angle(decStr+" degrees").deg ), unit=(u.deg, u.deg), frame='icrs')
pointing = SkyCoord(str(ra) + " " + str(dec), unit=(u.deg, u.deg), frame='icrs')
location = EarthLocation.from_geodetic(lat=lati * u.deg, lon=longit * u.deg, height=elevation)
atime = astropy.time.Time(getJulianTime(image_header), format='jd', scale='utc', location=location)
pointingAltAz = pointing.transform_to(AltAz(obstime=atime, location=location))
am = float(pointingAltAz.secz)
return am
# Convert time units to BJD_TDB if pre-reduced file not in proper units
def timeConvert(timeList, timeFormat, pDict, info_dict):
# Perform appropriate conversion for each time format if needed
if timeFormat == 'JD_UTC':
convertedTimes = utc_tdb.JDUTC_to_BJDTDB(timeList, ra=pDict['ra'], dec=pDict['dec'],
lat=info_dict['lat'], longi=info_dict['long'], alt=info_dict['elev'])
elif timeFormat == 'MJD_UTC':
convertedTimes = utc_tdb.JDUTC_to_BJDTDB(timeList + 2400000.5, ra=pDict['ra'], dec=pDict['dec'],
lat=info_dict['lat'], longi=info_dict['long'], alt=info_dict['elev'])
timeList = convertedTimes[0]
return timeList
# Convert magnitude units to flux if pre-reduced file not in flux already
def fluxConvert(fluxList, errorList, fluxFormat):
# If units already in flux, do nothing, perform appropriate conversions to flux otherwise
if fluxFormat == 'magnitude':
convertedPositiveErrors = 10. ** ((-1. * (fluxList + errorList)) / 2.5)
convertedNegativeErrors = 10. ** ((-1. * (fluxList - errorList)) / 2.5)
fluxList = 10. ** ((-1. * fluxList) / 2.5)
elif fluxFormat == 'millimagnitude':
convertedPositiveErrors = 10. ** ((-1. * ((fluxList + errorList) / 1000.)) / 2.5)
convertedNegativeErrors = 10. ** ((-1. * ((fluxList - errorList) / 1000.)) / 2.5)
fluxList = 10. ** (-1. * (fluxList / 1000.) / 2.5)
# Use distance from mean of upper/lower error bounds to calculate new sigma values
positiveErrorDistance = abs(convertedPositiveErrors - fluxList)
negativeErrorDistance = abs(convertedNegativeErrors - fluxList)
meanErrorList = (positiveErrorDistance * negativeErrorDistance) ** 0.5
return fluxList, meanErrorList
# Check for difference between NEA and initialization file
def check_parameters(init_parameters, parameters):
different = False
uncert = 1 / 36
for key, value in parameters.items():
if key in ['ra', 'dec'] and init_parameters[key]:
if not parameters[key] - uncert <= init_parameters[key] <= parameters[key] + uncert:
different = True
break
continue
if value != init_parameters[key]:
different = True
break
if different:
log_info("\nDifference(s) found between initialization file parameters and "
"those scraped by EXOTIC from the NASA Exoplanet Archive."
"\nWould you like:"
"\n (1) EXOTIC to adopt of all of your defined parameters or"
"\n (2) to review the ones scraped from the Archive that differ?")
opt = user_input("Enter 1 or 2: ", type_=int, values=[1, 2])
if opt == 2:
return True
else:
return False
# --------PLANETARY PARAMETERS UI------------------------------------------
# Get the user's confirmation of values that will later be used in lightcurve fit
def get_planetary_parameters(candplanetbool, userpdict, pdict=None):
log_info("*******************************************")
log_info("Planetary Parameters for Lightcurve Fitting")
# The order of planet_params list must match the pDict that is declared when scraping the NASA Exoplanet Archive
planet_params = ["Target Star RA in the form: HH:MM:SS (ignore the decimal values)",
"Target Star DEC in form: <sign>DD:MM:SS (ignore the decimal values and don't forget the '+' or '-' sign!)",
"Planet's Name",
"Host Star's Name",
"Orbital Period (days)",
"Orbital Period Uncertainty (days) \n(Keep in mind that 1.2e-34 is the same as 1.2 x 10^-34)",
"Published Mid-Transit Time (BJD_UTC)",
"Mid-Transit Time Uncertainty (BJD-UTC)",
"Ratio of Planet to Stellar Radius (Rp/Rs)",
"Ratio of Planet to Stellar Radius (Rp/Rs) Uncertainty",
"Ratio of Distance to Stellar Radius (a/Rs)",
"Ratio of Distance to Stellar Radius (a/Rs) Uncertainty",
"Orbital Inclination (deg)",
"Orbital Inclination (deg) Uncertainty",
"Orbital Eccentricity (0 if null)",
"Star Effective Temperature (K)",
"Star Effective Temperature Positive Uncertainty (K)",
"Star Effective Temperature Negative Uncertainty (K)",
"Star Metallicity ([FE/H])",
"Star Metallicity Positive Uncertainty ([FE/H])",
"Star Metallicity Negative Uncertainty ([FE/H])",
"Star Surface Gravity (log(g))",
"Star Surface Gravity Positive Uncertainty (log(g))",
"Star Surface Gravity Negative Uncertainty (log(g))"]
# Conversion between hours to degrees if user entered ra and dec
if userpdict['ra'] is None:
userpdict['ra'] = user_input(f"\nEnter the {planet_params[0]}: ", type_=str)
if userpdict['dec'] is None:
userpdict['dec'] = user_input(f"\nEnter the {planet_params[1]}: ", type_=str)
if type(userpdict['ra']) and type(userpdict['dec']) is str:
userpdict['ra'], userpdict['dec'] = radec_hours_to_degree(userpdict['ra'], userpdict['dec'])
radeclist = ['ra', 'dec']
if not candplanetbool:
for idx, item in enumerate(radeclist):
uncert = 20 / 3600
if pdict[item] - uncert <= userpdict[item] <= pdict[item] + uncert:
continue
else:
log_info(f"\n\nWarning: {pdict['pName']} initialization file's {planet_params[idx]} does not match "
"the value scraped by EXOTIC from the NASA Exoplanet Archive.\n", warn=True)
log_info(f"\tNASA Exoplanet Archive value (degrees): {pdict[item]}", warn=True)
log_info(f"\tInitialization file value (degrees): {userpdict[item]}", warn=True)
log_info("\nWould you like to:"
"\n (1) use NASA Exoplanet Archive value, "
"\n (2) use initialization file value, or "
"\n (3) enter in a new value.", warn=True)
option = user_input("Which option do you choose? (1/2/3): ", type_=int, values=[1, 2, 3])
if option == 1:
userpdict[item] = pdict[item]
elif option == 2:
continue
else:
userpdict['ra'] = user_input(f"Enter the {planet_params[0]}: ", type_=str)
userpdict['dec'] = user_input(f"Enter the {planet_params[1]}: ", type_=str)
break
if type(userpdict['ra']) and type(userpdict['dec']) is str:
userpdict['ra'], userpdict['dec'] = radec_hours_to_degree(userpdict['ra'], userpdict['dec'])
# Exoplanet confirmed in NASA Exoplanet Archive
if not candplanetbool:
log_info(f"*** Here are the values scraped from the NASA Exoplanet Archive for {pdict['pName']} that were not "
"set (or set to null) in your initialization file. ***")
for i, key in enumerate(userpdict):
if key in ('ra', 'dec'):
continue
if key in ('pName', 'sName'):
userpdict[key] = pdict[key]
# Initialization planetary parameters match NEA
if pdict[key] == userpdict[key]:
continue
# Initialization planetary parameters don't match NASA Exoplanet Archive
if userpdict[key] is not None:
log_info(f"\n\nWarning: {pdict['pName']} initialization file's {planet_params[i]} does not match "
"the value scraped by EXOTIC from the NASA Exoplanet Archive.\n", warn=True)
log_info(f"\tNASA Exoplanet Archive value: {pdict[key]}", warn=True)
log_info(f"\tInitialization file value: {userpdict[key]}", warn=True)
log_info("\nWould you like to: "
"\n (1) use NASA Exoplanet Archive value, "
"\n (2) use initialization file value, or "
"\n (3) enter in a new value.", warn=True)
option = user_input("Which option do you choose? (1/2/3): ", type_=int, values=[1, 2, 3])
if option == 1:
userpdict[key] = pdict[key]
elif option == 2:
continue
else:
userpdict[key] = user_input(f"Enter the {planet_params[i]}: ", type_=type(userpdict[key]))
# Did not use initialization file or null
else:
log_info(f"\n {pdict['pName']} {planet_params[i]}: {pdict[key]}")
agreement = user_input("Do you agree? (y/n): ", type_=str, values=['y', 'n'])
if agreement == 'y':
userpdict[key] = pdict[key]
else:
userpdict[key] = user_input(f"Enter the {planet_params[i]}: ", type_=type(pdict[key]))
# Exoplanet not confirmed in NASA Exoplanet Archive
else:
for i, key in enumerate(userpdict):
if key in ('ra', 'dec'):
continue
# Used initialization file and is not empty
if userpdict[key] is not None:
agreement = user_input(f"{planet_params[i]}: {userpdict[key]} \nDo you agree? (y/n): ",
type_=str, values=['y', 'n'])
if agreement == 'y':
continue
else:
userpdict[key] = user_input(f"Enter the {planet_params[i]}: ", type_=type(userpdict[key]))
# Did not use initialization file
else:
if key in ('pName', 'sName'):
userpdict[key] = user_input(f"\nEnter the {planet_params[i]}: ", type_=str)
else:
userpdict[key] = user_input(f"Enter the {planet_params[i]}: ", type_=float)
return userpdict
# Conversion of Right Ascension and Declination: hours -> degrees
def radec_hours_to_degree(ra, dec):
while True:
try:
ra = ra.replace(' ', '').replace(':', ' ')
dec = dec.replace(' ', '').replace(':', ' ')
c = SkyCoord(ra + ' ' + dec, unit=(u.hourangle, u.deg))
return c.ra.degree, c.dec.degree
except ValueError:
log_info("Error: The format entered for Right Ascension and/or Declination is not correct, "
"please try again.", error=True)
ra = input("Input the Right Ascension of target (HH:MM:SS): ")
dec = input("Input the Declination of target (<sign>DD:MM:SS): ")
class LimbDarkening:
def __init__(self, teff=None, teffpos=None, teffneg=None, met=None, metpos=None, metneg=None,
logg=None, loggpos=None, loggneg=None, wl_min=None, wl_max=None, filter_type=None):
self.priors = {'T*': teff, 'T*_uperr': teffpos, 'T*_lowerr': teffneg,
'FEH*': met, 'FEH*_uperr': metpos, 'FEH*_lowerr': metneg,
'LOGG*': logg, 'LOGG*_uperr': loggpos, 'LOGG*_lowerr': loggneg}
self.filter_type = filter_type
self.wl_min = wl_min
self.wl_max = wl_max
self.fwhm = fwhm
self.ld0 = self.ld1 = self.ld2 = self.ld3 = None
def nonlinear_ld(self):
if self.filter_type and not (self.wl_min or self.wl_max):
self._standard()
elif self.wl_min or self.wl_max:
self._custom()
else:
opt = user_input("\nWould you like EXOTIC to calculate your limb darkening parameters "
"with uncertainties? (y/n):", type_=str, values=['y', 'n'])
if opt == 'y':
opt = user_input("Please enter 1 to use a standard filter or 2 for a customized filter:",
type_=int, values=[1, 2])
if opt == 1:
self._standard()
elif opt == 2:
self._custom()
else:
self._user_entered()
return self.ld0, self.ld1, self.ld2, self.ld3, self.filter_type, self.wl_min * 1000, self.wl_max * 1000
def _standard_list(self):
log_info("\n\n***************************")
log_info("Limb Darkening Coefficients")
log_info("***************************")
log_info("\nThe standard bands that are available for limb darkening parameters (https://www.aavso.org/filters)"
"\nas well as filters for MObs and LCO (0.4m telescope) datasets:\n")
for key, value in self.fwhm.items():
log_info(f"\t{key[1]}: {key[0]} - ({value[0]:.2f}-{value[1]:.2f}) nm")
def _standard(self):
while True:
try:
if not self.filter_type:
self._standard_list()
self.filter_type = user_input("\nPlease enter in the filter type (EX: Johnson V, V, STB, RJ):",
type_=str)
for key, value in self.fwhm.items():
if self.filter_type in (key[0], key[1]) and self.filter_type != 'N/A':
self.filter_type = (key[0], key[1])
break
else:
raise KeyError
break
except KeyError:
log_info("\nError: The entered filter is not in the provided list of standard filters.", error=True)
self.filter_type = None
self.wl_min = self.fwhm[self.filter_type][0]
self.wl_max = self.fwhm[self.filter_type][1]
self.filter_type = self.filter_type[1]
self._calculate_ld()
def _custom(self):
self.filter_type = 'N/A'
if not self.wl_min:
self.wl_min = user_input("FWHM Minimum wavelength (nm):", type_=float)
if not self.wl_max:
self.wl_max = user_input("FWHM Maximum wavelength (nm):", type_=float)
self._calculate_ld()
def _user_entered(self):
self.filter_type = user_input("\nEnter in your filter name:", type_=str)
ld_0 = user_input("\nEnter in your first nonlinear term:", type_=float)
ld0_unc = user_input("Enter in your first nonlinear term uncertainty:", type_=float)
ld_1 = user_input("\nEnter in your second nonlinear term:", type_=float)
ld1_unc = user_input("Enter in your second nonlinear term uncertainty:", type_=float)
ld_2 = user_input("\nEnter in your third nonlinear term:", type_=float)
ld2_unc = user_input("Enter in your third nonlinear term uncertainty:", type_=float)
ld_3 = user_input("\nEnter in your fourth nonlinear term:", type_=float)
ld3_unc = user_input("Enter in your fourth nonlinear term uncertainty:", type_=float)
self.ld0, self.ld1, self.ld2, self.ld3 = (ld_0, ld0_unc), (ld_1, ld1_unc), (ld_2, ld2_unc), (ld_3, ld3_unc)
log.debug(f"Filter name: {self.filter_type}")
log.debug(f"User-defined nonlinear limb-darkening coefficients: {ld_0}+/-{ld0_unc}, {ld_1}+/-{ld1_unc}, "
f"{ld_2}+/-{ld2_unc}, {ld_3}+/-{ld3_unc}")
def _calculate_ld(self):
self.wl_min = self.wl_min / 1000
self.wl_max = self.wl_max / 1000
ld_params = createldgrid(np.array([self.wl_min]), np.array([self.wl_max]), self.priors)
self.ld0 = ld_params['LD'][0][0], ld_params['ERR'][0][0]
self.ld1 = ld_params['LD'][1][0], ld_params['ERR'][1][0]
self.ld2 = ld_params['LD'][2][0], ld_params['ERR'][2][0]
self.ld3 = ld_params['LD'][3][0], ld_params['ERR'][3][0]
log.debug("EXOTIC-calculated nonlinear limb-darkening coefficients: ")
log.debug(f"{ld_params['LD'][0][0]} +/- + {ld_params['ERR'][0][0]}")
log.debug(f"{ld_params['LD'][1][0]} +/- + {ld_params['ERR'][1][0]}")
log.debug(f"{ld_params['LD'][2][0]} +/- + {ld_params['ERR'][2][0]}")
log.debug(f"{ld_params['LD'][3][0]} +/- + {ld_params['ERR'][3][0]}")
def corruption_check(files):
valid_files = []
for file in files:
try:
hdul = fits.open(name=file, memmap=False, cache=False, lazy_load_hdus=False, ignore_missing_end=True)
valid_files.append(file)
except OSError as e:
log_info(f"Warning: corrupted file found and removed from reduction\n\t-File: {file}\n\t-Reason: {e}", warn=True)
finally:
if getattr(hdul, "close", None) and callable(hdul.close):
hdul.close()
del hdul
return valid_files
def check_wcs(fits_file, save_directory, plate_opt, rt=False):
wcs_file = None
if plate_opt == 'y' and not rt:
wcs_file = get_wcs(fits_file, save_directory)
if not wcs_file:
if search_wcs(fits_file).is_celestial:
log_info("Your FITS files have WCS (World Coordinate System) information in their headers. "
"EXOTIC will proceed to use these. "
"NOTE: If you do not trust your WCS coordinates, "
"please restart EXOTIC after enabling plate solutions via astrometry.net.")
wcs_file = fits_file
return wcs_file
def search_wcs(file):
with warnings.catch_warnings():
warnings.simplefilter('ignore', category=FITSFixedWarning)
header = fits.getheader(filename=file)
return WCS(header)
def get_wcs(file, directory=""):
log_info("\nGetting the plate solution for your imaging file to translate pixel coordinates on the sky. "
"\nPlease wait....")
animate_toggle(True)
wcs_obj = PlateSolution(file=file, directory=directory)
wcs_file = wcs_obj.plate_solution()
animate_toggle()
return wcs_file
# Getting the right ascension and declination for every pixel in imaging file if there is a plate solution
def get_radec(header):
wcs_header = WCS(header)
xaxis = np.arange(header['NAXIS1'])
yaxis = np.arange(header['NAXIS2'])
x, y = np.meshgrid(xaxis, yaxis)
return wcs_header.all_pix2world(x, y, 1)
def deg_to_pix(exp_ra, exp_dec, ra_list, dec_list):
dist = (ra_list - exp_ra) ** 2 + (dec_list - exp_dec) ** 2
return np.unravel_index(dist.argmin(), dist.shape)
# Check the ra and dec against the plate solution to see if the user entered in the correct values
def check_targetpixelwcs(pixx, pixy, expra, expdec, ralist, declist):
while True:
try:
uncert = 1 / 36
# Margins are within 100 arcseconds
if not (expra - uncert <= ralist[int(pixy)][int(pixx)] <= expra + uncert):
log_info("\nWarning: The X Pixel Coordinate entered does not match the target's Right Ascension.", warn=True)
raise ValueError
if not (expdec - uncert <= declist[int(pixy)][int(pixx)] <= expdec + uncert):
log_info("\nWarning: The Y Pixel Coordinate entered does not match the target's Declination.", warn=True)
raise ValueError
return pixx, pixy
except ValueError:
log_info(f"Your input pixel coordinates: [{pixx}, {pixy}]")
exotic_pixy, exotic_pixx = deg_to_pix(expra, expdec, ralist, declist)
log_info(f"EXOTIC's calculated pixel coordinates: [{exotic_pixx}, {exotic_pixy}]")
opt = user_input("Would you like to re-enter the pixel coordinates? (y/n): ", type_=str, values=['y', 'n'])
# User wants to change their coordinates
if opt == 'y':
searchopt = user_input(f"Here are the suggested pixel coordinates:"
f" X Pixel: {exotic_pixx}"
f" Y Pixel: {exotic_pixy}"
"\nWould you like to use these? (y/n): ",
type_=str, values=['y', 'n'])
# Use the coordinates found by code
if searchopt == 'y':
return exotic_pixx, exotic_pixy
# User enters their own coordinates to be re-checked
else:
pixx = user_input("Please re-enter the target star's X Pixel Coordinate: ", type_=int)
pixy = user_input("Please re-enter the target star's Y Pixel Coordinate: ", type_=int)
else:
# User does not want to change coordinates even though they don't match the expected ra and dec
return pixx, pixy
# Checks if comparison star is variable via querying SIMBAD
def variableStarCheck(ra, dec):
# Convert comparison star coordinates from pixel to WCS
sample = SkyCoord(ra * u.deg, dec * u.deg, frame='fk5')
radius = u.Quantity(20.0, u.arcsec)
# Query GAIA first to check for variability using the phot_variable_flag trait
gaia_result = gaia_query(sample, radius)
if not gaia_result:
log_info("Warning: Your comparison star cannot be resolved in the Gaia star database; "
"EXOTIC cannot check if it is variable or not. "
"\nEXOTIC will still include this star in the reduction. "
"\nPlease proceed with caution as we cannot check for stellar variability.\n", warn=True)
else:
# Individually go through the phot_variable_flag indicator for each star to see if variable or not
variableFlagList = gaia_result.columns["phot_variable_flag"]
constantCounter = 0
for currFlag in variableFlagList:
if currFlag == "VARIABLE":
return True
elif currFlag == "NOT_AVAILABLE":
continue
elif currFlag == "CONSTANT":
constantCounter += 1
if constantCounter == len(variableFlagList):
return False
# Query SIMBAD and search identifier result table to determine if comparison star is variable in any form
# This is a secondary check if GAIA query returns inconclusive results
star_name = simbad_query(sample)
if not star_name:
log_info("Warning: Your comparison star cannot be resolved in the SIMBAD star database; "
"EXOTIC cannot check if it is variable or not. "
"\nEXOTIC will still include this star in the reduction. "
"\nPlease proceed with caution as we cannot check for stellar variability.\n", warn=True)
return False
else:
identifiers = Simbad.query_objectids(star_name)
for currName in identifiers:
if "V*" in currName[0]:
return True
return False
@retry(stop=stop_after_delay(30))
def gaia_query(sample, radius):
try:
gaia_query = Gaia.cone_search(sample, radius)
return gaia_query.get_results()
except Exception:
return False
@retry(stop=stop_after_delay(30))
def simbad_query(sample):
try:
simbad_result = Simbad.query_region(sample, radius=20 * u.arcsec)
return simbad_result['MAIN_ID'][0].decode("utf-8")
except Exception:
return False
# Apply calibrations if applicable
def apply_cals(image_data, gen_dark, gen_bias, gen_flat, i):
if gen_dark.size != 0:
if i == 0:
log_info("Dark subtracting images.")
image_data = image_data - gen_dark
elif gen_bias.size != 0: # if a dark is not available, then at least subtract off the pedestal via the bias
if i == 0:
log_info("Bias-correcting images.")
image_data = image_data - gen_bias
else:
pass
if gen_flat.size != 0:
if i == 0:
log_info("Flattening images.")
gen_flat[gen_flat == 0] = 1
image_data = image_data / gen_flat
return image_data
# Aligns imaging data from .fits file to easily track the host and comparison star's positions
def transformation(image_data, file_name, roi=1):
# crop image to ROI
height = image_data.shape[1]
width = image_data.shape[2]
roix = slice(int(width * (0.5 - roi / 2)), int(width * (0.5 + roi / 2)))
roiy = slice(int(height * (0.5 - roi / 2)), int(height * (0.5 + roi / 2)))
# Find transformation from .FITS files and catch exceptions if not able to.
try:
results = aa.find_transform(image_data[1][roiy, roix], image_data[0][roiy, roix])
return results[0]
except Exception as ee:
log_info(ee)
ws = 5
# smooth image and try to align again
windows = view_as_windows(image_data[0], (ws,ws), step=1)
medimg = np.median(windows, axis=(2,3))
windows = view_as_windows(image_data[1], (ws,ws), step=1)
medimg1 = np.median(windows, axis=(2,3))
try:
results = aa.find_transform(medimg1[roiy, roix], medimg[roiy, roix])
return results[0]
except Exception as ee:
log_info(ee)
log_info(file_name)
for p in [99, 98, 95, 90]:
for it in [2, 1, 0]:
# create binary mask to align image
mask1 = image_data[1][roiy, roix] > np.percentile(image_data[1][roiy, roix], p)
mask1 = binary_erosion(mask1, iterations=it)
mask0 = image_data[0][roiy, roix] > np.percentile(image_data[0][roiy, roix], p)
mask0 = binary_erosion(mask0, iterations=it)
try:
results = aa.find_transform(mask1, mask0)
return results[0]
except Exception as ee:
log_info(ee)
log_info(f"Warning: Alignment failed - {file_name}", warn=True)
return SimilarityTransform(scale=1, rotation=0, translation=[0, 0])
def get_pixel_scale(wcs_header, header, pixel_init):
astrometry_scale = None
if wcs_header:
astrometry_scale = [key.value.split(' ') for key in wcs_header._cards if 'scale:' in str(key.value)]
if astrometry_scale:
image_scale_num = astrometry_scale[0][1]
image_scale_units = astrometry_scale[0][2]
image_scale = f"Image scale in {image_scale_units}: {image_scale_num}"
elif 'IM_SCALE' in header:
image_scale_num = header['IM_SCALE']
image_scale_units = header.comments['IM_SCALE']
image_scale = f"Image scale in {image_scale_units}: {image_scale_num}"
elif 'PIXSCALE' in header:
image_scale_num = header['PIXSCALE']
image_scale_units = header.comments['PIXSCALE']
image_scale = f"Image scale in {image_scale_units}: {image_scale_num}"
elif pixel_init:
image_scale = f"Image scale in arcsecs/pixel: {pixel_init}"
else:
log_info("Not able to find Image Scale in the Image Header.")
image_scale_num = user_input("Please enter Image Scale (arcsec/pixel): ", type_=float)
image_scale = f"Image scale in arcsecs/pixel: {image_scale_num}"
return image_scale
def exp_time_med(exptimes):
# exposure time
consistent_et = False
if len(exptimes) > 0:
consistent_et = all(elem == exptimes[0] for elem in exptimes)
exptimes = np.array(exptimes)
if consistent_et:
return exptimes[0]
else:
return np.median(exptimes)
# finds target in WCS image after applying proper motion correction from SIMBAD
def find_target(target, hdufile, verbose=False):
# query simbad to get proper motions
service = vo.dal.TAPService("http://simbad.u-strasbg.fr/simbad/sim-tap")
# http://simbad.u-strasbg.fr/simbad/tap/tapsearch.html
query = '''
SELECT basic.OID, ra, dec, main_id, pmra, pmdec
FROM basic JOIN ident ON oidref = oid
WHERE id = '{}';
'''.format(target)
result = service.search(query)
# TODO check that simbad returned a value
# set up astropy object
coord = SkyCoord(
ra=result['ra'][0] * u.deg,
dec=result['dec'][0] * u.deg,
distance=1 * u.pc,
pm_ra_cosdec=result['pmra'][0] * u.mas / u.yr,
pm_dec=result['pmdec'][0] * u.mas / u.yr,
frame="icrs",
obstime=astropy.time.Time("2000-1-1T00:00:00")
)
hdu = fits.open(hdufile)[0]
try:
dateobs = hdu.header['DATE_OBS']
except:
dateobs = hdu.header['DATE']
# ignore timezone
if len(dateobs.split('-')) == 4:
dateobs = '-'.join(dateobs.split('-')[:-1])
t = astropy.time.Time(dateobs, format='isot', scale='utc')
coordpm = coord.apply_space_motion(new_obstime=t)
# wcs coordinate translation
wcs = WCS(hdu.header)
pixcoord = wcs.wcs_world2pix([[coordpm.ra.value, coordpm.dec.value]], 0)
if verbose:
print("Simbad:", result)
print("\nObs Date:", t)
print("NEW:", coordpm.ra, coordpm.dec)
print("\nTarget Location:", np.round(pixcoord[0], 2))
return pixcoord[0]
def gaussian_psf(x, y, x0, y0, a, sigx, sigy, rot, b):
rx = (x - x0) * np.cos(rot) - (y - y0) * np.sin(rot)
ry = (x - x0) * np.sin(rot) + (y - y0) * np.cos(rot)
gausx = np.exp(-rx ** 2 / (2 * sigx ** 2))
gausy = np.exp(-ry ** 2 / (2 * sigy ** 2))
return a * gausx * gausy + b
def mesh_box(pos, box, maxx=0, maxy=0):
pos = [int(np.round(pos[0])), int(np.round(pos[1]))]
if maxx:
x = np.arange(max(0,pos[0] - box), min(maxx, pos[0] + box + 1))
else:
x = np.arange(max(0,pos[0] - box), pos[0] + box + 1)
if maxy:
y = np.arange(max(0,pos[1] - box), min(maxy, pos[1] + box + 1))
else:
y = np.arange(max(0,pos[1] - box), pos[1] + box + 1)
xv, yv = np.meshgrid(x, y)
return xv.astype(int), yv.astype(int)
# Method fits a 2D gaussian function that matches the star_psf to the star image and returns its pixel coordinates
def fit_centroid(data, pos, psf_function=gaussian_psf, box=15, weightedcenter=True):
# get sub field in image
xv, yv = mesh_box(pos, box, maxx=data.shape[1], maxy=data.shape[0])
subarray = data[yv, xv]
init = [np.nanmax(subarray) - np.nanmin(subarray), 1, 1, 0, np.nanmin(subarray)]
# compute flux weighted centroid in x and y
wx = np.sum(xv[0]*subarray.sum(0))/subarray.sum(0).sum()
wy = np.sum(yv[:,0]*subarray.sum(1))/subarray.sum(1).sum()
# lower bound: [xc, yc, amp, sigx, sigy, rotation, bg]
lo = [pos[0] - box * 0.5, pos[1] - box * 0.5, 0, 0.5, 0.5, -np.pi / 4, np.nanmin(subarray) - 1]
up = [pos[0] + box * 0.5, pos[1] + box * 0.5, 1e7, 20, 20, np.pi / 4, np.nanmax(subarray) + 1]
def fcn2min(pars):
model = psf_function(xv, yv, *pars)
return (subarray - model).flatten()
try:
res = least_squares(fcn2min, x0=[*pos, *init], bounds=[lo, up], jac='3-point', xtol=None, method='trf')
except:
log_info(f"Warning: Measured flux amplitude is really low---are you sure there is a star at {np.round(pos, 2)}?", warn=True)
res = least_squares(fcn2min, x0=[*pos, *init], jac='3-point', xtol=None, method='lm')
# override psf fit results with weighted centroid
if weightedcenter:
res.x[0] = wx
res.x[1] = wy
return res.x
# Method calculates the flux of the star (uses the skybg_phot method to do background sub)
def aperPhot(data, xc, yc, r=5, dr=5):
if dr > 0:
bgflux, sigmabg, Nbg = skybg_phot(data, xc, yc, r + 2, dr)
else:
bgflux, sigmabg, Nbg = 0, 0, 0
aperture = CircularAperture(positions=[(xc, yc)], r=r)
mask = aperture.to_mask(method='exact')[0]
data_cutout = mask.cutout(data)
aperture_sum = (mask.data * (data_cutout - bgflux)).sum()
return aperture_sum, bgflux
def skybg_phot(data, xc, yc, r=10, dr=5, ptol=99, debug=False):
# create a crude annulus to mask out bright background pixels
xv, yv = mesh_box([xc, yc], np.round(r + dr))
rv = ((xv - xc) ** 2 + (yv - yc) ** 2) ** 0.5
mask = (rv > r) & (rv < (r + dr))
try:
cutoff = np.nanpercentile(data[yv, xv][mask], ptol)
except IndexError:
log_info(f"Warning: IndexError, problem computing sky bg for {xc:.1f}, {yc:.1f}."
f"\nCheck if star is present or close to border.", warn=True)
# create pixel wise mask on entire image
x = np.arange(data.shape[1])
y = np.arange(data.shape[0])
xv, yv = np.meshgrid(x, y)
rv = ((xv - xc) ** 2 + (yv - yc) ** 2) ** 0.5
mask = (rv > r) & (rv < (r + dr))
cutoff = np.nanpercentile(data[yv, xv][mask], ptol)
dat = np.array(data[yv, xv], dtype=float)