Coverage for python/lsst/fgcmcal/fgcmOutputProducts.py: 15%

340 statements  

« prev     ^ index     » next       coverage.py v6.5.0, created at 2022-10-12 03:44 -0700

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 copy 

35 

36import numpy as np 

37import hpgeom as hpg 

38import esutil 

39from astropy import units 

40 

41import lsst.daf.base as dafBase 

42import lsst.pex.config as pexConfig 

43import lsst.pipe.base as pipeBase 

44from lsst.pipe.base import connectionTypes 

45from lsst.afw.image import TransmissionCurve 

46from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask 

47from lsst.meas.algorithms import ReferenceObjectLoader, LoadReferenceObjectsConfig 

48from lsst.pipe.tasks.photoCal import PhotoCalTask 

49import lsst.geom 

50import lsst.afw.image as afwImage 

51import lsst.afw.math as afwMath 

52import lsst.afw.table as afwTable 

53from lsst.meas.algorithms import DatasetConfig 

54from lsst.meas.algorithms.ingestIndexReferenceTask import addRefCatMetadata 

55 

56from .utilities import computeApproxPixelAreaFields 

57from .utilities import lookupStaticCalibrations 

58from .utilities import FGCM_ILLEGAL_VALUE 

59 

60import fgcm 

61 

62__all__ = ['FgcmOutputProductsConfig', 'FgcmOutputProductsTask'] 

63 

64 

65class FgcmOutputProductsConnections(pipeBase.PipelineTaskConnections, 

66 dimensions=("instrument",), 

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

68 camera = connectionTypes.PrerequisiteInput( 

69 doc="Camera instrument", 

70 name="camera", 

71 storageClass="Camera", 

72 dimensions=("instrument",), 

73 lookupFunction=lookupStaticCalibrations, 

74 isCalibration=True, 

75 ) 

76 

77 fgcmLookUpTable = connectionTypes.PrerequisiteInput( 

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

79 "chromatic corrections."), 

80 name="fgcmLookUpTable", 

81 storageClass="Catalog", 

82 dimensions=("instrument",), 

83 deferLoad=True, 

84 ) 

85 

86 fgcmVisitCatalog = connectionTypes.Input( 

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

88 name="fgcmVisitCatalog", 

89 storageClass="Catalog", 

90 dimensions=("instrument",), 

91 deferLoad=True, 

92 ) 

93 

94 fgcmStandardStars = connectionTypes.Input( 

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

96 name="fgcmStandardStars{cycleNumber}", 

97 storageClass="SimpleCatalog", 

98 dimensions=("instrument",), 

99 deferLoad=True, 

100 ) 

101 

102 fgcmZeropoints = connectionTypes.Input( 

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

104 name="fgcmZeropoints{cycleNumber}", 

105 storageClass="Catalog", 

106 dimensions=("instrument",), 

107 deferLoad=True, 

108 ) 

109 

110 fgcmAtmosphereParameters = connectionTypes.Input( 

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

112 name="fgcmAtmosphereParameters{cycleNumber}", 

113 storageClass="Catalog", 

114 dimensions=("instrument",), 

115 deferLoad=True, 

116 ) 

117 

118 refCat = connectionTypes.PrerequisiteInput( 

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

120 name="cal_ref_cat", 

121 storageClass="SimpleCatalog", 

122 dimensions=("skypix",), 

123 deferLoad=True, 

124 multiple=True, 

125 ) 

126 

127 fgcmPhotoCalib = connectionTypes.Output( 

128 doc=("Per-visit photometric calibrations derived from fgcm calibration. " 

129 "These catalogs use detector id for the id and are sorted for " 

130 "fast lookups of a detector."), 

131 name="fgcmPhotoCalibCatalog", 

132 storageClass="ExposureCatalog", 

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

134 multiple=True, 

135 ) 

136 

137 fgcmTransmissionAtmosphere = connectionTypes.Output( 

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

139 name="transmission_atmosphere_fgcm", 

140 storageClass="TransmissionCurve", 

141 dimensions=("instrument", 

142 "visit",), 

143 multiple=True, 

144 ) 

145 

146 fgcmOffsets = connectionTypes.Output( 

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

148 name="fgcmReferenceCalibrationOffsets", 

149 storageClass="Catalog", 

150 dimensions=("instrument",), 

151 multiple=False, 

152 ) 

153 

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

155 super().__init__(config=config) 

156 

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

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

159 

160 if not config.doReferenceCalibration: 

161 self.prerequisiteInputs.remove("refCat") 

162 if not config.doAtmosphereOutput: 

163 self.inputs.remove("fgcmAtmosphereParameters") 

164 if not config.doZeropointOutput: 

165 self.inputs.remove("fgcmZeropoints") 

166 if not config.doReferenceCalibration: 

167 self.outputs.remove("fgcmOffsets") 

168 

169 

170class FgcmOutputProductsConfig(pipeBase.PipelineTaskConfig, 

171 pipelineConnections=FgcmOutputProductsConnections): 

172 """Config for FgcmOutputProductsTask""" 

173 

174 cycleNumber = pexConfig.Field( 

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

176 dtype=int, 

177 default=0, 

178 deprecated=("This config is no longer used, and will be removed after v25. " 

179 "Please set config.connections.cycleNumber directly instead."), 

180 ) 

181 physicalFilterMap = pexConfig.DictField( 

182 doc="Mapping from 'physicalFilter' to band.", 

183 keytype=str, 

184 itemtype=str, 

185 default={}, 

186 ) 

187 # The following fields refer to calibrating from a reference 

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

189 doReferenceCalibration = pexConfig.Field( 

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

191 "This afterburner step is unnecessary if reference stars " 

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

193 dtype=bool, 

194 default=False, 

195 ) 

196 doRefcatOutput = pexConfig.Field( 

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

198 dtype=bool, 

199 default=False, 

200 deprecated="doRefcatOutput is no longer supported; this config will be removed after v24" 

201 ) 

202 doAtmosphereOutput = pexConfig.Field( 

203 doc="Output atmospheres in transmission_atmosphere_fgcm format", 

204 dtype=bool, 

205 default=True, 

206 ) 

207 doZeropointOutput = pexConfig.Field( 

208 doc="Output zeropoints in fgcm_photoCalib format", 

209 dtype=bool, 

210 default=True, 

211 ) 

212 doComposeWcsJacobian = pexConfig.Field( 

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

214 dtype=bool, 

215 default=True, 

216 ) 

217 doApplyMeanChromaticCorrection = pexConfig.Field( 

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

219 dtype=bool, 

220 default=True, 

221 ) 

222 refObjLoader = pexConfig.ConfigurableField( 

223 target=LoadIndexedReferenceObjectsTask, 

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

225 deprecated="refObjLoader is deprecated, and will be removed after v24", 

226 ) 

227 photoCal = pexConfig.ConfigurableField( 

228 target=PhotoCalTask, 

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

230 ) 

231 referencePixelizationNside = pexConfig.Field( 

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

233 dtype=int, 

234 default=64, 

235 ) 

236 referencePixelizationMinStars = pexConfig.Field( 

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

238 "to the specified reference catalog"), 

239 dtype=int, 

240 default=200, 

241 ) 

242 referenceMinMatch = pexConfig.Field( 

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

244 dtype=int, 

245 default=50, 

246 ) 

247 referencePixelizationNPixels = pexConfig.Field( 

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

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

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

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

252 dtype=int, 

253 default=100, 

254 ) 

255 datasetConfig = pexConfig.ConfigField( 

256 dtype=DatasetConfig, 

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

258 deprecated="The datasetConfig was only used for gen2; this config will be removed after v24.", 

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 

286 

287class FgcmOutputProductsTask(pipeBase.PipelineTask): 

288 """ 

289 Output products from FGCM global calibration. 

290 """ 

291 

292 ConfigClass = FgcmOutputProductsConfig 

293 _DefaultName = "fgcmOutputProducts" 

294 

295 def __init__(self, **kwargs): 

296 super().__init__(**kwargs) 

297 

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

299 handleDict = {} 

300 handleDict['camera'] = butlerQC.get(inputRefs.camera) 

301 handleDict['fgcmLookUpTable'] = butlerQC.get(inputRefs.fgcmLookUpTable) 

302 handleDict['fgcmVisitCatalog'] = butlerQC.get(inputRefs.fgcmVisitCatalog) 

303 handleDict['fgcmStandardStars'] = butlerQC.get(inputRefs.fgcmStandardStars) 

304 

305 if self.config.doZeropointOutput: 

306 handleDict['fgcmZeropoints'] = butlerQC.get(inputRefs.fgcmZeropoints) 

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

308 photoCalibRef for photoCalibRef in outputRefs.fgcmPhotoCalib} 

309 

310 if self.config.doAtmosphereOutput: 

311 handleDict['fgcmAtmosphereParameters'] = butlerQC.get(inputRefs.fgcmAtmosphereParameters) 

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

313 atmRef in outputRefs.fgcmTransmissionAtmosphere} 

314 

315 if self.config.doReferenceCalibration: 

316 refConfig = LoadReferenceObjectsConfig() 

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

318 for ref in inputRefs.refCat], 

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

320 name=self.config.connections.refCat, 

321 log=self.log, 

322 config=refConfig) 

323 else: 

324 self.refObjLoader = None 

325 

326 struct = self.run(handleDict, self.config.physicalFilterMap) 

327 

328 # Output the photoCalib exposure catalogs 

329 if struct.photoCalibCatalogs is not None: 

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

331 for visit, expCatalog in struct.photoCalibCatalogs: 

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

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

334 

335 # Output the atmospheres 

336 if struct.atmospheres is not None: 

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

338 for visit, atm in struct.atmospheres: 

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

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

341 

342 if self.config.doReferenceCalibration: 

343 # Turn offset into simple catalog for persistence if necessary 

344 schema = afwTable.Schema() 

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

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

347 offsetCat = afwTable.BaseCatalog(schema) 

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

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

350 

351 butlerQC.put(offsetCat, outputRefs.fgcmOffsets) 

352 

353 return 

354 

355 def run(self, handleDict, physicalFilterMap): 

356 """Run the output products task. 

357 

358 Parameters 

359 ---------- 

360 handleDict : `dict` 

361 All handles are `lsst.daf.butler.DeferredDatasetHandle` 

362 handle dictionary with keys: 

363 

364 ``"camera"`` 

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

366 ``"fgcmLookUpTable"`` 

367 handle for the FGCM look-up table. 

368 ``"fgcmVisitCatalog"`` 

369 handle for visit summary catalog. 

370 ``"fgcmStandardStars"`` 

371 handle for the output standard star catalog. 

372 ``"fgcmZeropoints"`` 

373 handle for the zeropoint data catalog. 

374 ``"fgcmAtmosphereParameters"`` 

375 handle for the atmosphere parameter catalog. 

376 ``"fgcmBuildStarsTableConfig"`` 

377 Config for `lsst.fgcmcal.fgcmBuildStarsTableTask`. 

378 physicalFilterMap : `dict` 

379 Dictionary of mappings from physical filter to FGCM band. 

380 

381 Returns 

382 ------- 

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

384 Output structure with keys: 

385 

386 offsets : `np.ndarray` 

387 Final reference offsets, per band. 

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

389 Generator that returns (visit, transmissionCurve) tuples. 

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

391 Generator that returns (visit, exposureCatalog) tuples. 

392 """ 

393 stdCat = handleDict['fgcmStandardStars'].get() 

394 md = stdCat.getMetadata() 

395 bands = md.getArray('BANDS') 

396 

397 if self.config.doReferenceCalibration: 

398 lutCat = handleDict['fgcmLookUpTable'].get() 

399 offsets = self._computeReferenceOffsets(stdCat, lutCat, physicalFilterMap, bands) 

400 else: 

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

402 

403 del stdCat 

404 

405 if self.config.doZeropointOutput: 

406 zptCat = handleDict['fgcmZeropoints'].get() 

407 visitCat = handleDict['fgcmVisitCatalog'].get() 

408 

409 pcgen = self._outputZeropoints(handleDict['camera'], zptCat, visitCat, offsets, bands, 

410 physicalFilterMap) 

411 else: 

412 pcgen = None 

413 

414 if self.config.doAtmosphereOutput: 

415 atmCat = handleDict['fgcmAtmosphereParameters'].get() 

416 atmgen = self._outputAtmospheres(handleDict, atmCat) 

417 else: 

418 atmgen = None 

419 

420 retStruct = pipeBase.Struct(offsets=offsets, 

421 atmospheres=atmgen) 

422 retStruct.photoCalibCatalogs = pcgen 

423 

424 return retStruct 

425 

426 def generateTractOutputProducts(self, handleDict, tract, 

427 visitCat, zptCat, atmCat, stdCat, 

428 fgcmBuildStarsConfig): 

429 """ 

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

431 

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

433 FgcmCalibrateTract. 

434 

435 Parameters 

436 ---------- 

437 handleDict : `dict` 

438 All handles are `lsst.daf.butler.DeferredDatasetHandle` 

439 handle dictionary with keys: 

440 

441 ``"camera"`` 

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

443 ``"fgcmLookUpTable"`` 

444 handle for the FGCM look-up table. 

445 tract : `int` 

446 Tract number 

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

448 FGCM visitCat from `FgcmBuildStarsTask` 

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

450 FGCM zeropoint catalog from `FgcmFitCycleTask` 

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

452 FGCM atmosphere parameter catalog from `FgcmFitCycleTask` 

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

454 FGCM standard star catalog from `FgcmFitCycleTask` 

455 fgcmBuildStarsConfig : `lsst.fgcmcal.FgcmBuildStarsConfig` 

456 Configuration object from `FgcmBuildStarsTask` 

457 

458 Returns 

459 ------- 

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

461 Output structure with keys: 

462 

463 offsets : `np.ndarray` 

464 Final reference offsets, per band. 

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

466 Generator that returns (visit, transmissionCurve) tuples. 

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

468 Generator that returns (visit, exposureCatalog) tuples. 

469 """ 

470 physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap 

471 

472 md = stdCat.getMetadata() 

473 bands = md.getArray('BANDS') 

474 

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

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

477 "in fgcmBuildStarsTask.") 

478 

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

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

481 

482 if self.config.doReferenceCalibration: 

483 lutCat = handleDict['fgcmLookUpTable'].get() 

484 offsets = self._computeReferenceOffsets(stdCat, lutCat, bands, physicalFilterMap) 

485 else: 

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

487 

488 if self.config.doZeropointOutput: 

489 pcgen = self._outputZeropoints(handleDict['camera'], zptCat, visitCat, offsets, bands, 

490 physicalFilterMap) 

491 else: 

492 pcgen = None 

493 

494 if self.config.doAtmosphereOutput: 

495 atmgen = self._outputAtmospheres(handleDict, atmCat) 

496 else: 

497 atmgen = None 

498 

499 retStruct = pipeBase.Struct(offsets=offsets, 

500 atmospheres=atmgen) 

501 retStruct.photoCalibCatalogs = pcgen 

502 

503 return retStruct 

504 

505 def _computeReferenceOffsets(self, stdCat, lutCat, physicalFilterMap, bands): 

506 """ 

507 Compute offsets relative to a reference catalog. 

508 

509 This method splits the star catalog into healpix pixels 

510 and computes the calibration transfer for a sample of 

511 these pixels to approximate the 'absolute' calibration 

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

513 absolute scale. 

514 

515 Parameters 

516 ---------- 

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

518 FGCM standard stars 

519 lutCat : `lsst.afw.table.SimpleCatalog` 

520 FGCM Look-up table 

521 physicalFilterMap : `dict` 

522 Dictionary of mappings from physical filter to FGCM band. 

523 bands : `list` [`str`] 

524 List of band names from FGCM output 

525 Returns 

526 ------- 

527 offsets : `numpy.array` of floats 

528 Per band zeropoint offsets 

529 """ 

530 

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

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

533 # calibration of each band. 

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

535 

536 goodStars = (minObs >= 1) 

537 stdCat = stdCat[goodStars] 

538 

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

540 (len(stdCat))) 

541 

542 # Associate each band with the appropriate physicalFilter and make 

543 # filterLabels 

544 filterLabels = [] 

545 

546 lutPhysicalFilters = lutCat[0]['physicalFilters'].split(',') 

547 lutStdPhysicalFilters = lutCat[0]['stdPhysicalFilters'].split(',') 

548 physicalFilterMapBands = list(physicalFilterMap.values()) 

549 physicalFilterMapFilters = list(physicalFilterMap.keys()) 

550 for band in bands: 

551 # Find a physical filter associated from the band by doing 

552 # a reverse lookup on the physicalFilterMap dict 

553 physicalFilterMapIndex = physicalFilterMapBands.index(band) 

554 physicalFilter = physicalFilterMapFilters[physicalFilterMapIndex] 

555 # Find the appropriate fgcm standard physicalFilter 

556 lutPhysicalFilterIndex = lutPhysicalFilters.index(physicalFilter) 

557 stdPhysicalFilter = lutStdPhysicalFilters[lutPhysicalFilterIndex] 

558 filterLabels.append(afwImage.FilterLabel(band=band, 

559 physical=stdPhysicalFilter)) 

560 

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

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

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

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

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

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

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

568 sourceMapper = afwTable.SchemaMapper(stdCat.schema) 

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

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

571 doc="instrumental flux (counts)") 

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

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

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

575 type='Flag', 

576 doc="bad flag") 

577 

578 # Split up the stars 

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

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

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

582 ipring = hpg.angle_to_pixel( 

583 self.config.referencePixelizationNside, 

584 stdCat['coord_ra'], 

585 stdCat['coord_dec'], 

586 degrees=False, 

587 ) 

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

589 

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

591 

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

593 (gdpix.size, 

594 self.config.referencePixelizationNside, 

595 self.config.referencePixelizationMinStars)) 

596 

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

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

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

600 else: 

601 # Sample out the pixels we want to use 

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

603 

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

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

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

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

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

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

610 

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

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

613 

614 refFluxFields = [None]*len(bands) 

615 

616 for p_index, pix in enumerate(gdpix): 

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

618 

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

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

621 # converts the index array into a boolean array 

622 selected[:] = False 

623 selected[i1a] = True 

624 

625 for b_index, filterLabel in enumerate(filterLabels): 

626 struct = self._computeOffsetOneBand(sourceMapper, badStarKey, b_index, 

627 filterLabel, stdCat, 

628 selected, refFluxFields) 

629 results['nstar'][p_index, b_index] = len(i1a) 

630 results['nmatch'][p_index, b_index] = len(struct.arrays.refMag) 

631 results['zp'][p_index, b_index] = struct.zp 

632 results['zpErr'][p_index, b_index] = struct.sigma 

633 

634 # And compute the summary statistics 

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

636 

637 for b_index, band in enumerate(bands): 

638 # make configurable 

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

640 offsets[b_index] = np.median(results['zp'][ok, b_index]) 

641 # use median absolute deviation to estimate Normal sigma 

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

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

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

645 band, offsets[b_index], madSigma) 

646 

647 return offsets 

648 

649 def _computeOffsetOneBand(self, sourceMapper, badStarKey, 

650 b_index, filterLabel, stdCat, selected, refFluxFields): 

651 """ 

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

653 stars for one pixel in one band 

654 

655 Parameters 

656 ---------- 

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

658 Mapper to go from stdCat to calibratable catalog 

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

660 Key for the field with bad stars 

661 b_index : `int` 

662 Index of the band in the star catalog 

663 filterLabel : `lsst.afw.image.FilterLabel` 

664 filterLabel with band and physical filter 

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

666 FGCM standard stars 

667 selected : `numpy.array(dtype=bool)` 

668 Boolean array of which stars are in the pixel 

669 refFluxFields : `list` 

670 List of names of flux fields for reference catalog 

671 """ 

672 

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

674 sourceCat.reserve(selected.sum()) 

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

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

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

678 * sourceCat['instFlux']) 

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

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

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

682 badStar = (stdCat['mag_std_noabs'][selected, b_index] > 90.0) 

683 for rec in sourceCat[badStar]: 

684 rec.set(badStarKey, True) 

685 

686 exposure = afwImage.ExposureF() 

687 exposure.setFilter(filterLabel) 

688 

689 if refFluxFields[b_index] is None: 

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

691 # to work around limitations of DirectMatch in PhotoCal 

692 ctr = stdCat[0].getCoord() 

693 rad = 0.05*lsst.geom.degrees 

694 refDataTest = self.refObjLoader.loadSkyCircle(ctr, rad, filterLabel.bandLabel) 

695 refFluxFields[b_index] = refDataTest.fluxField 

696 

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

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

699 calConfig.match.referenceSelection.signalToNoise.fluxField = refFluxFields[b_index] 

700 calConfig.match.referenceSelection.signalToNoise.errField = refFluxFields[b_index] + 'Err' 

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

702 config=calConfig, 

703 schema=sourceCat.getSchema()) 

704 

705 struct = calTask.run(exposure, sourceCat) 

706 

707 return struct 

708 

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

710 """ 

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

712 

713 Parameters 

714 ---------- 

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

716 SimpleCatalog as output by fgcmcal 

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

718 Zeropoint offsets to apply 

719 bands : `list` [`str`] 

720 List of band names from FGCM output 

721 

722 Returns 

723 ------- 

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

725 SimpleCatalog suitable for using as a reference catalog 

726 """ 

727 

728 sourceMapper = afwTable.SchemaMapper(fgcmStarCat.schema) 

729 minSchema = LoadIndexedReferenceObjectsTask.makeMinimalSchema(bands, 

730 addCentroid=False, 

731 addIsResolved=True, 

732 coordErrDim=0) 

733 sourceMapper.addMinimalSchema(minSchema) 

734 for band in bands: 

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

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

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

738 

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

740 formattedCat.reserve(len(fgcmStarCat)) 

741 formattedCat.extend(fgcmStarCat, mapper=sourceMapper) 

742 

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

744 

745 for b, band in enumerate(bands): 

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

747 # We want fluxes in nJy from calibrated AB magnitudes 

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

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

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

751 

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

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

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

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

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

757 

758 addRefCatMetadata(formattedCat) 

759 

760 return formattedCat 

761 

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

763 physicalFilterMap, tract=None): 

764 """Output the zeropoints in fgcm_photoCalib format. 

765 

766 Parameters 

767 ---------- 

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

769 Camera from the butler. 

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

771 FGCM zeropoint catalog from `FgcmFitCycleTask`. 

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

773 FGCM visitCat from `FgcmBuildStarsTask`. 

774 offsets : `numpy.array` 

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

776 bands : `list` [`str`] 

777 List of band names from FGCM output. 

778 physicalFilterMap : `dict` 

779 Dictionary of mappings from physical filter to FGCM band. 

780 tract: `int`, optional 

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

782 

783 Returns 

784 ------- 

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

786 Generator that returns (visit, exposureCatalog) tuples. 

787 """ 

788 # Select visit/ccds where we have a calibration 

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

790 # ccds. 

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

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

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

794 & (zptCat['fgcmZpt'] > FGCM_ILLEGAL_VALUE)) 

795 

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

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

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

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

800 for allBadVisit in allBadVisits: 

801 self.log.warning(f'No suitable photoCalib for visit {allBadVisit}') 

802 

803 # Get a mapping from filtername to the offsets 

804 offsetMapping = {} 

805 for f in physicalFilterMap: 

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

807 if physicalFilterMap[f] in bands: 

808 offsetMapping[f] = offsets[bands.index(physicalFilterMap[f])] 

809 

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

811 ccdMapping = {} 

812 for ccdIndex, detector in enumerate(camera): 

813 ccdMapping[detector.getId()] = ccdIndex 

814 

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

816 scalingMapping = {} 

817 for rec in visitCat: 

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

819 

820 if self.config.doComposeWcsJacobian: 

821 approxPixelAreaFields = computeApproxPixelAreaFields(camera) 

822 

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

824 lastVisit = -1 

825 zptVisitCatalog = None 

826 

827 metadata = dafBase.PropertyList() 

828 metadata.add("COMMENT", "Catalog id is detector id, sorted.") 

829 metadata.add("COMMENT", "Only detectors with data have entries.") 

830 

831 for rec in zptCat[selected]: 

832 # Retrieve overall scaling 

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

834 

835 # The postCalibrationOffset describe any zeropoint offsets 

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

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

838 # second part comes from the mean chromatic correction 

839 # (if configured). 

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

841 if self.config.doApplyMeanChromaticCorrection: 

842 postCalibrationOffset += rec['fgcmDeltaChrom'] 

843 

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

845 rec['fgcmfZptChebXyMax']) 

846 # Convert from FGCM AB to nJy 

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

848 rec['fgcmfZptChebXyMax'], 

849 offset=postCalibrationOffset, 

850 scaling=scaling) 

851 

852 if self.config.doComposeWcsJacobian: 

853 

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

855 fgcmSuperStarField, 

856 fgcmZptField]) 

857 else: 

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

859 # fgcmZptField 

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

861 

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

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

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

865 photoCalib = afwImage.PhotoCalib(calibrationMean=calibCenter, 

866 calibrationErr=calibErr, 

867 calibration=fgcmField, 

868 isConstant=False) 

869 

870 # Return full per-visit exposure catalogs 

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

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

873 # the ExposureCatalog 

874 if lastVisit > -1: 

875 # ensure that the detectors are in sorted order, for fast lookups 

876 zptVisitCatalog.sort() 

877 yield (int(lastVisit), zptVisitCatalog) 

878 else: 

879 # We need to create a new schema 

880 zptExpCatSchema = afwTable.ExposureTable.makeMinimalSchema() 

881 zptExpCatSchema.addField('visit', type='L', doc='Visit number') 

882 

883 # And start a new one 

884 zptVisitCatalog = afwTable.ExposureCatalog(zptExpCatSchema) 

885 zptVisitCatalog.setMetadata(metadata) 

886 

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

888 

889 catRecord = zptVisitCatalog.addNew() 

890 catRecord['id'] = int(rec['detector']) 

891 catRecord['visit'] = rec['visit'] 

892 catRecord.setPhotoCalib(photoCalib) 

893 

894 # Final output of last exposure catalog 

895 # ensure that the detectors are in sorted order, for fast lookups 

896 zptVisitCatalog.sort() 

897 yield (int(lastVisit), zptVisitCatalog) 

898 

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

900 """ 

901 Make a ChebyshevBoundedField from fgcm coefficients, with optional offset 

902 and scaling. 

903 

904 Parameters 

905 ---------- 

906 coefficients: `numpy.array` 

907 Flattened array of chebyshev coefficients 

908 xyMax: `list` of length 2 

909 Maximum x and y of the chebyshev bounding box 

910 offset: `float`, optional 

911 Absolute calibration offset. Default is 0.0 

912 scaling: `float`, optional 

913 Flat scaling value from fgcmBuildStars. Default is 1.0 

914 

915 Returns 

916 ------- 

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

918 """ 

919 

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

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

922 

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

924 lsst.geom.Point2I(*xyMax)) 

925 

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

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

928 

929 boundedField = afwMath.ChebyshevBoundedField(bbox, pars) 

930 

931 return boundedField 

932 

933 def _outputAtmospheres(self, handleDict, atmCat): 

934 """ 

935 Output the atmospheres. 

936 

937 Parameters 

938 ---------- 

939 handleDict : `dict` 

940 All data handles are `lsst.daf.butler.DeferredDatasetHandle` 

941 The handleDict has the follownig keys: 

942 

943 ``"fgcmLookUpTable"`` 

944 handle for the FGCM look-up table. 

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

946 FGCM atmosphere parameter catalog from fgcmFitCycleTask. 

947 

948 Returns 

949 ------- 

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

951 Generator that returns (visit, transmissionCurve) tuples. 

952 """ 

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

954 lutCat = handleDict['fgcmLookUpTable'].get() 

955 

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

957 elevation = lutCat[0]['elevation'] 

958 atmLambda = lutCat[0]['atmLambda'] 

959 lutCat = None 

960 

961 # Make the atmosphere table if possible 

962 try: 

963 atmTable = fgcm.FgcmAtmosphereTable.initWithTableName(atmosphereTableName) 

964 atmTable.loadTable() 

965 except IOError: 

966 atmTable = None 

967 

968 if atmTable is None: 

969 # Try to use MODTRAN instead 

970 try: 

971 modGen = fgcm.ModtranGenerator(elevation) 

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

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

974 except (ValueError, IOError) as e: 

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

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

977 

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

979 

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

981 if atmTable is not None: 

982 # Interpolate the atmosphere table 

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

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

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

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

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

988 zenith=zenith[i], 

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

990 atmCat[i]['lamStd']]) 

991 else: 

992 # Run modtran 

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

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

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

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

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

998 zenith=zenith[i], 

999 lambdaRange=lambdaRange, 

1000 lambdaStep=lambdaStep, 

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

1002 atmCat[i]['lamStd']]) 

1003 atmVals = modAtm['COMBINED'] 

1004 

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

1006 curve = TransmissionCurve.makeSpatiallyConstant(throughput=atmVals, 

1007 wavelengths=atmLambda, 

1008 throughputAtMin=atmVals[0], 

1009 throughputAtMax=atmVals[-1]) 

1010 

1011 yield (int(visit), curve)