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

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

172

173

174

175

176

177

178

179

180

181

182

183

184

185

186

187

188

189

190

191

192

193

194

195

196

197

198

199

200

201

202

203

204

205

206

207

208

209

210

211

212

213

214

215

216

217

218

219

220

221

222

223

224

225

226

227

228

229

230

231

232

233

234

235

236

237

238

239

240

241

242

243

244

245

246

247

248

249

250

251

252

253

254

255

256

257

258

259

260

261

262

263

264

265

266

267

268

269

270

271

272

273

274

275

276

277

278

279

280

281

282

283

284

285

286

287

288

289

290

291

292

293

294

295

296

297

298

299

300

301

302

303

304

305

306

307

308

309

310

311

312

313

314

315

316

317

318

319

320

321

322

323

324

325

326

327

328

329

330

331

332

333

334

335

336

337

338

339

340

341

342

343

344

345

346

347

348

349

350

351

352

353

354

355

356

357

358

359

360

361

362

363

364

365

366

367

368

369

370

371

372

373

374

375

376

377

378

379

380

381

382

383

384

385

386

387

388

389

390

391

392

393

394

395

396

397

398

399

400

401

402

403

404

405

406

407

408

409

410

411

412

413

414

415

416

417

418

419

420

421

422

423

424

425

426

427

428

429

430

431

432

433

434

435

436

437

438

439

440

441

442

443

444

445

446

447

448

449

450

451

452

453

454

455

456

457

458

459

460

461

462

463

464

465

466

467

468

469

470

471

472

473

474

475

476

477

478

479

480

481

482

483

484

485

486

487

488

489

490

491

492

493

494

495

496

497

498

499

500

501

502

503

504

505

506

507

508

509

510

511

512

513

514

515

516

517

518

519

520

521

522

523

524

525

526

527

528

529

530

531

532

533

534

535

536

537

538

539

540

541

542

543

544

545

546

547

548

549

550

551

552

553

554

555

556

557

558

559

560

561

562

563

564

565

566

567

568

569

570

571

572

573

574

575

576

577

578

579

580

581

582

583

584

585

586

587

588

589

590

591

592

593

594

595

596

597

598

599

600

601

602

603

604

605

606

607

608

609

610

611

612

613

614

615

616

617

618

619

620

621

622

623

624

625

626

627

628

629

630

631

632

633

634

635

636

637

638

639

640

641

642

643

644

645

646

647

648

649

650

651

652

653

654

655

656

657

658

659

660

661

662

663

664

665

666

667

668

669

670

671

672

673

674

675

676

677

678

679

680

681

682

683

684

685

686

687

688

689

690

691

692

693

694

695

696

697

698

699

700

701

702

703

704

705

706

707

708

709

710

711

712

713

714

715

716

717

718

719

720

721

722

723

724

725

726

727

728

729

730

731

732

733

734

735

736

737

738

739

740

741

742

743

744

745

746

747

748

749

750

751

752

753

754

755

756

757

758

759

760

761

762

763

764

765

766

767

768

769

770

771

772

773

774

775

776

777

778

779

780

781

782

783

784

785

786

787

788

789

790

791

792

793

794

795

796

797

798

799

800

801

802

803

804

805

806

807

808

809

810

811

812

813

814

815

816

817

818

819

820

821

822

823

824

825

826

827

828

829

830

831

832

833

834

835

836

837

838

839

840

841

842

843

844

845

846

847

848

849

850

851

852

853

854

855

856

857

858

859

860

861

862

863

864

865

866

867

868

869

870

871

872

873

874

875

876

877

878

879

880

881

882

883

884

885

886

887

888

889

890

891

892

893

894

895

896

897

898

899

900

901

902

903

904

905

906

907

908

909

910

911

912

913

914

915

916

917

918

919

920

921

922

923

924

925

926

927

928

929

930

931

932

933

934

935

936

937

938

939

940

941

942

943

944

945

946

947

948

949

950

951

952

953

954

955

956

957

958

959

960

961

962

""" 

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

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

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

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

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

additional methods to calculate band magnitudes using the LSST software stack 

after applying MW extinction: 

- calc_mags which use the magnitude calculations in LSST stack 

- extinction which use the extinction calculations in LSST stack 

 

""" 

from builtins import str 

import numpy as np 

 

from lsst.sims.photUtils.Sed import Sed 

from lsst.sims.catUtils.dust import EBVbase 

from lsst.sims.photUtils.BandpassDict import BandpassDict 

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

from lsst.sims.photUtils.PhotometricParameters import PhotometricParameters 

 

import sncosmo 

 

__all__ = ['SNObject'] 

 

 

_sn_ax_cache = None 

_sn_bx_cache = None 

_sn_ax_bx_wavelen = None 

 

class SNObject(sncosmo.Model): 

 

""" 

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

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

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

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

extinction. 

 

Parameters 

---------- 

ra : float 

ra of the SN in degrees 

dec : float 

dec of the SN in degrees 

 

 

Attributes 

---------- 

_ra : float or None 

ra of the SN in radians 

 

_dec : float or None 

dec of the SN in radians 

 

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

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

 

ebvofMW : float or None 

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

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

will be None, leading to exceptions in extinction calculation. 

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

quantities. 

rectifySED : Bool, True by Default 

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

instead of the negative value. 

 

 

Methods 

------- 

 

Examples 

-------- 

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

>>> SNObject._ra 

>>> 0.5235987755982988 

>>> SNObject._dec 

>>> 1.0471975511965976 

""" 

 

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

""" 

Instantiate object 

 

Parameters 

---------- 

ra : float 

ra of the SN in degrees 

dec : float 

dec of the SN in degrees 

 

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

source class to define the model 

""" 

 

dust = sncosmo.CCM89Dust() 

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

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

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

 

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

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

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

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

 

self.ModelSource = source 

self.set(mwebv=0.) 

 

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

# is not provided 

self._ra = None 

self._dec = None 

 

# ra, dec is input in degrees 

# If provided, set _ra, _dec in radians 

self._hascoords = True 

if dec is None: 

self._hascoords = False 

if ra is None: 

self._hascoords = False 

 

# Satisfied that coordinates provided are floats 

if self._hascoords: 

self.setCoords(ra, dec) 

 

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

# from DustMaps 

self.lsstmwebv = EBVbase() 

self.ebvofMW = None 

if self._hascoords: 

self.mwEBVfromMaps() 

 

 

# Behavior of model outside temporal range : 

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

# are set to 0. 

self._modelOutSideTemporalRange = 'zero' 

 

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

# attribute is set to True 

self.rectifySED = True 

return 

 

@property 

def SNstate(self): 

""" 

Dictionary summarizing the state of SNObject. Can be used to 

serialize to disk, and create SNObject from SNstate 

 

Returns : Dictionary with values of parameters of the model. 

""" 

statedict = dict() 

 

# SNCosmo Parameters 

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

for param_name in self.param_names: 

statedict[param_name] = self.get(param_name) 

 

# New Attributes 

# statedict['lsstmwebv'] = self.lsstmwebv 

statedict['_ra'] = self._ra 

statedict['_dec'] = self._dec 

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

 

return statedict 

 

@classmethod 

def fromSNState(cls, snState): 

""" 

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

 

Parameters 

---------- 

snState: Dictionary summarizing the state of SNObject 

 

Returns 

------- 

Instance of SNObject class with attributes set by snstate 

 

Example 

------- 

 

""" 

# Separate into SNCosmo parameters and SNObject parameters 

dust = sncosmo.CCM89Dust() 

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

effects=[dust, dust], 

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

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

 

sncosmoParams = cls.sncosmoParamDict(snState, sncosmoModel) 

 

# Now create the class 

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

 

# Set the SNObject coordinate properties 

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

setdec, setra = False, False 

if snState['_ra'] is not None: 

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

setra = True 

if snState['_dec'] is not None: 

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

setdec = True 

if setdec and setra: 

cls.setCoords(ra, dec) 

 

# Set the SNcosmo parameters 

cls.set(**sncosmoParams) 

 

# Set the ebvofMW by hand 

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

 

return cls 

 

@property 

def modelOutSideTemporalRange(self): 

""" 

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

definition. 

""" 

return self._modelOutSideTemporalRange 

 

@modelOutSideTemporalRange.setter 

def modelOutSideTemporalRange(self, value): 

if value != 'zero': 

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

return self._modelOutSideTemporalRange 

 

 

def equivalentSNCosmoModel(self): 

""" 

returns an SNCosmo Model which is equivalent to SNObject 

""" 

snState = self.SNstate 

dust = sncosmo.CCM89Dust() 

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

effects=[dust, dust], 

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

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

 

sncosmoParams = self.sncosmoParamDict(snState, sncosmoModel) 

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

sncosmoModel.set(**sncosmoParams) 

return sncosmoModel 

 

 

@staticmethod 

def equivsncosmoParamDict(SNstate, SNCosmoModel): 

""" 

return a dictionary that contains the parameters of SNCosmoModel 

that are contained in SNstate 

 

Parameters 

---------- 

SNstate : `SNObject.SNstate`, mandatory 

Dictionary defining the state of a SNObject 

SNCosmoModel : A `sncosmo.Model` instance, mandatory 

 

Returns 

------- 

sncosmoParams: Dictionary of sncosmo parameters 

 

""" 

sncosmoParams = dict() 

for param in SNstate: 

if param in SNCosmoModel.param_names: 

sncosmoParams[param] = SNstate[param] 

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

return sncosmoParams 

 

@staticmethod 

def sncosmoParamDict(SNstate, SNCosmoModel): 

""" 

return a dictionary that contains the parameters of SNCosmoModel 

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

equivalent SNCosmo model. 

 

Parameters 

---------- 

SNstate : `SNObject.SNstate`, mandatory 

Dictionary defining the state of a SNObject 

SNCosmoModel : A `sncosmo.Model` instance, mandatory 

 

Returns 

------- 

sncosmoParams: Dictionary of sncosmo parameters 

 

""" 

sncosmoParams = dict() 

for param in SNstate: 

if param in SNCosmoModel.param_names: 

sncosmoParams[param] = SNstate[param] 

return sncosmoParams 

 

 

 

def summary(self): 

''' 

summarizes the current state of the SNObject class in a returned 

string. 

 

Parameters 

---------- 

None 

 

Returns 

------- 

Summary State in string format 

 

Examples 

-------- 

>>> t = SNObject() 

>>> print t.summary() 

''' 

state = ' SNObject Summary \n' 

 

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

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

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

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

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

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

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

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

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

 

return state 

 

def setCoords(self, ra, dec): 

""" 

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

corresponding to the given values in degrees 

 

Parameters 

---------- 

ra: float, mandatory 

the ra in degrees 

dec: float, mandatory 

dec in degrees 

 

Returns 

------- 

None 

 

Examples 

-------- 

>>> t = SNObject() 

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

>>> t._ra 

>>> 0.5235987755982988 

>>> t._dec 

>>> 1.0471975511965976 

""" 

 

if ra is None or dec is None: 

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

'coordiantes?\n') 

self._ra = np.radians(ra) 

self._dec = np.radians(dec) 

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

 

self._hascoords = True 

 

return 

 

def set_MWebv(self, value): 

""" 

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

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

 

Parameters 

---------- 

value: float, mandatory 

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

applying extinction to the SNObject spectrum 

 

Returns 

------- 

None 

 

Examples 

-------- 

>>> t = SNObject() 

>>> t.set_MWebv(0.) 

>>> 0. 

""" 

self.ebvofMW = value 

return 

 

def mwEBVfromMaps(self): 

""" 

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

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

set this attribute to None. 

 

Parameters 

---------- 

None 

 

Returns 

------- 

None 

 

Examples 

-------- 

>>> t = SNObject() 

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

>>> t.mwEBVfromMaps() 

>>> t.ebvofMW 

>>> 0.977767825127 

 

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

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

will be set to None. 

 

""" 

 

if not self._hascoords: 

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

'is `None`') 

self.ebvofMW = self.lsstmwebv.calculateEbv( 

equatorialCoordinates=self.skycoord)[0] 

return 

 

 

def redshift(self, z, cosmo): 

""" 

Redshift the instance holding the intrinsic brightness of the object 

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

magnitude in rest frame. This requires knowing the cosmology 

 

 

Parameters 

---------- 

z : float, mandatory 

redshift at which the object must be placed. 

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

specifies the cosmology. 

Returns 

------- 

None, but it changes the instance 

 

""" 

import numbers 

 

# Check that the input redshift is a scalar 

try: 

assert isinstance(z, numbers.Number) 

except: 

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

'a scalar Numeric') 

 

# Ensure that the input redshift is greater than 0. 

try: 

assert z > 0. 

except: 

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

'should be greater than 0.') 

 

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

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

self.set(z=z) 

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

return 

 

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

applyExtinction=True): 

''' 

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

requested time and wavelengths with or without extinction from MW 

according to the SED extinction methods. The wavelengths may be 

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

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

for details on handling of exceptions. 

 

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

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

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

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

the range, and the model fluxes inside 

 

Parameters 

---------- 

time: float 

time of observation 

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

array containing wavelengths in nm 

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

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

Using the dict assumes that the wavelength sampling and range 

is the same for all elements of the dict. 

 

if provided, overrides wavelen input and the SED is 

obtained at the wavelength values native to bandpass 

object. 

 

 

Returns 

------- 

`sims_photutils.sed` object containing the wavelengths and SED 

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

 

 

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

will be raised. 

Examples 

-------- 

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

''' 

 

if wavelen is None and bandpass is None: 

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

bandpassobject must be provided') 

 

# if bandpassobject present, it overrides wavelen 

if bandpass is not None: 

if isinstance(bandpass, BandpassDict): 

firstfilter = bandpass.keys()[0] 

bp = bandpass[firstfilter] 

else: 

bp = bandpass 

# remember this is in nm 

wavelen = bp.wavelen 

 

flambda = np.zeros(len(wavelen)) 

 

 

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

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

# determined by self.modelOutSideTemporalRange 

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

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

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

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

# but helps us get past this stage. 

 

flambda = flambda * np.nan 

 

# Convert to Ang 

wave = wavelen * 10.0 

mask1 = wave >= self.minwave() 

mask2 = wave <= self.maxwave() 

mask = mask1 & mask2 

wave = wave[mask] 

 

# flux density dE/dlambda returned from SNCosmo in 

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

 

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

flambda[mask] = flambda[mask] * 10.0 

else: 

# use prescription for modelOutSideTemporalRange 

if self.modelOutSideTemporalRange != 'zero': 

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

# Else Do nothing as flambda is already 0. 

# This takes precedence over being outside wavelength range 

 

if self.rectifySED: 

# Note that this converts nans into 0. 

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

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

 

if not applyExtinction: 

return SEDfromSNcosmo 

 

# Apply LSST extinction 

global _sn_ax_cache 

global _sn_bx_cache 

global _sn_ax_bx_wavelen 

if _sn_ax_bx_wavelen is None \ 

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

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

 

ax, bx = SEDfromSNcosmo.setupCCM_ab() 

_sn_ax_cache = ax 

_sn_bx_cache = bx 

_sn_ax_bx_wavelen = np.copy(wavelen) 

else: 

ax = _sn_ax_cache 

bx = _sn_bx_cache 

 

if self.ebvofMW is None: 

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

' be set by hand using set_MWebv before this' 

'stage, or by using setcoords followed by' 

'mwEBVfromMaps\n') 

 

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

return SEDfromSNcosmo 

 

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

""" 

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

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

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

SALT model being used. 

 

Parameters 

---------- 

time : float, mandatory, 

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

of days. 

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

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

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

SALT model is sampled natively. 

Returns 

------- 

`numpy.ndarray` of dtype float. 

 

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

However, it may be different for the following reasons: 

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

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

is currently implemented, where outside the model range the value 

is set to 0. 

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

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

if `self.rectifySED = True` 

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

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

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

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

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

NaN). 

""" 

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

source = self.source 

 

# Set the default value of wavelength input 

if wavelen is None: 

# use native SALT grid in Ang 

wavelen = source._wave 

else: 

#input wavelen in nm, convert to Ang 

wavelen = wavelen.copy() 

wavelen *= 10.0 

 

flambda = np.zeros(len(wavelen)) 

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

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

# determined by self.modelOutSideTemporalRange 

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

if insidephaseRange: 

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

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

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

# but helps us get past this stage. 

 

flambda = flambda * np.nan 

 

mask1 = wavelen >= source.minwave() 

mask2 = wavelen <= source.maxwave() 

mask = mask1 & mask2 

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

wave = wavelen[mask] 

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

else: 

if self.modelOutSideTemporalRange == 'zero': 

# flambda was initialized as np.zeros before start of 

# conditional 

pass 

else: 

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

 

 

# rectify the flux 

if self.rectifySED: 

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

else: 

flux = flambda 

 

 

# convert per Ang to per nm 

flux *= 10.0 

# convert ang to nm 

wavelen = wavelen / 10. 

 

# If there is zero flux at 500nm, set 

# the flux in the slot closest to 500nm 

# equal to 0.01*minimum_non_zero_flux 

# (this is so SEDs used in PhoSim can have 

# finite normalization) 

if self.rectifySED: 

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

if flux[closest_to_500nm] == 0.0: 

non_zero_flux = np.where(flux>0.0) 

if len(non_zero_flux[0])>0: 

min_non_zero = np.min(flux[non_zero_flux]) 

flux[closest_to_500nm] = 0.01*min_non_zero 

 

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

# This has the cosmology built in. 

return sed 

 

def catsimBandFlux(self, time, bandpassobject): 

""" 

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

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

 

Parameters 

---------- 

time: mandatory, float 

MJD at which band fluxes are evaluated 

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

A particular bandpass which is an instantiation of one of 

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

Returns 

------- 

float value for flux in band in units of maggies 

 

Examples 

-------- 

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

>>> LSST_BandPass = BandpassDict.loadTotalBandpassesFromFiles() 

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

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

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

>>> 1.9856857972304903e-11 

 

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

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

""" 

# Speedup for cases outside temporal range of model 

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

return 0. 

SEDfromSNcosmo = self.SNObjectSED(time=time, 

bandpass=bandpassobject) 

return SEDfromSNcosmo.calcFlux(bandpass=bandpassobject) / 3631.0 

 

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

noNan=False): 

""" 

return the magnitude in the bandpass in the AB magnitude system 

 

Parameters 

---------- 

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

LSST Catsim bandpass instance for a particular bandpass 

time : mandatory, float 

MJD at which this is evaluated 

fluxinMaggies: float, defaults to None 

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

noNan : Bool, defaults to False 

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

associated with a flux of 0. 

Returns 

------- 

float value of band magnitude in AB system 

 

Examples 

-------- 

""" 

if fluxinMaggies is None: 

fluxinMaggies = self.catsimBandFlux(bandpassobject=bandpassobject, 

time=time) 

if noNan: 

if fluxinMaggies <= 0.: 

return 200.0 

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

return -2.5 * np.log10(fluxinMaggies) 

 

def catsimBandFluxError(self, time, bandpassobject, m5, 

fluxinMaggies=None, 

magnitude=None, 

photParams=None): 

""" 

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

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

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

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

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

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

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

than a np.nan. 

 

 

Parameters 

---------- 

time: mandatory, float 

MJD at which band fluxes are evaluated 

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

A particular bandpass which is an instantiation of one of 

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

m5 : float, mandatory 

fiveSigma Depth for the sky observation 

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

describes the hardware parameters of the Observing system 

magnitude : float, defaults to None 

AB magnitude of source in bandpass. 

fluxinMaggies : float, defaults to None 

flux in Maggies for source in bandpass 

Returns 

------- 

float 

 

Examples 

-------- 

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

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

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

""" 

if fluxinMaggies is None: 

fluxinMaggies = self.catsimBandFlux(time=time, 

bandpassobject=bandpassobject) 

if magnitude is None: 

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

bandpassobject=bandpassobject, noNan=True) 

else: 

mag = magnitude 

 

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

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

fluxinMaggies = 10.0**(-0.4 * mag) 

 

if photParams is None: 

photParams = PhotometricParameters() 

 

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

m5=m5, photParams=photParams) 

return fluxinMaggies / SNR 

 

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

magnitude=None): 

""" 

return the 68 percent uncertainty on the magnitude in the bandpass 

 

Parameters 

---------- 

time: mandatory, float 

MJD at which band fluxes are evaluated 

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

A particular bandpass which is an instantiation of one of 

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

m5 : 

photParams : 

magnitude : 

 

Returns 

------- 

float 

 

Examples 

-------- 

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

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

""" 

 

if magnitude is None: 

mag = self.catsimBandMag(time=time, 

bandpassobject=bandpassobject, 

noNan=True) 

else: 

mag = magnitude 

 

bandpass = bandpassobject 

 

if photParams is None: 

photParams = PhotometricParameters() 

 

magerr = calcMagError_m5(magnitude=mag, 

bandpass=bandpassobject, 

m5=m5, 

photParams=photParams) 

return magerr[0] 

 

 

def catsimManyBandFluxes(self, time, bandpassDict, 

observedBandPassInd=None): 

""" 

return the flux in the multiple bandpasses of a bandpassDict 

indicated by observedBandPassInd in units of maggies 

 

Parameters 

---------- 

time: mandatory, float 

MJD at which band fluxes are evaluated 

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

observedBandPassInd : optional, list of integers, defaults to None 

integer correspdonding to index of the bandpasses used in the 

observation in the ordered dict bandpassDict 

Returns 

------- 

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

 

Examples 

-------- 

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

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

""" 

SEDfromSNcosmo = self.SNObjectSED(time=time, 

bandpass=bandpassDict['u']) 

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

SEDfromSNcosmo.flambdaTofnu() 

f = SEDfromSNcosmo.manyFluxCalc(bandpassDict.phiArray, 

wavelen_step=wavelen_step, 

observedBandpassInd=observedBandPassInd) 

return f / 3631. 

 

def catsimManyBandMags(self, time, bandpassDict, 

observedBandPassInd=None): 

""" 

return the flux in the bandpass in units of the flux 

the AB magnitude reference spectrum would have in the 

same band. 

 

Parameters 

---------- 

time: mandatory, float 

MJD at which band fluxes are evaluated 

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

observedBandPassInd : optional, list of integers, defaults to None 

integer correspdonding to index of the bandpasses used in the 

observation in the ordered dict bandpassDict 

Returns 

------- 

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

 

Examples 

-------- 

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

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

""" 

f = self.catsimManyBandFluxes(time, 

bandpassDict, 

observedBandPassInd) 

 

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

return -2.5 * np.log10(f) 

 

 

def catsimManyBandADUs(self, time, bandpassDict, 

photParams=None, 

observedBandPassInds=None): 

""" 

time: float, mandatory 

MJD of the observation 

 

bandpassDict: mandatory, 

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

filters 

 

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

defaults to None 

Describes the observational parameters used in specifying the 

photometry of the ovservation 

observedBandPassInd: None 

Not used now 

""" 

SEDfromSNcosmo = self.SNObjectSED(time=time, 

bandpass=bandpassDict) 

 

bandpassNames = list(bandpassDict.keys()) 

adus = np.zeros(len(bandpassNames)) 

 

for i, filt in enumerate(bandpassNames): 

bandpass = bandpassDict[filt] 

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

 

return adus