-
Notifications
You must be signed in to change notification settings - Fork 1
/
Copy pathspec.py
1295 lines (972 loc) · 49.6 KB
/
spec.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
from db import *
# Core packages
from math import fsum
import scipy.constants as cte
from scipy.special import wofz,erf
from scipy.optimize import curve_fit
from scipy.signal import convolve
from scipy.interpolate import interp1d
# Astro-packages
from astropy.time import Time
from astropy.stats import sigma_clip
# Plot packages
import matplotlib.pyplot as plt
from matplotlib.ticker import AutoMinorLocator
plt.rc('xtick', direction='in', top='on')
#plt.rc('xtick.minor', visible=True)
plt.rc('ytick', direction='in', right='on')
#plt.rc('ytick.minor', visible=True)
class spec():
def __init__(self, spectrum, snr=0, rv0=0, offset=0, cut_edges=False, orig='IACOB', delimiter=' '):
'''
Parameters
----------
spectrum : str
Enter the input spectrum, either name(s) of the star(s), the fits files
separated by coma, a .txt/.lst file containing the filenames, or '*'
if you want to select all the fits files inside the working folder.
snr : str, optional
If 'best' as input, it finds the best SNR spectrum for the given name.
If 'bestHF' same as 'best' but prioritizing spectra from HERMES/FEROS.
If specified, and there are more than one spectra in the DB, it will fail.
rv0 : float, optional
Enter the radial velocity correction to apply to the spectrum in km/s.
offset : float, optional
Enter the offset in wavelength [A] of the spectrum to plot. Default is 0.
cut_edges : boolean, optional
If True, it cuts the edges of the spectrum where the flux is asymptotic.
Default is False.
orig : str, optional
If 'IACOB', it assumes that the spectrum comes from the IACOB database.
If 'txt', it assumes that the spectrum comes from a two-columns file with
wavelength and flux with no header in the file.
If 'syn', it assumes that the spectrum comes from a synthetic spectrum.
Default is IACOB.
delimiter : str, optional
Delimiter used in the txt file to separate lambda/wavelength and flux.
This is valid only if orig='txt'. Default is ' '.
'''
if type(spectrum) == list:
if len(spectrum) > 1:
print('Problem in spec(): More than one spectrum selected.\nExiting...')
return None
else: spectrum = spectrum[0]
if orig in ['IACOB','iacob']:
self.orig = 'IACOB'
elif orig in ['txt','TXT','ascii','ASCII']:
self.orig = 'ascii'
elif orig in ['syn','SYN','synthetic','SYNTHETIC']:
self.orig = 'synthetic'
if self.orig == 'IACOB':
self.fullpath = findstar(spectrum, snr=snr)
elif self.orig == 'ascii' or self.orig == 'synthetic':
self.fullpath = search(spectrum, datadir+'ASCII'+os.sep)
# So far the search function only returns the first match
if self.fullpath is None:
print('Problem in spec(): No spectrum found.\nExiting...')
return None
if self.orig == 'IACOB':
if len(self.fullpath) > 1:
print('Problem in spec(): More than one spectrum selected.\nExiting...')
return None
else:
self.fullpath = self.fullpath[0]
self.filename = self.fullpath.split(os.sep)[-1]
self.id_star = self.fullpath.split(os.sep)[-1].split('_')[0]
if self.orig == 'IACOB' or \
(self.orig in ['ascii','synthetic'] and any([i in self.filename for i in ['_V','_R']])):
self.resolution = int(re.split('(\d*\d+)',self.filename)[-2])
else:
self.resolution = 5000
print('Warning: Resolution not found in filename. Assuming R5000.')
self.offset = offset # Note, run self.waveflux to apply offset.
self.rv0 = rv0 # Note, run self.waveflux to apply the correction.
self.SpC = '' # Spectral classification
self.waveflux(cut_edges=cut_edges, delimiter=delimiter)
def get_spc(self):
'''
Function to retrieve the spectral classification of the star from Simbad
and add it to the class.
'''
try:
query = Simbad.query_object(self.id_star)
if query is None and 'HD' in self.id_star:
new_id_star = self.id_star.replace('HD', 'HD ')
query = Simbad.query_object(new_id_star)
self.SpC = query['SP_TYPE'][0]
#self.otypes = query['OTYPES'][0]
except:
print('Spectral classification could not be queried for: %s' % self.id_star)
self.SpC = ''
def waveflux(self, lwl=None, rwl=None, width=0, helcorr='hel', cut_edges=False, delimiter=' '):
'''
Function to load or update the wavelength and flux vectors and optionally apply
an offset or a radial velocity correction if they are different from 0 in the
class. It also adds the HJD, dlam to the class.
Parameters
----------
lwl : float, optional
Left wavelength limit of the spectra. Default is None (no change).
rwl : float, optional
Right wavelength limit of the spectra. Default is None (no change).
width : int, optional
Sets the width in [AA] where the line fits well in. Default is 10.
helcorr : str, optional
If 'hel' as input (default), it applies the heliocentric correction.
cut_edges : boolean, optional
If True, it cuts the edges of the spectrum where the flux is asymptotic.
Default is False.
delimiter : str, optional
Delimiter used in the txt file to separate lambda/wavelength and flux.
Default is ' ' (see ascii.read or db.findtable)
Note: All ascii files should be placed inside datadir/ASCII subfolder.
Returns
-------
In addition to update the class with new data, it returns the wavelength and flux
vectors together with the HJD.
'''
if self.orig == 'IACOB':
# Retrieve the key values fron the fits header
hdu = fits.open(self.fullpath) # Open the fits image file
hdu.verify('fix') # Fix possible issues with the keywords
header0 = hdu[0].header # Read header of primary extension
instrum = header0['INSTRUME'] # Instrument
lam0 = header0['CRVAL1'] # Get the wavelength of the first pixel
dlam = header0['CDELT1'] # Step of increase in wavelength
pix0 = header0['CRPIX1'] # Reference pixel (generally 1, FEROS -49)
spec_length = header0['NAXIS1'] # Length of the spectrum
# Alternatively use len(hdu[0].data[0]) (NOT/MERCATOR) or len(hdu[0].data)
# Correct Mercator CRVAL1 20101018-19:
if any(bad in self.fullpath for bad in ['_20101018_','_20101019_']) and lam0 == 3763.9375:
lam0 = 3763.61
# Assign barycentric rv correction at midpoint [km/s]
if 'I-VBAR' in header0:
self.vbar = header0['I-VBAR']
#self.vbar = header0['BVCOR'] # Specific keyword from MERCATOR
#self.vbar = header0['VHELIO'] # Specific keyword from NOT
else:
print('No helio/bary-centric correction applied to' + self.fullpath)
self.vbar = 0
# Assign Heliocentric Julian Date at midpoint
if 'I-HJD' in header0:
self.hjd = header0['I-HJD']
else:
self.hjd = Time(header0['DATE'], scale='utc').jd
# Assign SNR
if 'I-SNR' in header0:
self.snr = header0['I-SNR']
else:
self.snr = np.nan
# Assign SpC
if 'I-SPC' in header0:
self.SpC = header0['I-SPC']
# Make lists with wavelength and flux for each spectrum
# Those with log or from FEROS are already corrected from heliocentric velocity
wave = lam0 + dlam*(np.arange(spec_length) - pix0 + 1)
if '_log' in self.fullpath:
wave = np.exp(wave)
elif helcorr == 'hel' and not instrum == 'FEROS':
wave = wave*(1 + 1000*self.vbar/cte.c)
# If width is bigger than the spectrum, set it to the spectrum length
if width >= (wave[-1] - wave[0]):
width = wave[-1] - wave[0]
# Get the flux of the spectrum
try:
flux = hdu[0].data[0]
except:
flux = hdu[0].data
hdu.close()
elif self.orig == 'ascii' or self.orig == 'synthetic':
data = findtable(self.filename, path=self.fullpath.replace(self.filename,''), delimiter=delimiter, format='basic')
wl_keywords = ['wave','wavelength','lambda','lamb','ang','angstroms']
idx_wl = next((i for i, item in enumerate(data.colnames) if item.lower() in wl_keywords), None)
if idx_wl == None:
print('Problem in spec(): No wavelength column found in the first column of the ascii file. ')
return None
else:
wave = np.asarray(data[data.colnames[idx_wl]])
fx_keywords = ['flux','fluxes','norm_flux','flux_norm']
idx_fx = next((i for i, item in enumerate(data.colnames) if item.lower() in fx_keywords), None)
if idx_fx == None:
print('Problem in spec(): No flux column found in the second column of the ascii file.')
return None
else:
flux = data[data.colnames[idx_fx]]
dlam = (wave[-1]-wave[0])/(len(wave)-1)
self.vbar = 0
self.hjd = 0
if self.orig != 'synthetic':
self.get_spc()
# Apply the radial velocity correction (if provided)
wave = wave*(1 - 1000*self.rv0/cte.c)
# Apply the offset correction (if provided)
wave = wave - self.offset
# Cut the edges of the spectrum where the flux is asymptotic
if cut_edges == True:
if flux[1] < flux[0]:
l_cut = next(pix for pix in range(len(wave)) if flux[pix+1] > flux[pix])
else:
l_cut = next(pix for pix in range(len(wave)) if flux[pix+1] < flux[pix])
if flux[-2] < flux[-1]:
r_cut = next(pix for pix in reversed(range(len(wave))) if flux[pix-1] > flux[pix])
else:
r_cut = next(pix for pix in reversed(range(len(wave))) if flux[pix-1] < flux[pix])
wave = wave[l_cut:r_cut]
flux = flux[l_cut:r_cut]
# Cut the spectrum to the desired wavelength range given by lwl and rwl
if lwl != None and rwl != None:
if wave[0] > lwl+dlam or wave[-1] < rwl-dlam:
print('WARNING: Wavelenght limits outside spectrum wavelength range.')
flux = flux[(wave >= lwl-width/2) & (wave <= rwl+width/2)]
wave = wave[(wave >= lwl-width/2) & (wave <= rwl+width/2)]
if '_log' in self.fullpath and self.orig in ['IACOB','ascii']:
self.dlam = (wave[-1]-wave[0])/(len(wave)-1)
else:
self.dlam = dlam
self.wave_0 = wave
self.flux_0 = flux
self.wave = wave
self.flux = flux
hjd = self.hjd
return wave, flux, hjd
def fitline(self, line, width=15, tol=150., func='g', iter=3, fw3414=False, info=False,
outfit=False, plot=False):
'''
Function to fit a spectral line to a function. Different functions account for
different lines depending on their natural profile (e.g. metallic lines should be
fitted with either a gaussian, Voigt, rotational profiles). See fitline_readme.txt
and fitline_diagram.txt included with pyIACOB for more details.
Parameters
----------
line : float
Sets the central wavelength of the line to search and fit.
width : int, optional
Sets the width in [AA] where the line fits well in. Default is 15.
tol : int, optional
Sets the tolerance [km/s] to shifting the spectrum in order to fit the line.
func : str, optional
Choose the function to fit the line:
'g' Gaussian (default); 'l' Lorentzian; 'v' Voigt; 'r' Rotational.
'vr_H' Voigt profile with rotation for H lines.
'vr_Z' Voigt profile with rotation for metallic lines.
'vrg_H' Voigt profile with rotation and extra gaussian wings for H lines.
'vrg_Z' Voigt profile with rotation and extra gaussian wings for metallic lines.
iter : int, optional
Number of iterations to optimize window width. Default is 3.
fw3414 : boolean, optional
If 'True', it will output the FW difference at 3/4 - 1/4 of the line depth.
Default is 'False'.
info : boolean, optional
If 'True', it will print information for each line fitting.
Default is 'False'.
outfit : boolean, optional
If 'True', it will output the fitting as an array.
Default is 'False'.
plot : boolean, optional
If 'True', it will create and show the plots.
Returns
-------
Dictionary containing the results, fitting parameters and optionally the
original line, normalized flux and fitted flux.
Notes
-----
Emission lines are excluded, see "Filtering emission lines" section.
Returns: Parameters from the fitted line and a last value containing data
with input wavelength and flux, plus the flux normalized, the flux of the
fitted line, and the fitting parameters.
'''
fitsol = {'sol': 0, 'line': np.nan, 'RV_A': np.nan, 'RV_kms': np.nan, 'EW': np.nan,
'FWHM': np.nan, 'depth': np.nan, 'q_fit': np.nan, 'snr': np.nan}
#============================== Parameters =============================
# Catch line input containing more than one line
if type(line) == str and (',' in line or ' ' in line.strip()):
print('Problem in spec(): More than one line selected.\nExiting...')
return fitsol
line = float(line)
dlamb = line/self.resolution
# Maximum shift between the minimum of the fitted line and the tabulated value
tol_aa = float(tol)*(line)*1000/cte.c # Changes km/s to angstroms
# Maximum FWHM allowed (could be up to 18-20 for H lines in very cold stars)
FWHM_max = 17
#======== Determine the SNR based on the wavelength of the line ========
if self.orig != 'synthetic':
if line > 3000 and line < 4000:
snr_spec = self.snrcalc('uv')
elif line > 4000 and line < 5000:
snr_spec = self.snrcalc('b')
elif line > 5000 and line < 6000:
snr_spec = self.snrcalc('v')
elif line > 6000 and line < 7000:
snr_spec = self.snrcalc('r')
else:
snr_spec = np.nan
#=========== Dictionary and boundary limits for the functions ==========
fit_dic = {'g': ('A','lam0','sigma'),
'l': ('A','lam0','gamma','y'),
'v': ('A','lam0','sigma','gamma','y'),
'r': ('A','lam0','sigma','vsini'),
'vr': ('A','lam0','sigma','gamma','vsini','y'),
'vrg': ('A','lam0','sigma','gamma','vsini','A2','sigma2','y')}
# Fitting function: Gaussian | A,lam0,sig
if func == 'g':
fitfunc = f_gaussian1
bounds = ([-1,line-tol_aa,0],
[ 0,line+tol_aa,3]) # 3/6 for narrow/broad lines (default 6)
# Fitting function: Lorentzian | A,lam0,gamma,y
elif func == 'l':
fitfunc = f_lorentzian
bounds = ([-1,line-tol_aa, 0,1. ],
[ 0,line+tol_aa,10,1.01])
# Fitting function: Voigt profile | A,lam0,sigma,gamma,y
elif func == 'v':
fitfunc = f_voigt
bounds = ([-10,line-tol_aa,0. ,0. ,1. ], #'A' ~15 for As
[ 0,line+tol_aa,7.5,8.5,1.01])
# Fitting function: Rotational profile | A,lam0,sigma,vsini
elif func == 'r':
fitfunc = f_rot
bounds = ([.0,line-tol_aa,0. , 1],
[.3,line+tol_aa,2.5,410])
# Fitting function: Voigt x Rotational profile | A,lam0,sigma,gamma,vsini,y
elif func == 'vr_H':
fitfunc = f_voigtrot
bounds = ([-.3,line-tol_aa,0. ,0, 1,.0 ],
[ .0,line+tol_aa,4.5,8,410,.01])
elif func == 'vr_Z':
fitfunc = f_voigtrot
bounds = ([-.3,line-tol_aa,0. ,0, 1,.0 ],
[ .0,line+tol_aa,1.5,1,410,.01])
# Fitting function: Voigt x Rotational + Gaussian profile
# A,lam0,sigma,gamma,vsini,A2,sigma2,y
elif func == 'vrg_H':
fitfunc = f_vrg
bounds = ([-.5,line-tol_aa, 0, 0, 1,-.07,0,.0 ],
[ .0,line+tol_aa,10,10,410, .0 ,4,.01])
elif func == 'vrg_Z':
fitfunc = f_vrg
bounds = ([-.1,line-tol_aa,0. ,0, 1,-1.3,0,-.01], # y=-.01 = larger EWs
[ .0,line+tol_aa,1.5,1,410, 0. ,2, .01])
#============================= Line fitting ============================
i = 0; width_i = width
while i < iter:
# Extracting the window of the spectrum
window = (self.wave >= line-width_i/2) & (self.wave <= line+width_i/2)
if not any(window):
print('Line %sA not in spectra.\n' % line)
return fitsol
flux = self.flux[window]
wave = self.wave[window]
dlam_mean = (wave[-1]-wave[0])/(len(wave)-1)
# Auto-resampling
#if dlam_mean >= 0.025 and not 'log' in star.filename:
# factor = dlam_mean/0.025; star.resamp(factor)
# Find regions of the continuum to use during normalization (4 iter)
for j in range(4):
if j == 0:
mask_i = ~np.isnan(flux)
else:
mask_i = ~sigma_clip(flux/continuum_i, maxiters=None,
sigma_lower=1.4, sigma_upper=2.5, axis=-1).mask #1.4 before
c_fit = np.poly1d(np.polyfit(wave[mask_i], flux[mask_i], 1))
continuum_i = c_fit(wave)
# Final normalization of the iteration
flux_norm_i = flux / continuum_i
#========================= Fitting the line ========================
try:
popt_i = curve_fit(fitfunc, wave, flux_norm_i, bounds=bounds)[0]
flux_fit_i = fitfunc(wave, *popt_i)
# Calculate the empirical approximate FWHM
medval = (max(flux_fit_i) + min(flux_fit_i))/2
medpos = [np.where(flux_fit_i <= medval)[0][value] for value in (0,-1)]
FWHM = round(wave[medpos[1]] - wave[medpos[0]], 2)
# Checking step results
# The min FWHM will be defined by either 3 times the dlam, or 3/4 of the
# minimum theoretical FWHM given the input line and resolution
FWHM_min = np.max([3*dlam_mean, 3/4*dlamb])
if FWHM_min < FWHM < FWHM_max:
continuum = continuum_i
flux_norm = flux_norm_i
flux_fit = flux_fit_i
mask = mask_i
popt = popt_i
width = width_i
width_i = FWHM*7
i = i + 1
elif FWHM < FWHM_min:
print('WARNING: FWHM(%.1f) < minimum FWHM for %.3fA' % (FWHM,line))
break
elif FWHM > FWHM_max:
print('WARNING: FWHM(%.1f) > maximum FWHM for %.3fA ' % (FWHM,line))
break
except: break
#======================== Checking final results =======================
window = (self.wave >= line-width/2.) & (self.wave <= line+width/2.)
flux = self.flux[window]
wave = self.wave[window]
# If the line was not fitted
if i == 0:
if outfit == True:
fitsol['wave'] = wave
fitsol['flux_norm'] = flux_norm_i
fitsol['flux_fit'] = [np.nan]*len(wave)
if info is True:
print('Problem in spectrum %s' % self.filename)
print('Line %sA could not be fitted or does not exist.\n' % line)
return fitsol
# If the line was fitted
line_f = wave[flux_fit == min(flux_fit)][0]
# ...but the line is found outside the tolerance
if abs(line - line_f) > tol_aa:
if outfit == True:
fitsol['wave'] = wave
fitsol['flux_norm'] = flux_norm
fitsol['flux_fit'] = [np.nan]*len(wave)
if info is True:
print('Line %sA found outside tolerance.\n' % line)
return fitsol
RV_A = round((line_f - line), 3)
RV_kms = round(((line_f - line)/line)*cte.c/1000, 2) # max precision is 100 m/s
line_f = round(line_f, 3)
if info is True:
print('Line %sA found at %.3fA -> RV: %.1fkm/s\n' % (line,line_f,RV_kms))
#=========================== Calculate the EW ==========================
# stackoverflow.com/questions/34075111/calculate-equivalent-width-using-python-code
# When emission is considered abs should be removed and 1-flux_fit -> flux_fit-1
EW = 0.5*abs(fsum((wave[wl-1] - wave[wl])*((1 - flux_fit[wl-1]) +
(1 - flux_fit[wl])) for wl in range(1, len(flux_fit))))
EW = round(1000*EW)
#====================== Calculate the final FWHM =======================
medval = (max(flux_fit) + min(flux_fit))/2
medpos = [np.where(flux_fit <= medval)[0][value] for value in (0,-1)]
try:
l_val = np.interp(medval, [flux_fit[medpos[0]], flux_fit[medpos[0]-1]],
[wave[medpos[0]], wave[medpos[0]-1]])
except:
l_val = wave[medpos[0]]
try:
r_val = np.interp(medval, [flux_fit[medpos[1]], flux_fit[medpos[1]+1]],
[wave[medpos[1]], wave[medpos[1]+1]])
except:
r_val = wave[medpos[1]]
FWHM = round(r_val - l_val, 3)
#======= Calculate width difference at 3/4-1/4 of the line depth =======
if fw3414 is True:
try:
lowval = (max(flux_fit) + 3*min(flux_fit))/4
uppval = (3*max(flux_fit) + min(flux_fit))/4
medpos = []; FW34_14 = []
for par,val in zip(['FW14_Hb','FW34_Hb'],[lowval,uppval]):
medpos_i = [np.where(flux_fit <= val)[0][value] for value in (0,-1)]
medpos.append(medpos_i)
try: l_val = np.interp(val,[flux_fit[medpos_i[0]],flux_fit[medpos_i[0]-1]],
[wave[medpos_i[0]],wave[medpos_i[0]-1]])
except: l_val = wave[medpos_i[0]]
try: r_val = np.interp(val,[flux_fit[medpos_i[1]],flux_fit[medpos_i[1]+1]],
[wave[medpos_i[1]],wave[medpos_i[1]+1]])
except: r_val = wave[medpos_i[1]]
FW34_14.append(round(r_val-l_val, 3))
FW34_14 = FW34_14[1]-FW34_14[0]
except:
print('Problem calculating the FW at 1/4 and 3/4 of the Hb line.')
FW34_14 = np.nan
else:
FW34_14 = np.nan
#======================= Calculate the line depth ======================
depth = round(1 - min(flux_fit), 3)
#===================== Calculate the SNR continuum =====================
if self.orig == 'synthetic':
snr = np.nan
else:
sigma_cont = np.std(flux_norm[mask])
if sigma_cont <= 0.002:
print('SNR continuum for %.3fA is very high. Set to 1000' % line)
snr = 1000
else:
snr = int(1/sigma_cont)
# If the SNR measured on the continuum of the line is too different from
# the SNR measured on the wider region, the latter is used.
if snr == None or abs(snr_spec-snr)/snr_spec > 0.30:
snr = snr_spec
if snr_spec == None and snr == None:
snr = np.nan
#============================= Quality value ===========================
q_fit = 1/np.std(flux_norm[flux_fit<(1-0.2*depth)]/flux_fit[flux_fit<(1-0.2*depth)]) #simple
q_fit = round(q_fit, 3)
#================================ Plot =================================
if plot is True:
fig, ax = plt.subplots()
if fw3414 is True and FW34_14 != np.nan:
ax.plot([wave[medpos[0][0]],wave[medpos[0][1]]],[lowval,lowval],'gray',linestyle='--',lw=.5)
ax.plot([wave[medpos[1][0]],wave[medpos[1][1]]],[uppval,uppval],'gray',linestyle='--',lw=.5)
ax.plot(wave, flux, c='orange', lw=.5)
ax.plot(wave, continuum, c='r', lw=.5)
ax.plot(wave, flux_norm, c='b', lw=.5)
ax.plot(wave, flux_fit, c='g', lw=.5)
ax.plot(wave, np.where(mask==False, 1, np.nan) + 0.01, 'k', lw=.5)
ax.set_title('%s | %.2f | RV: %d | EW: %d | FWHM: %.2f | FW34-14: %.2f | SNR: %.f' %
(self.id_star,line_f,RV_kms,EW,FWHM,FW34_14,snr), fontsize=8)
ax.set_yticks([])
ax.set_xlabel('$\lambda$ $[\AA]$', size=13)
ax.set_ylabel('Normalized flux', size=13)
ax.tick_params(direction='in', top='on')
ax.figure.subplots_adjust(top=.9, bottom=.12, right=.88, left=.08)
plt.show(block=False)
#=======================================================================
#================= Packing the results in a dictionary =================
fitsol = {'sol':1, 'line':line_f, 'RV_A':RV_A, 'RV_kms':RV_kms,
'EW':EW, 'FWHM':FWHM, 'FW34_14':FW34_14, 'depth':depth,
'q_fit':q_fit, 'snr':snr}
for f_par,par in zip(fit_dic[func.split('_')[0]], popt):
fitsol[f_par] = round(par, 3)
if outfit == True:
fitsol['wave'] = wave
fitsol['flux_norm'] = flux_norm
fitsol['flux_fit'] = flux_fit
return fitsol
# Theorerical FWHM:
#if func == 'g': jFWHM = 2*np.sqrt(2*np.log(2))*popt[2]
#elif func == 'l': jFWHM = 2*abs(popt[2])
#elif func == 'v': jFWHM = 2*(.5346*popt[3]+np.sqrt(.2166*(popt[3]**2)+popt[2]**2))
#elif func == 'r': jFWHM = 1.7*popt[3]*line*1000/cte.c
#jFWHM = round(jFWHM, 3)
def snrcalc(self, zone='B'):
'''
Function to calculate the Signal to Noise Ratio in different regions of the
spectra if available.
Parameters
----------
zone : str, optional
Select the zone to calculate the spectra.
'b'/'B' -> 4000-5000 A
'v'/'V' -> 5000-6000 A
'r'/'R' -> 6000-7000 A
'all'/'ALL' -> 4000-7000 A
Returns
-------
Measured signal-to-noise ratio value.
'''
if self.orig == 'synthetic':
print('No SNR available for synthetic spectra.')
return np.nan
if zone in ['uv','UV']:
mask = (self.wave > 3000) & (self.wave < 4000)
elif zone in ['b','B']:
mask = (self.wave > 4000) & (self.wave < 5000)
elif zone in ['v','V']:
mask = (self.wave > 5000) & (self.wave < 6000)
elif zone in ['r','R']:
mask = (self.wave > 6000) & (self.wave < 7000)
elif zone in ['all','ALL']:
mask = (self.wave > 3000) & (self.wave < 7000)
if len(self.wave[mask]) == 0:
return np.nan
lambda0 = np.mean(self.wave[mask])
resol = 10000
sigma = lambda0/(2.35482*float(resol))
gauss = f_gaussian(np.arange(-5*sigma, 5*sigma, self.dlam), sigma)
kernel = gauss/np.trapz(gauss)
convoluted = 1 + convolve(self.flux[mask] - 1, kernel, mode='same')
flux_norm = self.flux[mask]/convoluted
snr_all = []
for gap in findlist('snr_gaps.txt'):
lwl,rwl = [float(i) for i in gap.split('-')]
flux_norm_i = flux_norm[(self.wave[mask] >= lwl) & (self.wave[mask] <= rwl)]
if len(flux_norm_i) == 0:
continue
sig_clip = 3
std = np.std(flux_norm_i)
flux_clean = np.where(abs(flux_norm_i - 1) > sig_clip*std, np.nan, flux_norm_i)
snr_all.append(1/np.nanstd(flux_clean))
self.snr = np.nanmean(snr_all)
if not np.isnan(np.nanmean(snr_all)):
return int(round(np.nanmean(snr_all)))
def cosmic(self, method='zscore',
dmin=0.05, zs_cut=5,
wl_split=100, ker_sig=2, ker_iter=3, sig_g=None,
protect_em_lines=True):
'''
Function to remove cosmic rays in the spectra by different approaches.
Parameters
----------
method : str, optional
Method for the cosmic ray removal strategy. Only zscore (def) or kernel.
dmin : int/float, optional
Minium distance between flux and cleaned flux to consider for replacement.
Default is 0.05.
zs_cut : int/float, optional
Threshold value used in the zscore method for finding rays. Default is 4.
Tip: For noisy spectra, rise this value up to 7-9.
wl_split : int/float, optional
In the 'kernel' method, wavelength size used to split the spectrum before
applying the cosmic removal. Default is 100.
ker_sig : float, optional
Sigma clipping value used to remove rays. Default is 2.
ker_iter : int, optional
Number of iterations of the sigma clipping to remove cosmic rays in the kernel method.
Default is 3.
sig_g : float, optional
Sigma of the gaussian function used to construct the kernel.
Default is the theoretical sigma based on wavelength and resolution.
protect_em_lines : boolean, optional
If True, some emission lines will be masked from cosmic rays removal.
Default is True.
Returns
-------
Nothing, but the flux is replaced and cleaned from rays.
'''
if self.orig == 'synthetic':
print('Cosmic rays removal not available for synthetic spectra.')
return None
if method == 'zscore':
# www.towardsdatascience.com/removing-spikes-from-raman-spectra-8a9fdda0ac22
# First we calculated (nabla)x(i):
delta_flux = np.diff(self.flux)
median_int = np.median(delta_flux)
mad_int = np.median([np.abs(delta_flux - median_int)])
modified_z_scores = 0.6745 * (delta_flux - median_int) / mad_int
# The multiplier 0.6745 is the 0.75th quartile of the standard normal
# distribution, to which the median absolute deviation converges to.
modified_z_scores = np.concatenate(([0], np.abs(modified_z_scores)))
flux_clean = np.where(modified_z_scores > zs_cut, np.nan, self.flux)
elif method == 'kernel':
wl_split = 100 # Range in angstroms in which the spectrum will be initially splitted
wl_range = self.wave[-1]-self.wave[0]
if wl_range < wl_split:
wl_split = wl_range
if 1 < wl_range/wl_split <= 2:
wl_split = wl_range/2 + 1
elif wl_range/wl_split > 2:
wl_split += (wl_range % wl_split) / int(wl_range/wl_split)
n = int(wl_range/wl_split)
flux_clean = []
for i in range(n):
mask = (self.wave >= self.wave[0]+i*wl_split) & (self.wave < self.wave[0]+(i+1)*wl_split)
if i == range(n)[-1]:
mask[-1] = True # To catch the last flux value from not using <= above
resolution = float(self.resolution) # resolution = 5000 # Empirical optimal value
dlam = self.dlam # dlam = 0.2564975 # Empirical optimal value
if sig_g is None:
lambda0 = np.mean(self.wave[mask])
sig_g = lambda0/(2.35482*resolution)
else:
sig_g = float(sig_g)
x = np.arange(-5*sig_g, 5*sig_g+dlam, dlam)
gauss = f_gaussian(x,sig_g)
kernel = gauss/np.trapz(gauss)
convoluted = 1 + convolve(self.flux[mask] - 1, kernel, mode='same')
flux_norm = self.flux[mask]/convoluted
for i in range(ker_iter):
std = np.nanstd(flux_norm)
flux_norm = np.where(abs(flux_norm - 1) > ker_sig*std, np.nan, flux_norm)
flux_clean = np.concatenate([flux_clean,np.where(np.isnan(flux_norm), np.nan, self.flux[mask])])
nans = np.isnan(flux_clean)
x = lambda z: z.nonzero()[0]
flux_clean[nans] = np.interp(x(nans), x(~nans), flux_clean[~nans])
# Recover the original flux in the regions of the spectrum where telluric lines are
flux_clean = np.where(
((self.wave > 3932.5) & (self.wave < 3934.5)) # Not sure if this is exactly telluric line
| ((self.wave > 5885.0) & (self.wave < 5900.0))
| ((self.wave > 6865.0) & (self.wave < 7035.0))
| ((self.wave > 7160.0) & (self.wave < 7340.0))
| ((self.wave > 7585.0) & (self.wave < 7700.0))
| ((self.wave > 8116.0) & (self.wave < 8380.0))
| ((self.wave > 8915.0) & (self.wave < 9220.0))
,self.flux, flux_clean)
# Recover the original flux in the regions of the spectrum where emission lines are
if protect_em_lines == True:
for wl_em in [3967.79,4958.911,5006.843,6300.304,6548.04,6583.46,6716.44,6730.82]:
mask = (self.wave > wl_em-0.8) & (self.wave < wl_em+0.8)
flux_clean[mask] = self.flux[mask]
# Recover the original flux if the difference is lower than dmin value
flux_clean = np.where(
(self.flux > flux_clean) & (self.flux - flux_clean > dmin), flux_clean, self.flux)
self.flux = flux_clean
return None
def cosmetic():
'''
IN DEVELOPMENT - Function to remove cosmetic defects recurrently found in the
spectra.
'''
# Correct from FEROS issues for window in
#[[4506,4507.5],[4693.7,4696],[4795.,4797.],[4900.7,4902.3],[6636.3,6638.3]]
return None
def degrade(self, resol, profile='g', vsini=None, vmac=None):
'''
Function to degrade a spectrum to a certain resolution by convolving it to a
gaussian (pure degradation) or to account for rotational+macroturbulence effect
for example if a synthetic spectrum is loaded.
Parameters
----------
resol : int/float, optional
Resolution of the gaussian profile used to degrade the spectrum.
profile : str
Use 'g' for gaussian profile convolution (Default).
Use 'rotmac' for rotational+macroturbulence profile convolution.
vsini : int/float, optiomal
Value of vsini. Only valid for rotational+macroturbulence profile.
vmac : int/float, optiomal
Value of vmac. Only valid for rotational+macroturbulence profile.
Returns
-------
Nothing, but the flux is replaced by the degraded one.
'''
lambda0 = np.nanmean(self.wave)
if profile == 'g' and (vsini==None and vmac==None):
sigma = lambda0/(2.35482*float(resol))
x = np.arange(-10*sigma, 10*sigma+self.dlam, self.dlam)
gauss = f_gaussian(x, sigma)
kernel = gauss/np.trapz(gauss)
self.resolution = resol
elif profile == 'rotmac' and vsini!=None and vmac!=None:
x = np.arange(-9, 9+self.dlam, self.dlam)
rotmac = f_rotmac(x, lambda0, vsini, vmac)
kernel = rotmac/np.trapz(rotmac)
# Remove the nans from the flux
mask = np.where(np.isnan(self.flux) == False)[0]
# Convolve the flux with the kernel
convoluted = 1 + convolve(self.flux[mask] - 1, kernel, mode='same')
# Replace the flux by the convoluted one recovering the original flux with the nans
self.flux[mask] = convoluted
def resamp(self, dlam, lwl=None, rwl=None, method='linear'):
'''
Function to resample a spectrum into a fixed delta-lambda and wavelength range.
Parameters
----------
dlam : float/int
New delta lambda to be used for the output spectra.
lwl : float/int, optional
Enter the forced initial wavelength to be used during interpolation.
If None, the original initial wavelength will be used.
rwl : float/int, optional
Enter the forced final wavelength to be used during interpolation.
If None, the original final wavelength will be used.
method : str, optional
Enter the interpolation method to be used. See doc for np.interp1d.
Default is 'linear'. See interp1d (keyword kind) for more options.
Returns
-------
Nothing, but the spectrum (wavelength,flux) is resampled.
'''
# Check that the input dlam is either a float or an integrer
if not isinstance(dlam, float) and not isinstance(dlam, int):
print('ERROR: The input delta lambda is not a float or an integrer.')
return None
self.dlam = dlam
if dlam > np.mean(self.wave)/self.resolution/3:
# It is divided by 3 to at least have 3 pixels in a gaussian
print('WARNING: The new delta lambda implies lossing information...')
if lwl is None or lwl < self.wave[0]:
lwl = self.wave[0]
if rwl is None or rwl > self.wave[-1]:
rwl = self.wave[-1]
# Interpolate the spectrum to the new delta lambda / step size