lsst.jointcal  21.0.0-8-ga5967ee+f8a261e37b
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
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 from lsst.meas.algorithms.sourceSelector import sourceSelectorRegistry
45 
46 from .dataIds import PerTractCcdDataIdContainer
47 
48 import lsst.jointcal
49 from lsst.jointcal import MinimizeResult
50 
51 __all__ = ["JointcalConfig", "JointcalRunner", "JointcalTask"]
52 
53 Photometry = collections.namedtuple('Photometry', ('fit', 'model'))
54 Astrometry = collections.namedtuple('Astrometry', ('fit', 'model', 'sky_to_tan_projection'))
55 
56 
57 # TODO: move this to MeasurementSet in lsst.verify per DM-12655.
58 def add_measurement(job, name, value):
59  meas = Measurement(job.metrics[name], value)
60  job.measurements.insert(meas)
61 
62 
63 class JointcalRunner(pipeBase.ButlerInitializedTaskRunner):
64  """Subclass of TaskRunner for jointcalTask
65 
66  jointcalTask.runDataRef() takes a number of arguments, one of which is a list of dataRefs
67  extracted from the command line (whereas most CmdLineTasks' runDataRef methods take
68  single dataRef, are are called repeatedly). This class transforms the processed
69  arguments generated by the ArgumentParser into the arguments expected by
70  Jointcal.runDataRef().
71 
72  See pipeBase.TaskRunner for more information.
73  """
74 
75  @staticmethod
76  def getTargetList(parsedCmd, **kwargs):
77  """
78  Return a list of tuples per tract, each containing (dataRefs, kwargs).
79 
80  Jointcal operates on lists of dataRefs simultaneously.
81  """
82  kwargs['profile_jointcal'] = parsedCmd.profile_jointcal
83  kwargs['butler'] = parsedCmd.butler
84 
85  # organize data IDs by tract
86  refListDict = {}
87  for ref in parsedCmd.id.refList:
88  refListDict.setdefault(ref.dataId["tract"], []).append(ref)
89  # we call runDataRef() once with each tract
90  result = [(refListDict[tract], kwargs) for tract in sorted(refListDict.keys())]
91  return result
92 
93  def __call__(self, args):
94  """
95  Parameters
96  ----------
97  args
98  Arguments for Task.runDataRef()
99 
100  Returns
101  -------
102  pipe.base.Struct
103  if self.doReturnResults is False:
104 
105  - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
106 
107  if self.doReturnResults is True:
108 
109  - ``result``: the result of calling jointcal.runDataRef()
110  - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
111  """
112  exitStatus = 0 # exit status for shell
113 
114  # NOTE: cannot call self.makeTask because that assumes args[0] is a single dataRef.
115  dataRefList, kwargs = args
116  butler = kwargs.pop('butler')
117  task = self.TaskClass(config=self.config, log=self.log, butler=butler)
118  result = None
119  try:
120  result = task.runDataRef(dataRefList, **kwargs)
121  exitStatus = result.exitStatus
122  job_path = butler.get('verify_job_filename')
123  result.job.write(job_path[0])
124  except Exception as e: # catch everything, sort it out later.
125  if self.doRaise:
126  raise e
127  else:
128  exitStatus = 1
129  eName = type(e).__name__
130  tract = dataRefList[0].dataId['tract']
131  task.log.fatal("Failed processing tract %s, %s: %s", tract, eName, e)
132 
133  # Put the butler back into kwargs for the other Tasks.
134  kwargs['butler'] = butler
135  if self.doReturnResults:
136  return pipeBase.Struct(result=result, exitStatus=exitStatus)
137  else:
138  return pipeBase.Struct(exitStatus=exitStatus)
139 
140 
141 class JointcalConfig(pexConfig.Config):
142  """Configuration for JointcalTask"""
143 
144  doAstrometry = pexConfig.Field(
145  doc="Fit astrometry and write the fitted result.",
146  dtype=bool,
147  default=True
148  )
149  doPhotometry = pexConfig.Field(
150  doc="Fit photometry and write the fitted result.",
151  dtype=bool,
152  default=True
153  )
154  coaddName = pexConfig.Field(
155  doc="Type of coadd, typically deep or goodSeeing",
156  dtype=str,
157  default="deep"
158  )
159  positionErrorPedestal = pexConfig.Field(
160  doc="Systematic term to apply to the measured position error (pixels)",
161  dtype=float,
162  default=0.02,
163  )
164  photometryErrorPedestal = pexConfig.Field(
165  doc="Systematic term to apply to the measured error on flux or magnitude as a "
166  "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
167  dtype=float,
168  default=0.0,
169  )
170  # TODO: DM-6885 matchCut should be an geom.Angle
171  matchCut = pexConfig.Field(
172  doc="Matching radius between fitted and reference stars (arcseconds)",
173  dtype=float,
174  default=3.0,
175  )
176  minMeasurements = pexConfig.Field(
177  doc="Minimum number of associated measured stars for a fitted star to be included in the fit",
178  dtype=int,
179  default=2,
180  )
181  minMeasuredStarsPerCcd = pexConfig.Field(
182  doc="Minimum number of measuredStars per ccdImage before printing warnings",
183  dtype=int,
184  default=100,
185  )
186  minRefStarsPerCcd = pexConfig.Field(
187  doc="Minimum number of measuredStars per ccdImage before printing warnings",
188  dtype=int,
189  default=30,
190  )
191  allowLineSearch = pexConfig.Field(
192  doc="Allow a line search during minimization, if it is reasonable for the model"
193  " (models with a significant non-linear component, e.g. constrainedPhotometry).",
194  dtype=bool,
195  default=False
196  )
197  astrometrySimpleOrder = pexConfig.Field(
198  doc="Polynomial order for fitting the simple astrometry model.",
199  dtype=int,
200  default=3,
201  )
202  astrometryChipOrder = pexConfig.Field(
203  doc="Order of the per-chip transform for the constrained astrometry model.",
204  dtype=int,
205  default=1,
206  )
207  astrometryVisitOrder = pexConfig.Field(
208  doc="Order of the per-visit transform for the constrained astrometry model.",
209  dtype=int,
210  default=5,
211  )
212  useInputWcs = pexConfig.Field(
213  doc="Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
214  dtype=bool,
215  default=True,
216  )
217  astrometryModel = pexConfig.ChoiceField(
218  doc="Type of model to fit to astrometry",
219  dtype=str,
220  default="constrained",
221  allowed={"simple": "One polynomial per ccd",
222  "constrained": "One polynomial per ccd, and one polynomial per visit"}
223  )
224  photometryModel = pexConfig.ChoiceField(
225  doc="Type of model to fit to photometry",
226  dtype=str,
227  default="constrainedMagnitude",
228  allowed={"simpleFlux": "One constant zeropoint per ccd and visit, fitting in flux space.",
229  "constrainedFlux": "Constrained zeropoint per ccd, and one polynomial per visit,"
230  " fitting in flux space.",
231  "simpleMagnitude": "One constant zeropoint per ccd and visit,"
232  " fitting in magnitude space.",
233  "constrainedMagnitude": "Constrained zeropoint per ccd, and one polynomial per visit,"
234  " fitting in magnitude space.",
235  }
236  )
237  applyColorTerms = pexConfig.Field(
238  doc="Apply photometric color terms to reference stars?"
239  "Requires that colorterms be set to a ColortermLibrary",
240  dtype=bool,
241  default=False
242  )
243  colorterms = pexConfig.ConfigField(
244  doc="Library of photometric reference catalog name to color term dict.",
245  dtype=ColortermLibrary,
246  )
247  photometryVisitOrder = pexConfig.Field(
248  doc="Order of the per-visit polynomial transform for the constrained photometry model.",
249  dtype=int,
250  default=7,
251  )
252  photometryDoRankUpdate = pexConfig.Field(
253  doc=("Do the rank update step during minimization. "
254  "Skipping this can help deal with models that are too non-linear."),
255  dtype=bool,
256  default=True,
257  )
258  astrometryDoRankUpdate = pexConfig.Field(
259  doc=("Do the rank update step during minimization (should not change the astrometry fit). "
260  "Skipping this can help deal with models that are too non-linear."),
261  dtype=bool,
262  default=True,
263  )
264  outlierRejectSigma = pexConfig.Field(
265  doc="How many sigma to reject outliers at during minimization.",
266  dtype=float,
267  default=5.0,
268  )
269  maxPhotometrySteps = pexConfig.Field(
270  doc="Maximum number of minimize iterations to take when fitting photometry.",
271  dtype=int,
272  default=20,
273  )
274  maxAstrometrySteps = pexConfig.Field(
275  doc="Maximum number of minimize iterations to take when fitting photometry.",
276  dtype=int,
277  default=20,
278  )
279  astrometryRefObjLoader = pexConfig.ConfigurableField(
280  target=LoadIndexedReferenceObjectsTask,
281  doc="Reference object loader for astrometric fit",
282  )
283  photometryRefObjLoader = pexConfig.ConfigurableField(
284  target=LoadIndexedReferenceObjectsTask,
285  doc="Reference object loader for photometric fit",
286  )
287  sourceSelector = sourceSelectorRegistry.makeField(
288  doc="How to select sources for cross-matching",
289  default="astrometry"
290  )
291  astrometryReferenceSelector = pexConfig.ConfigurableField(
292  target=ReferenceSourceSelectorTask,
293  doc="How to down-select the loaded astrometry reference catalog.",
294  )
295  photometryReferenceSelector = pexConfig.ConfigurableField(
296  target=ReferenceSourceSelectorTask,
297  doc="How to down-select the loaded photometry reference catalog.",
298  )
299  astrometryReferenceErr = pexConfig.Field(
300  doc=("Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
301  "If None, then raise an exception if the reference catalog is missing coordinate errors. "
302  "If specified, overrides any existing `coord_*Err` values."),
303  dtype=float,
304  default=None,
305  optional=True
306  )
307  writeInitMatrix = pexConfig.Field(
308  dtype=bool,
309  doc=("Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
310  "The output files will be of the form 'astrometry_preinit-mat.txt', in the current directory. "
311  "Note that these files are the dense versions of the matrix, and so may be very large."),
312  default=False
313  )
314  writeChi2FilesInitialFinal = pexConfig.Field(
315  dtype=bool,
316  doc="Write .csv files containing the contributions to chi2 for the initialization and final fit.",
317  default=False
318  )
319  writeChi2FilesOuterLoop = pexConfig.Field(
320  dtype=bool,
321  doc="Write .csv files containing the contributions to chi2 for the outer fit loop.",
322  default=False
323  )
324  writeInitialModel = pexConfig.Field(
325  dtype=bool,
326  doc=("Write the pre-initialization model to text files, for debugging."
327  " Output is written to `initial[Astro|Photo]metryModel.txt` in the current working directory."),
328  default=False
329  )
330  debugOutputPath = pexConfig.Field(
331  dtype=str,
332  default=".",
333  doc=("Path to write debug output files to. Used by "
334  "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
335  )
336  sourceFluxType = pexConfig.Field(
337  dtype=str,
338  doc="Source flux field to use in source selection and to get fluxes from the catalog.",
339  default='Calib'
340  )
341 
342  def validate(self):
343  super().validate()
344  if self.doPhotometrydoPhotometry and self.applyColorTermsapplyColorTerms and len(self.colortermscolorterms.data) == 0:
345  msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
346  raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
347  if self.doAstrometrydoAstrometry and not self.doPhotometrydoPhotometry and self.applyColorTermsapplyColorTerms:
348  msg = ("Only doing astrometry, but Colorterms are not applied for astrometry;"
349  "applyColorTerms=True will be ignored.")
350  lsst.log.warn(msg)
351 
352  def setDefaults(self):
353  # Use science source selector which can filter on extendedness, SNR, and whether blended
354  self.sourceSelectorsourceSelector.name = 'science'
355  # Use only stars because aperture fluxes of galaxies are biased and depend on seeing
356  self.sourceSelectorsourceSelector['science'].doUnresolved = True
357  # with dependable signal to noise ratio.
358  self.sourceSelectorsourceSelector['science'].doSignalToNoise = True
359  # Min SNR must be > 0 because jointcal cannot handle negative fluxes,
360  # and S/N > 10 to use sources that are not too faint, and thus better measured.
361  self.sourceSelectorsourceSelector['science'].signalToNoise.minimum = 10.
362  # Base SNR on CalibFlux because that is the flux jointcal that fits and must be positive
363  fluxField = f"slot_{self.sourceFluxType}Flux_instFlux"
364  self.sourceSelectorsourceSelector['science'].signalToNoise.fluxField = fluxField
365  self.sourceSelectorsourceSelector['science'].signalToNoise.errField = fluxField + "Err"
366  # Do not trust blended sources' aperture fluxes which also depend on seeing
367  self.sourceSelectorsourceSelector['science'].doIsolated = True
368  # Do not trust either flux or centroid measurements with flags,
369  # chosen from the usual QA flags for stars)
370  self.sourceSelectorsourceSelector['science'].doFlags = True
371  badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated',
372  'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag',
373  'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter']
374  self.sourceSelectorsourceSelector['science'].flags.bad = badFlags
375 
376  # Default to Gaia-DR2 (with proper motions) for astrometry and
377  # PS1-DR1 for photometry, with a reasonable initial filterMap.
378  self.astrometryRefObjLoaderastrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414"
379  self.astrometryRefObjLoaderastrometryRefObjLoader.requireProperMotion = True
380  self.astrometryRefObjLoaderastrometryRefObjLoader.anyFilterMapsToThis = 'phot_g_mean'
381  self.photometryRefObjLoaderphotometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110"
382 
383 
384 def writeModel(model, filename, log):
385  """Write model to outfile."""
386  with open(filename, "w") as file:
387  file.write(repr(model))
388  log.info("Wrote %s to file: %s", model, filename)
389 
390 
391 @dataclasses.dataclass
393  """The input data jointcal needs for each detector/visit."""
394  visit: int
395  """The visit identifier of this exposure."""
397  """The catalog derived from this exposure."""
398  visitInfo: lsst.afw.image.VisitInfo
399  """The VisitInfo of this exposure."""
401  """The detector of this exposure."""
402  photoCalib: lsst.afw.image.PhotoCalib
403  """The photometric calibration of this exposure."""
405  """The WCS of this exposure."""
406  bbox: lsst.geom.Box2I
407  """The bounding box of this exposure."""
409  """The filter of this exposure."""
410 
411 
412 class JointcalTask(pipeBase.PipelineTask, pipeBase.CmdLineTask):
413  """Astrometricly and photometricly calibrate across multiple visits of the
414  same field.
415 
416  Parameters
417  ----------
418  butler : `lsst.daf.persistence.Butler`
419  The butler is passed to the refObjLoader constructor in case it is
420  needed. Ignored if the refObjLoader argument provides a loader directly.
421  Used to initialize the astrometry and photometry refObjLoaders.
422  profile_jointcal : `bool`
423  Set to True to profile different stages of this jointcal run.
424  """
425 
426  ConfigClass = JointcalConfig
427  RunnerClass = JointcalRunner
428  _DefaultName = "jointcal"
429 
430  def __init__(self, butler=None, profile_jointcal=False, **kwargs):
431  super().__init__(**kwargs)
432  self.profile_jointcalprofile_jointcal = profile_jointcal
433  self.makeSubtask("sourceSelector")
434  if self.config.doAstrometry:
435  self.makeSubtask('astrometryRefObjLoader', butler=butler)
436  self.makeSubtask("astrometryReferenceSelector")
437  else:
438  self.astrometryRefObjLoaderastrometryRefObjLoader = None
439  if self.config.doPhotometry:
440  self.makeSubtask('photometryRefObjLoader', butler=butler)
441  self.makeSubtask("photometryReferenceSelector")
442  else:
443  self.photometryRefObjLoaderphotometryRefObjLoader = None
444 
445  # To hold various computed metrics for use by tests
446  self.jobjob = Job.load_metrics_package(subset='jointcal')
447 
448  # We don't currently need to persist the metadata.
449  # If we do in the future, we will have to add appropriate dataset templates
450  # to each obs package (the metadata template should look like `jointcal_wcs`).
451  def _getMetadataName(self):
452  return None
453 
454  @classmethod
455  def _makeArgumentParser(cls):
456  """Create an argument parser"""
457  parser = pipeBase.ArgumentParser(name=cls._DefaultName_DefaultName)
458  parser.add_argument("--profile_jointcal", default=False, action="store_true",
459  help="Profile steps of jointcal separately.")
460  parser.add_id_argument("--id", "calexp", help="data ID, e.g. --id visit=6789 ccd=0..9",
461  ContainerClass=PerTractCcdDataIdContainer)
462  return parser
463 
464  def _build_ccdImage(self, data, associations, jointcalControl):
465  """
466  Extract the necessary things from this dataRef to add a new ccdImage.
467 
468  Parameters
469  ----------
470  data : `JointcalInputData`
471  The loaded input data.
472  associations : `lsst.jointcal.Associations`
473  Object to add the info to, to construct a new CcdImage
474  jointcalControl : `jointcal.JointcalControl`
475  Control object for associations management
476 
477  Returns
478  ------
479  namedtuple
480  ``wcs``
481  The TAN WCS of this image, read from the calexp
482  (`lsst.afw.geom.SkyWcs`).
483  ``key``
484  A key to identify this dataRef by its visit and ccd ids
485  (`namedtuple`).
486  ``band``
487  This calexp's filter band (`str`) (used to e.g. load refcats)
488  """
489  goodSrc = self.sourceSelector.run(data.catalog)
490 
491  if len(goodSrc.sourceCat) == 0:
492  self.log.warn("No sources selected in visit %s ccd %s", data.visit, data.detector.getId())
493  else:
494  self.log.info("%d sources selected in visit %d ccd %d", len(goodSrc.sourceCat),
495  data.visit,
496  data.detector.getId())
497  associations.createCcdImage(goodSrc.sourceCat,
498  data.wcs,
499  data.visitInfo,
500  data.bbox,
501  data.filter.physicalLabel,
502  data.photoCalib,
503  data.detector,
504  data.visit,
505  data.detector.getId(),
506  jointcalControl)
507 
508  Result = collections.namedtuple('Result_from_build_CcdImage', ('wcs', 'key', 'band'))
509  Key = collections.namedtuple('Key', ('visit', 'ccd'))
510  return Result(data.wcs, Key(data.visit, data.detector.getId()), data.filter.bandLabel)
511 
512  def _readDataId(self, butler, dataId):
513  """Read all of the data for one dataId from the butler. (gen2 version)"""
514  # Not all instruments have `visit` in their dataIds.
515  if "visit" in dataId.keys():
516  visit = dataId["visit"]
517  else:
518  visit = butler.getButler().queryMetadata("calexp", ("visit"), butler.dataId)[0]
519 
520  catalog = butler.get('src',
521  flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS,
522  dataId=dataId)
523  return JointcalInputData(visit=visit,
524  catalog=catalog,
525  visitInfo=butler.get('calexp_visitInfo', dataId=dataId),
526  detector=butler.get('calexp_detector', dataId=dataId),
527  photoCalib=butler.get('calexp_photoCalib', dataId=dataId),
528  wcs=butler.get('calexp_wcs', dataId=dataId),
529  bbox=butler.get('calexp_bbox', dataId=dataId),
530  filter=butler.get('calexp_filterLabel', dataId=dataId))
531 
532  def loadData(self, dataRefs, associations, jointcalControl, profile_jointcal=False):
533  """Read the data that jointcal needs to run. (Gen2 version)"""
534  visit_ccd_to_dataRef = {}
535  oldWcsList = []
536  bands = []
537  load_cat_prof_file = 'jointcal_loadData.prof' if profile_jointcal else ''
538  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
539  # Need the bounding-box of the focal plane (the same for all visits) for photometry visit models
540  camera = dataRefs[0].get('camera', immediate=True)
541  self.focalPlaneBBoxfocalPlaneBBox = camera.getFpBBox()
542  for dataRef in dataRefs:
543  data = self._readDataId_readDataId(dataRef.getButler(), dataRef.dataId)
544  result = self._build_ccdImage_build_ccdImage(data, associations, jointcalControl)
545  oldWcsList.append(result.wcs)
546  visit_ccd_to_dataRef[result.key] = dataRef
547  bands.append(result.band)
548  bands = collections.Counter(bands)
549 
550  return oldWcsList, bands, visit_ccd_to_dataRef
551 
552  def _getDebugPath(self, filename):
553  """Constructs a path to filename using the configured debug path.
554  """
555  return os.path.join(self.config.debugOutputPath, filename)
556 
557  def _prep_sky(self, associations, bands):
558  """Prepare on-sky and other data that must be computed after data has
559  been read.
560  """
561  associations.computeCommonTangentPoint()
562 
563  boundingCircle = associations.computeBoundingCircle()
564  center = lsst.geom.SpherePoint(boundingCircle.getCenter())
565  radius = lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
566 
567  self.log.info(f"Data has center={center} with radius={radius.asDegrees()} degrees.")
568 
569  # Determine a default filter band associated with the catalog. See DM-9093
570  defaultBand = bands.most_common(1)[0][0]
571  self.log.debug("Using '%s' filter band for reference flux", defaultBand)
572 
573  return boundingCircle, center, radius, defaultBand
574 
575  @pipeBase.timeMethod
576  def runDataRef(self, dataRefs, profile_jointcal=False):
577  """
578  Jointly calibrate the astrometry and photometry across a set of images.
579 
580  NOTE: this is for gen2 middleware only.
581 
582  Parameters
583  ----------
584  dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef`
585  List of data references to the exposures to be fit.
586  profile_jointcal : `bool`
587  Profile the individual steps of jointcal.
588 
589  Returns
590  -------
591  result : `lsst.pipe.base.Struct`
592  Struct of metadata from the fit, containing:
593 
594  ``dataRefs``
595  The provided data references that were fit (with updated WCSs)
596  ``oldWcsList``
597  The original WCS from each dataRef
598  ``metrics``
599  Dictionary of internally-computed metrics for testing/validation.
600  """
601  if len(dataRefs) == 0:
602  raise ValueError('Need a non-empty list of data references!')
603 
604  exitStatus = 0 # exit status for shell
605 
606  sourceFluxField = "slot_%sFlux" % (self.config.sourceFluxType,)
607  jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
608  associations = lsst.jointcal.Associations()
609 
610  oldWcsList, bands, visit_ccd_to_dataRef = self.loadDataloadData(dataRefs,
611  associations,
612  jointcalControl,
613  profile_jointcal=profile_jointcal)
614 
615  boundingCircle, center, radius, defaultBand = self._prep_sky_prep_sky(associations, bands)
616  epoch = self._compute_proper_motion_epoch_compute_proper_motion_epoch(associations.getCcdImageList())
617 
618  tract = dataRefs[0].dataId['tract']
619 
620  if self.config.doAstrometry:
621  astrometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultBand, center, radius,
622  name="astrometry",
623  refObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
624  referenceSelector=self.astrometryReferenceSelector,
625  fit_function=self._fit_astrometry_fit_astrometry,
626  profile_jointcal=profile_jointcal,
627  tract=tract,
628  epoch=epoch)
629  self._write_astrometry_results_write_astrometry_results(associations, astrometry.model, visit_ccd_to_dataRef)
630  else:
631  astrometry = Astrometry(None, None, None)
632 
633  if self.config.doPhotometry:
634  photometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultBand, center, radius,
635  name="photometry",
636  refObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
637  referenceSelector=self.photometryReferenceSelector,
638  fit_function=self._fit_photometry_fit_photometry,
639  profile_jointcal=profile_jointcal,
640  tract=tract,
641  epoch=epoch,
642  reject_bad_fluxes=True)
643  self._write_photometry_results_write_photometry_results(associations, photometry.model, visit_ccd_to_dataRef)
644  else:
645  photometry = Photometry(None, None)
646 
647  return pipeBase.Struct(dataRefs=dataRefs,
648  oldWcsList=oldWcsList,
649  job=self.jobjob,
650  astrometryRefObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
651  photometryRefObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
652  defaultBand=defaultBand,
653  epoch=epoch,
654  exitStatus=exitStatus)
655 
656  def _get_refcat_coordinate_error_override(self, refCat, name):
657  """Check whether we should override the refcat coordinate errors, and
658  return the overridden error if necessary.
659 
660  Parameters
661  ----------
662  refCat : `lsst.afw.table.SimpleCatalog`
663  The reference catalog to check for a ``coord_raErr`` field.
664  name : `str`
665  Whether we are doing "astrometry" or "photometry".
666 
667  Returns
668  -------
669  refCoordErr : `float`
670  The refcat coordinate error to use, or NaN if we are not overriding
671  those fields.
672 
673  Raises
674  ------
675  lsst.pex.config.FieldValidationError
676  Raised if the refcat does not contain coordinate errors and
677  ``config.astrometryReferenceErr`` is not set.
678  """
679  # This value doesn't matter for photometry, so just set something to
680  # keep old refcats from causing problems.
681  if name.lower() == "photometry":
682  if 'coord_raErr' not in refCat.schema:
683  return 100
684  else:
685  return float('nan')
686 
687  if self.config.astrometryReferenceErr is None and 'coord_raErr' not in refCat.schema:
688  msg = ("Reference catalog does not contain coordinate errors, "
689  "and config.astrometryReferenceErr not supplied.")
690  raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
691  self.config,
692  msg)
693 
694  if self.config.astrometryReferenceErr is not None and 'coord_raErr' in refCat.schema:
695  self.log.warn("Overriding reference catalog coordinate errors with %f/coordinate [mas]",
696  self.config.astrometryReferenceErr)
697 
698  if self.config.astrometryReferenceErr is None:
699  return float('nan')
700  else:
701  return self.config.astrometryReferenceErr
702 
703  def _compute_proper_motion_epoch(self, ccdImageList):
704  """Return the proper motion correction epoch of the provided images.
705 
706  Parameters
707  ----------
708  ccdImageList : `list` [`lsst.jointcal.CcdImage`]
709  The images to compute the appropriate epoch for.
710 
711  Returns
712  -------
713  epoch : `astropy.time.Time`
714  The date to use for proper motion corrections.
715  """
716  mjds = [ccdImage.getMjd() for ccdImage in ccdImageList]
717  return astropy.time.Time(np.mean(mjds), format='mjd', scale="tai")
718 
719  def _do_load_refcat_and_fit(self, associations, defaultBand, center, radius,
720  tract="", profile_jointcal=False, match_cut=3.0,
721  reject_bad_fluxes=False, *,
722  name="", refObjLoader=None, referenceSelector=None,
723  fit_function=None, epoch=None):
724  """Load reference catalog, perform the fit, and return the result.
725 
726  Parameters
727  ----------
728  associations : `lsst.jointcal.Associations`
729  The star/reference star associations to fit.
730  defaultBand : `str`
731  filter to load from reference catalog.
732  center : `lsst.geom.SpherePoint`
733  ICRS center of field to load from reference catalog.
734  radius : `lsst.geom.Angle`
735  On-sky radius to load from reference catalog.
736  name : `str`
737  Name of thing being fit: "astrometry" or "photometry".
738  refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
739  Reference object loader to use to load a reference catalog.
740  referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
741  Selector to use to pick objects from the loaded reference catalog.
742  fit_function : callable
743  Function to call to perform fit (takes Associations object).
744  tract : `str`, optional
745  Name of tract currently being fit.
746  profile_jointcal : `bool`, optional
747  Separately profile the fitting step.
748  match_cut : `float`, optional
749  Radius in arcseconds to find cross-catalog matches to during
750  associations.associateCatalogs.
751  reject_bad_fluxes : `bool`, optional
752  Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
753  epoch : `astropy.time.Time`, optional
754  Epoch to which to correct refcat proper motion and parallax,
755  or `None` to not apply such corrections.
756 
757  Returns
758  -------
759  result : `Photometry` or `Astrometry`
760  Result of `fit_function()`
761  """
762  self.log.info("====== Now processing %s...", name)
763  # TODO: this should not print "trying to invert a singular transformation:"
764  # if it does that, something's not right about the WCS...
765  associations.associateCatalogs(match_cut)
766  add_measurement(self.jobjob, 'jointcal.associated_%s_fittedStars' % name,
767  associations.fittedStarListSize())
768 
769  applyColorterms = False if name.lower() == "astrometry" else self.config.applyColorTerms
770  refCat, fluxField = self._load_reference_catalog_load_reference_catalog(refObjLoader, referenceSelector,
771  center, radius, defaultBand,
772  applyColorterms=applyColorterms,
773  epoch=epoch)
774  refCoordErr = self._get_refcat_coordinate_error_override_get_refcat_coordinate_error_override(refCat, name)
775 
776  associations.collectRefStars(refCat,
777  self.config.matchCut*lsst.geom.arcseconds,
778  fluxField,
779  refCoordinateErr=refCoordErr,
780  rejectBadFluxes=reject_bad_fluxes)
781  add_measurement(self.jobjob, 'jointcal.collected_%s_refStars' % name,
782  associations.refStarListSize())
783 
784  associations.prepareFittedStars(self.config.minMeasurements)
785 
786  self._check_star_lists_check_star_lists(associations, name)
787  add_measurement(self.jobjob, 'jointcal.selected_%s_refStars' % name,
788  associations.nFittedStarsWithAssociatedRefStar())
789  add_measurement(self.jobjob, 'jointcal.selected_%s_fittedStars' % name,
790  associations.fittedStarListSize())
791  add_measurement(self.jobjob, 'jointcal.selected_%s_ccdImages' % name,
792  associations.nCcdImagesValidForFit())
793 
794  load_cat_prof_file = 'jointcal_fit_%s.prof'%name if profile_jointcal else ''
795  dataName = "{}_{}".format(tract, defaultBand)
796  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
797  result = fit_function(associations, dataName)
798  # TODO DM-12446: turn this into a "butler save" somehow.
799  # Save reference and measurement chi2 contributions for this data
800  if self.config.writeChi2FilesInitialFinal:
801  baseName = self._getDebugPath_getDebugPath(f"{name}_final_chi2-{dataName}")
802  result.fit.saveChi2Contributions(baseName+"{type}")
803  self.log.info("Wrote chi2 contributions files: %s", baseName)
804 
805  return result
806 
807  def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName,
808  applyColorterms=False, epoch=None):
809  """Load the necessary reference catalog sources, convert fluxes to
810  correct units, and apply color term corrections if requested.
811 
812  Parameters
813  ----------
814  refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
815  The reference catalog loader to use to get the data.
816  referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
817  Source selector to apply to loaded reference catalog.
818  center : `lsst.geom.SpherePoint`
819  The center around which to load sources.
820  radius : `lsst.geom.Angle`
821  The radius around ``center`` to load sources in.
822  filterName : `str`
823  The name of the camera filter to load fluxes for.
824  applyColorterms : `bool`
825  Apply colorterm corrections to the refcat for ``filterName``?
826  epoch : `astropy.time.Time`, optional
827  Epoch to which to correct refcat proper motion and parallax,
828  or `None` to not apply such corrections.
829 
830  Returns
831  -------
832  refCat : `lsst.afw.table.SimpleCatalog`
833  The loaded reference catalog.
834  fluxField : `str`
835  The name of the reference catalog flux field appropriate for ``filterName``.
836  """
837  skyCircle = refObjLoader.loadSkyCircle(center,
838  radius,
839  filterName,
840  epoch=epoch)
841 
842  selected = referenceSelector.run(skyCircle.refCat)
843  # Need memory contiguity to get reference filters as a vector.
844  if not selected.sourceCat.isContiguous():
845  refCat = selected.sourceCat.copy(deep=True)
846  else:
847  refCat = selected.sourceCat
848 
849  if applyColorterms:
850  refCatName = refObjLoader.ref_dataset_name
851  self.log.info("Applying color terms for filterName=%r reference catalog=%s",
852  filterName, refCatName)
853  colorterm = self.config.colorterms.getColorterm(
854  filterName=filterName, photoCatName=refCatName, doRaise=True)
855 
856  refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat, filterName)
857  refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
858  # TODO: I didn't want to use this, but I'll deal with it in DM-16903
859  refCat[skyCircle.fluxField+'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
860 
861  return refCat, skyCircle.fluxField
862 
863  def _check_star_lists(self, associations, name):
864  # TODO: these should be len(blah), but we need this properly wrapped first.
865  if associations.nCcdImagesValidForFit() == 0:
866  raise RuntimeError('No images in the ccdImageList!')
867  if associations.fittedStarListSize() == 0:
868  raise RuntimeError('No stars in the {} fittedStarList!'.format(name))
869  if associations.refStarListSize() == 0:
870  raise RuntimeError('No stars in the {} reference star list!'.format(name))
871 
872  def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
873  """Compute chi2, log it, validate the model, and return chi2.
874 
875  Parameters
876  ----------
877  associations : `lsst.jointcal.Associations`
878  The star/reference star associations to fit.
879  fit : `lsst.jointcal.FitterBase`
880  The fitter to use for minimization.
881  model : `lsst.jointcal.Model`
882  The model being fit.
883  chi2Label : `str`
884  Label to describe the chi2 (e.g. "Initialized", "Final").
885  writeChi2Name : `str`, optional
886  Filename prefix to write the chi2 contributions to.
887  Do not supply an extension: an appropriate one will be added.
888 
889  Returns
890  -------
891  chi2: `lsst.jointcal.Chi2Accumulator`
892  The chi2 object for the current fitter and model.
893 
894  Raises
895  ------
896  FloatingPointError
897  Raised if chi2 is infinite or NaN.
898  ValueError
899  Raised if the model is not valid.
900  """
901  if writeChi2Name is not None:
902  fullpath = self._getDebugPath_getDebugPath(writeChi2Name)
903  fit.saveChi2Contributions(fullpath+"{type}")
904  self.log.info("Wrote chi2 contributions files: %s", fullpath)
905 
906  chi2 = fit.computeChi2()
907  self.log.info("%s %s", chi2Label, chi2)
908  self._check_stars_check_stars(associations)
909  if not np.isfinite(chi2.chi2):
910  raise FloatingPointError(f'{chi2Label} chi2 is invalid: {chi2}')
911  if not model.validate(associations.getCcdImageList(), chi2.ndof):
912  raise ValueError("Model is not valid: check log messages for warnings.")
913  return chi2
914 
915  def _fit_photometry(self, associations, dataName=None):
916  """
917  Fit the photometric data.
918 
919  Parameters
920  ----------
921  associations : `lsst.jointcal.Associations`
922  The star/reference star associations to fit.
923  dataName : `str`
924  Name of the data being processed (e.g. "1234_HSC-Y"), for
925  identifying debugging files.
926 
927  Returns
928  -------
929  fit_result : `namedtuple`
930  fit : `lsst.jointcal.PhotometryFit`
931  The photometric fitter used to perform the fit.
932  model : `lsst.jointcal.PhotometryModel`
933  The photometric model that was fit.
934  """
935  self.log.info("=== Starting photometric fitting...")
936 
937  # TODO: should use pex.config.RegistryField here (see DM-9195)
938  if self.config.photometryModel == "constrainedFlux":
939  model = lsst.jointcal.ConstrainedFluxModel(associations.getCcdImageList(),
940  self.focalPlaneBBoxfocalPlaneBBox,
941  visitOrder=self.config.photometryVisitOrder,
942  errorPedestal=self.config.photometryErrorPedestal)
943  # potentially nonlinear problem, so we may need a line search to converge.
944  doLineSearch = self.config.allowLineSearch
945  elif self.config.photometryModel == "constrainedMagnitude":
946  model = lsst.jointcal.ConstrainedMagnitudeModel(associations.getCcdImageList(),
947  self.focalPlaneBBoxfocalPlaneBBox,
948  visitOrder=self.config.photometryVisitOrder,
949  errorPedestal=self.config.photometryErrorPedestal)
950  # potentially nonlinear problem, so we may need a line search to converge.
951  doLineSearch = self.config.allowLineSearch
952  elif self.config.photometryModel == "simpleFlux":
953  model = lsst.jointcal.SimpleFluxModel(associations.getCcdImageList(),
954  errorPedestal=self.config.photometryErrorPedestal)
955  doLineSearch = False # purely linear in model parameters, so no line search needed
956  elif self.config.photometryModel == "simpleMagnitude":
957  model = lsst.jointcal.SimpleMagnitudeModel(associations.getCcdImageList(),
958  errorPedestal=self.config.photometryErrorPedestal)
959  doLineSearch = False # purely linear in model parameters, so no line search needed
960 
961  fit = lsst.jointcal.PhotometryFit(associations, model)
962  # TODO DM-12446: turn this into a "butler save" somehow.
963  # Save reference and measurement chi2 contributions for this data
964  if self.config.writeChi2FilesInitialFinal:
965  baseName = f"photometry_initial_chi2-{dataName}"
966  else:
967  baseName = None
968  if self.config.writeInitialModel:
969  fullpath = self._getDebugPath_getDebugPath("initialPhotometryModel.txt")
970  writeModel(model, fullpath, self.log)
971  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialized", writeChi2Name=baseName)
972 
973  def getChi2Name(whatToFit):
974  if self.config.writeChi2FilesOuterLoop:
975  return f"photometry_init-%s_chi2-{dataName}" % whatToFit
976  else:
977  return None
978 
979  # The constrained model needs the visit transform fit first; the chip
980  # transform is initialized from the singleFrame PhotoCalib, so it's close.
981  dumpMatrixFile = self._getDebugPath_getDebugPath("photometry_preinit") if self.config.writeInitMatrix else ""
982  if self.config.photometryModel.startswith("constrained"):
983  # no line search: should be purely (or nearly) linear,
984  # and we want a large step size to initialize with.
985  fit.minimize("ModelVisit", dumpMatrixFile=dumpMatrixFile)
986  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize ModelVisit",
987  writeChi2Name=getChi2Name("ModelVisit"))
988  dumpMatrixFile = "" # so we don't redo the output on the next step
989 
990  fit.minimize("Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
991  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Model",
992  writeChi2Name=getChi2Name("Model"))
993 
994  fit.minimize("Fluxes") # no line search: always purely linear.
995  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Fluxes",
996  writeChi2Name=getChi2Name("Fluxes"))
997 
998  fit.minimize("Model Fluxes", doLineSearch=doLineSearch)
999  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize ModelFluxes",
1000  writeChi2Name=getChi2Name("ModelFluxes"))
1001 
1002  model.freezeErrorTransform()
1003  self.log.debug("Photometry error scales are frozen.")
1004 
1005  chi2 = self._iterate_fit_iterate_fit(associations,
1006  fit,
1007  self.config.maxPhotometrySteps,
1008  "photometry",
1009  "Model Fluxes",
1010  doRankUpdate=self.config.photometryDoRankUpdate,
1011  doLineSearch=doLineSearch,
1012  dataName=dataName)
1013 
1014  add_measurement(self.jobjob, 'jointcal.photometry_final_chi2', chi2.chi2)
1015  add_measurement(self.jobjob, 'jointcal.photometry_final_ndof', chi2.ndof)
1016  return Photometry(fit, model)
1017 
1018  def _fit_astrometry(self, associations, dataName=None):
1019  """
1020  Fit the astrometric data.
1021 
1022  Parameters
1023  ----------
1024  associations : `lsst.jointcal.Associations`
1025  The star/reference star associations to fit.
1026  dataName : `str`
1027  Name of the data being processed (e.g. "1234_HSC-Y"), for
1028  identifying debugging files.
1029 
1030  Returns
1031  -------
1032  fit_result : `namedtuple`
1033  fit : `lsst.jointcal.AstrometryFit`
1034  The astrometric fitter used to perform the fit.
1035  model : `lsst.jointcal.AstrometryModel`
1036  The astrometric model that was fit.
1037  sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1038  The model for the sky to tangent plane projection that was used in the fit.
1039  """
1040 
1041  self.log.info("=== Starting astrometric fitting...")
1042 
1043  associations.deprojectFittedStars()
1044 
1045  # NOTE: need to return sky_to_tan_projection so that it doesn't get garbage collected.
1046  # TODO: could we package sky_to_tan_projection and model together so we don't have to manage
1047  # them so carefully?
1048  sky_to_tan_projection = lsst.jointcal.OneTPPerVisitHandler(associations.getCcdImageList())
1049 
1050  if self.config.astrometryModel == "constrained":
1051  model = lsst.jointcal.ConstrainedAstrometryModel(associations.getCcdImageList(),
1052  sky_to_tan_projection,
1053  chipOrder=self.config.astrometryChipOrder,
1054  visitOrder=self.config.astrometryVisitOrder)
1055  elif self.config.astrometryModel == "simple":
1056  model = lsst.jointcal.SimpleAstrometryModel(associations.getCcdImageList(),
1057  sky_to_tan_projection,
1058  self.config.useInputWcs,
1059  nNotFit=0,
1060  order=self.config.astrometrySimpleOrder)
1061 
1062  fit = lsst.jointcal.AstrometryFit(associations, model, self.config.positionErrorPedestal)
1063  # TODO DM-12446: turn this into a "butler save" somehow.
1064  # Save reference and measurement chi2 contributions for this data
1065  if self.config.writeChi2FilesInitialFinal:
1066  baseName = f"astrometry_initial_chi2-{dataName}"
1067  else:
1068  baseName = None
1069  if self.config.writeInitialModel:
1070  fullpath = self._getDebugPath_getDebugPath("initialAstrometryModel.txt")
1071  writeModel(model, fullpath, self.log)
1072  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initial", writeChi2Name=baseName)
1073 
1074  def getChi2Name(whatToFit):
1075  if self.config.writeChi2FilesOuterLoop:
1076  return f"astrometry_init-%s_chi2-{dataName}" % whatToFit
1077  else:
1078  return None
1079 
1080  dumpMatrixFile = self._getDebugPath_getDebugPath("astrometry_preinit") if self.config.writeInitMatrix else ""
1081  # The constrained model needs the visit transform fit first; the chip
1082  # transform is initialized from the detector's cameraGeom, so it's close.
1083  if self.config.astrometryModel == "constrained":
1084  fit.minimize("DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1085  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize DistortionsVisit",
1086  writeChi2Name=getChi2Name("DistortionsVisit"))
1087  dumpMatrixFile = "" # so we don't redo the output on the next step
1088 
1089  fit.minimize("Distortions", dumpMatrixFile=dumpMatrixFile)
1090  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Distortions",
1091  writeChi2Name=getChi2Name("Distortions"))
1092 
1093  fit.minimize("Positions")
1094  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Positions",
1095  writeChi2Name=getChi2Name("Positions"))
1096 
1097  fit.minimize("Distortions Positions")
1098  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize DistortionsPositions",
1099  writeChi2Name=getChi2Name("DistortionsPositions"))
1100 
1101  chi2 = self._iterate_fit_iterate_fit(associations,
1102  fit,
1103  self.config.maxAstrometrySteps,
1104  "astrometry",
1105  "Distortions Positions",
1106  doRankUpdate=self.config.astrometryDoRankUpdate,
1107  dataName=dataName)
1108 
1109  add_measurement(self.jobjob, 'jointcal.astrometry_final_chi2', chi2.chi2)
1110  add_measurement(self.jobjob, 'jointcal.astrometry_final_ndof', chi2.ndof)
1111 
1112  return Astrometry(fit, model, sky_to_tan_projection)
1113 
1114  def _check_stars(self, associations):
1115  """Count measured and reference stars per ccd and warn/log them."""
1116  for ccdImage in associations.getCcdImageList():
1117  nMeasuredStars, nRefStars = ccdImage.countStars()
1118  self.log.debug("ccdImage %s has %s measured and %s reference stars",
1119  ccdImage.getName(), nMeasuredStars, nRefStars)
1120  if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1121  self.log.warn("ccdImage %s has only %s measuredStars (desired %s)",
1122  ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1123  if nRefStars < self.config.minRefStarsPerCcd:
1124  self.log.warn("ccdImage %s has only %s RefStars (desired %s)",
1125  ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1126 
1127  def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1128  dataName="",
1129  doRankUpdate=True,
1130  doLineSearch=False):
1131  """Run fitter.minimize up to max_steps times, returning the final chi2.
1132 
1133  Parameters
1134  ----------
1135  associations : `lsst.jointcal.Associations`
1136  The star/reference star associations to fit.
1137  fitter : `lsst.jointcal.FitterBase`
1138  The fitter to use for minimization.
1139  max_steps : `int`
1140  Maximum number of steps to run outlier rejection before declaring
1141  convergence failure.
1142  name : {'photometry' or 'astrometry'}
1143  What type of data are we fitting (for logs and debugging files).
1144  whatToFit : `str`
1145  Passed to ``fitter.minimize()`` to define the parameters to fit.
1146  dataName : `str`, optional
1147  Descriptive name for this dataset (e.g. tract and filter),
1148  for debugging.
1149  doRankUpdate : `bool`, optional
1150  Do an Eigen rank update during minimization, or recompute the full
1151  matrix and gradient?
1152  doLineSearch : `bool`, optional
1153  Do a line search for the optimum step during minimization?
1154 
1155  Returns
1156  -------
1157  chi2: `lsst.jointcal.Chi2Statistic`
1158  The final chi2 after the fit converges, or is forced to end.
1159 
1160  Raises
1161  ------
1162  FloatingPointError
1163  Raised if the fitter fails with a non-finite value.
1164  RuntimeError
1165  Raised if the fitter fails for some other reason;
1166  log messages will provide further details.
1167  """
1168  dumpMatrixFile = self._getDebugPath_getDebugPath(f"{name}_postinit") if self.config.writeInitMatrix else ""
1169  oldChi2 = lsst.jointcal.Chi2Statistic()
1170  oldChi2.chi2 = float("inf")
1171  for i in range(max_steps):
1172  if self.config.writeChi2FilesOuterLoop:
1173  writeChi2Name = f"{name}_iterate_{i}_chi2-{dataName}"
1174  else:
1175  writeChi2Name = None
1176  result = fitter.minimize(whatToFit,
1177  self.config.outlierRejectSigma,
1178  doRankUpdate=doRankUpdate,
1179  doLineSearch=doLineSearch,
1180  dumpMatrixFile=dumpMatrixFile)
1181  dumpMatrixFile = "" # clear it so we don't write the matrix again.
1182  chi2 = self._logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(),
1183  f"Fit iteration {i}", writeChi2Name=writeChi2Name)
1184 
1185  if result == MinimizeResult.Converged:
1186  if doRankUpdate:
1187  self.log.debug("fit has converged - no more outliers - redo minimization "
1188  "one more time in case we have lost accuracy in rank update.")
1189  # Redo minimization one more time in case we have lost accuracy in rank update
1190  result = fitter.minimize(whatToFit, self.config.outlierRejectSigma)
1191  chi2 = self._logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(), "Fit completed")
1192 
1193  # log a message for a large final chi2, TODO: DM-15247 for something better
1194  if chi2.chi2/chi2.ndof >= 4.0:
1195  self.log.error("Potentially bad fit: High chi-squared/ndof.")
1196 
1197  break
1198  elif result == MinimizeResult.Chi2Increased:
1199  self.log.warn("Still some outliers remaining but chi2 increased - retry")
1200  # Check whether the increase was large enough to cause trouble.
1201  chi2Ratio = chi2.chi2 / oldChi2.chi2
1202  if chi2Ratio > 1.5:
1203  self.log.warn('Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1204  chi2.chi2, oldChi2.chi2, chi2Ratio)
1205  # Based on a variety of HSC jointcal logs (see DM-25779), it
1206  # appears that chi2 increases more than a factor of ~2 always
1207  # result in the fit diverging rapidly and ending at chi2 > 1e10.
1208  # Using 10 as the "failure" threshold gives some room between
1209  # leaving a warning and bailing early.
1210  if chi2Ratio > 10:
1211  msg = ("Large chi2 increase between steps: fit likely cannot converge."
1212  " Try setting one or more of the `writeChi2*` config fields and looking"
1213  " at how individual star chi2-values evolve during the fit.")
1214  raise RuntimeError(msg)
1215  oldChi2 = chi2
1216  elif result == MinimizeResult.NonFinite:
1217  filename = self._getDebugPath_getDebugPath("{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1218  # TODO DM-12446: turn this into a "butler save" somehow.
1219  fitter.saveChi2Contributions(filename+"{type}")
1220  msg = "Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1221  raise FloatingPointError(msg.format(filename))
1222  elif result == MinimizeResult.Failed:
1223  raise RuntimeError("Chi2 minimization failure, cannot complete fit.")
1224  else:
1225  raise RuntimeError("Unxepected return code from minimize().")
1226  else:
1227  self.log.error("%s failed to converge after %d steps"%(name, max_steps))
1228 
1229  return chi2
1230 
1231  def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1232  """
1233  Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1234 
1235  Parameters
1236  ----------
1237  associations : `lsst.jointcal.Associations`
1238  The star/reference star associations to fit.
1239  model : `lsst.jointcal.AstrometryModel`
1240  The astrometric model that was fit.
1241  visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1242  Dict of ccdImage identifiers to dataRefs that were fit.
1243  """
1244 
1245  ccdImageList = associations.getCcdImageList()
1246  for ccdImage in ccdImageList:
1247  # TODO: there must be a better way to identify this ccdImage than a visit,ccd pair?
1248  ccd = ccdImage.ccdId
1249  visit = ccdImage.visit
1250  dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1251  self.log.info("Updating WCS for visit: %d, ccd: %d", visit, ccd)
1252  skyWcs = model.makeSkyWcs(ccdImage)
1253  try:
1254  dataRef.put(skyWcs, 'jointcal_wcs')
1255  except pexExceptions.Exception as e:
1256  self.log.fatal('Failed to write updated Wcs: %s', str(e))
1257  raise e
1258 
1259  def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1260  """
1261  Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1262 
1263  Parameters
1264  ----------
1265  associations : `lsst.jointcal.Associations`
1266  The star/reference star associations to fit.
1267  model : `lsst.jointcal.PhotometryModel`
1268  The photoometric model that was fit.
1269  visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1270  Dict of ccdImage identifiers to dataRefs that were fit.
1271  """
1272 
1273  ccdImageList = associations.getCcdImageList()
1274  for ccdImage in ccdImageList:
1275  # TODO: there must be a better way to identify this ccdImage than a visit,ccd pair?
1276  ccd = ccdImage.ccdId
1277  visit = ccdImage.visit
1278  dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1279  self.log.info("Updating PhotoCalib for visit: %d, ccd: %d", visit, ccd)
1280  photoCalib = model.toPhotoCalib(ccdImage)
1281  try:
1282  dataRef.put(photoCalib, 'jointcal_photoCalib')
1283  except pexExceptions.Exception as e:
1284  self.log.fatal('Failed to write updated PhotoCalib: %s', str(e))
1285  raise e
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:76
def _compute_proper_motion_epoch(self, ccdImageList)
Definition: jointcal.py:703
def _get_refcat_coordinate_error_override(self, refCat, name)
Definition: jointcal.py:656
def _getDebugPath(self, filename)
Definition: jointcal.py:552
def _check_star_lists(self, associations, name)
Definition: jointcal.py:863
def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1231
def _check_stars(self, associations)
Definition: jointcal.py:1114
def loadData(self, dataRefs, associations, jointcalControl, profile_jointcal=False)
Definition: jointcal.py:532
def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName, applyColorterms=False, epoch=None)
Definition: jointcal.py:808
def _readDataId(self, butler, dataId)
Definition: jointcal.py:512
def __init__(self, butler=None, profile_jointcal=False, **kwargs)
Definition: jointcal.py:430
def _build_ccdImage(self, data, associations, jointcalControl)
Definition: jointcal.py:464
def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", doRankUpdate=True, doLineSearch=False)
Definition: jointcal.py:1130
def _fit_photometry(self, associations, dataName=None)
Definition: jointcal.py:915
def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1259
def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None)
Definition: jointcal.py:872
def runDataRef(self, dataRefs, profile_jointcal=False)
Definition: jointcal.py:576
def _do_load_refcat_and_fit(self, associations, defaultBand, center, radius, tract="", profile_jointcal=False, match_cut=3.0, reject_bad_fluxes=False, *name="", refObjLoader=None, referenceSelector=None, fit_function=None, epoch=None)
Definition: jointcal.py:723
def _fit_astrometry(self, associations, dataName=None)
Definition: jointcal.py:1018
def _prep_sky(self, associations, bands)
Definition: jointcal.py:557
def add_measurement(job, name, value)
Definition: jointcal.py:58
def writeModel(model, filename, log)
Definition: jointcal.py:384