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 collections 

23import os 

24 

25import numpy as np 

26import astropy.units as u 

27 

28import lsst.geom 

29import lsst.utils 

30import lsst.pex.config as pexConfig 

31import lsst.pipe.base as pipeBase 

32from lsst.afw.image import fluxErrFromABMagErr 

33import lsst.pex.exceptions as pexExceptions 

34import lsst.afw.table 

35import lsst.log 

36import lsst.meas.algorithms 

37from lsst.pipe.tasks.colorterms import ColortermLibrary 

38from lsst.verify import Job, Measurement 

39 

40from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask, ReferenceSourceSelectorTask 

41from lsst.meas.algorithms.sourceSelector import sourceSelectorRegistry 

42 

43from .dataIds import PerTractCcdDataIdContainer 

44 

45import lsst.jointcal 

46from lsst.jointcal import MinimizeResult 

47 

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

49 

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

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

52 

53 

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

55def add_measurement(job, name, value): 

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

57 job.measurements.insert(meas) 

58 

59 

60class JointcalRunner(pipeBase.ButlerInitializedTaskRunner): 

61 """Subclass of TaskRunner for jointcalTask 

62 

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

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

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

66 arguments generated by the ArgumentParser into the arguments expected by 

67 Jointcal.runDataRef(). 

68 

69 See pipeBase.TaskRunner for more information. 

70 """ 

71 

72 @staticmethod 

73 def getTargetList(parsedCmd, **kwargs): 

74 """ 

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

76 

77 Jointcal operates on lists of dataRefs simultaneously. 

78 """ 

79 kwargs['profile_jointcal'] = parsedCmd.profile_jointcal 

80 kwargs['butler'] = parsedCmd.butler 

81 

82 # organize data IDs by tract 

83 refListDict = {} 

84 for ref in parsedCmd.id.refList: 

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

86 # we call runDataRef() once with each tract 

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

88 return result 

89 

90 def __call__(self, args): 

91 """ 

92 Parameters 

93 ---------- 

94 args 

95 Arguments for Task.runDataRef() 

96 

97 Returns 

98 ------- 

99 pipe.base.Struct 

100 if self.doReturnResults is False: 

101 

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

103 

104 if self.doReturnResults is True: 

105 

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

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

108 """ 

109 exitStatus = 0 # exit status for shell 

110 

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

112 dataRefList, kwargs = args 

113 butler = kwargs.pop('butler') 

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

115 result = None 

116 try: 

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

118 exitStatus = result.exitStatus 

119 job_path = butler.get('verify_job_filename') 

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

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

122 if self.doRaise: 

123 raise e 

124 else: 

125 exitStatus = 1 

126 eName = type(e).__name__ 

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

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

129 

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

131 kwargs['butler'] = butler 

132 if self.doReturnResults: 

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

134 else: 

135 return pipeBase.Struct(exitStatus=exitStatus) 

136 

137 

138class JointcalConfig(pexConfig.Config): 

139 """Configuration for JointcalTask""" 

140 

141 doAstrometry = pexConfig.Field( 

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

143 dtype=bool, 

144 default=True 

145 ) 

146 doPhotometry = pexConfig.Field( 

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

148 dtype=bool, 

149 default=True 

150 ) 

151 coaddName = pexConfig.Field( 

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

153 dtype=str, 

154 default="deep" 

155 ) 

156 positionErrorPedestal = pexConfig.Field( 

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

158 dtype=float, 

159 default=0.02, 

160 ) 

161 photometryErrorPedestal = pexConfig.Field( 

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

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

164 dtype=float, 

165 default=0.0, 

166 ) 

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

168 matchCut = pexConfig.Field( 

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

170 dtype=float, 

171 default=3.0, 

172 ) 

173 minMeasurements = pexConfig.Field( 

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

175 dtype=int, 

176 default=2, 

177 ) 

178 minMeasuredStarsPerCcd = pexConfig.Field( 

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

180 dtype=int, 

181 default=100, 

182 ) 

183 minRefStarsPerCcd = pexConfig.Field( 

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

185 dtype=int, 

186 default=30, 

187 ) 

188 allowLineSearch = pexConfig.Field( 

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

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

191 dtype=bool, 

192 default=False 

193 ) 

194 astrometrySimpleOrder = pexConfig.Field( 

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

196 dtype=int, 

197 default=3, 

198 ) 

199 astrometryChipOrder = pexConfig.Field( 

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

201 dtype=int, 

202 default=1, 

203 ) 

204 astrometryVisitOrder = pexConfig.Field( 

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

206 dtype=int, 

207 default=5, 

208 ) 

209 useInputWcs = pexConfig.Field( 

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

211 dtype=bool, 

212 default=True, 

213 ) 

214 astrometryModel = pexConfig.ChoiceField( 

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

216 dtype=str, 

217 default="constrained", 

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

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

220 ) 

221 photometryModel = pexConfig.ChoiceField( 

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

223 dtype=str, 

224 default="constrainedMagnitude", 

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

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

227 " fitting in flux space.", 

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

229 " fitting in magnitude space.", 

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

231 " fitting in magnitude space.", 

232 } 

233 ) 

234 applyColorTerms = pexConfig.Field( 

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

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

237 dtype=bool, 

238 default=False 

239 ) 

240 colorterms = pexConfig.ConfigField( 

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

242 dtype=ColortermLibrary, 

243 ) 

244 photometryVisitOrder = pexConfig.Field( 

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

246 dtype=int, 

247 default=7, 

248 ) 

249 photometryDoRankUpdate = pexConfig.Field( 

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

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

252 dtype=bool, 

253 default=True, 

254 ) 

255 astrometryDoRankUpdate = pexConfig.Field( 

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

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

258 dtype=bool, 

259 default=True, 

260 ) 

261 outlierRejectSigma = pexConfig.Field( 

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

263 dtype=float, 

264 default=5.0, 

265 ) 

266 maxPhotometrySteps = pexConfig.Field( 

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

268 dtype=int, 

269 default=20, 

270 ) 

271 maxAstrometrySteps = pexConfig.Field( 

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

273 dtype=int, 

274 default=20, 

275 ) 

276 astrometryRefObjLoader = pexConfig.ConfigurableField( 

277 target=LoadIndexedReferenceObjectsTask, 

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

279 ) 

280 photometryRefObjLoader = pexConfig.ConfigurableField( 

281 target=LoadIndexedReferenceObjectsTask, 

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

283 ) 

284 sourceSelector = sourceSelectorRegistry.makeField( 

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

286 default="astrometry" 

287 ) 

288 astrometryReferenceSelector = pexConfig.ConfigurableField( 

289 target=ReferenceSourceSelectorTask, 

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

291 ) 

292 photometryReferenceSelector = pexConfig.ConfigurableField( 

293 target=ReferenceSourceSelectorTask, 

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

295 ) 

296 astrometryReferenceErr = pexConfig.Field( 

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

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

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

300 dtype=float, 

301 default=None, 

302 optional=True 

303 ) 

304 writeInitMatrix = pexConfig.Field( 

305 dtype=bool, 

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

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

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

309 default=False 

310 ) 

311 writeChi2FilesInitialFinal = pexConfig.Field( 

312 dtype=bool, 

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

314 default=False 

315 ) 

316 writeChi2FilesOuterLoop = pexConfig.Field( 

317 dtype=bool, 

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

319 default=False 

320 ) 

321 writeInitialModel = pexConfig.Field( 

322 dtype=bool, 

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

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

325 default=False 

326 ) 

327 debugOutputPath = pexConfig.Field( 

328 dtype=str, 

329 default=".", 

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

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

332 ) 

333 sourceFluxType = pexConfig.Field( 

334 dtype=str, 

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

336 default='Calib' 

337 ) 

338 

339 def validate(self): 

340 super().validate() 

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

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

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

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

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

346 "applyColorTerms=True will be ignored.") 

347 lsst.log.warn(msg) 

348 

349 def setDefaults(self): 

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

351 self.sourceSelector.name = 'science' 

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

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

354 # with dependable signal to noise ratio. 

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

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

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

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

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

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

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

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

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

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

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

366 # chosen from the usual QA flags for stars) 

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

368 badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated', 

369 'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag', 

370 'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter'] 

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

372 

373 # Default to Gaia-DR2 for astrometry and PS1-DR1 for photometry, 

374 # with a reasonable initial filterMap. 

375 self.astrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414" 

376 self.astrometryRefObjLoader.filterMap = {'u': 'phot_g_mean', 

377 'g': 'phot_g_mean', 

378 'r': 'phot_g_mean', 

379 'i': 'phot_g_mean', 

380 'z': 'phot_g_mean', 

381 'y': 'phot_g_mean'} 

382 self.photometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110" 

383 

384 

385def writeModel(model, filename, log): 

386 """Write model to outfile.""" 

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

388 file.write(repr(model)) 

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

390 

391 

392class JointcalTask(pipeBase.CmdLineTask): 

393 """Jointly astrometrically and photometrically calibrate a group of images.""" 

394 

395 ConfigClass = JointcalConfig 

396 RunnerClass = JointcalRunner 

397 _DefaultName = "jointcal" 

398 

399 def __init__(self, butler=None, profile_jointcal=False, **kwargs): 

400 """ 

401 Instantiate a JointcalTask. 

402 

403 Parameters 

404 ---------- 

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

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

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

408 Used to initialize the astrometry and photometry refObjLoaders. 

409 profile_jointcal : `bool` 

410 Set to True to profile different stages of this jointcal run. 

411 """ 

412 pipeBase.CmdLineTask.__init__(self, **kwargs) 

413 self.profile_jointcal = profile_jointcal 

414 self.makeSubtask("sourceSelector") 

415 if self.config.doAstrometry: 

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

417 self.makeSubtask("astrometryReferenceSelector") 

418 else: 

419 self.astrometryRefObjLoader = None 

420 if self.config.doPhotometry: 

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

422 self.makeSubtask("photometryReferenceSelector") 

423 else: 

424 self.photometryRefObjLoader = None 

425 

426 # To hold various computed metrics for use by tests 

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

428 

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

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

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

432 def _getMetadataName(self): 

433 return None 

434 

435 @classmethod 

436 def _makeArgumentParser(cls): 

437 """Create an argument parser""" 

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

439 parser.add_argument("--profile_jointcal", default=False, action="store_true", 

440 help="Profile steps of jointcal separately.") 

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

442 ContainerClass=PerTractCcdDataIdContainer) 

443 return parser 

444 

445 def _build_ccdImage(self, dataRef, associations, jointcalControl): 

446 """ 

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

448 

449 Parameters 

450 ---------- 

451 dataRef : `lsst.daf.persistence.ButlerDataRef` 

452 DataRef to extract info from. 

453 associations : `lsst.jointcal.Associations` 

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

455 jointcalControl : `jointcal.JointcalControl` 

456 Control object for associations management 

457 

458 Returns 

459 ------ 

460 namedtuple 

461 ``wcs`` 

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

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

464 ``key`` 

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

466 (`namedtuple`). 

467 ``filter`` 

468 This calexp's filter (`str`). 

469 """ 

470 if "visit" in dataRef.dataId.keys(): 

471 visit = dataRef.dataId["visit"] 

472 else: 

473 visit = dataRef.getButler().queryMetadata("calexp", ("visit"), dataRef.dataId)[0] 

474 

475 src = dataRef.get("src", flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS, immediate=True) 

476 

477 visitInfo = dataRef.get('calexp_visitInfo') 

478 detector = dataRef.get('calexp_detector') 

479 ccdId = detector.getId() 

480 photoCalib = dataRef.get('calexp_photoCalib') 

481 tanWcs = dataRef.get('calexp_wcs') 

482 bbox = dataRef.get('calexp_bbox') 

483 filt = dataRef.get('calexp_filter') 

484 filterName = filt.getName() 

485 

486 goodSrc = self.sourceSelector.run(src) 

487 

488 if len(goodSrc.sourceCat) == 0: 

489 self.log.warn("No sources selected in visit %s ccd %s", visit, ccdId) 

490 else: 

491 self.log.info("%d sources selected in visit %d ccd %d", len(goodSrc.sourceCat), visit, ccdId) 

492 associations.createCcdImage(goodSrc.sourceCat, 

493 tanWcs, 

494 visitInfo, 

495 bbox, 

496 filterName, 

497 photoCalib, 

498 detector, 

499 visit, 

500 ccdId, 

501 jointcalControl) 

502 

503 Result = collections.namedtuple('Result_from_build_CcdImage', ('wcs', 'key', 'filter')) 

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

505 return Result(tanWcs, Key(visit, ccdId), filterName) 

506 

507 def _getDebugPath(self, filename): 

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

509 """ 

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

511 

512 @pipeBase.timeMethod 

513 def runDataRef(self, dataRefs, profile_jointcal=False): 

514 """ 

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

516 

517 Parameters 

518 ---------- 

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

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

521 profile_jointcal : `bool` 

522 Profile the individual steps of jointcal. 

523 

524 Returns 

525 ------- 

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

527 Struct of metadata from the fit, containing: 

528 

529 ``dataRefs`` 

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

531 ``oldWcsList`` 

532 The original WCS from each dataRef 

533 ``metrics`` 

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

535 """ 

536 if len(dataRefs) == 0: 

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

538 

539 exitStatus = 0 # exit status for shell 

540 

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

542 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField) 

543 associations = lsst.jointcal.Associations() 

544 

545 visit_ccd_to_dataRef = {} 

546 oldWcsList = [] 

547 filters = [] 

548 load_cat_prof_file = 'jointcal_build_ccdImage.prof' if profile_jointcal else '' 

549 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

550 # We need the bounding-box of the focal plane for photometry visit models. 

551 # NOTE: we only need to read it once, because its the same for all exposures of a camera. 

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

553 self.focalPlaneBBox = camera.getFpBBox() 

554 for ref in dataRefs: 

555 result = self._build_ccdImage(ref, associations, jointcalControl) 

556 oldWcsList.append(result.wcs) 

557 visit_ccd_to_dataRef[result.key] = ref 

558 filters.append(result.filter) 

559 filters = collections.Counter(filters) 

560 

561 associations.computeCommonTangentPoint() 

562 

563 boundingCircle = associations.computeBoundingCircle() 

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

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

566 

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

568 

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

570 defaultFilter = filters.most_common(1)[0][0] 

571 self.log.debug("Using %s band for reference flux", defaultFilter) 

572 

573 # TODO: need a better way to get the tract. 

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

575 

576 if self.config.doAstrometry: 

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

578 name="astrometry", 

579 refObjLoader=self.astrometryRefObjLoader, 

580 referenceSelector=self.astrometryReferenceSelector, 

581 fit_function=self._fit_astrometry, 

582 profile_jointcal=profile_jointcal, 

583 tract=tract) 

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

585 else: 

586 astrometry = Astrometry(None, None, None) 

587 

588 if self.config.doPhotometry: 

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

590 name="photometry", 

591 refObjLoader=self.photometryRefObjLoader, 

592 referenceSelector=self.photometryReferenceSelector, 

593 fit_function=self._fit_photometry, 

594 profile_jointcal=profile_jointcal, 

595 tract=tract, 

596 filters=filters, 

597 reject_bad_fluxes=True) 

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

599 else: 

600 photometry = Photometry(None, None) 

601 

602 return pipeBase.Struct(dataRefs=dataRefs, 

603 oldWcsList=oldWcsList, 

604 job=self.job, 

605 astrometryRefObjLoader=self.astrometryRefObjLoader, 

606 photometryRefObjLoader=self.photometryRefObjLoader, 

607 defaultFilter=defaultFilter, 

608 exitStatus=exitStatus) 

609 

610 def _get_refcat_coordinate_error_override(self, refCat, name): 

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

612 return the overridden error if necessary. 

613 

614 Parameters 

615 ---------- 

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

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

618 name : `str` 

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

620 

621 Returns 

622 ------- 

623 refCoordErr : `float` 

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

625 those fields. 

626 

627 Raises 

628 ------ 

629 lsst.pex.config.FieldValidationError 

630 Raised if the refcat does not contain coordinate errors and 

631 ``config.astrometryReferenceErr`` is not set. 

632 """ 

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

634 # keep old refcats from causing problems. 

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

636 if 'coord_raErr' not in refCat.schema: 

637 return 100 

638 else: 

639 return float('nan') 

640 

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

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

643 "and config.astrometryReferenceErr not supplied.") 

644 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr, 

645 self.config, 

646 msg) 

647 

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

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

650 self.config.astrometryReferenceErr) 

651 

652 if self.config.astrometryReferenceErr is None: 

653 return float('nan') 

654 else: 

655 return self.config.astrometryReferenceErr 

656 

657 def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius, 

658 filters=[], 

659 tract="", profile_jointcal=False, match_cut=3.0, 

660 reject_bad_fluxes=False, *, 

661 name="", refObjLoader=None, referenceSelector=None, 

662 fit_function=None): 

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

664 

665 Parameters 

666 ---------- 

667 associations : `lsst.jointcal.Associations` 

668 The star/reference star associations to fit. 

669 defaultFilter : `str` 

670 filter to load from reference catalog. 

671 center : `lsst.geom.SpherePoint` 

672 ICRS center of field to load from reference catalog. 

673 radius : `lsst.geom.Angle` 

674 On-sky radius to load from reference catalog. 

675 name : `str` 

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

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

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

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

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

681 fit_function : callable 

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

683 filters : `list` [`str`], optional 

684 List of filters to load from the reference catalog. 

685 tract : `str`, optional 

686 Name of tract currently being fit. 

687 profile_jointcal : `bool`, optional 

688 Separately profile the fitting step. 

689 match_cut : `float`, optional 

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

691 associations.associateCatalogs. 

692 reject_bad_fluxes : `bool`, optional 

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

694 

695 Returns 

696 ------- 

697 result : `Photometry` or `Astrometry` 

698 Result of `fit_function()` 

699 """ 

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

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

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

703 associations.associateCatalogs(match_cut) 

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

705 associations.fittedStarListSize()) 

706 

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

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

709 center, radius, defaultFilter, 

710 applyColorterms=applyColorterms) 

711 refCoordErr = self._get_refcat_coordinate_error_override(refCat, name) 

712 

713 associations.collectRefStars(refCat, 

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

715 fluxField, 

716 refCoordinateErr=refCoordErr, 

717 rejectBadFluxes=reject_bad_fluxes) 

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

719 associations.refStarListSize()) 

720 

721 associations.prepareFittedStars(self.config.minMeasurements) 

722 

723 self._check_star_lists(associations, name) 

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

725 associations.nFittedStarsWithAssociatedRefStar()) 

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

727 associations.fittedStarListSize()) 

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

729 associations.nCcdImagesValidForFit()) 

730 

731 load_cat_prof_file = 'jointcal_fit_%s.prof'%name if profile_jointcal else '' 

732 dataName = "{}_{}".format(tract, defaultFilter) 

733 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

734 result = fit_function(associations, dataName) 

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

736 # Save reference and measurement chi2 contributions for this data 

737 if self.config.writeChi2FilesInitialFinal: 

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

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

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

741 

742 return result 

743 

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

745 applyColorterms=False): 

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

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

748 

749 Parameters 

750 ---------- 

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

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

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

754 Source selector to apply to loaded reference catalog. 

755 center : `lsst.geom.SpherePoint` 

756 The center around which to load sources. 

757 radius : `lsst.geom.Angle` 

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

759 filterName : `str` 

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

761 applyColorterms : `bool` 

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

763 

764 Returns 

765 ------- 

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

767 The loaded reference catalog. 

768 fluxField : `str` 

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

770 """ 

771 skyCircle = refObjLoader.loadSkyCircle(center, 

772 radius, 

773 filterName) 

774 

775 selected = referenceSelector.run(skyCircle.refCat) 

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

777 if not selected.sourceCat.isContiguous(): 

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

779 else: 

780 refCat = selected.sourceCat 

781 

782 if applyColorterms: 

783 try: 

784 refCatName = refObjLoader.ref_dataset_name 

785 except AttributeError: 

786 # NOTE: we need this try:except: block in place until we've completely removed a.net support. 

787 raise RuntimeError("Cannot perform colorterm corrections with a.net refcats.") 

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

789 filterName, refCatName) 

790 colorterm = self.config.colorterms.getColorterm( 

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

792 

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

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

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

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

797 

798 return refCat, skyCircle.fluxField 

799 

800 def _check_star_lists(self, associations, name): 

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

802 if associations.nCcdImagesValidForFit() == 0: 

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

804 if associations.fittedStarListSize() == 0: 

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

806 if associations.refStarListSize() == 0: 

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

808 

809 def _logChi2AndValidate(self, associations, fit, model, chi2Label="Model", 

810 writeChi2Name=None): 

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

812 

813 Parameters 

814 ---------- 

815 associations : `lsst.jointcal.Associations` 

816 The star/reference star associations to fit. 

817 fit : `lsst.jointcal.FitterBase` 

818 The fitter to use for minimization. 

819 model : `lsst.jointcal.Model` 

820 The model being fit. 

821 chi2Label : str, optional 

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

823 writeChi2Name : `str`, optional 

824 Filename prefix to write the chi2 contributions to. 

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

826 

827 Returns 

828 ------- 

829 chi2: `lsst.jointcal.Chi2Accumulator` 

830 The chi2 object for the current fitter and model. 

831 

832 Raises 

833 ------ 

834 FloatingPointError 

835 Raised if chi2 is infinite or NaN. 

836 ValueError 

837 Raised if the model is not valid. 

838 """ 

839 if writeChi2Name is not None: 

840 fullpath = self._getDebugPath(writeChi2Name) 

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

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

843 

844 chi2 = fit.computeChi2() 

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

846 self._check_stars(associations) 

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

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

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

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

851 return chi2 

852 

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

854 """ 

855 Fit the photometric data. 

856 

857 Parameters 

858 ---------- 

859 associations : `lsst.jointcal.Associations` 

860 The star/reference star associations to fit. 

861 dataName : `str` 

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

863 identifying debugging files. 

864 

865 Returns 

866 ------- 

867 fit_result : `namedtuple` 

868 fit : `lsst.jointcal.PhotometryFit` 

869 The photometric fitter used to perform the fit. 

870 model : `lsst.jointcal.PhotometryModel` 

871 The photometric model that was fit. 

872 """ 

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

874 

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

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

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

878 self.focalPlaneBBox, 

879 visitOrder=self.config.photometryVisitOrder, 

880 errorPedestal=self.config.photometryErrorPedestal) 

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

882 doLineSearch = self.config.allowLineSearch 

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

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

885 self.focalPlaneBBox, 

886 visitOrder=self.config.photometryVisitOrder, 

887 errorPedestal=self.config.photometryErrorPedestal) 

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

889 doLineSearch = self.config.allowLineSearch 

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

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

892 errorPedestal=self.config.photometryErrorPedestal) 

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

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

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

896 errorPedestal=self.config.photometryErrorPedestal) 

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

898 

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

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

901 # Save reference and measurement chi2 contributions for this data 

902 if self.config.writeChi2FilesInitialFinal: 

903 baseName = f"photometry_initial_chi2-{dataName}" 

904 else: 

905 baseName = None 

906 if self.config.writeInitialModel: 

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

908 writeModel(model, fullpath, self.log) 

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

910 

911 def getChi2Name(whatToFit): 

912 if self.config.writeChi2FilesOuterLoop: 

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

914 else: 

915 return None 

916 

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

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

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

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

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

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

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

924 self._logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name("ModelVisit")) 

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

926 

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

928 self._logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name("Model")) 

929 

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

931 self._logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name("Fluxes")) 

932 

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

934 self._logChi2AndValidate(associations, fit, model, "Fit prepared", 

935 writeChi2Name=getChi2Name("ModelFluxes")) 

936 

937 model.freezeErrorTransform() 

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

939 

940 chi2 = self._iterate_fit(associations, 

941 fit, 

942 self.config.maxPhotometrySteps, 

943 "photometry", 

944 "Model Fluxes", 

945 doRankUpdate=self.config.photometryDoRankUpdate, 

946 doLineSearch=doLineSearch, 

947 dataName=dataName) 

948 

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

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

951 return Photometry(fit, model) 

952 

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

954 """ 

955 Fit the astrometric data. 

956 

957 Parameters 

958 ---------- 

959 associations : `lsst.jointcal.Associations` 

960 The star/reference star associations to fit. 

961 dataName : `str` 

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

963 identifying debugging files. 

964 

965 Returns 

966 ------- 

967 fit_result : `namedtuple` 

968 fit : `lsst.jointcal.AstrometryFit` 

969 The astrometric fitter used to perform the fit. 

970 model : `lsst.jointcal.AstrometryModel` 

971 The astrometric model that was fit. 

972 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler` 

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

974 """ 

975 

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

977 

978 associations.deprojectFittedStars() 

979 

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

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

982 # them so carefully? 

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

984 

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

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

987 sky_to_tan_projection, 

988 chipOrder=self.config.astrometryChipOrder, 

989 visitOrder=self.config.astrometryVisitOrder) 

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

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

992 sky_to_tan_projection, 

993 self.config.useInputWcs, 

994 nNotFit=0, 

995 order=self.config.astrometrySimpleOrder) 

996 

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

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

999 # Save reference and measurement chi2 contributions for this data 

1000 if self.config.writeChi2FilesInitialFinal: 

1001 baseName = f"astrometry_initial_chi2-{dataName}" 

1002 else: 

1003 baseName = None 

1004 if self.config.writeInitialModel: 

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

1006 writeModel(model, fullpath, self.log) 

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

1008 

1009 def getChi2Name(whatToFit): 

1010 if self.config.writeChi2FilesOuterLoop: 

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

1012 else: 

1013 return None 

1014 

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

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

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

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

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

1020 self._logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name("DistortionsVisit")) 

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

1022 

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

1024 self._logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name("Distortions")) 

1025 

1026 fit.minimize("Positions") 

1027 self._logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name("Positions")) 

1028 

1029 fit.minimize("Distortions Positions") 

1030 self._logChi2AndValidate(associations, fit, model, "Fit prepared", 

1031 writeChi2Name=getChi2Name("DistortionsPositions")) 

1032 

1033 chi2 = self._iterate_fit(associations, 

1034 fit, 

1035 self.config.maxAstrometrySteps, 

1036 "astrometry", 

1037 "Distortions Positions", 

1038 doRankUpdate=self.config.astrometryDoRankUpdate, 

1039 dataName=dataName) 

1040 

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

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

1043 

1044 return Astrometry(fit, model, sky_to_tan_projection) 

1045 

1046 def _check_stars(self, associations): 

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

1048 for ccdImage in associations.getCcdImageList(): 

1049 nMeasuredStars, nRefStars = ccdImage.countStars() 

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

1051 ccdImage.getName(), nMeasuredStars, nRefStars) 

1052 if nMeasuredStars < self.config.minMeasuredStarsPerCcd: 

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

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

1055 if nRefStars < self.config.minRefStarsPerCcd: 

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

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

1058 

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

1060 dataName="", 

1061 doRankUpdate=True, 

1062 doLineSearch=False): 

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

1064 

1065 Parameters 

1066 ---------- 

1067 associations : `lsst.jointcal.Associations` 

1068 The star/reference star associations to fit. 

1069 fitter : `lsst.jointcal.FitterBase` 

1070 The fitter to use for minimization. 

1071 max_steps : `int` 

1072 Maximum number of steps to run outlier rejection before declaring 

1073 convergence failure. 

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

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

1076 whatToFit : `str` 

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

1078 dataName : `str`, optional 

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

1080 for debugging. 

1081 doRankUpdate : `bool`, optional 

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

1083 matrix and gradient? 

1084 doLineSearch : `bool`, optional 

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

1086 

1087 Returns 

1088 ------- 

1089 chi2: `lsst.jointcal.Chi2Statistic` 

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

1091 

1092 Raises 

1093 ------ 

1094 FloatingPointError 

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

1096 RuntimeError 

1097 Raised if the fitter fails for some other reason; 

1098 log messages will provide further details. 

1099 """ 

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

1101 for i in range(max_steps): 

1102 if self.config.writeChi2FilesOuterLoop: 

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

1104 else: 

1105 writeChi2Name = None 

1106 result = fitter.minimize(whatToFit, 

1107 self.config.outlierRejectSigma, 

1108 doRankUpdate=doRankUpdate, 

1109 doLineSearch=doLineSearch, 

1110 dumpMatrixFile=dumpMatrixFile) 

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

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

1113 writeChi2Name=writeChi2Name) 

1114 

1115 if result == MinimizeResult.Converged: 

1116 if doRankUpdate: 

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

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

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

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

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

1122 

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

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

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

1126 

1127 break 

1128 elif result == MinimizeResult.Chi2Increased: 

1129 self.log.warn("still some outliers but chi2 increases - retry") 

1130 elif result == MinimizeResult.NonFinite: 

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

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

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

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

1135 raise FloatingPointError(msg.format(filename)) 

1136 elif result == MinimizeResult.Failed: 

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

1138 else: 

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

1140 else: 

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

1142 

1143 return chi2 

1144 

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

1146 """ 

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

1148 

1149 Parameters 

1150 ---------- 

1151 associations : `lsst.jointcal.Associations` 

1152 The star/reference star associations to fit. 

1153 model : `lsst.jointcal.AstrometryModel` 

1154 The astrometric model that was fit. 

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

1156 Dict of ccdImage identifiers to dataRefs that were fit. 

1157 """ 

1158 

1159 ccdImageList = associations.getCcdImageList() 

1160 for ccdImage in ccdImageList: 

1161 # TODO: there must be a better way to identify this ccdImage than a visit,ccd pair? 

1162 ccd = ccdImage.ccdId 

1163 visit = ccdImage.visit 

1164 dataRef = visit_ccd_to_dataRef[(visit, ccd)] 

1165 self.log.info("Updating WCS for visit: %d, ccd: %d", visit, ccd) 

1166 skyWcs = model.makeSkyWcs(ccdImage) 

1167 try: 

1168 dataRef.put(skyWcs, 'jointcal_wcs') 

1169 except pexExceptions.Exception as e: 

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

1171 raise e 

1172 

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

1174 """ 

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

1176 

1177 Parameters 

1178 ---------- 

1179 associations : `lsst.jointcal.Associations` 

1180 The star/reference star associations to fit. 

1181 model : `lsst.jointcal.PhotometryModel` 

1182 The photoometric model that was fit. 

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

1184 Dict of ccdImage identifiers to dataRefs that were fit. 

1185 """ 

1186 

1187 ccdImageList = associations.getCcdImageList() 

1188 for ccdImage in ccdImageList: 

1189 # TODO: there must be a better way to identify this ccdImage than a visit,ccd pair? 

1190 ccd = ccdImage.ccdId 

1191 visit = ccdImage.visit 

1192 dataRef = visit_ccd_to_dataRef[(visit, ccd)] 

1193 self.log.info("Updating PhotoCalib for visit: %d, ccd: %d", visit, ccd) 

1194 photoCalib = model.toPhotoCalib(ccdImage) 

1195 try: 

1196 dataRef.put(photoCalib, 'jointcal_photoCalib') 

1197 except pexExceptions.Exception as e: 

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

1199 raise e