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""" 

2Class describing the SN object itself. The SN object derives from 

3SNCosmo.Model and provides a sed of SNIa based on the SALT2 model-like 

4 model often called 'salt2-extended'. This model is extended to longer 

5 wavelength ranges compared to models in Guy10 or Betoule14, at the cost of 

6 larger model variance. SNObject has additional attributes such as ra, dec. and 

7additional methods to calculate band magnitudes using the LSST software stack 

8after applying MW extinction: 

9 - calc_mags which use the magnitude calculations in LSST stack 

10 - extinction which use the extinction calculations in LSST stack 

11 

12""" 

13from builtins import str 

14import numpy as np 

15 

16from lsst.sims.photUtils.Sed import Sed 

17from lsst.sims.catUtils.dust import EBVbase 

18from lsst.sims.photUtils.BandpassDict import BandpassDict 

19from lsst.sims.photUtils.SignalToNoise import calcSNR_m5, calcMagError_m5 

20from lsst.sims.photUtils.PhotometricParameters import PhotometricParameters 

21 

22import sncosmo 

23 

24__all__ = ['SNObject'] 

25 

26 

27_sn_ax_cache = None 

28_sn_bx_cache = None 

29_sn_ax_bx_wavelen = None 

30 

31class SNObject(sncosmo.Model): 

32 

33 """ 

34 Extension of the SNCosmo `TimeSeriesModel` to include more parameters and 

35 use methods in the catsim stack. We constrain ourselves to the use of a 

36 specific SALT model for the Supernova (Salt2-Extended), and set its MW 

37 extinction to be 0, since we will use the LSST software to calculate 

38 extinction. 

39 

40 Parameters 

41 ---------- 

42 ra : float 

43 ra of the SN in degrees 

44 dec : float 

45 dec of the SN in degrees 

46 

47 

48 Attributes 

49 ---------- 

50 _ra : float or None 

51 ra of the SN in radians 

52 

53 _dec : float or None 

54 dec of the SN in radians 

55 

56 skycoord : `np.ndarray' of size 2 or None 

57 np.array([[ra], [dec]]), which are in radians 

58 

59 ebvofMW : float or None 

60 mwebv value calculated from the self.skycoord if not None, or set to 

61 a value using self.set_MWebv. If neither of these are done, this value 

62 will be None, leading to exceptions in extinction calculation. 

63 Therefore, the value must be set explicitly to 0. to get unextincted 

64 quantities. 

65 rectifySED : Bool, True by Default 

66 if the SED is negative at the requested time and wavelength, return 0. 

67 instead of the negative value. 

68 

69 

70 Methods 

71 ------- 

72 

73 Examples 

74 -------- 

75 >>> SNObject = SNObject(ra=30., dec=60.) 

76 >>> SNObject._ra 

77 >>> 0.5235987755982988 

78 >>> SNObject._dec 

79 >>> 1.0471975511965976 

80 """ 

81 

82 def __init__(self, ra=None, dec=None, source='salt2-extended'): 

83 """ 

84 Instantiate object 

85 

86 Parameters 

87 ---------- 

88 ra : float 

89 ra of the SN in degrees 

90 dec : float 

91 dec of the SN in degrees 

92 

93 source : instance of `sncosmo.SALT2Source`, optional, defaults to using salt2-extended  

94 source class to define the model 

95 """ 

96 

97 dust = sncosmo.CCM89Dust() 

98 sncosmo.Model.__init__(self, source=source, effects=[dust, dust], 

99 effect_names=['host', 'mw'], 

100 effect_frames=['rest', 'obs']) 

101 

102 # Current implementation of Model has a default value of mwebv = 0. 

103 # ie. no extinction, but this is not part of the API, so should not 

104 # depend on it, set explicitly in order to unextincted SED from 

105 # SNCosmo. We will use catsim extinction from `lsst.sims.photUtils`. 

106 

107 self.ModelSource = source 

108 self.set(mwebv=0.) 

109 

110 # self._ra, self._dec is initialized as None for cases where ra, dec 

111 # is not provided 

112 self._ra = None 

113 self._dec = None 

114 

115 # ra, dec is input in degrees 

116 # If provided, set _ra, _dec in radians 

117 self._hascoords = True 

118 if dec is None: 

119 self._hascoords = False 

120 if ra is None: 

121 self._hascoords = False 

122 

123 # Satisfied that coordinates provided are floats 

124 if self._hascoords: 

125 self.setCoords(ra, dec) 

126 

127 # For values of ra, dec, the E(B-V) is calculated directly 

128 # from DustMaps 

129 self.lsstmwebv = EBVbase() 

130 self.ebvofMW = None 

131 if self._hascoords: 

132 self.mwEBVfromMaps() 

133 

134 

135 # Behavior of model outside temporal range : 

136 # if 'zero' then all fluxes outside the temporal range of the model 

137 # are set to 0. 

138 self._modelOutSideTemporalRange = 'zero' 

139 

140 # SED will be rectified to 0. for negative values of SED if this 

141 # attribute is set to True 

142 self.rectifySED = True 

143 return 

144 

145 @property 

146 def SNstate(self): 

147 """ 

148 Dictionary summarizing the state of SNObject. Can be used to 

149 serialize to disk, and create SNObject from SNstate 

150 

151 Returns : Dictionary with values of parameters of the model. 

152 """ 

153 statedict = dict() 

154 

155 # SNCosmo Parameters 

156 statedict['ModelSource'] = self.source.name 

157 for param_name in self.param_names: 

158 statedict[param_name] = self.get(param_name) 

159 

160 # New Attributes 

161 # statedict['lsstmwebv'] = self.lsstmwebv 

162 statedict['_ra'] = self._ra 

163 statedict['_dec'] = self._dec 

164 statedict['MWE(B-V)'] = self.ebvofMW 

165 

166 return statedict 

167 

168 @classmethod 

169 def fromSNState(cls, snState): 

170 """ 

171 creates an instance of SNObject with a state described by snstate. 

172 

173 Parameters 

174 ---------- 

175 snState: Dictionary summarizing the state of SNObject 

176 

177 Returns 

178 ------- 

179 Instance of SNObject class with attributes set by snstate 

180 

181 Example 

182 ------- 

183 

184 """ 

185 # Separate into SNCosmo parameters and SNObject parameters 

186 dust = sncosmo.CCM89Dust() 

187 sncosmoModel = sncosmo.Model(source=snState['ModelSource'], 

188 effects=[dust, dust], 

189 effect_names=['host', 'mw'], 

190 effect_frames=['rest', 'obs']) 

191 

192 sncosmoParams = cls.sncosmoParamDict(snState, sncosmoModel) 

193 

194 # Now create the class 

195 cls = SNObject(source=snState['ModelSource']) 

196 

197 # Set the SNObject coordinate properties 

198 # Have to be careful to not convert `None` type objects to degrees 

199 setdec, setra = False, False 

200 if snState['_ra'] is not None: 

201 ra = np.degrees(snState['_ra']) 

202 setra = True 

203 if snState['_dec'] is not None: 

204 dec = np.degrees(snState['_dec']) 

205 setdec = True 

206 if setdec and setra: 

207 cls.setCoords(ra, dec) 

208 

209 # Set the SNcosmo parameters 

210 cls.set(**sncosmoParams) 

211 

212 # Set the ebvofMW by hand 

213 cls.ebvofMW = snState['MWE(B-V)'] 

214 

215 return cls 

216 

217 @property 

218 def modelOutSideTemporalRange(self): 

219 """ 

220 Defines the behavior of the model when sampled at times beyond the model 

221 definition. 

222 """ 

223 return self._modelOutSideTemporalRange 

224 

225 @modelOutSideTemporalRange.setter 

226 def modelOutSideTemporalRange(self, value): 

227 if value != 'zero': 

228 raise ValueError('Model not implemented, defaulting to zero method\n') 

229 return self._modelOutSideTemporalRange 

230 

231 

232 def equivalentSNCosmoModel(self): 

233 """ 

234 returns an SNCosmo Model which is equivalent to SNObject 

235 """ 

236 snState = self.SNstate 

237 dust = sncosmo.CCM89Dust() 

238 sncosmoModel = sncosmo.Model(source=snState['ModelSource'], 

239 effects=[dust, dust], 

240 effect_names=['host', 'mw'], 

241 effect_frames=['rest', 'obs']) 

242 

243 sncosmoParams = self.sncosmoParamDict(snState, sncosmoModel) 

244 sncosmoParams['mwebv'] = snState['MWE(B-V)'] 

245 sncosmoModel.set(**sncosmoParams) 

246 return sncosmoModel 

247 

248 

249 @staticmethod 

250 def equivsncosmoParamDict(SNstate, SNCosmoModel): 

251 """ 

252 return a dictionary that contains the parameters of SNCosmoModel 

253 that are contained in SNstate 

254 

255 Parameters 

256 ---------- 

257 SNstate : `SNObject.SNstate`, mandatory 

258 Dictionary defining the state of a SNObject 

259 SNCosmoModel : A `sncosmo.Model` instance, mandatory 

260 

261 Returns 

262 ------- 

263 sncosmoParams: Dictionary of sncosmo parameters 

264 

265 """ 

266 sncosmoParams = dict() 

267 for param in SNstate: 

268 if param in SNCosmoModel.param_names: 

269 sncosmoParams[param] = SNstate[param] 

270 sncosmoParams['mwebv'] = SNstate['MWE(B-V)'] 

271 return sncosmoParams 

272 

273 @staticmethod 

274 def sncosmoParamDict(SNstate, SNCosmoModel): 

275 """ 

276 return a dictionary that contains the parameters of SNCosmoModel 

277 that are contained in SNstate. Note that this does not return the 

278 equivalent SNCosmo model. 

279 

280 Parameters 

281 ---------- 

282 SNstate : `SNObject.SNstate`, mandatory 

283 Dictionary defining the state of a SNObject 

284 SNCosmoModel : A `sncosmo.Model` instance, mandatory 

285 

286 Returns 

287 ------- 

288 sncosmoParams: Dictionary of sncosmo parameters 

289 

290 """ 

291 sncosmoParams = dict() 

292 for param in SNstate: 

293 if param in SNCosmoModel.param_names: 

294 sncosmoParams[param] = SNstate[param] 

295 return sncosmoParams 

296 

297 

298 

299 def summary(self): 

300 ''' 

301 summarizes the current state of the SNObject class in a returned 

302 string. 

303 

304 Parameters 

305 ---------- 

306 None 

307 

308 Returns 

309 ------- 

310 Summary State in string format 

311 

312 Examples 

313 -------- 

314 >>> t = SNObject() 

315 >>> print t.summary() 

316 ''' 

317 state = ' SNObject Summary \n' 

318 

319 state += 'Model = ' + '\n' 

320 state += 'z = ' + str(self.get('z')) + '\n' 

321 state += 'c = ' + str(self.get('c')) + '\n' 

322 state += 'x1 = ' + str(self.get('x1')) + '\n' 

323 state += 'x0 = ' + str(self.get('x0')) + '\n' 

324 state += 't0 = ' + str(self.get('t0')) + '\n' 

325 state += 'ra = ' + str(self._ra) + ' in radians \n' 

326 state += 'dec = ' + str(self._dec) + ' in radians \n' 

327 state += 'MW E(B-V) = ' + str(self.ebvofMW) + '\n' 

328 

329 return state 

330 

331 def setCoords(self, ra, dec): 

332 """ 

333 set the ra and dec coordinate of SNObject to values in radians 

334 corresponding to the given values in degrees 

335 

336 Parameters 

337 ---------- 

338 ra: float, mandatory 

339 the ra in degrees 

340 dec: float, mandatory 

341 dec in degrees 

342 

343 Returns 

344 ------- 

345 None 

346 

347 Examples 

348 -------- 

349 >>> t = SNObject() 

350 >>> t.setCoords(ra=30., dec=90.) 

351 >>> t._ra 

352 >>> 0.5235987755982988 

353 >>> t._dec 

354 >>> 1.0471975511965976 

355 """ 

356 

357 if ra is None or dec is None: 

358 raise ValueError('Why try to set coordinates without full' 

359 'coordiantes?\n') 

360 self._ra = np.radians(ra) 

361 self._dec = np.radians(dec) 

362 self.skycoord = np.array([[self._ra], [self._dec]]) 

363 

364 self._hascoords = True 

365 

366 return 

367 

368 def set_MWebv(self, value): 

369 """ 

370 if mwebv value is known, this can be used to set the attribute 

371 ebvofMW of the SNObject class to the value (float). 

372 

373 Parameters 

374 ---------- 

375 value: float, mandatory 

376 value of mw extinction parameter E(B-V) in mags to be used in 

377 applying extinction to the SNObject spectrum 

378 

379 Returns 

380 ------- 

381 None 

382 

383 Examples 

384 -------- 

385 >>> t = SNObject() 

386 >>> t.set_MWebv(0.) 

387 >>> 0. 

388 """ 

389 self.ebvofMW = value 

390 return 

391 

392 def mwEBVfromMaps(self): 

393 """ 

394 set the attribute ebvofMW of the class from the ra and dec 

395 of the SN. If the ra or dec attribute of the class is None, 

396 set this attribute to None. 

397 

398 Parameters 

399 ---------- 

400 None 

401 

402 Returns 

403 ------- 

404 None 

405 

406 Examples 

407 -------- 

408 >>> t = SNObject() 

409 >>> t.setCoords(ra=30., dec=60.) 

410 >>> t.mwEBVfromMaps() 

411 >>> t.ebvofMW 

412 >>> 0.977767825127 

413 

414 .. note:: This function must be run after the class has attributes ra 

415 and dec set. In case it is run before this, the mwebv value 

416 will be set to None. 

417 

418 """ 

419 

420 if not self._hascoords: 

421 raise ValueError('Cannot Calculate EBV from dust maps if ra or dec' 

422 'is `None`') 

423 self.ebvofMW = self.lsstmwebv.calculateEbv( 

424 equatorialCoordinates=self.skycoord)[0] 

425 return 

426 

427 

428 def redshift(self, z, cosmo): 

429 """ 

430 Redshift the instance holding the intrinsic brightness of the object 

431 fixed. By intrinsic brightness here, we mean the BessellB band asbolute 

432 magnitude in rest frame. This requires knowing the cosmology 

433 

434 

435 Parameters 

436 ---------- 

437 z : float, mandatory 

438 redshift at which the object must be placed. 

439 cosmo : instance of `astropy.cosmology` objects, mandatory 

440 specifies the cosmology. 

441 Returns 

442 ------- 

443 None, but it changes the instance 

444 

445 """ 

446 import numbers 

447 

448 # Check that the input redshift is a scalar 

449 try: 

450 assert isinstance(z, numbers.Number) 

451 except: 

452 raise TypeError('The argument z in method redshift should be' 

453 'a scalar Numeric') 

454 

455 # Ensure that the input redshift is greater than 0. 

456 try: 

457 assert z > 0. 

458 except: 

459 raise ValueError('The argument z in the method SNObject.redshift' 

460 'should be greater than 0.') 

461 

462 # Find the current value of the rest frame BessellB AB magnitude 

463 peakAbsMag = self.source_peakabsmag('BessellB', 'AB', cosmo=cosmo) 

464 self.set(z=z) 

465 self.set_source_peakabsmag(peakAbsMag, 'BessellB', 'AB', cosmo=cosmo) 

466 return 

467 

468 def SNObjectSED(self, time, wavelen=None, bandpass=None, 

469 applyExtinction=True): 

470 ''' 

471 return a `lsst.sims.photUtils.sed` object from the SN model at the 

472 requested time and wavelengths with or without extinction from MW 

473 according to the SED extinction methods. The wavelengths may be 

474 obtained from a `lsst.sims.Bandpass` object or a `lsst.sims.BandpassDict` 

475 object instead. (Currently, these have the same wavelengths). See notes 

476 for details on handling of exceptions. 

477 

478 If the sed is requested at times outside the validity range of the 

479 model, the flux density is returned as 0. If the time is within the 

480 range of validity of the model, but the wavelength range requested 

481 is outside the range, then the returned fluxes are np.nan outside 

482 the range, and the model fluxes inside 

483 

484 Parameters 

485 ---------- 

486 time: float 

487 time of observation 

488 wavelen: `np.ndarray` of floats, optional, defaults to None 

489 array containing wavelengths in nm 

490 bandpass: `lsst.sims.photUtils.Bandpass` object or 

491 `lsst.sims.photUtils.BandpassDict`, optional, defaults to `None`. 

492 Using the dict assumes that the wavelength sampling and range 

493 is the same for all elements of the dict. 

494 

495 if provided, overrides wavelen input and the SED is 

496 obtained at the wavelength values native to bandpass 

497 object. 

498 

499 

500 Returns 

501 ------- 

502 `sims_photutils.sed` object containing the wavelengths and SED 

503 values from the SN at time time in units of ergs/cm^2/sec/nm 

504 

505 

506 .. note: If both wavelen and bandpassobject are `None` then exception, 

507 will be raised. 

508 Examples 

509 -------- 

510 >>> sed = SN.SNObjectSED(time=0., wavelen=wavenm) 

511 ''' 

512 

513 if wavelen is None and bandpass is None: 

514 raise ValueError('A non None input to either wavelen or\ 

515 bandpassobject must be provided') 

516 

517 # if bandpassobject present, it overrides wavelen 

518 if bandpass is not None: 

519 if isinstance(bandpass, BandpassDict): 

520 firstfilter = bandpass.keys()[0] 

521 bp = bandpass[firstfilter] 

522 else: 

523 bp = bandpass 

524 # remember this is in nm 

525 wavelen = bp.wavelen 

526 

527 flambda = np.zeros(len(wavelen)) 

528 

529 

530 # self.mintime() and self.maxtime() are properties describing 

531 # the ranges of SNCosmo.Model in time. Behavior beyond this is  

532 # determined by self.modelOutSideTemporalRange 

533 if (time >= self.mintime()) and (time <= self.maxtime()): 

534 # If SNCosmo is requested a SED value beyond the wavelength range 

535 # of model it will crash. Try to prevent that by returning np.nan for 

536 # such wavelengths. This will still not help band flux calculations 

537 # but helps us get past this stage. 

538 

539 flambda = flambda * np.nan 

540 

541 # Convert to Ang 

542 wave = wavelen * 10.0 

543 mask1 = wave >= self.minwave() 

544 mask2 = wave <= self.maxwave() 

545 mask = mask1 & mask2 

546 wave = wave[mask] 

547 

548 # flux density dE/dlambda returned from SNCosmo in 

549 # ergs/cm^2/sec/Ang, convert to ergs/cm^2/sec/nm 

550 

551 flambda[mask] = self.flux(time=time, wave=wave) 

552 flambda[mask] = flambda[mask] * 10.0 

553 else: 

554 # use prescription for modelOutSideTemporalRange 

555 if self.modelOutSideTemporalRange != 'zero': 

556 raise NotImplementedError('Model not implemented, change to zero\n') 

557 # Else Do nothing as flambda is already 0. 

558 # This takes precedence over being outside wavelength range 

559 

560 if self.rectifySED: 

561 # Note that this converts nans into 0. 

562 flambda = np.where(flambda > 0., flambda, 0.) 

563 SEDfromSNcosmo = Sed(wavelen=wavelen, flambda=flambda) 

564 

565 if not applyExtinction: 

566 return SEDfromSNcosmo 

567 

568 # Apply LSST extinction 

569 global _sn_ax_cache 

570 global _sn_bx_cache 

571 global _sn_ax_bx_wavelen 

572 if _sn_ax_bx_wavelen is None \ 

573 or len(wavelen)!=len(_sn_ax_bx_wavelen) \ 

574 or (wavelen!=_sn_ax_bx_wavelen).any(): 

575 

576 ax, bx = SEDfromSNcosmo.setupCCM_ab() 

577 _sn_ax_cache = ax 

578 _sn_bx_cache = bx 

579 _sn_ax_bx_wavelen = np.copy(wavelen) 

580 else: 

581 ax = _sn_ax_cache 

582 bx = _sn_bx_cache 

583 

584 if self.ebvofMW is None: 

585 raise ValueError('ebvofMW attribute cannot be None Type and must' 

586 ' be set by hand using set_MWebv before this' 

587 'stage, or by using setcoords followed by' 

588 'mwEBVfromMaps\n') 

589 

590 SEDfromSNcosmo.addDust(a_x=ax, b_x=bx, ebv=self.ebvofMW) 

591 return SEDfromSNcosmo 

592 

593 def SNObjectSourceSED(self, time, wavelen=None): 

594 """ 

595 Return the rest Frame SED of SNObject at the phase corresponding to 

596 time, at rest frame wavelengths wavelen. If wavelen is None, 

597 then the SED is sampled at the rest frame wavelengths native to the 

598 SALT model being used. 

599 

600 Parameters 

601 ---------- 

602 time : float, mandatory, 

603 observer frame time at which the SED has been requested in units 

604 of days. 

605 wavelen : `np.ndarray`, optional, defaults to native SALT wavelengths 

606 array of wavelengths in the rest frame of the supernova in units 

607 of nm. If None, this defaults to the wavelengths at which the 

608 SALT model is sampled natively. 

609 Returns 

610 ------- 

611 `numpy.ndarray` of dtype float. 

612 

613 .. note: The result should usually match the SALT source spectrum. 

614 However, it may be different for the following reasons: 

615 1. If the time of observation is outside the model range, the values 

616 have to be inserted using additional models. Here only one model 

617 is currently implemented, where outside the model range the value 

618 is set to 0. 

619 2. If the wavelengths are beyond the range of the SALT model, the SED 

620 flambda values are set to `np.nan` and these are actually set to 0. 

621 if `self.rectifySED = True` 

622 3. If the `flambda` values of the SALT model are negative which happens 

623 in the less sampled phases of the model, these values are set to 0, 

624 if `self.rectifySED` = True. (Note: if `self.rectifySED` = True, then 

625 care will be taken to make sure that the flux at 500nm is not exactly 

626 zero, since that will cause PhoSim normalization of the SED to be 

627 NaN). 

628 """ 

629 phase = (time - self.get('t0')) / (1. + self.get('z')) 

630 source = self.source 

631 

632 # Set the default value of wavelength input 

633 if wavelen is None: 

634 # use native SALT grid in Ang 

635 wavelen = source._wave 

636 else: 

637 #input wavelen in nm, convert to Ang 

638 wavelen = wavelen.copy() 

639 wavelen *= 10.0 

640 

641 flambda = np.zeros(len(wavelen)) 

642 # self.mintime() and self.maxtime() are properties describing 

643 # the ranges of SNCosmo.Model in time. Behavior beyond this is  

644 # determined by self.modelOutSideTemporalRange 

645 insidephaseRange = (phase <= source.maxphase())and(phase >= source.minphase()) 

646 if insidephaseRange: 

647 # If SNCosmo is requested a SED value beyond the wavelength range 

648 # of model it will crash. Try to prevent that by returning np.nan for 

649 # such wavelengths. This will still not help band flux calculations 

650 # but helps us get past this stage. 

651 

652 flambda = flambda * np.nan 

653 

654 mask1 = wavelen >= source.minwave() 

655 mask2 = wavelen <= source.maxwave() 

656 mask = mask1 & mask2 

657 # Where we have to calculate fluxes because it is not `np.nan` 

658 wave = wavelen[mask] 

659 flambda[mask] = source.flux(phase, wave) 

660 else: 

661 if self.modelOutSideTemporalRange == 'zero': 

662 # flambda was initialized as np.zeros before start of 

663 # conditional 

664 pass 

665 else: 

666 raise NotImplementedError('Only modelOutSideTemporalRange=="zero" implemented') 

667 

668 

669 # rectify the flux 

670 if self.rectifySED: 

671 flux = np.where(flambda>0., flambda, 0.) 

672 else: 

673 flux = flambda 

674 

675 

676 # convert per Ang to per nm 

677 flux *= 10.0 

678 # convert ang to nm 

679 wavelen = wavelen / 10. 

680 

681 # If there is zero flux at 500nm, set 

682 # the flux in the slot closest to 500nm 

683 # equal to 0.01*minimum_non_zero_flux 

684 # (this is so SEDs used in PhoSim can have 

685 # finite normalization) 

686 if self.rectifySED: 

687 closest_to_500nm = np.argmin(np.abs(wavelen-500.0)) 

688 if flux[closest_to_500nm] == 0.0: 

689 non_zero_flux = np.where(flux>0.0) 

690 if len(non_zero_flux[0])>0: 

691 min_non_zero = np.min(flux[non_zero_flux]) 

692 flux[closest_to_500nm] = 0.01*min_non_zero 

693 

694 sed = Sed(wavelen=wavelen, flambda=flux) 

695 # This has the cosmology built in. 

696 return sed 

697 

698 def catsimBandFlux(self, time, bandpassobject): 

699 """ 

700 return the flux in the bandpass in units of maggies which is the flux 

701 the AB magnitude reference spectrum would have in the same band. 

702 

703 Parameters 

704 ---------- 

705 time: mandatory, float 

706 MJD at which band fluxes are evaluated 

707 bandpassobject: mandatory, `lsst.sims.photUtils.BandPass` object 

708 A particular bandpass which is an instantiation of one of 

709 (u, g, r, i, z, y) 

710 Returns 

711 ------- 

712 float value for flux in band in units of maggies 

713 

714 Examples 

715 -------- 

716 >>> bandpassnames = ['u', 'g', 'r', 'i', 'z', 'y'] 

717 >>> LSST_BandPass = BandpassDict.loadTotalBandpassesFromFiles() 

718 >>> SN = SNObject(ra=30., dec=-60.) 

719 >>> SN.set(z=0.96, t0=571181, x1=2.66, c=0.353, x0=1.796112e-06) 

720 >>> SN.catsimBandFlux(bandpassobject=LSST_BandPass['r'], time=571190.) 

721 >>> 1.9856857972304903e-11 

722 

723 .. note: If there is an unphysical value of sed in 

724 the wavelength range, it produces a flux of `np.nan` 

725 """ 

726 # Speedup for cases outside temporal range of model 

727 if time <= self.mintime() or time >= self.maxtime() : 

728 return 0. 

729 SEDfromSNcosmo = self.SNObjectSED(time=time, 

730 bandpass=bandpassobject) 

731 return SEDfromSNcosmo.calcFlux(bandpass=bandpassobject) / 3631.0 

732 

733 def catsimBandMag(self, bandpassobject, time, fluxinMaggies=None, 

734 noNan=False): 

735 """ 

736 return the magnitude in the bandpass in the AB magnitude system 

737 

738 Parameters 

739 ---------- 

740 bandpassobject : mandatory,`sims.photUtils.BandPass` instances 

741 LSST Catsim bandpass instance for a particular bandpass 

742 time : mandatory, float 

743 MJD at which this is evaluated 

744 fluxinMaggies: float, defaults to None 

745 provide the flux in maggies, if not provided, this will be evaluated 

746 noNan : Bool, defaults to False 

747 If True, an AB magnitude of 200.0 rather than nan values is 

748 associated with a flux of 0. 

749 Returns 

750 ------- 

751 float value of band magnitude in AB system 

752 

753 Examples 

754 -------- 

755 """ 

756 if fluxinMaggies is None: 

757 fluxinMaggies = self.catsimBandFlux(bandpassobject=bandpassobject, 

758 time=time) 

759 if noNan: 

760 if fluxinMaggies <= 0.: 

761 return 200.0 

762 with np.errstate(divide='ignore', invalid='ignore'): 

763 return -2.5 * np.log10(fluxinMaggies) 

764 

765 def catsimBandFluxError(self, time, bandpassobject, m5, 

766 fluxinMaggies=None, 

767 magnitude=None, 

768 photParams=None): 

769 """ 

770 return the flux uncertainty in the bandpass in units 'maggies' 

771 (the flux the AB magnitude reference spectrum would have in the 

772 same band.) for a source of given brightness. The source brightness 

773 may be calculated, but the need for calculation is overridden by a 

774 provided flux in bandpass (in units of maggies) which itself may be 

775 overridden by a provided magnitude. If the provided/calculated flux 

776 is 0. or negative the magnitude calculated is taken to be 200.0 rather 

777 than a np.nan. 

778 

779 

780 Parameters 

781 ---------- 

782 time: mandatory, float 

783 MJD at which band fluxes are evaluated 

784 bandpassobject: mandatory, `lsst.sims.photUtils.BandPass` object 

785 A particular bandpass which is an instantiation of one of 

786 (u, g, r, i, z, y) 

787 m5 : float, mandatory 

788 fiveSigma Depth for the sky observation 

789 photParams : instance of `sims.photUtils.PhotometricParameters`, defaults to `None`  

790 describes the hardware parameters of the Observing system 

791 magnitude : float, defaults to None 

792 AB magnitude of source in bandpass. 

793 fluxinMaggies : float, defaults to None 

794 flux in Maggies for source in bandpass 

795 Returns 

796 ------- 

797 float 

798 

799 Examples 

800 -------- 

801 .. note: If there is an unphysical value of sed the fluxinMaggies might 

802 be `np.nan`. The magnitude calculated from this is calculated using `noNan` 

803 and is therefore 200.0 rather than `np.nan`.  

804 """ 

805 if fluxinMaggies is None: 

806 fluxinMaggies = self.catsimBandFlux(time=time, 

807 bandpassobject=bandpassobject) 

808 if magnitude is None: 

809 mag = self.catsimBandMag(time=time, fluxinMaggies=fluxinMaggies, 

810 bandpassobject=bandpassobject, noNan=True) 

811 else: 

812 mag = magnitude 

813 

814 # recalculate fluxinMaggies as the previous one might have been `np.nan` 

815 # the noise is contaminated if this is `np.nan` 

816 fluxinMaggies = 10.0**(-0.4 * mag) 

817 

818 if photParams is None: 

819 photParams = PhotometricParameters() 

820 

821 SNR, gamma = calcSNR_m5(magnitude=mag, bandpass=bandpassobject, 

822 m5=m5, photParams=photParams) 

823 return fluxinMaggies / SNR 

824 

825 def catsimBandMagError(self, time, bandpassobject, m5, photParams=None, 

826 magnitude=None): 

827 """ 

828 return the 68 percent uncertainty on the magnitude in the bandpass 

829 

830 Parameters 

831 ---------- 

832 time: mandatory, float 

833 MJD at which band fluxes are evaluated 

834 bandpassobject: mandatory, `lsst.sims.photUtils.BandPass` object 

835 A particular bandpass which is an instantiation of one of 

836 (u, g, r, i, z, y) 

837 m5 : 

838 photParams : 

839 magnitude : 

840 

841 Returns 

842 ------- 

843 float 

844 

845 Examples 

846 -------- 

847 .. note: If there is an unphysical value of sed in 

848 the wavelength range, it produces a flux of `np.nan` 

849 """ 

850 

851 if magnitude is None: 

852 mag = self.catsimBandMag(time=time, 

853 bandpassobject=bandpassobject, 

854 noNan=True) 

855 else: 

856 mag = magnitude 

857 

858 bandpass = bandpassobject 

859 

860 if photParams is None: 

861 photParams = PhotometricParameters() 

862 

863 magerr = calcMagError_m5(magnitude=mag, 

864 bandpass=bandpassobject, 

865 m5=m5, 

866 photParams=photParams) 

867 return magerr[0] 

868 

869 

870 def catsimManyBandFluxes(self, time, bandpassDict, 

871 observedBandPassInd=None): 

872 """ 

873 return the flux in the multiple bandpasses of a bandpassDict 

874 indicated by observedBandPassInd in units of maggies 

875 

876 Parameters 

877 ---------- 

878 time: mandatory, float 

879 MJD at which band fluxes are evaluated 

880 bandpassDict: mandatory, `lsst.sims.photUtils.BandpassDict` instance 

881 observedBandPassInd : optional, list of integers, defaults to None 

882 integer correspdonding to index of the bandpasses used in the 

883 observation in the ordered dict bandpassDict 

884 Returns 

885 ------- 

886 `~numpy.ndarray` of length =len(observedBandPassInd) 

887 

888 Examples 

889 -------- 

890 .. note: If there is an unphysical value of sed in 

891 the wavelength range, it produces a flux of `np.nan` 

892 """ 

893 SEDfromSNcosmo = self.SNObjectSED(time=time, 

894 bandpass=bandpassDict['u']) 

895 wavelen_step = np.diff(SEDfromSNcosmo.wavelen)[0] 

896 SEDfromSNcosmo.flambdaTofnu() 

897 f = SEDfromSNcosmo.manyFluxCalc(bandpassDict.phiArray, 

898 wavelen_step=wavelen_step, 

899 observedBandpassInd=observedBandPassInd) 

900 return f / 3631. 

901 

902 def catsimManyBandMags(self, time, bandpassDict, 

903 observedBandPassInd=None): 

904 """ 

905 return the flux in the bandpass in units of the flux 

906 the AB magnitude reference spectrum would have in the 

907 same band. 

908 

909 Parameters 

910 ---------- 

911 time: mandatory, float 

912 MJD at which band fluxes are evaluated 

913 bandpassDict: mandatory, `lsst.sims.photUtils.BandpassDict` instance 

914 observedBandPassInd : optional, list of integers, defaults to None 

915 integer correspdonding to index of the bandpasses used in the 

916 observation in the ordered dict bandpassDict 

917 Returns 

918 ------- 

919 `~numpy.ndarray` of length =len(observedBandPassInd) 

920 

921 Examples 

922 -------- 

923 .. note: If there is an unphysical value of sed in 

924 the wavelength range, it produces a flux of `np.nan` 

925 """ 

926 f = self.catsimManyBandFluxes(time, 

927 bandpassDict, 

928 observedBandPassInd) 

929 

930 with np.errstate(invalid='ignore', divide='ignore'): 

931 return -2.5 * np.log10(f) 

932 

933 

934 def catsimManyBandADUs(self, time, bandpassDict, 

935 photParams=None, 

936 observedBandPassInds=None): 

937 """ 

938 time: float, mandatory 

939 MJD of the observation 

940 

941 bandpassDict: mandatory, 

942 Dictionary of instances of `sims.photUtils.Bandpass` for 

943 filters 

944 

945 photParams: Instance of `sims.photUtils.PhotometricParameters`, optional, 

946 defaults to None 

947 Describes the observational parameters used in specifying the 

948 photometry of the ovservation 

949 observedBandPassInd: None 

950 Not used now 

951 """ 

952 SEDfromSNcosmo = self.SNObjectSED(time=time, 

953 bandpass=bandpassDict) 

954 

955 bandpassNames = list(bandpassDict.keys()) 

956 adus = np.zeros(len(bandpassNames)) 

957 

958 for i, filt in enumerate(bandpassNames): 

959 bandpass = bandpassDict[filt] 

960 adus[i] = SEDfromSNcosmo.calcADU(bandpass, photParams=photParams) 

961 

962 return adus