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