lsst.jointcal  21.0.0-11-g25eff31+3a0149a11c
jointcal.py
Go to the documentation of this file.
1 # This file is part of jointcal.
2 #
3 # Developed for the LSST Data Management System.
4 # This product includes software developed by the LSST Project
5 # (https://www.lsst.org).
6 # See the COPYRIGHT file at the top-level directory of this distribution
7 # for details of code ownership.
8 #
9 # This program is free software: you can redistribute it and/or modify
10 # it under the terms of the GNU General Public License as published by
11 # the Free Software Foundation, either version 3 of the License, or
12 # (at your option) any later version.
13 #
14 # This program is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 # GNU General Public License for more details.
18 #
19 # You should have received a copy of the GNU General Public License
20 # along with this program. If not, see <https://www.gnu.org/licenses/>.
21 
22 import dataclasses
23 import collections
24 import os
25 
26 import astropy.time
27 import numpy as np
28 import astropy.units as u
29 
30 import lsst.geom
31 import lsst.utils
32 import lsst.pex.config as pexConfig
33 import lsst.pipe.base as pipeBase
34 from lsst.afw.image import fluxErrFromABMagErr
35 import lsst.pex.exceptions as pexExceptions
37 import lsst.afw.table
38 import lsst.log
39 from lsst.pipe.tasks.colorterms import ColortermLibrary
40 from lsst.verify import Job, Measurement
41 
42 from lsst.meas.algorithms import (LoadIndexedReferenceObjectsTask, ReferenceSourceSelectorTask,
43  ReferenceObjectLoader)
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 (gen2)
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['butler'] = parsedCmd.butler
83 
84  # organize data IDs by tract
85  refListDict = {}
86  for ref in parsedCmd.id.refList:
87  refListDict.setdefault(ref.dataId["tract"], []).append(ref)
88  # we call runDataRef() once with each tract
89  result = [(refListDict[tract], kwargs) for tract in sorted(refListDict.keys())]
90  return result
91 
92  def __call__(self, args):
93  """
94  Parameters
95  ----------
96  args
97  Arguments for Task.runDataRef()
98 
99  Returns
100  -------
101  pipe.base.Struct
102  if self.doReturnResults is False:
103 
104  - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
105 
106  if self.doReturnResults is True:
107 
108  - ``result``: the result of calling jointcal.runDataRef()
109  - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
110  """
111  exitStatus = 0 # exit status for shell
112 
113  # NOTE: cannot call self.makeTask because that assumes args[0] is a single dataRef.
114  dataRefList, kwargs = args
115  butler = kwargs.pop('butler')
116  task = self.TaskClass(config=self.config, log=self.log, butler=butler)
117  result = None
118  try:
119  result = task.runDataRef(dataRefList, **kwargs)
120  exitStatus = result.exitStatus
121  job_path = butler.get('verify_job_filename')
122  result.job.write(job_path[0])
123  except Exception as e: # catch everything, sort it out later.
124  if self.doRaise:
125  raise e
126  else:
127  exitStatus = 1
128  eName = type(e).__name__
129  tract = dataRefList[0].dataId['tract']
130  task.log.fatal("Failed processing tract %s, %s: %s", tract, eName, e)
131 
132  # Put the butler back into kwargs for the other Tasks.
133  kwargs['butler'] = butler
134  if self.doReturnResults:
135  return pipeBase.Struct(result=result, exitStatus=exitStatus)
136  else:
137  return pipeBase.Struct(exitStatus=exitStatus)
138 
139 
140 class JointcalTaskConnections(pipeBase.PipelineTaskConnections,
141  dimensions=("skymap", "tract", "instrument", "physical_filter")):
142  """Middleware input/output connections for jointcal data."""
143  inputCamera = pipeBase.connectionTypes.Input(
144  doc="The camera instrument that took these observations.",
145  name="camera",
146  storageClass="Camera",
147  dimensions=("instrument",),
148  isCalibration=True
149  )
150  inputSourceTableVisit = pipeBase.connectionTypes.Input(
151  doc="Source table in parquet format, per visit",
152  name="sourceTable_visit",
153  storageClass="DataFrame",
154  dimensions=("instrument", "visit"),
155  deferLoad=True,
156  multiple=True,
157  )
158  inputVisitSummary = pipeBase.connectionTypes.Input(
159  doc=("Per-visit consolidated exposure metadata built from calexps. "
160  "These catalogs use detector id for the id and must be sorted for "
161  "fast lookups of a detector."),
162  name="visitSummary",
163  storageClass="ExposureCatalog",
164  dimensions=("instrument", "visit"),
165  deferLoad=True,
166  multiple=True,
167  )
168  # TODO DM-28991: This does not load enough refcat data: the graph only to gets
169  # refcats that touch the tract, but we need to go larger because we're
170  # loading visits that may extend further. How to test that?
171  astrometryRefCat = pipeBase.connectionTypes.Input(
172  doc="The astrometry reference catalog to match to loaded input catalog sources.",
173  name="gaia_dr2_20200414",
174  storageClass="SimpleCatalog",
175  dimensions=("htm7",),
176  deferLoad=True,
177  multiple=True
178  )
179  photometryRefCat = pipeBase.connectionTypes.Input(
180  doc="The photometry reference catalog to match to loaded input catalog sources.",
181  name="ps1_pv3_3pi_20170110",
182  storageClass="SimpleCatalog",
183  dimensions=("htm7",),
184  deferLoad=True,
185  multiple=True
186  )
187 
188  outputWcs = pipeBase.connectionTypes.Output(
189  doc=("Per-tract, per-visit world coordinate systems derived from the fitted model."
190  " These catalogs only contain entries for detectors with an output, and use"
191  " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
192  name="jointcalSkyWcsCatalog",
193  storageClass="ExposureCatalog",
194  dimensions=("instrument", "visit", "skymap", "tract"),
195  multiple=True
196  )
197  outputPhotoCalib = pipeBase.connectionTypes.Output(
198  doc=("Per-tract, per-visit photometric calibrations derived from the fitted model."
199  " These catalogs only contain entries for detectors with an output, and use"
200  " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
201  name="jointcalPhotoCalibCatalog",
202  storageClass="ExposureCatalog",
203  dimensions=("instrument", "visit", "skymap", "tract"),
204  multiple=True
205  )
206 
207 
208 class JointcalConfig(pipeBase.PipelineTaskConfig,
209  pipelineConnections=JointcalTaskConnections):
210  """Configuration for JointcalTask"""
211 
212  doAstrometry = pexConfig.Field(
213  doc="Fit astrometry and write the fitted result.",
214  dtype=bool,
215  default=True
216  )
217  doPhotometry = pexConfig.Field(
218  doc="Fit photometry and write the fitted result.",
219  dtype=bool,
220  default=True
221  )
222  coaddName = pexConfig.Field(
223  doc="Type of coadd, typically deep or goodSeeing",
224  dtype=str,
225  default="deep"
226  )
227  # TODO DM-29008: Change this to "ApFlux_12_0" before gen2 removal.
228  sourceFluxType = pexConfig.Field(
229  dtype=str,
230  doc="Source flux field to use in source selection and to get fluxes from the catalog.",
231  default='Calib'
232  )
233  positionErrorPedestal = pexConfig.Field(
234  doc="Systematic term to apply to the measured position error (pixels)",
235  dtype=float,
236  default=0.02,
237  )
238  photometryErrorPedestal = pexConfig.Field(
239  doc="Systematic term to apply to the measured error on flux or magnitude as a "
240  "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
241  dtype=float,
242  default=0.0,
243  )
244  # TODO: DM-6885 matchCut should be an geom.Angle
245  matchCut = pexConfig.Field(
246  doc="Matching radius between fitted and reference stars (arcseconds)",
247  dtype=float,
248  default=3.0,
249  )
250  minMeasurements = pexConfig.Field(
251  doc="Minimum number of associated measured stars for a fitted star to be included in the fit",
252  dtype=int,
253  default=2,
254  )
255  minMeasuredStarsPerCcd = pexConfig.Field(
256  doc="Minimum number of measuredStars per ccdImage before printing warnings",
257  dtype=int,
258  default=100,
259  )
260  minRefStarsPerCcd = pexConfig.Field(
261  doc="Minimum number of measuredStars per ccdImage before printing warnings",
262  dtype=int,
263  default=30,
264  )
265  allowLineSearch = pexConfig.Field(
266  doc="Allow a line search during minimization, if it is reasonable for the model"
267  " (models with a significant non-linear component, e.g. constrainedPhotometry).",
268  dtype=bool,
269  default=False
270  )
271  astrometrySimpleOrder = pexConfig.Field(
272  doc="Polynomial order for fitting the simple astrometry model.",
273  dtype=int,
274  default=3,
275  )
276  astrometryChipOrder = pexConfig.Field(
277  doc="Order of the per-chip transform for the constrained astrometry model.",
278  dtype=int,
279  default=1,
280  )
281  astrometryVisitOrder = pexConfig.Field(
282  doc="Order of the per-visit transform for the constrained astrometry model.",
283  dtype=int,
284  default=5,
285  )
286  useInputWcs = pexConfig.Field(
287  doc="Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
288  dtype=bool,
289  default=True,
290  )
291  astrometryModel = pexConfig.ChoiceField(
292  doc="Type of model to fit to astrometry",
293  dtype=str,
294  default="constrained",
295  allowed={"simple": "One polynomial per ccd",
296  "constrained": "One polynomial per ccd, and one polynomial per visit"}
297  )
298  photometryModel = pexConfig.ChoiceField(
299  doc="Type of model to fit to photometry",
300  dtype=str,
301  default="constrainedMagnitude",
302  allowed={"simpleFlux": "One constant zeropoint per ccd and visit, fitting in flux space.",
303  "constrainedFlux": "Constrained zeropoint per ccd, and one polynomial per visit,"
304  " fitting in flux space.",
305  "simpleMagnitude": "One constant zeropoint per ccd and visit,"
306  " fitting in magnitude space.",
307  "constrainedMagnitude": "Constrained zeropoint per ccd, and one polynomial per visit,"
308  " fitting in magnitude space.",
309  }
310  )
311  applyColorTerms = pexConfig.Field(
312  doc="Apply photometric color terms to reference stars?"
313  "Requires that colorterms be set to a ColortermLibrary",
314  dtype=bool,
315  default=False
316  )
317  colorterms = pexConfig.ConfigField(
318  doc="Library of photometric reference catalog name to color term dict.",
319  dtype=ColortermLibrary,
320  )
321  photometryVisitOrder = pexConfig.Field(
322  doc="Order of the per-visit polynomial transform for the constrained photometry model.",
323  dtype=int,
324  default=7,
325  )
326  photometryDoRankUpdate = pexConfig.Field(
327  doc=("Do the rank update step during minimization. "
328  "Skipping this can help deal with models that are too non-linear."),
329  dtype=bool,
330  default=True,
331  )
332  astrometryDoRankUpdate = pexConfig.Field(
333  doc=("Do the rank update step during minimization (should not change the astrometry fit). "
334  "Skipping this can help deal with models that are too non-linear."),
335  dtype=bool,
336  default=True,
337  )
338  outlierRejectSigma = pexConfig.Field(
339  doc="How many sigma to reject outliers at during minimization.",
340  dtype=float,
341  default=5.0,
342  )
343  maxPhotometrySteps = pexConfig.Field(
344  doc="Maximum number of minimize iterations to take when fitting photometry.",
345  dtype=int,
346  default=20,
347  )
348  maxAstrometrySteps = pexConfig.Field(
349  doc="Maximum number of minimize iterations to take when fitting photometry.",
350  dtype=int,
351  default=20,
352  )
353  astrometryRefObjLoader = pexConfig.ConfigurableField(
354  target=LoadIndexedReferenceObjectsTask,
355  doc="Reference object loader for astrometric fit",
356  )
357  photometryRefObjLoader = pexConfig.ConfigurableField(
358  target=LoadIndexedReferenceObjectsTask,
359  doc="Reference object loader for photometric fit",
360  )
361  sourceSelector = sourceSelectorRegistry.makeField(
362  doc="How to select sources for cross-matching",
363  default="astrometry"
364  )
365  astrometryReferenceSelector = pexConfig.ConfigurableField(
366  target=ReferenceSourceSelectorTask,
367  doc="How to down-select the loaded astrometry reference catalog.",
368  )
369  photometryReferenceSelector = pexConfig.ConfigurableField(
370  target=ReferenceSourceSelectorTask,
371  doc="How to down-select the loaded photometry reference catalog.",
372  )
373  astrometryReferenceErr = pexConfig.Field(
374  doc=("Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
375  "If None, then raise an exception if the reference catalog is missing coordinate errors. "
376  "If specified, overrides any existing `coord_*Err` values."),
377  dtype=float,
378  default=None,
379  optional=True
380  )
381 
382  # configs for outputting debug information
383  writeInitMatrix = pexConfig.Field(
384  dtype=bool,
385  doc=("Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
386  "The output files will be of the form 'astrometry_preinit-mat.txt', in the current directory. "
387  "Note that these files are the dense versions of the matrix, and so may be very large."),
388  default=False
389  )
390  writeChi2FilesInitialFinal = pexConfig.Field(
391  dtype=bool,
392  doc="Write .csv files containing the contributions to chi2 for the initialization and final fit.",
393  default=False
394  )
395  writeChi2FilesOuterLoop = pexConfig.Field(
396  dtype=bool,
397  doc="Write .csv files containing the contributions to chi2 for the outer fit loop.",
398  default=False
399  )
400  writeInitialModel = pexConfig.Field(
401  dtype=bool,
402  doc=("Write the pre-initialization model to text files, for debugging."
403  " Output is written to `initial[Astro|Photo]metryModel.txt` in the current working directory."),
404  default=False
405  )
406  debugOutputPath = pexConfig.Field(
407  dtype=str,
408  default=".",
409  doc=("Path to write debug output files to. Used by "
410  "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
411  )
412  detailedProfile = pexConfig.Field(
413  dtype=bool,
414  default=False,
415  doc="Output separate profiling information for different parts of jointcal, e.g. data read, fitting"
416  )
417 
418  def validate(self):
419  super().validate()
420  if self.doPhotometrydoPhotometry and self.applyColorTermsapplyColorTerms and len(self.colortermscolorterms.data) == 0:
421  msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
422  raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
423  if self.doAstrometrydoAstrometry and not self.doPhotometrydoPhotometry and self.applyColorTermsapplyColorTerms:
424  msg = ("Only doing astrometry, but Colorterms are not applied for astrometry;"
425  "applyColorTerms=True will be ignored.")
426  lsst.log.warn(msg)
427 
428  def setDefaults(self):
429  # Use science source selector which can filter on extendedness, SNR, and whether blended
430  self.sourceSelectorsourceSelector.name = 'science'
431  # Use only stars because aperture fluxes of galaxies are biased and depend on seeing
432  self.sourceSelectorsourceSelector['science'].doUnresolved = True
433  # with dependable signal to noise ratio.
434  self.sourceSelectorsourceSelector['science'].doSignalToNoise = True
435  # Min SNR must be > 0 because jointcal cannot handle negative fluxes,
436  # and S/N > 10 to use sources that are not too faint, and thus better measured.
437  self.sourceSelectorsourceSelector['science'].signalToNoise.minimum = 10.
438  # Base SNR on CalibFlux because that is the flux jointcal that fits and must be positive
439  fluxField = f"slot_{self.sourceFluxType}Flux_instFlux"
440  self.sourceSelectorsourceSelector['science'].signalToNoise.fluxField = fluxField
441  self.sourceSelectorsourceSelector['science'].signalToNoise.errField = fluxField + "Err"
442  # Do not trust blended sources' aperture fluxes which also depend on seeing
443  self.sourceSelectorsourceSelector['science'].doIsolated = True
444  # Do not trust either flux or centroid measurements with flags,
445  # chosen from the usual QA flags for stars)
446  self.sourceSelectorsourceSelector['science'].doFlags = True
447  badFlags = ['base_PixelFlags_flag_edge', 'base_PixelFlags_flag_saturated',
448  'base_PixelFlags_flag_interpolatedCenter', 'base_SdssCentroid_flag',
449  'base_PsfFlux_flag', 'base_PixelFlags_flag_suspectCenter']
450  self.sourceSelectorsourceSelector['science'].flags.bad = badFlags
451 
452  # Default to Gaia-DR2 (with proper motions) for astrometry and
453  # PS1-DR1 for photometry, with a reasonable initial filterMap.
454  self.astrometryRefObjLoaderastrometryRefObjLoader.ref_dataset_name = "gaia_dr2_20200414"
455  self.astrometryRefObjLoaderastrometryRefObjLoader.requireProperMotion = True
456  self.astrometryRefObjLoaderastrometryRefObjLoader.anyFilterMapsToThis = 'phot_g_mean'
457  self.photometryRefObjLoaderphotometryRefObjLoader.ref_dataset_name = "ps1_pv3_3pi_20170110"
458 
459 
460 def writeModel(model, filename, log):
461  """Write model to outfile."""
462  with open(filename, "w") as file:
463  file.write(repr(model))
464  log.info("Wrote %s to file: %s", model, filename)
465 
466 
467 @dataclasses.dataclass
469  """The input data jointcal needs for each detector/visit."""
470  visit: int
471  """The visit identifier of this exposure."""
473  """The catalog derived from this exposure."""
474  visitInfo: lsst.afw.image.VisitInfo
475  """The VisitInfo of this exposure."""
477  """The detector of this exposure."""
478  photoCalib: lsst.afw.image.PhotoCalib
479  """The photometric calibration of this exposure."""
481  """The WCS of this exposure."""
482  bbox: lsst.geom.Box2I
483  """The bounding box of this exposure."""
485  """The filter of this exposure."""
486 
487 
488 class JointcalTask(pipeBase.PipelineTask, pipeBase.CmdLineTask):
489  """Astrometricly and photometricly calibrate across multiple visits of the
490  same field.
491 
492  Parameters
493  ----------
494  butler : `lsst.daf.persistence.Butler`
495  The butler is passed to the refObjLoader constructor in case it is
496  needed. Ignored if the refObjLoader argument provides a loader directly.
497  Used to initialize the astrometry and photometry refObjLoaders.
498  initInputs : `dict`, optional
499  Dictionary used to initialize PipelineTasks (empty for jointcal).
500  """
501 
502  ConfigClass = JointcalConfig
503  RunnerClass = JointcalRunner
504  _DefaultName = "jointcal"
505 
506  def __init__(self, butler=None, initInputs=None, **kwargs):
507  super().__init__(**kwargs)
508  self.makeSubtask("sourceSelector")
509  if self.config.doAstrometry:
510  if initInputs is None:
511  # gen3 middleware does refcat things internally (and will not have a butler here)
512  self.makeSubtask('astrometryRefObjLoader', butler=butler)
513  self.makeSubtask("astrometryReferenceSelector")
514  else:
515  self.astrometryRefObjLoaderastrometryRefObjLoader = None
516  if self.config.doPhotometry:
517  if initInputs is None:
518  # gen3 middleware does refcat things internally (and will not have a butler here)
519  self.makeSubtask('photometryRefObjLoader', butler=butler)
520  self.makeSubtask("photometryReferenceSelector")
521  else:
522  self.photometryRefObjLoaderphotometryRefObjLoader = None
523 
524  # To hold various computed metrics for use by tests
525  self.jobjob = Job.load_metrics_package(subset='jointcal')
526 
527  def runQuantum(self, butlerQC, inputRefs, outputRefs):
528  # We override runQuantum to set up the refObjLoaders and write the
529  # outputs to the correct refs.
530  inputs = butlerQC.get(inputRefs)
531  # We want the tract number for writing debug files
532  tract = butlerQC.quantum.dataId['tract']
533  if self.config.doAstrometry:
534  self.astrometryRefObjLoaderastrometryRefObjLoader = ReferenceObjectLoader(
535  dataIds=[ref.datasetRef.dataId
536  for ref in inputRefs.astrometryRefCat],
537  refCats=inputs.pop('astrometryRefCat'),
538  config=self.config.astrometryRefObjLoader,
539  log=self.log)
540  if self.config.doPhotometry:
541  self.photometryRefObjLoaderphotometryRefObjLoader = ReferenceObjectLoader(
542  dataIds=[ref.datasetRef.dataId
543  for ref in inputRefs.photometryRefCat],
544  refCats=inputs.pop('photometryRefCat'),
545  config=self.config.photometryRefObjLoader,
546  log=self.log)
547  outputs = self.runrun(**inputs, tract=tract)
548  if self.config.doAstrometry:
549  self._put_output_put_output(butlerQC, outputs.outputWcs, outputRefs.outputWcs,
550  inputs['inputCamera'], "setWcs")
551  if self.config.doPhotometry:
552  self._put_output_put_output(butlerQC, outputs.outputPhotoCalib, outputRefs.outputPhotoCalib,
553  inputs['inputCamera'], "setPhotoCalib")
554 
555  def _put_output(self, butlerQC, outputs, outputRefs, camera, setter):
556  """Persist the output datasets to their appropriate datarefs.
557 
558  Parameters
559  ----------
560  butlerQC : `ButlerQuantumContext`
561  A butler which is specialized to operate in the context of a
562  `lsst.daf.butler.Quantum`; This is the input to `runQuantum`.
563  outputs : `dict` [`tuple`, `lsst.afw.geom.SkyWcs`] or
564  `dict` [`tuple, `lsst.afw.image.PhotoCalib`]
565  The fitted objects to persist.
566  outputRefs : `list` [`OutputQuantizedConnection`]
567  The DatasetRefs to persist the data to.
568  camera : `lsst.afw.cameraGeom.Camera`
569  The camera for this instrument, to get detector ids from.
570  setter : `str`
571  The method to call on the ExposureCatalog to set each output.
572  """
574  schema.addField('visit', type='I', doc='Visit number')
575 
576  def new_catalog(visit, size):
577  """Return an catalog ready to be filled with appropriate output."""
578  catalog = lsst.afw.table.ExposureCatalog(schema)
579  catalog.resize(size)
580  catalog['visit'] = visit
581  metadata = lsst.daf.base.PropertyList()
582  metadata.add("COMMENT", "Catalog id is detector id, sorted.")
583  metadata.add("COMMENT", "Only detectors with data have entries.")
584  return catalog
585 
586  # count how many detectors have output for each visit
587  detectors_per_visit = collections.defaultdict(int)
588  for key in outputs:
589  # key is (visit, detector_id), and we only need visit here
590  detectors_per_visit[key[0]] += 1
591 
592  for ref in outputRefs:
593  visit = ref.dataId['visit']
594  catalog = new_catalog(visit, detectors_per_visit[visit])
595  # Iterate over every detector and skip the ones we don't have output for.
596  i = 0
597  for detector in camera:
598  detectorId = detector.getId()
599  key = (ref.dataId['visit'], detectorId)
600  if key not in outputs:
601  # skip detectors we don't have output for
602  self.log.debug("No %s output for detector %s in visit %s",
603  setter[3:], detectorId, visit)
604  continue
605 
606  catalog[i].setId(detectorId)
607  getattr(catalog[i], setter)(outputs[key])
608  i += 1
609 
610  catalog.sort() # ensure that the detectors are in sorted order, for fast lookups
611  butlerQC.put(catalog, ref)
612  self.log.info("Wrote %s detectors to %s", i, ref)
613 
614  def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None):
615  # Docstring inherited.
616 
617  # We take values out of the Parquet table, and put them in "flux_",
618  # and the config.sourceFluxType field is used during that extraction,
619  # so just use "flux" here.
620  sourceFluxField = "flux"
621  jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
622  associations = lsst.jointcal.Associations()
623  self.focalPlaneBBoxfocalPlaneBBox = inputCamera.getFpBBox()
624  oldWcsList, bands = self._load_data_load_data(inputSourceTableVisit,
625  inputVisitSummary,
626  associations,
627  jointcalControl)
628 
629  boundingCircle, center, radius, defaultFilter = self._prep_sky_prep_sky(associations, bands)
630  epoch = self._compute_proper_motion_epoch_compute_proper_motion_epoch(associations.getCcdImageList())
631 
632  if self.config.doAstrometry:
633  astrometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
634  name="astrometry",
635  refObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
636  referenceSelector=self.astrometryReferenceSelector,
637  fit_function=self._fit_astrometry_fit_astrometry,
638  tract=tract,
639  epoch=epoch)
640  astrometry_output = self._make_output_make_output(associations.getCcdImageList(),
641  astrometry.model,
642  "makeSkyWcs")
643  else:
644  astrometry_output = None
645 
646  if self.config.doPhotometry:
647  photometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultFilter, center, radius,
648  name="photometry",
649  refObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
650  referenceSelector=self.photometryReferenceSelector,
651  fit_function=self._fit_photometry_fit_photometry,
652  tract=tract,
653  epoch=epoch,
654  reject_bad_fluxes=True)
655  photometry_output = self._make_output_make_output(associations.getCcdImageList(),
656  photometry.model,
657  "toPhotoCalib")
658  else:
659  photometry_output = None
660 
661  return pipeBase.Struct(outputWcs=astrometry_output,
662  outputPhotoCalib=photometry_output,
663  job=self.jobjob,
664  astrometryRefObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
665  photometryRefObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader)
666 
667  def _make_schema_table(self):
668  """Return an afw SourceTable to use as a base for creating the
669  SourceCatalog to insert values from the dataFrame into.
670 
671  Returns
672  -------
673  table : `lsst.afw.table.SourceTable`
674  Table with schema and slots to use to make SourceCatalogs.
675  """
677  schema.addField("centroid_x", "D")
678  schema.addField("centroid_y", "D")
679  schema.addField("centroid_xErr", "F")
680  schema.addField("centroid_yErr", "F")
681  schema.addField("shape_xx", "D")
682  schema.addField("shape_yy", "D")
683  schema.addField("shape_xy", "D")
684  schema.addField("flux_instFlux", "D")
685  schema.addField("flux_instFluxErr", "D")
686  table = lsst.afw.table.SourceTable.make(schema)
687  table.defineCentroid("centroid")
688  table.defineShape("shape")
689  return table
690 
691  def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId):
692  """Return an afw SourceCatalog extracted from a visit-level dataframe,
693  limited to just one detector.
694 
695  Parameters
696  ----------
697  table : `lsst.afw.table.SourceTable`
698  Table factory to use to make the SourceCatalog that will be
699  populated with data from ``visitCatalog``.
700  visitCatalog : `pandas.DataFrame`
701  DataFrame to extract a detector catalog from.
702  detectorId : `int`
703  Numeric id of the detector to extract from ``visitCatalog``.
704 
705  Returns
706  -------
707  catalog : `lsst.afw.table.SourceCatalog`
708  Detector-level catalog extracted from ``visitCatalog``.
709  """
710  # map from dataFrame column to afw table column
711  mapping = {'sourceId': 'id',
712  'x': 'centroid_x',
713  'y': 'centroid_y',
714  'xErr': 'centroid_xErr',
715  'yErr': 'centroid_yErr',
716  'Ixx': 'shape_xx',
717  'Iyy': 'shape_yy',
718  'Ixy': 'shape_xy',
719  f'{self.config.sourceFluxType}_instFlux': 'flux_instFlux',
720  f'{self.config.sourceFluxType}_instFluxErr': 'flux_instFluxErr',
721  }
722  catalog = lsst.afw.table.SourceCatalog(table)
723  matched = visitCatalog['ccd'] == detectorId
724  catalog.resize(sum(matched))
725  view = visitCatalog.loc[matched]
726  for dfCol, afwCol in mapping.items():
727  catalog[afwCol] = view[dfCol]
728 
729  self.log.debug("%d sources selected in visit %d detector %d",
730  len(catalog),
731  view['visit'].iloc[0], # all visits in this catalog are the same, so take the first
732  detectorId)
733  return catalog
734 
735  def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations, jointcalControl):
736  """Read the data that jointcal needs to run. (Gen3 version)
737 
738  Modifies ``associations`` in-place with the loaded data.
739 
740  Parameters
741  ----------
742  inputSourceTableVisit : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
743  References to visit-level DataFrames to load the catalog data from.
744  inputVisitSummary : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
745  Visit-level exposure summary catalog with metadata.
746  associations : `lsst.jointcal.Associations`
747  Object to add the loaded data to by constructing new CcdImages.
748  jointcalControl : `jointcal.JointcalControl`
749  Control object for C++ associations management.
750 
751  Returns
752  -------
753  oldWcsList: `list` [`lsst.afw.geom.SkyWcs`]
754  The original WCS of the input data, to aid in writing tests.
755  bands : `list` [`str`]
756  The filter bands of each input dataset.
757  """
758  oldWcsList = []
759  bands = []
760  load_cat_prof_file = 'jointcal_load_data.prof' if self.config.detailedProfile else ''
761  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
762  table = self._make_schema_table_make_schema_table() # every detector catalog has the same layout
763  # No guarantee that the input is in the same order of visits, so we have to map one of them.
764  catalogMap = {ref.dataId['visit']: i for i, ref in enumerate(inputSourceTableVisit)}
765 
766  for visitSummaryRef in inputVisitSummary:
767  visitSummary = visitSummaryRef.get()
768  visitCatalog = inputSourceTableVisit[catalogMap[visitSummaryRef.dataId['visit']]].get()
769  selected = self.sourceSelector.run(visitCatalog)
770 
771  # Build a CcdImage for each detector in this visit.
772  detectors = {id: index for index, id in enumerate(visitSummary['id'])}
773  for id, index in detectors.items():
774  catalog = self._extract_detector_catalog_from_visit_catalog_extract_detector_catalog_from_visit_catalog(table, selected.sourceCat, id)
775  data = self._make_one_input_data_make_one_input_data(visitSummary[index], catalog)
776  result = self._build_ccdImage_build_ccdImage(data, associations, jointcalControl)
777  if result is not None:
778  oldWcsList.append(result.wcs)
779 
780  # A visit has only one band, so we can just use the first.
781  bands.append(visitSummary[0]['band'])
782  bands = collections.Counter(bands)
783 
784  return oldWcsList, bands
785 
786  def _make_one_input_data(self, visitRecord, catalog):
787  """Return a data structure for this detector+visit."""
788  return JointcalInputData(visit=visitRecord['visit'],
789  catalog=catalog,
790  visitInfo=visitRecord.getVisitInfo(),
791  detector=visitRecord.getDetector(),
792  photoCalib=visitRecord.getPhotoCalib(),
793  wcs=visitRecord.getWcs(),
794  bbox=visitRecord.getBBox(),
795  # ExposureRecord doesn't have a FilterLabel yet,
796  # so we have to make one.
797  filter=lsst.afw.image.FilterLabel(band=visitRecord['band'],
798  physical=visitRecord['physical_filter']))
799 
800  # We don't currently need to persist the metadata.
801  # If we do in the future, we will have to add appropriate dataset templates
802  # to each obs package (the metadata template should look like `jointcal_wcs`).
803  def _getMetadataName(self):
804  return None
805 
806  @classmethod
807  def _makeArgumentParser(cls):
808  """Create an argument parser"""
809  parser = pipeBase.ArgumentParser(name=cls._DefaultName_DefaultName)
810  parser.add_id_argument("--id", "calexp", help="data ID, e.g. --id visit=6789 ccd=0..9",
811  ContainerClass=PerTractCcdDataIdContainer)
812  return parser
813 
814  def _build_ccdImage(self, data, associations, jointcalControl):
815  """
816  Extract the necessary things from this dataRef to add a new ccdImage.
817 
818  Parameters
819  ----------
820  data : `JointcalInputData`
821  The loaded input data.
822  associations : `lsst.jointcal.Associations`
823  Object to add the info to, to construct a new CcdImage
824  jointcalControl : `jointcal.JointcalControl`
825  Control object for associations management
826 
827  Returns
828  ------
829  namedtuple or `None`
830  ``wcs``
831  The TAN WCS of this image, read from the calexp
832  (`lsst.afw.geom.SkyWcs`).
833  ``key``
834  A key to identify this dataRef by its visit and ccd ids
835  (`namedtuple`).
836  ``band``
837  This calexp's filter band (`str`) (used to e.g. load refcats)
838  `None`
839  if there are no sources in the loaded catalog.
840  """
841  if len(data.catalog) == 0:
842  self.log.warn("No sources selected in visit %s ccd %s", data.visit, data.detector.getId())
843  return None
844 
845  associations.createCcdImage(data.catalog,
846  data.wcs,
847  data.visitInfo,
848  data.bbox,
849  data.filter.physicalLabel,
850  data.photoCalib,
851  data.detector,
852  data.visit,
853  data.detector.getId(),
854  jointcalControl)
855 
856  Result = collections.namedtuple('Result_from_build_CcdImage', ('wcs', 'key', 'band'))
857  Key = collections.namedtuple('Key', ('visit', 'ccd'))
858  return Result(data.wcs, Key(data.visit, data.detector.getId()), data.filter.bandLabel)
859 
860  def _readDataId(self, butler, dataId):
861  """Read all of the data for one dataId from the butler. (gen2 version)"""
862  # Not all instruments have `visit` in their dataIds.
863  if "visit" in dataId.keys():
864  visit = dataId["visit"]
865  else:
866  visit = butler.getButler().queryMetadata("calexp", ("visit"), butler.dataId)[0]
867  detector = butler.get('calexp_detector', dataId=dataId)
868 
869  catalog = butler.get('src',
870  flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS,
871  dataId=dataId)
872  goodSrc = self.sourceSelector.run(catalog)
873  self.log.debug("%d sources selected in visit %d detector %d",
874  len(goodSrc.sourceCat),
875  visit,
876  detector.getId())
877  return JointcalInputData(visit=visit,
878  catalog=goodSrc.sourceCat,
879  visitInfo=butler.get('calexp_visitInfo', dataId=dataId),
880  detector=detector,
881  photoCalib=butler.get('calexp_photoCalib', dataId=dataId),
882  wcs=butler.get('calexp_wcs', dataId=dataId),
883  bbox=butler.get('calexp_bbox', dataId=dataId),
884  filter=butler.get('calexp_filterLabel', dataId=dataId))
885 
886  def loadData(self, dataRefs, associations, jointcalControl):
887  """Read the data that jointcal needs to run. (Gen2 version)"""
888  visit_ccd_to_dataRef = {}
889  oldWcsList = []
890  bands = []
891  load_cat_prof_file = 'jointcal_loadData.prof' if self.config.detailedProfile else ''
892  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
893  # Need the bounding-box of the focal plane (the same for all visits) for photometry visit models
894  camera = dataRefs[0].get('camera', immediate=True)
895  self.focalPlaneBBoxfocalPlaneBBox = camera.getFpBBox()
896  for dataRef in dataRefs:
897  data = self._readDataId_readDataId(dataRef.getButler(), dataRef.dataId)
898  result = self._build_ccdImage_build_ccdImage(data, associations, jointcalControl)
899  if result is None:
900  continue
901  oldWcsList.append(result.wcs)
902  visit_ccd_to_dataRef[result.key] = dataRef
903  bands.append(result.band)
904  if len(bands) == 0:
905  raise RuntimeError("No data to process: did source selector remove all sources?")
906  bands = collections.Counter(bands)
907 
908  return oldWcsList, bands, visit_ccd_to_dataRef
909 
910  def _getDebugPath(self, filename):
911  """Constructs a path to filename using the configured debug path.
912  """
913  return os.path.join(self.config.debugOutputPath, filename)
914 
915  def _prep_sky(self, associations, bands):
916  """Prepare on-sky and other data that must be computed after data has
917  been read.
918  """
919  associations.computeCommonTangentPoint()
920 
921  boundingCircle = associations.computeBoundingCircle()
922  center = lsst.geom.SpherePoint(boundingCircle.getCenter())
923  radius = lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
924 
925  self.log.info(f"Data has center={center} with radius={radius.asDegrees()} degrees.")
926 
927  # Determine a default filter band associated with the catalog. See DM-9093
928  defaultBand = bands.most_common(1)[0][0]
929  self.log.debug("Using '%s' filter band for reference flux", defaultBand)
930 
931  return boundingCircle, center, radius, defaultBand
932 
933  @pipeBase.timeMethod
934  def runDataRef(self, dataRefs):
935  """
936  Jointly calibrate the astrometry and photometry across a set of images.
937 
938  NOTE: this is for gen2 middleware only.
939 
940  Parameters
941  ----------
942  dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef`
943  List of data references to the exposures to be fit.
944 
945  Returns
946  -------
947  result : `lsst.pipe.base.Struct`
948  Struct of metadata from the fit, containing:
949 
950  ``dataRefs``
951  The provided data references that were fit (with updated WCSs)
952  ``oldWcsList``
953  The original WCS from each dataRef
954  ``metrics``
955  Dictionary of internally-computed metrics for testing/validation.
956  """
957  if len(dataRefs) == 0:
958  raise ValueError('Need a non-empty list of data references!')
959 
960  exitStatus = 0 # exit status for shell
961 
962  sourceFluxField = "slot_%sFlux" % (self.config.sourceFluxType,)
963  jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
964  associations = lsst.jointcal.Associations()
965 
966  oldWcsList, bands, visit_ccd_to_dataRef = self.loadDataloadData(dataRefs,
967  associations,
968  jointcalControl)
969 
970  boundingCircle, center, radius, defaultBand = self._prep_sky_prep_sky(associations, bands)
971  epoch = self._compute_proper_motion_epoch_compute_proper_motion_epoch(associations.getCcdImageList())
972 
973  tract = dataRefs[0].dataId['tract']
974 
975  if self.config.doAstrometry:
976  astrometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultBand, center, radius,
977  name="astrometry",
978  refObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
979  referenceSelector=self.astrometryReferenceSelector,
980  fit_function=self._fit_astrometry_fit_astrometry,
981  tract=tract,
982  epoch=epoch)
983  self._write_astrometry_results_write_astrometry_results(associations, astrometry.model, visit_ccd_to_dataRef)
984  else:
985  astrometry = Astrometry(None, None, None)
986 
987  if self.config.doPhotometry:
988  photometry = self._do_load_refcat_and_fit_do_load_refcat_and_fit(associations, defaultBand, center, radius,
989  name="photometry",
990  refObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
991  referenceSelector=self.photometryReferenceSelector,
992  fit_function=self._fit_photometry_fit_photometry,
993  tract=tract,
994  epoch=epoch,
995  reject_bad_fluxes=True)
996  self._write_photometry_results_write_photometry_results(associations, photometry.model, visit_ccd_to_dataRef)
997  else:
998  photometry = Photometry(None, None)
999 
1000  return pipeBase.Struct(dataRefs=dataRefs,
1001  oldWcsList=oldWcsList,
1002  job=self.jobjob,
1003  astrometryRefObjLoader=self.astrometryRefObjLoaderastrometryRefObjLoader,
1004  photometryRefObjLoader=self.photometryRefObjLoaderphotometryRefObjLoader,
1005  defaultBand=defaultBand,
1006  epoch=epoch,
1007  exitStatus=exitStatus)
1008 
1009  def _get_refcat_coordinate_error_override(self, refCat, name):
1010  """Check whether we should override the refcat coordinate errors, and
1011  return the overridden error if necessary.
1012 
1013  Parameters
1014  ----------
1015  refCat : `lsst.afw.table.SimpleCatalog`
1016  The reference catalog to check for a ``coord_raErr`` field.
1017  name : `str`
1018  Whether we are doing "astrometry" or "photometry".
1019 
1020  Returns
1021  -------
1022  refCoordErr : `float`
1023  The refcat coordinate error to use, or NaN if we are not overriding
1024  those fields.
1025 
1026  Raises
1027  ------
1028  lsst.pex.config.FieldValidationError
1029  Raised if the refcat does not contain coordinate errors and
1030  ``config.astrometryReferenceErr`` is not set.
1031  """
1032  # This value doesn't matter for photometry, so just set something to
1033  # keep old refcats from causing problems.
1034  if name.lower() == "photometry":
1035  if 'coord_raErr' not in refCat.schema:
1036  return 100
1037  else:
1038  return float('nan')
1039 
1040  if self.config.astrometryReferenceErr is None and 'coord_raErr' not in refCat.schema:
1041  msg = ("Reference catalog does not contain coordinate errors, "
1042  "and config.astrometryReferenceErr not supplied.")
1043  raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
1044  self.config,
1045  msg)
1046 
1047  if self.config.astrometryReferenceErr is not None and 'coord_raErr' in refCat.schema:
1048  self.log.warn("Overriding reference catalog coordinate errors with %f/coordinate [mas]",
1049  self.config.astrometryReferenceErr)
1050 
1051  if self.config.astrometryReferenceErr is None:
1052  return float('nan')
1053  else:
1054  return self.config.astrometryReferenceErr
1055 
1056  def _compute_proper_motion_epoch(self, ccdImageList):
1057  """Return the proper motion correction epoch of the provided images.
1058 
1059  Parameters
1060  ----------
1061  ccdImageList : `list` [`lsst.jointcal.CcdImage`]
1062  The images to compute the appropriate epoch for.
1063 
1064  Returns
1065  -------
1066  epoch : `astropy.time.Time`
1067  The date to use for proper motion corrections.
1068  """
1069  mjds = [ccdImage.getMjd() for ccdImage in ccdImageList]
1070  return astropy.time.Time(np.mean(mjds), format='mjd', scale="tai")
1071 
1072  def _do_load_refcat_and_fit(self, associations, defaultBand, center, radius,
1073  tract="", match_cut=3.0,
1074  reject_bad_fluxes=False, *,
1075  name="", refObjLoader=None, referenceSelector=None,
1076  fit_function=None, epoch=None):
1077  """Load reference catalog, perform the fit, and return the result.
1078 
1079  Parameters
1080  ----------
1081  associations : `lsst.jointcal.Associations`
1082  The star/reference star associations to fit.
1083  defaultBand : `str`
1084  filter to load from reference catalog.
1085  center : `lsst.geom.SpherePoint`
1086  ICRS center of field to load from reference catalog.
1087  radius : `lsst.geom.Angle`
1088  On-sky radius to load from reference catalog.
1089  name : `str`
1090  Name of thing being fit: "astrometry" or "photometry".
1091  refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
1092  Reference object loader to use to load a reference catalog.
1093  referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1094  Selector to use to pick objects from the loaded reference catalog.
1095  fit_function : callable
1096  Function to call to perform fit (takes Associations object).
1097  tract : `str`, optional
1098  Name of tract currently being fit.
1099  match_cut : `float`, optional
1100  Radius in arcseconds to find cross-catalog matches to during
1101  associations.associateCatalogs.
1102  reject_bad_fluxes : `bool`, optional
1103  Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
1104  epoch : `astropy.time.Time`, optional
1105  Epoch to which to correct refcat proper motion and parallax,
1106  or `None` to not apply such corrections.
1107 
1108  Returns
1109  -------
1110  result : `Photometry` or `Astrometry`
1111  Result of `fit_function()`
1112  """
1113  self.log.info("====== Now processing %s...", name)
1114  # TODO: this should not print "trying to invert a singular transformation:"
1115  # if it does that, something's not right about the WCS...
1116  associations.associateCatalogs(match_cut)
1117  add_measurement(self.jobjob, 'jointcal.associated_%s_fittedStars' % name,
1118  associations.fittedStarListSize())
1119 
1120  applyColorterms = False if name.lower() == "astrometry" else self.config.applyColorTerms
1121  refCat, fluxField = self._load_reference_catalog_load_reference_catalog(refObjLoader, referenceSelector,
1122  center, radius, defaultBand,
1123  applyColorterms=applyColorterms,
1124  epoch=epoch)
1125  refCoordErr = self._get_refcat_coordinate_error_override_get_refcat_coordinate_error_override(refCat, name)
1126 
1127  associations.collectRefStars(refCat,
1128  self.config.matchCut*lsst.geom.arcseconds,
1129  fluxField,
1130  refCoordinateErr=refCoordErr,
1131  rejectBadFluxes=reject_bad_fluxes)
1132  add_measurement(self.jobjob, 'jointcal.collected_%s_refStars' % name,
1133  associations.refStarListSize())
1134 
1135  associations.prepareFittedStars(self.config.minMeasurements)
1136 
1137  self._check_star_lists_check_star_lists(associations, name)
1138  add_measurement(self.jobjob, 'jointcal.selected_%s_refStars' % name,
1139  associations.nFittedStarsWithAssociatedRefStar())
1140  add_measurement(self.jobjob, 'jointcal.selected_%s_fittedStars' % name,
1141  associations.fittedStarListSize())
1142  add_measurement(self.jobjob, 'jointcal.selected_%s_ccdImages' % name,
1143  associations.nCcdImagesValidForFit())
1144 
1145  load_cat_prof_file = 'jointcal_fit_%s.prof'%name if self.config.detailedProfile else ''
1146  dataName = "{}_{}".format(tract, defaultBand)
1147  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
1148  result = fit_function(associations, dataName)
1149  # TODO DM-12446: turn this into a "butler save" somehow.
1150  # Save reference and measurement chi2 contributions for this data
1151  if self.config.writeChi2FilesInitialFinal:
1152  baseName = self._getDebugPath_getDebugPath(f"{name}_final_chi2-{dataName}")
1153  result.fit.saveChi2Contributions(baseName+"{type}")
1154  self.log.info("Wrote chi2 contributions files: %s", baseName)
1155 
1156  return result
1157 
1158  def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName,
1159  applyColorterms=False, epoch=None):
1160  """Load the necessary reference catalog sources, convert fluxes to
1161  correct units, and apply color term corrections if requested.
1162 
1163  Parameters
1164  ----------
1165  refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
1166  The reference catalog loader to use to get the data.
1167  referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1168  Source selector to apply to loaded reference catalog.
1169  center : `lsst.geom.SpherePoint`
1170  The center around which to load sources.
1171  radius : `lsst.geom.Angle`
1172  The radius around ``center`` to load sources in.
1173  filterName : `str`
1174  The name of the camera filter to load fluxes for.
1175  applyColorterms : `bool`
1176  Apply colorterm corrections to the refcat for ``filterName``?
1177  epoch : `astropy.time.Time`, optional
1178  Epoch to which to correct refcat proper motion and parallax,
1179  or `None` to not apply such corrections.
1180 
1181  Returns
1182  -------
1183  refCat : `lsst.afw.table.SimpleCatalog`
1184  The loaded reference catalog.
1185  fluxField : `str`
1186  The name of the reference catalog flux field appropriate for ``filterName``.
1187  """
1188  skyCircle = refObjLoader.loadSkyCircle(center,
1189  radius,
1190  filterName,
1191  epoch=epoch)
1192 
1193  selected = referenceSelector.run(skyCircle.refCat)
1194  # Need memory contiguity to get reference filters as a vector.
1195  if not selected.sourceCat.isContiguous():
1196  refCat = selected.sourceCat.copy(deep=True)
1197  else:
1198  refCat = selected.sourceCat
1199 
1200  if applyColorterms:
1201  refCatName = refObjLoader.ref_dataset_name
1202  self.log.info("Applying color terms for filterName=%r reference catalog=%s",
1203  filterName, refCatName)
1204  colorterm = self.config.colorterms.getColorterm(
1205  filterName=filterName, photoCatName=refCatName, doRaise=True)
1206 
1207  refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat, filterName)
1208  refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
1209  # TODO: I didn't want to use this, but I'll deal with it in DM-16903
1210  refCat[skyCircle.fluxField+'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
1211 
1212  return refCat, skyCircle.fluxField
1213 
1214  def _check_star_lists(self, associations, name):
1215  # TODO: these should be len(blah), but we need this properly wrapped first.
1216  if associations.nCcdImagesValidForFit() == 0:
1217  raise RuntimeError('No images in the ccdImageList!')
1218  if associations.fittedStarListSize() == 0:
1219  raise RuntimeError('No stars in the {} fittedStarList!'.format(name))
1220  if associations.refStarListSize() == 0:
1221  raise RuntimeError('No stars in the {} reference star list!'.format(name))
1222 
1223  def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
1224  """Compute chi2, log it, validate the model, and return chi2.
1225 
1226  Parameters
1227  ----------
1228  associations : `lsst.jointcal.Associations`
1229  The star/reference star associations to fit.
1230  fit : `lsst.jointcal.FitterBase`
1231  The fitter to use for minimization.
1232  model : `lsst.jointcal.Model`
1233  The model being fit.
1234  chi2Label : `str`
1235  Label to describe the chi2 (e.g. "Initialized", "Final").
1236  writeChi2Name : `str`, optional
1237  Filename prefix to write the chi2 contributions to.
1238  Do not supply an extension: an appropriate one will be added.
1239 
1240  Returns
1241  -------
1242  chi2: `lsst.jointcal.Chi2Accumulator`
1243  The chi2 object for the current fitter and model.
1244 
1245  Raises
1246  ------
1247  FloatingPointError
1248  Raised if chi2 is infinite or NaN.
1249  ValueError
1250  Raised if the model is not valid.
1251  """
1252  if writeChi2Name is not None:
1253  fullpath = self._getDebugPath_getDebugPath(writeChi2Name)
1254  fit.saveChi2Contributions(fullpath+"{type}")
1255  self.log.info("Wrote chi2 contributions files: %s", fullpath)
1256 
1257  chi2 = fit.computeChi2()
1258  self.log.info("%s %s", chi2Label, chi2)
1259  self._check_stars_check_stars(associations)
1260  if not np.isfinite(chi2.chi2):
1261  raise FloatingPointError(f'{chi2Label} chi2 is invalid: {chi2}')
1262  if not model.validate(associations.getCcdImageList(), chi2.ndof):
1263  raise ValueError("Model is not valid: check log messages for warnings.")
1264  return chi2
1265 
1266  def _fit_photometry(self, associations, dataName=None):
1267  """
1268  Fit the photometric data.
1269 
1270  Parameters
1271  ----------
1272  associations : `lsst.jointcal.Associations`
1273  The star/reference star associations to fit.
1274  dataName : `str`
1275  Name of the data being processed (e.g. "1234_HSC-Y"), for
1276  identifying debugging files.
1277 
1278  Returns
1279  -------
1280  fit_result : `namedtuple`
1281  fit : `lsst.jointcal.PhotometryFit`
1282  The photometric fitter used to perform the fit.
1283  model : `lsst.jointcal.PhotometryModel`
1284  The photometric model that was fit.
1285  """
1286  self.log.info("=== Starting photometric fitting...")
1287 
1288  # TODO: should use pex.config.RegistryField here (see DM-9195)
1289  if self.config.photometryModel == "constrainedFlux":
1290  model = lsst.jointcal.ConstrainedFluxModel(associations.getCcdImageList(),
1291  self.focalPlaneBBoxfocalPlaneBBox,
1292  visitOrder=self.config.photometryVisitOrder,
1293  errorPedestal=self.config.photometryErrorPedestal)
1294  # potentially nonlinear problem, so we may need a line search to converge.
1295  doLineSearch = self.config.allowLineSearch
1296  elif self.config.photometryModel == "constrainedMagnitude":
1297  model = lsst.jointcal.ConstrainedMagnitudeModel(associations.getCcdImageList(),
1298  self.focalPlaneBBoxfocalPlaneBBox,
1299  visitOrder=self.config.photometryVisitOrder,
1300  errorPedestal=self.config.photometryErrorPedestal)
1301  # potentially nonlinear problem, so we may need a line search to converge.
1302  doLineSearch = self.config.allowLineSearch
1303  elif self.config.photometryModel == "simpleFlux":
1304  model = lsst.jointcal.SimpleFluxModel(associations.getCcdImageList(),
1305  errorPedestal=self.config.photometryErrorPedestal)
1306  doLineSearch = False # purely linear in model parameters, so no line search needed
1307  elif self.config.photometryModel == "simpleMagnitude":
1308  model = lsst.jointcal.SimpleMagnitudeModel(associations.getCcdImageList(),
1309  errorPedestal=self.config.photometryErrorPedestal)
1310  doLineSearch = False # purely linear in model parameters, so no line search needed
1311 
1312  fit = lsst.jointcal.PhotometryFit(associations, model)
1313  # TODO DM-12446: turn this into a "butler save" somehow.
1314  # Save reference and measurement chi2 contributions for this data
1315  if self.config.writeChi2FilesInitialFinal:
1316  baseName = f"photometry_initial_chi2-{dataName}"
1317  else:
1318  baseName = None
1319  if self.config.writeInitialModel:
1320  fullpath = self._getDebugPath_getDebugPath("initialPhotometryModel.txt")
1321  writeModel(model, fullpath, self.log)
1322  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialized", writeChi2Name=baseName)
1323 
1324  def getChi2Name(whatToFit):
1325  if self.config.writeChi2FilesOuterLoop:
1326  return f"photometry_init-%s_chi2-{dataName}" % whatToFit
1327  else:
1328  return None
1329 
1330  # The constrained model needs the visit transform fit first; the chip
1331  # transform is initialized from the singleFrame PhotoCalib, so it's close.
1332  dumpMatrixFile = self._getDebugPath_getDebugPath("photometry_preinit") if self.config.writeInitMatrix else ""
1333  if self.config.photometryModel.startswith("constrained"):
1334  # no line search: should be purely (or nearly) linear,
1335  # and we want a large step size to initialize with.
1336  fit.minimize("ModelVisit", dumpMatrixFile=dumpMatrixFile)
1337  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize ModelVisit",
1338  writeChi2Name=getChi2Name("ModelVisit"))
1339  dumpMatrixFile = "" # so we don't redo the output on the next step
1340 
1341  fit.minimize("Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
1342  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Model",
1343  writeChi2Name=getChi2Name("Model"))
1344 
1345  fit.minimize("Fluxes") # no line search: always purely linear.
1346  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Fluxes",
1347  writeChi2Name=getChi2Name("Fluxes"))
1348 
1349  fit.minimize("Model Fluxes", doLineSearch=doLineSearch)
1350  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize ModelFluxes",
1351  writeChi2Name=getChi2Name("ModelFluxes"))
1352 
1353  model.freezeErrorTransform()
1354  self.log.debug("Photometry error scales are frozen.")
1355 
1356  chi2 = self._iterate_fit_iterate_fit(associations,
1357  fit,
1358  self.config.maxPhotometrySteps,
1359  "photometry",
1360  "Model Fluxes",
1361  doRankUpdate=self.config.photometryDoRankUpdate,
1362  doLineSearch=doLineSearch,
1363  dataName=dataName)
1364 
1365  add_measurement(self.jobjob, 'jointcal.photometry_final_chi2', chi2.chi2)
1366  add_measurement(self.jobjob, 'jointcal.photometry_final_ndof', chi2.ndof)
1367  return Photometry(fit, model)
1368 
1369  def _fit_astrometry(self, associations, dataName=None):
1370  """
1371  Fit the astrometric data.
1372 
1373  Parameters
1374  ----------
1375  associations : `lsst.jointcal.Associations`
1376  The star/reference star associations to fit.
1377  dataName : `str`
1378  Name of the data being processed (e.g. "1234_HSC-Y"), for
1379  identifying debugging files.
1380 
1381  Returns
1382  -------
1383  fit_result : `namedtuple`
1384  fit : `lsst.jointcal.AstrometryFit`
1385  The astrometric fitter used to perform the fit.
1386  model : `lsst.jointcal.AstrometryModel`
1387  The astrometric model that was fit.
1388  sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1389  The model for the sky to tangent plane projection that was used in the fit.
1390  """
1391 
1392  self.log.info("=== Starting astrometric fitting...")
1393 
1394  associations.deprojectFittedStars()
1395 
1396  # NOTE: need to return sky_to_tan_projection so that it doesn't get garbage collected.
1397  # TODO: could we package sky_to_tan_projection and model together so we don't have to manage
1398  # them so carefully?
1399  sky_to_tan_projection = lsst.jointcal.OneTPPerVisitHandler(associations.getCcdImageList())
1400 
1401  if self.config.astrometryModel == "constrained":
1402  model = lsst.jointcal.ConstrainedAstrometryModel(associations.getCcdImageList(),
1403  sky_to_tan_projection,
1404  chipOrder=self.config.astrometryChipOrder,
1405  visitOrder=self.config.astrometryVisitOrder)
1406  elif self.config.astrometryModel == "simple":
1407  model = lsst.jointcal.SimpleAstrometryModel(associations.getCcdImageList(),
1408  sky_to_tan_projection,
1409  self.config.useInputWcs,
1410  nNotFit=0,
1411  order=self.config.astrometrySimpleOrder)
1412 
1413  fit = lsst.jointcal.AstrometryFit(associations, model, self.config.positionErrorPedestal)
1414  # TODO DM-12446: turn this into a "butler save" somehow.
1415  # Save reference and measurement chi2 contributions for this data
1416  if self.config.writeChi2FilesInitialFinal:
1417  baseName = f"astrometry_initial_chi2-{dataName}"
1418  else:
1419  baseName = None
1420  if self.config.writeInitialModel:
1421  fullpath = self._getDebugPath_getDebugPath("initialAstrometryModel.txt")
1422  writeModel(model, fullpath, self.log)
1423  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initial", writeChi2Name=baseName)
1424 
1425  def getChi2Name(whatToFit):
1426  if self.config.writeChi2FilesOuterLoop:
1427  return f"astrometry_init-%s_chi2-{dataName}" % whatToFit
1428  else:
1429  return None
1430 
1431  dumpMatrixFile = self._getDebugPath_getDebugPath("astrometry_preinit") if self.config.writeInitMatrix else ""
1432  # The constrained model needs the visit transform fit first; the chip
1433  # transform is initialized from the detector's cameraGeom, so it's close.
1434  if self.config.astrometryModel == "constrained":
1435  fit.minimize("DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1436  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize DistortionsVisit",
1437  writeChi2Name=getChi2Name("DistortionsVisit"))
1438  dumpMatrixFile = "" # so we don't redo the output on the next step
1439 
1440  fit.minimize("Distortions", dumpMatrixFile=dumpMatrixFile)
1441  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Distortions",
1442  writeChi2Name=getChi2Name("Distortions"))
1443 
1444  fit.minimize("Positions")
1445  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize Positions",
1446  writeChi2Name=getChi2Name("Positions"))
1447 
1448  fit.minimize("Distortions Positions")
1449  self._logChi2AndValidate_logChi2AndValidate(associations, fit, model, "Initialize DistortionsPositions",
1450  writeChi2Name=getChi2Name("DistortionsPositions"))
1451 
1452  chi2 = self._iterate_fit_iterate_fit(associations,
1453  fit,
1454  self.config.maxAstrometrySteps,
1455  "astrometry",
1456  "Distortions Positions",
1457  doRankUpdate=self.config.astrometryDoRankUpdate,
1458  dataName=dataName)
1459 
1460  add_measurement(self.jobjob, 'jointcal.astrometry_final_chi2', chi2.chi2)
1461  add_measurement(self.jobjob, 'jointcal.astrometry_final_ndof', chi2.ndof)
1462 
1463  return Astrometry(fit, model, sky_to_tan_projection)
1464 
1465  def _check_stars(self, associations):
1466  """Count measured and reference stars per ccd and warn/log them."""
1467  for ccdImage in associations.getCcdImageList():
1468  nMeasuredStars, nRefStars = ccdImage.countStars()
1469  self.log.debug("ccdImage %s has %s measured and %s reference stars",
1470  ccdImage.getName(), nMeasuredStars, nRefStars)
1471  if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1472  self.log.warn("ccdImage %s has only %s measuredStars (desired %s)",
1473  ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1474  if nRefStars < self.config.minRefStarsPerCcd:
1475  self.log.warn("ccdImage %s has only %s RefStars (desired %s)",
1476  ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1477 
1478  def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1479  dataName="",
1480  doRankUpdate=True,
1481  doLineSearch=False):
1482  """Run fitter.minimize up to max_steps times, returning the final chi2.
1483 
1484  Parameters
1485  ----------
1486  associations : `lsst.jointcal.Associations`
1487  The star/reference star associations to fit.
1488  fitter : `lsst.jointcal.FitterBase`
1489  The fitter to use for minimization.
1490  max_steps : `int`
1491  Maximum number of steps to run outlier rejection before declaring
1492  convergence failure.
1493  name : {'photometry' or 'astrometry'}
1494  What type of data are we fitting (for logs and debugging files).
1495  whatToFit : `str`
1496  Passed to ``fitter.minimize()`` to define the parameters to fit.
1497  dataName : `str`, optional
1498  Descriptive name for this dataset (e.g. tract and filter),
1499  for debugging.
1500  doRankUpdate : `bool`, optional
1501  Do an Eigen rank update during minimization, or recompute the full
1502  matrix and gradient?
1503  doLineSearch : `bool`, optional
1504  Do a line search for the optimum step during minimization?
1505 
1506  Returns
1507  -------
1508  chi2: `lsst.jointcal.Chi2Statistic`
1509  The final chi2 after the fit converges, or is forced to end.
1510 
1511  Raises
1512  ------
1513  FloatingPointError
1514  Raised if the fitter fails with a non-finite value.
1515  RuntimeError
1516  Raised if the fitter fails for some other reason;
1517  log messages will provide further details.
1518  """
1519  dumpMatrixFile = self._getDebugPath_getDebugPath(f"{name}_postinit") if self.config.writeInitMatrix else ""
1520  oldChi2 = lsst.jointcal.Chi2Statistic()
1521  oldChi2.chi2 = float("inf")
1522  for i in range(max_steps):
1523  if self.config.writeChi2FilesOuterLoop:
1524  writeChi2Name = f"{name}_iterate_{i}_chi2-{dataName}"
1525  else:
1526  writeChi2Name = None
1527  result = fitter.minimize(whatToFit,
1528  self.config.outlierRejectSigma,
1529  doRankUpdate=doRankUpdate,
1530  doLineSearch=doLineSearch,
1531  dumpMatrixFile=dumpMatrixFile)
1532  dumpMatrixFile = "" # clear it so we don't write the matrix again.
1533  chi2 = self._logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(),
1534  f"Fit iteration {i}", writeChi2Name=writeChi2Name)
1535 
1536  if result == MinimizeResult.Converged:
1537  if doRankUpdate:
1538  self.log.debug("fit has converged - no more outliers - redo minimization "
1539  "one more time in case we have lost accuracy in rank update.")
1540  # Redo minimization one more time in case we have lost accuracy in rank update
1541  result = fitter.minimize(whatToFit, self.config.outlierRejectSigma)
1542  chi2 = self._logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(), "Fit completed")
1543 
1544  # log a message for a large final chi2, TODO: DM-15247 for something better
1545  if chi2.chi2/chi2.ndof >= 4.0:
1546  self.log.error("Potentially bad fit: High chi-squared/ndof.")
1547 
1548  break
1549  elif result == MinimizeResult.Chi2Increased:
1550  self.log.warn("Still some outliers remaining but chi2 increased - retry")
1551  # Check whether the increase was large enough to cause trouble.
1552  chi2Ratio = chi2.chi2 / oldChi2.chi2
1553  if chi2Ratio > 1.5:
1554  self.log.warn('Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1555  chi2.chi2, oldChi2.chi2, chi2Ratio)
1556  # Based on a variety of HSC jointcal logs (see DM-25779), it
1557  # appears that chi2 increases more than a factor of ~2 always
1558  # result in the fit diverging rapidly and ending at chi2 > 1e10.
1559  # Using 10 as the "failure" threshold gives some room between
1560  # leaving a warning and bailing early.
1561  if chi2Ratio > 10:
1562  msg = ("Large chi2 increase between steps: fit likely cannot converge."
1563  " Try setting one or more of the `writeChi2*` config fields and looking"
1564  " at how individual star chi2-values evolve during the fit.")
1565  raise RuntimeError(msg)
1566  oldChi2 = chi2
1567  elif result == MinimizeResult.NonFinite:
1568  filename = self._getDebugPath_getDebugPath("{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1569  # TODO DM-12446: turn this into a "butler save" somehow.
1570  fitter.saveChi2Contributions(filename+"{type}")
1571  msg = "Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1572  raise FloatingPointError(msg.format(filename))
1573  elif result == MinimizeResult.Failed:
1574  raise RuntimeError("Chi2 minimization failure, cannot complete fit.")
1575  else:
1576  raise RuntimeError("Unxepected return code from minimize().")
1577  else:
1578  self.log.error("%s failed to converge after %d steps"%(name, max_steps))
1579 
1580  return chi2
1581 
1582  def _make_output(self, ccdImageList, model, func):
1583  """Return the internal jointcal models converted to the afw
1584  structures that will be saved to disk.
1585 
1586  Parameters
1587  ----------
1588  ccdImageList : `lsst.jointcal.CcdImageList`
1589  The list of CcdImages to get the output for.
1590  model : `lsst.jointcal.AstrometryModel` or `lsst.jointcal.PhotometryModel`
1591  The internal jointcal model to convert for each `lsst.jointcal.CcdImage`.
1592  func : `str`
1593  The name of the function to call on ``model`` to get the converted
1594  structure. Must accept an `lsst.jointcal.CcdImage`.
1595 
1596  Returns
1597  -------
1598  output : `dict` [`tuple`, `lsst.jointcal.AstrometryModel`] or
1599  `dict` [`tuple`, `lsst.jointcal.PhotometryModel`]
1600  The data to be saved, keyed on (visit, detector).
1601  """
1602  output = {}
1603  for ccdImage in ccdImageList:
1604  ccd = ccdImage.ccdId
1605  visit = ccdImage.visit
1606  self.log.debug("%s for visit: %d, ccd: %d", func, visit, ccd)
1607  output[(visit, ccd)] = getattr(model, func)(ccdImage)
1608  return output
1609 
1610  def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1611  """
1612  Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1613 
1614  Parameters
1615  ----------
1616  associations : `lsst.jointcal.Associations`
1617  The star/reference star associations to fit.
1618  model : `lsst.jointcal.AstrometryModel`
1619  The astrometric model that was fit.
1620  visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1621  Dict of ccdImage identifiers to dataRefs that were fit.
1622  """
1623  ccdImageList = associations.getCcdImageList()
1624  output = self._make_output_make_output(ccdImageList, model, "makeSkyWcs")
1625  for key, skyWcs in output.items():
1626  dataRef = visit_ccd_to_dataRef[key]
1627  try:
1628  dataRef.put(skyWcs, 'jointcal_wcs')
1629  except pexExceptions.Exception as e:
1630  self.log.fatal('Failed to write updated Wcs: %s', str(e))
1631  raise e
1632 
1633  def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1634  """
1635  Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1636 
1637  Parameters
1638  ----------
1639  associations : `lsst.jointcal.Associations`
1640  The star/reference star associations to fit.
1641  model : `lsst.jointcal.PhotometryModel`
1642  The photoometric model that was fit.
1643  visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1644  Dict of ccdImage identifiers to dataRefs that were fit.
1645  """
1646 
1647  ccdImageList = associations.getCcdImageList()
1648  output = self._make_output_make_output(ccdImageList, model, "toPhotoCalib")
1649  for key, photoCalib in output.items():
1650  dataRef = visit_ccd_to_dataRef[key]
1651  try:
1652  dataRef.put(photoCalib, 'jointcal_photoCalib')
1653  except pexExceptions.Exception as e:
1654  self.log.fatal('Failed to write updated PhotoCalib: %s', str(e))
1655  raise e
static std::shared_ptr< SourceTable > make(Schema const &schema, std::shared_ptr< IdFactory > const &idFactory)
static Schema makeMinimalSchema()
The class that implements the relations between MeasuredStar and FittedStar.
Definition: Associations.h:54
Class that handles the astrometric least squares problem.
Definition: AstrometryFit.h:78
Simple structure to accumulate chi2 and ndof.
Definition: Chi2.h:52
A multi-component model, fitting mappings for sensors and visits simultaneously.
A projection handler in which all CCDs from the same visit have the same tangent point.
Class that handles the photometric least squares problem.
Definition: PhotometryFit.h:45
A model where there is one independent transform per CcdImage.
def getTargetList(parsedCmd, **kwargs)
Definition: jointcal.py:76
def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations, jointcalControl)
Definition: jointcal.py:735
def _compute_proper_motion_epoch(self, ccdImageList)
Definition: jointcal.py:1056
def _get_refcat_coordinate_error_override(self, refCat, name)
Definition: jointcal.py:1009
def _getDebugPath(self, filename)
Definition: jointcal.py:910
def _check_star_lists(self, associations, name)
Definition: jointcal.py:1214
def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1610
def _check_stars(self, associations)
Definition: jointcal.py:1465
def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName, applyColorterms=False, epoch=None)
Definition: jointcal.py:1159
def _readDataId(self, butler, dataId)
Definition: jointcal.py:860
def runQuantum(self, butlerQC, inputRefs, outputRefs)
Definition: jointcal.py:527
def _build_ccdImage(self, data, associations, jointcalControl)
Definition: jointcal.py:814
def runDataRef(self, dataRefs)
Definition: jointcal.py:934
def __init__(self, butler=None, initInputs=None, **kwargs)
Definition: jointcal.py:506
def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId)
Definition: jointcal.py:691
def _make_one_input_data(self, visitRecord, catalog)
Definition: jointcal.py:786
def _do_load_refcat_and_fit(self, associations, defaultBand, center, radius, tract="", match_cut=3.0, reject_bad_fluxes=False, *name="", refObjLoader=None, referenceSelector=None, fit_function=None, epoch=None)
Definition: jointcal.py:1076
def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", doRankUpdate=True, doLineSearch=False)
Definition: jointcal.py:1481
def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None)
Definition: jointcal.py:614
def _fit_photometry(self, associations, dataName=None)
Definition: jointcal.py:1266
def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef)
Definition: jointcal.py:1633
def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None)
Definition: jointcal.py:1223
def _make_output(self, ccdImageList, model, func)
Definition: jointcal.py:1582
def _put_output(self, butlerQC, outputs, outputRefs, camera, setter)
Definition: jointcal.py:555
def _fit_astrometry(self, associations, dataName=None)
Definition: jointcal.py:1369
def _prep_sky(self, associations, bands)
Definition: jointcal.py:915
def loadData(self, dataRefs, associations, jointcalControl)
Definition: jointcal.py:886
def add_measurement(job, name, value)
Definition: jointcal.py:58
def writeModel(model, filename, log)
Definition: jointcal.py:460