lsst.jointcal  14.0-14-g932474c+11
jointcal.py
Go to the documentation of this file.
1 # See COPYRIGHT file at the top of the source tree.
2 
3 from __future__ import division, absolute_import, print_function
4 from builtins import str
5 from builtins import range
6 
7 import collections
8 import numpy as np
9 
10 import lsst.utils
11 import lsst.pex.config as pexConfig
12 import lsst.pipe.base as pipeBase
13 import lsst.afw.image as afwImage
14 import lsst.afw.geom as afwGeom
15 import lsst.afw.coord as afwCoord
16 import lsst.pex.exceptions as pexExceptions
17 import lsst.afw.table
19 
20 from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask
21 from lsst.meas.algorithms.sourceSelector import sourceSelectorRegistry
22 
23 from .dataIds import PerTractCcdDataIdContainer
24 
25 import lsst.jointcal
26 from lsst.jointcal import MinimizeResult
27 
28 __all__ = ["JointcalConfig", "JointcalTask"]
29 
30 Photometry = collections.namedtuple('Photometry', ('fit', 'model'))
31 Astrometry = collections.namedtuple('Astrometry', ('fit', 'model', 'sky_to_tan_projection'))
32 
33 
34 class JointcalRunner(pipeBase.ButlerInitializedTaskRunner):
35  """Subclass of TaskRunner for jointcalTask
36 
37  jointcalTask.run() takes a number of arguments, one of which is a list of dataRefs
38  extracted from the command line (whereas most CmdLineTasks' run methods take
39  single dataRef, are are called repeatedly). This class transforms the processed
40  arguments generated by the ArgumentParser into the arguments expected by
41  Jointcal.run().
42 
43  See pipeBase.TaskRunner for more information.
44  """
45 
46  @staticmethod
47  def getTargetList(parsedCmd, **kwargs):
48  """
49  Return a list of tuples per tract, each containing (dataRefs, kwargs).
50 
51  Jointcal operates on lists of dataRefs simultaneously.
52  """
53  kwargs['profile_jointcal'] = parsedCmd.profile_jointcal
54  kwargs['butler'] = parsedCmd.butler
55 
56  # organize data IDs by tract
57  refListDict = {}
58  for ref in parsedCmd.id.refList:
59  refListDict.setdefault(ref.dataId["tract"], []).append(ref)
60  # we call run() once with each tract
61  result = [(refListDict[tract], kwargs) for tract in sorted(refListDict.keys())]
62  return result
63 
64  def __call__(self, args):
65  """
66  @param args Arguments for Task.run()
67 
68  @return
69  - A pipe.base.Struct containing these fields if self.doReturnResults is False:
70  - ``exitStatus`: 0 if the task completed successfully, 1 otherwise.
71  - A pipe.base.Struct containing these fields if self.doReturnResults is True:
72  - ``result``: the result of calling jointcal.run()
73  - ``exitStatus`: 0 if the task completed successfully, 1 otherwise.
74  """
75  exitStatus = 0 # exit status for shell
76 
77  # NOTE: cannot call self.makeTask because that assumes args[0] is a single dataRef.
78  dataRefList, kwargs = args
79  butler = kwargs.pop('butler')
80  task = self.TaskClass(config=self.config, log=self.log, butler=butler)
81  result = None
82  try:
83  result = task.run(dataRefList, **kwargs)
84  exitStatus = result.exitStatus
85  except Exception as e: # catch everything, sort it out later.
86  if self.doRaise:
87  raise e
88  else:
89  exitStatus = 1
90  eName = type(e).__name__
91  task.log.fatal("Failed on dataIds=%s: %s: %s", dataRefList, eName, e)
92 
93  if self.doReturnResults:
94  return pipeBase.Struct(result=result, exitStatus=exitStatus)
95  else:
96  return pipeBase.Struct(exitStatus=exitStatus)
97 
98 
99 class JointcalConfig(pexConfig.Config):
100  """Config for jointcalTask"""
101 
102  doAstrometry = pexConfig.Field(
103  doc="Fit astrometry and write the fitted result.",
104  dtype=bool,
105  default=True
106  )
107  doPhotometry = pexConfig.Field(
108  doc="Fit photometry and write the fitted result.",
109  dtype=bool,
110  default=True
111  )
112  coaddName = pexConfig.Field(
113  doc="Type of coadd, typically deep or goodSeeing",
114  dtype=str,
115  default="deep"
116  )
117  posError = pexConfig.Field(
118  doc="Constant term for error on position (in pixel unit)",
119  dtype=float,
120  default=0.02,
121  )
122  # TODO: DM-6885 matchCut should be an afw.geom.Angle
123  matchCut = pexConfig.Field(
124  doc="Matching radius between fitted and reference stars (arcseconds)",
125  dtype=float,
126  default=3.0,
127  )
128  minMeasurements = pexConfig.Field(
129  doc="Minimum number of associated measured stars for a fitted star to be included in the fit",
130  dtype=int,
131  default=2,
132  )
133  polyOrder = pexConfig.Field(
134  doc="Polynomial order for fitting distorsion",
135  dtype=int,
136  default=3,
137  )
138  astrometryModel = pexConfig.ChoiceField(
139  doc="Type of model to fit to astrometry",
140  dtype=str,
141  default="simplePoly",
142  allowed={"simplePoly": "One polynomial per ccd",
143  "constrainedPoly": "One polynomial per ccd, and one polynomial per visit"}
144  )
145  photometryModel = pexConfig.ChoiceField(
146  doc="Type of model to fit to photometry",
147  dtype=str,
148  default="simple",
149  allowed={"simple": "One constant zeropoint per ccd and visit",
150  "constrained": "Constrained zeropoint per ccd, and one polynomial per visit"}
151  )
152  photometryVisitDegree = pexConfig.Field(
153  doc="Degree of the per-visit polynomial transform for the constrained photometry model.",
154  dtype=int,
155  default=7,
156  )
157  astrometryRefObjLoader = pexConfig.ConfigurableField(
158  target=LoadIndexedReferenceObjectsTask,
159  doc="Reference object loader for astrometric fit",
160  )
161  photometryRefObjLoader = pexConfig.ConfigurableField(
162  target=LoadIndexedReferenceObjectsTask,
163  doc="Reference object loader for photometric fit",
164  )
165  sourceSelector = sourceSelectorRegistry.makeField(
166  doc="How to select sources for cross-matching",
167  default="astrometry"
168  )
169  writeChi2ContributionFiles = pexConfig.Field(
170  dtype=bool,
171  doc="Write initial/final fit files containing the contributions to chi2.",
172  default=False
173  )
174 
175  def setDefaults(self):
176  sourceSelector = self.sourceSelector["astrometry"]
177  sourceSelector.setDefaults()
178  # don't want to lose existing flags, just add to them.
179  sourceSelector.badFlags.extend(["slot_Shape_flag"])
180  # This should be used to set the FluxField value in jointcal::JointcalControl
181  sourceSelector.sourceFluxType = 'Calib'
182 
183 
184 class JointcalTask(pipeBase.CmdLineTask):
185  """Jointly astrometrically (photometrically later) calibrate a group of images."""
186 
187  ConfigClass = JointcalConfig
188  RunnerClass = JointcalRunner
189  _DefaultName = "jointcal"
190 
191  def __init__(self, butler=None, profile_jointcal=False, **kwargs):
192  """
193  Instantiate a JointcalTask.
194 
195  Parameters
196  ----------
197  butler : lsst.daf.persistence.Butler
198  The butler is passed to the refObjLoader constructor in case it is
199  needed. Ignored if the refObjLoader argument provides a loader directly.
200  Used to initialize the astrometry and photometry refObjLoaders.
201  profile_jointcal : bool
202  set to True to profile different stages of this jointcal run.
203  """
204  pipeBase.CmdLineTask.__init__(self, **kwargs)
205  self.profile_jointcal = profile_jointcal
206  self.makeSubtask("sourceSelector")
207  if self.config.doAstrometry:
208  self.makeSubtask('astrometryRefObjLoader', butler=butler)
209  if self.config.doPhotometry:
210  self.makeSubtask('photometryRefObjLoader', butler=butler)
211 
212  # To hold various computed metrics for use by tests
213  self.metrics = {}
214 
215  # We don't need to persist config and metadata at this stage.
216  # In this way, we don't need to put a specific entry in the camera mapper policy file
217  def _getConfigName(self):
218  return None
219 
220  def _getMetadataName(self):
221  return None
222 
223  @classmethod
224  def _makeArgumentParser(cls):
225  """Create an argument parser"""
226  parser = pipeBase.ArgumentParser(name=cls._DefaultName)
227  parser.add_argument("--profile_jointcal", default=False, action="store_true",
228  help="Profile steps of jointcal separately.")
229  parser.add_id_argument("--id", "calexp", help="data ID, e.g. --id visit=6789 ccd=0..9",
230  ContainerClass=PerTractCcdDataIdContainer)
231  return parser
232 
233  def _build_ccdImage(self, dataRef, associations, jointcalControl):
234  """
235  Extract the necessary things from this dataRef to add a new ccdImage.
236 
237  Parameters
238  ----------
239  dataRef : lsst.daf.persistence.ButlerDataRef
240  dataRef to extract info from.
241  associations : lsst.jointcal.Associations
242  object to add the info to, to construct a new CcdImage
243  jointcalControl : jointcal.JointcalControl
244  control object for associations management
245 
246  Returns
247  ------
248  namedtuple
249  wcs : lsst.afw.image.TanWcs
250  the TAN WCS of this image, read from the calexp
251  key : namedtuple
252  a key to identify this dataRef by its visit and ccd ids
253  filter : str
254  this calexp's filter
255  """
256  if "visit" in dataRef.dataId.keys():
257  visit = dataRef.dataId["visit"]
258  else:
259  visit = dataRef.getButler().queryMetadata("calexp", ("visit"), dataRef.dataId)[0]
260 
261  src = dataRef.get("src", flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS, immediate=True)
262 
263  visitInfo = dataRef.get('calexp_visitInfo')
264  detector = dataRef.get('calexp_detector')
265  ccdname = detector.getId()
266  calib = dataRef.get('calexp_calib')
267  tanWcs = dataRef.get('calexp_wcs')
268  bbox = dataRef.get('calexp_bbox')
269  filt = dataRef.get('calexp_filter')
270  filterName = filt.getName()
271  fluxMag0 = calib.getFluxMag0()
272  photoCalib = afwImage.PhotoCalib(1.0/fluxMag0[0], fluxMag0[1]/fluxMag0[0]**2, bbox)
273 
274  goodSrc = self.sourceSelector.selectSources(src)
275 
276  if len(goodSrc.sourceCat) == 0:
277  self.log.warn("No sources selected in visit %s ccd %s", visit, ccdname)
278  else:
279  self.log.info("%d sources selected in visit %d ccd %d", len(goodSrc.sourceCat), visit, ccdname)
280  associations.addImage(goodSrc.sourceCat, tanWcs, visitInfo, bbox, filterName, photoCalib, detector,
281  visit, ccdname, jointcalControl)
282 
283  Result = collections.namedtuple('Result_from_build_CcdImage', ('wcs', 'key', 'filter'))
284  Key = collections.namedtuple('Key', ('visit', 'ccd'))
285  return Result(tanWcs, Key(visit, ccdname), filterName)
286 
287  @pipeBase.timeMethod
288  def run(self, dataRefs, profile_jointcal=False):
289  """
290  Jointly calibrate the astrometry and photometry across a set of images.
291 
292  Parameters
293  ----------
294  dataRefs : list of lsst.daf.persistence.ButlerDataRef
295  List of data references to the exposures to be fit.
296  profile_jointcal : bool
297  Profile the individual steps of jointcal.
298 
299  Returns
300  -------
301  pipe.base.Struct
302  struct containing:
303  * dataRefs: the provided data references that were fit (with updated WCSs)
304  * oldWcsList: the original WCS from each dataRef
305  * metrics: dictionary of internally-computed metrics for testing/validation.
306  """
307  if len(dataRefs) == 0:
308  raise ValueError('Need a non-empty list of data references!')
309 
310  exitStatus = 0 # exit status for shell
311 
312  sourceFluxField = "slot_%sFlux" % (self.sourceSelector.config.sourceFluxType,)
313  jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
314  associations = lsst.jointcal.Associations()
315 
316  visit_ccd_to_dataRef = {}
317  oldWcsList = []
318  filters = []
319  load_cat_prof_file = 'jointcal_build_ccdImage.prof' if profile_jointcal else ''
320  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
321  # We need the bounding-box of the focal plane for photometry visit models.
322  # NOTE: we only need to read it once, because its the same for all exposures of a camera.
323  camera = dataRefs[0].get('camera', immediate=True)
324  self.focalPlaneBBox = camera.getFpBBox()
325  for ref in dataRefs:
326  result = self._build_ccdImage(ref, associations, jointcalControl)
327  oldWcsList.append(result.wcs)
328  visit_ccd_to_dataRef[result.key] = ref
329  filters.append(result.filter)
330  filters = collections.Counter(filters)
331 
332  centers = [ccdImage.getBoresightRaDec() for ccdImage in associations.getCcdImageList()]
333  commonTangentPoint = lsst.afw.coord.averageCoord(centers)
334  self.log.debug("Using common tangent point: %s", commonTangentPoint.getPosition())
335  associations.setCommonTangentPoint(commonTangentPoint.getPosition())
336 
337  # Use external reference catalogs handled by LSST stack mechanism
338  # Get the bounding box overlapping all associated images
339  # ==> This is probably a bad idea to do it this way <== To be improved
340  bbox = associations.getRaDecBBox()
341  center = afwCoord.Coord(bbox.getCenter(), afwGeom.degrees)
342  corner = afwCoord.Coord(bbox.getMax(), afwGeom.degrees)
343  radius = center.angularSeparation(corner).asRadians()
344 
345  # Get astrometry_net_data path
346  anDir = lsst.utils.getPackageDir('astrometry_net_data')
347  if anDir is None:
348  raise RuntimeError("astrometry_net_data is not setup")
349 
350  # Determine a default filter associated with the catalog. See DM-9093
351  defaultFilter = filters.most_common(1)[0][0]
352  self.log.debug("Using %s band for reference flux", defaultFilter)
353 
354  # TODO: need a better way to get the tract.
355  tract = dataRefs[0].dataId['tract']
356 
357  if self.config.doAstrometry:
358  astrometry = self._do_load_refcat_and_fit(associations, defaultFilter, center, radius,
359  name="Astrometry",
360  refObjLoader=self.astrometryRefObjLoader,
361  fit_function=self._fit_astrometry,
362  profile_jointcal=profile_jointcal,
363  tract=tract)
364  else:
365  astrometry = Astrometry(None, None, None)
366 
367  if self.config.doPhotometry:
368  photometry = self._do_load_refcat_and_fit(associations, defaultFilter, center, radius,
369  name="Photometry",
370  refObjLoader=self.photometryRefObjLoader,
371  fit_function=self._fit_photometry,
372  profile_jointcal=profile_jointcal,
373  tract=tract,
374  filters=filters,
375  reject_bad_fluxes=True)
376  else:
377  photometry = Photometry(None, None)
378 
379  load_cat_prof_file = 'jointcal_write_results.prof' if profile_jointcal else ''
380  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
381  self._write_results(associations, astrometry.model, photometry.model, visit_ccd_to_dataRef)
382 
383  return pipeBase.Struct(dataRefs=dataRefs,
384  oldWcsList=oldWcsList,
385  metrics=self.metrics,
386  exitStatus=exitStatus)
387 
388  def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
389  name="", refObjLoader=None, filters=[], fit_function=None,
390  tract=None, profile_jointcal=False, match_cut=3.0,
391  reject_bad_fluxes=False):
392  """Load reference catalog, perform the fit, and return the result.
393 
394  Parameters
395  ----------
396  associations : lsst.jointcal.Associations
397  The star/reference star associations to fit.
398  defaultFilter : str
399  filter to load from reference catalog.
400  center : lsst.afw.coord.Coord
401  Center of field to load from reference catalog.
402  radius : lsst.afw.geom.Angle
403  On-sky radius to load from reference catalog.
404  name : str
405  Name of thing being fit: "Astrometry" or "Photometry".
406  refObjLoader : lsst.meas.algorithms.LoadReferenceObjectsTask
407  Reference object loader to load from for fit.
408  filters : list of str, optional
409  List of filters to load from the reference catalog.
410  fit_function : function
411  function to call to perform fit (takes associations object).
412  tract : str
413  Name of tract currently being fit.
414  profile_jointcal : bool, optional
415  Separately profile the fitting step.
416  match_cut : float, optional
417  Radius in arcseconds to find cross-catalog matches to during
418  associations.associateCatalogs.
419  reject_bad_fluxes : bool, optional
420  Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
421 
422  Returns
423  -------
424  Result of `fit_function()`
425  """
426  self.log.info("====== Now processing %s...", name)
427  # TODO: this should not print "trying to invert a singular transformation:"
428  # if it does that, something's not right about the WCS...
429  associations.associateCatalogs(match_cut)
430  self.metrics['associated%sFittedStars' % name] = associations.fittedStarListSize()
431 
432  skyCircle = refObjLoader.loadSkyCircle(center,
433  afwGeom.Angle(radius, afwGeom.radians),
434  defaultFilter)
435 
436  # Need memory contiguity to get reference filters as a vector.
437  if not skyCircle.refCat.isContiguous():
438  refCat = skyCircle.refCat.copy(deep=True)
439  else:
440  refCat = skyCircle.refCat
441 
442  # load the reference catalog fluxes.
443  # TODO: Simon will file a ticket for making this better (and making it use the color terms)
444  refFluxes = {}
445  refFluxErrs = {}
446  for filt in filters:
447  filtKeys = lsst.meas.algorithms.getRefFluxKeys(refCat.schema, filt)
448  refFluxes[filt] = refCat.get(filtKeys[0])
449  refFluxErrs[filt] = refCat.get(filtKeys[1])
450 
451  associations.collectRefStars(refCat, self.config.matchCut*afwGeom.arcseconds,
452  skyCircle.fluxField, refFluxes, refFluxErrs, reject_bad_fluxes)
453  self.metrics['collected%sRefStars' % name] = associations.refStarListSize()
454 
455  associations.selectFittedStars(self.config.minMeasurements)
456  self._check_star_lists(associations, name)
457  self.metrics['selected%sRefStars' % name] = associations.refStarListSize()
458  self.metrics['selected%sFittedStars' % name] = associations.fittedStarListSize()
459  self.metrics['selected%sCcdImageList' % name] = associations.nCcdImagesValidForFit()
460 
461  load_cat_prof_file = 'jointcal_fit_%s.prof'%name if profile_jointcal else ''
462  dataName = "{}_{}".format(tract, defaultFilter)
463  with pipeBase.cmdLineTask.profile(load_cat_prof_file):
464  result = fit_function(associations, dataName)
465  # TODO DM-12446: turn this into a "butler save" somehow.
466  # Save reference and measurement chi2 contributions for this data
467  if self.config.writeChi2ContributionFiles:
468  baseName = "{}_final_chi2-{}.csv".format(name, dataName)
469  result.fit.saveChi2Contributions(baseName)
470 
471  return result
472 
473  def _check_star_lists(self, associations, name):
474  # TODO: these should be len(blah), but we need this properly wrapped first.
475  if associations.nCcdImagesValidForFit() == 0:
476  raise RuntimeError('No images in the ccdImageList!')
477  if associations.fittedStarListSize() == 0:
478  raise RuntimeError('No stars in the {} fittedStarList!'.format(name))
479  if associations.refStarListSize() == 0:
480  raise RuntimeError('No stars in the {} reference star list!'.format(name))
481 
482  def _fit_photometry(self, associations, dataName=None):
483  """
484  Fit the photometric data.
485 
486  Parameters
487  ----------
488  associations : lsst.jointcal.Associations
489  The star/reference star associations to fit.
490  dataName : str
491  Name of the data being processed (e.g. "1234_HSC-Y"), for
492  identifying debugging files.
493 
494  Returns
495  -------
496  namedtuple
497  fit : lsst.jointcal.PhotometryFit
498  The photometric fitter used to perform the fit.
499  model : lsst.jointcal.PhotometryModel
500  The photometric model that was fit.
501  """
502  self.log.info("=== Starting photometric fitting...")
503 
504  # TODO: should use pex.config.RegistryField here (see DM-9195)
505  if self.config.photometryModel == "constrained":
506  model = lsst.jointcal.ConstrainedPhotometryModel(associations.getCcdImageList(),
507  self.focalPlaneBBox,
508  visitDegree=self.config.photometryVisitDegree)
509  elif self.config.photometryModel == "simple":
510  model = lsst.jointcal.SimplePhotometryModel(associations.getCcdImageList())
511 
512  fit = lsst.jointcal.PhotometryFit(associations, model)
513  chi2 = fit.computeChi2()
514  # TODO DM-12446: turn this into a "butler save" somehow.
515  # Save reference and measurement chi2 contributions for this data
516  if self.config.writeChi2ContributionFiles:
517  baseName = "Photometry_initial_chi2-{}.csv".format(dataName)
518  fit.saveChi2Contributions(baseName)
519 
520  if not np.isfinite(chi2.chi2):
521  raise FloatingPointError('Initial chi2 is invalid: %s'%chi2)
522  self.log.info("Initialized: %s", str(chi2))
523  fit.minimize("Model")
524  chi2 = fit.computeChi2()
525  self.log.info(str(chi2))
526  fit.minimize("Fluxes")
527  chi2 = fit.computeChi2()
528  self.log.info(str(chi2))
529  fit.minimize("Model Fluxes")
530  chi2 = fit.computeChi2()
531  if not np.isfinite(chi2.chi2):
532  raise FloatingPointError('Pre-iteration chi2 is invalid: %s'%chi2)
533  self.log.info("Fit prepared with %s", str(chi2))
534 
535  chi2 = self._iterate_fit(fit, model, 20, "photometry", "Model Fluxes")
536 
537  self.metrics['photometryFinalChi2'] = chi2.chi2
538  self.metrics['photometryFinalNdof'] = chi2.ndof
539  return Photometry(fit, model)
540 
541  def _fit_astrometry(self, associations, dataName=None):
542  """
543  Fit the astrometric data.
544 
545  Parameters
546  ----------
547  associations : lsst.jointcal.Associations
548  The star/reference star associations to fit.
549  dataName : str
550  Name of the data being processed (e.g. "1234_HSC-Y"), for
551  identifying debugging files.
552 
553  Returns
554  -------
555  namedtuple
556  fit : lsst.jointcal.AstrometryFit
557  The astrometric fitter used to perform the fit.
558  model : lsst.jointcal.AstrometryModel
559  The astrometric model that was fit.
560  sky_to_tan_projection : lsst.jointcal.ProjectionHandler
561  The model for the sky to tangent plane projection that was used in the fit.
562  """
563 
564  self.log.info("=== Starting astrometric fitting...")
565 
566  associations.deprojectFittedStars()
567 
568  # NOTE: need to return sky_to_tan_projection so that it doesn't get garbage collected.
569  # TODO: could we package sky_to_tan_projection and model together so we don't have to manage
570  # them so carefully?
571  sky_to_tan_projection = lsst.jointcal.OneTPPerVisitHandler(associations.getCcdImageList())
572 
573  if self.config.astrometryModel == "constrainedPoly":
574  model = lsst.jointcal.ConstrainedPolyModel(associations.getCcdImageList(),
575  sky_to_tan_projection, True, 0)
576  elif self.config.astrometryModel == "simplePoly":
577  model = lsst.jointcal.SimplePolyModel(associations.getCcdImageList(),
578  sky_to_tan_projection,
579  True, 0, self.config.polyOrder)
580 
581  fit = lsst.jointcal.AstrometryFit(associations, model, self.config.posError)
582  chi2 = fit.computeChi2()
583  # TODO DM-12446: turn this into a "butler save" somehow.
584  # Save reference and measurement chi2 contributions for this data
585  if self.config.writeChi2ContributionFiles:
586  baseName = "Astrometry_initial_chi2-{}.csv".format(dataName)
587  fit.saveChi2Contributions(baseName)
588 
589  if not np.isfinite(chi2.chi2):
590  raise FloatingPointError('Initial chi2 is invalid: %s'%chi2)
591  self.log.info("Initialized: %s", str(chi2))
592  fit.minimize("Distortions")
593  chi2 = fit.computeChi2()
594  self.log.info(str(chi2))
595  fit.minimize("Positions")
596  chi2 = fit.computeChi2()
597  self.log.info(str(chi2))
598  fit.minimize("Distortions Positions")
599  chi2 = fit.computeChi2()
600  self.log.info(str(chi2))
601  if not np.isfinite(chi2.chi2):
602  raise FloatingPointError('Pre-iteration chi2 is invalid: %s'%chi2)
603  self.log.info("Fit prepared with %s", str(chi2))
604 
605  chi2 = self._iterate_fit(fit, model, 20, "astrometry", "Distortions Positions")
606 
607  self.metrics['astrometryFinalChi2'] = chi2.chi2
608  self.metrics['astrometryFinalNdof'] = chi2.ndof
609 
610  return Astrometry(fit, model, sky_to_tan_projection)
611 
612  def _iterate_fit(self, fit, model, max_steps, name, whatToFit):
613  """Run fit.minimize up to max_steps times, returning the final chi2."""
614 
615  for i in range(max_steps):
616  r = fit.minimize(whatToFit, 5) # outlier removal at 5 sigma.
617  chi2 = fit.computeChi2()
618  if not np.isfinite(chi2.chi2):
619  raise FloatingPointError('Fit iteration chi2 is invalid: %s'%chi2)
620  self.log.info(str(chi2))
621  if r == MinimizeResult.Converged:
622  self.log.debug("fit has converged - no more outliers - redo minimixation"
623  "one more time in case we have lost accuracy in rank update")
624  # Redo minimization one more time in case we have lost accuracy in rank update
625  r = fit.minimize(whatToFit, 5) # outliers removal at 5 sigma.
626  chi2 = fit.computeChi2()
627  self.log.info("Fit completed with: %s", str(chi2))
628  break
629  elif r == MinimizeResult.Failed:
630  self.log.warn("minimization failed")
631  break
632  elif r == MinimizeResult.Chi2Increased:
633  self.log.warn("still some ouliers but chi2 increases - retry")
634  else:
635  self.log.error("unxepected return code from minimize")
636  break
637  else:
638  self.log.error("%s failed to converge after %d steps"%(name, max_steps))
639 
640  return chi2
641 
642  def _write_results(self, associations, astrometry_model, photometry_model, visit_ccd_to_dataRef):
643  """
644  Write the fitted results (photometric and astrometric) to a new 'wcs' dataRef.
645 
646  Parameters
647  ----------
648  associations : lsst.jointcal.Associations
649  The star/reference star associations to fit.
650  astrometry_model : lsst.jointcal.AstrometryModel
651  The astrometric model that was fit.
652  photometry_model : lsst.jointcal.PhotometryModel
653  The photometric model that was fit.
654  visit_ccd_to_dataRef : dict of Key: lsst.daf.persistence.ButlerDataRef
655  dict of ccdImage identifiers to dataRefs that were fit
656  """
657 
658  ccdImageList = associations.getCcdImageList()
659  for ccdImage in ccdImageList:
660  # TODO: there must be a better way to identify this ccdImage than a visit,ccd pair?
661  ccd = ccdImage.ccdId
662  visit = ccdImage.visit
663  dataRef = visit_ccd_to_dataRef[(visit, ccd)]
664  exp = afwImage.ExposureI(0, 0)
665  if self.config.doAstrometry:
666  self.log.info("Updating WCS for visit: %d, ccd: %d", visit, ccd)
667  tanSip = astrometry_model.produceSipWcs(ccdImage)
668  tanWcs = lsst.jointcal.gtransfoToTanWcs(tanSip, ccdImage.imageFrame, False)
669  exp.setWcs(tanWcs)
670  try:
671  dataRef.put(exp, 'wcs')
672  except pexExceptions.Exception as e:
673  self.log.fatal('Failed to write updated Wcs: %s', str(e))
674  raise e
675  if self.config.doPhotometry:
676  self.log.info("Updating PhotoCalib for visit: %d, ccd: %d", visit, ccd)
677  photoCalib = photometry_model.toPhotoCalib(ccdImage)
678  try:
679  dataRef.put(photoCalib, 'photoCalib')
680  except pexExceptions.Exception as e:
681  self.log.fatal('Failed to write updated PhotoCalib: %s', str(e))
682  raise e
this is the model used to fit independent CCDs, meaning that there is no instrument model...
def _build_ccdImage(self, dataRef, associations, jointcalControl)
Definition: jointcal.py:233
def _fit_photometry(self, associations, dataName=None)
Definition: jointcal.py:482
def getTargetList(parsedCmd, kwargs)
Definition: jointcal.py:47
def _fit_astrometry(self, associations, dataName=None)
Definition: jointcal.py:541
def _write_results(self, associations, astrometry_model, photometry_model, visit_ccd_to_dataRef)
Definition: jointcal.py:642
def _iterate_fit(self, fit, model, max_steps, name, whatToFit)
Definition: jointcal.py:612
def _check_star_lists(self, associations, name)
Definition: jointcal.py:473
The class that implements the relations between MeasuredStar and FittedStar.
Definition: Associations.h:28
A projection handler in which all CCDs from the same visit have the same tangent point.
std::string getPackageDir(std::string const &packageName)
This is the model used to fit mappings as the combination of a transformation depending on the chip n...
Class that handles the photometric least squares problem.
Definition: PhotometryFit.h:20
boost::shared_ptr< lsst::afw::image::TanWcs > gtransfoToTanWcs(const lsst::jointcal::TanSipPix2RaDec wcsTransfo, const lsst::jointcal::Frame &ccdFrame, const bool noLowOrderSipTerms=false)
Transform the other way around.
Class that handles the astrometric least squares problem.
Definition: AstrometryFit.h:55
Photometry model with constraints, .
def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius, name="", refObjLoader=None, filters=[], fit_function=None, tract=None, profile_jointcal=False, match_cut=3.0, reject_bad_fluxes=False)
Definition: jointcal.py:391
def run(self, dataRefs, profile_jointcal=False)
Definition: jointcal.py:288
std::shared_ptr< Coord > averageCoord(std::vector< std::shared_ptr< Coord const >> const coords, CoordSystem system=UNKNOWN)
Photometric response model which has a single photometric factor per CcdImage.
def __init__(self, butler=None, profile_jointcal=False, kwargs)
Definition: jointcal.py:191