lsst.jointcal  21.0.0-5-gcff38f6+c276328c15
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.doPhotometry and self.applyColorTerms and len(self.colorterms.data) == 0:
345  msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
346  raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
347  if self.doAstrometry and not self.doPhotometry and self.applyColorTerms:
348  msg = ("Only doing astrometry, but Colorterms are not applied for astrometry;"
349  "applyColorTerms=True will be ignored.")
350  lsst.log.warn(msg)
351 
352  def setDefaults(self):
353  # Use science source selector which can filter on extendedness, SNR, and whether blended
354  self.sourceSelector.name = 'science'
355  # Use only stars because aperture fluxes of galaxies are biased and depend on seeing
356  self.sourceSelector['science'].doUnresolved = True
357  # with dependable signal to noise ratio.
358  self.sourceSelector['science'].doSignalToNoise = True
359  # Min SNR must be > 0 because jointcal cannot handle negative fluxes,
360  # and S/N > 10 to use sources that are not too faint, and thus better measured.
361  self.sourceSelector['science'].signalToNoise.minimum = 10.
362  # Base SNR on CalibFlux because that is the flux jointcal that fits and must be positive
363  fluxField = f"slot_{self.sourceFluxType}Flux_instFlux"
364  self.sourceSelector['science'].signalToNoise.fluxField = fluxField
365  self.sourceSelector['science'].signalToNoise.errField = fluxField + "Err"
366  # Do not trust blended sources' aperture fluxes which also depend on seeing
367  self.sourceSelector['science'].doIsolated = True
368  # Do not trust either flux or centroid measurements with flags,
369  # chosen from the usual QA flags for stars)
370  self.sourceSelector['science'].doFlags = True
371  badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated',
372  'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag',
373  'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter']
374  self.sourceSelector['science'].flags.bad = badFlags
375 
376  # Default to Gaia-DR2 (with proper motions) for astrometry and
377  # PS1-DR1 for photometry, with a reasonable initial filterMap.
378  self.astrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414"
379  self.astrometryRefObjLoader.requireProperMotion = True
380  self.astrometryRefObjLoader.anyFilterMapsToThis = 'phot_g_mean'
381  self.photometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110"
382 
383 
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."""
408  filter: lsst.afw.image.Filter
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_jointcal = profile_jointcal
433  self.makeSubtask("sourceSelector")
434  if self.config.doAstrometry:
435  self.makeSubtask('astrometryRefObjLoader', butler=butler)
436  self.makeSubtask("astrometryReferenceSelector")
437  else:
439  if self.config.doPhotometry:
440  self.makeSubtask('photometryRefObjLoader', butler=butler)
441  self.makeSubtask("photometryReferenceSelector")
442  else:
444 
445  # To hold various computed metrics for use by tests
446  self.job = Job.load_metrics_package(subset='jointcal')
447 
448  # We don't currently need to persist the metadata.
449  # If we do in the future, we will have to add appropriate dataset templates
450  # to each obs package (the metadata template should look like `jointcal_wcs`).
451  def _getMetadataName(self):
452  return None
453 
454  @classmethod
455  def _makeArgumentParser(cls):
456  """Create an argument parser"""
457  parser = pipeBase.ArgumentParser(name=cls._DefaultName)
458  parser.add_argument("--profile_jointcal", default=False, action="store_true",
459  help="Profile steps of jointcal separately.")
460  parser.add_id_argument("--id", "calexp", help="data ID, e.g. --id visit=6789 ccd=0..9",
461  ContainerClass=PerTractCcdDataIdContainer)
462  return parser
463 
464  def _build_ccdImage(self, data, associations, jointcalControl):
465  """
466  Extract the necessary things from this dataRef to add a new ccdImage.
467 
468  Parameters
469  ----------
470  data : `JointcalInputData`
471  The loaded input data.
472  associations : `lsst.jointcal.Associations`
473  Object to add the info to, to construct a new CcdImage
474  jointcalControl : `jointcal.JointcalControl`
475  Control object for associations management
476 
477  Returns
478  ------
479  namedtuple
480  ``wcs``
481  The TAN WCS of this image, read from the calexp
482  (`lsst.afw.geom.SkyWcs`).
483  ``key``
484  A key to identify this dataRef by its visit and ccd ids
485  (`namedtuple`).
486  ``filter``
487  This calexp's filter (`str`).
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.getName(),
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', 'filter'))
509  Key = collections.namedtuple('Key', ('visit', 'ccd'))
510  return Result(data.wcs, Key(data.visit, data.detector.getId()), data.filter.getName())
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_filter', 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  filters = []
537  load_cat_prof_file = 'jointcal_loadData.prof' if profile_jointcal else ''
538  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
539  # Need the bounding-box of the focal plane (the same for all visits) for photometry visit models
540  camera = dataRefs[0].get('camera', immediate=True)
541  self.focalPlaneBBox = camera.getFpBBox()
542  for dataRef in dataRefs:
543  data = self._readDataId(dataRef.getButler(), dataRef.dataId)
544  result = self._build_ccdImage(data, associations, jointcalControl)
545  oldWcsList.append(result.wcs)
546  visit_ccd_to_dataRef[result.key] = dataRef
547  filters.append(result.filter)
548  filters = collections.Counter(filters)
549 
550  return oldWcsList, filters, 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, filters):
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 associated with the catalog. See DM-9093
570  defaultFilter = filters.most_common(1)[0][0]
571  self.log.debug("Using %s band for reference flux", defaultFilter)
572 
573  return boundingCircle, center, radius, defaultFilter
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, filters, visit_ccd_to_dataRef = self.loadData(dataRefs,
611  associations,
612  jointcalControl,
613  profile_jointcal=profile_jointcal)
614 
615  boundingCircle, center, radius, defaultFilter = self._prep_sky(associations, filters)
616 
617  tract = dataRefs[0].dataId['tract']
618 
619  if self.config.doAstrometry:
620  astrometry = self._do_load_refcat_and_fit(associations, defaultFilter, center, radius,
621  name="astrometry",
622  refObjLoader=self.astrometryRefObjLoader,
623  referenceSelector=self.astrometryReferenceSelector,
624  fit_function=self._fit_astrometry,
625  profile_jointcal=profile_jointcal,
626  tract=tract)
627  self._write_astrometry_results(associations, astrometry.model, visit_ccd_to_dataRef)
628  else:
629  astrometry = Astrometry(None, None, None)
630 
631  if self.config.doPhotometry:
632  photometry = self._do_load_refcat_and_fit(associations, defaultFilter, center, radius,
633  name="photometry",
634  refObjLoader=self.photometryRefObjLoader,
635  referenceSelector=self.photometryReferenceSelector,
636  fit_function=self._fit_photometry,
637  profile_jointcal=profile_jointcal,
638  tract=tract,
639  filters=filters,
640  reject_bad_fluxes=True)
641  self._write_photometry_results(associations, photometry.model, visit_ccd_to_dataRef)
642  else:
643  photometry = Photometry(None, None)
644 
645  return pipeBase.Struct(dataRefs=dataRefs,
646  oldWcsList=oldWcsList,
647  job=self.job,
648  astrometryRefObjLoader=self.astrometryRefObjLoader,
649  photometryRefObjLoader=self.photometryRefObjLoader,
650  defaultFilter=defaultFilter,
651  exitStatus=exitStatus)
652 
653  def _get_refcat_coordinate_error_override(self, refCat, name):
654  """Check whether we should override the refcat coordinate errors, and
655  return the overridden error if necessary.
656 
657  Parameters
658  ----------
659  refCat : `lsst.afw.table.SimpleCatalog`
660  The reference catalog to check for a ``coord_raErr`` field.
661  name : `str`
662  Whether we are doing "astrometry" or "photometry".
663 
664  Returns
665  -------
666  refCoordErr : `float`
667  The refcat coordinate error to use, or NaN if we are not overriding
668  those fields.
669 
670  Raises
671  ------
672  lsst.pex.config.FieldValidationError
673  Raised if the refcat does not contain coordinate errors and
674  ``config.astrometryReferenceErr`` is not set.
675  """
676  # This value doesn't matter for photometry, so just set something to
677  # keep old refcats from causing problems.
678  if name.lower() == "photometry":
679  if 'coord_raErr' not in refCat.schema:
680  return 100
681  else:
682  return float('nan')
683 
684  if self.config.astrometryReferenceErr is None and 'coord_raErr' not in refCat.schema:
685  msg = ("Reference catalog does not contain coordinate errors, "
686  "and config.astrometryReferenceErr not supplied.")
687  raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
688  self.config,
689  msg)
690 
691  if self.config.astrometryReferenceErr is not None and 'coord_raErr' in refCat.schema:
692  self.log.warn("Overriding reference catalog coordinate errors with %f/coordinate [mas]",
693  self.config.astrometryReferenceErr)
694 
695  if self.config.astrometryReferenceErr is None:
696  return float('nan')
697  else:
698  return self.config.astrometryReferenceErr
699 
700  def _compute_proper_motion_epoch(self, ccdImageList):
701  """Return the proper motion correction epoch of the provided images.
702 
703  Parameters
704  ----------
705  ccdImageList : `list` [`lsst.jointcal.CcdImage`]
706  The images to compute the appropriate epoch for.
707 
708  Returns
709  -------
710  epoch : `astropy.time.Time`
711  The date to use for proper motion corrections.
712  """
713  mjds = [ccdImage.getMjd() for ccdImage in ccdImageList]
714  return astropy.time.Time(np.mean(mjds), format='mjd', scale="tai")
715 
716  def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
717  filters=[],
718  tract="", profile_jointcal=False, match_cut=3.0,
719  reject_bad_fluxes=False, *,
720  name="", refObjLoader=None, referenceSelector=None,
721  fit_function=None):
722  """Load reference catalog, perform the fit, and return the result.
723 
724  Parameters
725  ----------
726  associations : `lsst.jointcal.Associations`
727  The star/reference star associations to fit.
728  defaultFilter : `str`
729  filter to load from reference catalog.
730  center : `lsst.geom.SpherePoint`
731  ICRS center of field to load from reference catalog.
732  radius : `lsst.geom.Angle`
733  On-sky radius to load from reference catalog.
734  name : `str`
735  Name of thing being fit: "astrometry" or "photometry".
736  refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
737  Reference object loader to use to load a reference catalog.
738  referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
739  Selector to use to pick objects from the loaded reference catalog.
740  fit_function : callable
741  Function to call to perform fit (takes Associations object).
742  filters : `list` [`str`], optional
743  List of filters to load from the reference catalog.
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 
754  Returns
755  -------
756  result : `Photometry` or `Astrometry`
757  Result of `fit_function()`
758  """
759  self.log.info("====== Now processing %s...", name)
760  # TODO: this should not print "trying to invert a singular transformation:"
761  # if it does that, something's not right about the WCS...
762  associations.associateCatalogs(match_cut)
763  add_measurement(self.job, 'jointcal.associated_%s_fittedStars' % name,
764  associations.fittedStarListSize())
765 
766  applyColorterms = False if name.lower() == "astrometry" else self.config.applyColorTerms
767  epoch = self._compute_proper_motion_epoch(associations.getCcdImageList())
768  refCat, fluxField = self._load_reference_catalog(refObjLoader, referenceSelector,
769  center, radius, defaultFilter,
770  applyColorterms=applyColorterms,
771  epoch=epoch)
772  refCoordErr = self._get_refcat_coordinate_error_override(refCat, name)
773 
774  associations.collectRefStars(refCat,
775  self.config.matchCut*lsst.geom.arcseconds,
776  fluxField,
777  refCoordinateErr=refCoordErr,
778  rejectBadFluxes=reject_bad_fluxes)
779  add_measurement(self.job, 'jointcal.collected_%s_refStars' % name,
780  associations.refStarListSize())
781 
782  associations.prepareFittedStars(self.config.minMeasurements)
783 
784  self._check_star_lists(associations, name)
785  add_measurement(self.job, 'jointcal.selected_%s_refStars' % name,
786  associations.nFittedStarsWithAssociatedRefStar())
787  add_measurement(self.job, 'jointcal.selected_%s_fittedStars' % name,
788  associations.fittedStarListSize())
789  add_measurement(self.job, 'jointcal.selected_%s_ccdImages' % name,
790  associations.nCcdImagesValidForFit())
791 
792  load_cat_prof_file = 'jointcal_fit_%s.prof'%name if profile_jointcal else ''
793  dataName = "{}_{}".format(tract, defaultFilter)
794  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
795  result = fit_function(associations, dataName)
796  # TODO DM-12446: turn this into a "butler save" somehow.
797  # Save reference and measurement chi2 contributions for this data
798  if self.config.writeChi2FilesInitialFinal:
799  baseName = self._getDebugPath(f"{name}_final_chi2-{dataName}")
800  result.fit.saveChi2Contributions(baseName+"{type}")
801  self.log.info("Wrote chi2 contributions files: %s", baseName)
802 
803  return result
804 
805  def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName,
806  applyColorterms=False, epoch=None):
807  """Load the necessary reference catalog sources, convert fluxes to
808  correct units, and apply color term corrections if requested.
809 
810  Parameters
811  ----------
812  refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
813  The reference catalog loader to use to get the data.
814  referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
815  Source selector to apply to loaded reference catalog.
816  center : `lsst.geom.SpherePoint`
817  The center around which to load sources.
818  radius : `lsst.geom.Angle`
819  The radius around ``center`` to load sources in.
820  filterName : `str`
821  The name of the camera filter to load fluxes for.
822  applyColorterms : `bool`
823  Apply colorterm corrections to the refcat for ``filterName``?
824  epoch : `astropy.time.Time`, optional
825  Epoch to which to correct refcat proper motion and parallax,
826  or `None` to not apply such corrections.
827 
828  Returns
829  -------
830  refCat : `lsst.afw.table.SimpleCatalog`
831  The loaded reference catalog.
832  fluxField : `str`
833  The name of the reference catalog flux field appropriate for ``filterName``.
834  """
835  skyCircle = refObjLoader.loadSkyCircle(center,
836  radius,
837  filterName,
838  epoch=epoch)
839 
840  selected = referenceSelector.run(skyCircle.refCat)
841  # Need memory contiguity to get reference filters as a vector.
842  if not selected.sourceCat.isContiguous():
843  refCat = selected.sourceCat.copy(deep=True)
844  else:
845  refCat = selected.sourceCat
846 
847  if applyColorterms:
848  refCatName = refObjLoader.ref_dataset_name
849  self.log.info("Applying color terms for filterName=%r reference catalog=%s",
850  filterName, refCatName)
851  colorterm = self.config.colorterms.getColorterm(
852  filterName=filterName, photoCatName=refCatName, doRaise=True)
853 
854  refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat, filterName)
855  refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
856  # TODO: I didn't want to use this, but I'll deal with it in DM-16903
857  refCat[skyCircle.fluxField+'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
858 
859  return refCat, skyCircle.fluxField
860 
861  def _check_star_lists(self, associations, name):
862  # TODO: these should be len(blah), but we need this properly wrapped first.
863  if associations.nCcdImagesValidForFit() == 0:
864  raise RuntimeError('No images in the ccdImageList!')
865  if associations.fittedStarListSize() == 0:
866  raise RuntimeError('No stars in the {} fittedStarList!'.format(name))
867  if associations.refStarListSize() == 0:
868  raise RuntimeError('No stars in the {} reference star list!'.format(name))
869 
870  def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
871  """Compute chi2, log it, validate the model, and return chi2.
872 
873  Parameters
874  ----------
875  associations : `lsst.jointcal.Associations`
876  The star/reference star associations to fit.
877  fit : `lsst.jointcal.FitterBase`
878  The fitter to use for minimization.
879  model : `lsst.jointcal.Model`
880  The model being fit.
881  chi2Label : `str`
882  Label to describe the chi2 (e.g. "Initialized", "Final").
883  writeChi2Name : `str`, optional
884  Filename prefix to write the chi2 contributions to.
885  Do not supply an extension: an appropriate one will be added.
886 
887  Returns
888  -------
889  chi2: `lsst.jointcal.Chi2Accumulator`
890  The chi2 object for the current fitter and model.
891 
892  Raises
893  ------
894  FloatingPointError
895  Raised if chi2 is infinite or NaN.
896  ValueError
897  Raised if the model is not valid.
898  """
899  if writeChi2Name is not None:
900  fullpath = self._getDebugPath(writeChi2Name)
901  fit.saveChi2Contributions(fullpath+"{type}")
902  self.log.info("Wrote chi2 contributions files: %s", fullpath)
903 
904  chi2 = fit.computeChi2()
905  self.log.info("%s %s", chi2Label, chi2)
906  self._check_stars(associations)
907  if not np.isfinite(chi2.chi2):
908  raise FloatingPointError(f'{chi2Label} chi2 is invalid: {chi2}')
909  if not model.validate(associations.getCcdImageList(), chi2.ndof):
910  raise ValueError("Model is not valid: check log messages for warnings.")
911  return chi2
912 
913  def _fit_photometry(self, associations, dataName=None):
914  """
915  Fit the photometric data.
916 
917  Parameters
918  ----------
919  associations : `lsst.jointcal.Associations`
920  The star/reference star associations to fit.
921  dataName : `str`
922  Name of the data being processed (e.g. "1234_HSC-Y"), for
923  identifying debugging files.
924 
925  Returns
926  -------
927  fit_result : `namedtuple`
928  fit : `lsst.jointcal.PhotometryFit`
929  The photometric fitter used to perform the fit.
930  model : `lsst.jointcal.PhotometryModel`
931  The photometric model that was fit.
932  """
933  self.log.info("=== Starting photometric fitting...")
934 
935  # TODO: should use pex.config.RegistryField here (see DM-9195)
936  if self.config.photometryModel == "constrainedFlux":
937  model = lsst.jointcal.ConstrainedFluxModel(associations.getCcdImageList(),
938  self.focalPlaneBBox,
939  visitOrder=self.config.photometryVisitOrder,
940  errorPedestal=self.config.photometryErrorPedestal)
941  # potentially nonlinear problem, so we may need a line search to converge.
942  doLineSearch = self.config.allowLineSearch
943  elif self.config.photometryModel == "constrainedMagnitude":
944  model = lsst.jointcal.ConstrainedMagnitudeModel(associations.getCcdImageList(),
945  self.focalPlaneBBox,
946  visitOrder=self.config.photometryVisitOrder,
947  errorPedestal=self.config.photometryErrorPedestal)
948  # potentially nonlinear problem, so we may need a line search to converge.
949  doLineSearch = self.config.allowLineSearch
950  elif self.config.photometryModel == "simpleFlux":
951  model = lsst.jointcal.SimpleFluxModel(associations.getCcdImageList(),
952  errorPedestal=self.config.photometryErrorPedestal)
953  doLineSearch = False # purely linear in model parameters, so no line search needed
954  elif self.config.photometryModel == "simpleMagnitude":
955  model = lsst.jointcal.SimpleMagnitudeModel(associations.getCcdImageList(),
956  errorPedestal=self.config.photometryErrorPedestal)
957  doLineSearch = False # purely linear in model parameters, so no line search needed
958 
959  fit = lsst.jointcal.PhotometryFit(associations, model)
960  # TODO DM-12446: turn this into a "butler save" somehow.
961  # Save reference and measurement chi2 contributions for this data
962  if self.config.writeChi2FilesInitialFinal:
963  baseName = f"photometry_initial_chi2-{dataName}"
964  else:
965  baseName = None
966  if self.config.writeInitialModel:
967  fullpath = self._getDebugPath("initialPhotometryModel.txt")
968  writeModel(model, fullpath, self.log)
969  self._logChi2AndValidate(associations, fit, model, "Initialized", writeChi2Name=baseName)
970 
971  def getChi2Name(whatToFit):
972  if self.config.writeChi2FilesOuterLoop:
973  return f"photometry_init-%s_chi2-{dataName}" % whatToFit
974  else:
975  return None
976 
977  # The constrained model needs the visit transform fit first; the chip
978  # transform is initialized from the singleFrame PhotoCalib, so it's close.
979  dumpMatrixFile = self._getDebugPath("photometry_preinit") if self.config.writeInitMatrix else ""
980  if self.config.photometryModel.startswith("constrained"):
981  # no line search: should be purely (or nearly) linear,
982  # and we want a large step size to initialize with.
983  fit.minimize("ModelVisit", dumpMatrixFile=dumpMatrixFile)
984  self._logChi2AndValidate(associations, fit, model, "Initialize ModelVisit",
985  writeChi2Name=getChi2Name("ModelVisit"))
986  dumpMatrixFile = "" # so we don't redo the output on the next step
987 
988  fit.minimize("Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
989  self._logChi2AndValidate(associations, fit, model, "Initialize Model",
990  writeChi2Name=getChi2Name("Model"))
991 
992  fit.minimize("Fluxes") # no line search: always purely linear.
993  self._logChi2AndValidate(associations, fit, model, "Initialize Fluxes",
994  writeChi2Name=getChi2Name("Fluxes"))
995 
996  fit.minimize("Model Fluxes", doLineSearch=doLineSearch)
997  self._logChi2AndValidate(associations, fit, model, "Initialize ModelFluxes",
998  writeChi2Name=getChi2Name("ModelFluxes"))
999 
1000  model.freezeErrorTransform()
1001  self.log.debug("Photometry error scales are frozen.")
1002 
1003  chi2 = self._iterate_fit(associations,
1004  fit,
1005  self.config.maxPhotometrySteps,
1006  "photometry",
1007  "Model Fluxes",
1008  doRankUpdate=self.config.photometryDoRankUpdate,
1009  doLineSearch=doLineSearch,
1010  dataName=dataName)
1011 
1012  add_measurement(self.job, 'jointcal.photometry_final_chi2', chi2.chi2)
1013  add_measurement(self.job, 'jointcal.photometry_final_ndof', chi2.ndof)
1014  return Photometry(fit, model)
1015 
1016  def _fit_astrometry(self, associations, dataName=None):
1017  """
1018  Fit the astrometric data.
1019 
1020  Parameters
1021  ----------
1022  associations : `lsst.jointcal.Associations`
1023  The star/reference star associations to fit.
1024  dataName : `str`
1025  Name of the data being processed (e.g. "1234_HSC-Y"), for
1026  identifying debugging files.
1027 
1028  Returns
1029  -------
1030  fit_result : `namedtuple`
1031  fit : `lsst.jointcal.AstrometryFit`
1032  The astrometric fitter used to perform the fit.
1033  model : `lsst.jointcal.AstrometryModel`
1034  The astrometric model that was fit.
1035  sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1036  The model for the sky to tangent plane projection that was used in the fit.
1037  """
1038 
1039  self.log.info("=== Starting astrometric fitting...")
1040 
1041  associations.deprojectFittedStars()
1042 
1043  # NOTE: need to return sky_to_tan_projection so that it doesn't get garbage collected.
1044  # TODO: could we package sky_to_tan_projection and model together so we don't have to manage
1045  # them so carefully?
1046  sky_to_tan_projection = lsst.jointcal.OneTPPerVisitHandler(associations.getCcdImageList())
1047 
1048  if self.config.astrometryModel == "constrained":
1049  model = lsst.jointcal.ConstrainedAstrometryModel(associations.getCcdImageList(),
1050  sky_to_tan_projection,
1051  chipOrder=self.config.astrometryChipOrder,
1052  visitOrder=self.config.astrometryVisitOrder)
1053  elif self.config.astrometryModel == "simple":
1054  model = lsst.jointcal.SimpleAstrometryModel(associations.getCcdImageList(),
1055  sky_to_tan_projection,
1056  self.config.useInputWcs,
1057  nNotFit=0,
1058  order=self.config.astrometrySimpleOrder)
1059 
1060  fit = lsst.jointcal.AstrometryFit(associations, model, self.config.positionErrorPedestal)
1061  # TODO DM-12446: turn this into a "butler save" somehow.
1062  # Save reference and measurement chi2 contributions for this data
1063  if self.config.writeChi2FilesInitialFinal:
1064  baseName = f"astrometry_initial_chi2-{dataName}"
1065  else:
1066  baseName = None
1067  if self.config.writeInitialModel:
1068  fullpath = self._getDebugPath("initialAstrometryModel.txt")
1069  writeModel(model, fullpath, self.log)
1070  self._logChi2AndValidate(associations, fit, model, "Initial", writeChi2Name=baseName)
1071 
1072  def getChi2Name(whatToFit):
1073  if self.config.writeChi2FilesOuterLoop:
1074  return f"astrometry_init-%s_chi2-{dataName}" % whatToFit
1075  else:
1076  return None
1077 
1078  dumpMatrixFile = self._getDebugPath("astrometry_preinit") if self.config.writeInitMatrix else ""
1079  # The constrained model needs the visit transform fit first; the chip
1080  # transform is initialized from the detector's cameraGeom, so it's close.
1081  if self.config.astrometryModel == "constrained":
1082  fit.minimize("DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1083  self._logChi2AndValidate(associations, fit, model, "Initialize DistortionsVisit",
1084  writeChi2Name=getChi2Name("DistortionsVisit"))
1085  dumpMatrixFile = "" # so we don't redo the output on the next step
1086 
1087  fit.minimize("Distortions", dumpMatrixFile=dumpMatrixFile)
1088  self._logChi2AndValidate(associations, fit, model, "Initialize Distortions",
1089  writeChi2Name=getChi2Name("Distortions"))
1090 
1091  fit.minimize("Positions")
1092  self._logChi2AndValidate(associations, fit, model, "Initialize Positions",
1093  writeChi2Name=getChi2Name("Positions"))
1094 
1095  fit.minimize("Distortions Positions")
1096  self._logChi2AndValidate(associations, fit, model, "Initialize DistortionsPositions",
1097  writeChi2Name=getChi2Name("DistortionsPositions"))
1098 
1099  chi2 = self._iterate_fit(associations,
1100  fit,
1101  self.config.maxAstrometrySteps,
1102  "astrometry",
1103  "Distortions Positions",
1104  doRankUpdate=self.config.astrometryDoRankUpdate,
1105  dataName=dataName)
1106 
1107  add_measurement(self.job, 'jointcal.astrometry_final_chi2', chi2.chi2)
1108  add_measurement(self.job, 'jointcal.astrometry_final_ndof', chi2.ndof)
1109 
1110  return Astrometry(fit, model, sky_to_tan_projection)
1111 
1112  def _check_stars(self, associations):
1113  """Count measured and reference stars per ccd and warn/log them."""
1114  for ccdImage in associations.getCcdImageList():
1115  nMeasuredStars, nRefStars = ccdImage.countStars()
1116  self.log.debug("ccdImage %s has %s measured and %s reference stars",
1117  ccdImage.getName(), nMeasuredStars, nRefStars)
1118  if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1119  self.log.warn("ccdImage %s has only %s measuredStars (desired %s)",
1120  ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1121  if nRefStars < self.config.minRefStarsPerCcd:
1122  self.log.warn("ccdImage %s has only %s RefStars (desired %s)",
1123  ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1124 
1125  def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1126  dataName="",
1127  doRankUpdate=True,
1128  doLineSearch=False):
1129  """Run fitter.minimize up to max_steps times, returning the final chi2.
1130 
1131  Parameters
1132  ----------
1133  associations : `lsst.jointcal.Associations`
1134  The star/reference star associations to fit.
1135  fitter : `lsst.jointcal.FitterBase`
1136  The fitter to use for minimization.
1137  max_steps : `int`
1138  Maximum number of steps to run outlier rejection before declaring
1139  convergence failure.
1140  name : {'photometry' or 'astrometry'}
1141  What type of data are we fitting (for logs and debugging files).
1142  whatToFit : `str`
1143  Passed to ``fitter.minimize()`` to define the parameters to fit.
1144  dataName : `str`, optional
1145  Descriptive name for this dataset (e.g. tract and filter),
1146  for debugging.
1147  doRankUpdate : `bool`, optional
1148  Do an Eigen rank update during minimization, or recompute the full
1149  matrix and gradient?
1150  doLineSearch : `bool`, optional
1151  Do a line search for the optimum step during minimization?
1152 
1153  Returns
1154  -------
1155  chi2: `lsst.jointcal.Chi2Statistic`
1156  The final chi2 after the fit converges, or is forced to end.
1157 
1158  Raises
1159  ------
1160  FloatingPointError
1161  Raised if the fitter fails with a non-finite value.
1162  RuntimeError
1163  Raised if the fitter fails for some other reason;
1164  log messages will provide further details.
1165  """
1166  dumpMatrixFile = self._getDebugPath(f"{name}_postinit") if self.config.writeInitMatrix else ""
1167  oldChi2 = lsst.jointcal.Chi2Statistic()
1168  oldChi2.chi2 = float("inf")
1169  for i in range(max_steps):
1170  if self.config.writeChi2FilesOuterLoop:
1171  writeChi2Name = f"{name}_iterate_{i}_chi2-{dataName}"
1172  else:
1173  writeChi2Name = None
1174  result = fitter.minimize(whatToFit,
1175  self.config.outlierRejectSigma,
1176  doRankUpdate=doRankUpdate,
1177  doLineSearch=doLineSearch,
1178  dumpMatrixFile=dumpMatrixFile)
1179  dumpMatrixFile = "" # clear it so we don't write the matrix again.
1180  chi2 = self._logChi2AndValidate(associations, fitter, fitter.getModel(),
1181  f"Fit iteration {i}", writeChi2Name=writeChi2Name)
1182 
1183  if result == MinimizeResult.Converged:
1184  if doRankUpdate:
1185  self.log.debug("fit has converged - no more outliers - redo minimization "
1186  "one more time in case we have lost accuracy in rank update.")
1187  # Redo minimization one more time in case we have lost accuracy in rank update
1188  result = fitter.minimize(whatToFit, self.config.outlierRejectSigma)
1189  chi2 = self._logChi2AndValidate(associations, fitter, fitter.getModel(), "Fit completed")
1190 
1191  # log a message for a large final chi2, TODO: DM-15247 for something better
1192  if chi2.chi2/chi2.ndof >= 4.0:
1193  self.log.error("Potentially bad fit: High chi-squared/ndof.")
1194 
1195  break
1196  elif result == MinimizeResult.Chi2Increased:
1197  self.log.warn("Still some outliers remaining but chi2 increased - retry")
1198  # Check whether the increase was large enough to cause trouble.
1199  chi2Ratio = chi2.chi2 / oldChi2.chi2
1200  if chi2Ratio > 1.5:
1201  self.log.warn('Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1202  chi2.chi2, oldChi2.chi2, chi2Ratio)
1203  # Based on a variety of HSC jointcal logs (see DM-25779), it
1204  # appears that chi2 increases more than a factor of ~2 always
1205  # result in the fit diverging rapidly and ending at chi2 > 1e10.
1206  # Using 10 as the "failure" threshold gives some room between
1207  # leaving a warning and bailing early.
1208  if chi2Ratio > 10:
1209  msg = ("Large chi2 increase between steps: fit likely cannot converge."
1210  " Try setting one or more of the `writeChi2*` config fields and looking"
1211  " at how individual star chi2-values evolve during the fit.")
1212  raise RuntimeError(msg)
1213  oldChi2 = chi2
1214  elif result == MinimizeResult.NonFinite:
1215  filename = self._getDebugPath("{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1216  # TODO DM-12446: turn this into a "butler save" somehow.
1217  fitter.saveChi2Contributions(filename+"{type}")
1218  msg = "Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1219  raise FloatingPointError(msg.format(filename))
1220  elif result == MinimizeResult.Failed:
1221  raise RuntimeError("Chi2 minimization failure, cannot complete fit.")
1222  else:
1223  raise RuntimeError("Unxepected return code from minimize().")
1224  else:
1225  self.log.error("%s failed to converge after %d steps"%(name, max_steps))
1226 
1227  return chi2
1228 
1229  def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1230  """
1231  Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1232 
1233  Parameters
1234  ----------
1235  associations : `lsst.jointcal.Associations`
1236  The star/reference star associations to fit.
1237  model : `lsst.jointcal.AstrometryModel`
1238  The astrometric model that was fit.
1239  visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1240  Dict of ccdImage identifiers to dataRefs that were fit.
1241  """
1242 
1243  ccdImageList = associations.getCcdImageList()
1244  for ccdImage in ccdImageList:
1245  # TODO: there must be a better way to identify this ccdImage than a visit,ccd pair?
1246  ccd = ccdImage.ccdId
1247  visit = ccdImage.visit
1248  dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1249  self.log.info("Updating WCS for visit: %d, ccd: %d", visit, ccd)
1250  skyWcs = model.makeSkyWcs(ccdImage)
1251  try:
1252  dataRef.put(skyWcs, 'jointcal_wcs')
1253  except pexExceptions.Exception as e:
1254  self.log.fatal('Failed to write updated Wcs: %s', str(e))
1255  raise e
1256 
1257  def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1258  """
1259  Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1260 
1261  Parameters
1262  ----------
1263  associations : `lsst.jointcal.Associations`
1264  The star/reference star associations to fit.
1265  model : `lsst.jointcal.PhotometryModel`
1266  The photoometric model that was fit.
1267  visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1268  Dict of ccdImage identifiers to dataRefs that were fit.
1269  """
1270 
1271  ccdImageList = associations.getCcdImageList()
1272  for ccdImage in ccdImageList:
1273  # TODO: there must be a better way to identify this ccdImage than a visit,ccd pair?
1274  ccd = ccdImage.ccdId
1275  visit = ccdImage.visit
1276  dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1277  self.log.info("Updating PhotoCalib for visit: %d, ccd: %d", visit, ccd)
1278  photoCalib = model.toPhotoCalib(ccdImage)
1279  try:
1280  dataRef.put(photoCalib, 'jointcal_photoCalib')
1281  except pexExceptions.Exception as e:
1282  self.log.fatal('Failed to write updated PhotoCalib: %s', str(e))
1283  raise e
lsst::afw::image
lsst::jointcal.jointcal.JointcalTask._write_photometry_results
def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1257
lsst::jointcal.jointcal.JointcalTask.__init__
def __init__(self, butler=None, profile_jointcal=False, **kwargs)
Definition: jointcal.py:430
lsst::jointcal.jointcal.JointcalTask.astrometryRefObjLoader
astrometryRefObjLoader
Definition: jointcal.py:438
lsst::jointcal.jointcal.JointcalTask._getDebugPath
def _getDebugPath(self, filename)
Definition: jointcal.py:552
lsst::jointcal.jointcal.Astrometry
Astrometry
Definition: jointcal.py:54
lsst::meas::algorithms::sourceSelector
lsst::jointcal.jointcal.JointcalTask._readDataId
def _readDataId(self, butler, dataId)
Definition: jointcal.py:512
lsst::jointcal.jointcal.JointcalConfig
Definition: jointcal.py:141
lsst::jointcal.jointcal.JointcalTask._fit_photometry
def _fit_photometry(self, associations, dataName=None)
Definition: jointcal.py:913
lsst::jointcal.jointcal.add_measurement
def add_measurement(job, name, value)
Definition: jointcal.py:58
lsst::afw::image::Filter
lsst::jointcal.jointcal.JointcalRunner.getTargetList
def getTargetList(parsedCmd, **kwargs)
Definition: jointcal.py:76
lsst::jointcal.jointcal.JointcalTask._load_reference_catalog
def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName, applyColorterms=False, epoch=None)
Definition: jointcal.py:805
lsst::jointcal::SimpleMagnitudeModel
Definition: SimplePhotometryModel.h:123
lsst::jointcal.jointcal.JointcalTask._do_load_refcat_and_fit
def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius, filters=[], tract="", profile_jointcal=False, match_cut=3.0, reject_bad_fluxes=False, *name="", refObjLoader=None, referenceSelector=None, fit_function=None)
Definition: jointcal.py:716
lsst::jointcal.jointcal.JointcalTask.loadData
def loadData(self, dataRefs, associations, jointcalControl, profile_jointcal=False)
Definition: jointcal.py:532
lsst::jointcal.jointcal.JointcalTask.job
job
Definition: jointcal.py:446
lsst::jointcal.jointcal.Photometry
Photometry
Definition: jointcal.py:53
lsst::jointcal.jointcal.JointcalInputData
Definition: jointcal.py:392
lsst::jointcal.jointcal.JointcalTask._iterate_fit
def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", doRankUpdate=True, doLineSearch=False)
Definition: jointcal.py:1125
lsst::jointcal::ConstrainedMagnitudeModel
Definition: ConstrainedPhotometryModel.h:183
lsst::jointcal.jointcal.JointcalConfig.setDefaults
def setDefaults(self)
Definition: jointcal.py:352
lsst::jointcal::PhotometryFit
Class that handles the photometric least squares problem.
Definition: PhotometryFit.h:45
lsst::jointcal::ConstrainedFluxModel
Definition: ConstrainedPhotometryModel.h:137
lsst::jointcal.jointcal.JointcalTask._DefaultName
string _DefaultName
Definition: jointcal.py:428
lsst::jointcal::Chi2Statistic
Simple structure to accumulate chi2 and ndof.
Definition: Chi2.h:52
lsst::afw::image::VisitInfo
lsst::jointcal.jointcal.JointcalRunner.__call__
def __call__(self, args)
Definition: jointcal.py:93
lsst::jointcal.jointcal.JointcalTask._build_ccdImage
def _build_ccdImage(self, data, associations, jointcalControl)
Definition: jointcal.py:464
lsst::jointcal.jointcal.JointcalTask._write_astrometry_results
def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1229
lsst::jointcal.jointcal.JointcalTask.profile_jointcal
profile_jointcal
Definition: jointcal.py:432
lsst::jointcal.jointcal.JointcalTask._check_star_lists
def _check_star_lists(self, associations, name)
Definition: jointcal.py:861
lsst::jointcal::OneTPPerVisitHandler
A projection handler in which all CCDs from the same visit have the same tangent point.
Definition: ProjectionHandler.h:80
lsst::jointcal::SimpleAstrometryModel
A model where there is one independent transform per CcdImage.
Definition: SimpleAstrometryModel.h:62
lsst::jointcal.jointcal.JointcalTask
Definition: jointcal.py:412
lsst::jointcal.jointcal.JointcalTask.focalPlaneBBox
focalPlaneBBox
Definition: jointcal.py:541
lsst::afw::table
lsst::afw::image::PhotoCalib
lsst::jointcal::SimpleFluxModel
Definition: SimplePhotometryModel.h:86
lsst::utils
lsst::pex::config
lsst::jointcal
Definition: Associations.h:49
lsst::jointcal.jointcal.JointcalTask._check_stars
def _check_stars(self, associations)
Definition: jointcal.py:1112
lsst::jointcal.jointcal.JointcalTask._fit_astrometry
def _fit_astrometry(self, associations, dataName=None)
Definition: jointcal.py:1016
lsst::jointcal.jointcal.JointcalConfig.colorterms
colorterms
Definition: jointcal.py:243
lsst::jointcal::Associations
The class that implements the relations between MeasuredStar and FittedStar.
Definition: Associations.h:54
lsst::geom
lsst::jointcal.jointcal.JointcalTask.runDataRef
def runDataRef(self, dataRefs, profile_jointcal=False)
Definition: jointcal.py:576
lsst::afw::cameraGeom
lsst::jointcal::JointcalControl
Definition: JointcalControl.h:39
lsst::afw::geom::skyWcs
lsst::geom::Angle
lsst::geom::Box2I
lsst::jointcal.jointcal.JointcalTask._prep_sky
def _prep_sky(self, associations, filters)
Definition: jointcal.py:557
lsst::jointcal.jointcal.JointcalTask._compute_proper_motion_epoch
def _compute_proper_motion_epoch(self, ccdImageList)
Definition: jointcal.py:700
lsst::jointcal.jointcal.JointcalRunner
Definition: jointcal.py:63
lsst::jointcal.jointcal.JointcalTask.photometryRefObjLoader
photometryRefObjLoader
Definition: jointcal.py:443
lsst::jointcal::ConstrainedAstrometryModel
A multi-component model, fitting mappings for sensors and visits simultaneously.
Definition: ConstrainedAstrometryModel.h:62
lsst::jointcal::AstrometryFit
Class that handles the astrometric least squares problem.
Definition: AstrometryFit.h:78
lsst::pex::exceptions
lsst::geom::SpherePoint
lsst::jointcal.jointcal.JointcalTask._logChi2AndValidate
def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None)
Definition: jointcal.py:870
lsst::jointcal.jointcal.JointcalConfig.astrometryRefObjLoader
astrometryRefObjLoader
Definition: jointcal.py:279
lsst::jointcal.jointcal.JointcalConfig.doPhotometry
doPhotometry
Definition: jointcal.py:149
lsst::jointcal.jointcal.JointcalConfig.validate
def validate(self)
Definition: jointcal.py:342
lsst::jointcal.jointcal.JointcalConfig.photometryRefObjLoader
photometryRefObjLoader
Definition: jointcal.py:283
lsst::pipe::base
lsst::jointcal.jointcal.writeModel
def writeModel(model, filename, log)
Definition: jointcal.py:384
lsst::log
lsst::meas::algorithms
lsst::jointcal.jointcal.JointcalConfig.sourceSelector
sourceSelector
Definition: jointcal.py:287
lsst::jointcal.jointcal.JointcalConfig.doAstrometry
doAstrometry
Definition: jointcal.py:144
lsst::jointcal.jointcal.JointcalConfig.applyColorTerms
applyColorTerms
Definition: jointcal.py:237
lsst::jointcal.jointcal.JointcalTask._get_refcat_coordinate_error_override
def _get_refcat_coordinate_error_override(self, refCat, name)
Definition: jointcal.py:653
lsst::afw::table::SortedCatalogT
lsst::afw::cameraGeom::Detector