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