Coverage for python/lsst/fgcmcal/fgcmMakeLut.py: 22%

251 statements  

« prev     ^ index     » next       coverage.py v7.3.0, created at 2023-09-06 15:25 +0000

1# See COPYRIGHT file at the top of the source tree. 

2# 

3# This file is part of fgcmcal. 

4# 

5# Developed for the LSST Data Management System. 

6# This product includes software developed by the LSST Project 

7# (https://www.lsst.org). 

8# See the COPYRIGHT file at the top-level directory of this distribution 

9# for details of code ownership. 

10# 

11# This program is free software: you can redistribute it and/or modify 

12# it under the terms of the GNU General Public License as published by 

13# the Free Software Foundation, either version 3 of the License, or 

14# (at your option) any later version. 

15# 

16# This program is distributed in the hope that it will be useful, 

17# but WITHOUT ANY WARRANTY; without even the implied warranty of 

18# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 

19# GNU General Public License for more details. 

20# 

21# You should have received a copy of the GNU General Public License 

22# along with this program. If not, see <https://www.gnu.org/licenses/>. 

23"""Make a look-up-table (LUT) for FGCM calibration. 

24 

25This task computes a look-up-table for the range in expected atmosphere 

26variation and variation in instrumental throughput (as tracked by the 

27transmission_filter products). By pre-computing linearized integrals, 

28the FGCM fit is orders of magnitude faster for stars with a broad range 

29of colors and observing bands, yielding precision at the 1-2 mmag level. 

30 

31Computing a LUT requires running MODTRAN or with a pre-generated 

32atmosphere table packaged with fgcm. 

33""" 

34 

35import numpy as np 

36 

37import lsst.pex.config as pexConfig 

38import lsst.pipe.base as pipeBase 

39from lsst.pipe.base import connectionTypes 

40import lsst.afw.table as afwTable 

41import lsst.afw.cameraGeom as afwCameraGeom 

42from lsst.afw.image import TransmissionCurve 

43from .utilities import lookupStaticCalibrations 

44 

45import fgcm 

46 

47__all__ = ['FgcmMakeLutParametersConfig', 'FgcmMakeLutConfig', 'FgcmMakeLutTask'] 

48 

49 

50class FgcmMakeLutConnections(pipeBase.PipelineTaskConnections, 

51 dimensions=('instrument',), 

52 defaultTemplates={}): 

53 camera = connectionTypes.PrerequisiteInput( 

54 doc="Camera instrument", 

55 name="camera", 

56 storageClass="Camera", 

57 dimensions=("instrument",), 

58 isCalibration=True, 

59 ) 

60 

61 transmission_optics = connectionTypes.PrerequisiteInput( 

62 doc="Optics transmission curve information", 

63 name="transmission_optics", 

64 storageClass="TransmissionCurve", 

65 dimensions=("instrument",), 

66 isCalibration=True, 

67 deferLoad=True, 

68 ) 

69 

70 transmission_sensor = connectionTypes.PrerequisiteInput( 

71 doc="Sensor transmission curve information", 

72 name="transmission_sensor", 

73 storageClass="TransmissionCurve", 

74 dimensions=("instrument", "detector",), 

75 lookupFunction=lookupStaticCalibrations, 

76 isCalibration=True, 

77 deferLoad=True, 

78 multiple=True, 

79 ) 

80 

81 transmission_filter = connectionTypes.PrerequisiteInput( 

82 doc="Filter transmission curve information", 

83 name="transmission_filter", 

84 storageClass="TransmissionCurve", 

85 dimensions=("band", "instrument", "physical_filter",), 

86 lookupFunction=lookupStaticCalibrations, 

87 isCalibration=True, 

88 deferLoad=True, 

89 multiple=True, 

90 ) 

91 

92 fgcmLookUpTable = connectionTypes.Output( 

93 doc=("Atmosphere + instrument look-up-table for FGCM throughput and " 

94 "chromatic corrections."), 

95 name="fgcmLookUpTable", 

96 storageClass="Catalog", 

97 dimensions=("instrument",), 

98 ) 

99 

100 fgcmStandardAtmosphere = connectionTypes.Output( 

101 doc="Standard atmosphere used for FGCM calibration.", 

102 name="fgcm_standard_atmosphere", 

103 storageClass="TransmissionCurve", 

104 dimensions=("instrument",), 

105 ) 

106 

107 fgcmStandardPassbands = connectionTypes.Output( 

108 doc="Standard passbands used for FGCM calibration.", 

109 name="fgcm_standard_passband", 

110 storageClass="TransmissionCurve", 

111 dimensions=("instrument", "physical_filter"), 

112 multiple=True, 

113 ) 

114 

115 def __init__(self, *, config=None): 

116 if not config.doOpticsTransmission: 

117 del self.transmission_optics 

118 if not config.doSensorTransmission: 

119 del self.transmission_sensor 

120 

121 

122class FgcmMakeLutParametersConfig(pexConfig.Config): 

123 """Config for parameters if atmosphereTableName not available""" 

124 # TODO: When DM-16511 is done, it will be possible to get the 

125 # telescope elevation directly from the camera. 

126 elevation = pexConfig.Field( 

127 doc="Telescope elevation (m)", 

128 dtype=float, 

129 default=None, 

130 ) 

131 pmbRange = pexConfig.ListField( 

132 doc=("Barometric Pressure range (millibar) " 

133 "Recommended range depends on the site."), 

134 dtype=float, 

135 default=None, 

136 ) 

137 pmbSteps = pexConfig.Field( 

138 doc="Barometric Pressure number of steps", 

139 dtype=int, 

140 default=5, 

141 ) 

142 pwvRange = pexConfig.ListField( 

143 doc=("Precipitable Water Vapor range (mm) " 

144 "Recommended range depends on the site."), 

145 dtype=float, 

146 default=None, 

147 ) 

148 pwvSteps = pexConfig.Field( 

149 doc="Precipitable Water Vapor number of steps", 

150 dtype=int, 

151 default=15, 

152 ) 

153 o3Range = pexConfig.ListField( 

154 doc="Ozone range (dob)", 

155 dtype=float, 

156 default=[220.0, 310.0], 

157 ) 

158 o3Steps = pexConfig.Field( 

159 doc="Ozone number of steps", 

160 dtype=int, 

161 default=3, 

162 ) 

163 tauRange = pexConfig.ListField( 

164 doc="Aerosol Optical Depth range (unitless)", 

165 dtype=float, 

166 default=[0.002, 0.35], 

167 ) 

168 tauSteps = pexConfig.Field( 

169 doc="Aerosol Optical Depth number of steps", 

170 dtype=int, 

171 default=11, 

172 ) 

173 alphaRange = pexConfig.ListField( 

174 doc="Aerosol alpha range (unitless)", 

175 dtype=float, 

176 default=[0.0, 2.0], 

177 ) 

178 alphaSteps = pexConfig.Field( 

179 doc="Aerosol alpha number of steps", 

180 dtype=int, 

181 default=9, 

182 ) 

183 zenithRange = pexConfig.ListField( 

184 doc="Zenith angle range (degree)", 

185 dtype=float, 

186 default=[0.0, 70.0], 

187 ) 

188 zenithSteps = pexConfig.Field( 

189 doc="Zenith angle number of steps", 

190 dtype=int, 

191 default=21, 

192 ) 

193 # Note that the standard atmosphere parameters depend on the observatory 

194 # and elevation, and so these should be set on a per-camera basis. 

195 pmbStd = pexConfig.Field( 

196 doc=("Standard Atmosphere pressure (millibar); " 

197 "Recommended default depends on the site."), 

198 dtype=float, 

199 default=None, 

200 ) 

201 pwvStd = pexConfig.Field( 

202 doc=("Standard Atmosphere PWV (mm); " 

203 "Recommended default depends on the site."), 

204 dtype=float, 

205 default=None, 

206 ) 

207 o3Std = pexConfig.Field( 

208 doc="Standard Atmosphere O3 (dob)", 

209 dtype=float, 

210 default=263.0, 

211 ) 

212 tauStd = pexConfig.Field( 

213 doc="Standard Atmosphere aerosol optical depth", 

214 dtype=float, 

215 default=0.03, 

216 ) 

217 alphaStd = pexConfig.Field( 

218 doc="Standard Atmosphere aerosol alpha", 

219 dtype=float, 

220 default=1.0, 

221 ) 

222 airmassStd = pexConfig.Field( 

223 doc=("Standard Atmosphere airmass; " 

224 "Recommended default depends on the survey strategy."), 

225 dtype=float, 

226 default=None, 

227 ) 

228 lambdaNorm = pexConfig.Field( 

229 doc="Aerosol Optical Depth normalization wavelength (Angstrom)", 

230 dtype=float, 

231 default=7750.0, 

232 ) 

233 lambdaStep = pexConfig.Field( 

234 doc="Wavelength step for generating atmospheres (nm)", 

235 dtype=float, 

236 default=0.5, 

237 ) 

238 lambdaRange = pexConfig.ListField( 

239 doc="Wavelength range for LUT (Angstrom)", 

240 dtype=float, 

241 default=[3000.0, 11000.0], 

242 ) 

243 

244 

245class FgcmMakeLutConfig(pipeBase.PipelineTaskConfig, 

246 pipelineConnections=FgcmMakeLutConnections): 

247 """Config for FgcmMakeLutTask""" 

248 physicalFilters = pexConfig.ListField( 

249 doc="List of physicalFilter labels to generate look-up table.", 

250 dtype=str, 

251 default=[], 

252 ) 

253 stdPhysicalFilterOverrideMap = pexConfig.DictField( 

254 doc=("Override mapping from physical filter labels to 'standard' physical " 

255 "filter labels. The 'standard' physical filter defines the transmission " 

256 "curve that the FGCM standard bandpass will be based on. " 

257 "Any filter not listed here will be mapped to " 

258 "itself (e.g. g->g or HSC-G->HSC-G). Use this override for cross-" 

259 "filter calibration such as HSC-R->HSC-R2 and HSC-I->HSC-I2."), 

260 keytype=str, 

261 itemtype=str, 

262 default={}, 

263 ) 

264 atmosphereTableName = pexConfig.Field( 

265 doc="FGCM name or filename of precomputed atmospheres", 

266 dtype=str, 

267 default=None, 

268 optional=True, 

269 ) 

270 doOpticsTransmission = pexConfig.Field( 

271 doc="Include optics transmission?", 

272 dtype=bool, 

273 default=True, 

274 ) 

275 doSensorTransmission = pexConfig.Field( 

276 doc="Include sensor transmission?", 

277 dtype=bool, 

278 default=True, 

279 ) 

280 parameters = pexConfig.ConfigField( 

281 doc="Atmosphere parameters (required if no atmosphereTableName)", 

282 dtype=FgcmMakeLutParametersConfig, 

283 default=None, 

284 check=None) 

285 

286 def validate(self): 

287 """ 

288 Validate the config parameters. 

289 

290 This method behaves differently from the parent validate in the case 

291 that atmosphereTableName is set. In this case, the config values 

292 for standard values, step sizes, and ranges are loaded 

293 directly from the specified atmosphereTableName. 

294 """ 

295 # check that filterNames and stdFilterNames are okay 

296 self._fields['physicalFilters'].validate(self) 

297 self._fields['stdPhysicalFilterOverrideMap'].validate(self) 

298 

299 if self.atmosphereTableName is None: 

300 # Validate the parameters 

301 self._fields['parameters'].validate(self) 

302 

303 

304class FgcmMakeLutTask(pipeBase.PipelineTask): 

305 """ 

306 Make Look-Up Table for FGCM. 

307 

308 This task computes a look-up-table for the range in expected atmosphere 

309 variation and variation in instrumental throughput (as tracked by the 

310 transmission_filter products). By pre-computing linearized integrals, 

311 the FGCM fit is orders of magnitude faster for stars with a broad range 

312 of colors and observing bands, yielding precision at the 1-2 mmag level. 

313 

314 Computing a LUT requires running MODTRAN or with a pre-generated 

315 atmosphere table packaged with fgcm. 

316 """ 

317 

318 ConfigClass = FgcmMakeLutConfig 

319 _DefaultName = "fgcmMakeLut" 

320 

321 def __init__(self, initInputs=None, **kwargs): 

322 super().__init__(**kwargs) 

323 

324 def runQuantum(self, butlerQC, inputRefs, outputRefs): 

325 camera = butlerQC.get(inputRefs.camera) 

326 

327 if self.config.doOpticsTransmission: 

328 opticsHandle = butlerQC.get(inputRefs.transmission_optics) 

329 else: 

330 opticsHandle = None 

331 

332 if self.config.doSensorTransmission: 

333 sensorHandles = butlerQC.get(inputRefs.transmission_sensor) 

334 sensorHandleDict = {sensorHandle.dataId.byName()['detector']: sensorHandle for 

335 sensorHandle in sensorHandles} 

336 else: 

337 sensorHandleDict = {} 

338 

339 filterHandles = butlerQC.get(inputRefs.transmission_filter) 

340 filterHandleDict = {filterHandle.dataId['physical_filter']: filterHandle for 

341 filterHandle in filterHandles} 

342 

343 struct = self._fgcmMakeLut(camera, 

344 opticsHandle, 

345 sensorHandleDict, 

346 filterHandleDict) 

347 

348 butlerQC.put(struct.fgcmLookUpTable, outputRefs.fgcmLookUpTable) 

349 butlerQC.put(struct.fgcmStandardAtmosphere, outputRefs.fgcmStandardAtmosphere) 

350 

351 refDict = {passbandRef.dataId['physical_filter']: passbandRef for 

352 passbandRef in outputRefs.fgcmStandardPassbands} 

353 for physical_filter, passband in struct.fgcmStandardPassbands.items(): 

354 butlerQC.put(passband, refDict[physical_filter]) 

355 

356 def _fgcmMakeLut(self, camera, opticsHandle, sensorHandleDict, 

357 filterHandleDict): 

358 """ 

359 Make a FGCM Look-up Table 

360 

361 Parameters 

362 ---------- 

363 camera : `lsst.afw.cameraGeom.Camera` 

364 Camera from the butler. 

365 opticsHandle : `lsst.daf.butler.DeferredDatasetHandle` 

366 Reference to optics transmission curve. 

367 sensorHandleDict : `dict` of [`int`, `lsst.daf.butler.DeferredDatasetHandle`] 

368 Dictionary of references to sensor transmission curves. Key will 

369 be detector id. 

370 filterHandleDict : `dict` of [`str`, `lsst.daf.butler.DeferredDatasetHandle`] 

371 Dictionary of references to filter transmission curves. Key will 

372 be physical filter label. 

373 

374 Returns 

375 ------- 

376 retStruct : `lsst.pipe.base.Struct` 

377 Output structure with keys: 

378 

379 fgcmLookUpTable : `BaseCatalog` 

380 The FGCM look-up table. 

381 fgcmStandardAtmosphere : `lsst.afw.image.TransmissionCurve` 

382 Transmission curve for the FGCM standard atmosphere. 

383 fgcmStandardPassbands : `dict` [`str`, `lsst.afw.image.TransmissionCurve`] 

384 Dictionary of standard passbands, with the key as the 

385 physical filter name. 

386 """ 

387 # number of ccds from the length of the camera iterator 

388 nCcd = len(camera) 

389 self.log.info("Found %d ccds for look-up table" % (nCcd)) 

390 

391 # Load in optics, etc. 

392 self._loadThroughputs(camera, 

393 opticsHandle, 

394 sensorHandleDict, 

395 filterHandleDict) 

396 

397 lutConfig = self._createLutConfig(nCcd) 

398 

399 # make the lut object 

400 self.log.info("Making the LUT maker object") 

401 self.fgcmLutMaker = fgcm.FgcmLUTMaker(lutConfig) 

402 

403 # generate the throughput dictionary. 

404 

405 # these will be in Angstroms 

406 # note that lambdaStep is currently in nm, because of historical 

407 # reasons in the code. Convert to Angstroms here. 

408 throughputLambda = np.arange(self.fgcmLutMaker.lambdaRange[0], 

409 self.fgcmLutMaker.lambdaRange[1]+self.fgcmLutMaker.lambdaStep*10, 

410 self.fgcmLutMaker.lambdaStep*10.) 

411 

412 self.log.info("Built throughput lambda, %.1f-%.1f, step %.2f" % 

413 (throughputLambda[0], throughputLambda[-1], 

414 throughputLambda[1] - throughputLambda[0])) 

415 

416 throughputDict = {} 

417 for i, physicalFilter in enumerate(self.config.physicalFilters): 

418 tDict = {} 

419 tDict['LAMBDA'] = throughputLambda 

420 for ccdIndex, detector in enumerate(camera): 

421 tDict[ccdIndex] = self._getThroughputDetector(detector, physicalFilter, throughputLambda) 

422 throughputDict[physicalFilter] = tDict 

423 

424 # set the throughputs 

425 self.fgcmLutMaker.setThroughputs(throughputDict) 

426 

427 # make the LUT 

428 self.log.info("Making LUT") 

429 self.fgcmLutMaker.makeLUT() 

430 

431 # and save the LUT 

432 

433 # build the index values 

434 comma = ',' 

435 physicalFilterString = comma.join(self.config.physicalFilters) 

436 stdPhysicalFilterString = comma.join(self._getStdPhysicalFilterList()) 

437 

438 atmosphereTableName = 'NoTableWasUsed' 

439 if self.config.atmosphereTableName is not None: 

440 atmosphereTableName = self.config.atmosphereTableName 

441 

442 lutSchema = self._makeLutSchema(physicalFilterString, stdPhysicalFilterString, 

443 atmosphereTableName) 

444 

445 lutCat = self._makeLutCat(lutSchema, physicalFilterString, 

446 stdPhysicalFilterString, atmosphereTableName) 

447 

448 atmStd = TransmissionCurve.makeSpatiallyConstant( 

449 throughput=self.fgcmLutMaker.atmStdTrans.astype(np.float64), 

450 wavelengths=self.fgcmLutMaker.atmLambda.astype(np.float64), 

451 throughputAtMin=self.fgcmLutMaker.atmStdTrans[0], 

452 throughputAtMax=self.fgcmLutMaker.atmStdTrans[1], 

453 ) 

454 

455 fgcmStandardPassbands = {} 

456 for i, physical_filter in enumerate(self.fgcmLutMaker.filterNames): 

457 passband = self.fgcmLutMaker.throughputs[i]['THROUGHPUT_AVG']*self.fgcmLutMaker.atmStdTrans 

458 fgcmStandardPassbands[physical_filter] = TransmissionCurve.makeSpatiallyConstant( 

459 throughput=passband.astype(np.float64), 

460 wavelengths=self.fgcmLutMaker.atmLambda.astype(np.float64), 

461 throughputAtMin=passband[0], 

462 throughputAtMax=passband[-1], 

463 ) 

464 

465 retStruct = pipeBase.Struct( 

466 fgcmLookUpTable=lutCat, 

467 fgcmStandardAtmosphere=atmStd, 

468 fgcmStandardPassbands=fgcmStandardPassbands, 

469 ) 

470 

471 return retStruct 

472 

473 def _getStdPhysicalFilterList(self): 

474 """Get the standard physical filter lists from config.physicalFilters 

475 and config.stdPhysicalFilterOverrideMap 

476 

477 Returns 

478 ------- 

479 stdPhysicalFilters : `list` 

480 """ 

481 override = self.config.stdPhysicalFilterOverrideMap 

482 return [override.get(physicalFilter, physicalFilter) for 

483 physicalFilter in self.config.physicalFilters] 

484 

485 def _createLutConfig(self, nCcd): 

486 """ 

487 Create the fgcmLut config dictionary 

488 

489 Parameters 

490 ---------- 

491 nCcd: `int` 

492 Number of CCDs in the camera 

493 """ 

494 

495 # create the common stub of the lutConfig 

496 lutConfig = {} 

497 lutConfig['logger'] = self.log 

498 lutConfig['filterNames'] = self.config.physicalFilters 

499 lutConfig['stdFilterNames'] = self._getStdPhysicalFilterList() 

500 lutConfig['nCCD'] = nCcd 

501 

502 # atmosphereTable already validated if available 

503 if self.config.atmosphereTableName is not None: 

504 lutConfig['atmosphereTableName'] = self.config.atmosphereTableName 

505 else: 

506 # use the regular paramters (also validated if needed) 

507 lutConfig['elevation'] = self.config.parameters.elevation 

508 lutConfig['pmbRange'] = self.config.parameters.pmbRange 

509 lutConfig['pmbSteps'] = self.config.parameters.pmbSteps 

510 lutConfig['pwvRange'] = self.config.parameters.pwvRange 

511 lutConfig['pwvSteps'] = self.config.parameters.pwvSteps 

512 lutConfig['o3Range'] = self.config.parameters.o3Range 

513 lutConfig['o3Steps'] = self.config.parameters.o3Steps 

514 lutConfig['tauRange'] = self.config.parameters.tauRange 

515 lutConfig['tauSteps'] = self.config.parameters.tauSteps 

516 lutConfig['alphaRange'] = self.config.parameters.alphaRange 

517 lutConfig['alphaSteps'] = self.config.parameters.alphaSteps 

518 lutConfig['zenithRange'] = self.config.parameters.zenithRange 

519 lutConfig['zenithSteps'] = self.config.parameters.zenithSteps 

520 lutConfig['pmbStd'] = self.config.parameters.pmbStd 

521 lutConfig['pwvStd'] = self.config.parameters.pwvStd 

522 lutConfig['o3Std'] = self.config.parameters.o3Std 

523 lutConfig['tauStd'] = self.config.parameters.tauStd 

524 lutConfig['alphaStd'] = self.config.parameters.alphaStd 

525 lutConfig['airmassStd'] = self.config.parameters.airmassStd 

526 lutConfig['lambdaRange'] = self.config.parameters.lambdaRange 

527 lutConfig['lambdaStep'] = self.config.parameters.lambdaStep 

528 lutConfig['lambdaNorm'] = self.config.parameters.lambdaNorm 

529 

530 return lutConfig 

531 

532 def _loadThroughputs(self, camera, opticsHandle, sensorHandleDict, filterHandleDict): 

533 """Internal method to load throughput data for filters 

534 

535 Parameters 

536 ---------- 

537 camera: `lsst.afw.cameraGeom.Camera` 

538 Camera from the butler 

539 opticsHandle : `lsst.daf.butler.DeferredDatasetHandle` 

540 Reference to optics transmission curve. 

541 sensorHandleDict : `dict` of [`int`, `lsst.daf.butler.DeferredDatasetHandle`] 

542 Dictionary of references to sensor transmission curves. Key will 

543 be detector id. 

544 filterHandleDict : `dict` of [`str`, `lsst.daf.butler.DeferredDatasetHandle`] 

545 Dictionary of references to filter transmission curves. Key will 

546 be physical filter label. 

547 

548 Raises 

549 ------ 

550 ValueError : Raised if configured filter name does not match any of the 

551 available filter transmission curves. 

552 """ 

553 if self.config.doOpticsTransmission: 

554 self._opticsTransmission = opticsHandle.get() 

555 else: 

556 self._opticsTransmission = TransmissionCurve.makeSpatiallyConstant( 

557 throughput=np.ones(100), 

558 wavelengths=np.linspace( 

559 self.config.parameters.lambdaRange[0], 

560 self.config.parameters.lambdaRange[1], 

561 100, 

562 ), 

563 throughputAtMin=1.0, 

564 throughputAtMax=1.0, 

565 ) 

566 

567 self._sensorsTransmission = {} 

568 for detector in camera: 

569 if self.config.doSensorTransmission: 

570 self._sensorsTransmission[detector.getId()] = sensorHandleDict[detector.getId()].get() 

571 else: 

572 self._sensorsTransmission[detector.getId()] = TransmissionCurve.makeSpatiallyConstant( 

573 throughput=np.ones(100), 

574 wavelengths=np.linspace( 

575 self.config.parameters.lambdaRange[0], 

576 self.config.parameters.lambdaRange[1], 

577 100, 

578 ), 

579 throughputAtMin=1.0, 

580 throughputAtMax=1.0, 

581 ) 

582 

583 self._filtersTransmission = {} 

584 for physicalFilter in self.config.physicalFilters: 

585 self._filtersTransmission[physicalFilter] = filterHandleDict[physicalFilter].get() 

586 

587 def _getThroughputDetector(self, detector, physicalFilter, throughputLambda): 

588 """Internal method to get throughput for a detector. 

589 

590 Returns the throughput at the center of the detector for a given filter. 

591 

592 Parameters 

593 ---------- 

594 detector: `lsst.afw.cameraGeom._detector.Detector` 

595 Detector on camera 

596 physicalFilter: `str` 

597 Physical filter label 

598 throughputLambda: `np.array(dtype=np.float64)` 

599 Wavelength steps (Angstrom) 

600 

601 Returns 

602 ------- 

603 throughput: `np.array(dtype=np.float64)` 

604 Throughput (max 1.0) at throughputLambda 

605 """ 

606 

607 c = detector.getCenter(afwCameraGeom.FOCAL_PLANE) 

608 c.scale(1.0/detector.getPixelSize()[0]) # Assumes x and y pixel sizes in arcsec are the same 

609 

610 throughput = self._opticsTransmission.sampleAt(position=c, 

611 wavelengths=throughputLambda) 

612 

613 throughput *= self._sensorsTransmission[detector.getId()].sampleAt(position=c, 

614 wavelengths=throughputLambda) 

615 

616 throughput *= self._filtersTransmission[physicalFilter].sampleAt(position=c, 

617 wavelengths=throughputLambda) 

618 

619 # Clip the throughput from 0 to 1 

620 throughput = np.clip(throughput, 0.0, 1.0) 

621 

622 return throughput 

623 

624 def _makeLutSchema(self, physicalFilterString, stdPhysicalFilterString, 

625 atmosphereTableName): 

626 """ 

627 Make the LUT schema 

628 

629 Parameters 

630 ---------- 

631 physicalFilterString: `str` 

632 Combined string of all the physicalFilters 

633 stdPhysicalFilterString: `str` 

634 Combined string of all the standard physicalFilters 

635 atmosphereTableName: `str` 

636 Name of the atmosphere table used to generate LUT 

637 

638 Returns 

639 ------- 

640 lutSchema: `afwTable.schema` 

641 """ 

642 

643 lutSchema = afwTable.Schema() 

644 

645 lutSchema.addField('tablename', type=str, doc='Atmosphere table name', 

646 size=len(atmosphereTableName)) 

647 lutSchema.addField('elevation', type=float, doc="Telescope elevation used for LUT") 

648 lutSchema.addField('physicalFilters', type=str, doc='physicalFilters in LUT', 

649 size=len(physicalFilterString)) 

650 lutSchema.addField('stdPhysicalFilters', type=str, doc='Standard physicalFilters in LUT', 

651 size=len(stdPhysicalFilterString)) 

652 lutSchema.addField('pmb', type='ArrayD', doc='Barometric Pressure', 

653 size=self.fgcmLutMaker.pmb.size) 

654 lutSchema.addField('pmbFactor', type='ArrayD', doc='PMB scaling factor', 

655 size=self.fgcmLutMaker.pmb.size) 

656 lutSchema.addField('pmbElevation', type=np.float64, doc='PMB Scaling at elevation') 

657 lutSchema.addField('pwv', type='ArrayD', doc='Preciptable Water Vapor', 

658 size=self.fgcmLutMaker.pwv.size) 

659 lutSchema.addField('o3', type='ArrayD', doc='Ozone', 

660 size=self.fgcmLutMaker.o3.size) 

661 lutSchema.addField('tau', type='ArrayD', doc='Aerosol optical depth', 

662 size=self.fgcmLutMaker.tau.size) 

663 lutSchema.addField('lambdaNorm', type=np.float64, doc='AOD wavelength') 

664 lutSchema.addField('alpha', type='ArrayD', doc='Aerosol alpha', 

665 size=self.fgcmLutMaker.alpha.size) 

666 lutSchema.addField('zenith', type='ArrayD', doc='Zenith angle', 

667 size=self.fgcmLutMaker.zenith.size) 

668 lutSchema.addField('nCcd', type=np.int32, doc='Number of CCDs') 

669 

670 # and the standard values 

671 lutSchema.addField('pmbStd', type=np.float64, doc='PMB Standard') 

672 lutSchema.addField('pwvStd', type=np.float64, doc='PWV Standard') 

673 lutSchema.addField('o3Std', type=np.float64, doc='O3 Standard') 

674 lutSchema.addField('tauStd', type=np.float64, doc='Tau Standard') 

675 lutSchema.addField('alphaStd', type=np.float64, doc='Alpha Standard') 

676 lutSchema.addField('zenithStd', type=np.float64, doc='Zenith angle Standard') 

677 lutSchema.addField('lambdaRange', type='ArrayD', doc='Wavelength range', 

678 size=2) 

679 lutSchema.addField('lambdaStep', type=np.float64, doc='Wavelength step') 

680 lutSchema.addField('lambdaStd', type='ArrayD', doc='Standard Wavelength', 

681 size=len(self.fgcmLutMaker.filterNames)) 

682 lutSchema.addField('lambdaStdFilter', type='ArrayD', doc='Standard Wavelength (raw)', 

683 size=len(self.fgcmLutMaker.filterNames)) 

684 lutSchema.addField('i0Std', type='ArrayD', doc='I0 Standard', 

685 size=len(self.fgcmLutMaker.filterNames)) 

686 lutSchema.addField('i1Std', type='ArrayD', doc='I1 Standard', 

687 size=len(self.fgcmLutMaker.filterNames)) 

688 lutSchema.addField('i10Std', type='ArrayD', doc='I10 Standard', 

689 size=len(self.fgcmLutMaker.filterNames)) 

690 lutSchema.addField('i2Std', type='ArrayD', doc='I2 Standard', 

691 size=len(self.fgcmLutMaker.filterNames)) 

692 lutSchema.addField('lambdaB', type='ArrayD', doc='Wavelength for passband (no atm)', 

693 size=len(self.fgcmLutMaker.filterNames)) 

694 lutSchema.addField('atmLambda', type='ArrayD', doc='Atmosphere wavelengths (Angstrom)', 

695 size=self.fgcmLutMaker.atmLambda.size) 

696 lutSchema.addField('atmStdTrans', type='ArrayD', doc='Standard Atmosphere Throughput', 

697 size=self.fgcmLutMaker.atmStdTrans.size) 

698 

699 # and the look-up-tables 

700 lutSchema.addField('luttype', type=str, size=20, doc='Look-up table type') 

701 lutSchema.addField('lut', type='ArrayF', doc='Look-up table for luttype', 

702 size=self.fgcmLutMaker.lut['I0'].size) 

703 

704 return lutSchema 

705 

706 def _makeLutCat(self, lutSchema, physicalFilterString, stdPhysicalFilterString, 

707 atmosphereTableName): 

708 """ 

709 Make the LUT schema 

710 

711 Parameters 

712 ---------- 

713 lutSchema: `afwTable.schema` 

714 Lut catalog schema 

715 physicalFilterString: `str` 

716 Combined string of all the physicalFilters 

717 stdPhysicalFilterString: `str` 

718 Combined string of all the standard physicalFilters 

719 atmosphereTableName: `str` 

720 Name of the atmosphere table used to generate LUT 

721 

722 Returns 

723 ------- 

724 lutCat: `afwTable.BaseCatalog` 

725 Look-up table catalog for persistence. 

726 """ 

727 

728 # The somewhat strange format is to make sure that 

729 # the rows of the afwTable do not get too large 

730 # (see DM-11419) 

731 

732 lutCat = afwTable.BaseCatalog(lutSchema) 

733 lutCat.table.preallocate(14) 

734 

735 # first fill the first index 

736 rec = lutCat.addNew() 

737 

738 rec['tablename'] = atmosphereTableName 

739 rec['elevation'] = self.fgcmLutMaker.atmosphereTable.elevation 

740 rec['physicalFilters'] = physicalFilterString 

741 rec['stdPhysicalFilters'] = stdPhysicalFilterString 

742 rec['pmb'][:] = self.fgcmLutMaker.pmb 

743 rec['pmbFactor'][:] = self.fgcmLutMaker.pmbFactor 

744 rec['pmbElevation'] = self.fgcmLutMaker.pmbElevation 

745 rec['pwv'][:] = self.fgcmLutMaker.pwv 

746 rec['o3'][:] = self.fgcmLutMaker.o3 

747 rec['tau'][:] = self.fgcmLutMaker.tau 

748 rec['lambdaNorm'] = self.fgcmLutMaker.lambdaNorm 

749 rec['alpha'][:] = self.fgcmLutMaker.alpha 

750 rec['zenith'][:] = self.fgcmLutMaker.zenith 

751 rec['nCcd'] = self.fgcmLutMaker.nCCD 

752 

753 rec['pmbStd'] = self.fgcmLutMaker.pmbStd 

754 rec['pwvStd'] = self.fgcmLutMaker.pwvStd 

755 rec['o3Std'] = self.fgcmLutMaker.o3Std 

756 rec['tauStd'] = self.fgcmLutMaker.tauStd 

757 rec['alphaStd'] = self.fgcmLutMaker.alphaStd 

758 rec['zenithStd'] = self.fgcmLutMaker.zenithStd 

759 rec['lambdaRange'][:] = self.fgcmLutMaker.lambdaRange 

760 rec['lambdaStep'] = self.fgcmLutMaker.lambdaStep 

761 rec['lambdaStd'][:] = self.fgcmLutMaker.lambdaStd 

762 rec['lambdaStdFilter'][:] = self.fgcmLutMaker.lambdaStdFilter 

763 rec['i0Std'][:] = self.fgcmLutMaker.I0Std 

764 rec['i1Std'][:] = self.fgcmLutMaker.I1Std 

765 rec['i10Std'][:] = self.fgcmLutMaker.I10Std 

766 rec['i2Std'][:] = self.fgcmLutMaker.I2Std 

767 rec['lambdaB'][:] = self.fgcmLutMaker.lambdaB 

768 rec['atmLambda'][:] = self.fgcmLutMaker.atmLambda 

769 rec['atmStdTrans'][:] = self.fgcmLutMaker.atmStdTrans 

770 

771 rec['luttype'] = 'I0' 

772 rec['lut'][:] = self.fgcmLutMaker.lut['I0'].flatten() 

773 

774 # and add the rest 

775 rec = lutCat.addNew() 

776 rec['luttype'] = 'I1' 

777 rec['lut'][:] = self.fgcmLutMaker.lut['I1'].flatten() 

778 

779 derivTypes = ['D_PMB', 'D_LNPWV', 'D_O3', 'D_LNTAU', 'D_ALPHA', 'D_SECZENITH', 

780 'D_PMB_I1', 'D_LNPWV_I1', 'D_O3_I1', 'D_LNTAU_I1', 'D_ALPHA_I1', 

781 'D_SECZENITH_I1'] 

782 for derivType in derivTypes: 

783 rec = lutCat.addNew() 

784 rec['luttype'] = derivType 

785 rec['lut'][:] = self.fgcmLutMaker.lutDeriv[derivType].flatten() 

786 

787 return lutCat