lsst.jointcal 23.0.0+1690274f33
jointcal.py
Go to the documentation of this file.
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
37import lsst.afw.table
38import lsst.log
39from lsst.obs.base import Instrument
40from lsst.pipe.tasks.colorterms import ColortermLibrary
41from lsst.verify import Job, Measurement
42
43from lsst.meas.algorithms import (LoadIndexedReferenceObjectsTask, ReferenceSourceSelectorTask,
44 ReferenceObjectLoader)
45from lsst.meas.algorithms.sourceSelector import sourceSelectorRegistry
46
47from .dataIds import PerTractCcdDataIdContainer
48
49import lsst.jointcal
50from lsst.jointcal import MinimizeResult
51
52__all__ = ["JointcalConfig", "JointcalRunner", "JointcalTask"]
53
54Photometry = collections.namedtuple('Photometry', ('fit', 'model'))
55Astrometry = collections.namedtuple('Astrometry', ('fit', 'model', 'sky_to_tan_projection'))
56
57
58# TODO: move this to MeasurementSet in lsst.verify per DM-12655.
59def add_measurement(job, name, value):
60 meas = Measurement(job.metrics[name], value)
61 job.measurements.insert(meas)
62
63
64class JointcalRunner(pipeBase.ButlerInitializedTaskRunner):
65 """Subclass of TaskRunner for jointcalTask (gen2)
66
67 jointcalTask.runDataRef() takes a number of arguments, one of which is a list of dataRefs
68 extracted from the command line (whereas most CmdLineTasks' runDataRef methods take
69 single dataRef, are are called repeatedly). This class transforms the processed
70 arguments generated by the ArgumentParser into the arguments expected by
71 Jointcal.runDataRef().
72
73 See pipeBase.TaskRunner for more information.
74 """
75
76 @staticmethod
77 def getTargetList(parsedCmd, **kwargs):
78 """
79 Return a list of tuples per tract, each containing (dataRefs, kwargs).
80
81 Jointcal operates on lists of dataRefs simultaneously.
82 """
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
141def lookupStaticCalibrations(datasetType, registry, quantumDataId, collections):
142 """Lookup function that asserts/hopes that a static calibration dataset
143 exists in a particular collection, since this task can't provide a single
144 date/time to use to search for one properly.
145
146 This is mostly useful for the ``camera`` dataset, in cases where the task's
147 quantum dimensions do *not* include something temporal, like ``exposure``
148 or ``visit``.
149
150 Parameters
151 ----------
152 datasetType : `lsst.daf.butler.DatasetType`
153 Type of dataset being searched for.
154 registry : `lsst.daf.butler.Registry`
155 Data repository registry to search.
156 quantumDataId : `lsst.daf.butler.DataCoordinate`
157 Data ID of the quantum this camera should match.
158 collections : `Iterable` [ `str` ]
159 Collections that should be searched - but this lookup function works
160 by ignoring this in favor of a more-or-less hard-coded value.
161
162 Returns
163 -------
164 refs : `Iterator` [ `lsst.daf.butler.DatasetRef` ]
165 Iterator over dataset references; should have only one element.
166
167 Notes
168 -----
169 This implementation duplicates one in fgcmcal, and is at least quite
170 similar to another in cp_pipe. This duplicate has the most documentation.
171 Fixing this is DM-29661.
172 """
173 instrument = Instrument.fromName(quantumDataId["instrument"], registry)
174 unboundedCollection = instrument.makeUnboundedCalibrationRunName()
175 return registry.queryDatasets(datasetType,
176 dataId=quantumDataId,
177 collections=[unboundedCollection],
178 findFirst=True)
179
180
181def lookupVisitRefCats(datasetType, registry, quantumDataId, collections):
182 """Lookup function that finds all refcats for all visits that overlap a
183 tract, rather than just the refcats that directly overlap the tract.
184
185 Parameters
186 ----------
187 datasetType : `lsst.daf.butler.DatasetType`
188 Type of dataset being searched for.
189 registry : `lsst.daf.butler.Registry`
190 Data repository registry to search.
191 quantumDataId : `lsst.daf.butler.DataCoordinate`
192 Data ID of the quantum; expected to be something we can use as a
193 constraint to query for overlapping visits.
194 collections : `Iterable` [ `str` ]
195 Collections to search.
196
197 Returns
198 -------
199 refs : `Iterator` [ `lsst.daf.butler.DatasetRef` ]
200 Iterator over refcat references.
201 """
202 refs = set()
203 # Use .expanded() on the query methods below because we need data IDs with
204 # regions, both in the outer loop over visits (queryDatasets will expand
205 # any data ID we give it, but doing it up-front in bulk is much more
206 # efficient) and in the data IDs of the DatasetRefs this function yields
207 # (because the RefCatLoader relies on them to do some of its own
208 # filtering).
209 for visit_data_id in set(registry.queryDataIds("visit", dataId=quantumDataId).expanded()):
210 refs.update(
211 registry.queryDatasets(
212 datasetType,
213 collections=collections,
214 dataId=visit_data_id,
215 findFirst=True,
216 ).expanded()
217 )
218 yield from refs
219
220
221class JointcalTaskConnections(pipeBase.PipelineTaskConnections,
222 dimensions=("skymap", "tract", "instrument", "physical_filter")):
223 """Middleware input/output connections for jointcal data."""
224 inputCamera = pipeBase.connectionTypes.PrerequisiteInput(
225 doc="The camera instrument that took these observations.",
226 name="camera",
227 storageClass="Camera",
228 dimensions=("instrument",),
229 isCalibration=True,
230 lookupFunction=lookupStaticCalibrations,
231 )
232 inputSourceTableVisit = pipeBase.connectionTypes.Input(
233 doc="Source table in parquet format, per visit",
234 name="sourceTable_visit",
235 storageClass="DataFrame",
236 dimensions=("instrument", "visit"),
237 deferLoad=True,
238 multiple=True,
239 )
240 inputVisitSummary = pipeBase.connectionTypes.Input(
241 doc=("Per-visit consolidated exposure metadata built from calexps. "
242 "These catalogs use detector id for the id and must be sorted for "
243 "fast lookups of a detector."),
244 name="visitSummary",
245 storageClass="ExposureCatalog",
246 dimensions=("instrument", "visit"),
247 deferLoad=True,
248 multiple=True,
249 )
250 astrometryRefCat = pipeBase.connectionTypes.PrerequisiteInput(
251 doc="The astrometry reference catalog to match to loaded input catalog sources.",
252 name="gaia_dr2_20200414",
253 storageClass="SimpleCatalog",
254 dimensions=("skypix",),
255 deferLoad=True,
256 multiple=True,
257 lookupFunction=lookupVisitRefCats,
258 )
259 photometryRefCat = pipeBase.connectionTypes.PrerequisiteInput(
260 doc="The photometry reference catalog to match to loaded input catalog sources.",
261 name="ps1_pv3_3pi_20170110",
262 storageClass="SimpleCatalog",
263 dimensions=("skypix",),
264 deferLoad=True,
265 multiple=True,
266 lookupFunction=lookupVisitRefCats,
267 )
268
269 outputWcs = pipeBase.connectionTypes.Output(
270 doc=("Per-tract, per-visit world coordinate systems derived from the fitted model."
271 " These catalogs only contain entries for detectors with an output, and use"
272 " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
273 name="jointcalSkyWcsCatalog",
274 storageClass="ExposureCatalog",
275 dimensions=("instrument", "visit", "skymap", "tract"),
276 multiple=True
277 )
278 outputPhotoCalib = pipeBase.connectionTypes.Output(
279 doc=("Per-tract, per-visit photometric calibrations derived from the fitted model."
280 " These catalogs only contain entries for detectors with an output, and use"
281 " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
282 name="jointcalPhotoCalibCatalog",
283 storageClass="ExposureCatalog",
284 dimensions=("instrument", "visit", "skymap", "tract"),
285 multiple=True
286 )
287
288 def __init__(self, *, config=None):
289 super().__init__(config=config)
290 # When we are only doing one of astrometry or photometry, we don't
291 # need the reference catalog or produce the outputs for the other.
292 # This informs the middleware of that when the QuantumGraph is
293 # generated, so we don't block on getting something we won't need or
294 # create an expectation that downstream tasks will be able to consume
295 # something we won't produce.
296 if not config.doAstrometry:
297 self.prerequisiteInputs.remove("astrometryRefCat")
298 self.outputs.remove("outputWcs")
299 if not config.doPhotometry:
300 self.prerequisiteInputs.remove("photometryRefCat")
301 self.outputs.remove("outputPhotoCalib")
302
303
304class JointcalConfig(pipeBase.PipelineTaskConfig,
305 pipelineConnections=JointcalTaskConnections):
306 """Configuration for JointcalTask"""
307
308 doAstrometry = pexConfig.Field(
309 doc="Fit astrometry and write the fitted result.",
310 dtype=bool,
311 default=True
312 )
313 doPhotometry = pexConfig.Field(
314 doc="Fit photometry and write the fitted result.",
315 dtype=bool,
316 default=True
317 )
318 coaddName = pexConfig.Field(
319 doc="Type of coadd, typically deep or goodSeeing",
320 dtype=str,
321 default="deep"
322 )
323 # TODO DM-29008: Change this to "ApFlux_12_0" before gen2 removal.
324 sourceFluxType = pexConfig.Field(
325 dtype=str,
326 doc="Source flux field to use in source selection and to get fluxes from the catalog.",
327 default='Calib'
328 )
329 positionErrorPedestal = pexConfig.Field(
330 doc="Systematic term to apply to the measured position error (pixels)",
331 dtype=float,
332 default=0.02,
333 )
334 photometryErrorPedestal = pexConfig.Field(
335 doc="Systematic term to apply to the measured error on flux or magnitude as a "
336 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
337 dtype=float,
338 default=0.0,
339 )
340 # TODO: DM-6885 matchCut should be an geom.Angle
341 matchCut = pexConfig.Field(
342 doc="Matching radius between fitted and reference stars (arcseconds)",
343 dtype=float,
344 default=3.0,
345 )
346 minMeasurements = pexConfig.Field(
347 doc="Minimum number of associated measured stars for a fitted star to be included in the fit",
348 dtype=int,
349 default=2,
350 )
351 minMeasuredStarsPerCcd = pexConfig.Field(
352 doc="Minimum number of measuredStars per ccdImage before printing warnings",
353 dtype=int,
354 default=100,
355 )
356 minRefStarsPerCcd = pexConfig.Field(
357 doc="Minimum number of measuredStars per ccdImage before printing warnings",
358 dtype=int,
359 default=30,
360 )
361 allowLineSearch = pexConfig.Field(
362 doc="Allow a line search during minimization, if it is reasonable for the model"
363 " (models with a significant non-linear component, e.g. constrainedPhotometry).",
364 dtype=bool,
365 default=False
366 )
367 astrometrySimpleOrder = pexConfig.Field(
368 doc="Polynomial order for fitting the simple astrometry model.",
369 dtype=int,
370 default=3,
371 )
372 astrometryChipOrder = pexConfig.Field(
373 doc="Order of the per-chip transform for the constrained astrometry model.",
374 dtype=int,
375 default=1,
376 )
377 astrometryVisitOrder = pexConfig.Field(
378 doc="Order of the per-visit transform for the constrained astrometry model.",
379 dtype=int,
380 default=5,
381 )
382 useInputWcs = pexConfig.Field(
383 doc="Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
384 dtype=bool,
385 default=True,
386 )
387 astrometryModel = pexConfig.ChoiceField(
388 doc="Type of model to fit to astrometry",
389 dtype=str,
390 default="constrained",
391 allowed={"simple": "One polynomial per ccd",
392 "constrained": "One polynomial per ccd, and one polynomial per visit"}
393 )
394 photometryModel = pexConfig.ChoiceField(
395 doc="Type of model to fit to photometry",
396 dtype=str,
397 default="constrainedMagnitude",
398 allowed={"simpleFlux": "One constant zeropoint per ccd and visit, fitting in flux space.",
399 "constrainedFlux": "Constrained zeropoint per ccd, and one polynomial per visit,"
400 " fitting in flux space.",
401 "simpleMagnitude": "One constant zeropoint per ccd and visit,"
402 " fitting in magnitude space.",
403 "constrainedMagnitude": "Constrained zeropoint per ccd, and one polynomial per visit,"
404 " fitting in magnitude space.",
405 }
406 )
407 applyColorTerms = pexConfig.Field(
408 doc="Apply photometric color terms to reference stars?"
409 "Requires that colorterms be set to a ColortermLibrary",
410 dtype=bool,
411 default=False
412 )
413 colorterms = pexConfig.ConfigField(
414 doc="Library of photometric reference catalog name to color term dict.",
415 dtype=ColortermLibrary,
416 )
417 photometryVisitOrder = pexConfig.Field(
418 doc="Order of the per-visit polynomial transform for the constrained photometry model.",
419 dtype=int,
420 default=7,
421 )
422 photometryDoRankUpdate = pexConfig.Field(
423 doc=("Do the rank update step during minimization. "
424 "Skipping this can help deal with models that are too non-linear."),
425 dtype=bool,
426 default=True,
427 )
428 astrometryDoRankUpdate = pexConfig.Field(
429 doc=("Do the rank update step during minimization (should not change the astrometry fit). "
430 "Skipping this can help deal with models that are too non-linear."),
431 dtype=bool,
432 default=True,
433 )
434 outlierRejectSigma = pexConfig.Field(
435 doc="How many sigma to reject outliers at during minimization.",
436 dtype=float,
437 default=5.0,
438 )
439 astrometryOutlierRelativeTolerance = pexConfig.Field(
440 doc=("Convergence tolerance for outlier rejection threshold when fitting astrometry. Iterations will "
441 "stop when the fractional change in the chi2 cut level is below this value. If tolerance is set "
442 "to zero, iterations will continue until there are no more outliers. We suggest a value of 0.002"
443 "as a balance between a shorter minimization runtime and achieving a final fitted model that is"
444 "close to the solution found when removing all outliers."),
445 dtype=float,
446 default=0,
447 )
448 maxPhotometrySteps = pexConfig.Field(
449 doc="Maximum number of minimize iterations to take when fitting photometry.",
450 dtype=int,
451 default=20,
452 )
453 maxAstrometrySteps = pexConfig.Field(
454 doc="Maximum number of minimize iterations to take when fitting astrometry.",
455 dtype=int,
456 default=20,
457 )
458 astrometryRefObjLoader = pexConfig.ConfigurableField(
459 target=LoadIndexedReferenceObjectsTask,
460 doc="Reference object loader for astrometric fit",
461 )
462 photometryRefObjLoader = pexConfig.ConfigurableField(
463 target=LoadIndexedReferenceObjectsTask,
464 doc="Reference object loader for photometric fit",
465 )
466 sourceSelector = sourceSelectorRegistry.makeField(
467 doc="How to select sources for cross-matching",
468 default="astrometry"
469 )
470 astrometryReferenceSelector = pexConfig.ConfigurableField(
471 target=ReferenceSourceSelectorTask,
472 doc="How to down-select the loaded astrometry reference catalog.",
473 )
474 photometryReferenceSelector = pexConfig.ConfigurableField(
475 target=ReferenceSourceSelectorTask,
476 doc="How to down-select the loaded photometry reference catalog.",
477 )
478 astrometryReferenceErr = pexConfig.Field(
479 doc=("Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
480 "If None, then raise an exception if the reference catalog is missing coordinate errors. "
481 "If specified, overrides any existing `coord_*Err` values."),
482 dtype=float,
483 default=None,
484 optional=True
485 )
486
487 # configs for outputting debug information
488 writeInitMatrix = pexConfig.Field(
489 dtype=bool,
490 doc=("Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
491 "Output files will be written to `config.debugOutputPath` and will "
492 "be of the form 'astrometry_[pre|post]init-TRACT-FILTER-mat.txt'. "
493 "Note that these files are the dense versions of the matrix, and so may be very large."),
494 default=False
495 )
496 writeChi2FilesInitialFinal = pexConfig.Field(
497 dtype=bool,
498 doc=("Write .csv files containing the contributions to chi2 for the initialization and final fit. "
499 "Output files will be written to `config.debugOutputPath` and will "
500 "be of the form `astrometry_[initial|final]_chi2-TRACT-FILTER."),
501 default=False
502 )
503 writeChi2FilesOuterLoop = pexConfig.Field(
504 dtype=bool,
505 doc=("Write .csv files containing the contributions to chi2 for the outer fit loop. "
506 "Output files will be written to `config.debugOutputPath` and will "
507 "be of the form `astrometry_init-NN_chi2-TRACT-FILTER`."),
508 default=False
509 )
510 writeInitialModel = pexConfig.Field(
511 dtype=bool,
512 doc=("Write the pre-initialization model to text files, for debugging. "
513 "Output files will be written to `config.debugOutputPath` and will be "
514 "of the form `initial_astrometry_model-TRACT_FILTER.txt`."),
515 default=False
516 )
517 debugOutputPath = pexConfig.Field(
518 dtype=str,
519 default=".",
520 doc=("Path to write debug output files to. Used by "
521 "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
522 )
523 detailedProfile = pexConfig.Field(
524 dtype=bool,
525 default=False,
526 doc="Output separate profiling information for different parts of jointcal, e.g. data read, fitting"
527 )
528
529 def validate(self):
530 super().validate()
531 if self.doPhotometrydoPhotometry and self.applyColorTermsapplyColorTerms and len(self.colortermscolorterms.data) == 0:
532 msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
533 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
534 if self.doAstrometrydoAstrometry and not self.doPhotometrydoPhotometry and self.applyColorTermsapplyColorTerms:
535 msg = ("Only doing astrometry, but Colorterms are not applied for astrometry;"
536 "applyColorTerms=True will be ignored.")
537 lsst.log.warn(msg)
538
539 def setDefaults(self):
540 # Use science source selector which can filter on extendedness, SNR, and whether blended
541 self.sourceSelectorsourceSelector.name = 'science'
542 # Use only stars because aperture fluxes of galaxies are biased and depend on seeing
543 self.sourceSelectorsourceSelector['science'].doUnresolved = True
544 # with dependable signal to noise ratio.
545 self.sourceSelectorsourceSelector['science'].doSignalToNoise = True
546 # Min SNR must be > 0 because jointcal cannot handle negative fluxes,
547 # and S/N > 10 to use sources that are not too faint, and thus better measured.
548 self.sourceSelectorsourceSelector['science'].signalToNoise.minimum = 10.
549 # Base SNR on CalibFlux because that is the flux jointcal that fits and must be positive
550 fluxField = f"slot_{self.sourceFluxType}Flux_instFlux"
551 self.sourceSelectorsourceSelector['science'].signalToNoise.fluxField = fluxField
552 self.sourceSelectorsourceSelector['science'].signalToNoise.errField = fluxField + "Err"
553 # Do not trust blended sources' aperture fluxes which also depend on seeing
554 self.sourceSelectorsourceSelector['science'].doIsolated = True
555 # Do not trust either flux or centroid measurements with flags,
556 # chosen from the usual QA flags for stars)
557 self.sourceSelectorsourceSelector['science'].doFlags = True
558 badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated',
559 'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag',
560 'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter']
561 self.sourceSelectorsourceSelector['science'].flags.bad = badFlags
562
563 # Default to Gaia-DR2 (with proper motions) for astrometry and
564 # PS1-DR1 for photometry, with a reasonable initial filterMap.
565 self.astrometryRefObjLoaderastrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414"
566 self.astrometryRefObjLoaderastrometryRefObjLoader.requireProperMotion = True
567 self.astrometryRefObjLoaderastrometryRefObjLoader.anyFilterMapsToThis = 'phot_g_mean'
568 self.photometryRefObjLoaderphotometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110"
569
570
571def writeModel(model, filename, log):
572 """Write model to outfile."""
573 with open(filename, "w") as file:
574 file.write(repr(model))
575 log.info("Wrote %s to file: %s", model, filename)
576
577
578@dataclasses.dataclass
580 """The input data jointcal needs for each detector/visit."""
581 visit: int
582 """The visit identifier of this exposure."""
584 """The catalog derived from this exposure."""
585 visitInfo: lsst.afw.image.VisitInfo
586 """The VisitInfo of this exposure."""
588 """The detector of this exposure."""
589 photoCalib: lsst.afw.image.PhotoCalib
590 """The photometric calibration of this exposure."""
592 """The WCS of this exposure."""
593 bbox: lsst.geom.Box2I
594 """The bounding box of this exposure."""
596 """The filter of this exposure."""
597
598
599class JointcalTask(pipeBase.PipelineTask, pipeBase.CmdLineTask):
600 """Astrometricly and photometricly calibrate across multiple visits of the
601 same field.
602
603 Parameters
604 ----------
605 butler : `lsst.daf.persistence.Butler`
606 The butler is passed to the refObjLoader constructor in case it is
607 needed. Ignored if the refObjLoader argument provides a loader directly.
608 Used to initialize the astrometry and photometry refObjLoaders.
609 initInputs : `dict`, optional
610 Dictionary used to initialize PipelineTasks (empty for jointcal).
611 """
612
613 ConfigClass = JointcalConfig
614 RunnerClass = JointcalRunner
615 _DefaultName = "jointcal"
616
617 def __init__(self, butler=None, initInputs=None, **kwargs):
618 super().__init__(**kwargs)
619 self.makeSubtask("sourceSelector")
620 if self.config.doAstrometry:
621 if initInputs is None:
622 # gen3 middleware does refcat things internally (and will not have a butler here)
623 self.makeSubtask('astrometryRefObjLoader', butler=butler)
624 self.makeSubtask("astrometryReferenceSelector")
625 else:
626 self.astrometryRefObjLoaderastrometryRefObjLoader = None
627 if self.config.doPhotometry:
628 if initInputs is None:
629 # gen3 middleware does refcat things internally (and will not have a butler here)
630 self.makeSubtask('photometryRefObjLoader', butler=butler)
631 self.makeSubtask("photometryReferenceSelector")
632 else:
633 self.photometryRefObjLoaderphotometryRefObjLoader = None
634
635 # To hold various computed metrics for use by tests
636 self.jobjob = Job.load_metrics_package(subset='jointcal')
637
638 def runQuantum(self, butlerQC, inputRefs, outputRefs):
639 # We override runQuantum to set up the refObjLoaders and write the
640 # outputs to the correct refs.
641 inputs = butlerQC.get(inputRefs)
642 # We want the tract number for writing debug files
643 tract = butlerQC.quantum.dataId['tract']
644 if self.config.doAstrometry:
645 self.astrometryRefObjLoaderastrometryRefObjLoader = ReferenceObjectLoader(
646 dataIds=[ref.datasetRef.dataId
647 for ref in inputRefs.astrometryRefCat],
648 refCats=inputs.pop('astrometryRefCat'),
649 config=self.config.astrometryRefObjLoader,
650 log=self.log)
651 if self.config.doPhotometry:
652 self.photometryRefObjLoaderphotometryRefObjLoader = ReferenceObjectLoader(
653 dataIds=[ref.datasetRef.dataId
654 for ref in inputRefs.photometryRefCat],
655 refCats=inputs.pop('photometryRefCat'),
656 config=self.config.photometryRefObjLoader,
657 log=self.log)
658 outputs = self.runrun(**inputs, tract=tract)
659 if self.config.doAstrometry:
660 self._put_output_put_output(butlerQC, outputs.outputWcs, outputRefs.outputWcs,
661 inputs['inputCamera'], "setWcs")
662 if self.config.doPhotometry:
663 self._put_output_put_output(butlerQC, outputs.outputPhotoCalib, outputRefs.outputPhotoCalib,
664 inputs['inputCamera'], "setPhotoCalib")
665
666 def _put_output(self, butlerQC, outputs, outputRefs, camera, setter):
667 """Persist the output datasets to their appropriate datarefs.
668
669 Parameters
670 ----------
671 butlerQC : `ButlerQuantumContext`
672 A butler which is specialized to operate in the context of a
673 `lsst.daf.butler.Quantum`; This is the input to `runQuantum`.
674 outputs : `dict` [`tuple`, `lsst.afw.geom.SkyWcs`] or
675 `dict` [`tuple, `lsst.afw.image.PhotoCalib`]
676 The fitted objects to persist.
677 outputRefs : `list` [`OutputQuantizedConnection`]
678 The DatasetRefs to persist the data to.
680 The camera for this instrument, to get detector ids from.
681 setter : `str`
682 The method to call on the ExposureCatalog to set each output.
683 """
685 schema.addField('visit', type='I', doc='Visit number')
686
687 def new_catalog(visit, size):
688 """Return an catalog ready to be filled with appropriate output."""
689 catalog = lsst.afw.table.ExposureCatalog(schema)
690 catalog.resize(size)
691 catalog['visit'] = visit
692 metadata = lsst.daf.base.PropertyList()
693 metadata.add("COMMENT", "Catalog id is detector id, sorted.")
694 metadata.add("COMMENT", "Only detectors with data have entries.")
695 return catalog
696
697 # count how many detectors have output for each visit
698 detectors_per_visit = collections.defaultdict(int)
699 for key in outputs:
700 # key is (visit, detector_id), and we only need visit here
701 detectors_per_visit[key[0]] += 1
702
703 for ref in outputRefs:
704 visit = ref.dataId['visit']
705 catalog = new_catalog(visit, detectors_per_visit[visit])
706 # Iterate over every detector and skip the ones we don't have output for.
707 i = 0
708 for detector in camera:
709 detectorId = detector.getId()
710 key = (ref.dataId['visit'], detectorId)
711 if key not in outputs:
712 # skip detectors we don't have output for
713 self.log.debug("No %s output for detector %s in visit %s",
714 setter[3:], detectorId, visit)
715 continue
716
717 catalog[i].setId(detectorId)
718 getattr(catalog[i], setter)(outputs[key])
719 i += 1
720
721 catalog.sort() # ensure that the detectors are in sorted order, for fast lookups
722 butlerQC.put(catalog, ref)
723 self.log.info("Wrote %s detectors to %s", i, ref)
724
725 def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None):
726 # Docstring inherited.
727
728 # We take values out of the Parquet table, and put them in "flux_",
729 # and the config.sourceFluxType field is used during that extraction,
730 # so just use "flux" here.
731 sourceFluxField = "flux"
732 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
733 associations = lsst.jointcal.Associations()
734 self.focalPlaneBBoxfocalPlaneBBox = inputCamera.getFpBBox()
735 oldWcsList, bands = self._load_data_load_data(inputSourceTableVisit,
736 inputVisitSummary,
737 associations,
738 jointcalControl,
739 inputCamera)
740
741 boundingCircle, center, radius, defaultFilter, epoch = self._prep_sky_prep_sky(associations, bands)
742
743 if self.config.doAstrometry:
744 astrometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
745 name="astrometry",
746 refObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
747 referenceSelector=self.astrometryReferenceSelector,
748 fit_function=self._fit_astrometry_fit_astrometry,
749 tract=tract,
750 epoch=epoch)
751 astrometry_output = self._make_output_make_output(associations.getCcdImageList(),
752 astrometry.model,
753 "makeSkyWcs")
754 else:
755 astrometry_output = None
756
757 if self.config.doPhotometry:
758 photometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
759 name="photometry",
760 refObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
761 referenceSelector=self.photometryReferenceSelector,
762 fit_function=self._fit_photometry_fit_photometry,
763 tract=tract,
764 epoch=epoch,
765 reject_bad_fluxes=True)
766 photometry_output = self._make_output_make_output(associations.getCcdImageList(),
767 photometry.model,
768 "toPhotoCalib")
769 else:
770 photometry_output = None
771
772 return pipeBase.Struct(outputWcs=astrometry_output,
773 outputPhotoCalib=photometry_output,
774 job=self.jobjob,
775 astrometryRefObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
776 photometryRefObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader)
777
778 def _make_schema_table(self):
779 """Return an afw SourceTable to use as a base for creating the
780 SourceCatalog to insert values from the dataFrame into.
781
782 Returns
783 -------
785 Table with schema and slots to use to make SourceCatalogs.
786 """
788 schema.addField("centroid_x", "D")
789 schema.addField("centroid_y", "D")
790 schema.addField("centroid_xErr", "F")
791 schema.addField("centroid_yErr", "F")
792 schema.addField("shape_xx", "D")
793 schema.addField("shape_yy", "D")
794 schema.addField("shape_xy", "D")
795 schema.addField("flux_instFlux", "D")
796 schema.addField("flux_instFluxErr", "D")
797 table = lsst.afw.table.SourceTable.make(schema)
798 table.defineCentroid("centroid")
799 table.defineShape("shape")
800 return table
801
802 def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId,
803 detectorColumn, ixxColumns):
804 """Return an afw SourceCatalog extracted from a visit-level dataframe,
805 limited to just one detector.
806
807 Parameters
808 ----------
810 Table factory to use to make the SourceCatalog that will be
811 populated with data from ``visitCatalog``.
812 visitCatalog : `pandas.DataFrame`
813 DataFrame to extract a detector catalog from.
814 detectorId : `int`
815 Numeric id of the detector to extract from ``visitCatalog``.
816 detectorColumn : `str`
817 Name of the detector column in the catalog.
818 ixxColumns : `list` [`str`]
819 Names of the ixx/iyy/ixy columns in the catalog.
820
821 Returns
822 -------
824 Detector-level catalog extracted from ``visitCatalog``.
825 """
826 # map from dataFrame column to afw table column
827 mapping = {'x': 'centroid_x',
828 'y': 'centroid_y',
829 'xErr': 'centroid_xErr',
830 'yErr': 'centroid_yErr',
831 ixxColumns[0]: 'shape_xx',
832 ixxColumns[1]: 'shape_yy',
833 ixxColumns[2]: 'shape_xy',
834 f'{self.config.sourceFluxType}_instFlux': 'flux_instFlux',
835 f'{self.config.sourceFluxType}_instFluxErr': 'flux_instFluxErr',
836 }
837
838 catalog = lsst.afw.table.SourceCatalog(table)
839 matched = visitCatalog[detectorColumn] == detectorId
840 catalog.resize(sum(matched))
841 view = visitCatalog.loc[matched]
842 catalog['id'] = view.index
843 for dfCol, afwCol in mapping.items():
844 catalog[afwCol] = view[dfCol]
845
846 self.log.debug("%d sources selected in visit %d detector %d",
847 len(catalog),
848 view['visit'].iloc[0], # all visits in this catalog are the same, so take the first
849 detectorId)
850 return catalog
851
852 def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations,
853 jointcalControl, camera):
854 """Read the data that jointcal needs to run. (Gen3 version)
855
856 Modifies ``associations`` in-place with the loaded data.
857
858 Parameters
859 ----------
860 inputSourceTableVisit : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
861 References to visit-level DataFrames to load the catalog data from.
862 inputVisitSummary : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
863 Visit-level exposure summary catalog with metadata.
864 associations : `lsst.jointcal.Associations`
865 Object to add the loaded data to by constructing new CcdImages.
866 jointcalControl : `jointcal.JointcalControl`
867 Control object for C++ associations management.
869 Camera object for detector geometry.
870
871 Returns
872 -------
873 oldWcsList: `list` [`lsst.afw.geom.SkyWcs`]
874 The original WCS of the input data, to aid in writing tests.
875 bands : `list` [`str`]
876 The filter bands of each input dataset.
877 """
878 oldWcsList = []
879 filters = []
880 load_cat_prof_file = 'jointcal_load_data.prof' if self.config.detailedProfile else ''
881 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
882 table = self._make_schema_table_make_schema_table() # every detector catalog has the same layout
883 # No guarantee that the input is in the same order of visits, so we have to map one of them.
884 catalogMap = {ref.dataId['visit']: i for i, ref in enumerate(inputSourceTableVisit)}
885 detectorDict = {detector.getId(): detector for detector in camera}
886
887 columns = None
888
889 for visitSummaryRef in inputVisitSummary:
890 visitSummary = visitSummaryRef.get()
891
892 dataRef = inputSourceTableVisit[catalogMap[visitSummaryRef.dataId['visit']]]
893 if columns is None:
894 inColumns = dataRef.get(component='columns')
895 columns, detColumn, ixxColumns = self._get_sourceTable_visit_columns_get_sourceTable_visit_columns(inColumns)
896 visitCatalog = dataRef.get(parameters={'columns': columns})
897
898 selected = self.sourceSelector.run(visitCatalog)
899
900 # Build a CcdImage for each detector in this visit.
901 detectors = {id: index for index, id in enumerate(visitSummary['id'])}
902 for id, index in detectors.items():
903 catalog = self._extract_detector_catalog_from_visit_catalog_extract_detector_catalog_from_visit_catalog(table, selected.sourceCat, id,
904 detColumn, ixxColumns)
905 data = self._make_one_input_data_make_one_input_data(visitSummary[index], catalog, detectorDict)
906 result = self._build_ccdImage_build_ccdImage(data, associations, jointcalControl)
907 if result is not None:
908 oldWcsList.append(result.wcs)
909 # A visit has only one band, so we can just use the first.
910 filters.append(data.filter)
911 if len(filters) == 0:
912 raise RuntimeError("No data to process: did source selector remove all sources?")
913 filters = collections.Counter(filters)
914
915 return oldWcsList, filters
916
917 def _make_one_input_data(self, visitRecord, catalog, detectorDict):
918 """Return a data structure for this detector+visit."""
919 return JointcalInputData(visit=visitRecord['visit'],
920 catalog=catalog,
921 visitInfo=visitRecord.getVisitInfo(),
922 detector=detectorDict[visitRecord.getId()],
923 photoCalib=visitRecord.getPhotoCalib(),
924 wcs=visitRecord.getWcs(),
925 bbox=visitRecord.getBBox(),
926 # ExposureRecord doesn't have a FilterLabel yet,
927 # so we have to make one.
928 filter=lsst.afw.image.FilterLabel(band=visitRecord['band'],
929 physical=visitRecord['physical_filter']))
930
931 def _get_sourceTable_visit_columns(self, inColumns):
932 """
933 Get the sourceTable_visit columns from the config.
934
935 Parameters
936 ----------
937 inColumns : `list`
938 List of columns available in the sourceTable_visit
939
940 Returns
941 -------
942 columns : `list`
943 List of columns to read from sourceTable_visit.
944 detectorColumn : `str`
945 Name of the detector column.
946 ixxColumns : `list`
947 Name of the ixx/iyy/ixy columns.
948 """
949 if 'detector' in inColumns:
950 # Default name for Gen3.
951 detectorColumn = 'detector'
952 else:
953 # Default name for Gen2 and Gen2 conversions.
954 detectorColumn = 'ccd'
955
956 columns = ['visit', detectorColumn,
957 'sourceId', 'x', 'xErr', 'y', 'yErr',
958 self.config.sourceFluxType + '_instFlux', self.config.sourceFluxType + '_instFluxErr']
959
960 if 'ixx' in inColumns:
961 # New columns post-DM-31825
962 ixxColumns = ['ixx', 'iyy', 'ixy']
963 else:
964 # Old columns pre-DM-31825
965 ixxColumns = ['Ixx', 'Iyy', 'Ixy']
966 columns.extend(ixxColumns)
967
968 if self.sourceSelector.config.doFlags:
969 columns.extend(self.sourceSelector.config.flags.bad)
970 if self.sourceSelector.config.doUnresolved:
971 columns.append(self.sourceSelector.config.unresolved.name)
972 if self.sourceSelector.config.doIsolated:
973 columns.append(self.sourceSelector.config.isolated.parentName)
974 columns.append(self.sourceSelector.config.isolated.nChildName)
975
976 return columns, detectorColumn, ixxColumns
977
978 # We don't currently need to persist the metadata.
979 # If we do in the future, we will have to add appropriate dataset templates
980 # to each obs package (the metadata template should look like `jointcal_wcs`).
981 def _getMetadataName(self):
982 return None
983
984 @classmethod
985 def _makeArgumentParser(cls):
986 """Create an argument parser"""
987 parser = pipeBase.ArgumentParser(name=cls._DefaultName_DefaultName)
988 parser.add_id_argument("--id", "calexp", help="data ID, e.g. --id visit=6789 ccd=0..9",
989 ContainerClass=PerTractCcdDataIdContainer)
990 return parser
991
992 def _build_ccdImage(self, data, associations, jointcalControl):
993 """
994 Extract the necessary things from this catalog+metadata to add a new
995 ccdImage.
996
997 Parameters
998 ----------
999 data : `JointcalInputData`
1000 The loaded input data.
1001 associations : `lsst.jointcal.Associations`
1002 Object to add the info to, to construct a new CcdImage
1003 jointcalControl : `jointcal.JointcalControl`
1004 Control object for associations management
1005
1006 Returns
1007 ------
1008 namedtuple or `None`
1009 ``wcs``
1010 The TAN WCS of this image, read from the calexp
1012 ``key``
1013 A key to identify this dataRef by its visit and ccd ids
1014 (`namedtuple`).
1015 `None`
1016 if there are no sources in the loaded catalog.
1017 """
1018 if len(data.catalog) == 0:
1019 self.log.warn("No sources selected in visit %s ccd %s", data.visit, data.detector.getId())
1020 return None
1021
1022 associations.createCcdImage(data.catalog,
1023 data.wcs,
1024 data.visitInfo,
1025 data.bbox,
1026 data.filter.physicalLabel,
1027 data.photoCalib,
1028 data.detector,
1029 data.visit,
1030 data.detector.getId(),
1031 jointcalControl)
1032
1033 Result = collections.namedtuple('Result_from_build_CcdImage', ('wcs', 'key'))
1034 Key = collections.namedtuple('Key', ('visit', 'ccd'))
1035 return Result(data.wcs, Key(data.visit, data.detector.getId()))
1036
1037 def _readDataId(self, butler, dataId):
1038 """Read all of the data for one dataId from the butler. (gen2 version)"""
1039 # Not all instruments have `visit` in their dataIds.
1040 if "visit" in dataId.keys():
1041 visit = dataId["visit"]
1042 else:
1043 visit = butler.getButler().queryMetadata("calexp", ("visit"), butler.dataId)[0]
1044 detector = butler.get('calexp_detector', dataId=dataId)
1045
1046 catalog = butler.get('src',
1047 flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS,
1048 dataId=dataId)
1049 goodSrc = self.sourceSelector.run(catalog)
1050 self.log.debug("%d sources selected in visit %d detector %d",
1051 len(goodSrc.sourceCat),
1052 visit,
1053 detector.getId())
1054 return JointcalInputData(visit=visit,
1055 catalog=goodSrc.sourceCat,
1056 visitInfo=butler.get('calexp_visitInfo', dataId=dataId),
1057 detector=detector,
1058 photoCalib=butler.get('calexp_photoCalib', dataId=dataId),
1059 wcs=butler.get('calexp_wcs', dataId=dataId),
1060 bbox=butler.get('calexp_bbox', dataId=dataId),
1061 filter=butler.get('calexp_filterLabel', dataId=dataId))
1062
1063 def loadData(self, dataRefs, associations, jointcalControl):
1064 """Read the data that jointcal needs to run. (Gen2 version)"""
1065 visit_ccd_to_dataRef = {}
1066 oldWcsList = []
1067 filters = []
1068 load_cat_prof_file = 'jointcal_loadData.prof' if self.config.detailedProfile else ''
1069 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
1070 # Need the bounding-box of the focal plane (the same for all visits) for photometry visit models
1071 camera = dataRefs[0].get('camera', immediate=True)
1072 self.focalPlaneBBoxfocalPlaneBBox = camera.getFpBBox()
1073 for dataRef in dataRefs:
1074 data = self._readDataId_readDataId(dataRef.getButler(), dataRef.dataId)
1075 result = self._build_ccdImage_build_ccdImage(data, associations, jointcalControl)
1076 if result is None:
1077 continue
1078 oldWcsList.append(result.wcs)
1079 visit_ccd_to_dataRef[result.key] = dataRef
1080 filters.append(data.filter)
1081 if len(filters) == 0:
1082 raise RuntimeError("No data to process: did source selector remove all sources?")
1083 filters = collections.Counter(filters)
1084
1085 return oldWcsList, filters, visit_ccd_to_dataRef
1086
1087 def _getDebugPath(self, filename):
1088 """Constructs a path to filename using the configured debug path.
1089 """
1090 return os.path.join(self.config.debugOutputPath, filename)
1091
1092 def _prep_sky(self, associations, filters):
1093 """Prepare on-sky and other data that must be computed after data has
1094 been read.
1095 """
1096 associations.computeCommonTangentPoint()
1097
1098 boundingCircle = associations.computeBoundingCircle()
1099 center = lsst.geom.SpherePoint(boundingCircle.getCenter())
1100 radius = lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
1101
1102 self.log.info(f"Data has center={center} with radius={radius.asDegrees()} degrees.")
1103
1104 # Determine a default filter band associated with the catalog. See DM-9093
1105 defaultFilter = filters.most_common(1)[0][0]
1106 self.log.debug("Using '%s' filter for reference flux", defaultFilter.physicalLabel)
1107
1108 # compute and set the reference epoch of the observations, for proper motion corrections
1109 epoch = self._compute_proper_motion_epoch_compute_proper_motion_epoch(associations.getCcdImageList())
1110 associations.setEpoch(epoch.jyear)
1111
1112 return boundingCircle, center, radius, defaultFilter, epoch
1113
1114 @pipeBase.timeMethod
1115 def runDataRef(self, dataRefs):
1116 """
1117 Jointly calibrate the astrometry and photometry across a set of images.
1118
1119 NOTE: this is for gen2 middleware only.
1120
1121 Parameters
1122 ----------
1123 dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef`
1124 List of data references to the exposures to be fit.
1125
1126 Returns
1127 -------
1128 result : `lsst.pipe.base.Struct`
1129 Struct of metadata from the fit, containing:
1130
1131 ``dataRefs``
1132 The provided data references that were fit (with updated WCSs)
1133 ``oldWcsList``
1134 The original WCS from each dataRef
1135 ``metrics``
1136 Dictionary of internally-computed metrics for testing/validation.
1137 """
1138 if len(dataRefs) == 0:
1139 raise ValueError('Need a non-empty list of data references!')
1140
1141 exitStatus = 0 # exit status for shell
1142
1143 sourceFluxField = "slot_%sFlux" % (self.config.sourceFluxType,)
1144 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
1145 associations = lsst.jointcal.Associations()
1146
1147 oldWcsList, filters, visit_ccd_to_dataRef = self.loadDataloadData(dataRefs,
1148 associations,
1149 jointcalControl)
1150
1151 boundingCircle, center, radius, defaultFilter, epoch = self._prep_sky_prep_sky(associations, filters)
1152
1153 tract = dataRefs[0].dataId['tract']
1154
1155 if self.config.doAstrometry:
1156 astrometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
1157 name="astrometry",
1158 refObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
1159 referenceSelector=self.astrometryReferenceSelector,
1160 fit_function=self._fit_astrometry_fit_astrometry,
1161 tract=tract,
1162 epoch=epoch)
1163 self._write_astrometry_results_write_astrometry_results(associations, astrometry.model, visit_ccd_to_dataRef)
1164 else:
1165 astrometry = Astrometry(None, None, None)
1166
1167 if self.config.doPhotometry:
1168 photometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
1169 name="photometry",
1170 refObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
1171 referenceSelector=self.photometryReferenceSelector,
1172 fit_function=self._fit_photometry_fit_photometry,
1173 tract=tract,
1174 epoch=epoch,
1175 reject_bad_fluxes=True)
1176 self._write_photometry_results_write_photometry_results(associations, photometry.model, visit_ccd_to_dataRef)
1177 else:
1178 photometry = Photometry(None, None)
1179
1180 return pipeBase.Struct(dataRefs=dataRefs,
1181 oldWcsList=oldWcsList,
1182 job=self.jobjob,
1183 astrometryRefObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
1184 photometryRefObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
1185 defaultFilter=defaultFilter,
1186 epoch=epoch,
1187 exitStatus=exitStatus)
1188
1189 def _get_refcat_coordinate_error_override(self, refCat, name):
1190 """Check whether we should override the refcat coordinate errors, and
1191 return the overridden error if necessary.
1192
1193 Parameters
1194 ----------
1196 The reference catalog to check for a ``coord_raErr`` field.
1197 name : `str`
1198 Whether we are doing "astrometry" or "photometry".
1199
1200 Returns
1201 -------
1202 refCoordErr : `float`
1203 The refcat coordinate error to use, or NaN if we are not overriding
1204 those fields.
1205
1206 Raises
1207 ------
1208 lsst.pex.config.FieldValidationError
1209 Raised if the refcat does not contain coordinate errors and
1210 ``config.astrometryReferenceErr`` is not set.
1211 """
1212 # This value doesn't matter for photometry, so just set something to
1213 # keep old refcats from causing problems.
1214 if name.lower() == "photometry":
1215 if 'coord_raErr' not in refCat.schema:
1216 return 100
1217 else:
1218 return float('nan')
1219
1220 if self.config.astrometryReferenceErr is None and 'coord_raErr' not in refCat.schema:
1221 msg = ("Reference catalog does not contain coordinate errors, "
1222 "and config.astrometryReferenceErr not supplied.")
1223 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
1224 self.config,
1225 msg)
1226
1227 if self.config.astrometryReferenceErr is not None and 'coord_raErr' in refCat.schema:
1228 self.log.warn("Overriding reference catalog coordinate errors with %f/coordinate [mas]",
1229 self.config.astrometryReferenceErr)
1230
1231 if self.config.astrometryReferenceErr is None:
1232 return float('nan')
1233 else:
1234 return self.config.astrometryReferenceErr
1235
1236 def _compute_proper_motion_epoch(self, ccdImageList):
1237 """Return the proper motion correction epoch of the provided images.
1238
1239 Parameters
1240 ----------
1241 ccdImageList : `list` [`lsst.jointcal.CcdImage`]
1242 The images to compute the appropriate epoch for.
1243
1244 Returns
1245 -------
1246 epoch : `astropy.time.Time`
1247 The date to use for proper motion corrections.
1248 """
1249 return astropy.time.Time(np.mean([ccdImage.getEpoch() for ccdImage in ccdImageList]),
1250 format="jyear",
1251 scale="tai")
1252
1253 def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
1254 tract="", match_cut=3.0,
1255 reject_bad_fluxes=False, *,
1256 name="", refObjLoader=None, referenceSelector=None,
1257 fit_function=None, epoch=None):
1258 """Load reference catalog, perform the fit, and return the result.
1259
1260 Parameters
1261 ----------
1262 associations : `lsst.jointcal.Associations`
1263 The star/reference star associations to fit.
1264 defaultFilter : `lsst.afw.image.FilterLabel`
1265 filter to load from reference catalog.
1266 center : `lsst.geom.SpherePoint`
1267 ICRS center of field to load from reference catalog.
1268 radius : `lsst.geom.Angle`
1269 On-sky radius to load from reference catalog.
1270 name : `str`
1271 Name of thing being fit: "astrometry" or "photometry".
1272 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
1273 Reference object loader to use to load a reference catalog.
1274 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1275 Selector to use to pick objects from the loaded reference catalog.
1276 fit_function : callable
1277 Function to call to perform fit (takes Associations object).
1278 tract : `str`, optional
1279 Name of tract currently being fit.
1280 match_cut : `float`, optional
1281 Radius in arcseconds to find cross-catalog matches to during
1282 associations.associateCatalogs.
1283 reject_bad_fluxes : `bool`, optional
1284 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
1285 epoch : `astropy.time.Time`, optional
1286 Epoch to which to correct refcat proper motion and parallax,
1287 or `None` to not apply such corrections.
1288
1289 Returns
1290 -------
1291 result : `Photometry` or `Astrometry`
1292 Result of `fit_function()`
1293 """
1294 self.log.info("====== Now processing %s...", name)
1295 # TODO: this should not print "trying to invert a singular transformation:"
1296 # if it does that, something's not right about the WCS...
1297 associations.associateCatalogs(match_cut)
1298 add_measurement(self.jobjob, 'jointcal.associated_%s_fittedStars' % name,
1299 associations.fittedStarListSize())
1300
1301 applyColorterms = False if name.lower() == "astrometry" else self.config.applyColorTerms
1302 refCat, fluxField = self._load_reference_catalog_load_reference_catalog(refObjLoader, referenceSelector,
1303 center, radius, defaultFilter,
1304 applyColorterms=applyColorterms,
1305 epoch=epoch)
1306 refCoordErr = self._get_refcat_coordinate_error_override_get_refcat_coordinate_error_override(refCat, name)
1307
1308 associations.collectRefStars(refCat,
1309 self.config.matchCut*lsst.geom.arcseconds,
1310 fluxField,
1311 refCoordinateErr=refCoordErr,
1312 rejectBadFluxes=reject_bad_fluxes)
1313 add_measurement(self.jobjob, 'jointcal.collected_%s_refStars' % name,
1314 associations.refStarListSize())
1315
1316 associations.prepareFittedStars(self.config.minMeasurements)
1317
1318 self._check_star_lists_check_star_lists(associations, name)
1319 add_measurement(self.jobjob, 'jointcal.selected_%s_refStars' % name,
1320 associations.nFittedStarsWithAssociatedRefStar())
1321 add_measurement(self.jobjob, 'jointcal.selected_%s_fittedStars' % name,
1322 associations.fittedStarListSize())
1323 add_measurement(self.jobjob, 'jointcal.selected_%s_ccdImages' % name,
1324 associations.nCcdImagesValidForFit())
1325
1326 load_cat_prof_file = 'jointcal_fit_%s.prof'%name if self.config.detailedProfile else ''
1327 dataName = "{}_{}".format(tract, defaultFilter.physicalLabel)
1328 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
1329 result = fit_function(associations, dataName)
1330 # TODO DM-12446: turn this into a "butler save" somehow.
1331 # Save reference and measurement chi2 contributions for this data
1332 if self.config.writeChi2FilesInitialFinal:
1333 baseName = self._getDebugPath_getDebugPath(f"{name}_final_chi2-{dataName}")
1334 result.fit.saveChi2Contributions(baseName+"{type}")
1335 self.log.info("Wrote chi2 contributions files: %s", baseName)
1336
1337 return result
1338
1339 def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel,
1340 applyColorterms=False, epoch=None):
1341 """Load the necessary reference catalog sources, convert fluxes to
1342 correct units, and apply color term corrections if requested.
1343
1344 Parameters
1345 ----------
1346 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
1347 The reference catalog loader to use to get the data.
1348 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1349 Source selector to apply to loaded reference catalog.
1350 center : `lsst.geom.SpherePoint`
1351 The center around which to load sources.
1352 radius : `lsst.geom.Angle`
1353 The radius around ``center`` to load sources in.
1354 filterLabel : `lsst.afw.image.FilterLabel`
1355 The camera filter to load fluxes for.
1356 applyColorterms : `bool`
1357 Apply colorterm corrections to the refcat for ``filterName``?
1358 epoch : `astropy.time.Time`, optional
1359 Epoch to which to correct refcat proper motion and parallax,
1360 or `None` to not apply such corrections.
1361
1362 Returns
1363 -------
1365 The loaded reference catalog.
1366 fluxField : `str`
1367 The name of the reference catalog flux field appropriate for ``filterName``.
1368 """
1369 skyCircle = refObjLoader.loadSkyCircle(center,
1370 radius,
1371 filterLabel.bandLabel,
1372 epoch=epoch)
1373
1374 selected = referenceSelector.run(skyCircle.refCat)
1375 # Need memory contiguity to get reference filters as a vector.
1376 if not selected.sourceCat.isContiguous():
1377 refCat = selected.sourceCat.copy(deep=True)
1378 else:
1379 refCat = selected.sourceCat
1380
1381 if applyColorterms:
1382 refCatName = refObjLoader.ref_dataset_name
1383 self.log.info("Applying color terms for physical filter=%r reference catalog=%s",
1384 filterLabel.physicalLabel, refCatName)
1385 colorterm = self.config.colorterms.getColorterm(filterLabel.physicalLabel,
1386 refCatName,
1387 doRaise=True)
1388
1389 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat)
1390 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
1391 # TODO: I didn't want to use this, but I'll deal with it in DM-16903
1392 refCat[skyCircle.fluxField+'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
1393
1394 return refCat, skyCircle.fluxField
1395
1396 def _check_star_lists(self, associations, name):
1397 # TODO: these should be len(blah), but we need this properly wrapped first.
1398 if associations.nCcdImagesValidForFit() == 0:
1399 raise RuntimeError('No images in the ccdImageList!')
1400 if associations.fittedStarListSize() == 0:
1401 raise RuntimeError('No stars in the {} fittedStarList!'.format(name))
1402 if associations.refStarListSize() == 0:
1403 raise RuntimeError('No stars in the {} reference star list!'.format(name))
1404
1405 def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
1406 """Compute chi2, log it, validate the model, and return chi2.
1407
1408 Parameters
1409 ----------
1410 associations : `lsst.jointcal.Associations`
1411 The star/reference star associations to fit.
1413 The fitter to use for minimization.
1414 model : `lsst.jointcal.Model`
1415 The model being fit.
1416 chi2Label : `str`
1417 Label to describe the chi2 (e.g. "Initialized", "Final").
1418 writeChi2Name : `str`, optional
1419 Filename prefix to write the chi2 contributions to.
1420 Do not supply an extension: an appropriate one will be added.
1421
1422 Returns
1423 -------
1425 The chi2 object for the current fitter and model.
1426
1427 Raises
1428 ------
1429 FloatingPointError
1430 Raised if chi2 is infinite or NaN.
1431 ValueError
1432 Raised if the model is not valid.
1433 """
1434 if writeChi2Name is not None:
1435 fullpath = self._getDebugPath_getDebugPath(writeChi2Name)
1436 fit.saveChi2Contributions(fullpath+"{type}")
1437 self.log.info("Wrote chi2 contributions files: %s", fullpath)
1438
1439 chi2 = fit.computeChi2()
1440 self.log.info("%s %s", chi2Label, chi2)
1441 self._check_stars_check_stars(associations)
1442 if not np.isfinite(chi2.chi2):
1443 raise FloatingPointError(f'{chi2Label} chi2 is invalid: {chi2}')
1444 if not model.validate(associations.getCcdImageList(), chi2.ndof):
1445 raise ValueError("Model is not valid: check log messages for warnings.")
1446 return chi2
1447
1448 def _fit_photometry(self, associations, dataName=None):
1449 """
1450 Fit the photometric data.
1451
1452 Parameters
1453 ----------
1454 associations : `lsst.jointcal.Associations`
1455 The star/reference star associations to fit.
1456 dataName : `str`
1457 Name of the data being processed (e.g. "1234_HSC-Y"), for
1458 identifying debugging files.
1459
1460 Returns
1461 -------
1462 fit_result : `namedtuple`
1464 The photometric fitter used to perform the fit.
1466 The photometric model that was fit.
1467 """
1468 self.log.info("=== Starting photometric fitting...")
1469
1470 # TODO: should use pex.config.RegistryField here (see DM-9195)
1471 if self.config.photometryModel == "constrainedFlux":
1472 model = lsst.jointcal.ConstrainedFluxModel(associations.getCcdImageList(),
1473 self.focalPlaneBBoxfocalPlaneBBox,
1474 visitOrder=self.config.photometryVisitOrder,
1475 errorPedestal=self.config.photometryErrorPedestal)
1476 # potentially nonlinear problem, so we may need a line search to converge.
1477 doLineSearch = self.config.allowLineSearch
1478 elif self.config.photometryModel == "constrainedMagnitude":
1479 model = lsst.jointcal.ConstrainedMagnitudeModel(associations.getCcdImageList(),
1480 self.focalPlaneBBoxfocalPlaneBBox,
1481 visitOrder=self.config.photometryVisitOrder,
1482 errorPedestal=self.config.photometryErrorPedestal)
1483 # potentially nonlinear problem, so we may need a line search to converge.
1484 doLineSearch = self.config.allowLineSearch
1485 elif self.config.photometryModel == "simpleFlux":
1486 model = lsst.jointcal.SimpleFluxModel(associations.getCcdImageList(),
1487 errorPedestal=self.config.photometryErrorPedestal)
1488 doLineSearch = False # purely linear in model parameters, so no line search needed
1489 elif self.config.photometryModel == "simpleMagnitude":
1490 model = lsst.jointcal.SimpleMagnitudeModel(associations.getCcdImageList(),
1491 errorPedestal=self.config.photometryErrorPedestal)
1492 doLineSearch = False # purely linear in model parameters, so no line search needed
1493
1494 fit = lsst.jointcal.PhotometryFit(associations, model)
1495 # TODO DM-12446: turn this into a "butler save" somehow.
1496 # Save reference and measurement chi2 contributions for this data
1497 if self.config.writeChi2FilesInitialFinal:
1498 baseName = f"photometry_initial_chi2-{dataName}"
1499 else:
1500 baseName = None
1501 if self.config.writeInitialModel:
1502 fullpath = self._getDebugPath_getDebugPath(f"initial_photometry_model-{dataName}.txt")
1503 writeModel(model, fullpath, self.log)
1504 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialized", writeChi2Name=baseName)
1505
1506 def getChi2Name(whatToFit):
1507 if self.config.writeChi2FilesOuterLoop:
1508 return f"photometry_init-%s_chi2-{dataName}" % whatToFit
1509 else:
1510 return None
1511
1512 # The constrained model needs the visit transform fit first; the chip
1513 # transform is initialized from the singleFrame PhotoCalib, so it's close.
1514 if self.config.writeInitMatrix:
1515 dumpMatrixFile = self._getDebugPath_getDebugPath(f"photometry_preinit-{dataName}")
1516 else:
1517 dumpMatrixFile = ""
1518 if self.config.photometryModel.startswith("constrained"):
1519 # no line search: should be purely (or nearly) linear,
1520 # and we want a large step size to initialize with.
1521 fit.minimize("ModelVisit", dumpMatrixFile=dumpMatrixFile)
1522 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize ModelVisit",
1523 writeChi2Name=getChi2Name("ModelVisit"))
1524 dumpMatrixFile = "" # so we don't redo the output on the next step
1525
1526 fit.minimize("Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
1527 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Model",
1528 writeChi2Name=getChi2Name("Model"))
1529
1530 fit.minimize("Fluxes") # no line search: always purely linear.
1531 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Fluxes",
1532 writeChi2Name=getChi2Name("Fluxes"))
1533
1534 fit.minimize("Model Fluxes", doLineSearch=doLineSearch)
1535 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize ModelFluxes",
1536 writeChi2Name=getChi2Name("ModelFluxes"))
1537
1538 model.freezeErrorTransform()
1539 self.log.debug("Photometry error scales are frozen.")
1540
1541 chi2 = self._iterate_fit_iterate_fit(associations,
1542 fit,
1543 self.config.maxPhotometrySteps,
1544 "photometry",
1545 "Model Fluxes",
1546 doRankUpdate=self.config.photometryDoRankUpdate,
1547 doLineSearch=doLineSearch,
1548 dataName=dataName)
1549
1550 add_measurement(self.jobjob, 'jointcal.photometry_final_chi2', chi2.chi2)
1551 add_measurement(self.jobjob, 'jointcal.photometry_final_ndof', chi2.ndof)
1552 return Photometry(fit, model)
1553
1554 def _fit_astrometry(self, associations, dataName=None):
1555 """
1556 Fit the astrometric data.
1557
1558 Parameters
1559 ----------
1560 associations : `lsst.jointcal.Associations`
1561 The star/reference star associations to fit.
1562 dataName : `str`
1563 Name of the data being processed (e.g. "1234_HSC-Y"), for
1564 identifying debugging files.
1565
1566 Returns
1567 -------
1568 fit_result : `namedtuple`
1570 The astrometric fitter used to perform the fit.
1572 The astrometric model that was fit.
1573 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1574 The model for the sky to tangent plane projection that was used in the fit.
1575 """
1576
1577 self.log.info("=== Starting astrometric fitting...")
1578
1579 associations.deprojectFittedStars()
1580
1581 # NOTE: need to return sky_to_tan_projection so that it doesn't get garbage collected.
1582 # TODO: could we package sky_to_tan_projection and model together so we don't have to manage
1583 # them so carefully?
1584 sky_to_tan_projection = lsst.jointcal.OneTPPerVisitHandler(associations.getCcdImageList())
1585
1586 if self.config.astrometryModel == "constrained":
1587 model = lsst.jointcal.ConstrainedAstrometryModel(associations.getCcdImageList(),
1588 sky_to_tan_projection,
1589 chipOrder=self.config.astrometryChipOrder,
1590 visitOrder=self.config.astrometryVisitOrder)
1591 elif self.config.astrometryModel == "simple":
1592 model = lsst.jointcal.SimpleAstrometryModel(associations.getCcdImageList(),
1593 sky_to_tan_projection,
1594 self.config.useInputWcs,
1595 nNotFit=0,
1596 order=self.config.astrometrySimpleOrder)
1597
1598 fit = lsst.jointcal.AstrometryFit(associations, model, self.config.positionErrorPedestal)
1599 # TODO DM-12446: turn this into a "butler save" somehow.
1600 # Save reference and measurement chi2 contributions for this data
1601 if self.config.writeChi2FilesInitialFinal:
1602 baseName = f"astrometry_initial_chi2-{dataName}"
1603 else:
1604 baseName = None
1605 if self.config.writeInitialModel:
1606 fullpath = self._getDebugPath_getDebugPath(f"initial_astrometry_model-{dataName}.txt")
1607 writeModel(model, fullpath, self.log)
1608 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initial", writeChi2Name=baseName)
1609
1610 def getChi2Name(whatToFit):
1611 if self.config.writeChi2FilesOuterLoop:
1612 return f"astrometry_init-%s_chi2-{dataName}" % whatToFit
1613 else:
1614 return None
1615
1616 if self.config.writeInitMatrix:
1617 dumpMatrixFile = self._getDebugPath_getDebugPath(f"astrometry_preinit-{dataName}")
1618 else:
1619 dumpMatrixFile = ""
1620 # The constrained model needs the visit transform fit first; the chip
1621 # transform is initialized from the detector's cameraGeom, so it's close.
1622 if self.config.astrometryModel == "constrained":
1623 fit.minimize("DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1624 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize DistortionsVisit",
1625 writeChi2Name=getChi2Name("DistortionsVisit"))
1626 dumpMatrixFile = "" # so we don't redo the output on the next step
1627
1628 fit.minimize("Distortions", dumpMatrixFile=dumpMatrixFile)
1629 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Distortions",
1630 writeChi2Name=getChi2Name("Distortions"))
1631
1632 fit.minimize("Positions")
1633 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Positions",
1634 writeChi2Name=getChi2Name("Positions"))
1635
1636 fit.minimize("Distortions Positions")
1637 self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize DistortionsPositions",
1638 writeChi2Name=getChi2Name("DistortionsPositions"))
1639
1640 chi2 = self._iterate_fit_iterate_fit(associations,
1641 fit,
1642 self.config.maxAstrometrySteps,
1643 "astrometry",
1644 "Distortions Positions",
1645 sigmaRelativeTolerance=self.config.astrometryOutlierRelativeTolerance,
1646 doRankUpdate=self.config.astrometryDoRankUpdate,
1647 dataName=dataName)
1648
1649 add_measurement(self.jobjob, 'jointcal.astrometry_final_chi2', chi2.chi2)
1650 add_measurement(self.jobjob, 'jointcal.astrometry_final_ndof', chi2.ndof)
1651
1652 return Astrometry(fit, model, sky_to_tan_projection)
1653
1654 def _check_stars(self, associations):
1655 """Count measured and reference stars per ccd and warn/log them."""
1656 for ccdImage in associations.getCcdImageList():
1657 nMeasuredStars, nRefStars = ccdImage.countStars()
1658 self.log.debug("ccdImage %s has %s measured and %s reference stars",
1659 ccdImage.getName(), nMeasuredStars, nRefStars)
1660 if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1661 self.log.warn("ccdImage %s has only %s measuredStars (desired %s)",
1662 ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1663 if nRefStars < self.config.minRefStarsPerCcd:
1664 self.log.warn("ccdImage %s has only %s RefStars (desired %s)",
1665 ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1666
1667 def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1668 dataName="",
1669 sigmaRelativeTolerance=0,
1670 doRankUpdate=True,
1671 doLineSearch=False):
1672 """Run fitter.minimize up to max_steps times, returning the final chi2.
1673
1674 Parameters
1675 ----------
1676 associations : `lsst.jointcal.Associations`
1677 The star/reference star associations to fit.
1678 fitter : `lsst.jointcal.FitterBase`
1679 The fitter to use for minimization.
1680 max_steps : `int`
1681 Maximum number of steps to run outlier rejection before declaring
1682 convergence failure.
1683 name : {'photometry' or 'astrometry'}
1684 What type of data are we fitting (for logs and debugging files).
1685 whatToFit : `str`
1686 Passed to ``fitter.minimize()`` to define the parameters to fit.
1687 dataName : `str`, optional
1688 Descriptive name for this dataset (e.g. tract and filter),
1689 for debugging.
1690 sigmaRelativeTolerance : `float`, optional
1691 Convergence tolerance for the fractional change in the chi2 cut
1692 level for determining outliers. If set to zero, iterations will
1693 continue until there are no outliers.
1694 doRankUpdate : `bool`, optional
1695 Do an Eigen rank update during minimization, or recompute the full
1696 matrix and gradient?
1697 doLineSearch : `bool`, optional
1698 Do a line search for the optimum step during minimization?
1699
1700 Returns
1701 -------
1703 The final chi2 after the fit converges, or is forced to end.
1704
1705 Raises
1706 ------
1707 FloatingPointError
1708 Raised if the fitter fails with a non-finite value.
1709 RuntimeError
1710 Raised if the fitter fails for some other reason;
1711 log messages will provide further details.
1712 """
1713 if self.config.writeInitMatrix:
1714 dumpMatrixFile = self._getDebugPath_getDebugPath(f"{name}_postinit-{dataName}")
1715 else:
1716 dumpMatrixFile = ""
1717 oldChi2 = lsst.jointcal.Chi2Statistic()
1718 oldChi2.chi2 = float("inf")
1719 for i in range(max_steps):
1720 if self.config.writeChi2FilesOuterLoop:
1721 writeChi2Name = f"{name}_iterate_{i}_chi2-{dataName}"
1722 else:
1723 writeChi2Name = None
1724 result = fitter.minimize(whatToFit,
1725 self.config.outlierRejectSigma,
1726 sigmaRelativeTolerance=sigmaRelativeTolerance,
1727 doRankUpdate=doRankUpdate,
1728 doLineSearch=doLineSearch,
1729 dumpMatrixFile=dumpMatrixFile)
1730 dumpMatrixFile = "" # clear it so we don't write the matrix again.
1731 chi2 = self._logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(),
1732 f"Fit iteration {i}", writeChi2Name=writeChi2Name)
1733
1734 if result == MinimizeResult.Converged:
1735 if doRankUpdate:
1736 self.log.debug("fit has converged - no more outliers - redo minimization "
1737 "one more time in case we have lost accuracy in rank update.")
1738 # Redo minimization one more time in case we have lost accuracy in rank update
1739 result = fitter.minimize(whatToFit, self.config.outlierRejectSigma,
1740 sigmaRelativeTolerance=sigmaRelativeTolerance)
1741 chi2 = self._logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(), "Fit completed")
1742
1743 # log a message for a large final chi2, TODO: DM-15247 for something better
1744 if chi2.chi2/chi2.ndof >= 4.0:
1745 self.log.error("Potentially bad fit: High chi-squared/ndof.")
1746
1747 break
1748 elif result == MinimizeResult.Chi2Increased:
1749 self.log.warn("Still some outliers remaining but chi2 increased - retry")
1750 # Check whether the increase was large enough to cause trouble.
1751 chi2Ratio = chi2.chi2 / oldChi2.chi2
1752 if chi2Ratio > 1.5:
1753 self.log.warn('Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1754 chi2.chi2, oldChi2.chi2, chi2Ratio)
1755 # Based on a variety of HSC jointcal logs (see DM-25779), it
1756 # appears that chi2 increases more than a factor of ~2 always
1757 # result in the fit diverging rapidly and ending at chi2 > 1e10.
1758 # Using 10 as the "failure" threshold gives some room between
1759 # leaving a warning and bailing early.
1760 if chi2Ratio > 10:
1761 msg = ("Large chi2 increase between steps: fit likely cannot converge."
1762 " Try setting one or more of the `writeChi2*` config fields and looking"
1763 " at how individual star chi2-values evolve during the fit.")
1764 raise RuntimeError(msg)
1765 oldChi2 = chi2
1766 elif result == MinimizeResult.NonFinite:
1767 filename = self._getDebugPath_getDebugPath("{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1768 # TODO DM-12446: turn this into a "butler save" somehow.
1769 fitter.saveChi2Contributions(filename+"{type}")
1770 msg = "Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1771 raise FloatingPointError(msg.format(filename))
1772 elif result == MinimizeResult.Failed:
1773 raise RuntimeError("Chi2 minimization failure, cannot complete fit.")
1774 else:
1775 raise RuntimeError("Unxepected return code from minimize().")
1776 else:
1777 self.log.error("%s failed to converge after %d steps"%(name, max_steps))
1778
1779 return chi2
1780
1781 def _make_output(self, ccdImageList, model, func):
1782 """Return the internal jointcal models converted to the afw
1783 structures that will be saved to disk.
1784
1785 Parameters
1786 ----------
1787 ccdImageList : `lsst.jointcal.CcdImageList`
1788 The list of CcdImages to get the output for.
1790 The internal jointcal model to convert for each `lsst.jointcal.CcdImage`.
1791 func : `str`
1792 The name of the function to call on ``model`` to get the converted
1793 structure. Must accept an `lsst.jointcal.CcdImage`.
1794
1795 Returns
1796 -------
1797 output : `dict` [`tuple`, `lsst.jointcal.AstrometryModel`] or
1798 `dict` [`tuple`, `lsst.jointcal.PhotometryModel`]
1799 The data to be saved, keyed on (visit, detector).
1800 """
1801 output = {}
1802 for ccdImage in ccdImageList:
1803 ccd = ccdImage.ccdId
1804 visit = ccdImage.visit
1805 self.log.debug("%s for visit: %d, ccd: %d", func, visit, ccd)
1806 output[(visit, ccd)] = getattr(model, func)(ccdImage)
1807 return output
1808
1809 def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1810 """
1811 Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1812
1813 Parameters
1814 ----------
1815 associations : `lsst.jointcal.Associations`
1816 The star/reference star associations to fit.
1818 The astrometric model that was fit.
1819 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1820 Dict of ccdImage identifiers to dataRefs that were fit.
1821 """
1822 ccdImageList = associations.getCcdImageList()
1823 output = self._make_output_make_output(ccdImageList, model, "makeSkyWcs")
1824 for key, skyWcs in output.items():
1825 dataRef = visit_ccd_to_dataRef[key]
1826 try:
1827 dataRef.put(skyWcs, 'jointcal_wcs')
1828 except pexExceptions.Exception as e:
1829 self.log.fatal('Failed to write updated Wcs: %s', str(e))
1830 raise e
1831
1832 def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1833 """
1834 Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1835
1836 Parameters
1837 ----------
1838 associations : `lsst.jointcal.Associations`
1839 The star/reference star associations to fit.
1841 The photoometric model that was fit.
1842 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1843 Dict of ccdImage identifiers to dataRefs that were fit.
1844 """
1845
1846 ccdImageList = associations.getCcdImageList()
1847 output = self._make_output_make_output(ccdImageList, model, "toPhotoCalib")
1848 for key, photoCalib in output.items():
1849 dataRef = visit_ccd_to_dataRef[key]
1850 try:
1851 dataRef.put(photoCalib, 'jointcal_photoCalib')
1852 except pexExceptions.Exception as e:
1853 self.log.fatal('Failed to write updated PhotoCalib: %s', str(e))
1854 raise e
static std::shared_ptr< SourceTable > make(Schema const &schema, std::shared_ptr< IdFactory > const &idFactory)
static Schema makeMinimalSchema()
The class that implements the relations between MeasuredStar and FittedStar.
Definition: Associations.h:54
Class that handles the astrometric least squares problem.
Definition: AstrometryFit.h:78
Interface between AstrometryFit and the combinations of Mappings from pixels to some tangent plane (a...
Handler of an actual image from a single CCD.
Definition: CcdImage.h:64
Base class for Chi2Statistic and Chi2List, to allow addEntry inside Fitter for either class.
Definition: Chi2.h:44
Simple structure to accumulate chi2 and ndof.
Definition: Chi2.h:52
A multi-component model, fitting mappings for sensors and visits simultaneously.
Base class for fitters.
Definition: FitterBase.h:53
A projection handler in which all CCDs from the same visit have the same tangent point.
Class that handles the photometric least squares problem.
Definition: PhotometryFit.h:45
A model where there is one independent transform per CcdImage.
def getTargetList(parsedCmd, **kwargs)
Definition: jointcal.py:77
def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations, jointcalControl, camera)
Definition: jointcal.py:853
def _compute_proper_motion_epoch(self, ccdImageList)
Definition: jointcal.py:1236
def _get_refcat_coordinate_error_override(self, refCat, name)
Definition: jointcal.py:1189
def _getDebugPath(self, filename)
Definition: jointcal.py:1087
def _check_star_lists(self, associations, name)
Definition: jointcal.py:1396
def _make_one_input_data(self, visitRecord, catalog, detectorDict)
Definition: jointcal.py:917
def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", sigmaRelativeTolerance=0, doRankUpdate=True, doLineSearch=False)
Definition: jointcal.py:1671
def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1809
def _check_stars(self, associations)
Definition: jointcal.py:1654
def _prep_sky(self, associations, filters)
Definition: jointcal.py:1092
def _readDataId(self, butler, dataId)
Definition: jointcal.py:1037
def runQuantum(self, butlerQC, inputRefs, outputRefs)
Definition: jointcal.py:638
def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius, tract="", match_cut=3.0, reject_bad_fluxes=False, *name="", refObjLoader=None, referenceSelector=None, fit_function=None, epoch=None)
Definition: jointcal.py:1257
def _build_ccdImage(self, data, associations, jointcalControl)
Definition: jointcal.py:992
def runDataRef(self, dataRefs)
Definition: jointcal.py:1115
def __init__(self, butler=None, initInputs=None, **kwargs)
Definition: jointcal.py:617
def _get_sourceTable_visit_columns(self, inColumns)
Definition: jointcal.py:931
def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None)
Definition: jointcal.py:725
def _fit_photometry(self, associations, dataName=None)
Definition: jointcal.py:1448
def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1832
def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None)
Definition: jointcal.py:1405
def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId, detectorColumn, ixxColumns)
Definition: jointcal.py:803
def _make_output(self, ccdImageList, model, func)
Definition: jointcal.py:1781
def _put_output(self, butlerQC, outputs, outputRefs, camera, setter)
Definition: jointcal.py:666
def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel, applyColorterms=False, epoch=None)
Definition: jointcal.py:1340
def _fit_astrometry(self, associations, dataName=None)
Definition: jointcal.py:1554
def loadData(self, dataRefs, associations, jointcalControl)
Definition: jointcal.py:1063
def lookupStaticCalibrations(datasetType, registry, quantumDataId, collections)
Definition: jointcal.py:141
def add_measurement(job, name, value)
Definition: jointcal.py:59
def lookupVisitRefCats(datasetType, registry, quantumDataId, collections)
Definition: jointcal.py:181
def writeModel(model, filename, log)
Definition: jointcal.py:571
This is a virtual class that allows a lot of freedom in the choice of the projection from "Sky" (wher...