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