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# This file is part of jointcal. 

2# 

3# Developed for the LSST Data Management System. 

4# This product includes software developed by the LSST Project 

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

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

7# for details of code ownership. 

8# 

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

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

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

12# (at your option) any later version. 

13# 

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

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

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

17# GNU General Public License for more details. 

18# 

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

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

21 

22import dataclasses 

23import collections 

24import os 

25 

26import astropy.time 

27import numpy as np 

28import astropy.units as u 

29 

30import lsst.geom 

31import lsst.utils 

32import lsst.pex.config as pexConfig 

33import lsst.pipe.base as pipeBase 

34from lsst.afw.image import fluxErrFromABMagErr 

35import lsst.pex.exceptions as pexExceptions 

36import lsst.afw.cameraGeom 

37import lsst.afw.table 

38import lsst.log 

39from lsst.pipe.tasks.colorterms import ColortermLibrary 

40from lsst.verify import Job, Measurement 

41 

42from lsst.meas.algorithms import (LoadIndexedReferenceObjectsTask, ReferenceSourceSelectorTask, 

43 ReferenceObjectLoader) 

44from lsst.meas.algorithms.sourceSelector import sourceSelectorRegistry 

45 

46from .dataIds import PerTractCcdDataIdContainer 

47 

48import lsst.jointcal 

49from lsst.jointcal import MinimizeResult 

50 

51__all__ = ["JointcalConfig", "JointcalRunner", "JointcalTask"] 

52 

53Photometry = collections.namedtuple('Photometry', ('fit', 'model')) 

54Astrometry = collections.namedtuple('Astrometry', ('fit', 'model', 'sky_to_tan_projection')) 

55 

56 

57# TODO: move this to MeasurementSet in lsst.verify per DM-12655. 

58def add_measurement(job, name, value): 

59 meas = Measurement(job.metrics[name], value) 

60 job.measurements.insert(meas) 

61 

62 

63class JointcalRunner(pipeBase.ButlerInitializedTaskRunner): 

64 """Subclass of TaskRunner for jointcalTask (gen2) 

65 

66 jointcalTask.runDataRef() takes a number of arguments, one of which is a list of dataRefs 

67 extracted from the command line (whereas most CmdLineTasks' runDataRef methods take 

68 single dataRef, are are called repeatedly). This class transforms the processed 

69 arguments generated by the ArgumentParser into the arguments expected by 

70 Jointcal.runDataRef(). 

71 

72 See pipeBase.TaskRunner for more information. 

73 """ 

74 

75 @staticmethod 

76 def getTargetList(parsedCmd, **kwargs): 

77 """ 

78 Return a list of tuples per tract, each containing (dataRefs, kwargs). 

79 

80 Jointcal operates on lists of dataRefs simultaneously. 

81 """ 

82 kwargs['butler'] = parsedCmd.butler 

83 

84 # organize data IDs by tract 

85 refListDict = {} 

86 for ref in parsedCmd.id.refList: 

87 refListDict.setdefault(ref.dataId["tract"], []).append(ref) 

88 # we call runDataRef() once with each tract 

89 result = [(refListDict[tract], kwargs) for tract in sorted(refListDict.keys())] 

90 return result 

91 

92 def __call__(self, args): 

93 """ 

94 Parameters 

95 ---------- 

96 args 

97 Arguments for Task.runDataRef() 

98 

99 Returns 

100 ------- 

101 pipe.base.Struct 

102 if self.doReturnResults is False: 

103 

104 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise. 

105 

106 if self.doReturnResults is True: 

107 

108 - ``result``: the result of calling jointcal.runDataRef() 

109 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise. 

110 """ 

111 exitStatus = 0 # exit status for shell 

112 

113 # NOTE: cannot call self.makeTask because that assumes args[0] is a single dataRef. 

114 dataRefList, kwargs = args 

115 butler = kwargs.pop('butler') 

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

117 result = None 

118 try: 

119 result = task.runDataRef(dataRefList, **kwargs) 

120 exitStatus = result.exitStatus 

121 job_path = butler.get('verify_job_filename') 

122 result.job.write(job_path[0]) 

123 except Exception as e: # catch everything, sort it out later. 

124 if self.doRaise: 

125 raise e 

126 else: 

127 exitStatus = 1 

128 eName = type(e).__name__ 

129 tract = dataRefList[0].dataId['tract'] 

130 task.log.fatal("Failed processing tract %s, %s: %s", tract, eName, e) 

131 

132 # Put the butler back into kwargs for the other Tasks. 

133 kwargs['butler'] = butler 

134 if self.doReturnResults: 

135 return pipeBase.Struct(result=result, exitStatus=exitStatus) 

136 else: 

137 return pipeBase.Struct(exitStatus=exitStatus) 

138 

139 

140class JointcalTaskConnections(pipeBase.PipelineTaskConnections, 

141 dimensions=("skymap", "tract", "instrument", "physical_filter")): 

142 """Middleware input/output connections for jointcal data.""" 

143 inputCamera = pipeBase.connectionTypes.Input( 

144 doc="The camera instrument that took these observations.", 

145 name="camera", 

146 storageClass="Camera", 

147 dimensions=("instrument",), 

148 isCalibration=True 

149 ) 

150 inputSourceTableVisit = pipeBase.connectionTypes.Input( 

151 doc="Source table in parquet format, per visit", 

152 name="sourceTable_visit", 

153 storageClass="DataFrame", 

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

155 deferLoad=True, 

156 multiple=True, 

157 ) 

158 inputVisitSummary = pipeBase.connectionTypes.Input( 

159 doc=("Per-visit consolidated exposure metadata built from calexps. " 

160 "These catalogs use detector id for the id and must be sorted for " 

161 "fast lookups of a detector."), 

162 name="visitSummary", 

163 storageClass="ExposureCatalog", 

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

165 deferLoad=True, 

166 multiple=True, 

167 ) 

168 # TODO DM-28991: This does not load enough refcat data: the graph only to gets 

169 # refcats that touch the tract, but we need to go larger because we're 

170 # loading visits that may extend further. How to test that? 

171 astrometryRefCat = pipeBase.connectionTypes.Input( 

172 doc="The astrometry reference catalog to match to loaded input catalog sources.", 

173 name="gaia_dr2_20200414", 

174 storageClass="SimpleCatalog", 

175 dimensions=("htm7",), 

176 deferLoad=True, 

177 multiple=True 

178 ) 

179 photometryRefCat = pipeBase.connectionTypes.Input( 

180 doc="The photometry reference catalog to match to loaded input catalog sources.", 

181 name="ps1_pv3_3pi_20170110", 

182 storageClass="SimpleCatalog", 

183 dimensions=("htm7",), 

184 deferLoad=True, 

185 multiple=True 

186 ) 

187 

188 outputWcs = pipeBase.connectionTypes.Output( 

189 doc=("Per-tract, per-visit world coordinate systems derived from the fitted model." 

190 " These catalogs only contain entries for detectors with an output, and use" 

191 " the detector id for the catalog id, sorted on id for fast lookups of a detector."), 

192 name="jointcalSkyWcsCatalog", 

193 storageClass="ExposureCatalog", 

194 dimensions=("instrument", "visit", "skymap", "tract"), 

195 multiple=True 

196 ) 

197 outputPhotoCalib = pipeBase.connectionTypes.Output( 

198 doc=("Per-tract, per-visit photometric calibrations derived from the fitted model." 

199 " These catalogs only contain entries for detectors with an output, and use" 

200 " the detector id for the catalog id, sorted on id for fast lookups of a detector."), 

201 name="jointcalPhotoCalibCatalog", 

202 storageClass="ExposureCatalog", 

203 dimensions=("instrument", "visit", "skymap", "tract"), 

204 multiple=True 

205 ) 

206 

207 

208class JointcalConfig(pipeBase.PipelineTaskConfig, 

209 pipelineConnections=JointcalTaskConnections): 

210 """Configuration for JointcalTask""" 

211 

212 doAstrometry = pexConfig.Field( 

213 doc="Fit astrometry and write the fitted result.", 

214 dtype=bool, 

215 default=True 

216 ) 

217 doPhotometry = pexConfig.Field( 

218 doc="Fit photometry and write the fitted result.", 

219 dtype=bool, 

220 default=True 

221 ) 

222 coaddName = pexConfig.Field( 

223 doc="Type of coadd, typically deep or goodSeeing", 

224 dtype=str, 

225 default="deep" 

226 ) 

227 # TODO DM-29008: Change this to "ApFlux_12_0" before gen2 removal. 

228 sourceFluxType = pexConfig.Field( 

229 dtype=str, 

230 doc="Source flux field to use in source selection and to get fluxes from the catalog.", 

231 default='Calib' 

232 ) 

233 positionErrorPedestal = pexConfig.Field( 

234 doc="Systematic term to apply to the measured position error (pixels)", 

235 dtype=float, 

236 default=0.02, 

237 ) 

238 photometryErrorPedestal = pexConfig.Field( 

239 doc="Systematic term to apply to the measured error on flux or magnitude as a " 

240 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).", 

241 dtype=float, 

242 default=0.0, 

243 ) 

244 # TODO: DM-6885 matchCut should be an geom.Angle 

245 matchCut = pexConfig.Field( 

246 doc="Matching radius between fitted and reference stars (arcseconds)", 

247 dtype=float, 

248 default=3.0, 

249 ) 

250 minMeasurements = pexConfig.Field( 

251 doc="Minimum number of associated measured stars for a fitted star to be included in the fit", 

252 dtype=int, 

253 default=2, 

254 ) 

255 minMeasuredStarsPerCcd = pexConfig.Field( 

256 doc="Minimum number of measuredStars per ccdImage before printing warnings", 

257 dtype=int, 

258 default=100, 

259 ) 

260 minRefStarsPerCcd = pexConfig.Field( 

261 doc="Minimum number of measuredStars per ccdImage before printing warnings", 

262 dtype=int, 

263 default=30, 

264 ) 

265 allowLineSearch = pexConfig.Field( 

266 doc="Allow a line search during minimization, if it is reasonable for the model" 

267 " (models with a significant non-linear component, e.g. constrainedPhotometry).", 

268 dtype=bool, 

269 default=False 

270 ) 

271 astrometrySimpleOrder = pexConfig.Field( 

272 doc="Polynomial order for fitting the simple astrometry model.", 

273 dtype=int, 

274 default=3, 

275 ) 

276 astrometryChipOrder = pexConfig.Field( 

277 doc="Order of the per-chip transform for the constrained astrometry model.", 

278 dtype=int, 

279 default=1, 

280 ) 

281 astrometryVisitOrder = pexConfig.Field( 

282 doc="Order of the per-visit transform for the constrained astrometry model.", 

283 dtype=int, 

284 default=5, 

285 ) 

286 useInputWcs = pexConfig.Field( 

287 doc="Use the input calexp WCSs to initialize a SimpleAstrometryModel.", 

288 dtype=bool, 

289 default=True, 

290 ) 

291 astrometryModel = pexConfig.ChoiceField( 

292 doc="Type of model to fit to astrometry", 

293 dtype=str, 

294 default="constrained", 

295 allowed={"simple": "One polynomial per ccd", 

296 "constrained": "One polynomial per ccd, and one polynomial per visit"} 

297 ) 

298 photometryModel = pexConfig.ChoiceField( 

299 doc="Type of model to fit to photometry", 

300 dtype=str, 

301 default="constrainedMagnitude", 

302 allowed={"simpleFlux": "One constant zeropoint per ccd and visit, fitting in flux space.", 

303 "constrainedFlux": "Constrained zeropoint per ccd, and one polynomial per visit," 

304 " fitting in flux space.", 

305 "simpleMagnitude": "One constant zeropoint per ccd and visit," 

306 " fitting in magnitude space.", 

307 "constrainedMagnitude": "Constrained zeropoint per ccd, and one polynomial per visit," 

308 " fitting in magnitude space.", 

309 } 

310 ) 

311 applyColorTerms = pexConfig.Field( 

312 doc="Apply photometric color terms to reference stars?" 

313 "Requires that colorterms be set to a ColortermLibrary", 

314 dtype=bool, 

315 default=False 

316 ) 

317 colorterms = pexConfig.ConfigField( 

318 doc="Library of photometric reference catalog name to color term dict.", 

319 dtype=ColortermLibrary, 

320 ) 

321 photometryVisitOrder = pexConfig.Field( 

322 doc="Order of the per-visit polynomial transform for the constrained photometry model.", 

323 dtype=int, 

324 default=7, 

325 ) 

326 photometryDoRankUpdate = pexConfig.Field( 

327 doc=("Do the rank update step during minimization. " 

328 "Skipping this can help deal with models that are too non-linear."), 

329 dtype=bool, 

330 default=True, 

331 ) 

332 astrometryDoRankUpdate = pexConfig.Field( 

333 doc=("Do the rank update step during minimization (should not change the astrometry fit). " 

334 "Skipping this can help deal with models that are too non-linear."), 

335 dtype=bool, 

336 default=True, 

337 ) 

338 outlierRejectSigma = pexConfig.Field( 

339 doc="How many sigma to reject outliers at during minimization.", 

340 dtype=float, 

341 default=5.0, 

342 ) 

343 maxPhotometrySteps = pexConfig.Field( 

344 doc="Maximum number of minimize iterations to take when fitting photometry.", 

345 dtype=int, 

346 default=20, 

347 ) 

348 maxAstrometrySteps = pexConfig.Field( 

349 doc="Maximum number of minimize iterations to take when fitting photometry.", 

350 dtype=int, 

351 default=20, 

352 ) 

353 astrometryRefObjLoader = pexConfig.ConfigurableField( 

354 target=LoadIndexedReferenceObjectsTask, 

355 doc="Reference object loader for astrometric fit", 

356 ) 

357 photometryRefObjLoader = pexConfig.ConfigurableField( 

358 target=LoadIndexedReferenceObjectsTask, 

359 doc="Reference object loader for photometric fit", 

360 ) 

361 sourceSelector = sourceSelectorRegistry.makeField( 

362 doc="How to select sources for cross-matching", 

363 default="astrometry" 

364 ) 

365 astrometryReferenceSelector = pexConfig.ConfigurableField( 

366 target=ReferenceSourceSelectorTask, 

367 doc="How to down-select the loaded astrometry reference catalog.", 

368 ) 

369 photometryReferenceSelector = pexConfig.ConfigurableField( 

370 target=ReferenceSourceSelectorTask, 

371 doc="How to down-select the loaded photometry reference catalog.", 

372 ) 

373 astrometryReferenceErr = pexConfig.Field( 

374 doc=("Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. " 

375 "If None, then raise an exception if the reference catalog is missing coordinate errors. " 

376 "If specified, overrides any existing `coord_*Err` values."), 

377 dtype=float, 

378 default=None, 

379 optional=True 

380 ) 

381 

382 # configs for outputting debug information 

383 writeInitMatrix = pexConfig.Field( 

384 dtype=bool, 

385 doc=("Write the pre/post-initialization Hessian and gradient to text files, for debugging. " 

386 "The output files will be of the form 'astrometry_preinit-mat.txt', in the current directory. " 

387 "Note that these files are the dense versions of the matrix, and so may be very large."), 

388 default=False 

389 ) 

390 writeChi2FilesInitialFinal = pexConfig.Field( 

391 dtype=bool, 

392 doc="Write .csv files containing the contributions to chi2 for the initialization and final fit.", 

393 default=False 

394 ) 

395 writeChi2FilesOuterLoop = pexConfig.Field( 

396 dtype=bool, 

397 doc="Write .csv files containing the contributions to chi2 for the outer fit loop.", 

398 default=False 

399 ) 

400 writeInitialModel = pexConfig.Field( 

401 dtype=bool, 

402 doc=("Write the pre-initialization model to text files, for debugging." 

403 " Output is written to `initial[Astro|Photo]metryModel.txt` in the current working directory."), 

404 default=False 

405 ) 

406 debugOutputPath = pexConfig.Field( 

407 dtype=str, 

408 default=".", 

409 doc=("Path to write debug output files to. Used by " 

410 "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.") 

411 ) 

412 detailedProfile = pexConfig.Field( 

413 dtype=bool, 

414 default=False, 

415 doc="Output separate profiling information for different parts of jointcal, e.g. data read, fitting" 

416 ) 

417 

418 def validate(self): 

419 super().validate() 

420 if self.doPhotometry and self.applyColorTerms and len(self.colorterms.data) == 0: 

421 msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary." 

422 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg) 

423 if self.doAstrometry and not self.doPhotometry and self.applyColorTerms: 

424 msg = ("Only doing astrometry, but Colorterms are not applied for astrometry;" 

425 "applyColorTerms=True will be ignored.") 

426 lsst.log.warn(msg) 

427 

428 def setDefaults(self): 

429 # Use science source selector which can filter on extendedness, SNR, and whether blended 

430 self.sourceSelector.name = 'science' 

431 # Use only stars because aperture fluxes of galaxies are biased and depend on seeing 

432 self.sourceSelector['science'].doUnresolved = True 

433 # with dependable signal to noise ratio. 

434 self.sourceSelector['science'].doSignalToNoise = True 

435 # Min SNR must be > 0 because jointcal cannot handle negative fluxes, 

436 # and S/N > 10 to use sources that are not too faint, and thus better measured. 

437 self.sourceSelector['science'].signalToNoise.minimum = 10. 

438 # Base SNR on CalibFlux because that is the flux jointcal that fits and must be positive 

439 fluxField = f"slot_{self.sourceFluxType}Flux_instFlux" 

440 self.sourceSelector['science'].signalToNoise.fluxField = fluxField 

441 self.sourceSelector['science'].signalToNoise.errField = fluxField + "Err" 

442 # Do not trust blended sources' aperture fluxes which also depend on seeing 

443 self.sourceSelector['science'].doIsolated = True 

444 # Do not trust either flux or centroid measurements with flags, 

445 # chosen from the usual QA flags for stars) 

446 self.sourceSelector['science'].doFlags = True 

447 badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated', 

448 'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag', 

449 'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter'] 

450 self.sourceSelector['science'].flags.bad = badFlags 

451 

452 # Default to Gaia-DR2 (with proper motions) for astrometry and 

453 # PS1-DR1 for photometry, with a reasonable initial filterMap. 

454 self.astrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414" 

455 self.astrometryRefObjLoader.requireProperMotion = True 

456 self.astrometryRefObjLoader.anyFilterMapsToThis = 'phot_g_mean' 

457 self.photometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110" 

458 

459 

460def writeModel(model, filename, log): 

461 """Write model to outfile.""" 

462 with open(filename, "w") as file: 

463 file.write(repr(model)) 

464 log.info("Wrote %s to file: %s", model, filename) 

465 

466 

467@dataclasses.dataclass 

468class JointcalInputData: 

469 """The input data jointcal needs for each detector/visit.""" 

470 visit: int 

471 """The visit identifier of this exposure.""" 

472 catalog: lsst.afw.table.SourceCatalog 

473 """The catalog derived from this exposure.""" 

474 visitInfo: lsst.afw.image.VisitInfo 

475 """The VisitInfo of this exposure.""" 

476 detector: lsst.afw.cameraGeom.Detector 

477 """The detector of this exposure.""" 

478 photoCalib: lsst.afw.image.PhotoCalib 

479 """The photometric calibration of this exposure.""" 

480 wcs: lsst.afw.geom.skyWcs 

481 """The WCS of this exposure.""" 

482 bbox: lsst.geom.Box2I 

483 """The bounding box of this exposure.""" 

484 filter: lsst.afw.image.FilterLabel 

485 """The filter of this exposure.""" 

486 

487 

488class JointcalTask(pipeBase.PipelineTask, pipeBase.CmdLineTask): 

489 """Astrometricly and photometricly calibrate across multiple visits of the 

490 same field. 

491 

492 Parameters 

493 ---------- 

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

495 The butler is passed to the refObjLoader constructor in case it is 

496 needed. Ignored if the refObjLoader argument provides a loader directly. 

497 Used to initialize the astrometry and photometry refObjLoaders. 

498 initInputs : `dict`, optional 

499 Dictionary used to initialize PipelineTasks (empty for jointcal). 

500 """ 

501 

502 ConfigClass = JointcalConfig 

503 RunnerClass = JointcalRunner 

504 _DefaultName = "jointcal" 

505 

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

507 super().__init__(**kwargs) 

508 self.makeSubtask("sourceSelector") 

509 if self.config.doAstrometry: 

510 if initInputs is None: 

511 # gen3 middleware does refcat things internally (and will not have a butler here) 

512 self.makeSubtask('astrometryRefObjLoader', butler=butler) 

513 self.makeSubtask("astrometryReferenceSelector") 

514 else: 

515 self.astrometryRefObjLoader = None 

516 if self.config.doPhotometry: 

517 if initInputs is None: 

518 # gen3 middleware does refcat things internally (and will not have a butler here) 

519 self.makeSubtask('photometryRefObjLoader', butler=butler) 

520 self.makeSubtask("photometryReferenceSelector") 

521 else: 

522 self.photometryRefObjLoader = None 

523 

524 # To hold various computed metrics for use by tests 

525 self.job = Job.load_metrics_package(subset='jointcal') 

526 

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

528 # We override runQuantum to set up the refObjLoaders and write the 

529 # outputs to the correct refs. 

530 inputs = butlerQC.get(inputRefs) 

531 # We want the tract number for writing debug files 

532 tract = butlerQC.quantum.dataId['tract'] 

533 if self.config.doAstrometry: 

534 self.astrometryRefObjLoader = ReferenceObjectLoader( 

535 dataIds=[ref.datasetRef.dataId 

536 for ref in inputRefs.astrometryRefCat], 

537 refCats=inputs.pop('astrometryRefCat'), 

538 config=self.config.astrometryRefObjLoader, 

539 log=self.log) 

540 if self.config.doPhotometry: 

541 self.photometryRefObjLoader = ReferenceObjectLoader( 

542 dataIds=[ref.datasetRef.dataId 

543 for ref in inputRefs.photometryRefCat], 

544 refCats=inputs.pop('photometryRefCat'), 

545 config=self.config.photometryRefObjLoader, 

546 log=self.log) 

547 outputs = self.run(**inputs, tract=tract) 

548 if self.config.doAstrometry: 

549 self._put_output(butlerQC, outputs.outputWcs, outputRefs.outputWcs, 

550 inputs['inputCamera'], "setWcs") 

551 if self.config.doPhotometry: 

552 self._put_output(butlerQC, outputs.outputPhotoCalib, outputRefs.outputPhotoCalib, 

553 inputs['inputCamera'], "setPhotoCalib") 

554 

555 def _put_output(self, butlerQC, outputs, outputRefs, camera, setter): 

556 """Persist the output datasets to their appropriate datarefs. 

557 

558 Parameters 

559 ---------- 

560 butlerQC : `ButlerQuantumContext` 

561 A butler which is specialized to operate in the context of a 

562 `lsst.daf.butler.Quantum`; This is the input to `runQuantum`. 

563 outputs : `dict` [`tuple`, `lsst.afw.geom.SkyWcs`] or 

564 `dict` [`tuple, `lsst.afw.image.PhotoCalib`] 

565 The fitted objects to persist. 

566 outputRefs : `list` [`OutputQuantizedConnection`] 

567 The DatasetRefs to persist the data to. 

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

569 The camera for this instrument, to get detector ids from. 

570 setter : `str` 

571 The method to call on the ExposureCatalog to set each output. 

572 """ 

573 schema = lsst.afw.table.ExposureTable.makeMinimalSchema() 

574 schema.addField('visit', type='I', doc='Visit number') 

575 

576 def new_catalog(visit, size): 

577 """Return an catalog ready to be filled with appropriate output.""" 

578 catalog = lsst.afw.table.ExposureCatalog(schema) 

579 catalog.resize(size) 

580 catalog['visit'] = visit 

581 metadata = lsst.daf.base.PropertyList() 

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

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

584 return catalog 

585 

586 # count how many detectors have output for each visit 

587 detectors_per_visit = collections.defaultdict(int) 

588 for key in outputs: 

589 # key is (visit, detector_id), and we only need visit here 

590 detectors_per_visit[key[0]] += 1 

591 

592 for ref in outputRefs: 

593 visit = ref.dataId['visit'] 

594 catalog = new_catalog(visit, detectors_per_visit[visit]) 

595 # Iterate over every detector and skip the ones we don't have output for. 

596 i = 0 

597 for detector in camera: 

598 detectorId = detector.getId() 

599 key = (ref.dataId['visit'], detectorId) 

600 if key not in outputs: 

601 # skip detectors we don't have output for 

602 self.log.debug("No %s output for detector %s in visit %s", 

603 setter[3:], detectorId, visit) 

604 continue 

605 

606 catalog[i].setId(detectorId) 

607 getattr(catalog[i], setter)(outputs[key]) 

608 i += 1 

609 

610 catalog.sort() # ensure that the detectors are in sorted order, for fast lookups 

611 butlerQC.put(catalog, ref) 

612 self.log.info("Wrote %s detectors to %s", i, ref) 

613 

614 def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None): 

615 # Docstring inherited. 

616 

617 # We take values out of the Parquet table, and put them in "flux_", 

618 # and the config.sourceFluxType field is used during that extraction, 

619 # so just use "flux" here. 

620 sourceFluxField = "flux" 

621 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField) 

622 associations = lsst.jointcal.Associations() 

623 self.focalPlaneBBox = inputCamera.getFpBBox() 

624 oldWcsList, bands = self._load_data(inputSourceTableVisit, 

625 inputVisitSummary, 

626 associations, 

627 jointcalControl, 

628 inputCamera) 

629 

630 boundingCircle, center, radius, defaultFilter = self._prep_sky(associations, bands) 

631 epoch = self._compute_proper_motion_epoch(associations.getCcdImageList()) 

632 

633 if self.config.doAstrometry: 

634 astrometry = self._do_load_refcat_and_fit(associations, defaultFilter, center, radius, 

635 name="astrometry", 

636 refObjLoader=self.astrometryRefObjLoader, 

637 referenceSelector=self.astrometryReferenceSelector, 

638 fit_function=self._fit_astrometry, 

639 tract=tract, 

640 epoch=epoch) 

641 astrometry_output = self._make_output(associations.getCcdImageList(), 

642 astrometry.model, 

643 "makeSkyWcs") 

644 else: 

645 astrometry_output = None 

646 

647 if self.config.doPhotometry: 

648 photometry = self._do_load_refcat_and_fit(associations, defaultFilter, center, radius, 

649 name="photometry", 

650 refObjLoader=self.photometryRefObjLoader, 

651 referenceSelector=self.photometryReferenceSelector, 

652 fit_function=self._fit_photometry, 

653 tract=tract, 

654 epoch=epoch, 

655 reject_bad_fluxes=True) 

656 photometry_output = self._make_output(associations.getCcdImageList(), 

657 photometry.model, 

658 "toPhotoCalib") 

659 else: 

660 photometry_output = None 

661 

662 return pipeBase.Struct(outputWcs=astrometry_output, 

663 outputPhotoCalib=photometry_output, 

664 job=self.job, 

665 astrometryRefObjLoader=self.astrometryRefObjLoader, 

666 photometryRefObjLoader=self.photometryRefObjLoader) 

667 

668 def _make_schema_table(self): 

669 """Return an afw SourceTable to use as a base for creating the 

670 SourceCatalog to insert values from the dataFrame into. 

671 

672 Returns 

673 ------- 

674 table : `lsst.afw.table.SourceTable` 

675 Table with schema and slots to use to make SourceCatalogs. 

676 """ 

677 schema = lsst.afw.table.SourceTable.makeMinimalSchema() 

678 schema.addField("centroid_x", "D") 

679 schema.addField("centroid_y", "D") 

680 schema.addField("centroid_xErr", "F") 

681 schema.addField("centroid_yErr", "F") 

682 schema.addField("shape_xx", "D") 

683 schema.addField("shape_yy", "D") 

684 schema.addField("shape_xy", "D") 

685 schema.addField("flux_instFlux", "D") 

686 schema.addField("flux_instFluxErr", "D") 

687 table = lsst.afw.table.SourceTable.make(schema) 

688 table.defineCentroid("centroid") 

689 table.defineShape("shape") 

690 return table 

691 

692 def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId): 

693 """Return an afw SourceCatalog extracted from a visit-level dataframe, 

694 limited to just one detector. 

695 

696 Parameters 

697 ---------- 

698 table : `lsst.afw.table.SourceTable` 

699 Table factory to use to make the SourceCatalog that will be 

700 populated with data from ``visitCatalog``. 

701 visitCatalog : `pandas.DataFrame` 

702 DataFrame to extract a detector catalog from. 

703 detectorId : `int` 

704 Numeric id of the detector to extract from ``visitCatalog``. 

705 

706 Returns 

707 ------- 

708 catalog : `lsst.afw.table.SourceCatalog` 

709 Detector-level catalog extracted from ``visitCatalog``. 

710 """ 

711 # map from dataFrame column to afw table column 

712 mapping = {'sourceId': 'id', 

713 'x': 'centroid_x', 

714 'y': 'centroid_y', 

715 'xErr': 'centroid_xErr', 

716 'yErr': 'centroid_yErr', 

717 'Ixx': 'shape_xx', 

718 'Iyy': 'shape_yy', 

719 'Ixy': 'shape_xy', 

720 f'{self.config.sourceFluxType}_instFlux': 'flux_instFlux', 

721 f'{self.config.sourceFluxType}_instFluxErr': 'flux_instFluxErr', 

722 } 

723 catalog = lsst.afw.table.SourceCatalog(table) 

724 matched = visitCatalog['ccd'] == detectorId 

725 catalog.resize(sum(matched)) 

726 view = visitCatalog.loc[matched] 

727 for dfCol, afwCol in mapping.items(): 

728 catalog[afwCol] = view[dfCol] 

729 

730 self.log.debug("%d sources selected in visit %d detector %d", 

731 len(catalog), 

732 view['visit'].iloc[0], # all visits in this catalog are the same, so take the first 

733 detectorId) 

734 return catalog 

735 

736 def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations, 

737 jointcalControl, camera): 

738 """Read the data that jointcal needs to run. (Gen3 version) 

739 

740 Modifies ``associations`` in-place with the loaded data. 

741 

742 Parameters 

743 ---------- 

744 inputSourceTableVisit : `list` [`lsst.daf.butler.DeferredDatasetHandle`] 

745 References to visit-level DataFrames to load the catalog data from. 

746 inputVisitSummary : `list` [`lsst.daf.butler.DeferredDatasetHandle`] 

747 Visit-level exposure summary catalog with metadata. 

748 associations : `lsst.jointcal.Associations` 

749 Object to add the loaded data to by constructing new CcdImages. 

750 jointcalControl : `jointcal.JointcalControl` 

751 Control object for C++ associations management. 

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

753 Camera object for detector geometry. 

754 

755 Returns 

756 ------- 

757 oldWcsList: `list` [`lsst.afw.geom.SkyWcs`] 

758 The original WCS of the input data, to aid in writing tests. 

759 bands : `list` [`str`] 

760 The filter bands of each input dataset. 

761 """ 

762 oldWcsList = [] 

763 bands = [] 

764 load_cat_prof_file = 'jointcal_load_data.prof' if self.config.detailedProfile else '' 

765 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

766 table = self._make_schema_table() # every detector catalog has the same layout 

767 # No guarantee that the input is in the same order of visits, so we have to map one of them. 

768 catalogMap = {ref.dataId['visit']: i for i, ref in enumerate(inputSourceTableVisit)} 

769 detectorDict = {detector.getId(): detector for detector in camera} 

770 

771 for visitSummaryRef in inputVisitSummary: 

772 visitSummary = visitSummaryRef.get() 

773 visitCatalog = inputSourceTableVisit[catalogMap[visitSummaryRef.dataId['visit']]].get() 

774 selected = self.sourceSelector.run(visitCatalog) 

775 

776 # Build a CcdImage for each detector in this visit. 

777 detectors = {id: index for index, id in enumerate(visitSummary['id'])} 

778 for id, index in detectors.items(): 

779 catalog = self._extract_detector_catalog_from_visit_catalog(table, selected.sourceCat, id) 

780 data = self._make_one_input_data(visitSummary[index], catalog, detectorDict) 

781 result = self._build_ccdImage(data, associations, jointcalControl) 

782 if result is not None: 

783 oldWcsList.append(result.wcs) 

784 

785 # A visit has only one band, so we can just use the first. 

786 bands.append(visitSummary[0]['band']) 

787 bands = collections.Counter(bands) 

788 

789 return oldWcsList, bands 

790 

791 def _make_one_input_data(self, visitRecord, catalog, detectorDict): 

792 """Return a data structure for this detector+visit.""" 

793 return JointcalInputData(visit=visitRecord['visit'], 

794 catalog=catalog, 

795 visitInfo=visitRecord.getVisitInfo(), 

796 detector=detectorDict[visitRecord.getId()], 

797 photoCalib=visitRecord.getPhotoCalib(), 

798 wcs=visitRecord.getWcs(), 

799 bbox=visitRecord.getBBox(), 

800 # ExposureRecord doesn't have a FilterLabel yet, 

801 # so we have to make one. 

802 filter=lsst.afw.image.FilterLabel(band=visitRecord['band'], 

803 physical=visitRecord['physical_filter'])) 

804 

805 # We don't currently need to persist the metadata. 

806 # If we do in the future, we will have to add appropriate dataset templates 

807 # to each obs package (the metadata template should look like `jointcal_wcs`). 

808 def _getMetadataName(self): 

809 return None 

810 

811 @classmethod 

812 def _makeArgumentParser(cls): 

813 """Create an argument parser""" 

814 parser = pipeBase.ArgumentParser(name=cls._DefaultName) 

815 parser.add_id_argument("--id", "calexp", help="data ID, e.g. --id visit=6789 ccd=0..9", 

816 ContainerClass=PerTractCcdDataIdContainer) 

817 return parser 

818 

819 def _build_ccdImage(self, data, associations, jointcalControl): 

820 """ 

821 Extract the necessary things from this dataRef to add a new ccdImage. 

822 

823 Parameters 

824 ---------- 

825 data : `JointcalInputData` 

826 The loaded input data. 

827 associations : `lsst.jointcal.Associations` 

828 Object to add the info to, to construct a new CcdImage 

829 jointcalControl : `jointcal.JointcalControl` 

830 Control object for associations management 

831 

832 Returns 

833 ------ 

834 namedtuple or `None` 

835 ``wcs`` 

836 The TAN WCS of this image, read from the calexp 

837 (`lsst.afw.geom.SkyWcs`). 

838 ``key`` 

839 A key to identify this dataRef by its visit and ccd ids 

840 (`namedtuple`). 

841 ``band`` 

842 This calexp's filter band (`str`) (used to e.g. load refcats) 

843 `None` 

844 if there are no sources in the loaded catalog. 

845 """ 

846 if len(data.catalog) == 0: 

847 self.log.warn("No sources selected in visit %s ccd %s", data.visit, data.detector.getId()) 

848 return None 

849 

850 associations.createCcdImage(data.catalog, 

851 data.wcs, 

852 data.visitInfo, 

853 data.bbox, 

854 data.filter.physicalLabel, 

855 data.photoCalib, 

856 data.detector, 

857 data.visit, 

858 data.detector.getId(), 

859 jointcalControl) 

860 

861 Result = collections.namedtuple('Result_from_build_CcdImage', ('wcs', 'key', 'band')) 

862 Key = collections.namedtuple('Key', ('visit', 'ccd')) 

863 return Result(data.wcs, Key(data.visit, data.detector.getId()), data.filter.bandLabel) 

864 

865 def _readDataId(self, butler, dataId): 

866 """Read all of the data for one dataId from the butler. (gen2 version)""" 

867 # Not all instruments have `visit` in their dataIds. 

868 if "visit" in dataId.keys(): 

869 visit = dataId["visit"] 

870 else: 

871 visit = butler.getButler().queryMetadata("calexp", ("visit"), butler.dataId)[0] 

872 detector = butler.get('calexp_detector', dataId=dataId) 

873 

874 catalog = butler.get('src', 

875 flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS, 

876 dataId=dataId) 

877 goodSrc = self.sourceSelector.run(catalog) 

878 self.log.debug("%d sources selected in visit %d detector %d", 

879 len(goodSrc.sourceCat), 

880 visit, 

881 detector.getId()) 

882 return JointcalInputData(visit=visit, 

883 catalog=goodSrc.sourceCat, 

884 visitInfo=butler.get('calexp_visitInfo', dataId=dataId), 

885 detector=detector, 

886 photoCalib=butler.get('calexp_photoCalib', dataId=dataId), 

887 wcs=butler.get('calexp_wcs', dataId=dataId), 

888 bbox=butler.get('calexp_bbox', dataId=dataId), 

889 filter=butler.get('calexp_filterLabel', dataId=dataId)) 

890 

891 def loadData(self, dataRefs, associations, jointcalControl): 

892 """Read the data that jointcal needs to run. (Gen2 version)""" 

893 visit_ccd_to_dataRef = {} 

894 oldWcsList = [] 

895 bands = [] 

896 load_cat_prof_file = 'jointcal_loadData.prof' if self.config.detailedProfile else '' 

897 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

898 # Need the bounding-box of the focal plane (the same for all visits) for photometry visit models 

899 camera = dataRefs[0].get('camera', immediate=True) 

900 self.focalPlaneBBox = camera.getFpBBox() 

901 for dataRef in dataRefs: 

902 data = self._readDataId(dataRef.getButler(), dataRef.dataId) 

903 result = self._build_ccdImage(data, associations, jointcalControl) 

904 if result is None: 

905 continue 

906 oldWcsList.append(result.wcs) 

907 visit_ccd_to_dataRef[result.key] = dataRef 

908 bands.append(result.band) 

909 if len(bands) == 0: 

910 raise RuntimeError("No data to process: did source selector remove all sources?") 

911 bands = collections.Counter(bands) 

912 

913 return oldWcsList, bands, visit_ccd_to_dataRef 

914 

915 def _getDebugPath(self, filename): 

916 """Constructs a path to filename using the configured debug path. 

917 """ 

918 return os.path.join(self.config.debugOutputPath, filename) 

919 

920 def _prep_sky(self, associations, bands): 

921 """Prepare on-sky and other data that must be computed after data has 

922 been read. 

923 """ 

924 associations.computeCommonTangentPoint() 

925 

926 boundingCircle = associations.computeBoundingCircle() 

927 center = lsst.geom.SpherePoint(boundingCircle.getCenter()) 

928 radius = lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians) 

929 

930 self.log.info(f"Data has center={center} with radius={radius.asDegrees()} degrees.") 

931 

932 # Determine a default filter band associated with the catalog. See DM-9093 

933 defaultBand = bands.most_common(1)[0][0] 

934 self.log.debug("Using '%s' filter band for reference flux", defaultBand) 

935 

936 return boundingCircle, center, radius, defaultBand 

937 

938 @pipeBase.timeMethod 

939 def runDataRef(self, dataRefs): 

940 """ 

941 Jointly calibrate the astrometry and photometry across a set of images. 

942 

943 NOTE: this is for gen2 middleware only. 

944 

945 Parameters 

946 ---------- 

947 dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef` 

948 List of data references to the exposures to be fit. 

949 

950 Returns 

951 ------- 

952 result : `lsst.pipe.base.Struct` 

953 Struct of metadata from the fit, containing: 

954 

955 ``dataRefs`` 

956 The provided data references that were fit (with updated WCSs) 

957 ``oldWcsList`` 

958 The original WCS from each dataRef 

959 ``metrics`` 

960 Dictionary of internally-computed metrics for testing/validation. 

961 """ 

962 if len(dataRefs) == 0: 

963 raise ValueError('Need a non-empty list of data references!') 

964 

965 exitStatus = 0 # exit status for shell 

966 

967 sourceFluxField = "slot_%sFlux" % (self.config.sourceFluxType,) 

968 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField) 

969 associations = lsst.jointcal.Associations() 

970 

971 oldWcsList, bands, visit_ccd_to_dataRef = self.loadData(dataRefs, 

972 associations, 

973 jointcalControl) 

974 

975 boundingCircle, center, radius, defaultBand = self._prep_sky(associations, bands) 

976 epoch = self._compute_proper_motion_epoch(associations.getCcdImageList()) 

977 

978 tract = dataRefs[0].dataId['tract'] 

979 

980 if self.config.doAstrometry: 

981 astrometry = self._do_load_refcat_and_fit(associations, defaultBand, center, radius, 

982 name="astrometry", 

983 refObjLoader=self.astrometryRefObjLoader, 

984 referenceSelector=self.astrometryReferenceSelector, 

985 fit_function=self._fit_astrometry, 

986 tract=tract, 

987 epoch=epoch) 

988 self._write_astrometry_results(associations, astrometry.model, visit_ccd_to_dataRef) 

989 else: 

990 astrometry = Astrometry(None, None, None) 

991 

992 if self.config.doPhotometry: 

993 photometry = self._do_load_refcat_and_fit(associations, defaultBand, center, radius, 

994 name="photometry", 

995 refObjLoader=self.photometryRefObjLoader, 

996 referenceSelector=self.photometryReferenceSelector, 

997 fit_function=self._fit_photometry, 

998 tract=tract, 

999 epoch=epoch, 

1000 reject_bad_fluxes=True) 

1001 self._write_photometry_results(associations, photometry.model, visit_ccd_to_dataRef) 

1002 else: 

1003 photometry = Photometry(None, None) 

1004 

1005 return pipeBase.Struct(dataRefs=dataRefs, 

1006 oldWcsList=oldWcsList, 

1007 job=self.job, 

1008 astrometryRefObjLoader=self.astrometryRefObjLoader, 

1009 photometryRefObjLoader=self.photometryRefObjLoader, 

1010 defaultBand=defaultBand, 

1011 epoch=epoch, 

1012 exitStatus=exitStatus) 

1013 

1014 def _get_refcat_coordinate_error_override(self, refCat, name): 

1015 """Check whether we should override the refcat coordinate errors, and 

1016 return the overridden error if necessary. 

1017 

1018 Parameters 

1019 ---------- 

1020 refCat : `lsst.afw.table.SimpleCatalog` 

1021 The reference catalog to check for a ``coord_raErr`` field. 

1022 name : `str` 

1023 Whether we are doing "astrometry" or "photometry". 

1024 

1025 Returns 

1026 ------- 

1027 refCoordErr : `float` 

1028 The refcat coordinate error to use, or NaN if we are not overriding 

1029 those fields. 

1030 

1031 Raises 

1032 ------ 

1033 lsst.pex.config.FieldValidationError 

1034 Raised if the refcat does not contain coordinate errors and 

1035 ``config.astrometryReferenceErr`` is not set. 

1036 """ 

1037 # This value doesn't matter for photometry, so just set something to 

1038 # keep old refcats from causing problems. 

1039 if name.lower() == "photometry": 

1040 if 'coord_raErr' not in refCat.schema: 

1041 return 100 

1042 else: 

1043 return float('nan') 

1044 

1045 if self.config.astrometryReferenceErr is None and 'coord_raErr' not in refCat.schema: 

1046 msg = ("Reference catalog does not contain coordinate errors, " 

1047 "and config.astrometryReferenceErr not supplied.") 

1048 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr, 

1049 self.config, 

1050 msg) 

1051 

1052 if self.config.astrometryReferenceErr is not None and 'coord_raErr' in refCat.schema: 

1053 self.log.warn("Overriding reference catalog coordinate errors with %f/coordinate [mas]", 

1054 self.config.astrometryReferenceErr) 

1055 

1056 if self.config.astrometryReferenceErr is None: 

1057 return float('nan') 

1058 else: 

1059 return self.config.astrometryReferenceErr 

1060 

1061 def _compute_proper_motion_epoch(self, ccdImageList): 

1062 """Return the proper motion correction epoch of the provided images. 

1063 

1064 Parameters 

1065 ---------- 

1066 ccdImageList : `list` [`lsst.jointcal.CcdImage`] 

1067 The images to compute the appropriate epoch for. 

1068 

1069 Returns 

1070 ------- 

1071 epoch : `astropy.time.Time` 

1072 The date to use for proper motion corrections. 

1073 """ 

1074 mjds = [ccdImage.getMjd() for ccdImage in ccdImageList] 

1075 return astropy.time.Time(np.mean(mjds), format='mjd', scale="tai") 

1076 

1077 def _do_load_refcat_and_fit(self, associations, defaultBand, center, radius, 

1078 tract="", match_cut=3.0, 

1079 reject_bad_fluxes=False, *, 

1080 name="", refObjLoader=None, referenceSelector=None, 

1081 fit_function=None, epoch=None): 

1082 """Load reference catalog, perform the fit, and return the result. 

1083 

1084 Parameters 

1085 ---------- 

1086 associations : `lsst.jointcal.Associations` 

1087 The star/reference star associations to fit. 

1088 defaultBand : `str` 

1089 filter to load from reference catalog. 

1090 center : `lsst.geom.SpherePoint` 

1091 ICRS center of field to load from reference catalog. 

1092 radius : `lsst.geom.Angle` 

1093 On-sky radius to load from reference catalog. 

1094 name : `str` 

1095 Name of thing being fit: "astrometry" or "photometry". 

1096 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask` 

1097 Reference object loader to use to load a reference catalog. 

1098 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask` 

1099 Selector to use to pick objects from the loaded reference catalog. 

1100 fit_function : callable 

1101 Function to call to perform fit (takes Associations object). 

1102 tract : `str`, optional 

1103 Name of tract currently being fit. 

1104 match_cut : `float`, optional 

1105 Radius in arcseconds to find cross-catalog matches to during 

1106 associations.associateCatalogs. 

1107 reject_bad_fluxes : `bool`, optional 

1108 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr. 

1109 epoch : `astropy.time.Time`, optional 

1110 Epoch to which to correct refcat proper motion and parallax, 

1111 or `None` to not apply such corrections. 

1112 

1113 Returns 

1114 ------- 

1115 result : `Photometry` or `Astrometry` 

1116 Result of `fit_function()` 

1117 """ 

1118 self.log.info("====== Now processing %s...", name) 

1119 # TODO: this should not print "trying to invert a singular transformation:" 

1120 # if it does that, something's not right about the WCS... 

1121 associations.associateCatalogs(match_cut) 

1122 add_measurement(self.job, 'jointcal.associated_%s_fittedStars' % name, 

1123 associations.fittedStarListSize()) 

1124 

1125 applyColorterms = False if name.lower() == "astrometry" else self.config.applyColorTerms 

1126 refCat, fluxField = self._load_reference_catalog(refObjLoader, referenceSelector, 

1127 center, radius, defaultBand, 

1128 applyColorterms=applyColorterms, 

1129 epoch=epoch) 

1130 refCoordErr = self._get_refcat_coordinate_error_override(refCat, name) 

1131 

1132 associations.collectRefStars(refCat, 

1133 self.config.matchCut*lsst.geom.arcseconds, 

1134 fluxField, 

1135 refCoordinateErr=refCoordErr, 

1136 rejectBadFluxes=reject_bad_fluxes) 

1137 add_measurement(self.job, 'jointcal.collected_%s_refStars' % name, 

1138 associations.refStarListSize()) 

1139 

1140 associations.prepareFittedStars(self.config.minMeasurements) 

1141 

1142 self._check_star_lists(associations, name) 

1143 add_measurement(self.job, 'jointcal.selected_%s_refStars' % name, 

1144 associations.nFittedStarsWithAssociatedRefStar()) 

1145 add_measurement(self.job, 'jointcal.selected_%s_fittedStars' % name, 

1146 associations.fittedStarListSize()) 

1147 add_measurement(self.job, 'jointcal.selected_%s_ccdImages' % name, 

1148 associations.nCcdImagesValidForFit()) 

1149 

1150 load_cat_prof_file = 'jointcal_fit_%s.prof'%name if self.config.detailedProfile else '' 

1151 dataName = "{}_{}".format(tract, defaultBand) 

1152 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

1153 result = fit_function(associations, dataName) 

1154 # TODO DM-12446: turn this into a "butler save" somehow. 

1155 # Save reference and measurement chi2 contributions for this data 

1156 if self.config.writeChi2FilesInitialFinal: 

1157 baseName = self._getDebugPath(f"{name}_final_chi2-{dataName}") 

1158 result.fit.saveChi2Contributions(baseName+"{type}") 

1159 self.log.info("Wrote chi2 contributions files: %s", baseName) 

1160 

1161 return result 

1162 

1163 def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName, 

1164 applyColorterms=False, epoch=None): 

1165 """Load the necessary reference catalog sources, convert fluxes to 

1166 correct units, and apply color term corrections if requested. 

1167 

1168 Parameters 

1169 ---------- 

1170 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask` 

1171 The reference catalog loader to use to get the data. 

1172 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask` 

1173 Source selector to apply to loaded reference catalog. 

1174 center : `lsst.geom.SpherePoint` 

1175 The center around which to load sources. 

1176 radius : `lsst.geom.Angle` 

1177 The radius around ``center`` to load sources in. 

1178 filterName : `str` 

1179 The name of the camera filter to load fluxes for. 

1180 applyColorterms : `bool` 

1181 Apply colorterm corrections to the refcat for ``filterName``? 

1182 epoch : `astropy.time.Time`, optional 

1183 Epoch to which to correct refcat proper motion and parallax, 

1184 or `None` to not apply such corrections. 

1185 

1186 Returns 

1187 ------- 

1188 refCat : `lsst.afw.table.SimpleCatalog` 

1189 The loaded reference catalog. 

1190 fluxField : `str` 

1191 The name of the reference catalog flux field appropriate for ``filterName``. 

1192 """ 

1193 skyCircle = refObjLoader.loadSkyCircle(center, 

1194 radius, 

1195 filterName, 

1196 epoch=epoch) 

1197 

1198 selected = referenceSelector.run(skyCircle.refCat) 

1199 # Need memory contiguity to get reference filters as a vector. 

1200 if not selected.sourceCat.isContiguous(): 

1201 refCat = selected.sourceCat.copy(deep=True) 

1202 else: 

1203 refCat = selected.sourceCat 

1204 

1205 if applyColorterms: 

1206 refCatName = refObjLoader.ref_dataset_name 

1207 self.log.info("Applying color terms for filterName=%r reference catalog=%s", 

1208 filterName, refCatName) 

1209 colorterm = self.config.colorterms.getColorterm( 

1210 filterName=filterName, photoCatName=refCatName, doRaise=True) 

1211 

1212 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat, filterName) 

1213 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy) 

1214 # TODO: I didn't want to use this, but I'll deal with it in DM-16903 

1215 refCat[skyCircle.fluxField+'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9 

1216 

1217 return refCat, skyCircle.fluxField 

1218 

1219 def _check_star_lists(self, associations, name): 

1220 # TODO: these should be len(blah), but we need this properly wrapped first. 

1221 if associations.nCcdImagesValidForFit() == 0: 

1222 raise RuntimeError('No images in the ccdImageList!') 

1223 if associations.fittedStarListSize() == 0: 

1224 raise RuntimeError('No stars in the {} fittedStarList!'.format(name)) 

1225 if associations.refStarListSize() == 0: 

1226 raise RuntimeError('No stars in the {} reference star list!'.format(name)) 

1227 

1228 def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None): 

1229 """Compute chi2, log it, validate the model, and return chi2. 

1230 

1231 Parameters 

1232 ---------- 

1233 associations : `lsst.jointcal.Associations` 

1234 The star/reference star associations to fit. 

1235 fit : `lsst.jointcal.FitterBase` 

1236 The fitter to use for minimization. 

1237 model : `lsst.jointcal.Model` 

1238 The model being fit. 

1239 chi2Label : `str` 

1240 Label to describe the chi2 (e.g. "Initialized", "Final"). 

1241 writeChi2Name : `str`, optional 

1242 Filename prefix to write the chi2 contributions to. 

1243 Do not supply an extension: an appropriate one will be added. 

1244 

1245 Returns 

1246 ------- 

1247 chi2: `lsst.jointcal.Chi2Accumulator` 

1248 The chi2 object for the current fitter and model. 

1249 

1250 Raises 

1251 ------ 

1252 FloatingPointError 

1253 Raised if chi2 is infinite or NaN. 

1254 ValueError 

1255 Raised if the model is not valid. 

1256 """ 

1257 if writeChi2Name is not None: 

1258 fullpath = self._getDebugPath(writeChi2Name) 

1259 fit.saveChi2Contributions(fullpath+"{type}") 

1260 self.log.info("Wrote chi2 contributions files: %s", fullpath) 

1261 

1262 chi2 = fit.computeChi2() 

1263 self.log.info("%s %s", chi2Label, chi2) 

1264 self._check_stars(associations) 

1265 if not np.isfinite(chi2.chi2): 

1266 raise FloatingPointError(f'{chi2Label} chi2 is invalid: {chi2}') 

1267 if not model.validate(associations.getCcdImageList(), chi2.ndof): 

1268 raise ValueError("Model is not valid: check log messages for warnings.") 

1269 return chi2 

1270 

1271 def _fit_photometry(self, associations, dataName=None): 

1272 """ 

1273 Fit the photometric data. 

1274 

1275 Parameters 

1276 ---------- 

1277 associations : `lsst.jointcal.Associations` 

1278 The star/reference star associations to fit. 

1279 dataName : `str` 

1280 Name of the data being processed (e.g. "1234_HSC-Y"), for 

1281 identifying debugging files. 

1282 

1283 Returns 

1284 ------- 

1285 fit_result : `namedtuple` 

1286 fit : `lsst.jointcal.PhotometryFit` 

1287 The photometric fitter used to perform the fit. 

1288 model : `lsst.jointcal.PhotometryModel` 

1289 The photometric model that was fit. 

1290 """ 

1291 self.log.info("=== Starting photometric fitting...") 

1292 

1293 # TODO: should use pex.config.RegistryField here (see DM-9195) 

1294 if self.config.photometryModel == "constrainedFlux": 

1295 model = lsst.jointcal.ConstrainedFluxModel(associations.getCcdImageList(), 

1296 self.focalPlaneBBox, 

1297 visitOrder=self.config.photometryVisitOrder, 

1298 errorPedestal=self.config.photometryErrorPedestal) 

1299 # potentially nonlinear problem, so we may need a line search to converge. 

1300 doLineSearch = self.config.allowLineSearch 

1301 elif self.config.photometryModel == "constrainedMagnitude": 

1302 model = lsst.jointcal.ConstrainedMagnitudeModel(associations.getCcdImageList(), 

1303 self.focalPlaneBBox, 

1304 visitOrder=self.config.photometryVisitOrder, 

1305 errorPedestal=self.config.photometryErrorPedestal) 

1306 # potentially nonlinear problem, so we may need a line search to converge. 

1307 doLineSearch = self.config.allowLineSearch 

1308 elif self.config.photometryModel == "simpleFlux": 

1309 model = lsst.jointcal.SimpleFluxModel(associations.getCcdImageList(), 

1310 errorPedestal=self.config.photometryErrorPedestal) 

1311 doLineSearch = False # purely linear in model parameters, so no line search needed 

1312 elif self.config.photometryModel == "simpleMagnitude": 

1313 model = lsst.jointcal.SimpleMagnitudeModel(associations.getCcdImageList(), 

1314 errorPedestal=self.config.photometryErrorPedestal) 

1315 doLineSearch = False # purely linear in model parameters, so no line search needed 

1316 

1317 fit = lsst.jointcal.PhotometryFit(associations, model) 

1318 # TODO DM-12446: turn this into a "butler save" somehow. 

1319 # Save reference and measurement chi2 contributions for this data 

1320 if self.config.writeChi2FilesInitialFinal: 

1321 baseName = f"photometry_initial_chi2-{dataName}" 

1322 else: 

1323 baseName = None 

1324 if self.config.writeInitialModel: 

1325 fullpath = self._getDebugPath("initialPhotometryModel.txt") 

1326 writeModel(model, fullpath, self.log) 

1327 self._logChi2AndValidate(associations, fit, model, "Initialized", writeChi2Name=baseName) 

1328 

1329 def getChi2Name(whatToFit): 

1330 if self.config.writeChi2FilesOuterLoop: 

1331 return f"photometry_init-%s_chi2-{dataName}" % whatToFit 

1332 else: 

1333 return None 

1334 

1335 # The constrained model needs the visit transform fit first; the chip 

1336 # transform is initialized from the singleFrame PhotoCalib, so it's close. 

1337 dumpMatrixFile = self._getDebugPath("photometry_preinit") if self.config.writeInitMatrix else "" 

1338 if self.config.photometryModel.startswith("constrained"): 

1339 # no line search: should be purely (or nearly) linear, 

1340 # and we want a large step size to initialize with. 

1341 fit.minimize("ModelVisit", dumpMatrixFile=dumpMatrixFile) 

1342 self._logChi2AndValidate(associations, fit, model, "Initialize ModelVisit", 

1343 writeChi2Name=getChi2Name("ModelVisit")) 

1344 dumpMatrixFile = "" # so we don't redo the output on the next step 

1345 

1346 fit.minimize("Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile) 

1347 self._logChi2AndValidate(associations, fit, model, "Initialize Model", 

1348 writeChi2Name=getChi2Name("Model")) 

1349 

1350 fit.minimize("Fluxes") # no line search: always purely linear. 

1351 self._logChi2AndValidate(associations, fit, model, "Initialize Fluxes", 

1352 writeChi2Name=getChi2Name("Fluxes")) 

1353 

1354 fit.minimize("Model Fluxes", doLineSearch=doLineSearch) 

1355 self._logChi2AndValidate(associations, fit, model, "Initialize ModelFluxes", 

1356 writeChi2Name=getChi2Name("ModelFluxes")) 

1357 

1358 model.freezeErrorTransform() 

1359 self.log.debug("Photometry error scales are frozen.") 

1360 

1361 chi2 = self._iterate_fit(associations, 

1362 fit, 

1363 self.config.maxPhotometrySteps, 

1364 "photometry", 

1365 "Model Fluxes", 

1366 doRankUpdate=self.config.photometryDoRankUpdate, 

1367 doLineSearch=doLineSearch, 

1368 dataName=dataName) 

1369 

1370 add_measurement(self.job, 'jointcal.photometry_final_chi2', chi2.chi2) 

1371 add_measurement(self.job, 'jointcal.photometry_final_ndof', chi2.ndof) 

1372 return Photometry(fit, model) 

1373 

1374 def _fit_astrometry(self, associations, dataName=None): 

1375 """ 

1376 Fit the astrometric data. 

1377 

1378 Parameters 

1379 ---------- 

1380 associations : `lsst.jointcal.Associations` 

1381 The star/reference star associations to fit. 

1382 dataName : `str` 

1383 Name of the data being processed (e.g. "1234_HSC-Y"), for 

1384 identifying debugging files. 

1385 

1386 Returns 

1387 ------- 

1388 fit_result : `namedtuple` 

1389 fit : `lsst.jointcal.AstrometryFit` 

1390 The astrometric fitter used to perform the fit. 

1391 model : `lsst.jointcal.AstrometryModel` 

1392 The astrometric model that was fit. 

1393 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler` 

1394 The model for the sky to tangent plane projection that was used in the fit. 

1395 """ 

1396 

1397 self.log.info("=== Starting astrometric fitting...") 

1398 

1399 associations.deprojectFittedStars() 

1400 

1401 # NOTE: need to return sky_to_tan_projection so that it doesn't get garbage collected. 

1402 # TODO: could we package sky_to_tan_projection and model together so we don't have to manage 

1403 # them so carefully? 

1404 sky_to_tan_projection = lsst.jointcal.OneTPPerVisitHandler(associations.getCcdImageList()) 

1405 

1406 if self.config.astrometryModel == "constrained": 

1407 model = lsst.jointcal.ConstrainedAstrometryModel(associations.getCcdImageList(), 

1408 sky_to_tan_projection, 

1409 chipOrder=self.config.astrometryChipOrder, 

1410 visitOrder=self.config.astrometryVisitOrder) 

1411 elif self.config.astrometryModel == "simple": 

1412 model = lsst.jointcal.SimpleAstrometryModel(associations.getCcdImageList(), 

1413 sky_to_tan_projection, 

1414 self.config.useInputWcs, 

1415 nNotFit=0, 

1416 order=self.config.astrometrySimpleOrder) 

1417 

1418 fit = lsst.jointcal.AstrometryFit(associations, model, self.config.positionErrorPedestal) 

1419 # TODO DM-12446: turn this into a "butler save" somehow. 

1420 # Save reference and measurement chi2 contributions for this data 

1421 if self.config.writeChi2FilesInitialFinal: 

1422 baseName = f"astrometry_initial_chi2-{dataName}" 

1423 else: 

1424 baseName = None 

1425 if self.config.writeInitialModel: 

1426 fullpath = self._getDebugPath("initialAstrometryModel.txt") 

1427 writeModel(model, fullpath, self.log) 

1428 self._logChi2AndValidate(associations, fit, model, "Initial", writeChi2Name=baseName) 

1429 

1430 def getChi2Name(whatToFit): 

1431 if self.config.writeChi2FilesOuterLoop: 

1432 return f"astrometry_init-%s_chi2-{dataName}" % whatToFit 

1433 else: 

1434 return None 

1435 

1436 dumpMatrixFile = self._getDebugPath("astrometry_preinit") if self.config.writeInitMatrix else "" 

1437 # The constrained model needs the visit transform fit first; the chip 

1438 # transform is initialized from the detector's cameraGeom, so it's close. 

1439 if self.config.astrometryModel == "constrained": 

1440 fit.minimize("DistortionsVisit", dumpMatrixFile=dumpMatrixFile) 

1441 self._logChi2AndValidate(associations, fit, model, "Initialize DistortionsVisit", 

1442 writeChi2Name=getChi2Name("DistortionsVisit")) 

1443 dumpMatrixFile = "" # so we don't redo the output on the next step 

1444 

1445 fit.minimize("Distortions", dumpMatrixFile=dumpMatrixFile) 

1446 self._logChi2AndValidate(associations, fit, model, "Initialize Distortions", 

1447 writeChi2Name=getChi2Name("Distortions")) 

1448 

1449 fit.minimize("Positions") 

1450 self._logChi2AndValidate(associations, fit, model, "Initialize Positions", 

1451 writeChi2Name=getChi2Name("Positions")) 

1452 

1453 fit.minimize("Distortions Positions") 

1454 self._logChi2AndValidate(associations, fit, model, "Initialize DistortionsPositions", 

1455 writeChi2Name=getChi2Name("DistortionsPositions")) 

1456 

1457 chi2 = self._iterate_fit(associations, 

1458 fit, 

1459 self.config.maxAstrometrySteps, 

1460 "astrometry", 

1461 "Distortions Positions", 

1462 doRankUpdate=self.config.astrometryDoRankUpdate, 

1463 dataName=dataName) 

1464 

1465 add_measurement(self.job, 'jointcal.astrometry_final_chi2', chi2.chi2) 

1466 add_measurement(self.job, 'jointcal.astrometry_final_ndof', chi2.ndof) 

1467 

1468 return Astrometry(fit, model, sky_to_tan_projection) 

1469 

1470 def _check_stars(self, associations): 

1471 """Count measured and reference stars per ccd and warn/log them.""" 

1472 for ccdImage in associations.getCcdImageList(): 

1473 nMeasuredStars, nRefStars = ccdImage.countStars() 

1474 self.log.debug("ccdImage %s has %s measured and %s reference stars", 

1475 ccdImage.getName(), nMeasuredStars, nRefStars) 

1476 if nMeasuredStars < self.config.minMeasuredStarsPerCcd: 

1477 self.log.warn("ccdImage %s has only %s measuredStars (desired %s)", 

1478 ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd) 

1479 if nRefStars < self.config.minRefStarsPerCcd: 

1480 self.log.warn("ccdImage %s has only %s RefStars (desired %s)", 

1481 ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd) 

1482 

1483 def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, 

1484 dataName="", 

1485 doRankUpdate=True, 

1486 doLineSearch=False): 

1487 """Run fitter.minimize up to max_steps times, returning the final chi2. 

1488 

1489 Parameters 

1490 ---------- 

1491 associations : `lsst.jointcal.Associations` 

1492 The star/reference star associations to fit. 

1493 fitter : `lsst.jointcal.FitterBase` 

1494 The fitter to use for minimization. 

1495 max_steps : `int` 

1496 Maximum number of steps to run outlier rejection before declaring 

1497 convergence failure. 

1498 name : {'photometry' or 'astrometry'} 

1499 What type of data are we fitting (for logs and debugging files). 

1500 whatToFit : `str` 

1501 Passed to ``fitter.minimize()`` to define the parameters to fit. 

1502 dataName : `str`, optional 

1503 Descriptive name for this dataset (e.g. tract and filter), 

1504 for debugging. 

1505 doRankUpdate : `bool`, optional 

1506 Do an Eigen rank update during minimization, or recompute the full 

1507 matrix and gradient? 

1508 doLineSearch : `bool`, optional 

1509 Do a line search for the optimum step during minimization? 

1510 

1511 Returns 

1512 ------- 

1513 chi2: `lsst.jointcal.Chi2Statistic` 

1514 The final chi2 after the fit converges, or is forced to end. 

1515 

1516 Raises 

1517 ------ 

1518 FloatingPointError 

1519 Raised if the fitter fails with a non-finite value. 

1520 RuntimeError 

1521 Raised if the fitter fails for some other reason; 

1522 log messages will provide further details. 

1523 """ 

1524 dumpMatrixFile = self._getDebugPath(f"{name}_postinit") if self.config.writeInitMatrix else "" 

1525 oldChi2 = lsst.jointcal.Chi2Statistic() 

1526 oldChi2.chi2 = float("inf") 

1527 for i in range(max_steps): 

1528 if self.config.writeChi2FilesOuterLoop: 

1529 writeChi2Name = f"{name}_iterate_{i}_chi2-{dataName}" 

1530 else: 

1531 writeChi2Name = None 

1532 result = fitter.minimize(whatToFit, 

1533 self.config.outlierRejectSigma, 

1534 doRankUpdate=doRankUpdate, 

1535 doLineSearch=doLineSearch, 

1536 dumpMatrixFile=dumpMatrixFile) 

1537 dumpMatrixFile = "" # clear it so we don't write the matrix again. 

1538 chi2 = self._logChi2AndValidate(associations, fitter, fitter.getModel(), 

1539 f"Fit iteration {i}", writeChi2Name=writeChi2Name) 

1540 

1541 if result == MinimizeResult.Converged: 

1542 if doRankUpdate: 

1543 self.log.debug("fit has converged - no more outliers - redo minimization " 

1544 "one more time in case we have lost accuracy in rank update.") 

1545 # Redo minimization one more time in case we have lost accuracy in rank update 

1546 result = fitter.minimize(whatToFit, self.config.outlierRejectSigma) 

1547 chi2 = self._logChi2AndValidate(associations, fitter, fitter.getModel(), "Fit completed") 

1548 

1549 # log a message for a large final chi2, TODO: DM-15247 for something better 

1550 if chi2.chi2/chi2.ndof >= 4.0: 

1551 self.log.error("Potentially bad fit: High chi-squared/ndof.") 

1552 

1553 break 

1554 elif result == MinimizeResult.Chi2Increased: 

1555 self.log.warn("Still some outliers remaining but chi2 increased - retry") 

1556 # Check whether the increase was large enough to cause trouble. 

1557 chi2Ratio = chi2.chi2 / oldChi2.chi2 

1558 if chi2Ratio > 1.5: 

1559 self.log.warn('Significant chi2 increase by a factor of %.4g / %.4g = %.4g', 

1560 chi2.chi2, oldChi2.chi2, chi2Ratio) 

1561 # Based on a variety of HSC jointcal logs (see DM-25779), it 

1562 # appears that chi2 increases more than a factor of ~2 always 

1563 # result in the fit diverging rapidly and ending at chi2 > 1e10. 

1564 # Using 10 as the "failure" threshold gives some room between 

1565 # leaving a warning and bailing early. 

1566 if chi2Ratio > 10: 

1567 msg = ("Large chi2 increase between steps: fit likely cannot converge." 

1568 " Try setting one or more of the `writeChi2*` config fields and looking" 

1569 " at how individual star chi2-values evolve during the fit.") 

1570 raise RuntimeError(msg) 

1571 oldChi2 = chi2 

1572 elif result == MinimizeResult.NonFinite: 

1573 filename = self._getDebugPath("{}_failure-nonfinite_chi2-{}.csv".format(name, dataName)) 

1574 # TODO DM-12446: turn this into a "butler save" somehow. 

1575 fitter.saveChi2Contributions(filename+"{type}") 

1576 msg = "Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}" 

1577 raise FloatingPointError(msg.format(filename)) 

1578 elif result == MinimizeResult.Failed: 

1579 raise RuntimeError("Chi2 minimization failure, cannot complete fit.") 

1580 else: 

1581 raise RuntimeError("Unxepected return code from minimize().") 

1582 else: 

1583 self.log.error("%s failed to converge after %d steps"%(name, max_steps)) 

1584 

1585 return chi2 

1586 

1587 def _make_output(self, ccdImageList, model, func): 

1588 """Return the internal jointcal models converted to the afw 

1589 structures that will be saved to disk. 

1590 

1591 Parameters 

1592 ---------- 

1593 ccdImageList : `lsst.jointcal.CcdImageList` 

1594 The list of CcdImages to get the output for. 

1595 model : `lsst.jointcal.AstrometryModel` or `lsst.jointcal.PhotometryModel` 

1596 The internal jointcal model to convert for each `lsst.jointcal.CcdImage`. 

1597 func : `str` 

1598 The name of the function to call on ``model`` to get the converted 

1599 structure. Must accept an `lsst.jointcal.CcdImage`. 

1600 

1601 Returns 

1602 ------- 

1603 output : `dict` [`tuple`, `lsst.jointcal.AstrometryModel`] or 

1604 `dict` [`tuple`, `lsst.jointcal.PhotometryModel`] 

1605 The data to be saved, keyed on (visit, detector). 

1606 """ 

1607 output = {} 

1608 for ccdImage in ccdImageList: 

1609 ccd = ccdImage.ccdId 

1610 visit = ccdImage.visit 

1611 self.log.debug("%s for visit: %d, ccd: %d", func, visit, ccd) 

1612 output[(visit, ccd)] = getattr(model, func)(ccdImage) 

1613 return output 

1614 

1615 def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef): 

1616 """ 

1617 Write the fitted astrometric results to a new 'jointcal_wcs' dataRef. 

1618 

1619 Parameters 

1620 ---------- 

1621 associations : `lsst.jointcal.Associations` 

1622 The star/reference star associations to fit. 

1623 model : `lsst.jointcal.AstrometryModel` 

1624 The astrometric model that was fit. 

1625 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef` 

1626 Dict of ccdImage identifiers to dataRefs that were fit. 

1627 """ 

1628 ccdImageList = associations.getCcdImageList() 

1629 output = self._make_output(ccdImageList, model, "makeSkyWcs") 

1630 for key, skyWcs in output.items(): 

1631 dataRef = visit_ccd_to_dataRef[key] 

1632 try: 

1633 dataRef.put(skyWcs, 'jointcal_wcs') 

1634 except pexExceptions.Exception as e: 

1635 self.log.fatal('Failed to write updated Wcs: %s', str(e)) 

1636 raise e 

1637 

1638 def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef): 

1639 """ 

1640 Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef. 

1641 

1642 Parameters 

1643 ---------- 

1644 associations : `lsst.jointcal.Associations` 

1645 The star/reference star associations to fit. 

1646 model : `lsst.jointcal.PhotometryModel` 

1647 The photoometric model that was fit. 

1648 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef` 

1649 Dict of ccdImage identifiers to dataRefs that were fit. 

1650 """ 

1651 

1652 ccdImageList = associations.getCcdImageList() 

1653 output = self._make_output(ccdImageList, model, "toPhotoCalib") 

1654 for key, photoCalib in output.items(): 

1655 dataRef = visit_ccd_to_dataRef[key] 

1656 try: 

1657 dataRef.put(photoCalib, 'jointcal_photoCalib') 

1658 except pexExceptions.Exception as e: 

1659 self.log.fatal('Failed to write updated PhotoCalib: %s', str(e)) 

1660 raise e