Coverage for python/lsst/sims/maf/utils/snNSNUtils.py : 10%

Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
1from lsst.sims.photUtils import SignalToNoise
2from lsst.sims.photUtils import PhotometricParameters
3from lsst.sims.photUtils import Bandpass, Sed
5import numpy as np
6import matplotlib.pyplot as plt
7import math
8from scipy.constants import *
9from functools import wraps
10import os
11import h5py
12import multiprocessing
13from astropy.table import Table
14import pandas as pd
15from scipy import interpolate, integrate
16from scipy.interpolate import RegularGridInterpolator
17from astropy import (units as u, constants as const)
18from astropy.cosmology import FlatLambdaCDM
20STERADIAN2SQDEG = 180.**2 / np.pi**2
21# Mpc^3 -> Mpc^3/sr
22norm = 1. / (4. * np.pi)
25class LCfast:
26 """
27 class to simulate supernovae light curves in a fast way
28 The method relies on templates and broadcasting to increase speed
29 Parameters
30 ---------------
31 reference_lc:
32 x1: float
33 SN stretch
34 color: float
35 SN color
36 telescope: Telescope()
37 telescope for the study
38 mjdCol: str, opt
39 name of the MJD col in data to simulate (default: observationStartMJD)
40 RACol: str, opt
41 name of the RA col in data to simulate (default: fieldRA)
42 DecCol: str, opt
43 name of the Dec col in data to simulate (default: fieldDec)
44 filterCol: str, opt
45 name of the filter col in data to simulate (default: filter)
46 exptimeCol: str, opt
47 name of the exposure time col in data to simulate (default: visitExposureTime)
48 m5Col: str, opt
49 name of the fiveSigmaDepth col in data to simulate (default: fiveSigmaDepth)
50 seasonCol: str, opt
51 name of the season col in data to simulate (default: season)
52 snr_min: float, opt
53 minimal Signal-to-Noise Ratio to apply on LC points (default: 5)
54 """
56 def __init__(self, reference_lc, x1, color,
57 telescope, mjdCol='observationStartMJD',
58 RACol='fieldRA', DecCol='fieldDec',
59 filterCol='filter', exptimeCol='visitExposureTime',
60 m5Col='fiveSigmaDepth', seasonCol='season',
61 nexpCol='numExposures',
62 snr_min=5.):
64 # grab all vals
65 self.RACol = RACol
66 self.DecCol = DecCol
67 self.filterCol = filterCol
68 self.mjdCol = mjdCol
69 self.m5Col = m5Col
70 self.exptimeCol = exptimeCol
71 self.seasonCol = seasonCol
72 self.nexpCol = nexpCol
73 self.x1 = x1
74 self.color = color
76 # Loading reference file
77 self.reference_lc = reference_lc
79 self.telescope = telescope
81 # This cutoffs are used to select observations:
82 # phase = (mjd - DayMax)/(1.+z)
83 # selection: min_rf_phase < phase < max_rf_phase
84 # and blue_cutoff < mean_rest_frame < red_cutoff
85 # where mean_rest_frame = telescope.mean_wavelength/(1.+z)
86 self.blue_cutoff = 380.
87 self.red_cutoff = 800.
89 # SN parameters for Fisher matrix estimation
90 self.param_Fisher = ['x0', 'x1', 'daymax', 'color']
92 self.snr_min = snr_min
94 # getting the telescope zp
95 self.zp = {}
96 for b in 'ugrizy':
97 self.zp[b] = telescope.zp(b)
99 """
100 test = np.array(['u','g','g'])
101 index = np.argwhere(zp['band'] == test[:,None])
102 print(index)
103 print(zp['zp'][index][:,1])
104 print(toto)
105 """
107 def __call__(self, obs, gen_par=None, bands='grizy'):
108 """ Simulation of the light curve
110 Parameters
111 ----------------
112 obs: array
113 array of observations
114 gen_par: array, opt
115 simulation parameters (default: None)
116 bands: str, opt
117 filters to consider for simulation (default: grizy)
118 Returns
119 ------------
120 astropy table with:
121 columns: band, flux, fluxerr, snr_m5,flux_e,zp,zpsys,time
122 metadata : SNID,RA,Dec,DayMax,X1,Color,z
123 """
125 if len(obs) == 0:
126 return None
128 tab_tot = pd.DataFrame()
130 # multiprocessing here: one process (processBand) per band
132 for band in bands:
133 idx = obs[self.filterCol] == band
134 # print('multiproc',band,j,len(obs[idx]))
135 if len(obs[idx]) > 0:
136 res = self.processBand(obs[idx], band, gen_par)
137 tab_tot = tab_tot.append(res, ignore_index=True)
139 # return produced LC
140 return tab_tot
143 def processBand(self, sel_obs, band, gen_par, j=-1, output_q=None):
144 """ LC simulation of a set of obs corresponding to a band
145 The idea is to use python broadcasting so as to estimate
146 all the requested values (flux, flux error, Fisher components, ...)
147 in a single path (i.e no loop!)
148 Parameters
149 ---------------
150 sel_obs: array
151 array of observations
152 band: str
153 band of observations
154 gen_par: array
155 simulation parameters
156 j: int, opt
157 index for multiprocessing (default: -1)
158 output_q: multiprocessing.Queue(),opt
159 queue for multiprocessing (default: None)
160 Returns
161 -------
162 astropy table with fields corresponding to LC components
163 """
165 # method used for interpolation
166 method = 'linear'
167 interpType = 'regular'
169 # if there are no observations in this filter: return None
170 if len(sel_obs) == 0:
171 if output_q is not None:
172 output_q.put({j: None})
173 else:
174 return None
176 # Get the fluxes (from griddata reference)
178 # xi = MJD-T0
179 xi = sel_obs[self.mjdCol]-gen_par['daymax'][:, np.newaxis]
181 # yi = redshift simulated values
182 # requested to avoid interpolation problems near boundaries
183 yi = np.round(gen_par['z'], 4)
184 # yi = gen_par['z']
186 # p = phases of LC points = xi/(1.+z)
187 p = xi/(1.+yi[:, np.newaxis])
188 yi_arr = np.ones_like(p)*yi[:, np.newaxis]
190 if interpType == 'regular':
192 pts = (p, yi_arr)
193 fluxes_obs = self.reference_lc.flux[band](pts)
194 fluxes_obs_err = self.reference_lc.fluxerr[band](pts)
196 # Fisher components estimation
198 dFlux = {}
200 # loop on Fisher parameters
201 for val in self.param_Fisher:
202 dFlux[val] = self.reference_lc.param[band][val](pts)
203 # get the reference components
204 # z_c = self.reference_lc.lc_ref[band]['d'+val]
205 # get Fisher components from interpolation
206 # dFlux[val] = griddata((x, y), z_c, (p, yi_arr),
207 # method=method, fill_value=0.)
209 # replace crazy fluxes by dummy values
210 fluxes_obs[fluxes_obs <= 0.] = 1.e-10
211 fluxes_obs_err[fluxes_obs_err <= 0.] = 1.e-10
213 # Fisher matrix components estimation
214 # loop on SN parameters (x0,x1,color)
215 # estimate: dF/dxi*dF/dxj/sigma_flux**2
216 Derivative_for_Fisher = {}
217 for ia, vala in enumerate(self.param_Fisher):
218 for jb, valb in enumerate(self.param_Fisher):
219 if jb >= ia:
220 Derivative_for_Fisher[vala +
221 valb] = dFlux[vala] * dFlux[valb]
223 # remove LC points outside the restframe phase range
224 min_rf_phase = gen_par['min_rf_phase'][:, np.newaxis]
225 max_rf_phase = gen_par['max_rf_phase'][:, np.newaxis]
226 flag = (p >= min_rf_phase) & (p <= max_rf_phase)
228 # remove LC points outside the (blue-red) range
229 mean_restframe_wavelength = np.array(
230 [self.telescope.mean_wavelength[band]]*len(sel_obs))
231 mean_restframe_wavelength = np.tile(
232 mean_restframe_wavelength, (len(gen_par), 1))/(1.+gen_par['z'][:, np.newaxis])
233 flag &= (mean_restframe_wavelength > self.blue_cutoff) & (
234 mean_restframe_wavelength < self.red_cutoff)
236 flag_idx = np.argwhere(flag)
238 # Correct fluxes_err (m5 in generation probably different from m5 obs)
240 # gamma_obs = self.telescope.gamma(
241 # sel_obs[self.m5Col], [band]*len(sel_obs), sel_obs[self.exptimeCol])
243 gamma_obs = self.reference_lc.gamma[band](
244 (sel_obs[self.m5Col], sel_obs[self.exptimeCol]/sel_obs[self.nexpCol], sel_obs[self.nexpCol]))
246 mag_obs = -2.5*np.log10(fluxes_obs/3631.)
248 m5 = np.asarray([self.reference_lc.m5_ref[band]]*len(sel_obs))
250 gammaref = np.asarray([self.reference_lc.gamma_ref[band]]*len(sel_obs))
252 m5_tile = np.tile(m5, (len(p), 1))
254 srand_ref = self.srand(
255 np.tile(gammaref, (len(p), 1)), mag_obs, m5_tile)
257 srand_obs = self.srand(np.tile(gamma_obs, (len(p), 1)), mag_obs, np.tile(
258 sel_obs[self.m5Col], (len(p), 1)))
260 correct_m5 = srand_ref/srand_obs
262 """
263 print(band, gammaref, gamma_obs, m5,
264 sel_obs[self.m5Col], sel_obs[self.exptimeCol])
265 """
266 fluxes_obs_err = fluxes_obs_err/correct_m5
268 # now apply the flag to select LC points
269 fluxes = np.ma.array(fluxes_obs, mask=~flag)
270 fluxes_err = np.ma.array(fluxes_obs_err, mask=~flag)
271 phases = np.ma.array(p, mask=~flag)
272 snr_m5 = np.ma.array(fluxes_obs/fluxes_obs_err, mask=~flag)
274 nvals = len(phases)
276 obs_time = np.ma.array(
277 np.tile(sel_obs[self.mjdCol], (nvals, 1)), mask=~flag)
278 seasons = np.ma.array(
279 np.tile(sel_obs[self.seasonCol], (nvals, 1)), mask=~flag)
281 z_vals = gen_par['z'][flag_idx[:, 0]]
282 daymax_vals = gen_par['daymax'][flag_idx[:, 0]]
283 mag_obs = np.ma.array(mag_obs, mask=~flag)
284 Fisher_Mat = {}
285 for key, vals in Derivative_for_Fisher.items():
286 Fisher_Mat[key] = np.ma.array(vals, mask=~flag)
288 # Store in a panda dataframe
289 lc = pd.DataFrame()
291 ndata = len(fluxes_err[~fluxes_err.mask])
293 if ndata > 0:
295 lc['flux'] = fluxes[~fluxes.mask]
296 lc['fluxerr'] = fluxes_err[~fluxes_err.mask]
297 lc['phase'] = phases[~phases.mask]
298 lc['snr_m5'] = snr_m5[~snr_m5.mask]
299 lc['time'] = obs_time[~obs_time.mask]
300 lc['mag'] = mag_obs[~mag_obs.mask]
301 lc['band'] = ['LSST::'+band]*len(lc)
302 lc.loc[:, 'zp'] = self.zp[band]
303 lc['season'] = seasons[~seasons.mask]
304 lc['season'] = lc['season'].astype(int)
305 lc['z'] = z_vals
306 lc['daymax'] = daymax_vals
307 for key, vals in Fisher_Mat.items():
308 lc.loc[:, 'F_{}'.format(
309 key)] = vals[~vals.mask]/(lc['fluxerr'].values**2)
310 # lc.loc[:, 'F_{}'.format(key)] = 999.
311 lc.loc[:, 'x1'] = self.x1
312 lc.loc[:, 'color'] = self.color
314 lc.loc[:, 'n_aft'] = (np.sign(lc['phase']) == 1) & (
315 lc['snr_m5'] >= self.snr_min)
316 lc.loc[:, 'n_bef'] = (np.sign(lc['phase'])
317 == -1) & (lc['snr_m5'] >= self.snr_min)
319 lc.loc[:, 'n_phmin'] = (lc['phase'] <= -5.)
320 lc.loc[:, 'n_phmax'] = (lc['phase'] >= 20)
322 # transform boolean to int because of some problems in the sum()
324 for colname in ['n_aft', 'n_bef', 'n_phmin', 'n_phmax']:
325 lc.loc[:, colname] = lc[colname].astype(int)
327 """
328 idb = (lc['z'] > 0.65) & (lc['z'] < 0.9)
329 print(lc[idb][['z', 'ratio', 'm5', 'flux_e_sec', 'snr_m5']])
330 """
331 if output_q is not None:
332 output_q.put({j: lc})
333 else:
334 return lc
336 def srand(self, gamma, mag, m5):
337 """
338 Method to estimate :math:`srand=\sqrt((0.04-\gamma)*x+\gamma*x^2)`
339 with :math:`x = 10^{0.4*(m-m_5)}`
340 Parameters
341 ---------------
342 gamma: float
343 gamma value
344 mag: float
345 magnitude
346 m5: float
347 fiveSigmaDepth value
348 Returns
349 ----------
350 srand = np.sqrt((0.04-gamma)*x+gamma*x**2)
351 with x = 10**(0.4*(mag-m5))
352 """
354 x = 10**(0.4*(mag-m5))
355 return np.sqrt((0.04-gamma)*x+gamma*x**2)
358class Throughputs(object):
359 """ class to handle instrument throughput
360 Parameters
361 -------------
362 through_dir : str, opt
363 throughput directory
364 Default : LSST_THROUGHPUTS_BASELINE
365 atmos_dir : str, opt
366 directory of atmos files
367 Default : THROUGHPUTS_DIR
368 telescope_files : list(str),opt
369 list of of throughput files
370 Default : ['detector.dat', 'lens1.dat','lens2.dat',
371 'lens3.dat','m1.dat', 'm2.dat', 'm3.dat']
372 filterlist: list(str), opt
373 list of filters to consider
374 Default : 'ugrizy'
375 wave_min : float, opt
376 min wavelength for throughput
377 Default : 300
378 wave_max : float, opt
379 max wavelength for throughput
380 Default : 1150
381 atmos : bool, opt
382 to include atmosphere affects
383 Default : True
384 aerosol : bool, opt
385 to include aerosol effects
386 Default : True
387 Returns
388 ---------
389 Accessible throughputs (per band):
390 lsst_system: system throughput (lens+mirrors+filters)
391 lsst_atmos: lsst_system+atmosphere
392 lsst_atmos_aerosol: lsst_system+atmosphere+aerosol
393 """
395 def __init__(self, **kwargs):
397 params = {}
398 params['through_dir'] = 'LSST_THROUGHPUTS_BASELINE'
399 params['atmos_dir'] = 'THROUGHPUTS_DIR'
400 params['atmos'] = True
401 params['aerosol'] = True
402 params['telescope_files'] = ['detector.dat', 'lens1.dat',
403 'lens2.dat', 'lens3.dat',
404 'm1.dat', 'm2.dat', 'm3.dat']
405 params['filterlist'] = 'ugrizy'
406 params['wave_min'] = 300.
407 params['wave_max'] = 1150.
408 for par in ['through_dir', 'atmos_dir', 'atmos', 'aerosol',
409 'telescope_files', 'filterlist', 'wave_min', 'wave_max']:
410 if par in kwargs.keys():
411 params[par] = kwargs[par]
412 # params[par]=str(kwargs[par])
414 self.atmos = params['atmos']
415 self.throughputsDir = os.getenv(params['through_dir'])
416 if os.path.exists(os.path.join
417 (os.getenv(params['atmos_dir']), 'atmos')):
418 self.atmosDir = os.path.join(
419 os.getenv(params['atmos_dir']), 'atmos')
420 else:
421 self.atmosDir = os.getenv(params['atmos_dir'])
423 self.telescope_files = params['telescope_files']
424 self.filter_files = ['filter_'+f+'.dat' for f in params['filterlist']]
425 if 'filter_files' in kwargs.keys():
426 self.filter_files = kwargs['filter_files']
427 self.wave_min = params['wave_min']
428 self.wave_max = params['wave_max']
430 self.filterlist = params['filterlist']
431 self.filtercolors = {'u': 'b', 'g': 'c',
432 'r': 'g', 'i': 'y', 'z': 'r', 'y': 'm'}
434 self.lsst_std = {}
435 self.lsst_system = {}
436 self.mean_wavelength = {}
437 self.lsst_detector = {}
438 self.lsst_atmos = {}
439 self.lsst_atmos_aerosol = {}
440 self.airmass = -1.
441 self.aerosol_b = params['aerosol']
442 self.Load_System()
443 self.Load_DarkSky()
445 if params['atmos']:
446 self.Load_Atmosphere()
447 else:
448 for f in self.filterlist:
449 self.lsst_atmos[f] = self.lsst_system[f]
450 self.lsst_atmos_aerosol[f] = self.lsst_system[f]
452 # self.lsst_telescope={}
454 # self.Load_Telescope()
456 self.Mean_Wave()
458 @property
459 def system(self):
460 return self.lsst_system
462 @property
463 def telescope(self):
464 return self.lsst_telescope
466 @property
467 def atmosphere(self):
468 return self.lsst_atmos
470 @property
471 def aerosol(self):
472 return self.lsst_atmos_aerosol
474 def Load_System(self):
475 """ Load files required to estimate throughputs
476 """
478 for f in self.filterlist:
479 self.lsst_std[f] = Bandpass()
480 self.lsst_system[f] = Bandpass()
482 telfiles = ''
483 if len(self.telescope_files) > 0:
484 index = [i for i, x in enumerate(
485 self.filter_files) if f+'.dat' in x]
486 telfiles = self.telescope_files+[self.filter_files[index[0]]]
487 else:
488 telfiles = self.filter_files
489 self.lsst_system[f].readThroughputList(telfiles,
490 rootDir=self.throughputsDir,
491 wavelen_min=self.wave_min,
492 wavelen_max=self.wave_max)
493 """
494 def Load_Telescope(self):
495 for system in self.telescope_files+self.filter_files:
496 self.lsst_telescope[system] = Bandpass()
497 self.lsst_telescope[system].readThroughputList([system],
498 rootDir=self.throughputsDir,
499 wavelen_min=self.wave_min,
500 wavelen_max=self.wave_max)
501 """
503 def Load_DarkSky(self):
504 """ Load DarkSky
505 """
506 self.darksky = Sed()
507 self.darksky.readSED_flambda(os.path.join(
508 self.throughputsDir, 'darksky.dat'))
510 def Load_Atmosphere(self, airmass=1.2):
511 """ Load atmosphere files
512 and convolve with transmissions
513 Parameters
514 --------------
515 airmass : float,opt
516 airmass value
517 Default : 1.2
518 """
519 self.airmass = airmass
520 if self.airmass > 0.:
521 atmosphere = Bandpass()
522 path_atmos = os.path.join(
523 self.atmosDir, 'atmos_%d.dat' % (self.airmass*10))
524 if os.path.exists(path_atmos):
525 atmosphere.readThroughput(os.path.join(
526 self.atmosDir, 'atmos_%d.dat' % (self.airmass*10)))
527 else:
528 atmosphere.readThroughput(
529 os.path.join(self.atmosDir, 'atmos.dat'))
530 self.atmos = Bandpass(wavelen=atmosphere.wavelen, sb=atmosphere.sb)
532 for f in self.filterlist:
533 wavelen, sb = self.lsst_system[f].multiplyThroughputs(
534 atmosphere.wavelen, atmosphere.sb)
535 self.lsst_atmos[f] = Bandpass(wavelen=wavelen, sb=sb)
537 if self.aerosol_b:
538 atmosphere_aero = Bandpass()
539 atmosphere_aero.readThroughput(os.path.join(
540 self.atmosDir, 'atmos_%d_aerosol.dat' % (self.airmass*10)))
541 self.atmos_aerosol = Bandpass(
542 wavelen=atmosphere_aero.wavelen, sb=atmosphere_aero.sb)
544 for f in self.filterlist:
545 wavelen, sb = self.lsst_system[f].multiplyThroughputs(
546 atmosphere_aero.wavelen, atmosphere_aero.sb)
547 self.lsst_atmos_aerosol[f] = Bandpass(
548 wavelen=wavelen, sb=sb)
549 else:
550 for f in self.filterlist:
551 self.lsst_atmos[f] = self.lsst_system[f]
552 self.lsst_atmos_aerosol[f] = self.lsst_system[f]
554 def Plot_Throughputs(self):
555 """ Plot the throughputs
556 """
557 # colors=['b','g','r','m','c',[0.8,0,0]]
558 style = [',', ',', ',', ',']
559 for i, band in enumerate(self.filterlist):
561 plt.plot(self.lsst_system[band].wavelen,
562 self.lsst_system[band].sb,
563 linestyle='--', color=self.filtercolors[band],
564 label='%s - syst' % (band))
565 plt.plot(self.lsst_atmos[band].wavelen,
566 self.lsst_atmos[band].sb,
567 linestyle='-.', color=self.filtercolors[band],
568 label='%s - syst+atm' % (band))
569 if len(self.lsst_atmos_aerosol) > 0:
570 plt.plot(self.lsst_atmos_aerosol[band].wavelen,
571 self.lsst_atmos_aerosol[band].sb,
572 linestyle='-', color=self.filtercolors[band],
573 label='%s - syst+atm+aero' % (band))
575 plt.plot(self.atmos.wavelen, self.atmos.sb, 'k:',
576 label='X =%.1f atmos' % (self.airmass), linestyle='-')
577 if len(self.lsst_atmos_aerosol) > 0:
578 plt.plot(self.atmos_aerosol.wavelen, self.atmos_aerosol.sb, 'k:',
579 label='X =%.1f atm+aero' % (self.airmass), linestyle='--')
580 # plt.legend(loc=(0.85, 0.1), fontsize='smaller',
581 # fancybox=True, numpoints=1)
583 plt.legend(loc=(0.82, 0.1), fancybox=True, numpoints=1)
585 plt.xlabel('Wavelength (nm)')
586 plt.ylabel('Sb (0-1)')
587 plt.title('System throughput')
588 # plt.show()
590 def Plot_DarkSky(self):
591 """ Plot darksky
592 """
593 # self.Load_DarkSky()
594 plt.plot(self.darksky.wavelen,
595 self.darksky.flambda, 'k:', linestyle='-')
596 plt.xlabel('Wavelength (nm)')
597 plt.ylabel('flambda (ergs/cm$^2$/s/nm)')
598 plt.title('Dark Sky SED')
599 plt.show()
601 def Plot_Throughputs_Spectrum(self, wavelength, fluxes, z):
602 """ Plot throughput spectrum
603 """
604 fig, ax1 = plt.subplots()
605 style = [',', ',', ',', ',']
607 for i, band in enumerate(self.filterlist):
609 plt.plot(self.lsst_system[band].wavelen,
610 self.lsst_system[band].sb,
611 linestyle='-', color=self.filtercolors[band],
612 label='%s - system' % (band))
614 def Mean_Wave(self):
615 """ Estimate mean wave
616 """
617 for band in self.filterlist:
618 self.mean_wavelength[band] = np.sum(
619 self.lsst_atmos[band].wavelen*self.lsst_atmos[band].sb)\
620 / np.sum(self.lsst_atmos[band].sb)
622# decorator to access parameters of the class
625def get_val_decor(func):
626 @wraps(func)
627 def func_deco(theclass, what, xlist):
628 for x in xlist:
629 if x not in theclass.data[what].keys():
630 func(theclass, what, x)
631 return func_deco
634class Telescope(Throughputs):
635 """ Telescope class
636 inherits from Throughputs
637 estimate quantities defined in LSE-40
638 The following quantities are accessible:
639 mag_sky: sky magnitude
640 m5: 5-sigma depth
641 Sigmab: see eq. (36) of LSE-40
642 zp: see eq. (43) of LSE-40
643 counts_zp:
644 Skyb: see eq. (40) of LSE-40
645 flux_sky:
646 Parameters
647 -------------
648 through_dir : str, opt
649 throughput directory
650 Default : LSST_THROUGHPUTS_BASELINE
651 atmos_dir : str, opt
652 directory of atmos files
653 Default : THROUGHPUTS_DIR
654 telescope_files : list(str),opt
655 list of of throughput files
656 Default : ['detector.dat', 'lens1.dat','lens2.dat',
657 'lens3.dat','m1.dat', 'm2.dat', 'm3.dat']
658 filterlist: list(str), opt
659 list of filters to consider
660 Default : 'ugrizy'
661 wave_min : float, opt
662 min wavelength for throughput
663 Default : 300
664 wave_max : float, opt
665 max wavelength for throughput
666 Default : 1150
667 atmos : bool, opt
668 to include atmosphere affects
669 Default : True
670 aerosol : bool, opt
671 to include aerosol effects
672 Default : True
673 airmass : float, opt
674 airmass value
675 Default : 1.
676 Returns
677 ---------
678 Accessible throughputs (per band, from Throughput class):
679 lsst_system: system throughput (lens+mirrors+filters)
680 lsst_atmos: lsst_system+atmosphere
681 lsst_atmos_aerosol: lsst_system+atmosphere+aerosol
682 """
684 def __init__(self, name='unknown', airmass=1., **kwargs):
686 self.name = name
687 Throughputs.__init__(self, **kwargs)
689 params = ['mag_sky', 'm5', 'FWHMeff', 'Tb',
690 'Sigmab', 'zp', 'counts_zp', 'Skyb', 'flux_sky']
692 self.data = {}
693 for par in params:
694 self.data[par] = {}
696 self.data['FWHMeff'] = dict(
697 zip('ugrizy', [0.92, 0.87, 0.83, 0.80, 0.78, 0.76]))
699 # self.atmos = atmos
701 self.Load_Atmosphere(airmass)
703 @get_val_decor
704 def get(self, what, band):
705 """
706 Decorator to access quantities
707 Parameters
708 ---------------
709 what: str
710 parameter to estimate
711 band: str
712 filter
713 """
714 filter_trans = self.system[band]
715 wavelen_min, wavelen_max, wavelen_step = filter_trans.getWavelenLimits(
716 None, None, None)
718 bandpass = Bandpass(wavelen=filter_trans.wavelen, sb=filter_trans.sb)
720 flatSedb = Sed()
721 flatSedb.setFlatSED(wavelen_min, wavelen_max, wavelen_step)
722 flux0b = np.power(10., -0.4*self.mag_sky(band))
723 flatSedb.multiplyFluxNorm(flux0b)
724 photParams = PhotometricParameters(bandpass=band)
725 norm = photParams.platescale**2/2.*photParams.exptime/photParams.gain
726 trans = filter_trans
728 if self.atmos:
729 trans = self.atmosphere[band]
730 self.data['m5'][band] = SignalToNoise.calcM5(
731 flatSedb, trans, filter_trans,
732 photParams=photParams,
733 FWHMeff=self.FWHMeff(band))
734 adu_int = flatSedb.calcADU(bandpass=trans, photParams=photParams)
735 self.data['flux_sky'][band] = adu_int*norm
737 @get_val_decor
738 def get_inputs(self, what, band):
739 """
740 decorator to access Tb, Sigmab, mag_sky
741 Parameters
742 ---------------
743 what: str
744 parameter to estimate
745 band: str
746 filter
747 """
748 myup = self.Calc_Integ_Sed(self.darksky, self.system[band])
749 self.data['Tb'][band] = self.Calc_Integ(self.atmosphere[band])
750 self.data['Sigmab'][band] = self.Calc_Integ(self.system[band])
751 self.data['mag_sky'][band] = -2.5 * \
752 np.log10(myup/(3631.*self.Sigmab(band)))
754 @get_val_decor
755 def get_zp(self, what, band):
756 """
757 decorator get zero points
758 formula used here are extracted from LSE-40
759 Parameters
760 ---------------
761 what: str
762 parameter to estimate
763 band: str
764 filter
765 """
766 photParams = PhotometricParameters(bandpass=band)
767 Diameter = 2.*np.sqrt(photParams.effarea*1.e-4 /
768 np.pi) # diameter in meter
769 Cte = 3631.*np.pi*Diameter**2*2.*photParams.exptime/4/h/1.e36
771 self.data['Skyb'][band] = Cte*np.power(Diameter/6.5, 2.)\
772 * np.power(2.*photParams.exptime/30., 2.)\
773 * np.power(photParams.platescale, 2.)\
774 * 10.**0.4*(25.-self.mag_sky(band))\
775 * self.Sigmab(band)
777 Zb = 181.8*np.power(Diameter/6.5, 2.)*self.Tb(band)
778 mbZ = 25.+2.5*np.log10(Zb)
779 filtre_trans = self.system[band]
780 wavelen_min, wavelen_max, wavelen_step = filtre_trans.getWavelenLimits(
781 None, None, None)
782 bandpass = Bandpass(wavelen=filtre_trans.wavelen, sb=filtre_trans.sb)
783 flatSed = Sed()
784 flatSed.setFlatSED(wavelen_min, wavelen_max, wavelen_step)
785 flux0 = np.power(10., -0.4*mbZ)
786 flatSed.multiplyFluxNorm(flux0)
787 photParams = PhotometricParameters(bandpass=band)
788 # number of counts for exptime
789 counts = flatSed.calcADU(bandpass, photParams=photParams)
790 self.data['zp'][band] = mbZ
791 self.data['counts_zp'][band] = counts/2.*photParams.exptime
793 def return_value(self, what, band):
794 """
795 accessor
796 Parameters
797 ---------------
798 what: str
799 parameter to estimate
800 band: str
801 filter
802 """
803 if len(band) > 1:
804 return self.data[what]
805 else:
806 return self.data[what][band]
808 def m5(self, filtre):
809 """m5 accessor
810 """
811 self.get('m5', filtre)
812 return self.return_value('m5', filtre)
814 def Tb(self, filtre):
815 """Tb accessor
816 """
817 self.get_inputs('Tb', filtre)
818 return self.return_value('Tb', filtre)
820 def mag_sky(self, filtre):
821 """mag_sky accessor
822 """
823 self.get_inputs('mag_sky', filtre)
824 return self.return_value('mag_sky', filtre)
826 def Sigmab(self, filtre):
827 """
828 Sigmab accessor
829 Parameters
830 ----------------
831 band: str
832 filter
833 """
834 self.get_inputs('Sigmab', filtre)
835 return self.return_value('Sigmab', filtre)
837 def zp(self, filtre):
838 """
839 zp accessor
840 Parameters
841 ----------------
842 band: str
843 filter
844 """
845 self.get_zp('zp', filtre)
846 return self.return_value('zp', filtre)
848 def FWHMeff(self, filtre):
849 """
850 FWHMeff accessor
851 Parameters
852 ----------------
853 band: str
854 filter
855 """
856 return self.return_value('FWHMeff', filtre)
858 def Calc_Integ(self, bandpass):
859 """
860 integration over bandpass
861 Parameters
862 --------------
863 bandpass : float
864 Returns
865 ---------
866 integration
867 """
868 resu = 0.
869 dlam = 0
870 for i, wave in enumerate(bandpass.wavelen):
871 if i < len(bandpass.wavelen)-1:
872 dlam = bandpass.wavelen[i+1]-wave
873 resu += dlam*bandpass.sb[i]/wave
874 # resu+=dlam*bandpass.sb[i]
876 return resu
878 def Calc_Integ_Sed(self, sed, bandpass, wavelen=None, fnu=None):
879 """
880 SED integration
881 Parameters
882 --------------
883 sed : float
884 sed to integrate
885 bandpass : float
886 bandpass
887 wavelength : float, opt
888 wavelength values
889 Default : None
890 fnu : float, opt
891 fnu values
892 Default : None
893 Returns
894 ----------
895 integrated sed over the bandpass
896 """
897 use_self = sed._checkUseSelf(wavelen, fnu)
898 # Use self values if desired, otherwise use values passed to function.
899 if use_self:
900 # Calculate fnu if required.
901 if sed.fnu is None:
902 # If fnu not present, calculate. (does not regrid).
903 sed.flambdaTofnu()
904 wavelen = sed.wavelen
905 fnu = sed.fnu
906 # Make sure wavelen/fnu are on the same wavelength grid as bandpass.
907 wavelen, fnu = sed.resampleSED(
908 wavelen, fnu, wavelen_match=bandpass.wavelen)
910 # Calculate the number of photons.
911 nphoton = (fnu / wavelen * bandpass.sb).sum()
912 dlambda = wavelen[1] - wavelen[0]
913 return nphoton * dlambda
915 def flux_to_mag(self, flux, band, zp=None):
916 """
917 Flux to magnitude conversion
918 Parameters
919 --------------
920 flux : float
921 input fluxes
922 band : str
923 input band
924 zp : float, opt
925 zeropoints
926 Default : None
927 Returns
928 ---------
929 magnitudes
930 """
931 if zp is None:
932 zp = self.zero_points(band)
933 # print 'zp',zp,band
934 m = -2.5 * np.log10(flux) + zp
935 return m
937 def mag_to_flux(self, mag, band, zp=None):
938 """
939 Magnitude to flux conversion
940 Parameters
941 --------------
942 mag : float
943 input mags
944 band : str
945 input band
946 zp : float, opt
947 zeropoints
948 Default : None
949 Returns
950 ---------
951 fluxes
952 """
953 if zp is None:
954 zp = self.zero_points(band)
955 return np.power(10., -0.4 * (mag-zp))
957 def zero_points(self, band):
958 """
959 Zero points estimation
960 Parameters
961 --------------
962 band : list(str)
963 list of bands
964 Returns
965 ---------
966 array of zp
967 """
968 return np.asarray([self.zp[b] for b in band])
970 def mag_to_flux_e_sec(self, mag, band, exptime):
971 """
972 Mag to flux (in photoelec/sec) conversion
973 Parameters
974 --------------
975 mag : float
976 input magnitudes
977 band : str
978 input bands
979 exptime : float
980 input exposure times
981 Returns
982 ----------
983 counts : float
984 number of ADU counts
985 e_per_sec : float
986 flux in photoelectron per sec.
987 """
988 if not hasattr(mag, '__iter__'):
989 wavelen_min, wavelen_max, wavelen_step = self.atmosphere[band].getWavelenLimits(
990 None, None, None)
991 sed = Sed()
992 sed.setFlatSED()
993 flux0 = 3631.*10**(-0.4*mag) # flux in Jy
994 flux0 = sed.calcFluxNorm(mag, self.atmosphere[band])
995 sed.multiplyFluxNorm(flux0)
996 photParams = PhotometricParameters(nexp=exptime/15.)
997 counts = sed.calcADU(
998 bandpass=self.atmosphere[band], photParams=photParams)
999 e_per_sec = counts
1000 e_per_sec /= exptime/photParams.gain
1001 # print('hello',photParams.gain,exptime)
1002 return counts, e_per_sec
1003 else:
1004 return np.asarray([self.mag_to_flux_e_sec(m, b, expt) for m, b, expt in zip(mag, band, exptime)])
1006 def gamma(self, mag, band, exptime):
1007 """
1008 gamma parameter estimation
1009 cf eq(5) of the paper LSST : from science drivers to reference design and anticipated data products
1010 with sigma_rand = 0.2 and m=m5
1011 Parameters
1012 --------------
1013 mag : float
1014 magnitudes
1015 band : str
1016 band
1017 exptime : float
1018 exposure time
1019 Returns
1020 ----------
1021 gamma (float)
1022 """
1024 if not hasattr(mag, '__iter__'):
1025 photParams = PhotometricParameters(nexp=exptime/15.)
1026 counts, e_per_sec = self.mag_to_flux_e_sec(mag, band, exptime)
1027 return 0.04-1./(photParams.gain*counts)
1028 else:
1029 return np.asarray([self.gamma(m, b, e) for m, b, e in zip(mag, band, exptime)])
1032class Load_Reference:
1033 """
1034 class to load template files requested for LCFast
1035 These files should be stored in a reference_files directory
1037 Parameters
1038 ---------------
1039 server: str, opt
1040 where to get the files (default: https://me.lsst.eu/gris/DESC_SN_pipeline/Reference_Files)
1041 templateDir: str, opt
1042 where to put the files (default: reference_files)
1044 """
1046 def __init__(self, server='https://me.lsst.eu/gris/DESC_SN_pipeline',
1047 templateDir=None):
1049 if templateDir is None:
1050 sims_maf_contrib_dir = os.getenv("SIMS_MAF_CONTRIB_DIR")
1051 templateDir = os.path.join(sims_maf_contrib_dir, 'data/SNe_data')
1053 self.server = server
1054 # define instrument
1055 self.Instrument = {}
1056 self.Instrument['name'] = 'LSST' # name of the telescope (internal)
1057 # dir of throughput
1058 self.Instrument['throughput_dir'] = 'LSST_THROUGHPUTS_BASELINE'
1059 self.Instrument['atmos_dir'] = 'THROUGHPUTS_DIR' # dir of atmos
1060 self.Instrument['airmass'] = 1.2 # airmass value
1061 self.Instrument['atmos'] = True # atmos
1062 self.Instrument['aerosol'] = False # aerosol
1064 x1_colors = [(-2.0, 0.2), (0.0, 0.0)]
1066 lc_reference = {}
1068 # create this directory if it does not exist
1069 if not os.path.isdir(templateDir):
1070 os.system('mkdir {}'.format(templateDir))
1072 list_files = ['gamma.hdf5']
1073 for j in range(len(x1_colors)):
1074 x1 = x1_colors[j][0]
1075 color = x1_colors[j][1]
1076 fname = 'LC_{}_{}_380.0_800.0_ebvofMW_0.0_vstack.hdf5'.format(
1077 x1, color)
1078 list_files += [fname]
1080 self.check_grab(templateDir, list_files)
1082 # gamma_reference
1083 self.gamma_reference = '{}/gamma.hdf5'.format(templateDir)
1085 # print('Loading reference files')
1087 resultdict = {}
1089 for j in range(len(x1_colors)):
1090 x1 = x1_colors[j][0]
1091 color = x1_colors[j][1]
1092 fname = '{}/LC_{}_{}_380.0_800.0_ebvofMW_0.0_vstack.hdf5'.format(
1093 templateDir, x1, color)
1094 resultdict[j] = self.load(fname)
1096 for j in range(len(x1_colors)):
1097 if resultdict[j] is not None:
1098 lc_reference[x1_colors[j]] = resultdict[j]
1100 self.ref = lc_reference
1102 def load(self, fname):
1103 """
1104 Method to load reference files
1106 Parameters
1107 ---------------
1108 fname: str
1109 file name
1111 Returns
1112 -----------
1114 """
1115 lc_ref = GetReference(
1116 fname, self.gamma_reference, self.Instrument)
1118 return lc_ref
1120 def check_grab(self, templateDir, listfiles):
1121 """
1122 Method that check if files are on disk.
1123 If not: grab them from a server (self.server)
1125 Parameters
1126 ---------------
1127 templateDir: str
1128 directory where files are (or will be)
1129 listfiles: list(str)
1130 list of files that are (will be) in templateDir
1131 """
1133 for fi in listfiles:
1134 # check whether the file is available; if not-> get it!
1135 fname = '{}/{}'.format(templateDir, fi)
1136 if not os.path.isfile(fname):
1137 if 'gamma' in fname:
1138 fullname = '{}/reference_files/{}'.format(self.server, fi)
1139 else:
1140 fullname = '{}/Template_LC/{}'.format(self.server, fi)
1141 print('wget path:', fullname)
1142 cmd = 'wget --no-clobber --no-verbose {} --directory-prefix {}'.format(
1143 fullname, templateDir)
1144 os.system(cmd)
1147class GetReference:
1148 """
1149 Class to load reference data
1150 used for the fast SN simulator
1151 Parameters
1152 ----------------
1153 lcName: str
1154 name of the reference file to load (lc)
1155 gammaName: str
1156 name of the reference file to load (gamma)
1157 tel_par: dict
1158 telescope parameters
1159 param_Fisher : list(str),opt
1160 list of SN parameter for Fisher estimation to consider
1161 (default: ['x0', 'x1', 'color', 'daymax'])
1162 Returns
1163 -----------
1164 The following dict can be accessed:
1165 mag_to_flux_e_sec : Interp1D of mag to flux(e.sec-1) conversion
1166 flux : dict of RegularGridInterpolator of fluxes (key: filters, (x,y)=(phase, z), result=flux)
1167 fluxerr : dict of RegularGridInterpolator of flux errors (key: filters, (x,y)=(phase, z), result=fluxerr)
1168 param : dict of dict of RegularGridInterpolator of flux derivatives wrt SN parameters
1169 (key: filters plus param_Fisher parameters; (x,y)=(phase, z), result=flux derivatives)
1170 gamma : dict of RegularGridInterpolator of gamma values (key: filters)
1171 """""
1173 def __init__(self, lcName, gammaName, tel_par, param_Fisher=['x0', 'x1', 'color', 'daymax']):
1175 # Load the file - lc reference
1177 f = h5py.File(lcName, 'r')
1178 keys = list(f.keys())
1179 # lc_ref_tot = Table.read(filename, path=keys[0])
1180 lc_ref_tot = Table.from_pandas(pd.read_hdf(lcName))
1182 idx = lc_ref_tot['z'] > 0.005
1183 lc_ref_tot = np.copy(lc_ref_tot[idx])
1185 # telescope requested
1186 telescope = Telescope(name=tel_par['name'],
1187 throughput_dir=tel_par['throughput_dir'],
1188 atmos_dir=tel_par['atmos_dir'],
1189 atmos=tel_par['atmos'],
1190 aerosol=tel_par['aerosol'],
1191 airmass=tel_par['airmass'])
1193 # Load the file - gamma values
1194 if not os.path.exists(gammaName):
1195 print('gamma file {} does not exist')
1196 print('will generate it - few minutes')
1197 mag_range = np.arange(15., 38., 1.)
1198 exptimes = np.arange(1., 3000., 10.)
1199 Gamma('ugrizy', telescope, gammaName,
1200 mag_range=mag_range,
1201 exptimes=exptimes)
1202 print('end of gamma estimation')
1204 fgamma = h5py.File(gammaName, 'r')
1206 # Load references needed for the following
1207 self.lc_ref = {}
1208 self.gamma_ref = {}
1209 self.gamma = {}
1210 self.m5_ref = {}
1211 self.mag_to_flux_e_sec = {}
1213 self.flux = {}
1214 self.fluxerr = {}
1215 self.param = {}
1217 bands = np.unique(lc_ref_tot['band'])
1218 mag_range = np.arange(10., 38., 0.01)
1219 # exptimes = np.linspace(15.,30.,2)
1220 # exptimes = [15.,30.,60.,100.]
1222 # gammArray = self.loopGamma(bands, mag_range, exptimes,telescope)
1224 method = 'linear'
1226 # for each band: load data to be used for interpolation
1227 for band in bands:
1228 idx = lc_ref_tot['band'] == band
1229 lc_sel = Table(lc_ref_tot[idx])
1231 lc_sel['z'] = lc_sel['z'].data.round(decimals=2)
1232 lc_sel['phase'] = lc_sel['phase'].data.round(decimals=1)
1234 """
1235 select phases between -20 and 50 only
1236 """
1237 idx = lc_sel['phase'] < 50.
1238 idx &= lc_sel['phase'] > -20.
1239 lc_sel = lc_sel[idx]
1241 fluxes_e_sec = telescope.mag_to_flux_e_sec(
1242 mag_range, [band]*len(mag_range), [30]*len(mag_range))
1243 self.mag_to_flux_e_sec[band] = interpolate.interp1d(
1244 mag_range, fluxes_e_sec[:, 1], fill_value=0., bounds_error=False)
1246 # these reference data will be used for griddata interp.
1247 self.lc_ref[band] = lc_sel
1248 self.gamma_ref[band] = lc_sel['gamma'][0]
1249 self.m5_ref[band] = np.unique(lc_sel['m5'])[0]
1251 # Another interpolator, faster than griddata: regulargridinterpolator
1253 # Fluxes and errors
1254 zmin, zmax, zstep, nz = self.limVals(lc_sel, 'z')
1255 phamin, phamax, phastep, npha = self.limVals(lc_sel, 'phase')
1257 zstep = np.round(zstep, 1)
1258 phastep = np.round(phastep, 1)
1260 zv = np.linspace(zmin, zmax, nz)
1261 # zv = np.round(zv,2)
1262 # print(band,zv)
1263 phav = np.linspace(phamin, phamax, npha)
1265 print('Loading ', lcName, band, len(lc_sel), npha, nz)
1266 index = np.lexsort((lc_sel['z'], lc_sel['phase']))
1267 flux = np.reshape(lc_sel[index]['flux'], (npha, nz))
1268 fluxerr = np.reshape(lc_sel[index]['fluxerr'], (npha, nz))
1270 self.flux[band] = RegularGridInterpolator(
1271 (phav, zv), flux, method=method, bounds_error=False, fill_value=0.)
1272 self.fluxerr[band] = RegularGridInterpolator(
1273 (phav, zv), fluxerr, method=method, bounds_error=False, fill_value=0.)
1275 # Flux derivatives
1276 self.param[band] = {}
1277 for par in param_Fisher:
1278 valpar = np.reshape(
1279 lc_sel[index]['d{}'.format(par)], (npha, nz))
1280 self.param[band][par] = RegularGridInterpolator(
1281 (phav, zv), valpar, method=method, bounds_error=False, fill_value=0.)
1283 # gamma estimator
1285 rec = Table.read(gammaName, path='gamma_{}'.format(band))
1287 rec['mag'] = rec['mag'].data.round(decimals=4)
1288 rec['single_exptime'] = rec['single_exptime'].data.round(
1289 decimals=4)
1291 magmin, magmax, magstep, nmag = self.limVals(rec, 'mag')
1292 expmin, expmax, expstep, nexpo = self.limVals(
1293 rec, 'single_exptime')
1294 nexpmin, nexpmax, nexpstep, nnexp = self.limVals(rec, 'nexp')
1295 mag = np.linspace(magmin, magmax, nmag)
1296 exp = np.linspace(expmin, expmax, nexpo)
1297 nexp = np.linspace(nexpmin, nexpmax, nnexp)
1299 index = np.lexsort(
1300 (rec['nexp'], np.round(rec['single_exptime'], 4), rec['mag']))
1301 gammab = np.reshape(rec[index]['gamma'], (nmag, nexpo, nnexp))
1302 fluxb = np.reshape(rec[index]['flux_e_sec'], (nmag, nexpo, nnexp))
1303 self.gamma[band] = RegularGridInterpolator(
1304 (mag, exp, nexp), gammab, method='linear', bounds_error=False, fill_value=0.)
1305 """
1306 self.mag_to_flux[band] = RegularGridInterpolator(
1307 (mag, exp, nexp), fluxb, method='linear', bounds_error=False, fill_value=0.)
1310 print('hello', rec.columns)
1311 rec['mag'] = rec['mag'].data.round(decimals=4)
1312 rec['exptime'] = rec['exptime'].data.round(decimals=4)
1314 magmin, magmax, magstep, nmag = self.limVals(rec, 'mag')
1315 expmin, expmax, expstep, nexp = self.limVals(rec, 'exptime')
1316 mag = np.linspace(magmin, magmax, nmag)
1317 exp = np.linspace(expmin, expmax, nexp)
1319 index = np.lexsort((np.round(rec['exptime'], 4), rec['mag']))
1320 gammab = np.reshape(rec[index]['gamma'], (nmag, nexp))
1321 self.gamma[band] = RegularGridInterpolator(
1322 (mag, exp), gammab, method=method, bounds_error=False, fill_value=0.)
1323 """
1324 # print(band, gammab, mag, exp)
1326 def limVals(self, lc, field):
1327 """ Get unique values of a field in a table
1328 Parameters
1329 ----------
1330 lc: Table
1331 astropy Table (here probably a LC)
1332 field: str
1333 name of the field of interest
1334 Returns
1335 -------
1336 vmin: float
1337 min value of the field
1338 vmax: float
1339 max value of the field
1340 vstep: float
1341 step value for this field (median)
1342 nvals: int
1343 number of unique values
1344 """
1346 lc.sort(field)
1347 vals = np.unique(lc[field].data.round(decimals=4))
1348 # print(vals)
1349 vmin = np.min(vals)
1350 vmax = np.max(vals)
1351 vstep = np.median(vals[1:]-vals[:-1])
1353 return vmin, vmax, vstep, len(vals)
1355 def Read_Ref(self, fi, j=-1, output_q=None):
1356 """" Load the reference file and
1357 make a single astopy Table from a set of.
1358 Parameters
1359 ----------
1360 fi: str,
1361 name of the file to be loaded
1362 Returns
1363 -------
1364 tab_tot: astropy table
1365 single table = vstack of all the tables in fi.
1366 """
1368 tab_tot = Table()
1369 """
1370 keys=np.unique([int(z*100) for z in zvals])
1371 print(keys)
1372 """
1373 f = h5py.File(fi, 'r')
1374 keys = f.keys()
1375 zvals = np.arange(0.01, 0.9, 0.01)
1376 zvals_arr = np.array(zvals)
1378 for kk in keys:
1380 tab_b = Table.read(fi, path=kk)
1382 if tab_b is not None:
1383 tab_tot = vstack([tab_tot, tab_b], metadata_conflicts='silent')
1384 """
1385 diff = tab_b['z']-zvals_arr[:, np.newaxis]
1386 # flag = np.abs(diff)<1.e-3
1387 flag_idx = np.where(np.abs(diff) < 1.e-3)
1388 if len(flag_idx[1]) > 0:
1389 tab_tot = vstack([tab_tot, tab_b[flag_idx[1]]])
1390 """
1392 """
1393 print(flag,flag_idx[1])
1394 print('there man',tab_b[flag_idx[1]])
1395 mtile = np.tile(tab_b['z'],(len(zvals),1))
1396 # print('mtile',mtile*flag)
1398 masked_array = np.ma.array(mtile,mask=~flag)
1400 print('resu masked',masked_array,masked_array.shape)
1401 print('hhh',masked_array[~masked_array.mask])
1404 for val in zvals:
1405 print('hello',tab_b[['band','z','time']],'and',val)
1406 if np.abs(np.unique(tab_b['z'])-val)<0.01:
1407 # print('loading ref',np.unique(tab_b['z']))
1408 tab_tot=vstack([tab_tot,tab_b])
1409 break
1410 """
1411 if output_q is not None:
1412 output_q.put({j: tab_tot})
1413 else:
1414 return tab_tot
1416 def Read_Multiproc(self, tab):
1417 """
1418 Multiprocessing method to read references
1419 Parameters
1420 ---------------
1421 tab: astropy Table of data
1422 Returns
1423 -----------
1424 stacked astropy Table of data
1425 """
1426 # distrib=np.unique(tab['z'])
1427 nlc = len(tab)
1428 print('ici pal', nlc)
1429 # n_multi=8
1430 if nlc >= 8:
1431 n_multi = min(nlc, 8)
1432 nvals = nlc/n_multi
1433 batch = range(0, nlc, nvals)
1434 batch = np.append(batch, nlc)
1435 else:
1436 batch = range(0, nlc)
1438 # lc_ref_tot={}
1439 # print('there pal',batch)
1440 result_queue = multiprocessing.Queue()
1441 for i in range(len(batch)-1):
1443 ida = int(batch[i])
1444 idb = int(batch[i+1])
1446 p = multiprocessing.Process(
1447 name='Subprocess_main-'+str(i), target=self.Read_Ref, args=(tab[ida:idb], i, result_queue))
1448 p.start()
1450 resultdict = {}
1451 for j in range(len(batch)-1):
1452 resultdict.update(result_queue.get())
1454 for p in multiprocessing.active_children():
1455 p.join()
1457 tab_res = Table()
1458 for j in range(len(batch)-1):
1459 if resultdict[j] is not None:
1460 tab_res = vstack([tab_res, resultdict[j]])
1462 return tab_res
1465class SN_Rate:
1466 """
1467 Estimate production rates of typeIa SN
1468 Available rates: Ripoche, Perrett, Dilday
1469 Parameters
1470 ---------------
1471 rate : str,opt
1472 type of rate chosen (Ripoche, Perrett, Dilday) (default : Perrett)
1473 H0 : float, opt
1474 Hubble constant value :math:`H_{0}`(default : 70.)
1475 Om0 : float, opt
1476 matter density value :math:`\Omega_{0}` (default : 0.25)
1477 min_rf_phase : float, opt
1478 min rest-frame phase (default : -15.)
1479 max_rf_phase : float, opt
1480 max rest-frame phase (default : 30.)
1481 """
1483 def __init__(self, rate='Perrett', H0=70, Om0=0.25,
1484 min_rf_phase=-15., max_rf_phase=30.):
1486 self.astropy_cosmo = FlatLambdaCDM(H0=H0, Om0=Om0)
1487 self.rate = rate
1488 self.min_rf_phase = min_rf_phase
1489 self.max_rf_phase = max_rf_phase
1491 def __call__(self, zmin=0.1, zmax=0.2,
1492 dz=0.01, survey_area=9.6,
1493 bins=None, account_for_edges=False,
1494 duration=140., duration_z=None):
1495 """
1496 call method
1497 Parameters
1498 ----------------
1499 zmin : float, opt
1500 minimal redshift (default : 0.1)
1501 zmax : float,opt
1502 max redshift (default : 0.2)
1503 dz : float, opt
1504 redshift bin (default : 0.001)
1505 survey_area : float, opt
1506 area of the survey (:math:`deg^{2}`) (default : 9.6 :math:`deg^{2}`)
1507 bins : list(float), opt
1508 redshift bins (default : None)
1509 account_for_edges : bool
1510 to account for season edges. If true, duration of the survey will be reduced by (1+z)*(maf_rf_phase-min_rf_phase)/365.25 (default : False)
1511 duration : float, opt
1512 survey duration (in days) (default : 140 days)
1513 duration_z : list(float),opt
1514 survey duration (as a function of z) (default : None)
1515 Returns
1516 -----------
1517 Lists :
1518 zz : float
1519 redshift values
1520 rate : float
1521 production rate
1522 err_rate : float
1523 production rate error
1524 nsn : float
1525 number of SN
1526 err_nsn : float
1527 error on the number of SN
1528 """
1530 if bins is None:
1531 thebins = np.arange(zmin, zmax+dz, dz)
1532 zz = 0.5 * (thebins[1:] + thebins[:-1])
1533 else:
1534 zz = bins
1535 thebins = bins
1537 rate, err_rate = self.SNRate(zz)
1538 error_rel = err_rate/rate
1540 area = survey_area / STERADIAN2SQDEG
1541 # or area= self.survey_area/41253.
1543 dvol = norm*self.astropy_cosmo.comoving_volume(thebins).value
1545 dvol = dvol[1:] - dvol[:-1]
1547 if account_for_edges:
1548 margin = (1.+zz) * (self.max_rf_phase-self.min_rf_phase) / 365.25
1549 effective_duration = duration / 365.25 - margin
1550 effective_duration[effective_duration <= 0.] = 0.
1551 else:
1552 # duration in days!
1553 effective_duration = duration/365.25
1554 if duration_z is not None:
1555 effective_duration = duration_z(zz)/365.25
1557 normz = (1.+zz)
1558 nsn = rate * area * dvol * effective_duration / normz
1559 err_nsn = err_rate*area * dvol * effective_duration / normz
1561 return zz, rate, err_rate, nsn, err_nsn
1563 def RipocheRate(self, z):
1564 """The SNLS SNIa rate according to the (unpublished) Ripoche et al study.
1565 Parameters
1566 --------------
1567 z : float
1568 redshift
1569 Returns
1570 ----------
1571 rate : float
1572 error_rate : float
1573 """
1574 rate = 1.53e-4*0.343
1575 expn = 2.14
1576 my_z = np.copy(z)
1577 my_z[my_z > 1.] = 1.
1578 rate_sn = rate * np.power((1+my_z)/1.5, expn)
1579 return rate_sn, 0.2*rate_sn
1581 def PerrettRate(self, z):
1582 """The SNLS SNIa rate according to (Perrett et al, 201?)
1583 Parameters
1584 --------------
1585 z : float
1586 redshift
1587 Returns
1588 ----------
1589 rate : float
1590 error_rate : float
1591 """
1592 rate = 0.17E-4
1593 expn = 2.11
1594 err_rate = 0.03E-4
1595 err_expn = 0.28
1596 my_z = np.copy(z)
1597 rate_sn = rate * np.power(1+my_z, expn)
1598 err_rate_sn = np.power(1+my_z, 2.*expn)*np.power(err_rate, 2.)
1599 err_rate_sn += np.power(rate_sn*np.log(1+my_z)*err_expn, 2.)
1601 return rate_sn, np.power(err_rate_sn, 0.5)
1603 def DildayRate(self, z):
1604 """The Dilday rate according to
1605 Parameters
1606 --------------
1607 z : float
1608 redshift
1609 Returns
1610 ----------
1611 rate : float
1612 error_rate : float
1613 """
1615 rate = 2.6e-5
1616 expn = 1.5
1617 err_rate = 0.01
1618 err_expn = 0.6
1619 my_z = np.copy(z)
1620 my_z[my_z > 1.] = 1.
1621 rate_sn = rate * np.power(1+my_z, expn)
1622 err_rate_sn = rate_sn*np.log(1+my_z)*err_expn
1623 return rate_sn, err_rate_sn
1625 """
1626 def flat_rate(self, z):
1627 return 1., 0.1
1628 """
1630 def SNRate(self, z):
1631 """SN rate estimation
1632 Parameters
1633 --------------
1634 z : float
1635 redshift
1636 Returns
1637 ----------
1638 rate : float
1639 error_rate : float
1640 """
1641 if self.rate == 'Ripoche':
1642 return self.RipocheRate(z)
1643 if self.rate == 'Perrett':
1644 return self.PerrettRate(z)
1645 if self.rate == 'Dilday':
1646 return self.DildayRate(z)
1648 def PlotNSN(self, zmin=0.1, zmax=0.2,
1649 dz=0.01, survey_area=9.6,
1650 bins=None, account_for_edges=False,
1651 duration=140., duration_z=None, norm=False):
1652 """ Plot integrated number of supernovae as a function of redshift
1653 uses the __call__ function
1654 Parameters
1655 --------------
1656 zmin : float, opt
1657 minimal redshift (default : 0.1)
1658 zmax : float,opt
1659 max redshift (default : 0.2)
1660 dz : float, opt
1661 redshift bin (default : 0.001)
1662 survey_area : float, opt
1663 area of the survey (:math:`deg^{2}`) (default : 9.6 :math:`deg^{2}`)
1664 bins : list(float), opt
1665 redshift bins (default : None)
1666 account_for_edges : bool
1667 to account for season edges. If true, duration of the survey will be reduced by (1+z)*(maf_rf_phase-min_rf_phase)/365.25 (default : False)
1668 duration : float, opt
1669 survey duration (in days) (default : 140 days)
1670 duration_z : list(float),opt
1671 survey duration (as a function of z) (default : None)
1672 norm: bool, opt
1673 to normalise the results (default: False)
1674 """
1675 import pylab as plt
1677 zz, rate, err_rate, nsn, err_nsn = self.__call__(
1678 zmin=zmin, zmax=zmax, dz=dz, bins=bins,
1679 account_for_edges=account_for_edges,
1680 duration=duration, survey_area=survey_area)
1682 nsn_sum = np.cumsum(nsn)
1684 if norm is False:
1685 plt.errorbar(zz, nsn_sum, yerr=np.sqrt(np.cumsum(err_nsn**2)))
1686 else:
1687 plt.errorbar(zz, nsn_sum/nsn_sum[-1])
1688 plt.xlabel('z')
1689 plt.ylabel('N$_{SN}$ <')
1690 plt.grid()
1693class CovColor:
1694 """
1695 class to estimate CovColor from lc using Fisher matrix element
1696 Parameters
1697 ---------------
1698 lc: pandas df
1699 lc to process. Should contain the Fisher matrix components
1700 ie the sum of the derivative of the fluxes wrt SN parameters
1701 """
1703 def __init__(self, lc):
1705 self.Cov_colorcolor = self.varColor(lc)
1707 def varColor(self, lc):
1708 """
1709 Method to estimate the variance color from matrix element
1710 Parameters
1711 --------------
1712 lc: pandas df
1713 data to process containing the derivative of the flux with respect to SN parameters
1714 Returns
1715 ----------
1716 float: Cov_colorcolor
1717 """
1718 a1 = lc['F_x0x0']
1719 a2 = lc['F_x0x1']
1720 a3 = lc['F_x0daymax']
1721 a4 = lc['F_x0color']
1723 b1 = a2
1724 b2 = lc['F_x1x1']
1725 b3 = lc['F_x1daymax']
1726 b4 = lc['F_x1color']
1728 c1 = a3
1729 c2 = b3
1730 c3 = lc['F_daymaxdaymax']
1731 c4 = lc['F_daymaxcolor']
1733 d1 = a4
1734 d2 = b4
1735 d3 = c4
1736 d4 = lc['F_colorcolor']
1738 detM = a1*self.det(b2, b3, b4, c2, c3, c4, d2, d3, d4)
1739 detM -= b1*self.det(a2, a3, a4, c2, c3, c4, d2, d3, d4)
1740 detM += c1*self.det(a2, a3, a4, b2, b3, b4, d2, d3, d4)
1741 detM -= d1*self.det(a2, a3, a4, b2, b3, b4, c2, c3, c4)
1743 res = -a3*b2*c1+a2*b3*c1+a3*b1*c2-a1*b3*c2-a2*b1*c3+a1*b2*c3
1745 return res/detM
1747 def det(self, a1, a2, a3, b1, b2, b3, c1, c2, c3):
1748 """
1749 Method to estimate the det of a matrix from its values
1750 Parameters
1751 ---------------
1752 Values of the matrix
1753 ( a1 a2 a3)
1754 (b1 b2 b3)
1755 (c1 c2 c3)
1756 Returns
1757 -----------
1758 det value
1759 """
1760 resp = a1*b2*c3+b1*c2*a3+c1*a2*b3
1761 resm = a3*b2*c1+b3*c2*a1+c3*a2*b1
1763 return resp-resm