Hide keyboard shortcuts

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

1""" 

2Mixins for the InstanceCatalog class to provide SN catalogs in catsim. There 

3are three classes here: 

4 - SNFunctionality: provides common functions required by all SN instance 

5 catalogs. It does not make sense to instantiate this class, but rather 

6 it should be used as a mixin alongside another class. 

7 - SNIaCatalog: Dynamically created catalogs by sampling user specified 

8 distributions of SN parameters on the fly based on host galaxies in the 

9 catsim database. 

10 - FrozenSNCat: catalogs that are 'frozen' on the catsim database. For a 

11 user to use one of these catalogs, such a catalog would have to be 

12 uploaded to catsim. Examples of such catalogs that are in the catsim 

13 database are the tables `TwinkSN` and `TwinkSNKraken` 

14""" 

15from builtins import str 

16from builtins import range 

17import numpy as np 

18import numbers 

19 

20from lsst.sims.catalogs.definitions import InstanceCatalog 

21from lsst.sims.catalogs.decorators import compound 

22from lsst.sims.photUtils import (BandpassDict, Bandpass) 

23from lsst.sims.catUtils.mixins import CosmologyMixin 

24import lsst.sims.photUtils.PhotometricParameters as PhotometricParameters 

25from lsst.sims.catUtils.supernovae import SNObject 

26from lsst.sims.catUtils.supernovae import SNUniverse 

27from lsst.sims.catUtils.mixins import EBVmixin 

28from lsst.sims.utils import _galacticFromEquatorial 

29import astropy 

30 

31 

32__all__ = ['SNIaCatalog', 'SNFunctionality', 'FrozenSNCat'] 

33cosmo = CosmologyMixin() 

34 

35 

36class SNFunctionality(InstanceCatalog, EBVmixin, CosmologyMixin, SNUniverse): 

37 """ 

38 SNFunctionality is a mixin that provides functionality of getting fluxes 

39 and magnitudes for SN defined by parameters of `~sims_catUtils.SNObject` as 

40 defined in `~sims_catUtils/python/lsst/sims/catUtils/supernovae/SNObject` 

41 

42 

43 This class is not meant to be used by itself, as it does not have any way 

44 of obtaining its attributes, but as a mixin to classes like SNIaCatalog 

45 which define these attributes. 

46 """ 

47 

48 # Write the location of SED file (for example for PhoSim) 

49 writeSedFile = False 

50 # prefix to use for SED File name 

51 sn_sedfile_prefix = '' 

52 

53 # t_0, c, x_1, x_0 are parameters characterizing a SALT 

54 # based SN model as defined in sncosmo 

55 column_outputs = ['snid', 'snra', 'sndec', 'z', 't0', 'c', 'x1', 'x0'] 

56 

57 _lsstmwebv = None 

58 suppressHighzSN = True 

59 maxTimeSNVisible = 100. 

60 maxz = 1.2 

61 # Flux variables are convenient to display in exponential format to avoid 

62 # having them cut off 

63 variables = ['flux_u', 'flux_g', 'flux_r', 'flux_i', 'flux_z', 'flux_y'] 

64 variables += ['flux', 'flux_err', 'mag_err'] 

65 

66 override_formats = {'snra': '%8e', 'sndec': '%8e', 'c': '%8e', 

67 'x0': '%8e'} 

68 for var in variables: 

69 override_formats[var] = '%8e' 

70 # You can add parameters like fluxes and magnitudes by adding the following 

71 # variables to the list 

72 # 'flux_u', 'flux_g', 'flux_r', 'flux_i', 'flux_z', 'flux_y' , 

73 # 'mag_u', 'mag_g', 'mag_r', 'mag_i', 'mag_z', 'mag_y'] 

74 cannot_be_null = ['x0', 'z', 't0'] 

75 

76 _sn_object_cache = None 

77 

78 @astropy.utils.lazyproperty 

79 def mjdobs(self): 

80 ''' 

81 The time of observation for the catalog, which is set to be equal 

82 to obs_metadata.mjd 

83 ''' 

84 return self.obs_metadata.mjd.TAI 

85 

86 @astropy.utils.lazyproperty 

87 def badvalues(self): 

88 ''' 

89 The representation of bad values in this catalog is numpy.nan 

90 ''' 

91 return np.nan 

92 

93 @property 

94 def suppressDimSN(self): 

95 """ 

96 Boolean to decide whether to output observations of SN that are too dim 

97 should be represented in the catalog or not. By default set to True 

98 """ 

99 if not hasattr(self, '_suppressDimSN'): 

100 suppressDimSN_default = True 

101 self._suppressDimSN = suppressDimSN_default 

102 return self._suppressDimSN 

103 

104 @suppressDimSN.setter 

105 def suppressDimSN(self, suppressDimSN): 

106 """ 

107 set the value of suppressDimSN of the catalog Parameters 

108 Parameters 

109 ---------- 

110 supressDimSN : Boolean, mandatory 

111 Value to set suppressDimSN to 

112 """ 

113 self._suppressDimSN = suppressDimSN 

114 return self._suppressDimSN 

115 

116 @astropy.utils.lazyproperty 

117 def photometricparameters(self, expTime=15., nexp=2): 

118 lsstPhotometricParameters = PhotometricParameters(exptime=expTime, 

119 nexp=nexp) 

120 return lsstPhotometricParameters 

121 

122 @astropy.utils.lazyproperty 

123 def lsstBandpassDict(self): 

124 return BandpassDict.loadTotalBandpassesFromFiles() 

125 

126 @astropy.utils.lazyproperty 

127 def observedIndices(self): 

128 bandPassNames = self.obs_metadata.bandpass 

129 return [self.lsstBandpassDict.keys().index(x) for x in bandPassNames] 

130 

131 @compound('TsedFilepath', 'magNorm') 

132 def get_phosimVars(self): 

133 """ 

134 Obtain variables sedFilepath to be used to obtain unique filenames 

135 for each SED for phoSim and MagNorm which is also used. Note that aside 

136 from acting as a getter, this also writes spectra to  

137 `self.sn_sedfile_prefix`snid_mjd_band.dat for each observation of 

138 interest 

139 """ 

140 # construct the unique filename 

141 # method: snid_mjd(to 4 places of decimal)_bandpassname 

142 mjd = "_{:0.4f}_".format(self.mjdobs) 

143 mjd += self.obs_metadata.bandpass + '.dat' 

144 fnames = np.array([self.sn_sedfile_prefix + str(int(elem)) + mjd 

145 if isinstance(elem, numbers.Number) else 

146 self.sn_sedfile_prefix + str(elem) + mjd 

147 for elem in self.column_by_name('snid')], dtype='str') 

148 

149 c, x1, x0, t0, z = self.column_by_name('c'),\ 

150 self.column_by_name('x1'),\ 

151 self.column_by_name('x0'),\ 

152 self.column_by_name('t0'),\ 

153 self.column_by_name('redshift') 

154 

155 bp = Bandpass() 

156 bp.imsimBandpass() 

157 

158 magNorms = np.zeros(len(fnames)) 

159 

160 snobject = SNObject() 

161 snobject.rectifySED = True 

162 for i in range(len(self.column_by_name('snid'))): 

163 # if t0 is nan, this was set by the catalog for dim SN, or SN 

164 # outside redshift range, We will not provide a SED file for these 

165 if np.isnan(t0[i]): 

166 magNorms[i] = np.nan 

167 fnames[i] = None 

168 

169 else: 

170 snobject.set(c=c[i], x1=x1[i], x0=x0[i], t0=t0[i], 

171 z=z[i]) 

172 if snobject.modelOutSideTemporalRange == 'zero': 

173 if self.mjdobs > snobject.maxtime() or self.mjdobs < snobject.mintime(): 

174 magNorms[i] = np.nan 

175 fnames[i] = None 

176 continue 

177 

178 # SED in rest frame 

179 sed = snobject.SNObjectSourceSED(time=self.mjdobs) 

180 try: 

181 magNorms[i] = sed.calcMag(bandpass=bp) 

182 except: 

183 # sed.flambda = 1.0e-20 

184 magNorms[i] = 1000. # sed.calcMag(bandpass=bp) 

185 

186 if self.writeSedFile: 

187 sed.writeSED(fnames[i]) 

188 

189 return (fnames, magNorms) 

190 

191 def get_snid(self): 

192 # Not necessarily unique if the same galaxy hosts two SN 

193 # Use refIdCol to access the relevant id column of the dbobj 

194 # Should revert to galTileID for galaxyTiled catalogDBObj and 

195 # id for galaxyObj catalogDBObj 

196 # (email from Scott) 

197 return self.column_by_name(self.refIdCol) 

198 

199 def load_SNsed(self): 

200 """ 

201 returns a list of SN seds in `lsst.sims.photUtils.Sed` observed within 

202 the spatio-temporal range specified by obs_metadata 

203 

204 """ 

205 c, x1, x0, t0, _z, ra, dec = self.column_by_name('c'),\ 

206 self.column_by_name('x1'),\ 

207 self.column_by_name('x0'),\ 

208 self.column_by_name('t0'),\ 

209 self.column_by_name('redshift'),\ 

210 self.column_by_name('raJ2000'),\ 

211 self.column_by_name('decJ2000') 

212 

213 SNobject = SNObject() 

214 

215 raDeg = np.degrees(ra) 

216 decDeg = np.degrees(dec) 

217 

218 sedlist = [] 

219 for i in range(self.numobjs): 

220 SNobject.set(z=_z[i], c=c[i], x1=x1[i], t0=t0[i], x0=x0[i]) 

221 SNobject.setCoords(ra=raDeg[i], dec=decDeg[i]) 

222 SNobject.mwEBVfromMaps() 

223 sed = SNobject.SNObjectSED(time=self.mjdobs, 

224 bandpass=self.lsstBandpassDict, 

225 applyExitinction=True) 

226 sedlist.append(sed) 

227 

228 return sedlist 

229 

230 @property 

231 def numobjs(self): 

232 return len(self.column_by_name('id')) 

233 

234 def get_time(self): 

235 """ 

236 mjd at SALT2 'peak' 

237 """ 

238 return np.repeat(self.mjdobs, self.numobjs) 

239 

240 def get_band(self): 

241 bandname = self.obs_metadata.bandpass 

242 return np.repeat(bandname, self.numobjs) 

243 

244 @compound('flux', 'mag', 'flux_err', 'mag_err', 'adu') 

245 def get_snbrightness(self): 

246 """ 

247 getters for brightness related parameters of sn 

248 """ 

249 if self._sn_object_cache is None or len(self._sn_object_cache) > 1000000: 

250 self._sn_object_cache = {} 

251 

252 c, x1, x0, t0, _z, ra, dec = self.column_by_name('c'),\ 

253 self.column_by_name('x1'),\ 

254 self.column_by_name('x0'),\ 

255 self.column_by_name('t0'),\ 

256 self.column_by_name('redshift'),\ 

257 self.column_by_name('raJ2000'),\ 

258 self.column_by_name('decJ2000') 

259 

260 raDeg = np.degrees(ra) 

261 decDeg = np.degrees(dec) 

262 

263 ebv = self.column_by_name('EBV') 

264 id_list = self.column_by_name('snid') 

265 

266 bandname = self.obs_metadata.bandpass 

267 if isinstance(bandname, list): 

268 raise ValueError('bandname expected to be string, but is list\n') 

269 bandpass = self.lsstBandpassDict[bandname] 

270 

271 # Initialize return array so that it contains the values you would get 

272 # if you passed through a t0=self.badvalues supernova 

273 vals = np.array([[0.0]*len(t0), [np.inf]*len(t0), 

274 [np.nan]*len(t0), [np.inf]*len(t0), 

275 [0.0]*len(t0)]).transpose() 

276 

277 for i in np.where(np.logical_and(np.isfinite(t0), 

278 np.abs(self.mjdobs - t0) < self.maxTimeSNVisible))[0]: 

279 

280 if id_list[i] in self._sn_object_cache: 

281 SNobject = self._sn_object_cache[id_list[i]] 

282 else: 

283 SNobject = SNObject() 

284 SNobject.set(z=_z[i], c=c[i], x1=x1[i], t0=t0[i], x0=x0[i]) 

285 SNobject.setCoords(ra=raDeg[i], dec=decDeg[i]) 

286 SNobject.set_MWebv(ebv[i]) 

287 self._sn_object_cache[id_list[i]] = SNobject 

288 

289 if self.mjdobs <= SNobject.maxtime() and self.mjdobs >= SNobject.mintime(): 

290 

291 # Calculate fluxes 

292 fluxinMaggies = SNobject.catsimBandFlux(time=self.mjdobs, 

293 bandpassobject=bandpass) 

294 mag = SNobject.catsimBandMag(time=self.mjdobs, 

295 fluxinMaggies=fluxinMaggies, 

296 bandpassobject=bandpass) 

297 vals[i, 0] = fluxinMaggies 

298 vals[i, 1] = mag 

299 flux_err = SNobject.catsimBandFluxError(time=self.mjdobs, 

300 bandpassobject=bandpass, 

301 m5=self.obs_metadata.m5[ 

302 bandname], 

303 photParams=self.photometricparameters, 

304 fluxinMaggies=fluxinMaggies, 

305 magnitude=mag) 

306 

307 mag_err = SNobject.catsimBandMagError(time=self.mjdobs, 

308 bandpassobject=bandpass, 

309 m5=self.obs_metadata.m5[ 

310 bandname], 

311 photParams=self.photometricparameters, 

312 magnitude=mag) 

313 sed = SNobject.SNObjectSED(time=self.mjdobs, 

314 bandpass=self.lsstBandpassDict, 

315 applyExtinction=True) 

316 adu = sed.calcADU(bandpass, photParams=self.photometricparameters) 

317 vals[i, 2] = flux_err 

318 vals[i, 3] = mag_err 

319 vals[i, 4] = adu 

320 

321 return (vals[:, 0], vals[:, 1], vals[:, 2], vals[:, 3], vals[:, 4]) 

322 

323 @compound('flux_u', 'flux_g', 'flux_r', 'flux_i', 'flux_z', 'flux_y', 

324 'mag_u', 'mag_g', 'mag_r', 'mag_i', 'mag_z', 'mag_y', 

325 'adu_u', 'adu_g', 'adu_r', 'adu_i', 'adu_z', 'adu_y', 'mwebv') 

326 def get_snfluxes(self): 

327 

328 c, x1, x0, t0, _z, ra, dec = self.column_by_name('c'),\ 

329 self.column_by_name('x1'),\ 

330 self.column_by_name('x0'),\ 

331 self.column_by_name('t0'),\ 

332 self.column_by_name('redshift'),\ 

333 self.column_by_name('raJ2000'),\ 

334 self.column_by_name('decJ2000') 

335 

336 raDeg = np.degrees(ra) 

337 decDeg = np.degrees(dec) 

338 

339 snobject = SNObject() 

340 # Initialize return array 

341 vals = np.zeros(shape=(self.numobjs, 19)) 

342 for i, _ in enumerate(vals): 

343 snobject.set(z=_z[i], c=c[i], x1=x1[i], t0=t0[i], x0=x0[i]) 

344 snobject.setCoords(ra=raDeg[i], dec=decDeg[i]) 

345 snobject.mwEBVfromMaps() 

346 # Calculate fluxes 

347 vals[i, :6] = snobject.catsimManyBandFluxes(time=self.mjdobs, 

348 bandpassDict=self.lsstBandpassDict, 

349 observedBandPassInd=None) 

350 # Calculate magnitudes 

351 vals[i, 6:12] = snobject.catsimManyBandMags(time=self.mjdobs, 

352 bandpassDict=self.lsstBandpassDict, 

353 observedBandPassInd=None) 

354 

355 vals[i, 12:18] = snobject.catsimManyBandADUs(time=self.mjdobs, 

356 bandpassDict=self.lsstBandpassDict, 

357 photParams=self.photometricparameters) 

358 vals[i, 18] = snobject.ebvofMW 

359 return (vals[:, 0], vals[:, 1], vals[:, 2], vals[:, 3], 

360 vals[:, 4], vals[:, 5], vals[:, 6], vals[:, 7], 

361 vals[:, 8], vals[:, 9], vals[:, 10], vals[:, 11], 

362 vals[:, 12], vals[:, 13], vals[:, 14], vals[:, 15], 

363 vals[:, 16], vals[:, 17], vals[:, 18]) 

364 

365 #def get_EBV(self): 

366 # return self.column_by_name('EBV') 

367 

368 

369class SNIaCatalog (SNFunctionality): 

370 

371 """ 

372 `lsst.sims.catalogs.measures.instance.InstanceCatalog` class with SN 

373 characterized by the following attributes 

374 

375 Attributes 

376 ---------- 

377 column_outputs : 

378 suppressHighzSN : 

379 maxTimeSNVisible : 

380 maxz : 

381 variables : 

382 override_formats : 

383 cannot_be_null : 

384 mjdobs : 

385 badvalues position : 

386 3-tuple of floats (ra, dec, redshift), velocity : 3 tuple of floats 

387 velocity wrt host galaxy in Km/s, the supernova model (eg. SALT2) 

388 and parameters of the supernova model that predict the SED. 

389 """ 

390 

391 @compound('snra', 'sndec', 'z', 'vra', 'vdec', 'vr') 

392 def get_angularCoordinates(self): 

393 ''' 

394 Obtain the coordinates and velocity of the SN from the host galaxy 

395 

396 Returns 

397 ------- 

398 `np.ndarray` of coordinara, dec, z, vra, vdec, and vr 

399 

400 ''' 

401 hostra, hostdec, hostz = self.column_by_name('raJ2000'),\ 

402 self.column_by_name('decJ2000'),\ 

403 self.column_by_name('redshift') 

404 snra, sndec, snz, snvra, snvdec, snvr = self.SNCoordinatesFromHost( 

405 hostra, hostdec, hostz) 

406 

407 return ([snra, sndec, snz, snvra, snvdec, snvr]) 

408 

409 @compound('glon', 'glat') 

410 def get_galacticCoords(self): 

411 return _galacticFromEquatorial(self.column_by_name('snra'), self.column_by_name('sndec')) 

412 

413 @compound('c', 'x1', 'x0', 't0') 

414 def get_snparams(self): 

415 hostz, hostid, hostmu = self.column_by_name('redshift'),\ 

416 self.column_by_name('snid'),\ 

417 self.column_by_name('cosmologicalDistanceModulus') 

418 

419 vals = self.SNparamDistFromHost(hostz, hostid, hostmu) 

420 return (vals[:, 0], vals[:, 1], vals[:, 2], vals[:, 3]) 

421 

422 

423class FrozenSNCat(SNFunctionality): 

424 

425 """ 

426 `lsst.sims.catalogs.measures.instance.InstanceCatalog` class with SN 

427 characterized by the following attributes 

428 

429 Attributes 

430 ---------- 

431 column_outputs : 

432 suppressHighzSN : 

433 maxTimeSNVisible : 

434 maxz : 

435 variables : 

436 override_formats : 

437 cannot_be_null : 

438 mjdobs : 

439 badvalues position : 

440 3-tuple of floats (ra, dec, redshift), velocity : 3 tuple of floats 

441 velocity wrt host galaxy in Km/s, the supernova model (eg. SALT2) 

442 and parameters of the supernova model that predict the SED. 

443 """ 

444 

445 surveyStartDate = 59580.0 # For Kraken_1042 / Minion_1016 

446 

447 @compound('snra', 'sndec', 'z', 'vra', 'vdec', 'vr') 

448 def get_angularCoordinates(self): 

449 ''' 

450 Obtain the coordinates and velocity of the SN from the host galaxy 

451 

452 Returns 

453 ------- 

454 `np.ndarray` of coordinara, dec, z, vra, vdec, and vr 

455 

456 ''' 

457 snra, sndec, snz = self.column_by_name('raJ2000'),\ 

458 self.column_by_name('decJ2000'),\ 

459 self.column_by_name('Tredshift') 

460 snvra = np.zeros(self.numobjs) 

461 snvdec = np.zeros(self.numobjs) 

462 snvr = np.zeros(self.numobjs) 

463 

464 return (snra, sndec, snz, snvra, snvdec, snvr) 

465 

466 @compound('c', 'x1', 'x0', 't0') 

467 def get_snparams(self): 

468 """ 

469 getter for SN parameters (SALT2) 

470 """ 

471 

472 c, x1, x0 = self.column_by_name('Tc'), \ 

473 self.column_by_name('Tx1'),\ 

474 self.column_by_name('Tx0') 

475 t0 = self.column_by_name('Tt0') + self.surveyStartDate 

476 if self.suppressDimSN: 

477 t0 = np.where(np.abs(t0 - self.mjdobs) > self.maxTimeSNVisible, 

478 self.badvalues, t0) 

479 

480 return (c, x1, x0, t0) 

481