-
Notifications
You must be signed in to change notification settings - Fork 0
/
post_calcos.py
1640 lines (1511 loc) · 92 KB
/
post_calcos.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
'''
This is the FaintCOS v1.2 post_calcos script
Purpose:
1. Few percent-level accurate dark current subtraction of a collection of COS FUV exposures
(Worseck et al. 2016, ApJ, 825, 144; Makan et al. 2021, ApJ, 912, 38)
2. Estimation of scattered light from geocoronal Lyman alpha for G140L (Worseck et al. 2016)
3. Co-addition and rebinning of exposures from multiple COS FUV data sets
(several visits, both FUV M gratings, several central wavelengths)
4. Calculation of proper Poisson flux error bars (e.g. Worseck et al. 2016)
5. Generation of 1D spectrum with HLSP-compliant metadata (as of Apr 09, 2021,
https://outerspace.stsci.edu/display/MASTDOCS/HLSP+Contributor+Guide)
Input:
1. Path to science data (single target, either COS G140L *or* G130M/G160M gratings)
2. Optional path to target-specific configuration file faintcos_config.py,
otherwise take the default in the FaintCOS installation directory
Output:
1. MEF file with two FITS binary tables, extension 1: co-added and rebinned COS FUV spectrum,
extension 2: provenance metadata
2. If HLSP output is desired (see config file) the primary output spectrum is copied to an
hlsp* file following the HLSP naming convention
Note: The choice of the primary output spectrum depends on the number of data sets.
For a target with a single data set take the combined spectrum rebinned by an integer number of pixels.
For a target with multiple data sets, the individual data sets are rebinned onto a common wavelength
grid with constant dispersion, by default assuming that target variability and flux calibration
errors are negligible. For multiple data sets the rebinned spectra of individual data sets
are optional (see config file).
Authors: Kirill Makan, Gabor Worseck
'''
version = '1.2'
import os
import sys
import math
from astropy.io import fits
from astropy.table import Table, vstack, hstack, Column
from astropy.stats import poisson_conf_interval
from astropy.time import Time
from datetime import datetime, timezone
from scipy.interpolate import interp1d
from scipy.stats import linregress
import numpy as np
from shutil import copy2
def calc_conf_lim_feldman_cousins(counts, bkg):
'''
Purpose:
Calculation of two-sided equal-tailed 68.26% confidence (1 sigma) statistical uncertainty for Poisson
counts accounting for known background (frequentist method: Feldman & Cousins 1998, Phys.Rev. D., 57, 3873).
The confidence interval is for the signal, i.e. counts minus background. If counts < background
compute a one-sided uncertainty from the one-sided 1 sigma upper limit (84.13% c.l.).
Feldman & Cousins call this a sensitivity limit (their Section 6).
Parameters:
1. counts: numpy array with (integer) Poisson counts (signal + background)
2. bkg: numpy array with (non-integer) background
Returns:
cnt_err_down, cnt_err_up: two numpy arrays with lower and upper statistical uncertainty,
the sum of both corresponds to a 68.26% confidence level (1 sigma)
'''
# in order to install CustomConfLim, run 'python setup.py build' and
# 'python setup.py install' in the CustomConfLimits folder
import CustomConfLim
cnt_err_up = np.zeros(shape=len(counts), dtype=np.float32)
cnt_err_down = np.zeros(shape=len(counts), dtype=np.float32)
print("Calculating stat. flux errors (frequentist method: Feldman & Cousins 1998):")
for i in range(len(counts)):
sys.stdout.write("|")
for p in range(50):
if i > p*len(counts)/50.:
sys.stdout.write("=")
else:
sys.stdout.write(" ")
sys.stdout.write("| " + \
str(round(100.*float(i)/len(counts))) + \
" %\r")
sys.stdout.flush()
if counts[i] > 150:
limits = (math.sqrt(counts[i]), \
math.sqrt(counts[i]))
else:
if bkg[i] > counts[i]:
limits = CustomConfLim.bootstrap_bkg_conf_lim(counts[i], \
bkg[i], 10000)
elif bkg[i] > 0.0:
limits = \
CustomConfLim.feldman_cousins_conf_lim(round(counts[i]), \
bkg[i])
else:
limits = (0.0, 0.0)
cnt_err_down[i] = limits[0]
cnt_err_up[i] = limits[1]
if cnt_err_down[i] < 0:
cnt_err_down[i] = 0.0
if cnt_err_up[i] < 0:
cnt_err_up[i] = 0.0
print("\n")
return cnt_err_down, cnt_err_up
def calc_conf_lim_kraft(counts, bkg):
'''
Purpose:
Calculation of two-sided minimal 68.26% integrated posterior density (1 sigma) statistical uncertainty
for Poisson counts accounting for known background (Bayesian method: Kraft et al. 1991, ApJ, 374, 344;
implemented in astropy.stats module). The credible interval is for the signal, i.e. counts minus background,
and corresponds to the posterior maximum. The statistical uncertainty is calculated from the minimal
68.27% credible interval (as in Kraft et al., hopefully).
WARNING: Most astronomers prefer/assume equal-tailed credible intervals. The difference is subtle,
but publications should be specific.
Parameters:
1. counts: numpy array with (integer) Poisson counts (signal + background)
2. bkg: numpy array with (non-integer) background
Returns:
cnt_err_down, cnt_err_up: two numpy arrays with lower and upper statistical uncertainty,
the sum of both corresponds to a 68.26% credibility (1 sigma)
'''
cnt_err_up = np.zeros(shape=len(counts), dtype=np.float32)
cnt_err_down = np.zeros(shape=len(counts), dtype=np.float32)
print("Calculating stat. flux errors (Bayesian method: Kraft et al. 1991):")
for i in range(len(counts)):
sys.stdout.write("|")
for p in range(50):
if i > p*len(counts)/50.:
sys.stdout.write("=")
else:
sys.stdout.write(" ")
sys.stdout.write("| " + \
str(round(100.*float(i)/len(counts))) + \
" %\r")
sys.stdout.flush()
if counts[i] > 99:
cnt_err_up[i] = math.sqrt(counts[i] - bkg[i])
cnt_err_down[i] = cnt_err_up[i]
else:
try:
limits = poisson_conf_interval(counts[i], \
background=bkg[i], \
confidence_level=0.6827, \
interval='kraft-burrows-nousek')
if counts[i] > bkg[i]:
cnt_err_down[i] = counts[i] - bkg[i] - limits[0]
else:
cnt_err_down[i] = 0.0
cnt_err_up[i] = limits[1] - (counts[i] - bkg[i])
except ValueError:
cnt_err_down[i] = 0.0
cnt_err_up[i] = 0.0
if cnt_err_down[i] < 0:
cnt_err_down[i] = 0.0
if cnt_err_up[i] < 0:
cnt_err_up[i] = 0.0
print("\n")
return cnt_err_down, cnt_err_up
def calc_lya_scatter_model(gcounts, wave_arr, dq):
'''
Purpose:
Calculation of scattered geocoronal Lyman alpha emission in COS G140L spectra
using the simple model (Gaussian profile scaled to maximum Lyman alpha counts)
from Worseck et al. 2016, ApJ, 825, 144.
The maximum of the geocoronal Lyman alpha counts is estimated from the mean
at 1215.37-1215.97A (broad unresolved profile).
Parameters:
1. gcounts: numpy array with 1D extracted count spectrum
2. wave_arr: numpy array with wavelength grid in Angstroem
Returns:
1. Blya: numpy array with estimated scattered light profile given in counts
2. Blya_err_up: numpy array with upper error of Blya
3. Blya_err_down: numpy array with lower error of Blya
The sum of both errors gives the 1 sigma (68% confidence level) uncertainty
of the scattered light given the scarce data (estimated by parametric bootstrap
in Worseck et al. 2016). The error is likely underestimated, because it assumes
a rigid Gaussian shape for the scattered light profile. When calculating the
total background uncertainty, this error must be treated as a systematic error.
'''
a = 1.7348e-5
lam_0 = 1254.6
b = 100.9
# Relative upper error of the Lya contamination
rel_err = np.array([[ 801., 853., 905., 957., 1009., 1061., 1113.,\
1165., 1217., 1269., 1321., 1373., 1425., 1477.,\
1529., 1581., 1633., 1685., 1737., 1789., 1841.,\
1893., 1945., 1997.],\
[2.976524 ,2.029496 ,1.3314583 ,0.8864973 ,\
0.5695007 ,0.34901363 ,0.1948071 ,0.09911883,\
0.05120203 ,0.04276433 ,0.04487024 ,0.05724196,\
0.11888563 ,0.23237737 ,0.40143612 ,0.6552495 ,\
1.0261706 ,1.4999093 ,2.2647116 ,3.3137312 ,\
4.849577 ,7.3380733 ,11.06501 ,16.91639 ],\
[0.7585422 ,0.67552364 ,0.57972205 ,0.4754494 ,\
0.37029144 ,0.2608441 ,0.169427 ,0.09658019,\
0.05248015 ,0.04616276 ,0.04640767 ,0.06427228,\
0.11368275 ,0.19347219 ,0.28874132 ,0.4030617 ,\
0.51270425 ,0.613789 ,0.70610535 ,0.78464127,\
0.84194565 ,0.8906022 ,0.92548954 ,0.9521525 ]])
err_up_sys = np.interp(wave_arr, rel_err[0], rel_err[1])
err_down_sys = np.interp(wave_arr, rel_err[0], rel_err[2])
CountsInLyaRegion = gcounts[(wave_arr <= 1215.97) & (wave_arr >= 1215.37) & (dq == 0)]
if len(CountsInLyaRegion)<3:
print('Warning: Grid wire on Lyman alpha. Peak determined from 5A region.')
CountsInLyaRegion = gcounts[(wave_arr <= 1218.17) & (wave_arr >= 1213.17) & (dq == 0)]
Clya = np.percentile(CountsInLyaRegion,75)
else:
Clya = np.mean(CountsInLyaRegion)
Blya = a*Clya*np.exp(-np.power(wave_arr - lam_0, 2.)/(2.*b**2))
Clya_err_stat = math.sqrt(Clya)/math.sqrt(len(CountsInLyaRegion))/Clya
err_up = np.sqrt(err_up_sys**2 + Clya_err_stat**2)
err_down = np.sqrt(err_down_sys**2 + Clya_err_stat**2)
return Blya, err_up*Blya, err_down*Blya
def createCDRPrimaryHeader(hdul, wave):
'''
Purpose:
Create custom primary header for FITS file with 1D spectrum of single subexposure,
i.e. a former CALCOS _x1d file after including the improved dark current estimate
and scattered light (for G140L). The FITS header is almost HLSP-compliant and
includes the relevant metadata for a COS FUV 1D spectrum.
Parameters:
1. hdul: HDU list from the CALCOS _x1d file
2. wave: numpy array with wavelengths of the 1D spectrum in Angstroem
Returns:
hdu: custom header for primary HDU of the cdr file
'''
hdu = fits.PrimaryHDU()
hdu.header['DATE'] = (datetime.utcnow().isoformat(timespec='seconds'),\
"File creation date")
hdu.header['FILETYPE'] = (hdul[0].header['FILETYPE'],
'type of data found in data file')
hdu.header['TELESCOP'] = (hdul[0].header['TELESCOP'],
'telescope used to acquire data')
hdu.header['INSTRUME'] = (hdul[0].header['INSTRUME'],
'identifier for instrument used to acquire data')
hdu.header['EQUINOX'] = (hdul[0].header['EQUINOX'],
'equinox of celestial coord. system')
hdu.header['TARGNAME'] = (hdul[0].header['TARGNAME'],
'proposer\'s target name')
hdu.header['RA_TARG'] = (hdul[0].header['RA_TARG'],
'[deg] right ascention of the target')
hdu.header['DEC_TARG'] = (hdul[0].header['DEC_TARG'],
'[deg] declination of the target')
hdu.header['PROPOSID'] = (hdul[0].header['PROPOSID'],
'PEP proposal identifier')
hdu.header['OPUS_VER'] = (hdul[0].header['OPUS_VER'],
'data processing software system version')
hdu.header['CSYS_VER'] = (hdul[0].header['CSYS_VER'],
'calibration software system version id')
hdu.header['CAL_VER'] = (hdul[0].header['CAL_VER'],
'CALCOS code version')
hdu.header['FCOS_VER'] = (version,
'FaintCOS code version')
hdu.header['DETECTOR'] = (hdul[0].header['DETECTOR'],
'FUV or NUV')
hdu.header['SEGMENT'] = (hdul[0].header['SEGMENT'],
'FUV detector segment name (FUVA, FUVB or BOTH)')
hdu.header['LIFE_ADJ'] = (hdul[0].header['LIFE_ADJ'],
'Life Time Adjustment Position')
hdu.header['APERTURE'] = (hdul[0].header['APERTURE'],
'aperture name')
hdu.header['OPT_ELEM'] = (hdul[0].header['OPT_ELEM'],
'optical element in use')
hdu.header['CENWAVE'] = (hdul[0].header['CENWAVE'],
'[Angstrom] grating central wavelength')
hdu.header['BANDWID'] = (max(wave)-min(wave),
'[Angstrom] bandwidth of the data')
res = cos_res[(cos_res['LP'] == hdu.header['LIFE_ADJ']) & \
(cos_res['OPT_ELEM'] == hdu.header['OPT_ELEM']) &\
(cos_res['CENWAVE'] == hdu.header['CENWAVE'])]['R'][0]
hdu.header['SPECRES'] = (res,
'approx. resolving power at CENWAVE')
hdu.header['CENTRWV'] = ((max(wave)+min(wave))/2.0,
'[Angstrom] central wavelength of the data')
hdu.header['MINWAVE'] = (min(wave),
'[Angstrom] minimum wavelength in spectrum')
hdu.header['MAXWAVE'] = (max(wave),
'[Angstrom] maximum wavelength in spectrum')
hdu.header['TIMESYS'] = ('UTC ', 'time scale of time-related keywords')
hdu.header['DATE-BEG'] = (hdul[1].header['DATE-OBS']+'T'+hdul[1].header['TIME-OBS'], 'exposure start time (ISO-8601 DateTime)')
t_end = Time(hdul[1].header['EXPEND'], format='mjd', scale='utc', precision=0)
hdu.header['DATE-END'] = (t_end.isot, 'exposure end time (ISO-8601 DateTime)')
hdu.header['MJD-BEG'] = (hdul[1].header['EXPSTART'], '[d] exposure start time (Modified Julian Date)')
hdu.header['MJD-END'] = (hdul[1].header['EXPEND'], '[d] exposure end time (Modified Julian Date)')
hdu.header['EXPTIME'] = (hdul[1].header['EXPTIME'], '[d] exposure duration (seconds)--calculated')
hdu.header['ASN_ID'] = (hdul[0].header['ASN_ID'],
'unique identifier assigned to association')
return hdu
def createDatasetPrimaryHeader(headers, wave):
'''
Purpose:
Create custom primary header for FITS file with coadded calibrated spectrum of COS data set,
i.e. from a collection of exposures with the same grating setting from one _asn association file.
The FITS header includes the relevant metadata for a COS FUV 1D spectrum in HLSP-compliant keywords.
HLSP keywords are added if HLSP output is desired.
Parameters:
1. headers: array with primary HDUs of individual exposures from _asn association file
2. wave: numpy array with wavelengths of the coadded 1D spectrum in Angstroem
Returns:
hdu: custom header for primary HDU of the coadded spectrum from a single data set.
This is written into the _dataset_sum file (unbinned spectrum) and the
_bin file (spectrum binned by given number of original pixels)
'''
hdu = fits.PrimaryHDU()
hdu.header['DATE'] = (datetime.utcnow().isoformat(timespec='seconds'),\
"File creation date")
hdu.header['FILETYPE'] = (headers[0]['FILETYPE'],
'type of data found in data file')
hdu.header['OBSERVAT'] = (headers[0]['TELESCOP'],
'observatory used to obtain observation')
hdu.header['TELESCOP'] = (headers[0]['TELESCOP'],
'telescope used to acquire data')
hdu.header['INSTRUME'] = (headers[0]['INSTRUME'],
'identifier for instrument used to acquire data')
hdu.header['EQUINOX'] = (headers[0]['EQUINOX'],
'equinox of celestial coord. system')
hdu.header['TARGNAME'] = (headers[0]['TARGNAME'],
'proposer\'s target name')
hdu.header['RA_TARG'] = (headers[0]['RA_TARG'],
'[deg] right ascention of the target')
hdu.header['DEC_TARG'] = (headers[0]['DEC_TARG'],
'[deg] declination of the target')
hdu.header['PROPOSID'] = (headers[0]['PROPOSID'],
'PEP proposal identifier')
hdu.header['OPUS_VER'] = (headers[0]['OPUS_VER'],
'data processing software system version')
hdu.header['CSYS_VER'] = (headers[0]['CSYS_VER'],
'calibration software system version id')
hdu.header['CAL_VER'] = (headers[0]['CAL_VER'],
'CALCOS code version')
hdu.header['FCOS_VER'] = (version,
'FaintCOS code version')
hdu.header['DETECTOR'] = (headers[0]['DETECTOR'],
'FUV or NUV')
hdu.header['SEGMENT'] = (headers[0]['SEGMENT'],
'FUV detector segment name (FUVA, FUVB or BOTH)')
hdu.header['LIFE_ADJ'] = (headers[0]['LIFE_ADJ'],
'Life Time Adjustment Position')
hdu.header['APERTURE'] = (headers[0]['APERTURE'],
'aperture name')
hdu.header['DISPERSR'] = (headers[0]['OPT_ELEM'],
'name of dispersive element used')
hdu.header['CENWAVE'] = (headers[0]['CENWAVE'],
'[Angstrom] grating central wavelength')
hdu.header['BANDWID'] = (max(wave)-min(wave),
'[Angstrom] bandwidth of the data')
res = cos_res[(cos_res['LP'] == headers[0]['LIFE_ADJ']) & \
(cos_res['OPT_ELEM'] == headers[0]['OPT_ELEM']) &\
(cos_res['CENWAVE'] == headers[0]['CENWAVE'])]['R'][0]
hdu.header['SPECRES'] = (res,
'approx. resolving power at CENWAVE')
hdu.header['CENTRWV'] = ((max(wave)+min(wave))/2.0,
'[Angstrom] central wavelength of the data')
hdu.header['MINWAVE'] = (min(wave),
'[Angstrom] minimum wavelength in spectrum')
hdu.header['MAXWAVE'] = (max(wave),
'[Angstrom] maximum wavelength in spectrum')
hdu.header['BINNING'] = (BIN_PX, '[pixel] binning factor w.r.t. orig. pixel size')
hdu.header['BUNIT'] = ('erg s^-1 cm^-2 A^-1', 'brightness unit')
date_beg = []
mjd_beg = []
mjd_end = []
exptime = 0
for hdr in headers:
if hdr['DATE-BEG'] not in date_beg:
date_beg.append(hdr['DATE-BEG'])
mjd_beg.append(hdr['MJD-BEG'])
mjd_end.append(hdr['MJD-END'])
exptime += hdr['EXPTIME']
sorted_dates = [x for _,x in sorted(zip(mjd_beg, date_beg))]
hdu.header['TIMESYS'] = (headers[0]['TIMESYS'], 'time scale of time-related keywords')
hdu.header['DATE-BEG'] = (sorted_dates.pop(0), 'ISO-8601 DateTime of first exposure start')
t_end = Time(max(mjd_end), format='mjd', scale='utc', precision=0)
hdu.header['DATE-END'] = (t_end.isot, 'ISO-8601 DateTime of last exposure end')
hdu.header['MJD-BEG'] = (min(mjd_beg), '[d] Mod. Julian Date of first exposure start')
hdu.header['MJD-END'] = (max(mjd_end), '[d] Mod. Julian Date of last exposure end')
hdu.header['MJD-MID'] = (0.5*(max(mjd_end)+min(mjd_beg)), '[d] MJD mid-time -- likely meaningless')
hdu.header['XPOSURE'] = (exptime, '[s] observation duration --calculated')
hdu.header['FILE_ID'] = (headers[0]['ASN_ID'], 'Dataset identifier')
if HLSP_write:
hdu.header['HLSPTARG'] = (headers[0]['TARGNAME'], 'HLSP target designation')
hdu.header['HLSPID'] = (HLSP_id, 'HLSP identifier (acronym)')
hdu.header['HLSPNAME'] = (HLSP_name, 'title for HLSP project')
hdu.header['HLSPLEAD'] = (HLSP_lead, 'full name of HLSP project lead')
hdu.header['HLSPVER'] = (HLSP_ver, 'version identifier for HLSP product')
hdu.header['DOI'] = (HLSP_doi, 'HLSP Digital Object Identifier')
hdu.header['REFERENC'] = (HLSP_referenc, 'bibliogr. identifier (ADS bibcode)')
hdu.header['LICENSE'] = ('CC BY 4.0', 'license for use of data')
hdu.header['LICENURL'] = ('https://creativecommons.org/licenses/by/4.0/', 'Data license URL')
return hdu
def createCoAddPrimaryHeader(headers, wave, exptime):
'''
Purpose:
Create custom primary header for FITS file with coadded and rebinned spectrum from several COS data sets,
i.e. from a collection of exposures from several _asn association files for the same target
The FITS header includes the relevant metadata for a COS FUV 1D spectrum in HLSP-compliant keywords.
HLSP keywords are added if HLSP output is desired.
Parameters:
1. headers: array with primary HDUs of individual exposures
2. wave: numpy array with wavelengths of the coadded 1D spectrum in Angstroem
3. exptime: numpy array with pixel exposure time in seconds
Returns:
hdu: custom header for primary HDU of the coadded and rebinned spectrum (_spectrum file).
'''
unique_pid = []
unique_lifeadj = []
unique_aper = []
unique_disp = []
unique_cenwave = []
for hdr in headers:
if hdr['PROPOSID'] not in unique_pid:
unique_pid.append(hdr['PROPOSID'])
if hdr['LIFE_ADJ'] not in unique_lifeadj:
unique_lifeadj.append(hdr['LIFE_ADJ'])
if hdr['APERTURE'] not in unique_aper:
unique_aper.append(hdr['APERTURE'])
if hdr['OPT_ELEM'] not in unique_disp:
unique_disp.append(hdr['OPT_ELEM'])
if hdr['CENWAVE'] not in unique_cenwave:
unique_cenwave.append(hdr['CENWAVE'])
hdu = fits.PrimaryHDU()
hdu.header['DATE'] = (datetime.utcnow().isoformat(timespec='seconds'),\
"File creation date")
hdu.header['FILETYPE'] = (headers[0]['FILETYPE'],
'type of data found in data file')
hdu.header['OBSERVAT'] = (headers[0]['TELESCOP'],
'observatory used to obtain observation')
hdu.header['TELESCOP'] = (headers[0]['TELESCOP'],
'telescope used to acquire data')
hdu.header['INSTRUME'] = (headers[0]['INSTRUME'],
'identifier for instrument used to acquire data')
hdu.header['EQUINOX'] = (headers[0]['EQUINOX'],
'equinox of celestial coord. system')
hdu.header['TARGNAME'] = (headers[0]['TARGNAME'],
'proposer\'s target name')
hdu.header['RA_TARG'] = (headers[0]['RA_TARG'],
'[deg] right ascention of the target')
hdu.header['DEC_TARG'] = (headers[0]['DEC_TARG'],
'[deg] declination of the target')
if len(unique_pid) == 1:
hdu.header['PROPOSID'] = (unique_pid[0], 'PEP proposal identifier')
else:
hdu.header['PROPOSID'] = ('MULTI', 'PEP proposal identifier')
hdu.header['OPUS_VER'] = (headers[0]['OPUS_VER'],
'data processing software system version')
hdu.header['CSYS_VER'] = (headers[0]['CSYS_VER'],
'calibration software system version id')
hdu.header['CAL_VER'] = (headers[0]['CAL_VER'],
'CALCOS code version')
hdu.header['FCOS_VER'] = (version,
'FaintCOS code version')
hdu.header['DETECTOR'] = (headers[0]['DETECTOR'],
'FUV or NUV')
segments = [hdr['SEGMENT'] for hdr in headers]
segm = ""
if 'FUVA' in segments and 'FUVB' in segments:
segm = 'BOTH'
else:
segm = segments[0]
hdu.header['SEGMENT'] = (segm, 'FUV detector segment name (FUVA, FUVB or BOTH)')
if len(unique_lifeadj) == 1:
hdu.header['LIFE_ADJ'] = (unique_lifeadj[0], 'Life Time Adjustment Position')
else:
hdu.header['LIFE_ADJ'] = ('MULTI', 'Life Time Adjustment Position')
if len(unique_aper) == 1:
hdu.header['APERTURE'] = (unique_aper[0], 'aperture name')
else:
hdu.header['APERTURE'] = ('MULTI', 'aperture name')
if len(unique_disp) == 1:
hdu.header['DISPERSR'] = (unique_disp[0], 'name of dispersive element used')
else:
hdu.header['DISPERSR'] = ('MULTI', 'name of dispersive element used')
if len(unique_cenwave) == 1:
hdu.header['CENWAVE'] = (unique_cenwave[0], '[Angstrom] grating central wavelength')
else:
hdu.header['CENWAVE'] = ('MULTI', 'grating central wavelength')
hdu.header['BANDWID'] = (max(wave)-min(wave),
'[Angstrom] bandwidth of the data')
hdu.header['SPECRES'] = (min([hdr['SPECRES'] for hdr in headers]),
'min. res. power at CENWAVE from all data sets')
hdu.header['CENTRWV'] = ((max(wave)+min(wave))/2.0,
'[Angstrom] central wavelength of the data')
hdu.header['MINWAVE'] = (min(wave),
'[Angstrom] minimum wavelength in spectrum')
hdu.header['MAXWAVE'] = (max(wave),
'[Angstrom] maximum wavelength in spectrum')
hdu.header['BINNING'] = (BIN_SIZE, '[Angstrom] bin size of wavelength axis')
hdu.header['BUNIT'] = ('erg s^-1 cm^-2 A^-1', 'brightness unit')
hdu.header['TIMESYS'] = (headers[0]['TIMESYS'], 'time scale of time-related keywords')
date_beg = []
mjd_beg = []
mjd_end = []
tottime = 0
for hdr in headers:
if hdr['DATE-BEG'] not in date_beg:
date_beg.append(hdr['DATE-BEG'])
mjd_beg.append(hdr['MJD-BEG'])
mjd_end.append(hdr['MJD-END'])
tottime += hdr['EXPTIME']
sorted_dates = [x for _,x in sorted(zip(mjd_beg, date_beg))]
hdu.header['DATE-BEG'] = (sorted_dates.pop(0), 'ISO-8601 DateTime of first exposure start')
t_end = Time(max(mjd_end), format='mjd', scale='utc', precision=0)
hdu.header['DATE-END'] = (t_end.isot, 'ISO-8601 DateTime of last exposure end')
hdu.header['MJD-BEG'] = (min(mjd_beg), '[d] Mod. Julian Date of first exposure start')
hdu.header['MJD-END'] = (max(mjd_end), '[d] Mod. Julian Date of last exposure end')
hdu.header['MJD-MID'] = (0.5*(max(mjd_end)+min(mjd_beg)), '[d] MJD mid-time -- likely meaningless')
hdu.header['XPOSURE'] = (round(tottime,3), '[s] total exposure time of data sets')
hdu.header['MXPOSURE'] = (round(float(max(exptime)),3), '[s] max. exposure time in spectral range')
if HLSP_write:
hdu.header['HLSPTARG'] = (headers[0]['TARGNAME'], 'HLSP target designation')
hdu.header['HLSPID'] = (HLSP_id, 'HLSP identifier (acronym)')
hdu.header['HLSPNAME'] = (HLSP_name, 'title for HLSP project')
hdu.header['HLSPLEAD'] = (HLSP_lead, 'full name of HLSP project lead')
hdu.header['HLSPVER'] = (HLSP_ver, 'version identifier for HLSP product')
hdu.header['DOI'] = (HLSP_doi, 'HSLP Digital Object Identifier')
hdu.header['REFERENC'] = (HLSP_referenc, 'bibliogr. identifier (ADS bibcode)')
hdu.header['LICENSE'] = ('CC BY 4.0', 'license for use of data')
hdu.header['LICENURL'] = ('https://creativecommons.org/licenses/by/4.0/', 'Data license URL')
return hdu
if __name__ == "__main__":
print("FaintCOS v"+version, flush=True)
path_sci = "."
# Get paths to science data and FaintCOS config file from command line.
# If 2nd path to config file is given then frontload it into the search path for imports to force local import.
# Then import FaintCOS config file
if len(sys.argv) == 3:
path_sci = sys.argv[1]
sys.path.insert(0, sys.argv[2])
else:
if len(sys.argv) == 2:
path_sci = sys.argv[1]
from faintcos_config import *
# find all corrtag files in the science directory
path_corrtag = [f for f in os.listdir(path_sci) if "corrtag" in f]
if path_sci != ".":
path_corrtag = [path_sci + s for s in path_corrtag]
else:
path_sci = ""
# load all relevant metadata from corrtags and datasets into astropy tables and sort them
# dataset table is built from unique dataset IDs in corrtag files
datasets = Table(names=("PROPOSID", "FILE_ID", "TARGNAME", "DATE-OBS", "DATE-BEG", "DATE-END", "MJD-BEG", "MJD-END", \
"LIFE_ADJ", "APERTURE", "DISPERSR", "CENWAVE", "XPOSURE", "EXPA", "EXPB"),\
dtype=('i4', 'S10', 'S25', 'S10', 'S19', 'S19', 'f8', 'f8', 'i4', 'S4', 'S5', 'i4', 'f8', 'i4', 'i4'))
corrtags = Table(names=("TARGNAME", "FILE_ID", "CORRTAG_FILE", "DATE-BEG",\
"DISPERSR", "SEGMENT", "CENWAVE",\
"EXPTIME"), \
dtype=('S25', 'S10', 'S30', 'S19', 'S5' , 'S5', 'i4', 'f8'))
datasets['XPOSURE'].format = '4.3f'
corrtags['EXPTIME'].format = '4.3f'
for f in path_corrtag:
h0 = fits.open(f)[0].header
h1 = fits.open(f)[1].header
if h1['EXPFLAG'] != 'NORMAL':
continue
corrtags.add_row([h0['TARGNAME'], h0['ASN_ID'], f.split('/')[-1], \
h1['DATE-OBS'] + "T" + h1['TIME-OBS'],\
h0['OPT_ELEM'], h0['SEGMENT'], \
h0['CENWAVE'], h1['EXPTIME']])
corrtags.sort(["TARGNAME", 'FILE_ID', "DATE-BEG", 'SEGMENT'])
unique_datasets = np.unique(np.array(corrtags['FILE_ID']))
for asn_id in unique_datasets:
visit = corrtags[corrtags['FILE_ID'] == asn_id]
hdul = fits.open(path_sci + visit[len(visit)-1]['CORRTAG_FILE'])
mjd_end = hdul[1].header['EXPEND']
t_end = Time(mjd_end, format='mjd', scale='utc', precision=0)
date_end = t_end.isot
hdul.close()
hdul = fits.open(path_sci + visit[0]['CORRTAG_FILE'])
propid = hdul[0].header['PROPOSID']
target = hdul[0].header['TARGNAME']
date_obs = hdul[1].header['DATE-OBS']
date_beg = hdul[1].header['DATE-OBS']+"T"+hdul[1].header['TIME-OBS']
mjd_beg = hdul[1].header['EXPSTART']
life_adj = hdul[0].header['LIFE_ADJ']
aper = hdul[0].header['APERTURE']
hdul.close()
exp_time_a = np.array(visit[visit['SEGMENT'] == 'FUVA']['EXPTIME'])
exp_time_b = np.array(visit[visit['SEGMENT'] == 'FUVB']['EXPTIME'])
exp_time = max(np.array([np.sum(exp_time_a), np.sum(exp_time_b)]))
num_fuva = len(np.array(visit[visit['SEGMENT'] == 'FUVA']))
num_fuvb = len(np.array(visit[visit['SEGMENT'] == 'FUVB']))
cenwave = visit[0]['CENWAVE']
opt_elem = visit[0]['DISPERSR']
datasets.add_row([propid, asn_id, target, date_obs, date_beg, date_end, mjd_beg, mjd_end, \
life_adj, aper, opt_elem, cenwave, exp_time, num_fuva, num_fuvb])
datasets.sort(['PROPOSID','FILE_ID','DATE-OBS'])
# reduce individual exposures as in Worseck et al. 2016: For every exposure estimate the dark current in
# science extraction aperture from dark frames with similar pulse height distribution (sensitive to environmental conditions)
# as in offset calibration windows in science exposure _corrtag file. For G140L exposures also estimate the scattered light
if REDUCE_EXPOSURES:
# Get directory with CALCOS calibration files
try:
path_ref = os.environ['lref']
except KeyError:
print("ERROR: lref is not defined!")
sys.exit()
# Get directory with dark frames (_corrtag files, ideally reduced with the same CALCOS version and calibration files)
try:
path_dark = os.environ['ldark']
except KeyError:
print("ERROR: ldark is not defined!")
sys.exit()
# load all dark frames in the 'ldark' directory
print("Loading dark frames...", end=" ", flush=True)
path_darkframes = [f for f in os.listdir(path_dark) if "corrtag" in f]
path_darkframes = [path_dark + s for s in path_darkframes]
dark_file = []
dark_expstart = []
dark_segment = []
dark_voltage = []
for d in path_darkframes:
dark_tmp = fits.open(d)
dark_file.append(dark_tmp[0].header['FILENAME'])
dark_expstart.append(dark_tmp[1].header['EXPSTART'])
dark_segment.append(dark_tmp[0].header['SEGMENT'])
dark_voltage.append(dark_tmp[1].header["HVLEVEL" + \
dark_segment[-1].split("FUV")[1]])
dark_tmp.close()
darkframes = Table([dark_file, path_darkframes, \
dark_expstart, dark_segment, dark_voltage],\
names = ("FILE", "PATH", "EXPSTART", \
"SEGMENT", "VOLTAGE"))
print("OK")
print(str(len(darkframes['FILE'])) + " dark frames have been found!")
# print sorted lists of corrtag files and datasets in terminal
print("Valid 'corrtag' files in the input directory:")
corrtags['FILE_ID','TARGNAME','CORRTAG_FILE','CENWAVE','EXPTIME'].pprint()
print("##################################################################")
print("\n")
print("Valid datasets in the input directory:")
datasets['FILE_ID','TARGNAME','DATE-OBS','XPOSURE','DISPERSR','CENWAVE','EXPA','EXPB'].pprint()
print("###################################################################")
print("\n")
a = input("Do you wish to proceed? (y/n)")
if a != 'y' and a != 'Y':
print("Canceled by user!")
sys.exit()
# determine background (dark current and scattered light) for every corrtag
corr_files = [path_sci + s for s in list(corrtags['CORRTAG_FILE'])]
for c in corr_files:
corr_hdul = fits.open(c)
corr_data = Table(corr_hdul[1].data)
corr_prefix = corr_hdul[0].header['FILENAME'].split('_')[0]
print("Working on " + str(c.split("/")[-1]))
segm = corr_hdul[0].header['SEGMENT']
opt_elem = corr_hdul[0].header['OPT_ELEM']
cenwave = corr_hdul[0].header['CENTRWV']
voltage = corr_hdul[1].header["HVLEVEL" + segm.split("FUV")[1]]
exp_start = corr_hdul[1].header['EXPSTART']
# open xtractab to get apertures
xtractab = Table(fits.getdata(path_ref + \
corr_hdul[0].header['XTRACTAB'].split('$')[-1]))
xtractab = xtractab[(xtractab['SEGMENT'] == segm) &\
(xtractab['OPT_ELEM'] == opt_elem) &\
(xtractab['CENWAVE'] == cenwave) &\
(xtractab['APERTURE'] == 'PSA')]
ap_spec = [float(xtractab['B_SPEC']) - float(xtractab['HEIGHT'])/2.,\
float(xtractab['B_SPEC']) + float(xtractab['HEIGHT'])/2.]
if 'B_HGT1' in xtractab.columns:
ap_bkg = [float(xtractab['B_BKG1']) - float(xtractab['B_HGT1'])/2.,\
float(xtractab['B_BKG1']) + float(xtractab['B_HGT1'])/2.,\
float(xtractab['B_BKG2']) - float(xtractab['B_HGT2'])/2.,\
float(xtractab['B_BKG2']) + float(xtractab['B_HGT2'])/2.]
else:
ap_bkg = [float(xtractab['B_BKG1']) - float(xtractab['BHEIGHT'])/2.,\
float(xtractab['B_BKG1']) + float(xtractab['BHEIGHT'])/2.,\
float(xtractab['B_BKG2']) - float(xtractab['BHEIGHT'])/2.,\
float(xtractab['B_BKG2']) + float(xtractab['BHEIGHT'])/2.]
# open pulse height calibration file to get pulse height limits
# these limits should be identical to those given in the FaintCOS config file
pha_file = path_ref + corr_hdul[0].header['PHATAB'].split('$')[-1]
pha_limits = Table(fits.open(pha_file)[1].data)
pha_max = pha_limits[(pha_limits['OPT_ELEM'] == opt_elem) & \
(pha_limits['SEGMENT'] == segm)]['ULT']
pha_min = pha_limits[(pha_limits['OPT_ELEM'] == opt_elem) & \
(pha_limits['SEGMENT'] == segm)]['LLT']
# open the uncalibrated 1D spectrum produced by CALCOS (_x1d)
path_x1d = path_sci + corr_prefix + "_x1d.fits"
hdul_x1d = fits.open(path_x1d)
data_x1d = Table(hdul_x1d[1].data)
data_x1d = data_x1d[data_x1d['SEGMENT'] == segm]
# select only valid dark frames for this corrtag (taken within defined time frame with the same detector segment at the same detector voltage)
# if there are not enough matching dark frames it is very likely that the detector voltages do not match, so you have no choice as to accept a
# systematic error in the dark current (either over- or underestimate depending on the calibration windows and the gain sag state of the detector)
val_darks = []
sel_darks = darkframes[(darkframes['SEGMENT'] == segm) & \
(darkframes['VOLTAGE'] == voltage) & \
(darkframes['EXPSTART'] > exp_start - DARK_EXPSTART_INTERVAL) & \
(darkframes['EXPSTART'] < exp_start + DARK_EXPSTART_INTERVAL)]
for d in sel_darks['PATH']:
val_darks.append(fits.open(d))
if len(val_darks) >= MIN_DARKS:
print(str(len(val_darks)) + " valid dark frames.")
else:
print("Not enough dark frames for this exposure taken at voltage level {} !".format(voltage))
print("This likely means that the science voltage level has not been monitored by STScI.")
print("A mismatch in the science vs. dark frame voltage will lead to ~10% systematic errors in the estimated dark current (Makan et al. 2021, ApJ, 912, 38).")
a = input("Do you wish to include dark frames taken at all voltage levels? (y/n)")
if a == 'y' or a == 'Y':
val_darks = []
sel_darks = darkframes[(darkframes['SEGMENT'] == segm) & \
(darkframes['EXPSTART'] > exp_start - DARK_EXPSTART_INTERVAL) & \
(darkframes['EXPSTART'] < exp_start + DARK_EXPSTART_INTERVAL)]
for d in sel_darks['PATH']:
val_darks.append(fits.open(d))
if len(val_darks) >= MIN_DARKS:
print(str(len(val_darks)) + " valid dark frames.")
else:
print("Canceled by user!")
sys.exit()
# find the dispersion function for the Doppler shift xdopp(wavelength), this defines the calibration windows
wl = np.array(corr_data[(corr_data['WAVELENGTH'] > 1) & \
(corr_data['YFULL'] > ap_spec[0]) & \
(corr_data['YFULL'] < ap_spec[1]) & \
(corr_data['DQ'] == 0)]['WAVELENGTH'])
xdopp = np.array(corr_data[(corr_data['WAVELENGTH'] > 1) & \
(corr_data['YFULL'] > ap_spec[0]) & \
(corr_data['YFULL'] < ap_spec[1]) & \
(corr_data['DQ'] == 0)]['XDOPP'])
xdopp_disp = linregress(wl, xdopp)[0:2]
# find bad regions (with extended geocoronal emission) in xdopp coordinates
if opt_elem == 'G140L':
BAD_REGIONS = BAD_REGIONS_G140L
else:
BAD_REGIONS = BAD_REGIONS_G130M
bad_reg_darks = []
for reg in BAD_REGIONS:
bad_reg_darks.append(np.array(reg)*xdopp_disp[0] + xdopp_disp[1])
# remove bad regions from corrtag
for reg in BAD_REGIONS:
corr_data = corr_data[(corr_data['WAVELENGTH'] < reg[0]) | \
(corr_data['WAVELENGTH'] > reg[1])]
# remove counts outside of detector's active area
corr_data = corr_data[(corr_data['XDOPP'] < 15000) & \
(corr_data['XDOPP'] > 1500)]
# produce a cumulative pulse height distribution for the corrtag
# select background calibration windows in the corrtag data
corr_bkg1 = corr_data[(corr_data['YFULL'] < ap_bkg[1]) & \
(corr_data['YFULL'] > ap_bkg[0]) & \
((corr_data['DQ'] == 0) | \
(corr_data['DQ'] == 512) | \
(corr_data['DQ'] == 8192) | \
(corr_data['DQ'] == 8704))]
corr_bkg2 = corr_data[(corr_data['YFULL'] < ap_bkg[3]) & \
(corr_data['YFULL'] > ap_bkg[2]) & \
((corr_data['DQ'] == 0) | \
(corr_data['DQ'] == 512) | \
(corr_data['DQ'] == 8192) | \
(corr_data['DQ'] == 8704))]
corr_bkg = np.append(np.array(corr_bkg1['PHA']), \
np.array(corr_bkg2['PHA']))
# cumulative pulse height distribution
unique, counts = np.unique(corr_bkg, return_counts=True)
corr_pha_dist = np.zeros(shape=32, dtype=np.int32)
for i in range(len(unique)):
corr_pha_dist[unique[i]] = counts[i]
corr_data_cumsum = np.cumsum(corr_pha_dist)
corr_max_counts = corr_data_cumsum[-1]
corr_data_cumsum = corr_data_cumsum / corr_max_counts
# produce cumulative pulse height distribution for every dark frame
# and compare it to the pulse height distribution of the corrtag
# with a Kolmogorov–Smirnov test
KS_values = np.zeros(shape=(len(val_darks)), dtype=np.float32)
dark_max_counts = np.zeros(shape=(len(val_darks)), dtype=np.float32)
d = 0
for i in range(len(val_darks)):
dark_data = Table(val_darks[i][1].data)
# remove bad regions from the dark frame
for reg in bad_reg_darks:
dark_data = dark_data[(dark_data['XDOPP'] < reg[0]) | \
(dark_data['XDOPP'] > reg[1])]
# remove counts outside of detector's active area
dark_data = dark_data[(dark_data['XDOPP'] < 15000) & \
(dark_data['XDOPP'] > 1500)]
dark_bkg1 = dark_data[(dark_data['YFULL'] < ap_bkg[1]) & \
(dark_data['YFULL'] > ap_bkg[0]) & \
((dark_data['DQ'] == 0) | \
(dark_data['DQ'] == 8192))]
dark_bkg2 = dark_data[(dark_data['YFULL'] < ap_bkg[3]) & \
(dark_data['YFULL'] > ap_bkg[2]) & \
((dark_data['DQ'] == 0) | \
(dark_data['DQ'] == 8192))]
dark_bkg = np.append(np.array(dark_bkg1['PHA']), \
np.array(dark_bkg2['PHA']))
# cumulative pulse height distribution
unique, counts = np.unique(dark_bkg, return_counts=True)
dark_pha_dist = np.zeros(shape=32, dtype=np.int32)
for j in range(len(unique)):
dark_pha_dist[unique[j]] = counts[j]
dark_data_cumsum = np.cumsum(dark_pha_dist)
dark_max_counts[i] = dark_data_cumsum[-1]
dark_data_cumsum = dark_data_cumsum / dark_max_counts[i]
KS_values[i] = max(abs(dark_data_cumsum - corr_data_cumsum))
# sort K-S statistics
dark_sorted_KS = KS_values.argsort()
KS_test = KS_THRESHOLD
ind_threshold = 0
while(True):
for di in range(len(dark_sorted_KS)):
if KS_values[dark_sorted_KS[di]] < KS_test:
ind_threshold = di
if ind_threshold + 1 < MIN_DARKS:
KS_test = KS_test + KS_STEP
else:
break
number_of_darks = ind_threshold + 1
dark_total_counts = 0
# save prefixes of the used darks
dark_prefixes = val_darks[dark_sorted_KS[0]][0].\
header['FILENAME'].split("_")[0]
for i in range(1, number_of_darks):
dark_prefixes = dark_prefixes + ", " + \
val_darks[dark_sorted_KS[i]][0].\
header['FILENAME'].split("_")[0]
# calculate total counts of the used darks
for i in range(number_of_darks):
dark_total_counts = dark_total_counts + \
dark_max_counts[dark_sorted_KS[i]]
print("Number of used darks: " + str(number_of_darks))
print("KS threshold: " + str(KS_test))
# calculate scaling factor between the corrtag and combined darks
scaling_factor = corr_max_counts / dark_total_counts
scaling_factor_err = math.sqrt(scaling_factor * \
(1. + 1. / dark_total_counts) / dark_total_counts)
print("Average scaling factor: " + str(scaling_factor*number_of_darks))
# coadd best dark frames
darks_combined = Table(val_darks[dark_sorted_KS[0]][1].data)
for i in range(1, number_of_darks):
darks_combined = vstack([darks_combined, \
Table(val_darks[dark_sorted_KS[i]][1].data)])
# in the combined dark extract the 2d dark current spectrum within the science aperture using the science PHA limits
dark_psa = darks_combined[(darks_combined['YFULL'] >= ap_spec[0]) &\
(darks_combined['YFULL'] <= ap_spec[1]) & \
(darks_combined['PHA'] >= pha_min) & \
(darks_combined['PHA'] <= pha_max)]
dark_psa_xfull = np.array(dark_psa['XFULL'])
# shift xfull according to XSHIFT (darks are taken at nominal FPPOS whereas science data use multiple offsets)
xshift = corr_hdul[1].header['SHIFT1' + segm.split("FUV")[1]]
dark_psa_xshift = dark_psa_xfull - xshift
# collapse the 2d spectrum to 1d
dark_hist, tmp = np.histogram(dark_psa_xshift, 16384, range=(0, 16384))
dark_hist_mean = np.zeros(shape=16384, dtype=np.float32)
dark_hist_error = np.zeros(shape=16384, dtype=np.float32)
science_dq = np.array(data_x1d['DQ'])[0]
dark_hist_masked = np.ma.masked_where(science_dq != 0, dark_hist)
unmasked_hist_indices = np.ma.masked_where(science_dq != 0, \
np.arange(len(dark_hist_masked))).compressed()
half_width = int((BKG_AV - 1) / 2)
first_mean_value = unmasked_hist_indices[0] + half_width
# shift first value, if dq == 0
last_mean_value = unmasked_hist_indices[-1] - half_width
# Estimate smoothed dark current with running average between first and last valid points
for n in range(first_mean_value, last_mean_value + 1):
win = dark_hist_masked[n - half_width:n + half_width].compressed()
if len(win) > 0:
dark_hist_mean[n] = np.sum(win) / len(win)
dark_hist_error[n] = np.sqrt(np.sum(win)) / len(win)
else:
dark_hist_mean[n] = 0.0
dark_hist_error[n] = 0.0
# extrapolate first and last values to the left and right edges
for n in range(half_width):
dark_hist_mean[first_mean_value - half_width + n] = \
dark_hist_mean[first_mean_value]
dark_hist_error[first_mean_value - half_width + n] = \
dark_hist_error[first_mean_value]
dark_hist_mean[last_mean_value + n] = \
dark_hist_mean[last_mean_value]
dark_hist_error[last_mean_value + n] = \
dark_hist_error[last_mean_value]
# rescale to the corrtag file
dark_hist_mean = scaling_factor * dark_hist_mean
# error propagation due to the scaling factor
dark_hist_error = np.sqrt(np.power(dark_hist_mean, 2) * \
scaling_factor_err * scaling_factor_err + \
scaling_factor * scaling_factor * \
np.power(dark_hist_error, 2))
# Back out the flux calibration curve used by CALCOS
calib = np.divide(np.array(data_x1d['NET'])[0], \
np.array(data_x1d['FLUX'])[0], \
out=np.full(shape=(16384), fill_value=np.nan), \
where=np.array(data_x1d['FLUX'])[0] != 0)
wave = np.array(data_x1d['WAVELENGTH'])[0]
wave_inter = np.delete(wave, np.where(np.isnan(calib)))
calib_inter = np.delete(calib, np.where(np.isnan(calib)))
interp_calib = interp1d(wave_inter, calib_inter, \
kind='quadratic', \
fill_value=0, \
bounds_error=False)
calib = interp_calib(wave)
# Reset CALCOS data flags and weights: Good data (DQ=0) and low-response regions (DQ=1024) are included in coadd,
# all other regions (including grid wires with DQ=4) are excluded from coadd.
# Experience shows that grid wires cannot be flatfielded out in Poisson-limited data!
dq = np.array(data_x1d['DQ'])[0]
dq_wgt = np.zeros(shape=len(dq), dtype=np.int32)
dq[dq == 1024] = 2
for q in range(len(dq_wgt)):
if dq[q] == 0 or dq[q] == 2:
dq_wgt[q] = 1
else:
dq_wgt[q] = 0
# Back out 1D low-order flatfield calibration applied by CALCOS
gross = np.array(data_x1d['GROSS'][0])
net = np.array(data_x1d['NET'][0])
used_px = np.where((gross > 0) & (net > 0))[0]
valid_gross = gross[used_px]
valid_net = net[used_px]
valid_dq = dq[used_px]
valid_wave = wave[used_px]
flt_raw = valid_gross/valid_net
interp_flt = interp1d(valid_wave[np.where(valid_dq == 0)[0]], \