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", "physical_filter"), 

155 deferLoad=True, 

156 multiple=True, 

157 ) 

158 inputVisitSummary = pipeBase.connectionTypes.Input( 

159 doc="Visit summary tables built from the calexps for these observations.", 

160 name="visitSummary", 

161 storageClass="ExposureCatalog", 

162 dimensions=("instrument", "visit", "physical_filter"), 

163 deferLoad=True, 

164 multiple=True, 

165 ) 

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

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

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

169 astrometryRefCat = pipeBase.connectionTypes.Input( 

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

171 name="gaia_dr2_20200414", 

172 storageClass="SimpleCatalog", 

173 dimensions=("htm7",), 

174 deferLoad=True, 

175 multiple=True 

176 ) 

177 photometryRefCat = pipeBase.connectionTypes.Input( 

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

179 name="ps1_pv3_3pi_20170110", 

180 storageClass="SimpleCatalog", 

181 dimensions=("htm7",), 

182 deferLoad=True, 

183 multiple=True 

184 ) 

185 

186 outputWcs = pipeBase.connectionTypes.Output( 

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

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

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

190 name="jointcalSkyWcsCatalog", 

191 storageClass="ExposureCatalog", 

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

193 multiple=True 

194 ) 

195 outputPhotoCalib = pipeBase.connectionTypes.Output( 

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

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

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

199 name="jointcalPhotoCalibCatalog", 

200 storageClass="ExposureCatalog", 

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

202 multiple=True 

203 ) 

204 

205 

206class JointcalConfig(pipeBase.PipelineTaskConfig, 

207 pipelineConnections=JointcalTaskConnections): 

208 """Configuration for JointcalTask""" 

209 

210 doAstrometry = pexConfig.Field( 

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

212 dtype=bool, 

213 default=True 

214 ) 

215 doPhotometry = pexConfig.Field( 

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

217 dtype=bool, 

218 default=True 

219 ) 

220 coaddName = pexConfig.Field( 

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

222 dtype=str, 

223 default="deep" 

224 ) 

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

226 sourceFluxType = pexConfig.Field( 

227 dtype=str, 

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

229 default='Calib' 

230 ) 

231 positionErrorPedestal = pexConfig.Field( 

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

233 dtype=float, 

234 default=0.02, 

235 ) 

236 photometryErrorPedestal = pexConfig.Field( 

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

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

239 dtype=float, 

240 default=0.0, 

241 ) 

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

243 matchCut = pexConfig.Field( 

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

245 dtype=float, 

246 default=3.0, 

247 ) 

248 minMeasurements = pexConfig.Field( 

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

250 dtype=int, 

251 default=2, 

252 ) 

253 minMeasuredStarsPerCcd = pexConfig.Field( 

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

255 dtype=int, 

256 default=100, 

257 ) 

258 minRefStarsPerCcd = pexConfig.Field( 

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

260 dtype=int, 

261 default=30, 

262 ) 

263 allowLineSearch = pexConfig.Field( 

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

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

266 dtype=bool, 

267 default=False 

268 ) 

269 astrometrySimpleOrder = pexConfig.Field( 

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

271 dtype=int, 

272 default=3, 

273 ) 

274 astrometryChipOrder = pexConfig.Field( 

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

276 dtype=int, 

277 default=1, 

278 ) 

279 astrometryVisitOrder = pexConfig.Field( 

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

281 dtype=int, 

282 default=5, 

283 ) 

284 useInputWcs = pexConfig.Field( 

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

286 dtype=bool, 

287 default=True, 

288 ) 

289 astrometryModel = pexConfig.ChoiceField( 

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

291 dtype=str, 

292 default="constrained", 

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

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

295 ) 

296 photometryModel = pexConfig.ChoiceField( 

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

298 dtype=str, 

299 default="constrainedMagnitude", 

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

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

302 " fitting in flux space.", 

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

304 " fitting in magnitude space.", 

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

306 " fitting in magnitude space.", 

307 } 

308 ) 

309 applyColorTerms = pexConfig.Field( 

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

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

312 dtype=bool, 

313 default=False 

314 ) 

315 colorterms = pexConfig.ConfigField( 

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

317 dtype=ColortermLibrary, 

318 ) 

319 photometryVisitOrder = pexConfig.Field( 

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

321 dtype=int, 

322 default=7, 

323 ) 

324 photometryDoRankUpdate = pexConfig.Field( 

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

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

327 dtype=bool, 

328 default=True, 

329 ) 

330 astrometryDoRankUpdate = pexConfig.Field( 

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

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

333 dtype=bool, 

334 default=True, 

335 ) 

336 outlierRejectSigma = pexConfig.Field( 

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

338 dtype=float, 

339 default=5.0, 

340 ) 

341 maxPhotometrySteps = pexConfig.Field( 

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

343 dtype=int, 

344 default=20, 

345 ) 

346 maxAstrometrySteps = pexConfig.Field( 

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

348 dtype=int, 

349 default=20, 

350 ) 

351 astrometryRefObjLoader = pexConfig.ConfigurableField( 

352 target=LoadIndexedReferenceObjectsTask, 

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

354 ) 

355 photometryRefObjLoader = pexConfig.ConfigurableField( 

356 target=LoadIndexedReferenceObjectsTask, 

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

358 ) 

359 sourceSelector = sourceSelectorRegistry.makeField( 

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

361 default="astrometry" 

362 ) 

363 astrometryReferenceSelector = pexConfig.ConfigurableField( 

364 target=ReferenceSourceSelectorTask, 

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

366 ) 

367 photometryReferenceSelector = pexConfig.ConfigurableField( 

368 target=ReferenceSourceSelectorTask, 

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

370 ) 

371 astrometryReferenceErr = pexConfig.Field( 

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

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

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

375 dtype=float, 

376 default=None, 

377 optional=True 

378 ) 

379 

380 # configs for outputting debug information 

381 writeInitMatrix = pexConfig.Field( 

382 dtype=bool, 

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

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

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

386 default=False 

387 ) 

388 writeChi2FilesInitialFinal = pexConfig.Field( 

389 dtype=bool, 

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

391 default=False 

392 ) 

393 writeChi2FilesOuterLoop = pexConfig.Field( 

394 dtype=bool, 

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

396 default=False 

397 ) 

398 writeInitialModel = pexConfig.Field( 

399 dtype=bool, 

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

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

402 default=False 

403 ) 

404 debugOutputPath = pexConfig.Field( 

405 dtype=str, 

406 default=".", 

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

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

409 ) 

410 detailedProfile = pexConfig.Field( 

411 dtype=bool, 

412 default=False, 

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

414 ) 

415 

416 def validate(self): 

417 super().validate() 

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

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

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

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

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

423 "applyColorTerms=True will be ignored.") 

424 lsst.log.warn(msg) 

425 

426 def setDefaults(self): 

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

428 self.sourceSelector.name = 'science' 

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

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

431 # with dependable signal to noise ratio. 

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

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

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

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

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

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

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

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

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

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

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

443 # chosen from the usual QA flags for stars) 

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

445 badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated', 

446 'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag', 

447 'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter'] 

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

449 

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

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

452 self.astrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414" 

453 self.astrometryRefObjLoader.requireProperMotion = True 

454 self.astrometryRefObjLoader.anyFilterMapsToThis = 'phot_g_mean' 

455 self.photometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110" 

456 

457 

458def writeModel(model, filename, log): 

459 """Write model to outfile.""" 

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

461 file.write(repr(model)) 

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

463 

464 

465@dataclasses.dataclass 

466class JointcalInputData: 

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

468 visit: int 

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

470 catalog: lsst.afw.table.SourceCatalog 

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

472 visitInfo: lsst.afw.image.VisitInfo 

473 """The VisitInfo of this exposure.""" 

474 detector: lsst.afw.cameraGeom.Detector 

475 """The detector of this exposure.""" 

476 photoCalib: lsst.afw.image.PhotoCalib 

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

478 wcs: lsst.afw.geom.skyWcs 

479 """The WCS of this exposure.""" 

480 bbox: lsst.geom.Box2I 

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

482 filter: lsst.afw.image.FilterLabel 

483 """The filter of this exposure.""" 

484 

485 

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

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

488 same field. 

489 

490 Parameters 

491 ---------- 

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

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

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

495 Used to initialize the astrometry and photometry refObjLoaders. 

496 initInputs : `dict`, optional 

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

498 """ 

499 

500 ConfigClass = JointcalConfig 

501 RunnerClass = JointcalRunner 

502 _DefaultName = "jointcal" 

503 

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

505 super().__init__(**kwargs) 

506 self.makeSubtask("sourceSelector") 

507 if self.config.doAstrometry: 

508 if initInputs is None: 

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

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

511 self.makeSubtask("astrometryReferenceSelector") 

512 else: 

513 self.astrometryRefObjLoader = None 

514 if self.config.doPhotometry: 

515 if initInputs is None: 

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

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

518 self.makeSubtask("photometryReferenceSelector") 

519 else: 

520 self.photometryRefObjLoader = None 

521 

522 # To hold various computed metrics for use by tests 

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

524 

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

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

527 # outputs to the correct refs. 

528 inputs = butlerQC.get(inputRefs) 

529 # We want the tract number for writing debug files 

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

531 if self.config.doAstrometry: 

532 self.astrometryRefObjLoader = ReferenceObjectLoader( 

533 dataIds=[ref.datasetRef.dataId 

534 for ref in inputRefs.astrometryRefCat], 

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

536 config=self.config.astrometryRefObjLoader, 

537 log=self.log) 

538 if self.config.doPhotometry: 

539 self.photometryRefObjLoader = ReferenceObjectLoader( 

540 dataIds=[ref.datasetRef.dataId 

541 for ref in inputRefs.photometryRefCat], 

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

543 config=self.config.photometryRefObjLoader, 

544 log=self.log) 

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

546 if self.config.doAstrometry: 

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

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

549 if self.config.doPhotometry: 

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

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

552 

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

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

555 

556 Parameters 

557 ---------- 

558 butlerQC : `ButlerQuantumContext` 

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

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

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

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

563 The fitted objects to persist. 

564 outputRefs : `list` [`OutputQuantizedConnection`] 

565 The DatasetRefs to persist the data to. 

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

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

568 setter : `str` 

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

570 """ 

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

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

573 

574 def new_catalog(visit, size): 

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

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

577 catalog.resize(size) 

578 catalog['visit'] = visit 

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

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

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

582 return catalog 

583 

584 # count how many detectors have output for each visit 

585 detectors_per_visit = collections.defaultdict(int) 

586 for key in outputs: 

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

588 detectors_per_visit[key[0]] += 1 

589 

590 for ref in outputRefs: 

591 visit = ref.dataId['visit'] 

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

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

594 i = 0 

595 for detector in camera: 

596 detectorId = detector.getId() 

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

598 if key not in outputs: 

599 # skip detectors we don't have output for 

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

601 setter[3:], detectorId, visit) 

602 continue 

603 

604 catalog[i].setId(detectorId) 

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

606 i += 1 

607 

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

609 butlerQC.put(catalog, ref) 

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

611 

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

613 # Docstring inherited. 

614 

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

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

617 # so just use "flux" here. 

618 sourceFluxField = "flux" 

619 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField) 

620 associations = lsst.jointcal.Associations() 

621 self.focalPlaneBBox = inputCamera.getFpBBox() 

622 oldWcsList, bands = self._load_data(inputSourceTableVisit, 

623 inputVisitSummary, 

624 associations, 

625 jointcalControl) 

626 

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

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

629 

630 if self.config.doAstrometry: 

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

632 name="astrometry", 

633 refObjLoader=self.astrometryRefObjLoader, 

634 referenceSelector=self.astrometryReferenceSelector, 

635 fit_function=self._fit_astrometry, 

636 tract=tract, 

637 epoch=epoch) 

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

639 astrometry.model, 

640 "makeSkyWcs") 

641 else: 

642 astrometry_output = None 

643 

644 if self.config.doPhotometry: 

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

646 name="photometry", 

647 refObjLoader=self.photometryRefObjLoader, 

648 referenceSelector=self.photometryReferenceSelector, 

649 fit_function=self._fit_photometry, 

650 tract=tract, 

651 epoch=epoch, 

652 reject_bad_fluxes=True) 

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

654 photometry.model, 

655 "toPhotoCalib") 

656 else: 

657 photometry_output = None 

658 

659 return pipeBase.Struct(outputWcs=astrometry_output, 

660 outputPhotoCalib=photometry_output, 

661 job=self.job, 

662 astrometryRefObjLoader=self.astrometryRefObjLoader, 

663 photometryRefObjLoader=self.photometryRefObjLoader) 

664 

665 def _make_schema_table(self): 

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

667 SourceCatalog to insert values from the dataFrame into. 

668 

669 Returns 

670 ------- 

671 table : `lsst.afw.table.SourceTable` 

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

673 """ 

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

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

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

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

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

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

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

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

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

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

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

685 table.defineCentroid("centroid") 

686 table.defineShape("shape") 

687 return table 

688 

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

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

691 limited to just one detector. 

692 

693 Parameters 

694 ---------- 

695 table : `lsst.afw.table.SourceTable` 

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

697 populated with data from ``visitCatalog``. 

698 visitCatalog : `pandas.DataFrame` 

699 DataFrame to extract a detector catalog from. 

700 detectorId : `int` 

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

702 

703 Returns 

704 ------- 

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

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

707 """ 

708 # map from dataFrame column to afw table column 

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

710 'x': 'centroid_x', 

711 'y': 'centroid_y', 

712 'xErr': 'centroid_xErr', 

713 'yErr': 'centroid_yErr', 

714 'Ixx': 'shape_xx', 

715 'Iyy': 'shape_yy', 

716 'Ixy': 'shape_xy', 

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

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

719 } 

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

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

722 catalog.resize(sum(matched)) 

723 view = visitCatalog.loc[matched] 

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

725 catalog[afwCol] = view[dfCol] 

726 

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

728 len(catalog), 

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

730 detectorId) 

731 return catalog 

732 

733 def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations, jointcalControl): 

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

735 

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

737 

738 Parameters 

739 ---------- 

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

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

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

743 Visit-level exposure summary catalog with metadata. 

744 associations : `lsst.jointcal.Associations` 

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

746 jointcalControl : `jointcal.JointcalControl` 

747 Control object for C++ associations management. 

748 

749 Returns 

750 ------- 

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

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

753 bands : `list` [`str`] 

754 The filter bands of each input dataset. 

755 """ 

756 oldWcsList = [] 

757 bands = [] 

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

759 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

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

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

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

763 

764 for visitSummaryRef in inputVisitSummary: 

765 visitSummary = visitSummaryRef.get() 

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

767 selected = self.sourceSelector.run(visitCatalog) 

768 

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

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

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

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

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

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

775 if result is not None: 

776 oldWcsList.append(result.wcs) 

777 

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

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

780 bands = collections.Counter(bands) 

781 

782 return oldWcsList, bands 

783 

784 def _make_one_input_data(self, visitRecord, catalog): 

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

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

787 catalog=catalog, 

788 visitInfo=visitRecord.getVisitInfo(), 

789 detector=visitRecord.getDetector(), 

790 photoCalib=visitRecord.getPhotoCalib(), 

791 wcs=visitRecord.getWcs(), 

792 bbox=visitRecord.getBBox(), 

793 # ExposureRecord doesn't have a FilterLabel yet, 

794 # so we have to make one. 

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

796 physical=visitRecord['physical_filter'])) 

797 

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

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

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

801 def _getMetadataName(self): 

802 return None 

803 

804 @classmethod 

805 def _makeArgumentParser(cls): 

806 """Create an argument parser""" 

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

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

809 ContainerClass=PerTractCcdDataIdContainer) 

810 return parser 

811 

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

813 """ 

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

815 

816 Parameters 

817 ---------- 

818 data : `JointcalInputData` 

819 The loaded input data. 

820 associations : `lsst.jointcal.Associations` 

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

822 jointcalControl : `jointcal.JointcalControl` 

823 Control object for associations management 

824 

825 Returns 

826 ------ 

827 namedtuple or `None` 

828 ``wcs`` 

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

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

831 ``key`` 

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

833 (`namedtuple`). 

834 ``band`` 

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

836 `None` 

837 if there are no sources in the loaded catalog. 

838 """ 

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

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

841 return None 

842 

843 associations.createCcdImage(data.catalog, 

844 data.wcs, 

845 data.visitInfo, 

846 data.bbox, 

847 data.filter.physicalLabel, 

848 data.photoCalib, 

849 data.detector, 

850 data.visit, 

851 data.detector.getId(), 

852 jointcalControl) 

853 

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

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

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

857 

858 def _readDataId(self, butler, dataId): 

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

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

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

862 visit = dataId["visit"] 

863 else: 

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

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

866 

867 catalog = butler.get('src', 

868 flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS, 

869 dataId=dataId) 

870 goodSrc = self.sourceSelector.run(catalog) 

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

872 len(goodSrc.sourceCat), 

873 visit, 

874 detector.getId()) 

875 return JointcalInputData(visit=visit, 

876 catalog=goodSrc.sourceCat, 

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

878 detector=detector, 

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

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

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

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

883 

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

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

886 visit_ccd_to_dataRef = {} 

887 oldWcsList = [] 

888 bands = [] 

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

890 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

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

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

893 self.focalPlaneBBox = camera.getFpBBox() 

894 for dataRef in dataRefs: 

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

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

897 if result is None: 

898 continue 

899 oldWcsList.append(result.wcs) 

900 visit_ccd_to_dataRef[result.key] = dataRef 

901 bands.append(result.band) 

902 if len(bands) == 0: 

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

904 bands = collections.Counter(bands) 

905 

906 return oldWcsList, bands, visit_ccd_to_dataRef 

907 

908 def _getDebugPath(self, filename): 

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

910 """ 

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

912 

913 def _prep_sky(self, associations, bands): 

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

915 been read. 

916 """ 

917 associations.computeCommonTangentPoint() 

918 

919 boundingCircle = associations.computeBoundingCircle() 

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

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

922 

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

924 

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

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

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

928 

929 return boundingCircle, center, radius, defaultBand 

930 

931 @pipeBase.timeMethod 

932 def runDataRef(self, dataRefs): 

933 """ 

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

935 

936 NOTE: this is for gen2 middleware only. 

937 

938 Parameters 

939 ---------- 

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

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

942 

943 Returns 

944 ------- 

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

946 Struct of metadata from the fit, containing: 

947 

948 ``dataRefs`` 

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

950 ``oldWcsList`` 

951 The original WCS from each dataRef 

952 ``metrics`` 

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

954 """ 

955 if len(dataRefs) == 0: 

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

957 

958 exitStatus = 0 # exit status for shell 

959 

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

961 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField) 

962 associations = lsst.jointcal.Associations() 

963 

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

965 associations, 

966 jointcalControl) 

967 

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

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

970 

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

972 

973 if self.config.doAstrometry: 

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

975 name="astrometry", 

976 refObjLoader=self.astrometryRefObjLoader, 

977 referenceSelector=self.astrometryReferenceSelector, 

978 fit_function=self._fit_astrometry, 

979 tract=tract, 

980 epoch=epoch) 

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

982 else: 

983 astrometry = Astrometry(None, None, None) 

984 

985 if self.config.doPhotometry: 

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

987 name="photometry", 

988 refObjLoader=self.photometryRefObjLoader, 

989 referenceSelector=self.photometryReferenceSelector, 

990 fit_function=self._fit_photometry, 

991 tract=tract, 

992 epoch=epoch, 

993 reject_bad_fluxes=True) 

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

995 else: 

996 photometry = Photometry(None, None) 

997 

998 return pipeBase.Struct(dataRefs=dataRefs, 

999 oldWcsList=oldWcsList, 

1000 job=self.job, 

1001 astrometryRefObjLoader=self.astrometryRefObjLoader, 

1002 photometryRefObjLoader=self.photometryRefObjLoader, 

1003 defaultBand=defaultBand, 

1004 epoch=epoch, 

1005 exitStatus=exitStatus) 

1006 

1007 def _get_refcat_coordinate_error_override(self, refCat, name): 

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

1009 return the overridden error if necessary. 

1010 

1011 Parameters 

1012 ---------- 

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

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

1015 name : `str` 

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

1017 

1018 Returns 

1019 ------- 

1020 refCoordErr : `float` 

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

1022 those fields. 

1023 

1024 Raises 

1025 ------ 

1026 lsst.pex.config.FieldValidationError 

1027 Raised if the refcat does not contain coordinate errors and 

1028 ``config.astrometryReferenceErr`` is not set. 

1029 """ 

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

1031 # keep old refcats from causing problems. 

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

1033 if 'coord_raErr' not in refCat.schema: 

1034 return 100 

1035 else: 

1036 return float('nan') 

1037 

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

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

1040 "and config.astrometryReferenceErr not supplied.") 

1041 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr, 

1042 self.config, 

1043 msg) 

1044 

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

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

1047 self.config.astrometryReferenceErr) 

1048 

1049 if self.config.astrometryReferenceErr is None: 

1050 return float('nan') 

1051 else: 

1052 return self.config.astrometryReferenceErr 

1053 

1054 def _compute_proper_motion_epoch(self, ccdImageList): 

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

1056 

1057 Parameters 

1058 ---------- 

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

1060 The images to compute the appropriate epoch for. 

1061 

1062 Returns 

1063 ------- 

1064 epoch : `astropy.time.Time` 

1065 The date to use for proper motion corrections. 

1066 """ 

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

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

1069 

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

1071 tract="", match_cut=3.0, 

1072 reject_bad_fluxes=False, *, 

1073 name="", refObjLoader=None, referenceSelector=None, 

1074 fit_function=None, epoch=None): 

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

1076 

1077 Parameters 

1078 ---------- 

1079 associations : `lsst.jointcal.Associations` 

1080 The star/reference star associations to fit. 

1081 defaultBand : `str` 

1082 filter to load from reference catalog. 

1083 center : `lsst.geom.SpherePoint` 

1084 ICRS center of field to load from reference catalog. 

1085 radius : `lsst.geom.Angle` 

1086 On-sky radius to load from reference catalog. 

1087 name : `str` 

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

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

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

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

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

1093 fit_function : callable 

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

1095 tract : `str`, optional 

1096 Name of tract currently being fit. 

1097 match_cut : `float`, optional 

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

1099 associations.associateCatalogs. 

1100 reject_bad_fluxes : `bool`, optional 

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

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

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

1104 or `None` to not apply such corrections. 

1105 

1106 Returns 

1107 ------- 

1108 result : `Photometry` or `Astrometry` 

1109 Result of `fit_function()` 

1110 """ 

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

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

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

1114 associations.associateCatalogs(match_cut) 

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

1116 associations.fittedStarListSize()) 

1117 

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

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

1120 center, radius, defaultBand, 

1121 applyColorterms=applyColorterms, 

1122 epoch=epoch) 

1123 refCoordErr = self._get_refcat_coordinate_error_override(refCat, name) 

1124 

1125 associations.collectRefStars(refCat, 

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

1127 fluxField, 

1128 refCoordinateErr=refCoordErr, 

1129 rejectBadFluxes=reject_bad_fluxes) 

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

1131 associations.refStarListSize()) 

1132 

1133 associations.prepareFittedStars(self.config.minMeasurements) 

1134 

1135 self._check_star_lists(associations, name) 

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

1137 associations.nFittedStarsWithAssociatedRefStar()) 

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

1139 associations.fittedStarListSize()) 

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

1141 associations.nCcdImagesValidForFit()) 

1142 

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

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

1145 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

1146 result = fit_function(associations, dataName) 

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

1148 # Save reference and measurement chi2 contributions for this data 

1149 if self.config.writeChi2FilesInitialFinal: 

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

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

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

1153 

1154 return result 

1155 

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

1157 applyColorterms=False, epoch=None): 

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

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

1160 

1161 Parameters 

1162 ---------- 

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

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

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

1166 Source selector to apply to loaded reference catalog. 

1167 center : `lsst.geom.SpherePoint` 

1168 The center around which to load sources. 

1169 radius : `lsst.geom.Angle` 

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

1171 filterName : `str` 

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

1173 applyColorterms : `bool` 

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

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

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

1177 or `None` to not apply such corrections. 

1178 

1179 Returns 

1180 ------- 

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

1182 The loaded reference catalog. 

1183 fluxField : `str` 

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

1185 """ 

1186 skyCircle = refObjLoader.loadSkyCircle(center, 

1187 radius, 

1188 filterName, 

1189 epoch=epoch) 

1190 

1191 selected = referenceSelector.run(skyCircle.refCat) 

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

1193 if not selected.sourceCat.isContiguous(): 

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

1195 else: 

1196 refCat = selected.sourceCat 

1197 

1198 if applyColorterms: 

1199 refCatName = refObjLoader.ref_dataset_name 

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

1201 filterName, refCatName) 

1202 colorterm = self.config.colorterms.getColorterm( 

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

1204 

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

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

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

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

1209 

1210 return refCat, skyCircle.fluxField 

1211 

1212 def _check_star_lists(self, associations, name): 

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

1214 if associations.nCcdImagesValidForFit() == 0: 

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

1216 if associations.fittedStarListSize() == 0: 

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

1218 if associations.refStarListSize() == 0: 

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

1220 

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

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

1223 

1224 Parameters 

1225 ---------- 

1226 associations : `lsst.jointcal.Associations` 

1227 The star/reference star associations to fit. 

1228 fit : `lsst.jointcal.FitterBase` 

1229 The fitter to use for minimization. 

1230 model : `lsst.jointcal.Model` 

1231 The model being fit. 

1232 chi2Label : `str` 

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

1234 writeChi2Name : `str`, optional 

1235 Filename prefix to write the chi2 contributions to. 

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

1237 

1238 Returns 

1239 ------- 

1240 chi2: `lsst.jointcal.Chi2Accumulator` 

1241 The chi2 object for the current fitter and model. 

1242 

1243 Raises 

1244 ------ 

1245 FloatingPointError 

1246 Raised if chi2 is infinite or NaN. 

1247 ValueError 

1248 Raised if the model is not valid. 

1249 """ 

1250 if writeChi2Name is not None: 

1251 fullpath = self._getDebugPath(writeChi2Name) 

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

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

1254 

1255 chi2 = fit.computeChi2() 

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

1257 self._check_stars(associations) 

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

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

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

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

1262 return chi2 

1263 

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

1265 """ 

1266 Fit the photometric data. 

1267 

1268 Parameters 

1269 ---------- 

1270 associations : `lsst.jointcal.Associations` 

1271 The star/reference star associations to fit. 

1272 dataName : `str` 

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

1274 identifying debugging files. 

1275 

1276 Returns 

1277 ------- 

1278 fit_result : `namedtuple` 

1279 fit : `lsst.jointcal.PhotometryFit` 

1280 The photometric fitter used to perform the fit. 

1281 model : `lsst.jointcal.PhotometryModel` 

1282 The photometric model that was fit. 

1283 """ 

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

1285 

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

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

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

1289 self.focalPlaneBBox, 

1290 visitOrder=self.config.photometryVisitOrder, 

1291 errorPedestal=self.config.photometryErrorPedestal) 

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

1293 doLineSearch = self.config.allowLineSearch 

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

1295 model = lsst.jointcal.ConstrainedMagnitudeModel(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 == "simpleFlux": 

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

1303 errorPedestal=self.config.photometryErrorPedestal) 

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

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

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

1307 errorPedestal=self.config.photometryErrorPedestal) 

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

1309 

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

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

1312 # Save reference and measurement chi2 contributions for this data 

1313 if self.config.writeChi2FilesInitialFinal: 

1314 baseName = f"photometry_initial_chi2-{dataName}" 

1315 else: 

1316 baseName = None 

1317 if self.config.writeInitialModel: 

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

1319 writeModel(model, fullpath, self.log) 

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

1321 

1322 def getChi2Name(whatToFit): 

1323 if self.config.writeChi2FilesOuterLoop: 

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

1325 else: 

1326 return None 

1327 

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

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

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

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

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

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

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

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

1336 writeChi2Name=getChi2Name("ModelVisit")) 

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

1338 

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

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

1341 writeChi2Name=getChi2Name("Model")) 

1342 

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

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

1345 writeChi2Name=getChi2Name("Fluxes")) 

1346 

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

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

1349 writeChi2Name=getChi2Name("ModelFluxes")) 

1350 

1351 model.freezeErrorTransform() 

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

1353 

1354 chi2 = self._iterate_fit(associations, 

1355 fit, 

1356 self.config.maxPhotometrySteps, 

1357 "photometry", 

1358 "Model Fluxes", 

1359 doRankUpdate=self.config.photometryDoRankUpdate, 

1360 doLineSearch=doLineSearch, 

1361 dataName=dataName) 

1362 

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

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

1365 return Photometry(fit, model) 

1366 

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

1368 """ 

1369 Fit the astrometric data. 

1370 

1371 Parameters 

1372 ---------- 

1373 associations : `lsst.jointcal.Associations` 

1374 The star/reference star associations to fit. 

1375 dataName : `str` 

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

1377 identifying debugging files. 

1378 

1379 Returns 

1380 ------- 

1381 fit_result : `namedtuple` 

1382 fit : `lsst.jointcal.AstrometryFit` 

1383 The astrometric fitter used to perform the fit. 

1384 model : `lsst.jointcal.AstrometryModel` 

1385 The astrometric model that was fit. 

1386 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler` 

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

1388 """ 

1389 

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

1391 

1392 associations.deprojectFittedStars() 

1393 

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

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

1396 # them so carefully? 

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

1398 

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

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

1401 sky_to_tan_projection, 

1402 chipOrder=self.config.astrometryChipOrder, 

1403 visitOrder=self.config.astrometryVisitOrder) 

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

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

1406 sky_to_tan_projection, 

1407 self.config.useInputWcs, 

1408 nNotFit=0, 

1409 order=self.config.astrometrySimpleOrder) 

1410 

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

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

1413 # Save reference and measurement chi2 contributions for this data 

1414 if self.config.writeChi2FilesInitialFinal: 

1415 baseName = f"astrometry_initial_chi2-{dataName}" 

1416 else: 

1417 baseName = None 

1418 if self.config.writeInitialModel: 

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

1420 writeModel(model, fullpath, self.log) 

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

1422 

1423 def getChi2Name(whatToFit): 

1424 if self.config.writeChi2FilesOuterLoop: 

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

1426 else: 

1427 return None 

1428 

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

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

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

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

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

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

1435 writeChi2Name=getChi2Name("DistortionsVisit")) 

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

1437 

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

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

1440 writeChi2Name=getChi2Name("Distortions")) 

1441 

1442 fit.minimize("Positions") 

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

1444 writeChi2Name=getChi2Name("Positions")) 

1445 

1446 fit.minimize("Distortions Positions") 

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

1448 writeChi2Name=getChi2Name("DistortionsPositions")) 

1449 

1450 chi2 = self._iterate_fit(associations, 

1451 fit, 

1452 self.config.maxAstrometrySteps, 

1453 "astrometry", 

1454 "Distortions Positions", 

1455 doRankUpdate=self.config.astrometryDoRankUpdate, 

1456 dataName=dataName) 

1457 

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

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

1460 

1461 return Astrometry(fit, model, sky_to_tan_projection) 

1462 

1463 def _check_stars(self, associations): 

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

1465 for ccdImage in associations.getCcdImageList(): 

1466 nMeasuredStars, nRefStars = ccdImage.countStars() 

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

1468 ccdImage.getName(), nMeasuredStars, nRefStars) 

1469 if nMeasuredStars < self.config.minMeasuredStarsPerCcd: 

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

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

1472 if nRefStars < self.config.minRefStarsPerCcd: 

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

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

1475 

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

1477 dataName="", 

1478 doRankUpdate=True, 

1479 doLineSearch=False): 

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

1481 

1482 Parameters 

1483 ---------- 

1484 associations : `lsst.jointcal.Associations` 

1485 The star/reference star associations to fit. 

1486 fitter : `lsst.jointcal.FitterBase` 

1487 The fitter to use for minimization. 

1488 max_steps : `int` 

1489 Maximum number of steps to run outlier rejection before declaring 

1490 convergence failure. 

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

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

1493 whatToFit : `str` 

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

1495 dataName : `str`, optional 

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

1497 for debugging. 

1498 doRankUpdate : `bool`, optional 

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

1500 matrix and gradient? 

1501 doLineSearch : `bool`, optional 

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

1503 

1504 Returns 

1505 ------- 

1506 chi2: `lsst.jointcal.Chi2Statistic` 

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

1508 

1509 Raises 

1510 ------ 

1511 FloatingPointError 

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

1513 RuntimeError 

1514 Raised if the fitter fails for some other reason; 

1515 log messages will provide further details. 

1516 """ 

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

1518 oldChi2 = lsst.jointcal.Chi2Statistic() 

1519 oldChi2.chi2 = float("inf") 

1520 for i in range(max_steps): 

1521 if self.config.writeChi2FilesOuterLoop: 

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

1523 else: 

1524 writeChi2Name = None 

1525 result = fitter.minimize(whatToFit, 

1526 self.config.outlierRejectSigma, 

1527 doRankUpdate=doRankUpdate, 

1528 doLineSearch=doLineSearch, 

1529 dumpMatrixFile=dumpMatrixFile) 

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

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

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

1533 

1534 if result == MinimizeResult.Converged: 

1535 if doRankUpdate: 

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

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

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

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

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

1541 

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

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

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

1545 

1546 break 

1547 elif result == MinimizeResult.Chi2Increased: 

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

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

1550 chi2Ratio = chi2.chi2 / oldChi2.chi2 

1551 if chi2Ratio > 1.5: 

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

1553 chi2.chi2, oldChi2.chi2, chi2Ratio) 

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

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

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

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

1558 # leaving a warning and bailing early. 

1559 if chi2Ratio > 10: 

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

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

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

1563 raise RuntimeError(msg) 

1564 oldChi2 = chi2 

1565 elif result == MinimizeResult.NonFinite: 

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

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

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

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

1570 raise FloatingPointError(msg.format(filename)) 

1571 elif result == MinimizeResult.Failed: 

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

1573 else: 

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

1575 else: 

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

1577 

1578 return chi2 

1579 

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

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

1582 structures that will be saved to disk. 

1583 

1584 Parameters 

1585 ---------- 

1586 ccdImageList : `lsst.jointcal.CcdImageList` 

1587 The list of CcdImages to get the output for. 

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

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

1590 func : `str` 

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

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

1593 

1594 Returns 

1595 ------- 

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

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

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

1599 """ 

1600 output = {} 

1601 for ccdImage in ccdImageList: 

1602 ccd = ccdImage.ccdId 

1603 visit = ccdImage.visit 

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

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

1606 return output 

1607 

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

1609 """ 

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

1611 

1612 Parameters 

1613 ---------- 

1614 associations : `lsst.jointcal.Associations` 

1615 The star/reference star associations to fit. 

1616 model : `lsst.jointcal.AstrometryModel` 

1617 The astrometric model that was fit. 

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

1619 Dict of ccdImage identifiers to dataRefs that were fit. 

1620 """ 

1621 ccdImageList = associations.getCcdImageList() 

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

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

1624 dataRef = visit_ccd_to_dataRef[key] 

1625 try: 

1626 dataRef.put(skyWcs, 'jointcal_wcs') 

1627 except pexExceptions.Exception as e: 

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

1629 raise e 

1630 

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

1632 """ 

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

1634 

1635 Parameters 

1636 ---------- 

1637 associations : `lsst.jointcal.Associations` 

1638 The star/reference star associations to fit. 

1639 model : `lsst.jointcal.PhotometryModel` 

1640 The photoometric model that was fit. 

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

1642 Dict of ccdImage identifiers to dataRefs that were fit. 

1643 """ 

1644 

1645 ccdImageList = associations.getCcdImageList() 

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

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

1648 dataRef = visit_ccd_to_dataRef[key] 

1649 try: 

1650 dataRef.put(photoCalib, 'jointcal_photoCalib') 

1651 except pexExceptions.Exception as e: 

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

1653 raise e