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

from builtins import zip 

from builtins import object 

import numpy as np 

import ephem 

from lsst.sims.utils import haversine, _raDecFromAltAz, _altAzPaFromRaDec, Site, ObservationMetaData, calcLmstLast 

import warnings 

from lsst.sims.skybrightness.utils import wrapRA, mjd2djd 

from .interpComponents import ScatteredStar, Airglow, LowerAtm, UpperAtm, MergedSpec, TwilightInterp, MoonInterp, ZodiacalInterp 

from lsst.sims.photUtils import Sed 

 

 

__all__ = ['justReturn', 'stupidFast_RaDec2AltAz', 'stupidFast_altAz2RaDec', 'SkyModel'] 

 

 

def justReturn(inval): 

""" 

Really, just return the input. 

 

Parameters 

---------- 

input : anything 

 

Returns 

------- 

input : anything 

Just return whatever you sent in. 

""" 

return inval 

 

 

def inrange(inval, minimum=-1., maximum=1.): 

""" 

Make sure values are within min/max 

""" 

inval = np.array(inval) 

below = np.where(inval < minimum) 

inval[below] = minimum 

above = np.where(inval > maximum) 

inval[above] = maximum 

return inval 

 

 

def stupidFast_RaDec2AltAz(ra, dec, lat, lon, mjd): 

""" 

Convert Ra,Dec to Altitude and Azimuth. 

 

Coordinate transformation is killing performance. Just use simple equations to speed it up 

and ignore abberation, precesion, nutation, nutrition, etc. 

 

Parameters 

---------- 

ra : array_like 

RA, in radians. 

dec : array_like 

Dec, in radians. Must be same length as `ra`. 

lat : float 

Latitude of the observatory in radians. 

lon : float 

Longitude of the observatory in radians. 

mjd : float 

Modified Julian Date. 

 

Returns 

------- 

alt : numpy.array 

Altitude, same length as `ra` and `dec`. Radians. 

az : numpy.array 

Azimuth, same length as `ra` and `dec`. Radians. 

""" 

lmst, last = calcLmstLast(mjd, lon) 

lmst = lmst/12.*np.pi # convert to rad 

ha = lmst-ra 

sindec = np.sin(dec) 

sinlat = np.sin(lat) 

coslat = np.cos(lat) 

sinalt = sindec*sinlat+np.cos(dec)*coslat*np.cos(ha) 

sinalt = inrange(sinalt) 

alt = np.arcsin(sinalt) 

cosaz = (sindec-np.sin(alt)*sinlat)/(np.cos(alt)*coslat) 

cosaz = inrange(cosaz) 

az = np.arccos(cosaz) 

signflip = np.where(np.sin(ha) > 0) 

az[signflip] = 2.*np.pi-az[signflip] 

return alt, az 

 

 

def stupidFast_altAz2RaDec(alt, az, lat, lon, mjd): 

""" 

Convert alt, az to RA, Dec without taking into account abberation, precesion, diffraction, ect. 

 

Parameters 

---------- 

alt : numpy.array 

Altitude, same length as `ra` and `dec`. Radians. 

az : numpy.array 

Azimuth, same length as `ra` and `dec`. Must be same length as `alt`. Radians. 

lat : float 

Latitude of the observatory in radians. 

lon : float 

Longitude of the observatory in radians. 

mjd : float 

Modified Julian Date. 

 

Returns 

------- 

ra : array_like 

RA, in radians. 

dec : array_like 

Dec, in radians. 

""" 

lmst, last = calcLmstLast(mjd, lon) 

lmst = lmst/12.*np.pi # convert to rad 

sindec = np.sin(lat)*np.sin(alt) + np.cos(lat)*np.cos(alt)*np.cos(az) 

sindec = inrange(sindec) 

dec = np.arcsin(sindec) 

ha = np.arctan2(-np.sin(az)*np.cos(alt), -np.cos(az)*np.sin(lat)*np.cos(alt)+np.sin(alt)*np.cos(lat)) 

ra = (lmst-ha) 

raneg = np.where(ra < 0) 

ra[raneg] = ra[raneg] + 2.*np.pi 

return ra, dec 

 

 

def calcAzRelMoon(azs, moonAz): 

azRelMoon = wrapRA(azs - moonAz) 

125 ↛ 129line 125 didn't jump to line 129, because the condition on line 125 was never false if isinstance(azs, np.ndarray): 

over = np.where(azRelMoon > np.pi) 

azRelMoon[over] = 2. * np.pi - azRelMoon[over] 

else: 

if azRelMoon > np.pi: 

azRelMoon = 2.0 * np.pi - azRelMoon 

return azRelMoon 

 

 

class SkyModel(object): 

 

def __init__(self, observatory='LSST', 

twilight=True, zodiacal=True, moon=True, 

airglow=True, lowerAtm=False, upperAtm=False, scatteredStar=False, 

mergedSpec=True, mags=False, preciseAltAz=False, airmass_limit=2.5): 

""" 

Instatiate the SkyModel. This loads all the required template spectra/magnitudes 

that will be used for interpolation. 

 

Parameters 

---------- 

Observatory : Site object 

object with attributes lat, lon, elev. But default loads LSST. 

 

twilight : bool (True) 

Include twilight component (True) 

zodiacal : bool (True) 

Include zodiacal light component (True) 

moon : bool (True) 

Include scattered moonlight component (True) 

airglow : bool (True) 

Include airglow component 

lowerAtm : bool (False) 

Include lower atmosphere component. This component is part of `mergedSpec`. 

upperAtm : bool (False) 

Include upper atmosphere component. This component is part of `mergedSpec`. 

scatteredStar : bool (False) 

Include scattered starlight component. This component is part of `mergedSpec`. 

mergedSpec : bool (True) 

Compute the lowerAtm, upperAtm, and scatteredStar simultaneously since they are all 

functions of only airmass. 

mags : bool (False) 

By default, the sky model computes a 17,001 element spectrum. If `mags` is True, 

the model will return the LSST ugrizy magnitudes (in that order). 

preciseAltAz : bool (False) 

If False, use the fast alt, az to ra, dec coordinate 

transformations that do not take abberation, diffraction, etc 

into account. Results in errors up to ~1.5 degrees, 

but an order of magnitude faster than coordinate transforms in sims_utils. 

airmass_limit : float (2.5) 

Most of the models are only accurate to airmass 2.5. If set higher, airmass values 

higher than 2.5 are set to 2.5. 

""" 

 

self.moon = moon 

self.lowerAtm = lowerAtm 

self.twilight = twilight 

self.zodiacal = zodiacal 

self.upperAtm = upperAtm 

self.airglow = airglow 

self.scatteredStar = scatteredStar 

self.mergedSpec = mergedSpec 

self.mags = mags 

self.preciseAltAz = preciseAltAz 

 

# set this as a way to track if coords have been set 

self.azs = None 

 

# Airmass limit. 

self.airmassLimit = airmass_limit 

 

if self.mags: 

self.npix = 6 

else: 

self.npix = 17001 

 

self.components = {'moon': self.moon, 'lowerAtm': self.lowerAtm, 'twilight': self.twilight, 

'upperAtm': self.upperAtm, 'airglow': self.airglow, 'zodiacal': self.zodiacal, 

'scatteredStar': self.scatteredStar, 'mergedSpec': self.mergedSpec} 

 

# Check that the merged component isn't being run with other components 

mergedComps = [self.lowerAtm, self.upperAtm, self.scatteredStar] 

for comp in mergedComps: 

208 ↛ 209line 208 didn't jump to line 209, because the condition on line 208 was never true if comp & self.mergedSpec: 

warnings.warn("Adding component multiple times to the final output spectra.") 

 

interpolators = {'scatteredStar': ScatteredStar, 'airglow': Airglow, 'lowerAtm': LowerAtm, 

'upperAtm': UpperAtm, 'mergedSpec': MergedSpec, 'moon': MoonInterp, 

'zodiacal': ZodiacalInterp, 'twilight': TwilightInterp} 

 

# Load up the interpolation objects for each component 

self.interpObjs = {} 

for key in self.components: 

if self.components[key]: 

self.interpObjs[key] = interpolators[key](mags=self.mags) 

 

# Set up a pyephem observatory object 

222 ↛ 223line 222 didn't jump to line 223, because the condition on line 222 was never true if hasattr(observatory, 'latitude_rad') & hasattr(observatory, 'longitude_rad') & hasattr(observatory, 'height'): 

self.telescope = observatory 

self.Observatory = ephem.Observer() 

self.Observatory.lat = self.telescope.latitude_rad 

self.Observatory.lon = self.telescope.longitude_rad 

self.Observatory.elevation = self.telescope.height 

228 ↛ 235line 228 didn't jump to line 235, because the condition on line 228 was never false elif observatory == 'LSST': 

self.telescope = Site('LSST') 

self.Observatory = ephem.Observer() 

self.Observatory.lat = self.telescope.latitude_rad 

self.Observatory.lon = self.telescope.longitude_rad 

self.Observatory.elevation = self.telescope.height 

else: 

self.Observatory = observatory 

 

# Note that observing conditions have not been set 

self.paramsSet = False 

 

def setComponents(self, twilight=True, zodiacal=True, moon=True, 

airglow=True, lowerAtm=False, upperAtm=False, scatteredStar=False, 

mergedSpec=True): 

""" 

Convience function for turning on/off different sky components. 

""" 

self.moon = moon 

self.lowerAtm = lowerAtm 

self.twilight = twilight 

self.zodiacal = zodiacal 

self.upperAtm = upperAtm 

self.airglow = airglow 

self.scatteredStar = scatteredStar 

self.mergedSpec = mergedSpec 

 

def _initPoints(self): 

""" 

Set up an array for all the interpolation points 

""" 

 

names = ['airmass', 'nightTimes', 'alt', 'az', 'azRelMoon', 'moonSunSep', 'moonAltitude', 

'altEclip', 'azEclipRelSun', 'sunAlt', 'azRelSun', 'solarFlux'] 

types = [float]*len(names) 

self.points = np.zeros(self.npts, list(zip(names, types))) 

 

def setRaDecMjd(self, lon, lat, mjd, degrees=False, azAlt=False, solarFlux=130., 

filterNames=['u', 'g', 'r', 'i', 'z', 'y']): 

""" 

Set the sky parameters by computing the sky conditions on a given MJD and sky location. 

 

 

 

lon: Longitude-like (RA or Azimuth). Can be single number, list, or numpy array 

lat: Latitude-like (Dec or Altitude) 

mjd: Modified Julian Date for the calculation. Must be single number. 

degrees: (False) Assumes lon and lat are radians unless degrees=True 

azAlt: (False) Assume lon, lat are RA, Dec unless azAlt=True 

solarFlux: solar flux in SFU Between 50 and 310. Default=130. 1 SFU=10^4 Jy. 

filterNames: list of fitlers to return magnitudes for (if initialized with mags=True). 

""" 

self.filterNames = filterNames 

if self.mags: 

self.npix = len(self.filterNames) 

# Wrap in array just in case single points were passed 

if np.size(lon) == 1: 

lon = np.array([lon]).ravel() 

lat = np.array([lat]).ravel() 

else: 

lon = np.array(lon) 

lat = np.array(lat) 

if degrees: 

self.ra = np.radians(lon) 

self.dec = np.radians(lat) 

else: 

self.ra = lon 

self.dec = lat 

296 ↛ 297line 296 didn't jump to line 297, because the condition on line 296 was never true if np.size(mjd) > 1: 

raise ValueError('mjd must be single value.') 

self.mjd = mjd 

if azAlt: 

self.azs = self.ra.copy() 

self.alts = self.dec.copy() 

302 ↛ 303line 302 didn't jump to line 303, because the condition on line 302 was never true if self.preciseAltAz: 

self.ra, self.dec = _raDecFromAltAz(self.alts, self.azs, 

ObservationMetaData(mjd=self.mjd, site=self.telescope)) 

else: 

self.ra, self.dec = stupidFast_altAz2RaDec(self.alts, self.azs, 

self.telescope.latitude_rad, 

self.telescope.longitude_rad, mjd) 

else: 

310 ↛ 311line 310 didn't jump to line 311, because the condition on line 310 was never true if self.preciseAltAz: 

self.alts, self.azs, pa = _altAzPaFromRaDec(self.ra, self.dec, 

ObservationMetaData(mjd=self.mjd, 

site=self.telescope)) 

else: 

self.alts, self.azs = stupidFast_RaDec2AltAz(self.ra, self.dec, 

self.telescope.latitude_rad, 

self.telescope.longitude_rad, mjd) 

 

self.npts = self.ra.size 

self._initPoints() 

 

self.solarFlux = solarFlux 

self.points['solarFlux'] = self.solarFlux 

 

self._setupPointGrid() 

 

self.paramsSet = True 

 

# Assume large airmasses are the same as airmass=2.5 

to_fudge = np.where((self.points['airmass'] > 2.5) & (self.points['airmass'] < self.airmassLimit)) 

self.points['airmass'][to_fudge] = 2.499 

self.points['alt'][to_fudge] = np.pi/2-np.arccos(1./self.airmassLimit) 

 

# Interpolate the templates to the set paramters 

self._interpSky() 

 

def setRaDecAltAzMjd(self, ra, dec, alt, az, mjd, degrees=False, solarFlux=130., 

filterNames=['u', 'g', 'r', 'i', 'z', 'y']): 

""" 

Set the sky parameters by computing the sky conditions on a given MJD and sky location. 

 

Use if you already have alt az coordinates so you can skip the coordinate conversion. 

""" 

self.filterNames = filterNames 

345 ↛ 348line 345 didn't jump to line 348, because the condition on line 345 was never false if self.mags: 

self.npix = len(self.filterNames) 

# Wrap in array just in case single points were passed 

348 ↛ 349line 348 didn't jump to line 349, because the condition on line 348 was never true if not type(ra).__module__ == np.__name__: 

if np.size(ra) == 1: 

ra = np.array([ra]).ravel() 

dec = np.array([dec]).ravel() 

alt = np.array(alt).ravel() 

az = np.array(az).ravel() 

else: 

ra = np.array(ra) 

dec = np.array(dec) 

alt = np.array(alt) 

az = np.array(az) 

359 ↛ 360line 359 didn't jump to line 360, because the condition on line 359 was never true if degrees: 

self.ra = np.radians(ra) 

self.dec = np.radians(dec) 

self.alts = np.radians(alt) 

self.azs = np.radians(az) 

else: 

self.ra = ra 

self.dec = dec 

self.azs = az 

self.alts = alt 

369 ↛ 370line 369 didn't jump to line 370, because the condition on line 369 was never true if np.size(mjd) > 1: 

raise ValueError('mjd must be single value.') 

self.mjd = mjd 

 

self.npts = self.ra.size 

self._initPoints() 

 

self.solarFlux = solarFlux 

self.points['solarFlux'] = self.solarFlux 

 

self._setupPointGrid() 

 

self.paramsSet = True 

 

# Assume large airmasses are the same as airmass=2.5 

to_fudge = np.where((self.points['airmass'] > 2.5) & (self.points['airmass'] < self.airmassLimit)) 

self.points['airmass'][to_fudge] = 2.5 

self.points['alt'][to_fudge] = np.pi/2-np.arccos(1./self.airmassLimit) 

 

# Interpolate the templates to the set paramters 

self._interpSky() 

 

def getComputedVals(self): 

""" 

Return the intermediate values that are caluculated by setRaDecMjd and used for interpolation. 

All of these values are also accesible as class atributes, this is a convience method to grab them 

all at once and document the formats. 

 

Returns 

------- 

out : dict 

Dictionary of all the intermediate calculated values that may be of use outside 

(the key:values in the output dict) 

ra : numpy.array 

RA of the interpolation points (radians) 

dec : np.array 

Dec of the interpolation points (radians) 

alts : np.array 

Altitude (radians) 

azs : np.array 

Azimuth of interpolation points (radians) 

airmass : np.array 

Airmass values for each point, computed via 1./np.cos(np.pi/2.-self.alts). 

solarFlux : float 

The solar flux used (SFU). 

sunAz : float 

Azimuth of the sun (radians) 

sunAlt : float 

Altitude of the sun (radians) 

sunRA : float 

RA of the sun (radians) 

sunDec : float 

Dec of the sun (radians) 

azRelSun : np.array 

Azimuth of each point relative to the sun (0=same direction as sun) (radians) 

moonAz : float 

Azimuth of the moon (radians) 

moonAlt : float 

Altitude of the moon (radians) 

moonRA : float 

RA of the moon (radians) 

moonDec : float 

Dec of the moon (radians). Note, if you want distances 

moonPhase : float 

Phase of the moon (0-100) 

moonSunSep : float 

Seperation of moon and sun (degrees) 

azRelMoon : np.array 

Azimuth of each point relative to teh moon 

eclipLon : np.array 

Ecliptic longitude (radians) of each point 

eclipLat : np.array 

Ecliptic latitude (radians) of each point 

sunEclipLon: np.array 

Ecliptic longitude (radians) of each point with the sun at longitude zero 

 

Note that since the alt and az can be calculated using the fast approximation, if one wants 

to compute the distance between the the points and the sun or moon, it is probably better to 

use the ra,dec positions rather than the alt,az positions. 

""" 

 

result = {} 

attributes = ['ra', 'dec', 'alts', 'azs', 'airmass', 'solarFlux', 'moonPhase', 

'moonAz', 'moonAlt', 'sunAlt', 'sunAz', 'azRelSun', 'moonSunSep', 

'azRelMoon', 'eclipLon', 'eclipLat', 'moonRA', 'moonDec', 'sunRA', 

'sunDec', 'sunEclipLon'] 

 

for attribute in attributes: 

457 ↛ 460line 457 didn't jump to line 460, because the condition on line 457 was never false if hasattr(self, attribute): 

result[attribute] = getattr(self, attribute) 

else: 

result[attribute] = None 

 

return result 

 

def _setupPointGrid(self): 

""" 

Setup the points for the interpolation functions. 

""" 

# Switch to Dublin Julian Date for pyephem 

self.Observatory.date = mjd2djd(self.mjd) 

 

sun = ephem.Sun() 

sun.compute(self.Observatory) 

self.sunAlt = sun.alt 

self.sunAz = sun.az 

self.sunRA = sun.ra 

self.sunDec = sun.dec 

 

# Compute airmass the same way as ESO model 

self.airmass = 1./np.cos(np.pi/2.-self.alts) 

 

self.points['airmass'] = self.airmass 

self.points['nightTimes'] = 2 

self.points['alt'] = self.alts 

self.points['az'] = self.azs 

 

if self.twilight: 

self.points['sunAlt'] = self.sunAlt 

self.azRelSun = wrapRA(self.azs - self.sunAz) 

self.points['azRelSun'] = self.azRelSun 

 

if self.moon: 

moon = ephem.Moon() 

moon.compute(self.Observatory) 

self.moonPhase = moon.phase 

self.moonAlt = moon.alt 

self.moonAz = moon.az 

self.moonRA = moon.ra 

self.moonDec = moon.dec 

# Calc azimuth relative to moon 

self.azRelMoon = calcAzRelMoon(self.azs, self.moonAz) 

self.moonTargSep = haversine(self.azs, self.alts, self.moonAz, self.moonAlt) 

self.points['moonAltitude'] += np.degrees(self.moonAlt) 

self.points['azRelMoon'] += self.azRelMoon 

self.moonSunSep = self.moonPhase/100.*180. 

self.points['moonSunSep'] += self.moonSunSep 

 

if self.zodiacal: 

self.eclipLon = np.zeros(self.npts) 

self.eclipLat = np.zeros(self.npts) 

 

for i, temp in enumerate(self.ra): 

eclip = ephem.Ecliptic(ephem.Equatorial(self.ra[i], self.dec[i], epoch='2000')) 

self.eclipLon[i] += eclip.lon 

self.eclipLat[i] += eclip.lat 

# Subtract off the sun ecliptic longitude 

sunEclip = ephem.Ecliptic(sun) 

self.sunEclipLon = sunEclip.lon 

self.points['altEclip'] += self.eclipLat 

self.points['azEclipRelSun'] += wrapRA(self.eclipLon - self.sunEclipLon) 

 

self.mask = np.where((self.airmass > self.airmassLimit) | (self.airmass < 1.))[0] 

self.goodPix = np.where((self.airmass <= self.airmassLimit) & (self.airmass >= 1.))[0] 

 

def setParams(self, airmass=1., azs=90., alts=None, moonPhase=31.67, moonAlt=45., 

moonAz=0., sunAlt=-12., sunAz=0., sunEclipLon=0., 

eclipLon=135., eclipLat=90., degrees=True, solarFlux=130., 

filterNames=['u', 'g', 'r', 'i', 'z', 'y']): 

""" 

Set parameters manually. 

Note, you can put in unphysical combinations of paramters if you want to 

(e.g., put a full moon at zenith at sunset). 

if the alts kwarg is set it will override the airmass kwarg. 

MoonPhase is percent of moon illuminated (0-100) 

""" 

 

# Convert all values to radians for internal use. 

self.filterNames = filterNames 

538 ↛ 539line 538 didn't jump to line 539, because the condition on line 538 was never true if self.mags: 

self.npix = len(self.filterNames) 

if degrees: 

convertFunc = np.radians 

else: 

convertFunc = justReturn 

 

self.solarFlux = solarFlux 

self.sunAlt = convertFunc(sunAlt) 

self.moonPhase = moonPhase 

self.moonAlt = convertFunc(moonAlt) 

self.moonAz = convertFunc(moonAz) 

self.eclipLon = convertFunc(eclipLon) 

self.eclipLat = convertFunc(eclipLat) 

self.sunEclipLon = convertFunc(sunEclipLon) 

self.azs = convertFunc(azs) 

554 ↛ 558line 554 didn't jump to line 558, because the condition on line 554 was never false if alts is not None: 

self.airmass = 1./np.cos(np.pi/2.-convertFunc(alts)) 

self.alts = convertFunc(alts) 

else: 

self.airmass = airmass 

self.alts = np.pi/2.-np.arccos(1./airmass) 

self.moonTargSep = haversine(self.azs, self.alts, moonAz, self.moonAlt) 

self.npts = np.size(self.airmass) 

self._initPoints() 

 

self.points['airmass'] = self.airmass 

self.points['nightTimes'] = 2 

self.points['alt'] = self.alts 

self.points['az'] = self.azs 

self.azRelMoon = calcAzRelMoon(self.azs, self.moonAz) 

self.points['moonAltitude'] += np.degrees(self.moonAlt) 

self.points['azRelMoon'] = self.azRelMoon 

self.points['moonSunSep'] += self.moonPhase/100.*180. 

 

self.eclipLon = convertFunc(eclipLon) 

self.eclipLat = convertFunc(eclipLat) 

 

self.sunEclipLon = convertFunc(sunEclipLon) 

self.points['altEclip'] += self.eclipLat 

self.points['azEclipRelSun'] += wrapRA(self.eclipLon - self.sunEclipLon) 

 

self.sunAz = convertFunc(sunAz) 

self.points['sunAlt'] = self.sunAlt 

self.points['azRelSun'] = wrapRA(self.azs - self.sunAz) 

self.points['solarFlux'] = solarFlux 

 

self.paramsSet = True 

 

self.mask = np.where((self.airmass > self.airmassLimit) | (self.airmass < 1.))[0] 

self.goodPix = np.where((self.airmass <= self.airmassLimit) & (self.airmass >= 1.))[0] 

# Interpolate the templates to the set paramters 

self._interpSky() 

 

def _interpSky(self): 

""" 

Interpolate the template spectra to the set RA, Dec and MJD. 

 

the results are stored as attributes of the class: 

.wave = the wavelength in nm 

.spec = array of spectra with units of ergs/s/cm^2/nm 

""" 

 

601 ↛ 602line 601 didn't jump to line 602, because the condition on line 601 was never true if not self.paramsSet: 

raise ValueError( 

'No parameters have been set. Must run setRaDecMjd or setParams before running interpSky.') 

 

# set up array to hold the resulting spectra for each ra, dec point. 

self.spec = np.zeros((self.npts, self.npix), dtype=float) 

 

# Rebuild the components dict so things can be turned on/off 

self.components = {'moon': self.moon, 'lowerAtm': self.lowerAtm, 'twilight': self.twilight, 

'upperAtm': self.upperAtm, 'airglow': self.airglow, 'zodiacal': self.zodiacal, 

'scatteredStar': self.scatteredStar, 'mergedSpec': self.mergedSpec} 

 

# Loop over each component and add it to the result. 

mask = np.ones(self.npts) 

for key in self.components: 

if self.components[key]: 

result = self.interpObjs[key](self.points[self.goodPix], filterNames=self.filterNames) 

# Make sure the component has something 

619 ↛ 620line 619 didn't jump to line 620, because the condition on line 619 was never true if np.size(result['spec']) == 0: 

self.spec[self.mask, :] = np.nan 

return 

if np.max(result['spec']) > 0: 

mask[np.where(np.sum(result['spec'], axis=1) == 0)] = 0 

self.spec[self.goodPix] += result['spec'] 

if not hasattr(self, 'wave'): 

self.wave = result['wave'] 

else: 

628 ↛ 629line 628 didn't jump to line 629, because the condition on line 628 was never true if not np.allclose(result['wave'], self.wave, rtol=1e-5, atol=1e-5): 

warnings.warn('Wavelength arrays of components do not match.') 

630 ↛ 632line 630 didn't jump to line 632, because the condition on line 630 was never false if self.airmassLimit <= 2.5: 

self.spec[np.where(mask == 0), :] = 0 

self.spec[self.mask, :] = np.nan 

 

def returnWaveSpec(self): 

""" 

Return the wavelength and spectra. 

Wavelenth in nm 

spectra is flambda in ergs/cm^2/s/nm 

""" 

640 ↛ 641line 640 didn't jump to line 641, because the condition on line 640 was never true if self.azs is None: 

raise ValueError('No coordinates set. Use setRaDecMjd, setRaDecAltAzMjd, or setParams methods before calling returnWaveSpec.') 

642 ↛ 643line 642 didn't jump to line 643, because the condition on line 642 was never true if self.mags: 

raise ValueError('SkyModel set to interpolate magnitudes. Initialize object with mags=False') 

# Mask out high airmass points 

# self.spec[self.mask] *= 0 

return self.wave, self.spec 

 

def returnMags(self, bandpasses=None): 

""" 

Convert the computed spectra to a magnitude using the supplied bandpass, 

or, if self.mags=True, return the mags in the LSST filters 

 

If mags=True when initialized, return mags returns an structured array with 

dtype names u,g,r,i,z,y. 

 

bandpasses: optional dictionary with bandpass name keys and bandpass object values. 

 

""" 

659 ↛ 660line 659 didn't jump to line 660, because the condition on line 659 was never true if self.azs is None: 

raise ValueError('No coordinates set. Use setRaDecMjd, setRaDecAltAzMjd, or setParams methods before calling returnMags.') 

 

if self.mags: 

663 ↛ 664line 663 didn't jump to line 664, because the condition on line 663 was never true if bandpasses: 

warnings.warn('Ignoring set bandpasses and returning LSST ugrizy.') 

mags = -2.5*np.log10(self.spec)+np.log10(3631.) 

# Mask out high airmass 

mags[self.mask] *= np.nan 

mags = mags.swapaxes(0, 1) 

magsBack = {} 

for i, f in enumerate(self.filterNames): 

magsBack[f] = mags[i] 

else: 

magsBack = {} 

for key in bandpasses: 

mags = np.zeros(self.npts, dtype=float)-666 

tempSed = Sed() 

isThrough = np.where(bandpasses[key].sb > 0) 

minWave = bandpasses[key].wavelen[isThrough].min() 

maxWave = bandpasses[key].wavelen[isThrough].max() 

inBand = np.where((self.wave >= minWave) & (self.wave <= maxWave)) 

for i, ra in enumerate(self.ra): 

# Check that there is flux in the band, otherwise calcMag fails 

683 ↛ 681line 683 didn't jump to line 681, because the condition on line 683 was never false if np.max(self.spec[i, inBand]) > 0: 

tempSed.setSED(self.wave, flambda=self.spec[i, :]) 

mags[i] = tempSed.calcMag(bandpasses[key]) 

# Mask out high airmass 

mags[self.mask] *= np.nan 

magsBack[key] = mags 

return magsBack