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 

39import lsst.meas.algorithms 

40from lsst.pipe.tasks.colorterms import ColortermLibrary 

41from lsst.verify import Job, Measurement 

42 

43from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask, ReferenceSourceSelectorTask 

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 

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['profile_jointcal'] = parsedCmd.profile_jointcal 

83 kwargs['butler'] = parsedCmd.butler 

84 

85 # organize data IDs by tract 

86 refListDict = {} 

87 for ref in parsedCmd.id.refList: 

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

89 # we call runDataRef() once with each tract 

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

91 return result 

92 

93 def __call__(self, args): 

94 """ 

95 Parameters 

96 ---------- 

97 args 

98 Arguments for Task.runDataRef() 

99 

100 Returns 

101 ------- 

102 pipe.base.Struct 

103 if self.doReturnResults is False: 

104 

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

106 

107 if self.doReturnResults is True: 

108 

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

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

111 """ 

112 exitStatus = 0 # exit status for shell 

113 

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

115 dataRefList, kwargs = args 

116 butler = kwargs.pop('butler') 

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

118 result = None 

119 try: 

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

121 exitStatus = result.exitStatus 

122 job_path = butler.get('verify_job_filename') 

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

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

125 if self.doRaise: 

126 raise e 

127 else: 

128 exitStatus = 1 

129 eName = type(e).__name__ 

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

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

132 

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

134 kwargs['butler'] = butler 

135 if self.doReturnResults: 

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

137 else: 

138 return pipeBase.Struct(exitStatus=exitStatus) 

139 

140 

141class JointcalConfig(pexConfig.Config): 

142 """Configuration for JointcalTask""" 

143 

144 doAstrometry = pexConfig.Field( 

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

146 dtype=bool, 

147 default=True 

148 ) 

149 doPhotometry = pexConfig.Field( 

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

151 dtype=bool, 

152 default=True 

153 ) 

154 coaddName = pexConfig.Field( 

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

156 dtype=str, 

157 default="deep" 

158 ) 

159 positionErrorPedestal = pexConfig.Field( 

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

161 dtype=float, 

162 default=0.02, 

163 ) 

164 photometryErrorPedestal = pexConfig.Field( 

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

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

167 dtype=float, 

168 default=0.0, 

169 ) 

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

171 matchCut = pexConfig.Field( 

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

173 dtype=float, 

174 default=3.0, 

175 ) 

176 minMeasurements = pexConfig.Field( 

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

178 dtype=int, 

179 default=2, 

180 ) 

181 minMeasuredStarsPerCcd = pexConfig.Field( 

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

183 dtype=int, 

184 default=100, 

185 ) 

186 minRefStarsPerCcd = pexConfig.Field( 

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

188 dtype=int, 

189 default=30, 

190 ) 

191 allowLineSearch = pexConfig.Field( 

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

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

194 dtype=bool, 

195 default=False 

196 ) 

197 astrometrySimpleOrder = pexConfig.Field( 

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

199 dtype=int, 

200 default=3, 

201 ) 

202 astrometryChipOrder = pexConfig.Field( 

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

204 dtype=int, 

205 default=1, 

206 ) 

207 astrometryVisitOrder = pexConfig.Field( 

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

209 dtype=int, 

210 default=5, 

211 ) 

212 useInputWcs = pexConfig.Field( 

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

214 dtype=bool, 

215 default=True, 

216 ) 

217 astrometryModel = pexConfig.ChoiceField( 

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

219 dtype=str, 

220 default="constrained", 

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

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

223 ) 

224 photometryModel = pexConfig.ChoiceField( 

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

226 dtype=str, 

227 default="constrainedMagnitude", 

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

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

230 " fitting in flux space.", 

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

232 " fitting in magnitude space.", 

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

234 " fitting in magnitude space.", 

235 } 

236 ) 

237 applyColorTerms = pexConfig.Field( 

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

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

240 dtype=bool, 

241 default=False 

242 ) 

243 colorterms = pexConfig.ConfigField( 

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

245 dtype=ColortermLibrary, 

246 ) 

247 photometryVisitOrder = pexConfig.Field( 

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

249 dtype=int, 

250 default=7, 

251 ) 

252 photometryDoRankUpdate = pexConfig.Field( 

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

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

255 dtype=bool, 

256 default=True, 

257 ) 

258 astrometryDoRankUpdate = pexConfig.Field( 

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

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

261 dtype=bool, 

262 default=True, 

263 ) 

264 outlierRejectSigma = pexConfig.Field( 

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

266 dtype=float, 

267 default=5.0, 

268 ) 

269 maxPhotometrySteps = pexConfig.Field( 

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

271 dtype=int, 

272 default=20, 

273 ) 

274 maxAstrometrySteps = pexConfig.Field( 

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

276 dtype=int, 

277 default=20, 

278 ) 

279 astrometryRefObjLoader = pexConfig.ConfigurableField( 

280 target=LoadIndexedReferenceObjectsTask, 

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

282 ) 

283 photometryRefObjLoader = pexConfig.ConfigurableField( 

284 target=LoadIndexedReferenceObjectsTask, 

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

286 ) 

287 sourceSelector = sourceSelectorRegistry.makeField( 

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

289 default="astrometry" 

290 ) 

291 astrometryReferenceSelector = pexConfig.ConfigurableField( 

292 target=ReferenceSourceSelectorTask, 

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

294 ) 

295 photometryReferenceSelector = pexConfig.ConfigurableField( 

296 target=ReferenceSourceSelectorTask, 

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

298 ) 

299 astrometryReferenceErr = pexConfig.Field( 

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

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

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

303 dtype=float, 

304 default=None, 

305 optional=True 

306 ) 

307 writeInitMatrix = pexConfig.Field( 

308 dtype=bool, 

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

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

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

312 default=False 

313 ) 

314 writeChi2FilesInitialFinal = pexConfig.Field( 

315 dtype=bool, 

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

317 default=False 

318 ) 

319 writeChi2FilesOuterLoop = pexConfig.Field( 

320 dtype=bool, 

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

322 default=False 

323 ) 

324 writeInitialModel = pexConfig.Field( 

325 dtype=bool, 

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

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

328 default=False 

329 ) 

330 debugOutputPath = pexConfig.Field( 

331 dtype=str, 

332 default=".", 

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

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

335 ) 

336 sourceFluxType = pexConfig.Field( 

337 dtype=str, 

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

339 default='Calib' 

340 ) 

341 

342 def validate(self): 

343 super().validate() 

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

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

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

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

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

349 "applyColorTerms=True will be ignored.") 

350 lsst.log.warn(msg) 

351 

352 def setDefaults(self): 

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

354 self.sourceSelector.name = 'science' 

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

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

357 # with dependable signal to noise ratio. 

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

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

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

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

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

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

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

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

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

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

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

369 # chosen from the usual QA flags for stars) 

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

371 badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated', 

372 'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag', 

373 'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter'] 

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

375 

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

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

378 self.astrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414" 

379 self.astrometryRefObjLoader.requireProperMotion = True 

380 self.astrometryRefObjLoader.anyFilterMapsToThis = 'phot_g_mean' 

381 self.photometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110" 

382 

383 

384def writeModel(model, filename, log): 

385 """Write model to outfile.""" 

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

387 file.write(repr(model)) 

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

389 

390 

391@dataclasses.dataclass 

392class JointcalInputData: 

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

394 visit: int 

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

396 catalog: lsst.afw.table.SourceCatalog 

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

398 visitInfo: lsst.afw.image.VisitInfo 

399 """The VisitInfo of this exposure.""" 

400 detector: lsst.afw.cameraGeom.Detector 

401 """The detector of this exposure.""" 

402 photoCalib: lsst.afw.image.PhotoCalib 

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

404 wcs: lsst.afw.geom.skyWcs 

405 """The WCS of this exposure.""" 

406 bbox: lsst.geom.Box2I 

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

408 filter: lsst.afw.image.FilterLabel 

409 """The filter of this exposure.""" 

410 

411 

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

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

414 same field. 

415 

416 Parameters 

417 ---------- 

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

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

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

421 Used to initialize the astrometry and photometry refObjLoaders. 

422 profile_jointcal : `bool` 

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

424 """ 

425 

426 ConfigClass = JointcalConfig 

427 RunnerClass = JointcalRunner 

428 _DefaultName = "jointcal" 

429 

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

431 super().__init__(**kwargs) 

432 self.profile_jointcal = profile_jointcal 

433 self.makeSubtask("sourceSelector") 

434 if self.config.doAstrometry: 

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

436 self.makeSubtask("astrometryReferenceSelector") 

437 else: 

438 self.astrometryRefObjLoader = None 

439 if self.config.doPhotometry: 

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

441 self.makeSubtask("photometryReferenceSelector") 

442 else: 

443 self.photometryRefObjLoader = None 

444 

445 # To hold various computed metrics for use by tests 

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

447 

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

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

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

451 def _getMetadataName(self): 

452 return None 

453 

454 @classmethod 

455 def _makeArgumentParser(cls): 

456 """Create an argument parser""" 

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

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

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

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

461 ContainerClass=PerTractCcdDataIdContainer) 

462 return parser 

463 

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

465 """ 

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

467 

468 Parameters 

469 ---------- 

470 data : `JointcalInputData` 

471 The loaded input data. 

472 associations : `lsst.jointcal.Associations` 

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

474 jointcalControl : `jointcal.JointcalControl` 

475 Control object for associations management 

476 

477 Returns 

478 ------ 

479 namedtuple 

480 ``wcs`` 

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

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

483 ``key`` 

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

485 (`namedtuple`). 

486 ``band`` 

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

488 """ 

489 goodSrc = self.sourceSelector.run(data.catalog) 

490 

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

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

493 else: 

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

495 data.visit, 

496 data.detector.getId()) 

497 associations.createCcdImage(goodSrc.sourceCat, 

498 data.wcs, 

499 data.visitInfo, 

500 data.bbox, 

501 data.filter.physicalLabel, 

502 data.photoCalib, 

503 data.detector, 

504 data.visit, 

505 data.detector.getId(), 

506 jointcalControl) 

507 

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

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

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

511 

512 def _readDataId(self, butler, dataId): 

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

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

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

516 visit = dataId["visit"] 

517 else: 

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

519 

520 catalog = butler.get('src', 

521 flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS, 

522 dataId=dataId) 

523 return JointcalInputData(visit=visit, 

524 catalog=catalog, 

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

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

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

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

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

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

531 

532 def loadData(self, dataRefs, associations, jointcalControl, profile_jointcal=False): 

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

534 visit_ccd_to_dataRef = {} 

535 oldWcsList = [] 

536 bands = [] 

537 load_cat_prof_file = 'jointcal_loadData.prof' if profile_jointcal else '' 

538 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

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

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

541 self.focalPlaneBBox = camera.getFpBBox() 

542 for dataRef in dataRefs: 

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

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

545 oldWcsList.append(result.wcs) 

546 visit_ccd_to_dataRef[result.key] = dataRef 

547 bands.append(result.band) 

548 bands = collections.Counter(bands) 

549 

550 return oldWcsList, bands, visit_ccd_to_dataRef 

551 

552 def _getDebugPath(self, filename): 

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

554 """ 

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

556 

557 def _prep_sky(self, associations, bands): 

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

559 been read. 

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 band associated with the catalog. See DM-9093 

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

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

572 

573 return boundingCircle, center, radius, defaultBand 

574 

575 @pipeBase.timeMethod 

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

577 """ 

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

579 

580 NOTE: this is for gen2 middleware only. 

581 

582 Parameters 

583 ---------- 

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

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

586 profile_jointcal : `bool` 

587 Profile the individual steps of jointcal. 

588 

589 Returns 

590 ------- 

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

592 Struct of metadata from the fit, containing: 

593 

594 ``dataRefs`` 

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

596 ``oldWcsList`` 

597 The original WCS from each dataRef 

598 ``metrics`` 

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

600 """ 

601 if len(dataRefs) == 0: 

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

603 

604 exitStatus = 0 # exit status for shell 

605 

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

607 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField) 

608 associations = lsst.jointcal.Associations() 

609 

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

611 associations, 

612 jointcalControl, 

613 profile_jointcal=profile_jointcal) 

614 

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

616 

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

618 

619 if self.config.doAstrometry: 

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

621 name="astrometry", 

622 refObjLoader=self.astrometryRefObjLoader, 

623 referenceSelector=self.astrometryReferenceSelector, 

624 fit_function=self._fit_astrometry, 

625 profile_jointcal=profile_jointcal, 

626 tract=tract) 

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

628 else: 

629 astrometry = Astrometry(None, None, None) 

630 

631 if self.config.doPhotometry: 

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

633 name="photometry", 

634 refObjLoader=self.photometryRefObjLoader, 

635 referenceSelector=self.photometryReferenceSelector, 

636 fit_function=self._fit_photometry, 

637 profile_jointcal=profile_jointcal, 

638 tract=tract, 

639 reject_bad_fluxes=True) 

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

641 else: 

642 photometry = Photometry(None, None) 

643 

644 return pipeBase.Struct(dataRefs=dataRefs, 

645 oldWcsList=oldWcsList, 

646 job=self.job, 

647 astrometryRefObjLoader=self.astrometryRefObjLoader, 

648 photometryRefObjLoader=self.photometryRefObjLoader, 

649 defaultBand=defaultBand, 

650 exitStatus=exitStatus) 

651 

652 def _get_refcat_coordinate_error_override(self, refCat, name): 

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

654 return the overridden error if necessary. 

655 

656 Parameters 

657 ---------- 

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

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

660 name : `str` 

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

662 

663 Returns 

664 ------- 

665 refCoordErr : `float` 

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

667 those fields. 

668 

669 Raises 

670 ------ 

671 lsst.pex.config.FieldValidationError 

672 Raised if the refcat does not contain coordinate errors and 

673 ``config.astrometryReferenceErr`` is not set. 

674 """ 

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

676 # keep old refcats from causing problems. 

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

678 if 'coord_raErr' not in refCat.schema: 

679 return 100 

680 else: 

681 return float('nan') 

682 

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

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

685 "and config.astrometryReferenceErr not supplied.") 

686 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr, 

687 self.config, 

688 msg) 

689 

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

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

692 self.config.astrometryReferenceErr) 

693 

694 if self.config.astrometryReferenceErr is None: 

695 return float('nan') 

696 else: 

697 return self.config.astrometryReferenceErr 

698 

699 def _compute_proper_motion_epoch(self, ccdImageList): 

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

701 

702 Parameters 

703 ---------- 

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

705 The images to compute the appropriate epoch for. 

706 

707 Returns 

708 ------- 

709 epoch : `astropy.time.Time` 

710 The date to use for proper motion corrections. 

711 """ 

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

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

714 

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

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

717 reject_bad_fluxes=False, *, 

718 name="", refObjLoader=None, referenceSelector=None, 

719 fit_function=None): 

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

721 

722 Parameters 

723 ---------- 

724 associations : `lsst.jointcal.Associations` 

725 The star/reference star associations to fit. 

726 defaultBand : `str` 

727 filter to load from reference catalog. 

728 center : `lsst.geom.SpherePoint` 

729 ICRS center of field to load from reference catalog. 

730 radius : `lsst.geom.Angle` 

731 On-sky radius to load from reference catalog. 

732 name : `str` 

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

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

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

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

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

738 fit_function : callable 

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

740 tract : `str`, optional 

741 Name of tract currently being fit. 

742 profile_jointcal : `bool`, optional 

743 Separately profile the fitting step. 

744 match_cut : `float`, optional 

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

746 associations.associateCatalogs. 

747 reject_bad_fluxes : `bool`, optional 

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

749 

750 Returns 

751 ------- 

752 result : `Photometry` or `Astrometry` 

753 Result of `fit_function()` 

754 """ 

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

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

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

758 associations.associateCatalogs(match_cut) 

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

760 associations.fittedStarListSize()) 

761 

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

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

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

765 center, radius, defaultBand, 

766 applyColorterms=applyColorterms, 

767 epoch=epoch) 

768 refCoordErr = self._get_refcat_coordinate_error_override(refCat, name) 

769 

770 associations.collectRefStars(refCat, 

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

772 fluxField, 

773 refCoordinateErr=refCoordErr, 

774 rejectBadFluxes=reject_bad_fluxes) 

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

776 associations.refStarListSize()) 

777 

778 associations.prepareFittedStars(self.config.minMeasurements) 

779 

780 self._check_star_lists(associations, name) 

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

782 associations.nFittedStarsWithAssociatedRefStar()) 

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

784 associations.fittedStarListSize()) 

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

786 associations.nCcdImagesValidForFit()) 

787 

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

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

790 with pipeBase.cmdLineTask.profile(load_cat_prof_file): 

791 result = fit_function(associations, dataName) 

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

793 # Save reference and measurement chi2 contributions for this data 

794 if self.config.writeChi2FilesInitialFinal: 

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

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

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

798 

799 return result 

800 

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

802 applyColorterms=False, epoch=None): 

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

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

805 

806 Parameters 

807 ---------- 

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

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

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

811 Source selector to apply to loaded reference catalog. 

812 center : `lsst.geom.SpherePoint` 

813 The center around which to load sources. 

814 radius : `lsst.geom.Angle` 

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

816 filterName : `str` 

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

818 applyColorterms : `bool` 

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

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

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

822 or `None` to not apply such corrections. 

823 

824 Returns 

825 ------- 

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

827 The loaded reference catalog. 

828 fluxField : `str` 

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

830 """ 

831 skyCircle = refObjLoader.loadSkyCircle(center, 

832 radius, 

833 filterName, 

834 epoch=epoch) 

835 

836 selected = referenceSelector.run(skyCircle.refCat) 

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

838 if not selected.sourceCat.isContiguous(): 

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

840 else: 

841 refCat = selected.sourceCat 

842 

843 if applyColorterms: 

844 refCatName = refObjLoader.ref_dataset_name 

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

846 filterName, refCatName) 

847 colorterm = self.config.colorterms.getColorterm( 

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

849 

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

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

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

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

854 

855 return refCat, skyCircle.fluxField 

856 

857 def _check_star_lists(self, associations, name): 

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

859 if associations.nCcdImagesValidForFit() == 0: 

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

861 if associations.fittedStarListSize() == 0: 

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

863 if associations.refStarListSize() == 0: 

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

865 

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

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

868 

869 Parameters 

870 ---------- 

871 associations : `lsst.jointcal.Associations` 

872 The star/reference star associations to fit. 

873 fit : `lsst.jointcal.FitterBase` 

874 The fitter to use for minimization. 

875 model : `lsst.jointcal.Model` 

876 The model being fit. 

877 chi2Label : `str` 

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

879 writeChi2Name : `str`, optional 

880 Filename prefix to write the chi2 contributions to. 

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

882 

883 Returns 

884 ------- 

885 chi2: `lsst.jointcal.Chi2Accumulator` 

886 The chi2 object for the current fitter and model. 

887 

888 Raises 

889 ------ 

890 FloatingPointError 

891 Raised if chi2 is infinite or NaN. 

892 ValueError 

893 Raised if the model is not valid. 

894 """ 

895 if writeChi2Name is not None: 

896 fullpath = self._getDebugPath(writeChi2Name) 

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

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

899 

900 chi2 = fit.computeChi2() 

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

902 self._check_stars(associations) 

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

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

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

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

907 return chi2 

908 

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

910 """ 

911 Fit the photometric data. 

912 

913 Parameters 

914 ---------- 

915 associations : `lsst.jointcal.Associations` 

916 The star/reference star associations to fit. 

917 dataName : `str` 

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

919 identifying debugging files. 

920 

921 Returns 

922 ------- 

923 fit_result : `namedtuple` 

924 fit : `lsst.jointcal.PhotometryFit` 

925 The photometric fitter used to perform the fit. 

926 model : `lsst.jointcal.PhotometryModel` 

927 The photometric model that was fit. 

928 """ 

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

930 

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

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

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

934 self.focalPlaneBBox, 

935 visitOrder=self.config.photometryVisitOrder, 

936 errorPedestal=self.config.photometryErrorPedestal) 

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

938 doLineSearch = self.config.allowLineSearch 

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

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

941 self.focalPlaneBBox, 

942 visitOrder=self.config.photometryVisitOrder, 

943 errorPedestal=self.config.photometryErrorPedestal) 

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

945 doLineSearch = self.config.allowLineSearch 

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

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

948 errorPedestal=self.config.photometryErrorPedestal) 

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

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

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

952 errorPedestal=self.config.photometryErrorPedestal) 

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

954 

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

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

957 # Save reference and measurement chi2 contributions for this data 

958 if self.config.writeChi2FilesInitialFinal: 

959 baseName = f"photometry_initial_chi2-{dataName}" 

960 else: 

961 baseName = None 

962 if self.config.writeInitialModel: 

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

964 writeModel(model, fullpath, self.log) 

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

966 

967 def getChi2Name(whatToFit): 

968 if self.config.writeChi2FilesOuterLoop: 

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

970 else: 

971 return None 

972 

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

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

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

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

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

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

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

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

981 writeChi2Name=getChi2Name("ModelVisit")) 

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

983 

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

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

986 writeChi2Name=getChi2Name("Model")) 

987 

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

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

990 writeChi2Name=getChi2Name("Fluxes")) 

991 

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

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

994 writeChi2Name=getChi2Name("ModelFluxes")) 

995 

996 model.freezeErrorTransform() 

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

998 

999 chi2 = self._iterate_fit(associations, 

1000 fit, 

1001 self.config.maxPhotometrySteps, 

1002 "photometry", 

1003 "Model Fluxes", 

1004 doRankUpdate=self.config.photometryDoRankUpdate, 

1005 doLineSearch=doLineSearch, 

1006 dataName=dataName) 

1007 

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

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

1010 return Photometry(fit, model) 

1011 

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

1013 """ 

1014 Fit the astrometric data. 

1015 

1016 Parameters 

1017 ---------- 

1018 associations : `lsst.jointcal.Associations` 

1019 The star/reference star associations to fit. 

1020 dataName : `str` 

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

1022 identifying debugging files. 

1023 

1024 Returns 

1025 ------- 

1026 fit_result : `namedtuple` 

1027 fit : `lsst.jointcal.AstrometryFit` 

1028 The astrometric fitter used to perform the fit. 

1029 model : `lsst.jointcal.AstrometryModel` 

1030 The astrometric model that was fit. 

1031 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler` 

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

1033 """ 

1034 

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

1036 

1037 associations.deprojectFittedStars() 

1038 

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

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

1041 # them so carefully? 

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

1043 

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

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

1046 sky_to_tan_projection, 

1047 chipOrder=self.config.astrometryChipOrder, 

1048 visitOrder=self.config.astrometryVisitOrder) 

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

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

1051 sky_to_tan_projection, 

1052 self.config.useInputWcs, 

1053 nNotFit=0, 

1054 order=self.config.astrometrySimpleOrder) 

1055 

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

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

1058 # Save reference and measurement chi2 contributions for this data 

1059 if self.config.writeChi2FilesInitialFinal: 

1060 baseName = f"astrometry_initial_chi2-{dataName}" 

1061 else: 

1062 baseName = None 

1063 if self.config.writeInitialModel: 

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

1065 writeModel(model, fullpath, self.log) 

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

1067 

1068 def getChi2Name(whatToFit): 

1069 if self.config.writeChi2FilesOuterLoop: 

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

1071 else: 

1072 return None 

1073 

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

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

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

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

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

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

1080 writeChi2Name=getChi2Name("DistortionsVisit")) 

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

1082 

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

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

1085 writeChi2Name=getChi2Name("Distortions")) 

1086 

1087 fit.minimize("Positions") 

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

1089 writeChi2Name=getChi2Name("Positions")) 

1090 

1091 fit.minimize("Distortions Positions") 

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

1093 writeChi2Name=getChi2Name("DistortionsPositions")) 

1094 

1095 chi2 = self._iterate_fit(associations, 

1096 fit, 

1097 self.config.maxAstrometrySteps, 

1098 "astrometry", 

1099 "Distortions Positions", 

1100 doRankUpdate=self.config.astrometryDoRankUpdate, 

1101 dataName=dataName) 

1102 

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

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

1105 

1106 return Astrometry(fit, model, sky_to_tan_projection) 

1107 

1108 def _check_stars(self, associations): 

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

1110 for ccdImage in associations.getCcdImageList(): 

1111 nMeasuredStars, nRefStars = ccdImage.countStars() 

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

1113 ccdImage.getName(), nMeasuredStars, nRefStars) 

1114 if nMeasuredStars < self.config.minMeasuredStarsPerCcd: 

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

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

1117 if nRefStars < self.config.minRefStarsPerCcd: 

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

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

1120 

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

1122 dataName="", 

1123 doRankUpdate=True, 

1124 doLineSearch=False): 

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

1126 

1127 Parameters 

1128 ---------- 

1129 associations : `lsst.jointcal.Associations` 

1130 The star/reference star associations to fit. 

1131 fitter : `lsst.jointcal.FitterBase` 

1132 The fitter to use for minimization. 

1133 max_steps : `int` 

1134 Maximum number of steps to run outlier rejection before declaring 

1135 convergence failure. 

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

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

1138 whatToFit : `str` 

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

1140 dataName : `str`, optional 

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

1142 for debugging. 

1143 doRankUpdate : `bool`, optional 

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

1145 matrix and gradient? 

1146 doLineSearch : `bool`, optional 

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

1148 

1149 Returns 

1150 ------- 

1151 chi2: `lsst.jointcal.Chi2Statistic` 

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

1153 

1154 Raises 

1155 ------ 

1156 FloatingPointError 

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

1158 RuntimeError 

1159 Raised if the fitter fails for some other reason; 

1160 log messages will provide further details. 

1161 """ 

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

1163 oldChi2 = lsst.jointcal.Chi2Statistic() 

1164 oldChi2.chi2 = float("inf") 

1165 for i in range(max_steps): 

1166 if self.config.writeChi2FilesOuterLoop: 

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

1168 else: 

1169 writeChi2Name = None 

1170 result = fitter.minimize(whatToFit, 

1171 self.config.outlierRejectSigma, 

1172 doRankUpdate=doRankUpdate, 

1173 doLineSearch=doLineSearch, 

1174 dumpMatrixFile=dumpMatrixFile) 

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

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

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

1178 

1179 if result == MinimizeResult.Converged: 

1180 if doRankUpdate: 

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

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

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

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

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

1186 

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

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

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

1190 

1191 break 

1192 elif result == MinimizeResult.Chi2Increased: 

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

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

1195 chi2Ratio = chi2.chi2 / oldChi2.chi2 

1196 if chi2Ratio > 1.5: 

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

1198 chi2.chi2, oldChi2.chi2, chi2Ratio) 

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

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

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

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

1203 # leaving a warning and bailing early. 

1204 if chi2Ratio > 10: 

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

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

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

1208 raise RuntimeError(msg) 

1209 oldChi2 = chi2 

1210 elif result == MinimizeResult.NonFinite: 

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

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

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

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

1215 raise FloatingPointError(msg.format(filename)) 

1216 elif result == MinimizeResult.Failed: 

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

1218 else: 

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

1220 else: 

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

1222 

1223 return chi2 

1224 

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

1226 """ 

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

1228 

1229 Parameters 

1230 ---------- 

1231 associations : `lsst.jointcal.Associations` 

1232 The star/reference star associations to fit. 

1233 model : `lsst.jointcal.AstrometryModel` 

1234 The astrometric model that was fit. 

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

1236 Dict of ccdImage identifiers to dataRefs that were fit. 

1237 """ 

1238 

1239 ccdImageList = associations.getCcdImageList() 

1240 for ccdImage in ccdImageList: 

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

1242 ccd = ccdImage.ccdId 

1243 visit = ccdImage.visit 

1244 dataRef = visit_ccd_to_dataRef[(visit, ccd)] 

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

1246 skyWcs = model.makeSkyWcs(ccdImage) 

1247 try: 

1248 dataRef.put(skyWcs, 'jointcal_wcs') 

1249 except pexExceptions.Exception as e: 

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

1251 raise e 

1252 

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

1254 """ 

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

1256 

1257 Parameters 

1258 ---------- 

1259 associations : `lsst.jointcal.Associations` 

1260 The star/reference star associations to fit. 

1261 model : `lsst.jointcal.PhotometryModel` 

1262 The photoometric model that was fit. 

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

1264 Dict of ccdImage identifiers to dataRefs that were fit. 

1265 """ 

1266 

1267 ccdImageList = associations.getCcdImageList() 

1268 for ccdImage in ccdImageList: 

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

1270 ccd = ccdImage.ccdId 

1271 visit = ccdImage.visit 

1272 dataRef = visit_ccd_to_dataRef[(visit, ccd)] 

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

1274 photoCalib = model.toPhotoCalib(ccdImage) 

1275 try: 

1276 dataRef.put(photoCalib, 'jointcal_photoCalib') 

1277 except pexExceptions.Exception as e: 

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

1279 raise e