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# 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 the final fgcmcal output products. 

24 

25This task takes the final output from fgcmFitCycle and produces the following 

26outputs for use in the DM stack: the FGCM standard stars in a reference 

27catalog format; the model atmospheres in "transmission_atmosphere_fgcm" 

28format; and the zeropoints in "fgcm_photoCalib" format. Optionally, the 

29task can transfer the 'absolute' calibration from a reference catalog 

30to put the fgcm standard stars in units of Jansky. This is accomplished 

31by matching stars in a sample of healpix pixels, and applying the median 

32offset per band. 

33""" 

34import sys 

35import traceback 

36import copy 

37 

38import numpy as np 

39import healpy as hp 

40import esutil 

41from astropy import units 

42 

43import lsst.pex.config as pexConfig 

44import lsst.pipe.base as pipeBase 

45from lsst.pipe.base import connectionTypes 

46from lsst.afw.image import TransmissionCurve 

47from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask 

48from lsst.meas.algorithms import ReferenceObjectLoader 

49from lsst.pipe.tasks.photoCal import PhotoCalTask 

50import lsst.geom 

51import lsst.afw.image as afwImage 

52import lsst.afw.math as afwMath 

53import lsst.afw.table as afwTable 

54from lsst.meas.algorithms import IndexerRegistry 

55from lsst.meas.algorithms import DatasetConfig 

56from lsst.meas.algorithms.ingestIndexReferenceTask import addRefCatMetadata 

57 

58from .utilities import computeApproxPixelAreaFields 

59from .utilities import lookupStaticCalibrations 

60 

61import fgcm 

62 

63__all__ = ['FgcmOutputProductsConfig', 'FgcmOutputProductsTask', 'FgcmOutputProductsRunner'] 

64 

65 

66class FgcmOutputProductsConnections(pipeBase.PipelineTaskConnections, 

67 dimensions=("instrument",), 

68 defaultTemplates={"cycleNumber": "0"}): 

69 camera = connectionTypes.PrerequisiteInput( 

70 doc="Camera instrument", 

71 name="camera", 

72 storageClass="Camera", 

73 dimensions=("instrument",), 

74 lookupFunction=lookupStaticCalibrations, 

75 isCalibration=True, 

76 ) 

77 

78 fgcmLookUpTable = connectionTypes.PrerequisiteInput( 

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

80 "chromatic corrections."), 

81 name="fgcmLookUpTable", 

82 storageClass="Catalog", 

83 dimensions=("instrument",), 

84 deferLoad=True, 

85 ) 

86 

87 fgcmVisitCatalog = connectionTypes.PrerequisiteInput( 

88 doc="Catalog of visit information for fgcm", 

89 name="fgcmVisitCatalog", 

90 storageClass="Catalog", 

91 dimensions=("instrument",), 

92 deferLoad=True, 

93 ) 

94 

95 fgcmStandardStars = connectionTypes.PrerequisiteInput( 

96 doc="Catalog of standard star data from fgcm fit", 

97 name="fgcmStandardStars{cycleNumber}", 

98 storageClass="SimpleCatalog", 

99 dimensions=("instrument",), 

100 deferLoad=True, 

101 ) 

102 

103 fgcmZeropoints = connectionTypes.PrerequisiteInput( 

104 doc="Catalog of zeropoints from fgcm fit", 

105 name="fgcmZeropoints{cycleNumber}", 

106 storageClass="Catalog", 

107 dimensions=("instrument",), 

108 deferLoad=True, 

109 ) 

110 

111 fgcmAtmosphereParameters = connectionTypes.PrerequisiteInput( 

112 doc="Catalog of atmosphere parameters from fgcm fit", 

113 name="fgcmAtmosphereParameters{cycleNumber}", 

114 storageClass="Catalog", 

115 dimensions=("instrument",), 

116 deferLoad=True, 

117 ) 

118 

119 refCat = connectionTypes.PrerequisiteInput( 

120 doc="Reference catalog to use for photometric calibration", 

121 name="cal_ref_cat", 

122 storageClass="SimpleCatalog", 

123 dimensions=("skypix",), 

124 deferLoad=True, 

125 multiple=True, 

126 ) 

127 

128 fgcmBuildStarsTableConfig = connectionTypes.PrerequisiteInput( 

129 doc="Config used to build FGCM input stars", 

130 name="fgcmBuildStarsTable_config", 

131 storageClass="Config", 

132 ) 

133 

134 fgcmPhotoCalib = connectionTypes.Output( 

135 doc="Per-visit photoCalib exposure catalogs produced from fgcm calibration", 

136 name="fgcmPhotoCalibCatalog", 

137 storageClass="ExposureCatalog", 

138 dimensions=("instrument", "visit",), 

139 multiple=True, 

140 ) 

141 

142 fgcmTransmissionAtmosphere = connectionTypes.Output( 

143 doc="Per-visit atmosphere transmission files produced from fgcm calibration", 

144 name="transmission_atmosphere_fgcm", 

145 storageClass="TransmissionCurve", 

146 dimensions=("instrument", 

147 "visit",), 

148 multiple=True, 

149 ) 

150 

151 fgcmOffsets = connectionTypes.Output( 

152 doc="Per-band offsets computed from doReferenceCalibration", 

153 name="fgcmReferenceCalibrationOffsets", 

154 storageClass="Catalog", 

155 dimensions=("instrument",), 

156 multiple=False, 

157 ) 

158 

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

160 super().__init__(config=config) 

161 

162 if str(int(config.connections.cycleNumber)) != config.connections.cycleNumber: 

163 raise ValueError("cycleNumber must be of integer format") 

164 if config.connections.refCat != config.refObjLoader.ref_dataset_name: 

165 raise ValueError("connections.refCat must be the same as refObjLoader.ref_dataset_name") 

166 

167 if config.doRefcatOutput: 

168 raise ValueError("FgcmOutputProductsTask (Gen3) does not support doRefcatOutput") 

169 

170 if not config.doReferenceCalibration: 

171 self.prerequisiteInputs.remove("refCat") 

172 if not config.doAtmosphereOutput: 

173 self.prerequisiteInputs.remove("fgcmAtmosphereParameters") 

174 if not config.doZeropointOutput: 

175 self.prerequisiteInputs.remove("fgcmZeropoints") 

176 if not config.doReferenceCalibration: 

177 self.outputs.remove("fgcmOffsets") 

178 

179 

180class FgcmOutputProductsConfig(pipeBase.PipelineTaskConfig, 

181 pipelineConnections=FgcmOutputProductsConnections): 

182 """Config for FgcmOutputProductsTask""" 

183 

184 cycleNumber = pexConfig.Field( 

185 doc="Final fit cycle from FGCM fit", 

186 dtype=int, 

187 default=None, 

188 ) 

189 

190 # The following fields refer to calibrating from a reference 

191 # catalog, but in the future this might need to be expanded 

192 doReferenceCalibration = pexConfig.Field( 

193 doc=("Transfer 'absolute' calibration from reference catalog? " 

194 "This afterburner step is unnecessary if reference stars " 

195 "were used in the full fit in FgcmFitCycleTask."), 

196 dtype=bool, 

197 default=False, 

198 ) 

199 doRefcatOutput = pexConfig.Field( 

200 doc="Output standard stars in reference catalog format", 

201 dtype=bool, 

202 default=True, 

203 ) 

204 doAtmosphereOutput = pexConfig.Field( 

205 doc="Output atmospheres in transmission_atmosphere_fgcm format", 

206 dtype=bool, 

207 default=True, 

208 ) 

209 doZeropointOutput = pexConfig.Field( 

210 doc="Output zeropoints in fgcm_photoCalib format", 

211 dtype=bool, 

212 default=True, 

213 ) 

214 doComposeWcsJacobian = pexConfig.Field( 

215 doc="Compose Jacobian of WCS with fgcm calibration for output photoCalib?", 

216 dtype=bool, 

217 default=True, 

218 ) 

219 doApplyMeanChromaticCorrection = pexConfig.Field( 

220 doc="Apply the mean chromatic correction to the zeropoints?", 

221 dtype=bool, 

222 default=True, 

223 ) 

224 refObjLoader = pexConfig.ConfigurableField( 

225 target=LoadIndexedReferenceObjectsTask, 

226 doc="reference object loader for 'absolute' photometric calibration", 

227 ) 

228 photoCal = pexConfig.ConfigurableField( 

229 target=PhotoCalTask, 

230 doc="task to perform 'absolute' calibration", 

231 ) 

232 referencePixelizationNside = pexConfig.Field( 

233 doc="Healpix nside to pixelize catalog to compare to reference catalog", 

234 dtype=int, 

235 default=64, 

236 ) 

237 referencePixelizationMinStars = pexConfig.Field( 

238 doc=("Minimum number of stars per healpix pixel to select for comparison" 

239 "to the specified reference catalog"), 

240 dtype=int, 

241 default=200, 

242 ) 

243 referenceMinMatch = pexConfig.Field( 

244 doc="Minimum number of stars matched to reference catalog to be used in statistics", 

245 dtype=int, 

246 default=50, 

247 ) 

248 referencePixelizationNPixels = pexConfig.Field( 

249 doc=("Number of healpix pixels to sample to do comparison. " 

250 "Doing too many will take a long time and not yield any more " 

251 "precise results because the final number is the median offset " 

252 "(per band) from the set of pixels."), 

253 dtype=int, 

254 default=100, 

255 ) 

256 datasetConfig = pexConfig.ConfigField( 

257 dtype=DatasetConfig, 

258 doc="Configuration for writing/reading ingested catalog", 

259 ) 

260 

261 def setDefaults(self): 

262 pexConfig.Config.setDefaults(self) 

263 

264 # In order to transfer the "absolute" calibration from a reference 

265 # catalog to the relatively calibrated FGCM standard stars (one number 

266 # per band), we use the PhotoCalTask to match stars in a sample of healpix 

267 # pixels. These basic settings ensure that only well-measured, good stars 

268 # from the source and reference catalogs are used for the matching. 

269 

270 # applyColorTerms needs to be False if doReferenceCalibration is False, 

271 # as is the new default after DM-16702 

272 self.photoCal.applyColorTerms = False 

273 self.photoCal.fluxField = 'instFlux' 

274 self.photoCal.magErrFloor = 0.003 

275 self.photoCal.match.referenceSelection.doSignalToNoise = True 

276 self.photoCal.match.referenceSelection.signalToNoise.minimum = 10.0 

277 self.photoCal.match.sourceSelection.doSignalToNoise = True 

278 self.photoCal.match.sourceSelection.signalToNoise.minimum = 10.0 

279 self.photoCal.match.sourceSelection.signalToNoise.fluxField = 'instFlux' 

280 self.photoCal.match.sourceSelection.signalToNoise.errField = 'instFluxErr' 

281 self.photoCal.match.sourceSelection.doFlags = True 

282 self.photoCal.match.sourceSelection.flags.good = [] 

283 self.photoCal.match.sourceSelection.flags.bad = ['flag_badStar'] 

284 self.photoCal.match.sourceSelection.doUnresolved = False 

285 self.datasetConfig.ref_dataset_name = 'fgcm_stars' 

286 self.datasetConfig.format_version = 1 

287 

288 def validate(self): 

289 super().validate() 

290 

291 # Force the connections to conform with cycleNumber 

292 self.connections.cycleNumber = str(self.cycleNumber) 

293 

294 

295class FgcmOutputProductsRunner(pipeBase.ButlerInitializedTaskRunner): 

296 """Subclass of TaskRunner for fgcmOutputProductsTask 

297 

298 fgcmOutputProductsTask.run() takes one argument, the butler, and 

299 does not run on any data in the repository. 

300 This runner does not use any parallelization. 

301 """ 

302 

303 @staticmethod 

304 def getTargetList(parsedCmd): 

305 """ 

306 Return a list with one element, the butler. 

307 """ 

308 return [parsedCmd.butler] 

309 

310 def __call__(self, butler): 

311 """ 

312 Parameters 

313 ---------- 

314 butler: `lsst.daf.persistence.Butler` 

315 

316 Returns 

317 ------- 

318 exitStatus: `list` with `pipeBase.Struct` 

319 exitStatus (0: success; 1: failure) 

320 if self.doReturnResults also 

321 results (`np.array` with absolute zeropoint offsets) 

322 """ 

323 task = self.TaskClass(butler=butler, config=self.config, log=self.log) 

324 

325 exitStatus = 0 

326 if self.doRaise: 

327 results = task.runDataRef(butler) 

328 else: 

329 try: 

330 results = task.runDataRef(butler) 

331 except Exception as e: 

332 exitStatus = 1 

333 task.log.fatal("Failed: %s" % e) 

334 if not isinstance(e, pipeBase.TaskError): 

335 traceback.print_exc(file=sys.stderr) 

336 

337 task.writeMetadata(butler) 

338 

339 if self.doReturnResults: 

340 # The results here are the zeropoint offsets for each band 

341 return [pipeBase.Struct(exitStatus=exitStatus, 

342 results=results)] 

343 else: 

344 return [pipeBase.Struct(exitStatus=exitStatus)] 

345 

346 def run(self, parsedCmd): 

347 """ 

348 Run the task, with no multiprocessing 

349 

350 Parameters 

351 ---------- 

352 parsedCmd: `lsst.pipe.base.ArgumentParser` parsed command line 

353 """ 

354 

355 resultList = [] 

356 

357 if self.precall(parsedCmd): 

358 targetList = self.getTargetList(parsedCmd) 

359 # make sure that we only get 1 

360 resultList = self(targetList[0]) 

361 

362 return resultList 

363 

364 

365class FgcmOutputProductsTask(pipeBase.PipelineTask, pipeBase.CmdLineTask): 

366 """ 

367 Output products from FGCM global calibration. 

368 """ 

369 

370 ConfigClass = FgcmOutputProductsConfig 

371 RunnerClass = FgcmOutputProductsRunner 

372 _DefaultName = "fgcmOutputProducts" 

373 

374 def __init__(self, butler=None, **kwargs): 

375 super().__init__(**kwargs) 

376 

377 # no saving of metadata for now 

378 def _getMetadataName(self): 

379 return None 

380 

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

382 dataRefDict = {} 

383 dataRefDict['camera'] = butlerQC.get(inputRefs.camera) 

384 dataRefDict['fgcmLookUpTable'] = butlerQC.get(inputRefs.fgcmLookUpTable) 

385 dataRefDict['fgcmVisitCatalog'] = butlerQC.get(inputRefs.fgcmVisitCatalog) 

386 dataRefDict['fgcmStandardStars'] = butlerQC.get(inputRefs.fgcmStandardStars) 

387 

388 if self.config.doZeropointOutput: 

389 dataRefDict['fgcmZeropoints'] = butlerQC.get(inputRefs.fgcmZeropoints) 

390 photoCalibRefDict = {photoCalibRef.dataId.byName()['visit']: 

391 photoCalibRef for photoCalibRef in outputRefs.fgcmPhotoCalib} 

392 

393 if self.config.doAtmosphereOutput: 

394 dataRefDict['fgcmAtmosphereParameters'] = butlerQC.get(inputRefs.fgcmAtmosphereParameters) 

395 atmRefDict = {atmRef.dataId.byName()['visit']: atmRef for 

396 atmRef in outputRefs.fgcmTransmissionAtmosphere} 

397 

398 if self.config.doReferenceCalibration: 

399 refConfig = self.config.refObjLoader 

400 self.refObjLoader = ReferenceObjectLoader(dataIds=[ref.datasetRef.dataId 

401 for ref in inputRefs.refCat], 

402 refCats=butlerQC.get(inputRefs.refCat), 

403 config=refConfig, 

404 log=self.log) 

405 else: 

406 self.refObjLoader = None 

407 

408 dataRefDict['fgcmBuildStarsTableConfig'] = butlerQC.get(inputRefs.fgcmBuildStarsTableConfig) 

409 

410 fgcmBuildStarsConfig = butlerQC.get(inputRefs.fgcmBuildStarsTableConfig) 

411 filterMap = fgcmBuildStarsConfig.filterMap 

412 

413 if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian: 

414 raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied " 

415 "in fgcmBuildStarsTask.") 

416 if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian: 

417 self.log.warn("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.") 

418 

419 struct = self.run(dataRefDict, filterMap, returnCatalogs=True) 

420 

421 # Output the photoCalib exposure catalogs 

422 if struct.photoCalibCatalogs is not None: 

423 self.log.info("Outputting photoCalib catalogs.") 

424 for visit, expCatalog in struct.photoCalibCatalogs: 

425 butlerQC.put(expCatalog, photoCalibRefDict[visit]) 

426 self.log.info("Done outputting photoCalib catalogs.") 

427 

428 # Output the atmospheres 

429 if struct.atmospheres is not None: 

430 self.log.info("Outputting atmosphere transmission files.") 

431 for visit, atm in struct.atmospheres: 

432 butlerQC.put(atm, atmRefDict[visit]) 

433 self.log.info("Done outputting atmosphere files.") 

434 

435 if self.config.doReferenceCalibration: 

436 # Turn offset into simple catalog for persistence if necessary 

437 schema = afwTable.Schema() 

438 schema.addField('offset', type=np.float64, 

439 doc="Post-process calibration offset (mag)") 

440 offsetCat = afwTable.BaseCatalog(schema) 

441 offsetCat.resize(len(struct.offsets)) 

442 offsetCat['offset'][:] = struct.offsets 

443 

444 butlerQC.put(offsetCat, outputRefs.fgcmOffsets) 

445 

446 return 

447 

448 @pipeBase.timeMethod 

449 def runDataRef(self, butler): 

450 """ 

451 Make FGCM output products for use in the stack 

452 

453 Parameters 

454 ---------- 

455 butler: `lsst.daf.persistence.Butler` 

456 cycleNumber: `int` 

457 Final fit cycle number, override config. 

458 

459 Returns 

460 ------- 

461 offsets: `lsst.pipe.base.Struct` 

462 A structure with array of zeropoint offsets 

463 

464 Raises 

465 ------ 

466 RuntimeError: 

467 Raised if any one of the following is true: 

468 

469 - butler cannot find "fgcmBuildStars_config" or 

470 "fgcmBuildStarsTable_config". 

471 - butler cannot find "fgcmFitCycle_config". 

472 - "fgcmFitCycle_config" does not refer to 

473 `self.config.cycleNumber`. 

474 - butler cannot find "fgcmAtmosphereParameters" and 

475 `self.config.doAtmosphereOutput` is `True`. 

476 - butler cannot find "fgcmStandardStars" and 

477 `self.config.doReferenceCalibration` is `True` or 

478 `self.config.doRefcatOutput` is `True`. 

479 - butler cannot find "fgcmZeropoints" and 

480 `self.config.doZeropointOutput` is `True`. 

481 """ 

482 if self.config.doReferenceCalibration: 

483 # We need the ref obj loader to get the flux field 

484 self.makeSubtask("refObjLoader", butler=butler) 

485 

486 # Check to make sure that the fgcmBuildStars config exists, to retrieve 

487 # the visit and ccd dataset tags 

488 if not butler.datasetExists('fgcmBuildStarsTable_config') and \ 

489 not butler.datasetExists('fgcmBuildStars_config'): 

490 raise RuntimeError("Cannot find fgcmBuildStarsTable_config or fgcmBuildStars_config, " 

491 "which is prereq for fgcmOutputProducts") 

492 

493 if butler.datasetExists('fgcmBuildStarsTable_config'): 

494 fgcmBuildStarsConfig = butler.get('fgcmBuildStarsTable_config') 

495 else: 

496 fgcmBuildStarsConfig = butler.get('fgcmBuildStars_config') 

497 visitDataRefName = fgcmBuildStarsConfig.visitDataRefName 

498 ccdDataRefName = fgcmBuildStarsConfig.ccdDataRefName 

499 filterMap = fgcmBuildStarsConfig.filterMap 

500 

501 if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian: 

502 raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied " 

503 "in fgcmBuildStarsTask.") 

504 

505 if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian: 

506 self.log.warn("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.") 

507 

508 # And make sure that the atmosphere was output properly 

509 if (self.config.doAtmosphereOutput 

510 and not butler.datasetExists('fgcmAtmosphereParameters', fgcmcycle=self.config.cycleNumber)): 

511 raise RuntimeError(f"Atmosphere parameters are missing for cycle {self.config.cycleNumber}.") 

512 

513 if not butler.datasetExists('fgcmStandardStars', 

514 fgcmcycle=self.config.cycleNumber): 

515 raise RuntimeError("Standard stars are missing for cycle %d." % 

516 (self.config.cycleNumber)) 

517 

518 if (self.config.doZeropointOutput 

519 and (not butler.datasetExists('fgcmZeropoints', fgcmcycle=self.config.cycleNumber))): 

520 raise RuntimeError("Zeropoints are missing for cycle %d." % 

521 (self.config.cycleNumber)) 

522 

523 dataRefDict = {} 

524 # This is the _actual_ camera 

525 dataRefDict['camera'] = butler.get('camera') 

526 dataRefDict['fgcmLookUpTable'] = butler.dataRef('fgcmLookUpTable') 

527 dataRefDict['fgcmVisitCatalog'] = butler.dataRef('fgcmVisitCatalog') 

528 dataRefDict['fgcmStandardStars'] = butler.dataRef('fgcmStandardStars', 

529 fgcmcycle=self.config.cycleNumber) 

530 

531 if self.config.doZeropointOutput: 

532 dataRefDict['fgcmZeropoints'] = butler.dataRef('fgcmZeropoints', 

533 fgcmcycle=self.config.cycleNumber) 

534 if self.config.doAtmosphereOutput: 

535 dataRefDict['fgcmAtmosphereParameters'] = butler.dataRef('fgcmAtmosphereParameters', 

536 fgcmcycle=self.config.cycleNumber) 

537 

538 struct = self.run(dataRefDict, filterMap, butler=butler, returnCatalogs=False) 

539 

540 if struct.photoCalibs is not None: 

541 self.log.info("Outputting photoCalib files.") 

542 

543 filterMapping = {} 

544 for visit, detector, filtername, photoCalib in struct.photoCalibs: 

545 if filtername not in filterMapping: 

546 # We need to find the mapping from encoded filter to dataid filter, 

547 # and this trick allows us to do that. 

548 dataId = {visitDataRefName: visit, 

549 ccdDataRefName: detector} 

550 dataRef = butler.dataRef('raw', dataId=dataId) 

551 filterMapping[filtername] = dataRef.dataId['filter'] 

552 

553 butler.put(photoCalib, 'fgcm_photoCalib', 

554 dataId={visitDataRefName: visit, 

555 ccdDataRefName: detector, 

556 'filter': filterMapping[filtername]}) 

557 

558 self.log.info("Done outputting photoCalib files.") 

559 

560 if struct.atmospheres is not None: 

561 self.log.info("Outputting atmosphere transmission files.") 

562 for visit, atm in struct.atmospheres: 

563 butler.put(atm, "transmission_atmosphere_fgcm", 

564 dataId={visitDataRefName: visit}) 

565 self.log.info("Done outputting atmosphere transmissions.") 

566 

567 return pipeBase.Struct(offsets=struct.offsets) 

568 

569 def run(self, dataRefDict, filterMap, returnCatalogs=True, butler=None): 

570 """Run the output products task. 

571 

572 Parameters 

573 ---------- 

574 dataRefDict : `dict` 

575 All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or 

576 `lsst.daf.butler.DeferredDatasetHandle` (gen3) 

577 dataRef dictionary with keys: 

578 

579 ``"camera"`` 

580 Camera object (`lsst.afw.cameraGeom.Camera`) 

581 ``"fgcmLookUpTable"`` 

582 dataRef for the FGCM look-up table. 

583 ``"fgcmVisitCatalog"`` 

584 dataRef for visit summary catalog. 

585 ``"fgcmStandardStars"`` 

586 dataRef for the output standard star catalog. 

587 ``"fgcmZeropoints"`` 

588 dataRef for the zeropoint data catalog. 

589 ``"fgcmAtmosphereParameters"`` 

590 dataRef for the atmosphere parameter catalog. 

591 ``"fgcmBuildStarsTableConfig"`` 

592 Config for `lsst.fgcmcal.fgcmBuildStarsTableTask`. 

593 filterMap : `dict` 

594 Dictionary of mappings from filter to FGCM band. 

595 returnCatalogs : `bool`, optional 

596 Return photoCalibs as per-visit exposure catalogs. 

597 butler : `lsst.daf.persistence.Butler`, optional 

598 Gen2 butler used for reference star outputs 

599 

600 Returns 

601 ------- 

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

603 Output structure with keys: 

604 

605 offsets : `np.ndarray` 

606 Final reference offsets, per band. 

607 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)] 

608 Generator that returns (visit, transmissionCurve) tuples. 

609 photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)] 

610 Generator that returns (visit, ccd, filtername, photoCalib) tuples. 

611 (returned if returnCatalogs is False). 

612 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)] 

613 Generator that returns (visit, exposureCatalog) tuples. 

614 (returned if returnCatalogs is True). 

615 """ 

616 stdCat = dataRefDict['fgcmStandardStars'].get() 

617 md = stdCat.getMetadata() 

618 bands = md.getArray('BANDS') 

619 

620 if self.config.doReferenceCalibration: 

621 offsets = self._computeReferenceOffsets(stdCat, bands) 

622 else: 

623 offsets = np.zeros(len(bands)) 

624 

625 # This is Gen2 only, and requires the butler. 

626 if self.config.doRefcatOutput and butler is not None: 

627 self._outputStandardStars(butler, stdCat, offsets, bands, self.config.datasetConfig) 

628 

629 del stdCat 

630 

631 if self.config.doZeropointOutput: 

632 zptCat = dataRefDict['fgcmZeropoints'].get() 

633 visitCat = dataRefDict['fgcmVisitCatalog'].get() 

634 

635 pcgen = self._outputZeropoints(dataRefDict['camera'], zptCat, visitCat, offsets, bands, 

636 filterMap, returnCatalogs=returnCatalogs) 

637 else: 

638 pcgen = None 

639 

640 if self.config.doAtmosphereOutput: 

641 atmCat = dataRefDict['fgcmAtmosphereParameters'].get() 

642 atmgen = self._outputAtmospheres(dataRefDict, atmCat) 

643 else: 

644 atmgen = None 

645 

646 retStruct = pipeBase.Struct(offsets=offsets, 

647 atmospheres=atmgen) 

648 if returnCatalogs: 

649 retStruct.photoCalibCatalogs = pcgen 

650 else: 

651 retStruct.photoCalibs = pcgen 

652 

653 return retStruct 

654 

655 def generateTractOutputProducts(self, dataRefDict, tract, 

656 visitCat, zptCat, atmCat, stdCat, 

657 fgcmBuildStarsConfig, 

658 returnCatalogs=True, 

659 butler=None): 

660 """ 

661 Generate the output products for a given tract, as specified in the config. 

662 

663 This method is here to have an alternate entry-point for 

664 FgcmCalibrateTract. 

665 

666 Parameters 

667 ---------- 

668 dataRefDict : `dict` 

669 All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or 

670 `lsst.daf.butler.DeferredDatasetHandle` (gen3) 

671 dataRef dictionary with keys: 

672 

673 ``"camera"`` 

674 Camera object (`lsst.afw.cameraGeom.Camera`) 

675 ``"fgcmLookUpTable"`` 

676 dataRef for the FGCM look-up table. 

677 tract : `int` 

678 Tract number 

679 visitCat : `lsst.afw.table.BaseCatalog` 

680 FGCM visitCat from `FgcmBuildStarsTask` 

681 zptCat : `lsst.afw.table.BaseCatalog` 

682 FGCM zeropoint catalog from `FgcmFitCycleTask` 

683 atmCat : `lsst.afw.table.BaseCatalog` 

684 FGCM atmosphere parameter catalog from `FgcmFitCycleTask` 

685 stdCat : `lsst.afw.table.SimpleCatalog` 

686 FGCM standard star catalog from `FgcmFitCycleTask` 

687 fgcmBuildStarsConfig : `lsst.fgcmcal.FgcmBuildStarsConfig` 

688 Configuration object from `FgcmBuildStarsTask` 

689 returnCatalogs : `bool`, optional 

690 Return photoCalibs as per-visit exposure catalogs. 

691 butler: `lsst.daf.persistence.Butler`, optional 

692 Gen2 butler used for reference star outputs 

693 

694 Returns 

695 ------- 

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

697 Output structure with keys: 

698 

699 offsets : `np.ndarray` 

700 Final reference offsets, per band. 

701 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)] 

702 Generator that returns (visit, transmissionCurve) tuples. 

703 photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)] 

704 Generator that returns (visit, ccd, filtername, photoCalib) tuples. 

705 (returned if returnCatalogs is False). 

706 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)] 

707 Generator that returns (visit, exposureCatalog) tuples. 

708 (returned if returnCatalogs is True). 

709 """ 

710 filterMap = fgcmBuildStarsConfig.filterMap 

711 

712 md = stdCat.getMetadata() 

713 bands = md.getArray('BANDS') 

714 

715 if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian: 

716 raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied " 

717 "in fgcmBuildStarsTask.") 

718 

719 if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian: 

720 self.log.warn("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.") 

721 

722 if self.config.doReferenceCalibration: 

723 offsets = self._computeReferenceOffsets(stdCat, bands) 

724 else: 

725 offsets = np.zeros(len(bands)) 

726 

727 if self.config.doRefcatOutput and butler is not None: 

728 # Create a special config that has the tract number in it 

729 datasetConfig = copy.copy(self.config.datasetConfig) 

730 datasetConfig.ref_dataset_name = '%s_%d' % (self.config.datasetConfig.ref_dataset_name, 

731 tract) 

732 self._outputStandardStars(butler, stdCat, offsets, bands, datasetConfig) 

733 

734 if self.config.doZeropointOutput: 

735 pcgen = self._outputZeropoints(dataRefDict['camera'], zptCat, visitCat, offsets, bands, 

736 filterMap, returnCatalogs=returnCatalogs) 

737 else: 

738 pcgen = None 

739 

740 if self.config.doAtmosphereOutput: 

741 atmgen = self._outputAtmospheres(dataRefDict, atmCat) 

742 else: 

743 atmgen = None 

744 

745 retStruct = pipeBase.Struct(offsets=offsets, 

746 atmospheres=atmgen) 

747 if returnCatalogs: 

748 retStruct.photoCalibCatalogs = pcgen 

749 else: 

750 retStruct.photoCalibs = pcgen 

751 

752 return retStruct 

753 

754 def _computeReferenceOffsets(self, stdCat, bands): 

755 """ 

756 Compute offsets relative to a reference catalog. 

757 

758 This method splits the star catalog into healpix pixels 

759 and computes the calibration transfer for a sample of 

760 these pixels to approximate the 'absolute' calibration 

761 values (on for each band) to apply to transfer the 

762 absolute scale. 

763 

764 Parameters 

765 ---------- 

766 stdCat : `lsst.afw.table.SimpleCatalog` 

767 FGCM standard stars 

768 bands : `list` [`str`] 

769 List of band names from FGCM output 

770 Returns 

771 ------- 

772 offsets : `numpy.array` of floats 

773 Per band zeropoint offsets 

774 """ 

775 

776 # Only use stars that are observed in all the bands that were actually used 

777 # This will ensure that we use the same healpix pixels for the absolute 

778 # calibration of each band. 

779 minObs = stdCat['ngood'].min(axis=1) 

780 

781 goodStars = (minObs >= 1) 

782 stdCat = stdCat[goodStars] 

783 

784 self.log.info("Found %d stars with at least 1 good observation in each band" % 

785 (len(stdCat))) 

786 

787 # We have to make a table for each pixel with flux/fluxErr 

788 # This is a temporary table generated for input to the photoCal task. 

789 # These fluxes are not instFlux (they are top-of-the-atmosphere approximate and 

790 # have had chromatic corrections applied to get to the standard system 

791 # specified by the atmosphere/instrumental parameters), nor are they 

792 # in Jansky (since they don't have a proper absolute calibration: the overall 

793 # zeropoint is estimated from the telescope size, etc.) 

794 sourceMapper = afwTable.SchemaMapper(stdCat.schema) 

795 sourceMapper.addMinimalSchema(afwTable.SimpleTable.makeMinimalSchema()) 

796 sourceMapper.editOutputSchema().addField('instFlux', type=np.float64, 

797 doc="instrumental flux (counts)") 

798 sourceMapper.editOutputSchema().addField('instFluxErr', type=np.float64, 

799 doc="instrumental flux error (counts)") 

800 badStarKey = sourceMapper.editOutputSchema().addField('flag_badStar', 

801 type='Flag', 

802 doc="bad flag") 

803 

804 # Split up the stars 

805 # Note that there is an assumption here that the ra/dec coords stored 

806 # on-disk are in radians, and therefore that starObs['coord_ra'] / 

807 # starObs['coord_dec'] return radians when used as an array of numpy float64s. 

808 theta = np.pi/2. - stdCat['coord_dec'] 

809 phi = stdCat['coord_ra'] 

810 

811 ipring = hp.ang2pix(self.config.referencePixelizationNside, theta, phi) 

812 h, rev = esutil.stat.histogram(ipring, rev=True) 

813 

814 gdpix, = np.where(h >= self.config.referencePixelizationMinStars) 

815 

816 self.log.info("Found %d pixels (nside=%d) with at least %d good stars" % 

817 (gdpix.size, 

818 self.config.referencePixelizationNside, 

819 self.config.referencePixelizationMinStars)) 

820 

821 if gdpix.size < self.config.referencePixelizationNPixels: 

822 self.log.warn("Found fewer good pixels (%d) than preferred in configuration (%d)" % 

823 (gdpix.size, self.config.referencePixelizationNPixels)) 

824 else: 

825 # Sample out the pixels we want to use 

826 gdpix = np.random.choice(gdpix, size=self.config.referencePixelizationNPixels, replace=False) 

827 

828 results = np.zeros(gdpix.size, dtype=[('hpix', 'i4'), 

829 ('nstar', 'i4', len(bands)), 

830 ('nmatch', 'i4', len(bands)), 

831 ('zp', 'f4', len(bands)), 

832 ('zpErr', 'f4', len(bands))]) 

833 results['hpix'] = ipring[rev[rev[gdpix]]] 

834 

835 # We need a boolean index to deal with catalogs... 

836 selected = np.zeros(len(stdCat), dtype=np.bool) 

837 

838 refFluxFields = [None]*len(bands) 

839 

840 for p, pix in enumerate(gdpix): 

841 i1a = rev[rev[pix]: rev[pix + 1]] 

842 

843 # the stdCat afwTable can only be indexed with boolean arrays, 

844 # and not numpy index arrays (see DM-16497). This little trick 

845 # converts the index array into a boolean array 

846 selected[:] = False 

847 selected[i1a] = True 

848 

849 for b, band in enumerate(bands): 

850 

851 struct = self._computeOffsetOneBand(sourceMapper, badStarKey, b, band, stdCat, 

852 selected, refFluxFields) 

853 results['nstar'][p, b] = len(i1a) 

854 results['nmatch'][p, b] = len(struct.arrays.refMag) 

855 results['zp'][p, b] = struct.zp 

856 results['zpErr'][p, b] = struct.sigma 

857 

858 # And compute the summary statistics 

859 offsets = np.zeros(len(bands)) 

860 

861 for b, band in enumerate(bands): 

862 # make configurable 

863 ok, = np.where(results['nmatch'][:, b] >= self.config.referenceMinMatch) 

864 offsets[b] = np.median(results['zp'][ok, b]) 

865 # use median absolute deviation to estimate Normal sigma 

866 # see https://en.wikipedia.org/wiki/Median_absolute_deviation 

867 madSigma = 1.4826*np.median(np.abs(results['zp'][ok, b] - offsets[b])) 

868 self.log.info("Reference catalog offset for %s band: %.12f +/- %.12f" % 

869 (band, offsets[b], madSigma)) 

870 

871 return offsets 

872 

873 def _computeOffsetOneBand(self, sourceMapper, badStarKey, 

874 b, band, stdCat, selected, refFluxFields): 

875 """ 

876 Compute the zeropoint offset between the fgcm stdCat and the reference 

877 stars for one pixel in one band 

878 

879 Parameters 

880 ---------- 

881 sourceMapper: `lsst.afw.table.SchemaMapper` 

882 Mapper to go from stdCat to calibratable catalog 

883 badStarKey: `lsst.afw.table.Key` 

884 Key for the field with bad stars 

885 b: `int` 

886 Index of the band in the star catalog 

887 band: `str` 

888 Name of band for reference catalog 

889 stdCat: `lsst.afw.table.SimpleCatalog` 

890 FGCM standard stars 

891 selected: `numpy.array(dtype=np.bool)` 

892 Boolean array of which stars are in the pixel 

893 refFluxFields: `list` 

894 List of names of flux fields for reference catalog 

895 """ 

896 

897 sourceCat = afwTable.SimpleCatalog(sourceMapper.getOutputSchema()) 

898 sourceCat.reserve(selected.sum()) 

899 sourceCat.extend(stdCat[selected], mapper=sourceMapper) 

900 sourceCat['instFlux'] = 10.**(stdCat['mag_std_noabs'][selected, b]/(-2.5)) 

901 sourceCat['instFluxErr'] = (np.log(10.)/2.5)*(stdCat['magErr_std'][selected, b] 

902 * sourceCat['instFlux']) 

903 # Make sure we only use stars that have valid measurements 

904 # (This is perhaps redundant with requirements above that the 

905 # stars be observed in all bands, but it can't hurt) 

906 badStar = (stdCat['mag_std_noabs'][selected, b] > 90.0) 

907 for rec in sourceCat[badStar]: 

908 rec.set(badStarKey, True) 

909 

910 exposure = afwImage.ExposureF() 

911 exposure.setFilter(afwImage.Filter(band)) 

912 

913 if refFluxFields[b] is None: 

914 # Need to find the flux field in the reference catalog 

915 # to work around limitations of DirectMatch in PhotoCal 

916 ctr = stdCat[0].getCoord() 

917 rad = 0.05*lsst.geom.degrees 

918 refDataTest = self.refObjLoader.loadSkyCircle(ctr, rad, band) 

919 refFluxFields[b] = refDataTest.fluxField 

920 

921 # Make a copy of the config so that we can modify it 

922 calConfig = copy.copy(self.config.photoCal.value) 

923 calConfig.match.referenceSelection.signalToNoise.fluxField = refFluxFields[b] 

924 calConfig.match.referenceSelection.signalToNoise.errField = refFluxFields[b] + 'Err' 

925 calTask = self.config.photoCal.target(refObjLoader=self.refObjLoader, 

926 config=calConfig, 

927 schema=sourceCat.getSchema()) 

928 

929 struct = calTask.run(exposure, sourceCat) 

930 

931 return struct 

932 

933 def _outputStandardStars(self, butler, stdCat, offsets, bands, datasetConfig): 

934 """ 

935 Output standard stars in indexed reference catalog format. 

936 This is not currently supported in Gen3. 

937 

938 Parameters 

939 ---------- 

940 butler : `lsst.daf.persistence.Butler` 

941 stdCat : `lsst.afw.table.SimpleCatalog` 

942 FGCM standard star catalog from fgcmFitCycleTask 

943 offsets : `numpy.array` of floats 

944 Per band zeropoint offsets 

945 bands : `list` [`str`] 

946 List of band names from FGCM output 

947 datasetConfig : `lsst.meas.algorithms.DatasetConfig` 

948 Config for reference dataset 

949 """ 

950 

951 self.log.info("Outputting standard stars to %s" % (datasetConfig.ref_dataset_name)) 

952 

953 indexer = IndexerRegistry[self.config.datasetConfig.indexer.name]( 

954 self.config.datasetConfig.indexer.active) 

955 

956 # We determine the conversion from the native units (typically radians) to 

957 # degrees for the first star. This allows us to treat coord_ra/coord_dec as 

958 # numpy arrays rather than Angles, which would we approximately 600x slower. 

959 # TODO: Fix this after DM-16524 (HtmIndexer.indexPoints should take coords 

960 # (as Angles) for input 

961 conv = stdCat[0]['coord_ra'].asDegrees()/float(stdCat[0]['coord_ra']) 

962 indices = np.array(indexer.indexPoints(stdCat['coord_ra']*conv, 

963 stdCat['coord_dec']*conv)) 

964 

965 formattedCat = self._formatCatalog(stdCat, offsets, bands) 

966 

967 # Write the master schema 

968 dataId = indexer.makeDataId('master_schema', 

969 datasetConfig.ref_dataset_name) 

970 masterCat = afwTable.SimpleCatalog(formattedCat.schema) 

971 addRefCatMetadata(masterCat) 

972 butler.put(masterCat, 'ref_cat', dataId=dataId) 

973 

974 # Break up the pixels using a histogram 

975 h, rev = esutil.stat.histogram(indices, rev=True) 

976 gd, = np.where(h > 0) 

977 selected = np.zeros(len(formattedCat), dtype=np.bool) 

978 for i in gd: 

979 i1a = rev[rev[i]: rev[i + 1]] 

980 

981 # the formattedCat afwTable can only be indexed with boolean arrays, 

982 # and not numpy index arrays (see DM-16497). This little trick 

983 # converts the index array into a boolean array 

984 selected[:] = False 

985 selected[i1a] = True 

986 

987 # Write the individual pixel 

988 dataId = indexer.makeDataId(indices[i1a[0]], 

989 datasetConfig.ref_dataset_name) 

990 butler.put(formattedCat[selected], 'ref_cat', dataId=dataId) 

991 

992 # And save the dataset configuration 

993 dataId = indexer.makeDataId(None, datasetConfig.ref_dataset_name) 

994 butler.put(datasetConfig, 'ref_cat_config', dataId=dataId) 

995 

996 self.log.info("Done outputting standard stars.") 

997 

998 def _formatCatalog(self, fgcmStarCat, offsets, bands): 

999 """ 

1000 Turn an FGCM-formatted star catalog, applying zeropoint offsets. 

1001 

1002 Parameters 

1003 ---------- 

1004 fgcmStarCat : `lsst.afw.Table.SimpleCatalog` 

1005 SimpleCatalog as output by fgcmcal 

1006 offsets : `list` with len(self.bands) entries 

1007 Zeropoint offsets to apply 

1008 bands : `list` [`str`] 

1009 List of band names from FGCM output 

1010 

1011 Returns 

1012 ------- 

1013 formattedCat: `lsst.afw.table.SimpleCatalog` 

1014 SimpleCatalog suitable for using as a reference catalog 

1015 """ 

1016 

1017 sourceMapper = afwTable.SchemaMapper(fgcmStarCat.schema) 

1018 minSchema = LoadIndexedReferenceObjectsTask.makeMinimalSchema(bands, 

1019 addCentroid=False, 

1020 addIsResolved=True, 

1021 coordErrDim=0) 

1022 sourceMapper.addMinimalSchema(minSchema) 

1023 for band in bands: 

1024 sourceMapper.editOutputSchema().addField('%s_nGood' % (band), type=np.int32) 

1025 sourceMapper.editOutputSchema().addField('%s_nTotal' % (band), type=np.int32) 

1026 sourceMapper.editOutputSchema().addField('%s_nPsfCandidate' % (band), type=np.int32) 

1027 

1028 formattedCat = afwTable.SimpleCatalog(sourceMapper.getOutputSchema()) 

1029 formattedCat.reserve(len(fgcmStarCat)) 

1030 formattedCat.extend(fgcmStarCat, mapper=sourceMapper) 

1031 

1032 # Note that we don't have to set `resolved` because the default is False 

1033 

1034 for b, band in enumerate(bands): 

1035 mag = fgcmStarCat['mag_std_noabs'][:, b].astype(np.float64) + offsets[b] 

1036 # We want fluxes in nJy from calibrated AB magnitudes 

1037 # (after applying offset). Updated after RFC-549 and RFC-575. 

1038 flux = (mag*units.ABmag).to_value(units.nJy) 

1039 fluxErr = (np.log(10.)/2.5)*flux*fgcmStarCat['magErr_std'][:, b].astype(np.float64) 

1040 

1041 formattedCat['%s_flux' % (band)][:] = flux 

1042 formattedCat['%s_fluxErr' % (band)][:] = fluxErr 

1043 formattedCat['%s_nGood' % (band)][:] = fgcmStarCat['ngood'][:, b] 

1044 formattedCat['%s_nTotal' % (band)][:] = fgcmStarCat['ntotal'][:, b] 

1045 formattedCat['%s_nPsfCandidate' % (band)][:] = fgcmStarCat['npsfcand'][:, b] 

1046 

1047 addRefCatMetadata(formattedCat) 

1048 

1049 return formattedCat 

1050 

1051 def _outputZeropoints(self, camera, zptCat, visitCat, offsets, bands, 

1052 filterMap, returnCatalogs=True, 

1053 tract=None): 

1054 """Output the zeropoints in fgcm_photoCalib format. 

1055 

1056 Parameters 

1057 ---------- 

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

1059 Camera from the butler. 

1060 zptCat : `lsst.afw.table.BaseCatalog` 

1061 FGCM zeropoint catalog from `FgcmFitCycleTask`. 

1062 visitCat : `lsst.afw.table.BaseCatalog` 

1063 FGCM visitCat from `FgcmBuildStarsTask`. 

1064 offsets : `numpy.array` 

1065 Float array of absolute calibration offsets, one for each filter. 

1066 bands : `list` [`str`] 

1067 List of band names from FGCM output. 

1068 filterMap : `dict` 

1069 Dictionary of mappings from filter to FGCM band. 

1070 returnCatalogs : `bool`, optional 

1071 Return photoCalibs as per-visit exposure catalogs. 

1072 tract: `int`, optional 

1073 Tract number to output. Default is None (global calibration) 

1074 

1075 Returns 

1076 ------- 

1077 photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)] 

1078 Generator that returns (visit, ccd, filtername, photoCalib) tuples. 

1079 (returned if returnCatalogs is False). 

1080 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)] 

1081 Generator that returns (visit, exposureCatalog) tuples. 

1082 (returned if returnCatalogs is True). 

1083 """ 

1084 # Select visit/ccds where we have a calibration 

1085 # This includes ccds where we were able to interpolate from neighboring 

1086 # ccds. 

1087 cannot_compute = fgcm.fgcmUtilities.zpFlagDict['CANNOT_COMPUTE_ZEROPOINT'] 

1088 selected = (((zptCat['fgcmFlag'] & cannot_compute) == 0) 

1089 & (zptCat['fgcmZptVar'] > 0.0)) 

1090 

1091 # Log warnings for any visit which has no valid zeropoints 

1092 badVisits = np.unique(zptCat['visit'][~selected]) 

1093 goodVisits = np.unique(zptCat['visit'][selected]) 

1094 allBadVisits = badVisits[~np.isin(badVisits, goodVisits)] 

1095 for allBadVisit in allBadVisits: 

1096 self.log.warn(f'No suitable photoCalib for visit {allBadVisit}') 

1097 

1098 # Get a mapping from filtername to the offsets 

1099 offsetMapping = {} 

1100 for f in filterMap: 

1101 # Not every filter in the map will necesarily have a band. 

1102 if filterMap[f] in bands: 

1103 offsetMapping[f] = offsets[bands.index(filterMap[f])] 

1104 

1105 # Get a mapping from "ccd" to the ccd index used for the scaling 

1106 ccdMapping = {} 

1107 for ccdIndex, detector in enumerate(camera): 

1108 ccdMapping[detector.getId()] = ccdIndex 

1109 

1110 # And a mapping to get the flat-field scaling values 

1111 scalingMapping = {} 

1112 for rec in visitCat: 

1113 scalingMapping[rec['visit']] = rec['scaling'] 

1114 

1115 if self.config.doComposeWcsJacobian: 

1116 approxPixelAreaFields = computeApproxPixelAreaFields(camera) 

1117 

1118 # The zptCat is sorted by visit, which is useful 

1119 lastVisit = -1 

1120 zptCounter = 0 

1121 zptVisitCatalog = None 

1122 for rec in zptCat[selected]: 

1123 

1124 # Retrieve overall scaling 

1125 scaling = scalingMapping[rec['visit']][ccdMapping[rec['detector']]] 

1126 

1127 # The postCalibrationOffset describe any zeropoint offsets 

1128 # to apply after the fgcm calibration. The first part comes 

1129 # from the reference catalog match (used in testing). The 

1130 # second part comes from the mean chromatic correction 

1131 # (if configured). 

1132 postCalibrationOffset = offsetMapping[rec['filtername']] 

1133 if self.config.doApplyMeanChromaticCorrection: 

1134 postCalibrationOffset += rec['fgcmDeltaChrom'] 

1135 

1136 fgcmSuperStarField = self._getChebyshevBoundedField(rec['fgcmfZptSstarCheb'], 

1137 rec['fgcmfZptChebXyMax']) 

1138 # Convert from FGCM AB to nJy 

1139 fgcmZptField = self._getChebyshevBoundedField((rec['fgcmfZptCheb']*units.AB).to_value(units.nJy), 

1140 rec['fgcmfZptChebXyMax'], 

1141 offset=postCalibrationOffset, 

1142 scaling=scaling) 

1143 

1144 if self.config.doComposeWcsJacobian: 

1145 

1146 fgcmField = afwMath.ProductBoundedField([approxPixelAreaFields[rec['detector']], 

1147 fgcmSuperStarField, 

1148 fgcmZptField]) 

1149 else: 

1150 # The photoCalib is just the product of the fgcmSuperStarField and the 

1151 # fgcmZptField 

1152 fgcmField = afwMath.ProductBoundedField([fgcmSuperStarField, fgcmZptField]) 

1153 

1154 # The "mean" calibration will be set to the center of the ccd for reference 

1155 calibCenter = fgcmField.evaluate(fgcmField.getBBox().getCenter()) 

1156 calibErr = (np.log(10.0)/2.5)*calibCenter*np.sqrt(rec['fgcmZptVar']) 

1157 photoCalib = afwImage.PhotoCalib(calibrationMean=calibCenter, 

1158 calibrationErr=calibErr, 

1159 calibration=fgcmField, 

1160 isConstant=False) 

1161 

1162 if not returnCatalogs: 

1163 # Return individual photoCalibs 

1164 yield (int(rec['visit']), int(rec['detector']), rec['filtername'], photoCalib) 

1165 else: 

1166 # Return full per-visit exposure catalogs 

1167 if rec['visit'] != lastVisit: 

1168 # This is a new visit. If the last visit was not -1, yield 

1169 # the ExposureCatalog 

1170 if lastVisit > -1: 

1171 yield (int(lastVisit), zptVisitCatalog) 

1172 else: 

1173 # We need to create a new schema 

1174 zptExpCatSchema = afwTable.ExposureTable.makeMinimalSchema() 

1175 zptExpCatSchema.addField('visit', type='I', doc='Visit number') 

1176 zptExpCatSchema.addField('detector_id', type='I', doc='Detector number') 

1177 

1178 # And start a new one 

1179 zptVisitCatalog = afwTable.ExposureCatalog(zptExpCatSchema) 

1180 zptVisitCatalog.resize(len(camera)) 

1181 

1182 # Reset the counter 

1183 zptCounter = 0 

1184 

1185 lastVisit = int(rec['visit']) 

1186 

1187 zptVisitCatalog[zptCounter].setPhotoCalib(photoCalib) 

1188 zptVisitCatalog[zptCounter]['visit'] = int(rec['visit']) 

1189 zptVisitCatalog[zptCounter]['detector_id'] = int(rec['detector']) 

1190 

1191 zptCounter += 1 

1192 

1193 # Final output of last exposure catalog 

1194 if returnCatalogs: 

1195 yield (int(lastVisit), zptVisitCatalog) 

1196 

1197 def _getChebyshevBoundedField(self, coefficients, xyMax, offset=0.0, scaling=1.0): 

1198 """ 

1199 Make a ChebyshevBoundedField from fgcm coefficients, with optional offset 

1200 and scaling. 

1201 

1202 Parameters 

1203 ---------- 

1204 coefficients: `numpy.array` 

1205 Flattened array of chebyshev coefficients 

1206 xyMax: `list` of length 2 

1207 Maximum x and y of the chebyshev bounding box 

1208 offset: `float`, optional 

1209 Absolute calibration offset. Default is 0.0 

1210 scaling: `float`, optional 

1211 Flat scaling value from fgcmBuildStars. Default is 1.0 

1212 

1213 Returns 

1214 ------- 

1215 boundedField: `lsst.afw.math.ChebyshevBoundedField` 

1216 """ 

1217 

1218 orderPlus1 = int(np.sqrt(coefficients.size)) 

1219 pars = np.zeros((orderPlus1, orderPlus1)) 

1220 

1221 bbox = lsst.geom.Box2I(lsst.geom.Point2I(0.0, 0.0), 

1222 lsst.geom.Point2I(*xyMax)) 

1223 

1224 pars[:, :] = (coefficients.reshape(orderPlus1, orderPlus1) 

1225 * (10.**(offset/-2.5))*scaling) 

1226 

1227 boundedField = afwMath.ChebyshevBoundedField(bbox, pars) 

1228 

1229 return boundedField 

1230 

1231 def _outputAtmospheres(self, dataRefDict, atmCat): 

1232 """ 

1233 Output the atmospheres. 

1234 

1235 Parameters 

1236 ---------- 

1237 dataRefDict : `dict` 

1238 All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or 

1239 `lsst.daf.butler.DeferredDatasetHandle` (gen3) 

1240 dataRef dictionary with keys: 

1241 

1242 ``"fgcmLookUpTable"`` 

1243 dataRef for the FGCM look-up table. 

1244 atmCat : `lsst.afw.table.BaseCatalog` 

1245 FGCM atmosphere parameter catalog from fgcmFitCycleTask. 

1246 

1247 Returns 

1248 ------- 

1249 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)] 

1250 Generator that returns (visit, transmissionCurve) tuples. 

1251 """ 

1252 # First, we need to grab the look-up table and key info 

1253 lutCat = dataRefDict['fgcmLookUpTable'].get() 

1254 

1255 atmosphereTableName = lutCat[0]['tablename'] 

1256 elevation = lutCat[0]['elevation'] 

1257 atmLambda = lutCat[0]['atmLambda'] 

1258 lutCat = None 

1259 

1260 # Make the atmosphere table if possible 

1261 try: 

1262 atmTable = fgcm.FgcmAtmosphereTable.initWithTableName(atmosphereTableName) 

1263 atmTable.loadTable() 

1264 except IOError: 

1265 atmTable = None 

1266 

1267 if atmTable is None: 

1268 # Try to use MODTRAN instead 

1269 try: 

1270 modGen = fgcm.ModtranGenerator(elevation) 

1271 lambdaRange = np.array([atmLambda[0], atmLambda[-1]])/10. 

1272 lambdaStep = (atmLambda[1] - atmLambda[0])/10. 

1273 except (ValueError, IOError) as e: 

1274 raise RuntimeError("FGCM look-up-table generated with modtran, " 

1275 "but modtran not configured to run.") from e 

1276 

1277 zenith = np.degrees(np.arccos(1./atmCat['secZenith'])) 

1278 

1279 for i, visit in enumerate(atmCat['visit']): 

1280 if atmTable is not None: 

1281 # Interpolate the atmosphere table 

1282 atmVals = atmTable.interpolateAtmosphere(pmb=atmCat[i]['pmb'], 

1283 pwv=atmCat[i]['pwv'], 

1284 o3=atmCat[i]['o3'], 

1285 tau=atmCat[i]['tau'], 

1286 alpha=atmCat[i]['alpha'], 

1287 zenith=zenith[i], 

1288 ctranslamstd=[atmCat[i]['cTrans'], 

1289 atmCat[i]['lamStd']]) 

1290 else: 

1291 # Run modtran 

1292 modAtm = modGen(pmb=atmCat[i]['pmb'], 

1293 pwv=atmCat[i]['pwv'], 

1294 o3=atmCat[i]['o3'], 

1295 tau=atmCat[i]['tau'], 

1296 alpha=atmCat[i]['alpha'], 

1297 zenith=zenith[i], 

1298 lambdaRange=lambdaRange, 

1299 lambdaStep=lambdaStep, 

1300 ctranslamstd=[atmCat[i]['cTrans'], 

1301 atmCat[i]['lamStd']]) 

1302 atmVals = modAtm['COMBINED'] 

1303 

1304 # Now need to create something to persist... 

1305 curve = TransmissionCurve.makeSpatiallyConstant(throughput=atmVals, 

1306 wavelengths=atmLambda, 

1307 throughputAtMin=atmVals[0], 

1308 throughputAtMax=atmVals[-1]) 

1309 

1310 yield (int(visit), curve)