lsst.jointcal  21.0.0-15-ge3888ae+eb634caa15
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 
22 import dataclasses
23 import collections
24 import os
25 
26 import astropy.time
27 import numpy as np
28 import astropy.units as u
29 
30 import lsst.geom
31 import lsst.utils
32 import lsst.pex.config as pexConfig
33 import lsst.pipe.base as pipeBase
34 from lsst.afw.image import fluxErrFromABMagErr
35 import lsst.pex.exceptions as pexExceptions
37 import lsst.afw.table
38 import lsst.log
39 from lsst.obs.base import Instrument
40 from lsst.pipe.tasks.colorterms import ColortermLibrary
41 from lsst.verify import Job, Measurement
42 
43 from lsst.meas.algorithms import (LoadIndexedReferenceObjectsTask, ReferenceSourceSelectorTask,
44  ReferenceObjectLoader)
45 from lsst.meas.algorithms.sourceSelector import sourceSelectorRegistry
46 
47 from .dataIds import PerTractCcdDataIdContainer
48 
49 import lsst.jointcal
50 from lsst.jointcal import MinimizeResult
51 
52 __all__ = ["JointcalConfig", "JointcalRunner", "JointcalTask"]
53 
54 Photometry = collections.namedtuple('Photometry', ('fit', 'model'))
55 Astrometry = collections.namedtuple('Astrometry', ('fit', 'model', 'sky_to_tan_projection'))
56 
57 
58 # TODO: move this to MeasurementSet in lsst.verify per DM-12655.
59 def add_measurement(job, name, value):
60  meas = Measurement(job.metrics[name], value)
61  job.measurements.insert(meas)
62 
63 
64 class 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 
141 def 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 
181 def 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 
221 class 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 
304 class 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  maxPhotometrySteps = pexConfig.Field(
440  doc="Maximum number of minimize iterations to take when fitting photometry.",
441  dtype=int,
442  default=20,
443  )
444  maxAstrometrySteps = pexConfig.Field(
445  doc="Maximum number of minimize iterations to take when fitting photometry.",
446  dtype=int,
447  default=20,
448  )
449  astrometryRefObjLoader = pexConfig.ConfigurableField(
450  target=LoadIndexedReferenceObjectsTask,
451  doc="Reference object loader for astrometric fit",
452  )
453  photometryRefObjLoader = pexConfig.ConfigurableField(
454  target=LoadIndexedReferenceObjectsTask,
455  doc="Reference object loader for photometric fit",
456  )
457  sourceSelector = sourceSelectorRegistry.makeField(
458  doc="How to select sources for cross-matching",
459  default="astrometry"
460  )
461  astrometryReferenceSelector = pexConfig.ConfigurableField(
462  target=ReferenceSourceSelectorTask,
463  doc="How to down-select the loaded astrometry reference catalog.",
464  )
465  photometryReferenceSelector = pexConfig.ConfigurableField(
466  target=ReferenceSourceSelectorTask,
467  doc="How to down-select the loaded photometry reference catalog.",
468  )
469  astrometryReferenceErr = pexConfig.Field(
470  doc=("Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
471  "If None, then raise an exception if the reference catalog is missing coordinate errors. "
472  "If specified, overrides any existing `coord_*Err` values."),
473  dtype=float,
474  default=None,
475  optional=True
476  )
477 
478  # configs for outputting debug information
479  writeInitMatrix = pexConfig.Field(
480  dtype=bool,
481  doc=("Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
482  "The output files will be of the form 'astrometry_preinit-mat.txt', in the current directory. "
483  "Note that these files are the dense versions of the matrix, and so may be very large."),
484  default=False
485  )
486  writeChi2FilesInitialFinal = pexConfig.Field(
487  dtype=bool,
488  doc="Write .csv files containing the contributions to chi2 for the initialization and final fit.",
489  default=False
490  )
491  writeChi2FilesOuterLoop = pexConfig.Field(
492  dtype=bool,
493  doc="Write .csv files containing the contributions to chi2 for the outer fit loop.",
494  default=False
495  )
496  writeInitialModel = pexConfig.Field(
497  dtype=bool,
498  doc=("Write the pre-initialization model to text files, for debugging."
499  " Output is written to `initial[Astro|Photo]metryModel.txt` in the current working directory."),
500  default=False
501  )
502  debugOutputPath = pexConfig.Field(
503  dtype=str,
504  default=".",
505  doc=("Path to write debug output files to. Used by "
506  "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
507  )
508  detailedProfile = pexConfig.Field(
509  dtype=bool,
510  default=False,
511  doc="Output separate profiling information for different parts of jointcal, e.g. data read, fitting"
512  )
513 
514  def validate(self):
515  super().validate()
516  if self.doPhotometrydoPhotometry and self.applyColorTermsapplyColorTerms and len(self.colortermscolorterms.data) == 0:
517  msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
518  raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
519  if self.doAstrometrydoAstrometry and not self.doPhotometrydoPhotometry and self.applyColorTermsapplyColorTerms:
520  msg = ("Only doing astrometry, but Colorterms are not applied for astrometry;"
521  "applyColorTerms=True will be ignored.")
522  lsst.log.warn(msg)
523 
524  def setDefaults(self):
525  # Use science source selector which can filter on extendedness, SNR, and whether blended
526  self.sourceSelectorsourceSelector.name = 'science'
527  # Use only stars because aperture fluxes of galaxies are biased and depend on seeing
528  self.sourceSelectorsourceSelector['science'].doUnresolved = True
529  # with dependable signal to noise ratio.
530  self.sourceSelectorsourceSelector['science'].doSignalToNoise = True
531  # Min SNR must be > 0 because jointcal cannot handle negative fluxes,
532  # and S/N > 10 to use sources that are not too faint, and thus better measured.
533  self.sourceSelectorsourceSelector['science'].signalToNoise.minimum = 10.
534  # Base SNR on CalibFlux because that is the flux jointcal that fits and must be positive
535  fluxField = f"slot_{self.sourceFluxType}Flux_instFlux"
536  self.sourceSelectorsourceSelector['science'].signalToNoise.fluxField = fluxField
537  self.sourceSelectorsourceSelector['science'].signalToNoise.errField = fluxField + "Err"
538  # Do not trust blended sources' aperture fluxes which also depend on seeing
539  self.sourceSelectorsourceSelector['science'].doIsolated = True
540  # Do not trust either flux or centroid measurements with flags,
541  # chosen from the usual QA flags for stars)
542  self.sourceSelectorsourceSelector['science'].doFlags = True
543  badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated',
544  'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag',
545  'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter']
546  self.sourceSelectorsourceSelector['science'].flags.bad = badFlags
547 
548  # Default to Gaia-DR2 (with proper motions) for astrometry and
549  # PS1-DR1 for photometry, with a reasonable initial filterMap.
550  self.astrometryRefObjLoaderastrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414"
551  self.astrometryRefObjLoaderastrometryRefObjLoader.requireProperMotion = True
552  self.astrometryRefObjLoaderastrometryRefObjLoader.anyFilterMapsToThis = 'phot_g_mean'
553  self.photometryRefObjLoaderphotometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110"
554 
555 
556 def writeModel(model, filename, log):
557  """Write model to outfile."""
558  with open(filename, "w") as file:
559  file.write(repr(model))
560  log.info("Wrote %s to file: %s", model, filename)
561 
562 
563 @dataclasses.dataclass
565  """The input data jointcal needs for each detector/visit."""
566  visit: int
567  """The visit identifier of this exposure."""
569  """The catalog derived from this exposure."""
570  visitInfo: lsst.afw.image.VisitInfo
571  """The VisitInfo of this exposure."""
573  """The detector of this exposure."""
574  photoCalib: lsst.afw.image.PhotoCalib
575  """The photometric calibration of this exposure."""
577  """The WCS of this exposure."""
578  bbox: lsst.geom.Box2I
579  """The bounding box of this exposure."""
581  """The filter of this exposure."""
582 
583 
584 class JointcalTask(pipeBase.PipelineTask, pipeBase.CmdLineTask):
585  """Astrometricly and photometricly calibrate across multiple visits of the
586  same field.
587 
588  Parameters
589  ----------
590  butler : `lsst.daf.persistence.Butler`
591  The butler is passed to the refObjLoader constructor in case it is
592  needed. Ignored if the refObjLoader argument provides a loader directly.
593  Used to initialize the astrometry and photometry refObjLoaders.
594  initInputs : `dict`, optional
595  Dictionary used to initialize PipelineTasks (empty for jointcal).
596  """
597 
598  ConfigClass = JointcalConfig
599  RunnerClass = JointcalRunner
600  _DefaultName = "jointcal"
601 
602  def __init__(self, butler=None, initInputs=None, **kwargs):
603  super().__init__(**kwargs)
604  self.makeSubtask("sourceSelector")
605  if self.config.doAstrometry:
606  if initInputs is None:
607  # gen3 middleware does refcat things internally (and will not have a butler here)
608  self.makeSubtask('astrometryRefObjLoader', butler=butler)
609  self.makeSubtask("astrometryReferenceSelector")
610  else:
611  self.astrometryRefObjLoaderastrometryRefObjLoader = None
612  if self.config.doPhotometry:
613  if initInputs is None:
614  # gen3 middleware does refcat things internally (and will not have a butler here)
615  self.makeSubtask('photometryRefObjLoader', butler=butler)
616  self.makeSubtask("photometryReferenceSelector")
617  else:
618  self.photometryRefObjLoaderphotometryRefObjLoader = None
619 
620  # To hold various computed metrics for use by tests
621  self.jobjob = Job.load_metrics_package(subset='jointcal')
622 
623  def runQuantum(self, butlerQC, inputRefs, outputRefs):
624  # We override runQuantum to set up the refObjLoaders and write the
625  # outputs to the correct refs.
626  inputs = butlerQC.get(inputRefs)
627  # We want the tract number for writing debug files
628  tract = butlerQC.quantum.dataId['tract']
629  if self.config.doAstrometry:
630  self.astrometryRefObjLoaderastrometryRefObjLoader = ReferenceObjectLoader(
631  dataIds=[ref.datasetRef.dataId
632  for ref in inputRefs.astrometryRefCat],
633  refCats=inputs.pop('astrometryRefCat'),
634  config=self.config.astrometryRefObjLoader,
635  log=self.log)
636  if self.config.doPhotometry:
637  self.photometryRefObjLoaderphotometryRefObjLoader = ReferenceObjectLoader(
638  dataIds=[ref.datasetRef.dataId
639  for ref in inputRefs.photometryRefCat],
640  refCats=inputs.pop('photometryRefCat'),
641  config=self.config.photometryRefObjLoader,
642  log=self.log)
643  outputs = self.runrun(**inputs, tract=tract)
644  if self.config.doAstrometry:
645  self._put_output_put_output(butlerQC, outputs.outputWcs, outputRefs.outputWcs,
646  inputs['inputCamera'], "setWcs")
647  if self.config.doPhotometry:
648  self._put_output_put_output(butlerQC, outputs.outputPhotoCalib, outputRefs.outputPhotoCalib,
649  inputs['inputCamera'], "setPhotoCalib")
650 
651  def _put_output(self, butlerQC, outputs, outputRefs, camera, setter):
652  """Persist the output datasets to their appropriate datarefs.
653 
654  Parameters
655  ----------
656  butlerQC : `ButlerQuantumContext`
657  A butler which is specialized to operate in the context of a
658  `lsst.daf.butler.Quantum`; This is the input to `runQuantum`.
659  outputs : `dict` [`tuple`, `lsst.afw.geom.SkyWcs`] or
660  `dict` [`tuple, `lsst.afw.image.PhotoCalib`]
661  The fitted objects to persist.
662  outputRefs : `list` [`OutputQuantizedConnection`]
663  The DatasetRefs to persist the data to.
664  camera : `lsst.afw.cameraGeom.Camera`
665  The camera for this instrument, to get detector ids from.
666  setter : `str`
667  The method to call on the ExposureCatalog to set each output.
668  """
670  schema.addField('visit', type='I', doc='Visit number')
671 
672  def new_catalog(visit, size):
673  """Return an catalog ready to be filled with appropriate output."""
674  catalog = lsst.afw.table.ExposureCatalog(schema)
675  catalog.resize(size)
676  catalog['visit'] = visit
677  metadata = lsst.daf.base.PropertyList()
678  metadata.add("COMMENT", "Catalog id is detector id, sorted.")
679  metadata.add("COMMENT", "Only detectors with data have entries.")
680  return catalog
681 
682  # count how many detectors have output for each visit
683  detectors_per_visit = collections.defaultdict(int)
684  for key in outputs:
685  # key is (visit, detector_id), and we only need visit here
686  detectors_per_visit[key[0]] += 1
687 
688  for ref in outputRefs:
689  visit = ref.dataId['visit']
690  catalog = new_catalog(visit, detectors_per_visit[visit])
691  # Iterate over every detector and skip the ones we don't have output for.
692  i = 0
693  for detector in camera:
694  detectorId = detector.getId()
695  key = (ref.dataId['visit'], detectorId)
696  if key not in outputs:
697  # skip detectors we don't have output for
698  self.log.debug("No %s output for detector %s in visit %s",
699  setter[3:], detectorId, visit)
700  continue
701 
702  catalog[i].setId(detectorId)
703  getattr(catalog[i], setter)(outputs[key])
704  i += 1
705 
706  catalog.sort() # ensure that the detectors are in sorted order, for fast lookups
707  butlerQC.put(catalog, ref)
708  self.log.info("Wrote %s detectors to %s", i, ref)
709 
710  def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None):
711  # Docstring inherited.
712 
713  # We take values out of the Parquet table, and put them in "flux_",
714  # and the config.sourceFluxType field is used during that extraction,
715  # so just use "flux" here.
716  sourceFluxField = "flux"
717  jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
718  associations = lsst.jointcal.Associations()
719  self.focalPlaneBBoxfocalPlaneBBox = inputCamera.getFpBBox()
720  oldWcsList, bands = self._load_data_load_data(inputSourceTableVisit,
721  inputVisitSummary,
722  associations,
723  jointcalControl,
724  inputCamera)
725 
726  boundingCircle, center, radius, defaultFilter = self._prep_sky_prep_sky(associations, bands)
727  epoch = self._compute_proper_motion_epoch_compute_proper_motion_epoch(associations.getCcdImageList())
728 
729  if self.config.doAstrometry:
730  astrometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
731  name="astrometry",
732  refObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
733  referenceSelector=self.astrometryReferenceSelector,
734  fit_function=self._fit_astrometry_fit_astrometry,
735  tract=tract,
736  epoch=epoch)
737  astrometry_output = self._make_output_make_output(associations.getCcdImageList(),
738  astrometry.model,
739  "makeSkyWcs")
740  else:
741  astrometry_output = None
742 
743  if self.config.doPhotometry:
744  photometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
745  name="photometry",
746  refObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
747  referenceSelector=self.photometryReferenceSelector,
748  fit_function=self._fit_photometry_fit_photometry,
749  tract=tract,
750  epoch=epoch,
751  reject_bad_fluxes=True)
752  photometry_output = self._make_output_make_output(associations.getCcdImageList(),
753  photometry.model,
754  "toPhotoCalib")
755  else:
756  photometry_output = None
757 
758  return pipeBase.Struct(outputWcs=astrometry_output,
759  outputPhotoCalib=photometry_output,
760  job=self.jobjob,
761  astrometryRefObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
762  photometryRefObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader)
763 
764  def _make_schema_table(self):
765  """Return an afw SourceTable to use as a base for creating the
766  SourceCatalog to insert values from the dataFrame into.
767 
768  Returns
769  -------
770  table : `lsst.afw.table.SourceTable`
771  Table with schema and slots to use to make SourceCatalogs.
772  """
774  schema.addField("centroid_x", "D")
775  schema.addField("centroid_y", "D")
776  schema.addField("centroid_xErr", "F")
777  schema.addField("centroid_yErr", "F")
778  schema.addField("shape_xx", "D")
779  schema.addField("shape_yy", "D")
780  schema.addField("shape_xy", "D")
781  schema.addField("flux_instFlux", "D")
782  schema.addField("flux_instFluxErr", "D")
783  table = lsst.afw.table.SourceTable.make(schema)
784  table.defineCentroid("centroid")
785  table.defineShape("shape")
786  return table
787 
788  def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId):
789  """Return an afw SourceCatalog extracted from a visit-level dataframe,
790  limited to just one detector.
791 
792  Parameters
793  ----------
794  table : `lsst.afw.table.SourceTable`
795  Table factory to use to make the SourceCatalog that will be
796  populated with data from ``visitCatalog``.
797  visitCatalog : `pandas.DataFrame`
798  DataFrame to extract a detector catalog from.
799  detectorId : `int`
800  Numeric id of the detector to extract from ``visitCatalog``.
801 
802  Returns
803  -------
804  catalog : `lsst.afw.table.SourceCatalog`
805  Detector-level catalog extracted from ``visitCatalog``.
806  """
807  # map from dataFrame column to afw table column
808  mapping = {'sourceId': 'id',
809  'x': 'centroid_x',
810  'y': 'centroid_y',
811  'xErr': 'centroid_xErr',
812  'yErr': 'centroid_yErr',
813  'Ixx': 'shape_xx',
814  'Iyy': 'shape_yy',
815  'Ixy': 'shape_xy',
816  f'{self.config.sourceFluxType}_instFlux': 'flux_instFlux',
817  f'{self.config.sourceFluxType}_instFluxErr': 'flux_instFluxErr',
818  }
819  # If the DataFrame we're reading was generated by a task running with
820  # Gen2 middleware, the column for the detector will be "ccd" for at
821  # least HSC (who knows what it might be in general!); that would be
822  # true even if the data repo is later converted to Gen3, because that
823  # doesn't touch the files themselves. In Gen3, the column for the
824  # detector will always be "detector".
825  detector_column = "detector" if "detector" in visitCatalog.columns else "ccd"
826  catalog = lsst.afw.table.SourceCatalog(table)
827  matched = visitCatalog[detector_column] == detectorId
828  catalog.resize(sum(matched))
829  view = visitCatalog.loc[matched]
830  for dfCol, afwCol in mapping.items():
831  catalog[afwCol] = view[dfCol]
832 
833  self.log.debug("%d sources selected in visit %d detector %d",
834  len(catalog),
835  view['visit'].iloc[0], # all visits in this catalog are the same, so take the first
836  detectorId)
837  return catalog
838 
839  def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations,
840  jointcalControl, camera):
841  """Read the data that jointcal needs to run. (Gen3 version)
842 
843  Modifies ``associations`` in-place with the loaded data.
844 
845  Parameters
846  ----------
847  inputSourceTableVisit : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
848  References to visit-level DataFrames to load the catalog data from.
849  inputVisitSummary : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
850  Visit-level exposure summary catalog with metadata.
851  associations : `lsst.jointcal.Associations`
852  Object to add the loaded data to by constructing new CcdImages.
853  jointcalControl : `jointcal.JointcalControl`
854  Control object for C++ associations management.
855  camera : `lsst.afw.cameraGeom.Camera`
856  Camera object for detector geometry.
857 
858  Returns
859  -------
860  oldWcsList: `list` [`lsst.afw.geom.SkyWcs`]
861  The original WCS of the input data, to aid in writing tests.
862  bands : `list` [`str`]
863  The filter bands of each input dataset.
864  """
865  oldWcsList = []
866  filters = []
867  load_cat_prof_file = 'jointcal_load_data.prof' if self.config.detailedProfile else ''
868  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
869  table = self._make_schema_table_make_schema_table() # every detector catalog has the same layout
870  # No guarantee that the input is in the same order of visits, so we have to map one of them.
871  catalogMap = {ref.dataId['visit']: i for i, ref in enumerate(inputSourceTableVisit)}
872  detectorDict = {detector.getId(): detector for detector in camera}
873 
874  for visitSummaryRef in inputVisitSummary:
875  visitSummary = visitSummaryRef.get()
876  visitCatalog = inputSourceTableVisit[catalogMap[visitSummaryRef.dataId['visit']]].get()
877  selected = self.sourceSelector.run(visitCatalog)
878 
879  # Build a CcdImage for each detector in this visit.
880  detectors = {id: index for index, id in enumerate(visitSummary['id'])}
881  for id, index in detectors.items():
882  catalog = self._extract_detector_catalog_from_visit_catalog_extract_detector_catalog_from_visit_catalog(table, selected.sourceCat, id)
883  data = self._make_one_input_data_make_one_input_data(visitSummary[index], catalog, detectorDict)
884  result = self._build_ccdImage_build_ccdImage(data, associations, jointcalControl)
885  if result is not None:
886  oldWcsList.append(result.wcs)
887  # A visit has only one band, so we can just use the first.
888  filters.append(data.filter)
889  if len(filters) == 0:
890  raise RuntimeError("No data to process: did source selector remove all sources?")
891  filters = collections.Counter(filters)
892 
893  return oldWcsList, filters
894 
895  def _make_one_input_data(self, visitRecord, catalog, detectorDict):
896  """Return a data structure for this detector+visit."""
897  return JointcalInputData(visit=visitRecord['visit'],
898  catalog=catalog,
899  visitInfo=visitRecord.getVisitInfo(),
900  detector=detectorDict[visitRecord.getId()],
901  photoCalib=visitRecord.getPhotoCalib(),
902  wcs=visitRecord.getWcs(),
903  bbox=visitRecord.getBBox(),
904  # ExposureRecord doesn't have a FilterLabel yet,
905  # so we have to make one.
906  filter=lsst.afw.image.FilterLabel(band=visitRecord['band'],
907  physical=visitRecord['physical_filter']))
908 
909  # We don't currently need to persist the metadata.
910  # If we do in the future, we will have to add appropriate dataset templates
911  # to each obs package (the metadata template should look like `jointcal_wcs`).
912  def _getMetadataName(self):
913  return None
914 
915  @classmethod
916  def _makeArgumentParser(cls):
917  """Create an argument parser"""
918  parser = pipeBase.ArgumentParser(name=cls._DefaultName_DefaultName)
919  parser.add_id_argument("--id", "calexp", help="data ID, e.g. --id visit=6789 ccd=0..9",
920  ContainerClass=PerTractCcdDataIdContainer)
921  return parser
922 
923  def _build_ccdImage(self, data, associations, jointcalControl):
924  """
925  Extract the necessary things from this catalog+metadata to add a new
926  ccdImage.
927 
928  Parameters
929  ----------
930  data : `JointcalInputData`
931  The loaded input data.
932  associations : `lsst.jointcal.Associations`
933  Object to add the info to, to construct a new CcdImage
934  jointcalControl : `jointcal.JointcalControl`
935  Control object for associations management
936 
937  Returns
938  ------
939  namedtuple or `None`
940  ``wcs``
941  The TAN WCS of this image, read from the calexp
942  (`lsst.afw.geom.SkyWcs`).
943  ``key``
944  A key to identify this dataRef by its visit and ccd ids
945  (`namedtuple`).
946  `None`
947  if there are no sources in the loaded catalog.
948  """
949  if len(data.catalog) == 0:
950  self.log.warn("No sources selected in visit %s ccd %s", data.visit, data.detector.getId())
951  return None
952 
953  associations.createCcdImage(data.catalog,
954  data.wcs,
955  data.visitInfo,
956  data.bbox,
957  data.filter.physicalLabel,
958  data.photoCalib,
959  data.detector,
960  data.visit,
961  data.detector.getId(),
962  jointcalControl)
963 
964  Result = collections.namedtuple('Result_from_build_CcdImage', ('wcs', 'key'))
965  Key = collections.namedtuple('Key', ('visit', 'ccd'))
966  return Result(data.wcs, Key(data.visit, data.detector.getId()))
967 
968  def _readDataId(self, butler, dataId):
969  """Read all of the data for one dataId from the butler. (gen2 version)"""
970  # Not all instruments have `visit` in their dataIds.
971  if "visit" in dataId.keys():
972  visit = dataId["visit"]
973  else:
974  visit = butler.getButler().queryMetadata("calexp", ("visit"), butler.dataId)[0]
975  detector = butler.get('calexp_detector', dataId=dataId)
976 
977  catalog = butler.get('src',
978  flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS,
979  dataId=dataId)
980  goodSrc = self.sourceSelector.run(catalog)
981  self.log.debug("%d sources selected in visit %d detector %d",
982  len(goodSrc.sourceCat),
983  visit,
984  detector.getId())
985  return JointcalInputData(visit=visit,
986  catalog=goodSrc.sourceCat,
987  visitInfo=butler.get('calexp_visitInfo', dataId=dataId),
988  detector=detector,
989  photoCalib=butler.get('calexp_photoCalib', dataId=dataId),
990  wcs=butler.get('calexp_wcs', dataId=dataId),
991  bbox=butler.get('calexp_bbox', dataId=dataId),
992  filter=butler.get('calexp_filterLabel', dataId=dataId))
993 
994  def loadData(self, dataRefs, associations, jointcalControl):
995  """Read the data that jointcal needs to run. (Gen2 version)"""
996  visit_ccd_to_dataRef = {}
997  oldWcsList = []
998  filters = []
999  load_cat_prof_file = 'jointcal_loadData.prof' if self.config.detailedProfile else ''
1000  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
1001  # Need the bounding-box of the focal plane (the same for all visits) for photometry visit models
1002  camera = dataRefs[0].get('camera', immediate=True)
1003  self.focalPlaneBBoxfocalPlaneBBox = camera.getFpBBox()
1004  for dataRef in dataRefs:
1005  data = self._readDataId_readDataId(dataRef.getButler(), dataRef.dataId)
1006  result = self._build_ccdImage_build_ccdImage(data, associations, jointcalControl)
1007  if result is None:
1008  continue
1009  oldWcsList.append(result.wcs)
1010  visit_ccd_to_dataRef[result.key] = dataRef
1011  filters.append(data.filter)
1012  if len(filters) == 0:
1013  raise RuntimeError("No data to process: did source selector remove all sources?")
1014  filters = collections.Counter(filters)
1015 
1016  return oldWcsList, filters, visit_ccd_to_dataRef
1017 
1018  def _getDebugPath(self, filename):
1019  """Constructs a path to filename using the configured debug path.
1020  """
1021  return os.path.join(self.config.debugOutputPath, filename)
1022 
1023  def _prep_sky(self, associations, filters):
1024  """Prepare on-sky and other data that must be computed after data has
1025  been read.
1026  """
1027  associations.computeCommonTangentPoint()
1028 
1029  boundingCircle = associations.computeBoundingCircle()
1030  center = lsst.geom.SpherePoint(boundingCircle.getCenter())
1031  radius = lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
1032 
1033  self.log.info(f"Data has center={center} with radius={radius.asDegrees()} degrees.")
1034 
1035  # Determine a default filter band associated with the catalog. See DM-9093
1036  defaultFilter = filters.most_common(1)[0][0]
1037  self.log.debug("Using '%s' filter for reference flux", defaultFilter.physicalLabel)
1038 
1039  return boundingCircle, center, radius, defaultFilter
1040 
1041  @pipeBase.timeMethod
1042  def runDataRef(self, dataRefs):
1043  """
1044  Jointly calibrate the astrometry and photometry across a set of images.
1045 
1046  NOTE: this is for gen2 middleware only.
1047 
1048  Parameters
1049  ----------
1050  dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef`
1051  List of data references to the exposures to be fit.
1052 
1053  Returns
1054  -------
1055  result : `lsst.pipe.base.Struct`
1056  Struct of metadata from the fit, containing:
1057 
1058  ``dataRefs``
1059  The provided data references that were fit (with updated WCSs)
1060  ``oldWcsList``
1061  The original WCS from each dataRef
1062  ``metrics``
1063  Dictionary of internally-computed metrics for testing/validation.
1064  """
1065  if len(dataRefs) == 0:
1066  raise ValueError('Need a non-empty list of data references!')
1067 
1068  exitStatus = 0 # exit status for shell
1069 
1070  sourceFluxField = "slot_%sFlux" % (self.config.sourceFluxType,)
1071  jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
1072  associations = lsst.jointcal.Associations()
1073 
1074  oldWcsList, filters, visit_ccd_to_dataRef = self.loadDataloadData(dataRefs,
1075  associations,
1076  jointcalControl)
1077 
1078  boundingCircle, center, radius, defaultFilter = self._prep_sky_prep_sky(associations, filters)
1079  epoch = self._compute_proper_motion_epoch_compute_proper_motion_epoch(associations.getCcdImageList())
1080 
1081  tract = dataRefs[0].dataId['tract']
1082 
1083  if self.config.doAstrometry:
1084  astrometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
1085  name="astrometry",
1086  refObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
1087  referenceSelector=self.astrometryReferenceSelector,
1088  fit_function=self._fit_astrometry_fit_astrometry,
1089  tract=tract,
1090  epoch=epoch)
1091  self._write_astrometry_results_write_astrometry_results(associations, astrometry.model, visit_ccd_to_dataRef)
1092  else:
1093  astrometry = Astrometry(None, None, None)
1094 
1095  if self.config.doPhotometry:
1096  photometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
1097  name="photometry",
1098  refObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
1099  referenceSelector=self.photometryReferenceSelector,
1100  fit_function=self._fit_photometry_fit_photometry,
1101  tract=tract,
1102  epoch=epoch,
1103  reject_bad_fluxes=True)
1104  self._write_photometry_results_write_photometry_results(associations, photometry.model, visit_ccd_to_dataRef)
1105  else:
1106  photometry = Photometry(None, None)
1107 
1108  return pipeBase.Struct(dataRefs=dataRefs,
1109  oldWcsList=oldWcsList,
1110  job=self.jobjob,
1111  astrometryRefObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
1112  photometryRefObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
1113  defaultFilter=defaultFilter,
1114  epoch=epoch,
1115  exitStatus=exitStatus)
1116 
1117  def _get_refcat_coordinate_error_override(self, refCat, name):
1118  """Check whether we should override the refcat coordinate errors, and
1119  return the overridden error if necessary.
1120 
1121  Parameters
1122  ----------
1123  refCat : `lsst.afw.table.SimpleCatalog`
1124  The reference catalog to check for a ``coord_raErr`` field.
1125  name : `str`
1126  Whether we are doing "astrometry" or "photometry".
1127 
1128  Returns
1129  -------
1130  refCoordErr : `float`
1131  The refcat coordinate error to use, or NaN if we are not overriding
1132  those fields.
1133 
1134  Raises
1135  ------
1136  lsst.pex.config.FieldValidationError
1137  Raised if the refcat does not contain coordinate errors and
1138  ``config.astrometryReferenceErr`` is not set.
1139  """
1140  # This value doesn't matter for photometry, so just set something to
1141  # keep old refcats from causing problems.
1142  if name.lower() == "photometry":
1143  if 'coord_raErr' not in refCat.schema:
1144  return 100
1145  else:
1146  return float('nan')
1147 
1148  if self.config.astrometryReferenceErr is None and 'coord_raErr' not in refCat.schema:
1149  msg = ("Reference catalog does not contain coordinate errors, "
1150  "and config.astrometryReferenceErr not supplied.")
1151  raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
1152  self.config,
1153  msg)
1154 
1155  if self.config.astrometryReferenceErr is not None and 'coord_raErr' in refCat.schema:
1156  self.log.warn("Overriding reference catalog coordinate errors with %f/coordinate [mas]",
1157  self.config.astrometryReferenceErr)
1158 
1159  if self.config.astrometryReferenceErr is None:
1160  return float('nan')
1161  else:
1162  return self.config.astrometryReferenceErr
1163 
1164  def _compute_proper_motion_epoch(self, ccdImageList):
1165  """Return the proper motion correction epoch of the provided images.
1166 
1167  Parameters
1168  ----------
1169  ccdImageList : `list` [`lsst.jointcal.CcdImage`]
1170  The images to compute the appropriate epoch for.
1171 
1172  Returns
1173  -------
1174  epoch : `astropy.time.Time`
1175  The date to use for proper motion corrections.
1176  """
1177  mjds = [ccdImage.getMjd() for ccdImage in ccdImageList]
1178  return astropy.time.Time(np.mean(mjds), format='mjd', scale="tai")
1179 
1180  def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
1181  tract="", match_cut=3.0,
1182  reject_bad_fluxes=False, *,
1183  name="", refObjLoader=None, referenceSelector=None,
1184  fit_function=None, epoch=None):
1185  """Load reference catalog, perform the fit, and return the result.
1186 
1187  Parameters
1188  ----------
1189  associations : `lsst.jointcal.Associations`
1190  The star/reference star associations to fit.
1191  defaultFilter : `lsst.afw.image.FilterLabel`
1192  filter to load from reference catalog.
1193  center : `lsst.geom.SpherePoint`
1194  ICRS center of field to load from reference catalog.
1195  radius : `lsst.geom.Angle`
1196  On-sky radius to load from reference catalog.
1197  name : `str`
1198  Name of thing being fit: "astrometry" or "photometry".
1199  refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
1200  Reference object loader to use to load a reference catalog.
1201  referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1202  Selector to use to pick objects from the loaded reference catalog.
1203  fit_function : callable
1204  Function to call to perform fit (takes Associations object).
1205  tract : `str`, optional
1206  Name of tract currently being fit.
1207  match_cut : `float`, optional
1208  Radius in arcseconds to find cross-catalog matches to during
1209  associations.associateCatalogs.
1210  reject_bad_fluxes : `bool`, optional
1211  Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
1212  epoch : `astropy.time.Time`, optional
1213  Epoch to which to correct refcat proper motion and parallax,
1214  or `None` to not apply such corrections.
1215 
1216  Returns
1217  -------
1218  result : `Photometry` or `Astrometry`
1219  Result of `fit_function()`
1220  """
1221  self.log.info("====== Now processing %s...", name)
1222  # TODO: this should not print "trying to invert a singular transformation:"
1223  # if it does that, something's not right about the WCS...
1224  associations.associateCatalogs(match_cut)
1225  add_measurement(self.jobjob, 'jointcal.associated_%s_fittedStars' % name,
1226  associations.fittedStarListSize())
1227 
1228  applyColorterms = False if name.lower() == "astrometry" else self.config.applyColorTerms
1229  refCat, fluxField = self._load_reference_catalog_load_reference_catalog(refObjLoader, referenceSelector,
1230  center, radius, defaultFilter,
1231  applyColorterms=applyColorterms,
1232  epoch=epoch)
1233  refCoordErr = self._get_refcat_coordinate_error_override_get_refcat_coordinate_error_override(refCat, name)
1234 
1235  associations.collectRefStars(refCat,
1236  self.config.matchCut*lsst.geom.arcseconds,
1237  fluxField,
1238  refCoordinateErr=refCoordErr,
1239  rejectBadFluxes=reject_bad_fluxes)
1240  add_measurement(self.jobjob, 'jointcal.collected_%s_refStars' % name,
1241  associations.refStarListSize())
1242 
1243  associations.prepareFittedStars(self.config.minMeasurements)
1244 
1245  self._check_star_lists_check_star_lists(associations, name)
1246  add_measurement(self.jobjob, 'jointcal.selected_%s_refStars' % name,
1247  associations.nFittedStarsWithAssociatedRefStar())
1248  add_measurement(self.jobjob, 'jointcal.selected_%s_fittedStars' % name,
1249  associations.fittedStarListSize())
1250  add_measurement(self.jobjob, 'jointcal.selected_%s_ccdImages' % name,
1251  associations.nCcdImagesValidForFit())
1252 
1253  load_cat_prof_file = 'jointcal_fit_%s.prof'%name if self.config.detailedProfile else ''
1254  dataName = "{}_{}".format(tract, defaultFilter.bandLabel)
1255  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
1256  result = fit_function(associations, dataName)
1257  # TODO DM-12446: turn this into a "butler save" somehow.
1258  # Save reference and measurement chi2 contributions for this data
1259  if self.config.writeChi2FilesInitialFinal:
1260  baseName = self._getDebugPath_getDebugPath(f"{name}_final_chi2-{dataName}")
1261  result.fit.saveChi2Contributions(baseName+"{type}")
1262  self.log.info("Wrote chi2 contributions files: %s", baseName)
1263 
1264  return result
1265 
1266  def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel,
1267  applyColorterms=False, epoch=None):
1268  """Load the necessary reference catalog sources, convert fluxes to
1269  correct units, and apply color term corrections if requested.
1270 
1271  Parameters
1272  ----------
1273  refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
1274  The reference catalog loader to use to get the data.
1275  referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1276  Source selector to apply to loaded reference catalog.
1277  center : `lsst.geom.SpherePoint`
1278  The center around which to load sources.
1279  radius : `lsst.geom.Angle`
1280  The radius around ``center`` to load sources in.
1281  filterLabel : `lsst.afw.image.FilterLabel`
1282  The camera filter to load fluxes for.
1283  applyColorterms : `bool`
1284  Apply colorterm corrections to the refcat for ``filterName``?
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  refCat : `lsst.afw.table.SimpleCatalog`
1292  The loaded reference catalog.
1293  fluxField : `str`
1294  The name of the reference catalog flux field appropriate for ``filterName``.
1295  """
1296  skyCircle = refObjLoader.loadSkyCircle(center,
1297  radius,
1298  filterLabel.bandLabel,
1299  epoch=epoch)
1300 
1301  selected = referenceSelector.run(skyCircle.refCat)
1302  # Need memory contiguity to get reference filters as a vector.
1303  if not selected.sourceCat.isContiguous():
1304  refCat = selected.sourceCat.copy(deep=True)
1305  else:
1306  refCat = selected.sourceCat
1307 
1308  if applyColorterms:
1309  refCatName = refObjLoader.ref_dataset_name
1310  self.log.info("Applying color terms for physical filter=%r reference catalog=%s",
1311  filterLabel.physicalLabel, refCatName)
1312  colorterm = self.config.colorterms.getColorterm(filterLabel.physicalLabel,
1313  refCatName,
1314  doRaise=True)
1315 
1316  refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat)
1317  refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
1318  # TODO: I didn't want to use this, but I'll deal with it in DM-16903
1319  refCat[skyCircle.fluxField+'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
1320 
1321  return refCat, skyCircle.fluxField
1322 
1323  def _check_star_lists(self, associations, name):
1324  # TODO: these should be len(blah), but we need this properly wrapped first.
1325  if associations.nCcdImagesValidForFit() == 0:
1326  raise RuntimeError('No images in the ccdImageList!')
1327  if associations.fittedStarListSize() == 0:
1328  raise RuntimeError('No stars in the {} fittedStarList!'.format(name))
1329  if associations.refStarListSize() == 0:
1330  raise RuntimeError('No stars in the {} reference star list!'.format(name))
1331 
1332  def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
1333  """Compute chi2, log it, validate the model, and return chi2.
1334 
1335  Parameters
1336  ----------
1337  associations : `lsst.jointcal.Associations`
1338  The star/reference star associations to fit.
1339  fit : `lsst.jointcal.FitterBase`
1340  The fitter to use for minimization.
1341  model : `lsst.jointcal.Model`
1342  The model being fit.
1343  chi2Label : `str`
1344  Label to describe the chi2 (e.g. "Initialized", "Final").
1345  writeChi2Name : `str`, optional
1346  Filename prefix to write the chi2 contributions to.
1347  Do not supply an extension: an appropriate one will be added.
1348 
1349  Returns
1350  -------
1351  chi2: `lsst.jointcal.Chi2Accumulator`
1352  The chi2 object for the current fitter and model.
1353 
1354  Raises
1355  ------
1356  FloatingPointError
1357  Raised if chi2 is infinite or NaN.
1358  ValueError
1359  Raised if the model is not valid.
1360  """
1361  if writeChi2Name is not None:
1362  fullpath = self._getDebugPath_getDebugPath(writeChi2Name)
1363  fit.saveChi2Contributions(fullpath+"{type}")
1364  self.log.info("Wrote chi2 contributions files: %s", fullpath)
1365 
1366  chi2 = fit.computeChi2()
1367  self.log.info("%s %s", chi2Label, chi2)
1368  self._check_stars_check_stars(associations)
1369  if not np.isfinite(chi2.chi2):
1370  raise FloatingPointError(f'{chi2Label} chi2 is invalid: {chi2}')
1371  if not model.validate(associations.getCcdImageList(), chi2.ndof):
1372  raise ValueError("Model is not valid: check log messages for warnings.")
1373  return chi2
1374 
1375  def _fit_photometry(self, associations, dataName=None):
1376  """
1377  Fit the photometric data.
1378 
1379  Parameters
1380  ----------
1381  associations : `lsst.jointcal.Associations`
1382  The star/reference star associations to fit.
1383  dataName : `str`
1384  Name of the data being processed (e.g. "1234_HSC-Y"), for
1385  identifying debugging files.
1386 
1387  Returns
1388  -------
1389  fit_result : `namedtuple`
1390  fit : `lsst.jointcal.PhotometryFit`
1391  The photometric fitter used to perform the fit.
1392  model : `lsst.jointcal.PhotometryModel`
1393  The photometric model that was fit.
1394  """
1395  self.log.info("=== Starting photometric fitting...")
1396 
1397  # TODO: should use pex.config.RegistryField here (see DM-9195)
1398  if self.config.photometryModel == "constrainedFlux":
1399  model = lsst.jointcal.ConstrainedFluxModel(associations.getCcdImageList(),
1400  self.focalPlaneBBoxfocalPlaneBBox,
1401  visitOrder=self.config.photometryVisitOrder,
1402  errorPedestal=self.config.photometryErrorPedestal)
1403  # potentially nonlinear problem, so we may need a line search to converge.
1404  doLineSearch = self.config.allowLineSearch
1405  elif self.config.photometryModel == "constrainedMagnitude":
1406  model = lsst.jointcal.ConstrainedMagnitudeModel(associations.getCcdImageList(),
1407  self.focalPlaneBBoxfocalPlaneBBox,
1408  visitOrder=self.config.photometryVisitOrder,
1409  errorPedestal=self.config.photometryErrorPedestal)
1410  # potentially nonlinear problem, so we may need a line search to converge.
1411  doLineSearch = self.config.allowLineSearch
1412  elif self.config.photometryModel == "simpleFlux":
1413  model = lsst.jointcal.SimpleFluxModel(associations.getCcdImageList(),
1414  errorPedestal=self.config.photometryErrorPedestal)
1415  doLineSearch = False # purely linear in model parameters, so no line search needed
1416  elif self.config.photometryModel == "simpleMagnitude":
1417  model = lsst.jointcal.SimpleMagnitudeModel(associations.getCcdImageList(),
1418  errorPedestal=self.config.photometryErrorPedestal)
1419  doLineSearch = False # purely linear in model parameters, so no line search needed
1420 
1421  fit = lsst.jointcal.PhotometryFit(associations, model)
1422  # TODO DM-12446: turn this into a "butler save" somehow.
1423  # Save reference and measurement chi2 contributions for this data
1424  if self.config.writeChi2FilesInitialFinal:
1425  baseName = f"photometry_initial_chi2-{dataName}"
1426  else:
1427  baseName = None
1428  if self.config.writeInitialModel:
1429  fullpath = self._getDebugPath_getDebugPath("initialPhotometryModel.txt")
1430  writeModel(model, fullpath, self.log)
1431  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialized", writeChi2Name=baseName)
1432 
1433  def getChi2Name(whatToFit):
1434  if self.config.writeChi2FilesOuterLoop:
1435  return f"photometry_init-%s_chi2-{dataName}" % whatToFit
1436  else:
1437  return None
1438 
1439  # The constrained model needs the visit transform fit first; the chip
1440  # transform is initialized from the singleFrame PhotoCalib, so it's close.
1441  dumpMatrixFile = self._getDebugPath_getDebugPath("photometry_preinit") if self.config.writeInitMatrix else ""
1442  if self.config.photometryModel.startswith("constrained"):
1443  # no line search: should be purely (or nearly) linear,
1444  # and we want a large step size to initialize with.
1445  fit.minimize("ModelVisit", dumpMatrixFile=dumpMatrixFile)
1446  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize ModelVisit",
1447  writeChi2Name=getChi2Name("ModelVisit"))
1448  dumpMatrixFile = "" # so we don't redo the output on the next step
1449 
1450  fit.minimize("Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
1451  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Model",
1452  writeChi2Name=getChi2Name("Model"))
1453 
1454  fit.minimize("Fluxes") # no line search: always purely linear.
1455  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Fluxes",
1456  writeChi2Name=getChi2Name("Fluxes"))
1457 
1458  fit.minimize("Model Fluxes", doLineSearch=doLineSearch)
1459  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize ModelFluxes",
1460  writeChi2Name=getChi2Name("ModelFluxes"))
1461 
1462  model.freezeErrorTransform()
1463  self.log.debug("Photometry error scales are frozen.")
1464 
1465  chi2 = self._iterate_fit_iterate_fit(associations,
1466  fit,
1467  self.config.maxPhotometrySteps,
1468  "photometry",
1469  "Model Fluxes",
1470  doRankUpdate=self.config.photometryDoRankUpdate,
1471  doLineSearch=doLineSearch,
1472  dataName=dataName)
1473 
1474  add_measurement(self.jobjob, 'jointcal.photometry_final_chi2', chi2.chi2)
1475  add_measurement(self.jobjob, 'jointcal.photometry_final_ndof', chi2.ndof)
1476  return Photometry(fit, model)
1477 
1478  def _fit_astrometry(self, associations, dataName=None):
1479  """
1480  Fit the astrometric data.
1481 
1482  Parameters
1483  ----------
1484  associations : `lsst.jointcal.Associations`
1485  The star/reference star associations to fit.
1486  dataName : `str`
1487  Name of the data being processed (e.g. "1234_HSC-Y"), for
1488  identifying debugging files.
1489 
1490  Returns
1491  -------
1492  fit_result : `namedtuple`
1493  fit : `lsst.jointcal.AstrometryFit`
1494  The astrometric fitter used to perform the fit.
1495  model : `lsst.jointcal.AstrometryModel`
1496  The astrometric model that was fit.
1497  sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1498  The model for the sky to tangent plane projection that was used in the fit.
1499  """
1500 
1501  self.log.info("=== Starting astrometric fitting...")
1502 
1503  associations.deprojectFittedStars()
1504 
1505  # NOTE: need to return sky_to_tan_projection so that it doesn't get garbage collected.
1506  # TODO: could we package sky_to_tan_projection and model together so we don't have to manage
1507  # them so carefully?
1508  sky_to_tan_projection = lsst.jointcal.OneTPPerVisitHandler(associations.getCcdImageList())
1509 
1510  if self.config.astrometryModel == "constrained":
1511  model = lsst.jointcal.ConstrainedAstrometryModel(associations.getCcdImageList(),
1512  sky_to_tan_projection,
1513  chipOrder=self.config.astrometryChipOrder,
1514  visitOrder=self.config.astrometryVisitOrder)
1515  elif self.config.astrometryModel == "simple":
1516  model = lsst.jointcal.SimpleAstrometryModel(associations.getCcdImageList(),
1517  sky_to_tan_projection,
1518  self.config.useInputWcs,
1519  nNotFit=0,
1520  order=self.config.astrometrySimpleOrder)
1521 
1522  fit = lsst.jointcal.AstrometryFit(associations, model, self.config.positionErrorPedestal)
1523  # TODO DM-12446: turn this into a "butler save" somehow.
1524  # Save reference and measurement chi2 contributions for this data
1525  if self.config.writeChi2FilesInitialFinal:
1526  baseName = f"astrometry_initial_chi2-{dataName}"
1527  else:
1528  baseName = None
1529  if self.config.writeInitialModel:
1530  fullpath = self._getDebugPath_getDebugPath("initialAstrometryModel.txt")
1531  writeModel(model, fullpath, self.log)
1532  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initial", writeChi2Name=baseName)
1533 
1534  def getChi2Name(whatToFit):
1535  if self.config.writeChi2FilesOuterLoop:
1536  return f"astrometry_init-%s_chi2-{dataName}" % whatToFit
1537  else:
1538  return None
1539 
1540  dumpMatrixFile = self._getDebugPath_getDebugPath("astrometry_preinit") if self.config.writeInitMatrix else ""
1541  # The constrained model needs the visit transform fit first; the chip
1542  # transform is initialized from the detector's cameraGeom, so it's close.
1543  if self.config.astrometryModel == "constrained":
1544  fit.minimize("DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1545  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize DistortionsVisit",
1546  writeChi2Name=getChi2Name("DistortionsVisit"))
1547  dumpMatrixFile = "" # so we don't redo the output on the next step
1548 
1549  fit.minimize("Distortions", dumpMatrixFile=dumpMatrixFile)
1550  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Distortions",
1551  writeChi2Name=getChi2Name("Distortions"))
1552 
1553  fit.minimize("Positions")
1554  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Positions",
1555  writeChi2Name=getChi2Name("Positions"))
1556 
1557  fit.minimize("Distortions Positions")
1558  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize DistortionsPositions",
1559  writeChi2Name=getChi2Name("DistortionsPositions"))
1560 
1561  chi2 = self._iterate_fit_iterate_fit(associations,
1562  fit,
1563  self.config.maxAstrometrySteps,
1564  "astrometry",
1565  "Distortions Positions",
1566  doRankUpdate=self.config.astrometryDoRankUpdate,
1567  dataName=dataName)
1568 
1569  add_measurement(self.jobjob, 'jointcal.astrometry_final_chi2', chi2.chi2)
1570  add_measurement(self.jobjob, 'jointcal.astrometry_final_ndof', chi2.ndof)
1571 
1572  return Astrometry(fit, model, sky_to_tan_projection)
1573 
1574  def _check_stars(self, associations):
1575  """Count measured and reference stars per ccd and warn/log them."""
1576  for ccdImage in associations.getCcdImageList():
1577  nMeasuredStars, nRefStars = ccdImage.countStars()
1578  self.log.debug("ccdImage %s has %s measured and %s reference stars",
1579  ccdImage.getName(), nMeasuredStars, nRefStars)
1580  if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1581  self.log.warn("ccdImage %s has only %s measuredStars (desired %s)",
1582  ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1583  if nRefStars < self.config.minRefStarsPerCcd:
1584  self.log.warn("ccdImage %s has only %s RefStars (desired %s)",
1585  ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1586 
1587  def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1588  dataName="",
1589  doRankUpdate=True,
1590  doLineSearch=False):
1591  """Run fitter.minimize up to max_steps times, returning the final chi2.
1592 
1593  Parameters
1594  ----------
1595  associations : `lsst.jointcal.Associations`
1596  The star/reference star associations to fit.
1597  fitter : `lsst.jointcal.FitterBase`
1598  The fitter to use for minimization.
1599  max_steps : `int`
1600  Maximum number of steps to run outlier rejection before declaring
1601  convergence failure.
1602  name : {'photometry' or 'astrometry'}
1603  What type of data are we fitting (for logs and debugging files).
1604  whatToFit : `str`
1605  Passed to ``fitter.minimize()`` to define the parameters to fit.
1606  dataName : `str`, optional
1607  Descriptive name for this dataset (e.g. tract and filter),
1608  for debugging.
1609  doRankUpdate : `bool`, optional
1610  Do an Eigen rank update during minimization, or recompute the full
1611  matrix and gradient?
1612  doLineSearch : `bool`, optional
1613  Do a line search for the optimum step during minimization?
1614 
1615  Returns
1616  -------
1617  chi2: `lsst.jointcal.Chi2Statistic`
1618  The final chi2 after the fit converges, or is forced to end.
1619 
1620  Raises
1621  ------
1622  FloatingPointError
1623  Raised if the fitter fails with a non-finite value.
1624  RuntimeError
1625  Raised if the fitter fails for some other reason;
1626  log messages will provide further details.
1627  """
1628  dumpMatrixFile = self._getDebugPath_getDebugPath(f"{name}_postinit") if self.config.writeInitMatrix else ""
1629  oldChi2 = lsst.jointcal.Chi2Statistic()
1630  oldChi2.chi2 = float("inf")
1631  for i in range(max_steps):
1632  if self.config.writeChi2FilesOuterLoop:
1633  writeChi2Name = f"{name}_iterate_{i}_chi2-{dataName}"
1634  else:
1635  writeChi2Name = None
1636  result = fitter.minimize(whatToFit,
1637  self.config.outlierRejectSigma,
1638  doRankUpdate=doRankUpdate,
1639  doLineSearch=doLineSearch,
1640  dumpMatrixFile=dumpMatrixFile)
1641  dumpMatrixFile = "" # clear it so we don't write the matrix again.
1642  chi2 = self._logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(),
1643  f"Fit iteration {i}", writeChi2Name=writeChi2Name)
1644 
1645  if result == MinimizeResult.Converged:
1646  if doRankUpdate:
1647  self.log.debug("fit has converged - no more outliers - redo minimization "
1648  "one more time in case we have lost accuracy in rank update.")
1649  # Redo minimization one more time in case we have lost accuracy in rank update
1650  result = fitter.minimize(whatToFit, self.config.outlierRejectSigma)
1651  chi2 = self._logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(), "Fit completed")
1652 
1653  # log a message for a large final chi2, TODO: DM-15247 for something better
1654  if chi2.chi2/chi2.ndof >= 4.0:
1655  self.log.error("Potentially bad fit: High chi-squared/ndof.")
1656 
1657  break
1658  elif result == MinimizeResult.Chi2Increased:
1659  self.log.warn("Still some outliers remaining but chi2 increased - retry")
1660  # Check whether the increase was large enough to cause trouble.
1661  chi2Ratio = chi2.chi2 / oldChi2.chi2
1662  if chi2Ratio > 1.5:
1663  self.log.warn('Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1664  chi2.chi2, oldChi2.chi2, chi2Ratio)
1665  # Based on a variety of HSC jointcal logs (see DM-25779), it
1666  # appears that chi2 increases more than a factor of ~2 always
1667  # result in the fit diverging rapidly and ending at chi2 > 1e10.
1668  # Using 10 as the "failure" threshold gives some room between
1669  # leaving a warning and bailing early.
1670  if chi2Ratio > 10:
1671  msg = ("Large chi2 increase between steps: fit likely cannot converge."
1672  " Try setting one or more of the `writeChi2*` config fields and looking"
1673  " at how individual star chi2-values evolve during the fit.")
1674  raise RuntimeError(msg)
1675  oldChi2 = chi2
1676  elif result == MinimizeResult.NonFinite:
1677  filename = self._getDebugPath_getDebugPath("{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1678  # TODO DM-12446: turn this into a "butler save" somehow.
1679  fitter.saveChi2Contributions(filename+"{type}")
1680  msg = "Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1681  raise FloatingPointError(msg.format(filename))
1682  elif result == MinimizeResult.Failed:
1683  raise RuntimeError("Chi2 minimization failure, cannot complete fit.")
1684  else:
1685  raise RuntimeError("Unxepected return code from minimize().")
1686  else:
1687  self.log.error("%s failed to converge after %d steps"%(name, max_steps))
1688 
1689  return chi2
1690 
1691  def _make_output(self, ccdImageList, model, func):
1692  """Return the internal jointcal models converted to the afw
1693  structures that will be saved to disk.
1694 
1695  Parameters
1696  ----------
1697  ccdImageList : `lsst.jointcal.CcdImageList`
1698  The list of CcdImages to get the output for.
1699  model : `lsst.jointcal.AstrometryModel` or `lsst.jointcal.PhotometryModel`
1700  The internal jointcal model to convert for each `lsst.jointcal.CcdImage`.
1701  func : `str`
1702  The name of the function to call on ``model`` to get the converted
1703  structure. Must accept an `lsst.jointcal.CcdImage`.
1704 
1705  Returns
1706  -------
1707  output : `dict` [`tuple`, `lsst.jointcal.AstrometryModel`] or
1708  `dict` [`tuple`, `lsst.jointcal.PhotometryModel`]
1709  The data to be saved, keyed on (visit, detector).
1710  """
1711  output = {}
1712  for ccdImage in ccdImageList:
1713  ccd = ccdImage.ccdId
1714  visit = ccdImage.visit
1715  self.log.debug("%s for visit: %d, ccd: %d", func, visit, ccd)
1716  output[(visit, ccd)] = getattr(model, func)(ccdImage)
1717  return output
1718 
1719  def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1720  """
1721  Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1722 
1723  Parameters
1724  ----------
1725  associations : `lsst.jointcal.Associations`
1726  The star/reference star associations to fit.
1727  model : `lsst.jointcal.AstrometryModel`
1728  The astrometric model that was fit.
1729  visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1730  Dict of ccdImage identifiers to dataRefs that were fit.
1731  """
1732  ccdImageList = associations.getCcdImageList()
1733  output = self._make_output_make_output(ccdImageList, model, "makeSkyWcs")
1734  for key, skyWcs in output.items():
1735  dataRef = visit_ccd_to_dataRef[key]
1736  try:
1737  dataRef.put(skyWcs, 'jointcal_wcs')
1738  except pexExceptions.Exception as e:
1739  self.log.fatal('Failed to write updated Wcs: %s', str(e))
1740  raise e
1741 
1742  def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1743  """
1744  Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1745 
1746  Parameters
1747  ----------
1748  associations : `lsst.jointcal.Associations`
1749  The star/reference star associations to fit.
1750  model : `lsst.jointcal.PhotometryModel`
1751  The photoometric model that was fit.
1752  visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1753  Dict of ccdImage identifiers to dataRefs that were fit.
1754  """
1755 
1756  ccdImageList = associations.getCcdImageList()
1757  output = self._make_output_make_output(ccdImageList, model, "toPhotoCalib")
1758  for key, photoCalib in output.items():
1759  dataRef = visit_ccd_to_dataRef[key]
1760  try:
1761  dataRef.put(photoCalib, 'jointcal_photoCalib')
1762  except pexExceptions.Exception as e:
1763  self.log.fatal('Failed to write updated PhotoCalib: %s', str(e))
1764  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
Simple structure to accumulate chi2 and ndof.
Definition: Chi2.h:52
A multi-component model, fitting mappings for sensors and visits simultaneously.
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:840
def _compute_proper_motion_epoch(self, ccdImageList)
Definition: jointcal.py:1164
def _get_refcat_coordinate_error_override(self, refCat, name)
Definition: jointcal.py:1117
def _getDebugPath(self, filename)
Definition: jointcal.py:1018
def _check_star_lists(self, associations, name)
Definition: jointcal.py:1323
def _make_one_input_data(self, visitRecord, catalog, detectorDict)
Definition: jointcal.py:895
def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1719
def _check_stars(self, associations)
Definition: jointcal.py:1574
def _prep_sky(self, associations, filters)
Definition: jointcal.py:1023
def _readDataId(self, butler, dataId)
Definition: jointcal.py:968
def runQuantum(self, butlerQC, inputRefs, outputRefs)
Definition: jointcal.py:623
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:1184
def _build_ccdImage(self, data, associations, jointcalControl)
Definition: jointcal.py:923
def runDataRef(self, dataRefs)
Definition: jointcal.py:1042
def __init__(self, butler=None, initInputs=None, **kwargs)
Definition: jointcal.py:602
def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId)
Definition: jointcal.py:788
def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", doRankUpdate=True, doLineSearch=False)
Definition: jointcal.py:1590
def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None)
Definition: jointcal.py:710
def _fit_photometry(self, associations, dataName=None)
Definition: jointcal.py:1375
def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1742
def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None)
Definition: jointcal.py:1332
def _make_output(self, ccdImageList, model, func)
Definition: jointcal.py:1691
def _put_output(self, butlerQC, outputs, outputRefs, camera, setter)
Definition: jointcal.py:651
def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel, applyColorterms=False, epoch=None)
Definition: jointcal.py:1267
def _fit_astrometry(self, associations, dataName=None)
Definition: jointcal.py:1478
def loadData(self, dataRefs, associations, jointcalControl)
Definition: jointcal.py:994
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:556