Coverage for python/lsst/sims/maf/metrics/snNSNMetric.py : 6%

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
1import numpy as np
2from lsst.sims.maf.metrics import BaseMetric
3from lsst.sims.maf.utils.snNSNUtils import Load_Reference, Telescope, LCfast
4from lsst.sims.maf.utils.snNSNUtils import SN_Rate, CovColor
5import pandas as pd
6import numpy.lib.recfunctions as rf
7import time
8from scipy.interpolate import interp1d
9import numpy.lib.recfunctions as nlr
12class SNNSNMetric(BaseMetric):
13 """
14 Estimate (nSN,zlim) of type Ia supernovae.
16 Parameters
17 ---------------
18 metricName : str, opt
19 metric name (default : SNSNRMetric)
20 mjdCol : str, opt
21 mjd column name (default : observationStartMJD)
22 RACol : str,opt
23 Right Ascension column name (default : fieldRA)
24 DecCol : str,opt
25 Declinaison column name (default : fieldDec)
26 filterCol : str,opt
27 filter column name (default: filter)
28 m5Col : str, opt
29 five-sigma depth column name (default : fiveSigmaDepth)
30 exptimeCol : str,opt
31 exposure time column name (default : visitExposureTime)
32 nightCol : str,opt
33 night column name (default : night)
34 obsidCol : str,opt
35 observation id column name (default : observationId)
36 nexpCol : str,opt
37 number of exposure column name (default : numExposures)
38 vistimeCol : str,opt
39 visit time column name (default : visitTime)
40 season : list,opt
41 list of seasons to process (float)(default: -1 = all seasons)
42 zmin : float,opt
43 min redshift for the study (default: 0.0)
44 zmax : float,opt
45 max redshift for the study (default: 1.2)
46 pixArea: float, opt
47 pixel area (default: 9.6)
48 verbose: bool,opt
49 verbose mode (default: False)
50 ploteffi: bool,opt
51 to plot observing efficiencies vs z (default: False)
52 n_bef: int, opt
53 number of LC points LC before T0 (default:5)
54 n_aft: int, opt
55 number of LC points after T0 (default: 10)
56 snr_min: float, opt
57 minimal SNR of LC points (default: 5.0)
58 n_phase_min: int, opt
59 number of LC points with phase<= -5(default:1)
60 n_phase_max: int, opt
61 number of LC points with phase>= 20 (default: 1)
62 """
64 def __init__(self, metricName='SNNSNMetric',
65 mjdCol='observationStartMJD', RACol='fieldRA', DecCol='fieldDec',
66 filterCol='filter', m5Col='fiveSigmaDepth', exptimeCol='visitExposureTime',
67 nightCol='night', obsidCol='observationId', nexpCol='numExposures',
68 vistimeCol='visitTime', season=[-1], zmin=0.0, zmax=1.2,
69 pixArea=9.6, verbose=False, ploteffi=False,
70 n_bef=4, n_aft=10, snr_min=5., n_phase_min=1,
71 n_phase_max=1, templateDir=None, zlim_coeff=-1., **kwargs):
73 self.mjdCol = mjdCol
74 self.m5Col = m5Col
75 self.filterCol = filterCol
76 self.RACol = RACol
77 self.DecCol = DecCol
78 self.exptimeCol = exptimeCol
79 self.seasonCol = 'season'
80 self.nightCol = nightCol
81 self.obsidCol = obsidCol
82 self.nexpCol = nexpCol
83 self.vistimeCol = vistimeCol
84 self.pixArea = pixArea
85 self.zlim_coeff = zlim_coeff
87 cols = [self.nightCol, self.m5Col, self.filterCol, self.mjdCol, self.obsidCol,
88 self.nexpCol, self.vistimeCol, self.exptimeCol, self.seasonCol]
90 super(SNNSNMetric, self).__init__(
91 col=cols, metricDtype='object', metricName=metricName, **kwargs)
93 self.season = season
94 # LC selection parameters
95 self.n_bef = n_bef # nb points before peak
96 self.n_aft = n_aft # nb points after peak
97 self.snr_min = snr_min # SNR cut for points before/after peak
98 self.n_phase_min = n_phase_min # nb of point with phase <=-5
99 self.n_phase_max = n_phase_max # nb of points with phase >=20
101 # loading reference LC files
102 lc_reference = Load_Reference(templateDir=templateDir).ref
104 self.lcFast = {}
105 telescope = Telescope(airmass=1.2)
106 for key, vals in lc_reference.items():
107 self.lcFast[key] = LCfast(vals, key[0], key[1], telescope,
108 self.mjdCol, self.RACol, self.DecCol,
109 self.filterCol, self.exptimeCol,
110 self.m5Col, self.seasonCol, self.nexpCol,
111 self.snr_min)
113 # loading parameters
114 self.zmin = zmin # zmin for the study
115 self.zmax = zmax # zmax for the study
116 self.zStep = 0.05 # zstep
117 self.daymaxStep = 3. # daymax step
118 self.min_rf_phase = -20. # min ref phase for LC points selection
119 self.max_rf_phase = 40. # max ref phase for LC points selection
121 self.min_rf_phase_qual = -15. # min ref phase for bounds effects
122 self.max_rf_phase_qual = 25. # max ref phase for bounds effects
124 # snrate
125 self.rateSN = SN_Rate(
126 min_rf_phase=self.min_rf_phase_qual, max_rf_phase=self.max_rf_phase_qual)
128 # verbose mode - useful for debug and code performance estimation
129 self.verbose = False
130 self.ploteffi = False
132 # supernovae parameters
133 self.params = ['x0', 'x1', 'daymax', 'color']
135 # r = [(-1.0, -1.0)]
136 self.bad = np.rec.fromrecords([(-1.0, -1.0)], names=['nSN', 'zlim'])
137 # self.bad = {'nSN': -1.0, 'zlim': -1.0}
139 def run(self, dataSlice, slicePoint=None):
140 """
141 run method of the metric
143 Parameters
144 ---------------
145 dataSlice: array
146 data to process
148 """
149 idarray = None
150 healpixID = -1
151 if slicePoint is not None:
152 if 'nside' in slicePoint.keys():
153 import healpy as hp
154 self.pixArea = hp.nside2pixarea(
155 slicePoint['nside'], degrees=True)
156 r = []
157 names = []
159 healpixID = hp.ang2pix(
160 slicePoint['nside'], np.rad2deg(slicePoint['ra']), np.rad2deg(slicePoint['dec']), nest=True, lonlat=True)
161 for kk, vv in slicePoint.items():
162 r.append(vv)
163 names.append(kk)
164 idarray = np.rec.fromrecords([r], names=names)
165 else:
166 idarray = np.rec.fromrecords([0., 0.], names=['RA', 'Dec'])
168 # Two things to do: concatenate data (per band, night) and estimate seasons
169 dataSlice = rf.drop_fields(dataSlice, ['season'])
171 dataSlice = self.coadd(pd.DataFrame(dataSlice))
173 dataSlice = self.getseason(dataSlice)
175 # get the seasons
176 seasons = self.season
178 # if seasons = -1: process the seasons seen in data
179 if self.season == [-1]:
180 seasons = np.unique(dataSlice[self.seasonCol])
182 # get redshift range for processing
183 zRange = list(np.arange(self.zmin, self.zmax, self.zStep))
184 if zRange[0] < 1.e-6:
185 zRange[0] = 0.01
187 self.zRange = np.unique(zRange)
188 # season infos
189 dfa = pd.DataFrame(np.copy(dataSlice))
190 season_info = dfa.groupby(['season']).apply(
191 lambda x: self.seasonInfo(x)).reset_index()
193 #print('season info', season_info)
194 # select seasons of at least 30 days
195 idx = season_info['season_length'] >= 60.
196 season_info = season_info[idx]
198 # check wether requested seasons can be processed
199 test_season = season_info[season_info['season'].isin(seasons)]
200 # if len(test_season) == 0:
202 if test_season.empty:
203 return nlr.merge_arrays([idarray, self.bad], flatten=True)
204 else:
205 seasons = test_season['season']
206 # print('test_seas', seasons)
207 # print('hh', season_info)
209 # get season length depending on the redshift
210 dur_z = season_info.groupby(['season']).apply(
211 lambda x: self.duration_z(x)).reset_index()
213 # remove dur_z with negative season lengths
214 idx = dur_z['season_length'] >= 10.
215 dur_z = dur_z[idx]
217 # generating simulation parameters
218 gen_par = dur_z.groupby(['z', 'season']).apply(
219 lambda x: self.calcDaymax(x)).reset_index()
221 if gen_par.empty:
222 return nlr.merge_arrays([idarray, self.bad], flatten=True)
223 resdf = pd.DataFrame()
225 for seas in seasons:
226 vara_df = self.run_season(
227 dataSlice, [seas], gen_par, dur_z)
228 #print('res', seas, vara_df)
229 if vara_df is not None:
230 resdf = pd.concat((resdf, vara_df))
232 # final result: median zlim for a faint sn
233 # and nsn_med for z<zlim
235 if resdf.empty:
236 return nlr.merge_arrays([idarray, self.bad], flatten=True)
238 resdf = resdf.round({'zlim': 3, 'nsn_med': 3})
239 x1_ref = -2.0
240 color_ref = 0.2
242 idx = np.abs(resdf['x1']-x1_ref) < 1.e-5
243 idx &= np.abs(resdf['color']-color_ref) < 1.e-5
244 idx &= resdf['zlim'] > 0
246 if not resdf[idx].empty:
247 zlim = resdf[idx]['zlim'].median()
248 nSN = resdf[idx]['nsn_med'].sum()
250 resd = np.rec.fromrecords([(nSN, zlim, healpixID)], names=[
251 'nSN', 'zlim', 'healpixID'])
253 res = nlr.merge_arrays([idarray, resd], flatten=True)
255 else:
257 res = nlr.merge_arrays([idarray, self.bad], flatten=True)
259 return res
261 def reducenSN(self, metricVal):
263 # At each slicepoint, return the sum nSN value.
265 return np.sum(metricVal['nSN'])
267 def reducezlim(self, metricVal):
269 # At each slicepoint, return the median zlim
271 return np.median(metricVal['zlim'])
273 def coadd(self, data):
274 """
275 Method to coadd data per band and per night
277 Parameters
278 ---------------
279 data: pandas df of observations
281 Returns
282 -----------
283 coadded data (pandas df)
285 """
287 keygroup = [self.filterCol, self.nightCol]
289 data.sort_values(by=keygroup, ascending=[
290 True, True], inplace=True)
292 coadd_df = data.groupby(keygroup).agg({self.nexpCol: ['sum'],
293 self.vistimeCol: ['sum'],
294 self.exptimeCol: ['sum'],
295 self.mjdCol: ['mean'],
296 self.RACol: ['mean'],
297 self.DecCol: ['mean'],
298 self.m5Col: ['mean']}).reset_index()
300 coadd_df.columns = [self.filterCol, self.nightCol, self.nexpCol,
301 self.vistimeCol, self.exptimeCol, self.mjdCol,
302 self.RACol, self.DecCol, self.m5Col]
304 coadd_df.loc[:, self.m5Col] += 1.25 * \
305 np.log10(coadd_df[self.vistimeCol]/30.)
307 coadd_df.sort_values(by=[self.filterCol, self.nightCol], ascending=[
308 True, True], inplace=True)
310 return coadd_df.to_records(index=False)
312 def getseason(self, obs, season_gap=80., mjdCol='observationStartMJD'):
313 """
314 Method to estimate seasons
315 Parameters
316 --------------
317 obs: numpy array
318 array of observations
319 season_gap: float, opt
320 minimal gap required to define a season (default: 80 days)
321 mjdCol: str, opt
322 col name for MJD infos (default: observationStartMJD)
324 Returns
325 ----------
326 original numpy array with seasonnumber appended
327 """
329 # check wether season has already been estimated
331 if 'season' in obs.dtype.names:
332 return obs
334 obs.sort(order=mjdCol)
336 seasoncalc = np.ones(obs.size, dtype=int)
338 if len(obs) > 1:
339 diff = np.diff(obs[mjdCol])
340 flag = np.where(diff > season_gap)[0]
342 if len(flag) > 0:
343 for i, indx in enumerate(flag):
344 seasoncalc[indx+1:] = i+2
346 obs = rf.append_fields(obs, 'season', seasoncalc)
348 return obs
350 def seasonInfo(self, grp):
351 """
352 Method to estimate seasonal info (cadence, season length, ...)
353 Parameters
354 --------------
355 grp: pandas df group
356 Returns
357 ---------
358 pandas df with the cfollowing cols:
359 """
360 df = pd.DataFrame([len(grp)], columns=['Nvisits'])
361 df['MJD_min'] = grp[self.mjdCol].min()
362 df['MJD_max'] = grp[self.mjdCol].max()
363 df['season_length'] = df['MJD_max']-df['MJD_min']
364 df['cadence'] = 0.
366 if len(grp) > 5:
367 to = grp.groupby(['night'])[self.mjdCol].median().sort_values()
368 df['cadence'] = np.mean(to.diff())
370 return df
372 def duration_z(self, grp):
373 """
374 Method to estimate the season length vs redshift
375 This is necessary to take into account boundary effects
376 when estimating the number of SN that can be detected
377 daymin, daymax = min and max MJD of a season
378 T0_min(z) = daymin-(1+z)*min_rf_phase_qual
379 T0_max(z) = daymax-(1+z)*max_rf_phase_qual
380 season_length(z) = T0_max(z)-T0_min(z)
381 Parameters
382 --------------
383 grp: pandas df group
384 data to process: season infos
385 Returns
386 ----------
387 pandas df with season_length, z, T0_min and T0_max cols
388 """
390 daymin = grp['MJD_min'].values
391 daymax = grp['MJD_max'].values
392 dur_z = pd.DataFrame(self.zRange, columns=['z'])
393 dur_z['T0_min'] = daymin-(1.+dur_z['z'])*self.min_rf_phase_qual
394 dur_z['T0_max'] = daymax-(1.+dur_z['z'])*self.max_rf_phase_qual
395 dur_z['season_length'] = dur_z['T0_max']-dur_z['T0_min']
397 return dur_z
399 def calcDaymax(self, grp):
400 """
401 Method to estimate T0 (daymax) values for simulation.
402 Parameters
403 --------------
404 grp: group (pandas df sense)
405 group of data to process with the following cols:
406 T0_min: T0 min value (per season)
407 T0_max: T0 max value (per season)
408 Returns
409 ----------
410 pandas df with daymax, min_rf_phase, max_rf_phase values
411 """
413 T0_max = grp['T0_max'].values
414 T0_min = grp['T0_min'].values
415 num = (T0_max-T0_min)/self.daymaxStep
416 if T0_max-T0_min > 10:
417 df = pd.DataFrame(np.linspace(
418 T0_min, T0_max, int(num)), columns=['daymax'])
419 else:
420 df = pd.DataFrame([-1], columns=['daymax'])
422 df['min_rf_phase'] = self.min_rf_phase_qual
423 df['max_rf_phase'] = self.max_rf_phase_qual
425 return df
427 def run_season(self, dataSlice, season, gen_par, dura_z):
428 """
429 Method to run on seasons
430 Parameters
431 --------------
432 dataSlice: numpy array, opt
433 data to process (scheduler simulations)
434 seasons: list(int)
435 list of seasons to process
436 Returns
437 ---------
438 effi_seasondf: pandas df
439 efficiency curves
440 zlimsdf: pandas df
441 redshift limits and number of supernovae
442 """
444 time_ref = time.time()
446 if self.verbose:
447 print('#### Processing season', season)
449 groupnames = ['season', 'x1', 'color']
451 gen_p = gen_par[gen_par['season'].isin(season)]
453 if gen_p.empty:
454 if self.verbose:
455 print('No generator parameter found')
456 return None
457 dur_z = dura_z[dura_z['season'].isin(season)]
458 obs = pd.DataFrame(np.copy(dataSlice))
459 obs = obs[obs['season'].isin(season)]
461 obs = obs.sort_values(by=['night'])
462 #print('data here', obs.columns)
463 #print(obs[['night', 'filter', 'observationStartMJD', 'fieldRA', 'fieldDec']])
464 """
465 import matplotlib.pyplot as plt
466 plt.plot(dataSlice['fieldRA'], dataSlice['fieldDec'], 'ko')
467 print('data', len(dataSlice))
468 plt.show()
469 """
470 # simulate supernovae and lc
471 if self.verbose:
472 print("SN generation")
473 print(season, obs)
474 sn = self.genSN(obs.to_records(
475 index=False), gen_p.to_records(index=False))
477 if self.verbose:
478 idx = np.abs(sn['x1']+2) < 1.e-5
479 idx &= np.abs(sn['z']-0.2) < 1.e-5
480 sel = sn[idx]
481 sel = sel.sort_values(by=['z', 'daymax'])
483 print('sn and lc', len(sn), sel.columns,
484 sel[['x1', 'color', 'z', 'daymax', 'Cov_colorcolor', 'n_bef', 'n_aft']])
486 # from these supernovae: estimate observation efficiency vs z
487 effi_seasondf = self.effidf(sn)
489 # zlims can only be estimated if efficiencies are ok
490 idx = effi_seasondf['z'] <= 0.2
491 x1ref = -2.0
492 colorref = 0.2
493 idx &= np.abs(effi_seasondf['x1']-x1ref) < 1.e-5
494 idx &= np.abs(effi_seasondf['color']-colorref) < 1.e-5
495 sel = effi_seasondf[idx]
497 if np.mean(sel['effi']) > 0.02:
498 # estimate zlims
499 zlimsdf = self.zlims(effi_seasondf, dur_z, groupnames)
501 # estimate number of medium supernovae
502 zlimsdf['nsn_med'], zlimsdf['var_nsn_med'] = zlimsdf.apply(lambda x: self.nsn_typedf(
503 x, 0.0, 0.0, effi_seasondf, dur_z), axis=1, result_type='expand').T.values
504 else:
505 return None
507 if self.verbose:
508 print('#### SEASON processed', time.time()-time_ref,
509 season)
511 return zlimsdf
513 def genSN(self, obs, gen_par):
514 """
515 Method to simulate LC and supernovae
516 Parameters
517 ---------------
518 obs: numpy array
519 array of observations(from scheduler)
520 gen_par: numpy array
521 array of parameters for simulation
522 """
524 time_ref = time.time()
525 # LC estimation
527 sn_tot = pd.DataFrame()
528 lc_tot = pd.DataFrame()
529 for key, vals in self.lcFast.items():
530 time_refs = time.time()
531 gen_par_cp = np.copy(gen_par)
532 if key == (-2.0, 0.2):
533 idx = gen_par_cp['z'] < 0.9
534 gen_par_cp = gen_par_cp[idx]
535 lc = vals(obs, gen_par_cp, bands='grizy')
536 if self.verbose:
537 print('End of simulation', key, time.time()-time_refs)
539 if self.verbose:
540 print('End of simulation after concat',
541 key, time.time()-time_refs)
543 # estimate SN
545 sn = pd.DataFrame()
546 if len(lc) > 0:
547 sn = self.process(pd.DataFrame(lc))
549 if self.verbose:
550 print('End of supernova', time.time()-time_refs)
552 if not sn.empty:
553 sn_tot = pd.concat([sn_tot, pd.DataFrame(sn)], sort=False)
555 if self.verbose:
556 print('End of supernova - all', time.time()-time_ref)
558 return sn_tot
560 def process(self, tab):
561 """
562 Method to process LC: sigma_color estimation and LC selection
563 Parameters
564 --------------
565 tab: pandas df of LC points with the following cols:
566 flux: flux
567 fluxerr: flux error
568 phase: phase
569 snr_m5: Signal-to-Noise Ratio
570 time: time(MJD)
571 mag: magnitude
572 m5: five-sigma depth
573 magerr: magnitude error
574 exposuretime: exposure time
575 band: filter
576 zp: zero-point
577 season: season number
578 healpixID: pixel ID
579 pixRA: pixel RA
580 pixDec: pixel Dec
581 z: redshift
582 daymax: T0
583 flux_e_sec: flux(in photoelec/sec)
584 flux_5: 5-sigma flux(in photoelec/sec)
585 F_x0x0, ...F_colorcolor: Fisher matrix elements
586 x1: x1 SN
587 color: color SN
588 n_aft: number of LC points before daymax
589 n_bef: number of LC points after daymax
590 n_phmin: number of LC points with a phase < -5
591 n_phmax: number of LC points with a phase > 20
592 Returns
593 ----------
594 """
595 # now groupby
596 tab = tab.round({'daymax': 3,
597 'z': 3, 'x1': 2, 'color': 2})
598 groups = tab.groupby(
599 ['daymax', 'season', 'z', 'x1', 'color'])
601 tosum = []
602 for ia, vala in enumerate(self.params):
603 for jb, valb in enumerate(self.params):
604 if jb >= ia:
605 tosum.append('F_'+vala+valb)
606 tosum += ['n_aft', 'n_bef', 'n_phmin', 'n_phmax']
607 # apply the sum on the group
608 sums = groups[tosum].sum().reset_index()
610 # select LC according to the number of points bef/aft peak
611 idx = sums['n_aft'] >= self.n_aft
612 idx &= sums['n_bef'] >= self.n_bef
613 idx &= sums['n_phmin'] >= self.n_phase_min
614 idx &= sums['n_phmax'] >= self.n_phase_max
616 if self.verbose:
617 print('selection parameters', self.n_bef,
618 self.n_aft, self.n_phase_min, self.n_phase_max)
619 finalsn = pd.DataFrame()
620 goodsn = pd.DataFrame(sums.loc[idx])
622 # estimate the color for SN that passed the selection cuts
623 if len(goodsn) > 0:
624 goodsn.loc[:, 'Cov_colorcolor'] = CovColor(goodsn).Cov_colorcolor
625 finalsn = pd.concat([finalsn, goodsn], sort=False)
627 badsn = pd.DataFrame(sums.loc[~idx])
629 # Supernovae that did not pass the cut have a sigma_color=10
630 if len(badsn) > 0:
631 badsn.loc[:, 'Cov_colorcolor'] = 100.
632 finalsn = pd.concat([finalsn, badsn], sort=False)
634 return finalsn
636 def effidf(self, sn_tot, color_cut=0.04):
637 """
638 Method estimating efficiency vs z for a sigma_color cut
639 Parameters
640 ---------------
641 sn_tot: pandas df
642 data used to estimate efficiencies
643 color_cut: float, opt
644 color selection cut(default: 0.04)
645 Returns
646 ----------
647 effi: pandas df with the following cols:
648 season: season
649 pixRA: RA of the pixel
650 pixDec: Dec of the pixel
651 healpixID: pixel ID
652 x1: SN stretch
653 color: SN color
654 z: redshift
655 effi: efficiency
656 effi_err: efficiency error(binomial)
657 """
659 sndf = pd.DataFrame(sn_tot)
661 listNames = ['season', 'x1', 'color']
662 groups = sndf.groupby(listNames)
664 # estimating efficiencies
665 effi = groups[['Cov_colorcolor', 'z']].apply(
666 lambda x: self.effiObsdf(x, color_cut)).reset_index(level=list(range(len(listNames))))
668 # this is to plot efficiencies and also sigma_color vs z
669 if self.ploteffi:
670 import matplotlib.pylab as plt
671 fig, ax = plt.subplots()
672 figb, axb = plt.subplots()
674 self.plot(ax, effi, 'effi', 'effi_err',
675 'Observing Efficiencies', ls='-')
676 sndf['sigma_color'] = np.sqrt(sndf['Cov_colorcolor'])
677 self.plot(axb, sndf, 'sigma_color', None, '$\sigma_{color}$')
678 # get efficiencies vs z
680 plt.show()
682 return effi
684 def plot(self, ax, effi, vary, erry=None, legy='', ls='None'):
685 """
686 Simple method to plot vs z
687 Parameters
688 --------------
689 ax:
690 axis where to plot
691 effi: pandas df
692 data to plot
693 vary: str
694 variable(column of effi) to plot
695 erry: str, opt
696 error on y-axis(default: None)
697 legy: str, opt
698 y-axis legend(default: '')
699 """
700 grb = effi.groupby(['x1', 'color'])
701 yerr = None
702 for key, grp in grb:
703 x1 = grp['x1'].unique()[0]
704 color = grp['color'].unique()[0]
705 if erry is not None:
706 yerr = grp[erry]
707 ax.errorbar(grp['z'], grp[vary], yerr=yerr,
708 marker='o', label='(x1,color)=({},{})'.format(x1, color), lineStyle=ls)
710 ftsize = 15
711 ax.set_xlabel('z', fontsize=ftsize)
712 ax.set_ylabel(legy, fontsize=ftsize)
713 ax.xaxis.set_tick_params(labelsize=ftsize)
714 ax.yaxis.set_tick_params(labelsize=ftsize)
715 ax.legend(fontsize=ftsize)
717 def effiObsdf(self, data, color_cut=0.04):
718 """
719 Method to estimate observing efficiencies for supernovae
720 Parameters
721 --------------
722 data: pandas df - grp
723 data to process
724 Returns
725 ----------
726 pandas df with the following cols:
727 - cols used to make the group
728 - effi, effi_err: observing efficiency and associated error
729 """
731 # reference df to estimate efficiencies
732 df = data.loc[lambda dfa: np.sqrt(dfa['Cov_colorcolor']) < 100000., :]
734 # selection on sigma_c<= 0.04
735 df_sel = df.loc[lambda dfa: np.sqrt(
736 dfa['Cov_colorcolor']) <= color_cut, :]
738 # make groups (with z)
739 group = df.groupby('z')
740 group_sel = df_sel.groupby('z')
742 # Take the ratio to get efficiencies
743 rb = (group_sel.size()/group.size())
744 err = np.sqrt(rb*(1.-rb)/group.size())
745 var = rb*(1.-rb)*group.size()
747 rb = rb.array
748 err = err.array
749 var = var.array
751 rb[np.isnan(rb)] = 0.
752 err[np.isnan(err)] = 0.
753 var[np.isnan(var)] = 0.
755 return pd.DataFrame({group.keys: list(group.groups.keys()),
756 'effi': rb,
757 'effi_err': err,
758 'effi_var': var})
760 def zlims(self, effi_seasondf, dur_z, groupnames):
761 """
762 Method to estimate redshift limits
763 Parameters
764 --------------
765 effi_seasondf: pandas df
766 season: season
767 pixRA: RA of the pixel
768 pixDec: Dec of the pixel
769 healpixID: pixel ID
770 x1: SN stretch
771 color: SN color
772 z: redshift
773 effi: efficiency
774 effi_err: efficiency error (binomial)
775 dur_z: pandas df with the following cols:
776 season: season
777 z: redshift
778 T0_min: min daymax
779 T0_max: max daymax
780 season_length: season length
781 groupnames: list(str)
782 list of columns to use to define the groups
783 Returns
784 ----------
785 pandas df with the following cols: pixRA: RA of the pixel
786 pixDec: Dec of the pixel
787 healpixID: pixel ID
788 season: season number
789 x1: SN stretch
790 color: SN color
791 zlim: redshift limit
792 """
794 res = effi_seasondf.groupby(groupnames).apply(
795 lambda x: self.zlimdf(x, dur_z)).reset_index(level=list(range(len(groupnames))))
797 return res
799 def zlimdf(self, grp, duration_z):
800 """
801 Method to estimate redshift limits
802 Parameters
803 --------------
804 grp: pandas df group
805 efficiencies to estimate redshift limits;
806 columns:
807 season: season
808 pixRA: RA of the pixel
809 pixDec: Dec of the pixel
810 healpixID: pixel ID
811 x1: SN stretch
812 color: SN color
813 z: redshift
814 effi: efficiency
815 effi_err: efficiency error (binomial)
816 duration_z: pandas df with the following cols:
817 season: season
818 z: redshift
819 T0_min: min daymax
820 T0_max: max daymax
821 season_length: season length
822 Returns
823 ----------
824 pandas df with the following cols:
825 zlimit: redshift limit
826 """
828 zlimit = 0.0
830 # z range for the study
831 zplot = np.arange(self.zmin, self.zmax, 0.01)
833 # print(grp['z'], grp['effi'])
835 if len(grp['z']) <= 3:
836 return pd.DataFrame({'zlim': [zlimit]})
837 # 'status': [int(status)]})
838 # interpolate efficiencies vs z
839 effiInterp = interp1d(
840 grp['z'], grp['effi'], kind='linear', bounds_error=False, fill_value=0.)
842 if self.zlim_coeff < 0.:
843 # in that case zlim is estimated from efficiencies
844 # first step: identify redshift domain with efficiency decrease
845 zlimit = self.zlim_from_effi(effiInterp, zplot)
846 #status = self.status['ok']
848 else:
849 zlimit = self.zlim_from_cumul(
850 grp, duration_z, effiInterp, zplot)
852 return pd.DataFrame({'zlim': [zlimit]})
853 # 'status': [int(status)]})
855 def zlim_from_cumul(self, grp, duration_z, effiInterp, zplot, rate='cte'):
856 """
857 Method to estimate the redshift limit from the cumulative
858 The redshift limit is estimated to be the z value corresponding to:
859 frac(NSN(z<zlimit))=zlimi_coeff
861 Parameters
862 ---------------
863 grp: pandas group
864 data to process
865 duration_z: array
866 duration as a function of the redshift
867 effiInterp: interp1d
868 interpolator for efficiencies
869 zplot: interp1d
870 interpolator for redshift values
871 rate: str, opt
872 rate to estimate the number of SN to estimate zlimit
873 rate = cte: rate independent of z
874 rate = SN_rate: rate from SN_Rate class
876 Returns
877 ----------
878 zlimit: float
879 the redshift limit
880 """
882 if rate == 'SN_rate':
883 # get rate
884 season = np.median(grp['season'])
885 idx = duration_z['season'] == season
886 seas_duration_z = duration_z[idx]
888 durinterp_z = interp1d(
889 seas_duration_z['z'], seas_duration_z['season_length'], bounds_error=False, fill_value=0.)
891 # estimate the rates and nsn vs z
892 zz, rate, err_rate, nsn, err_nsn = self.rateSN(zmin=self.zmin,
893 zmax=self.zmax,
894 duration_z=durinterp_z,
895 survey_area=self.pixArea)
897 # rate interpolation
898 rateInterp = interp1d(zz, nsn, kind='linear',
899 bounds_error=False, fill_value=0)
900 else:
901 # this is for a rate z-independent
902 nsn = np.ones(len(zplot))
903 rateInterp = interp1d(zplot, nsn, kind='linear',
904 bounds_error=False, fill_value=0)
906 nsn_cum = np.cumsum(effiInterp(zplot)*rateInterp(zplot))
908 if nsn_cum[-1] >= 1.e-5:
909 nsn_cum_norm = nsn_cum/nsn_cum[-1] # normalize
910 zlim = interp1d(nsn_cum_norm, zplot)
911 zlimit = zlim(self.zlim_coeff).item()
913 if self.ploteffi:
914 self.plot_NSN_cumul(grp, nsn_cum_norm, zplot)
915 else:
916 zlimit = 0.
918 return zlimit
920 def plot_NSN_cumul(self, grp, nsn_cum_norm, zplot):
921 """
922 Method to plot the NSN cumulative vs redshift
923 Parameters
924 --------------
925 grp: pandas group
926 data to process
927 """
929 import matplotlib.pylab as plt
930 fig, ax = plt.subplots()
931 x1 = grp['x1'].unique()[0]
932 color = grp['color'].unique()[0]
934 ax.plot(zplot, nsn_cum_norm,
935 label='(x1,color)=({},{})'.format(x1, color))
937 ftsize = 15
938 ax.set_ylabel('NSN ($z<$)', fontsize=ftsize)
939 ax.set_xlabel('z', fontsize=ftsize)
940 ax.xaxis.set_tick_params(labelsize=ftsize)
941 ax.yaxis.set_tick_params(labelsize=ftsize)
942 ax.set_xlim((0.0, 0.8))
943 ax.set_ylim((0.0, 1.05))
944 ax.plot([0., 1.2], [0.95, 0.95], ls='--', color='k')
945 plt.legend(fontsize=ftsize)
946 plt.show()
948 def zlim_from_effi(self, effiInterp, zplot):
949 """
950 Method to estimate the redshift limit from efficiency curves
951 The redshift limit is defined here as the redshift value beyond
952 which efficiency decreases up to zero.
953 Parameters
954 ---------------
955 effiInterp: interpolator
956 use to get efficiencies
957 zplot: numpy array
958 redshift values
959 Returns
960 -----------
961 zlimit: float
962 the redshift limit
963 """
965 # get efficiencies
966 effis = effiInterp(zplot)
967 # select data with efficiency decrease
968 idx = np.where(np.diff(effis) < -0.005)
969 z_effi = np.array(zplot[idx], dtype={
970 'names': ['z'], 'formats': [np.float]})
971 # from this make some "z-periods" to avoid accidental zdecrease at low z
972 z_gap = 0.05
973 seasoncalc = np.ones(z_effi.size, dtype=int)
974 diffz = np.diff(z_effi['z'])
975 flag = np.where(diffz > z_gap)[0]
977 if len(flag) > 0:
978 for i, indx in enumerate(flag):
979 seasoncalc[indx+1:] = i+2
980 z_effi = rf.append_fields(z_effi, 'season', seasoncalc)
982 # now take the highest season (end of the efficiency curve)
983 idd = z_effi['season'] == np.max(z_effi['season'])
984 zlimit = np.min(z_effi[idd]['z'])
986 return zlimit
988 def zlimdf_deprecated(self, grp, duration_z):
989 """
990 Method to estimate redshift limits
991 Parameters
992 --------------
993 grp: pandas df group
994 efficiencies to estimate redshift limits;
995 columns:
996 season: season
997 pixRA: RA of the pixel
998 pixDec: Dec of the pixel
999 healpixID: pixel ID
1000 x1: SN stretch
1001 color: SN color
1002 z: redshift
1003 effi: efficiency
1004 effi_err: efficiency error (binomial)
1005 duration_z: pandas df with the following cols:
1006 season: season
1007 z: redshift
1008 T0_min: min daymax
1009 T0_max: max daymax
1010 season_length: season length
1011 Returns
1012 ----------
1013 pandas df with the following cols:
1014 zlimit: redshift limit
1015 """
1017 zlimit = 0.0
1019 # z range for the study
1020 zplot = list(np.arange(self.zmin, self.zmax, 0.001))
1022 # print(grp['z'], grp['effi'])
1024 if len(grp['z']) <= 3:
1025 return pd.DataFrame({'zlim': [zlimit]})
1026 # 'status': [int(status)]})
1027 # interpolate efficiencies vs z
1028 effiInterp = interp1d(
1029 grp['z'], grp['effi'], kind='linear', bounds_error=False, fill_value=0.)
1031 # get rate
1032 season = np.median(grp['season'])
1033 idx = duration_z['season'] == season
1034 seas_duration_z = duration_z[idx]
1036 # print('hhh1', seas_duration_z)
1037 durinterp_z = interp1d(
1038 seas_duration_z['z'], seas_duration_z['season_length'], bounds_error=False, fill_value=0.)
1040 # estimate the rates and nsn vs z
1041 zz, rate, err_rate, nsn, err_nsn = self.rateSN(zmin=self.zmin,
1042 zmax=self.zmax,
1043 duration_z=durinterp_z,
1044 survey_area=self.pixArea)
1046 # rate interpolation
1047 rateInterp = interp1d(zz, nsn, kind='linear',
1048 bounds_error=False, fill_value=0)
1050 # estimate the cumulated number of SN vs z
1051 nsn_cum = np.cumsum(effiInterp(zplot)*rateInterp(zplot))
1053 if nsn_cum[-1] >= 1.e-5:
1054 nsn_cum_norm = nsn_cum/nsn_cum[-1] # normalize
1055 zlim = interp1d(nsn_cum_norm, zplot)
1056 zlimit = zlim(0.95).item()
1057 # status = self.status['ok']
1059 if self.ploteffi:
1060 import matplotlib.pylab as plt
1061 fig, ax = plt.subplots()
1062 x1 = grp['x1'].unique()[0]
1063 color = grp['color'].unique()[0]
1065 ax.plot(zplot, nsn_cum_norm,
1066 label='(x1,color)=({},{})'.format(x1, color))
1068 ftsize = 15
1069 ax.set_ylabel('NSN ($z<$)', fontsize=ftsize)
1070 ax.set_xlabel('z', fontsize=ftsize)
1071 ax.xaxis.set_tick_params(labelsize=ftsize)
1072 ax.yaxis.set_tick_params(labelsize=ftsize)
1073 ax.set_xlim((0.0, 1.2))
1074 ax.set_ylim((0.0, 1.05))
1075 ax.plot([0., 1.2], [0.95, 0.95], ls='--', color='k')
1076 plt.legend(fontsize=ftsize)
1077 plt.show()
1079 return pd.DataFrame({'zlim': [zlimit]})
1080 # 'status': [int(status)]})
1082 def nsn_typedf(self, grp, x1, color, effi_tot, duration_z, search=True):
1083 """
1084 Method to estimate the number of supernovae for a given type of SN
1085 Parameters
1086 --------------
1087 grp: pandas series with the following infos:
1088 pixRA: pixelRA
1089 pixDec: pixel Dec
1090 healpixID: pixel ID
1091 season: season
1092 x1: SN stretch
1093 color: SN color
1094 zlim: redshift limit
1095 x1, color: SN params to estimate the number
1096 effi_tot: pandas df with columns:
1097 season: season
1098 pixRA: RA of the pixel
1099 pixDec: Dec of the pixel
1100 healpixID: pixel ID
1101 x1: SN stretch
1102 color: SN color
1103 z: redshift
1104 effi: efficiency
1105 effi_err: efficiency error (binomial)
1106 duration_z: pandas df with the following cols:
1107 season: season
1108 z: redshift
1109 T0_min: min daymax
1110 T0_max: max daymax
1111 season_length: season length
1112 Returns
1113 ----------
1114 nsn: float
1115 number of supernovae
1116 """
1118 # get rate
1119 season = np.median(grp['season'])
1120 idx = duration_z['season'] == season
1121 seas_duration_z = duration_z[idx]
1123 # print('hhh2', seas_duration_z)
1124 durinterp_z = interp1d(
1125 seas_duration_z['z'], seas_duration_z['season_length'], bounds_error=False, fill_value=0.)
1127 if search:
1128 effisel = effi_tot.loc[lambda dfa: (
1129 dfa['x1'] == x1) & (dfa['color'] == color), :]
1130 else:
1131 effisel = effi_tot
1133 nsn, var_nsn = self.nsn(effisel, grp['zlim'], durinterp_z)
1135 return (nsn, var_nsn)
1137 def nsn(self, effi, zlim, duration_z):
1138 """
1139 Method to estimate the number of supernovae
1140 Parameters
1141 --------------
1142 effi: pandas df grp of efficiencies
1143 season: season
1144 pixRA: RA of the pixel
1145 pixDec: Dec of the pixel
1146 healpixID: pixel ID
1147 x1: SN stretch
1148 color: SN color
1149 z: redshift
1150 effi: efficiency
1151 effi_err: efficiency error (binomial)
1152 zlim: float
1153 redshift limit value
1154 duration_z: pandas df with the following cols:
1155 season: season
1156 z: redshift
1157 T0_min: min daymax
1158 T0_max: max daymax
1159 season_length: season length
1160 Returns
1161 ----------
1162 nsn, var_nsn : float
1163 number of supernovae (and variance) with z<zlim
1164 """
1166 if zlim < 1.e-3:
1167 return (-1.0, -1.0)
1169 dz = 0.001
1170 zplot = list(np.arange(self.zmin, self.zmax, dz))
1171 # interpolate efficiencies vs z
1172 effiInterp = interp1d(
1173 effi['z'], effi['effi'], kind='linear', bounds_error=False, fill_value=0.)
1174 # estimate the cumulated number of SN vs z
1175 zz, rate, err_rate, nsn, err_nsn = self.rateSN(zmin=self.zmin,
1176 zmax=self.zmax,
1177 dz=dz,
1178 duration_z=duration_z,
1179 survey_area=self.pixArea)
1181 nsn_cum = np.cumsum(effiInterp(zplot)*nsn)
1183 nsn_interp = interp1d(zplot, nsn_cum)
1185 nsn = nsn_interp(zlim).item()
1187 return [nsn, 0.0]