27 import astropy.units
as u
39 from lsst.verify
import Job, Measurement
44 from .dataIds
import PerTractCcdDataIdContainer
49 __all__ = [
"JointcalConfig",
"JointcalRunner",
"JointcalTask"]
51 Photometry = collections.namedtuple(
'Photometry', (
'fit',
'model'))
52 Astrometry = collections.namedtuple(
'Astrometry', (
'fit',
'model',
'sky_to_tan_projection'))
57 meas = Measurement(job.metrics[name], value)
58 job.measurements.insert(meas)
62 """Subclass of TaskRunner for jointcalTask
64 jointcalTask.runDataRef() takes a number of arguments, one of which is a list of dataRefs
65 extracted from the command line (whereas most CmdLineTasks' runDataRef methods take
66 single dataRef, are are called repeatedly). This class transforms the processed
67 arguments generated by the ArgumentParser into the arguments expected by
68 Jointcal.runDataRef().
70 See pipeBase.TaskRunner for more information.
76 Return a list of tuples per tract, each containing (dataRefs, kwargs).
78 Jointcal operates on lists of dataRefs simultaneously.
80 kwargs[
'profile_jointcal'] = parsedCmd.profile_jointcal
81 kwargs[
'butler'] = parsedCmd.butler
85 for ref
in parsedCmd.id.refList:
86 refListDict.setdefault(ref.dataId[
"tract"], []).
append(ref)
88 result = [(refListDict[tract], kwargs)
for tract
in sorted(refListDict.keys())]
96 Arguments for Task.runDataRef()
101 if self.doReturnResults is False:
103 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
105 if self.doReturnResults is True:
107 - ``result``: the result of calling jointcal.runDataRef()
108 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
113 dataRefList, kwargs = args
114 butler = kwargs.pop(
'butler')
115 task = self.TaskClass(config=self.config, log=self.log, butler=butler)
118 result = task.runDataRef(dataRefList, **kwargs)
119 exitStatus = result.exitStatus
120 job_path = butler.get(
'verify_job_filename')
121 result.job.write(job_path[0])
122 except Exception
as e:
127 eName =
type(e).__name__
128 tract = dataRefList[0].dataId[
'tract']
129 task.log.fatal(
"Failed processing tract %s, %s: %s", tract, eName, e)
132 kwargs[
'butler'] = butler
133 if self.doReturnResults:
134 return pipeBase.Struct(result=result, exitStatus=exitStatus)
136 return pipeBase.Struct(exitStatus=exitStatus)
140 """Configuration for JointcalTask"""
142 doAstrometry = pexConfig.Field(
143 doc=
"Fit astrometry and write the fitted result.",
147 doPhotometry = pexConfig.Field(
148 doc=
"Fit photometry and write the fitted result.",
152 coaddName = pexConfig.Field(
153 doc=
"Type of coadd, typically deep or goodSeeing",
157 positionErrorPedestal = pexConfig.Field(
158 doc=
"Systematic term to apply to the measured position error (pixels)",
162 photometryErrorPedestal = pexConfig.Field(
163 doc=
"Systematic term to apply to the measured error on flux or magnitude as a "
164 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
169 matchCut = pexConfig.Field(
170 doc=
"Matching radius between fitted and reference stars (arcseconds)",
174 minMeasurements = pexConfig.Field(
175 doc=
"Minimum number of associated measured stars for a fitted star to be included in the fit",
179 minMeasuredStarsPerCcd = pexConfig.Field(
180 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
184 minRefStarsPerCcd = pexConfig.Field(
185 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
189 allowLineSearch = pexConfig.Field(
190 doc=
"Allow a line search during minimization, if it is reasonable for the model"
191 " (models with a significant non-linear component, e.g. constrainedPhotometry).",
195 astrometrySimpleOrder = pexConfig.Field(
196 doc=
"Polynomial order for fitting the simple astrometry model.",
200 astrometryChipOrder = pexConfig.Field(
201 doc=
"Order of the per-chip transform for the constrained astrometry model.",
205 astrometryVisitOrder = pexConfig.Field(
206 doc=
"Order of the per-visit transform for the constrained astrometry model.",
210 useInputWcs = pexConfig.Field(
211 doc=
"Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
215 astrometryModel = pexConfig.ChoiceField(
216 doc=
"Type of model to fit to astrometry",
218 default=
"constrained",
219 allowed={
"simple":
"One polynomial per ccd",
220 "constrained":
"One polynomial per ccd, and one polynomial per visit"}
222 photometryModel = pexConfig.ChoiceField(
223 doc=
"Type of model to fit to photometry",
225 default=
"constrainedMagnitude",
226 allowed={
"simpleFlux":
"One constant zeropoint per ccd and visit, fitting in flux space.",
227 "constrainedFlux":
"Constrained zeropoint per ccd, and one polynomial per visit,"
228 " fitting in flux space.",
229 "simpleMagnitude":
"One constant zeropoint per ccd and visit,"
230 " fitting in magnitude space.",
231 "constrainedMagnitude":
"Constrained zeropoint per ccd, and one polynomial per visit,"
232 " fitting in magnitude space.",
235 applyColorTerms = pexConfig.Field(
236 doc=
"Apply photometric color terms to reference stars?"
237 "Requires that colorterms be set to a ColortermLibrary",
241 colorterms = pexConfig.ConfigField(
242 doc=
"Library of photometric reference catalog name to color term dict.",
243 dtype=ColortermLibrary,
245 photometryVisitOrder = pexConfig.Field(
246 doc=
"Order of the per-visit polynomial transform for the constrained photometry model.",
250 photometryDoRankUpdate = pexConfig.Field(
251 doc=(
"Do the rank update step during minimization. "
252 "Skipping this can help deal with models that are too non-linear."),
256 astrometryDoRankUpdate = pexConfig.Field(
257 doc=(
"Do the rank update step during minimization (should not change the astrometry fit). "
258 "Skipping this can help deal with models that are too non-linear."),
262 outlierRejectSigma = pexConfig.Field(
263 doc=
"How many sigma to reject outliers at during minimization.",
267 maxPhotometrySteps = pexConfig.Field(
268 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
272 maxAstrometrySteps = pexConfig.Field(
273 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
277 astrometryRefObjLoader = pexConfig.ConfigurableField(
278 target=LoadIndexedReferenceObjectsTask,
279 doc=
"Reference object loader for astrometric fit",
281 photometryRefObjLoader = pexConfig.ConfigurableField(
282 target=LoadIndexedReferenceObjectsTask,
283 doc=
"Reference object loader for photometric fit",
285 sourceSelector = sourceSelectorRegistry.makeField(
286 doc=
"How to select sources for cross-matching",
289 astrometryReferenceSelector = pexConfig.ConfigurableField(
290 target=ReferenceSourceSelectorTask,
291 doc=
"How to down-select the loaded astrometry reference catalog.",
293 photometryReferenceSelector = pexConfig.ConfigurableField(
294 target=ReferenceSourceSelectorTask,
295 doc=
"How to down-select the loaded photometry reference catalog.",
297 astrometryReferenceErr = pexConfig.Field(
298 doc=(
"Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
299 "If None, then raise an exception if the reference catalog is missing coordinate errors. "
300 "If specified, overrides any existing `coord_*Err` values."),
305 writeInitMatrix = pexConfig.Field(
307 doc=(
"Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
308 "The output files will be of the form 'astrometry_preinit-mat.txt', in the current directory. "
309 "Note that these files are the dense versions of the matrix, and so may be very large."),
312 writeChi2FilesInitialFinal = pexConfig.Field(
314 doc=
"Write .csv files containing the contributions to chi2 for the initialization and final fit.",
317 writeChi2FilesOuterLoop = pexConfig.Field(
319 doc=
"Write .csv files containing the contributions to chi2 for the outer fit loop.",
322 writeInitialModel = pexConfig.Field(
324 doc=(
"Write the pre-initialization model to text files, for debugging."
325 " Output is written to `initial[Astro|Photo]metryModel.txt` in the current working directory."),
328 debugOutputPath = pexConfig.Field(
331 doc=(
"Path to write debug output files to. Used by "
332 "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
334 sourceFluxType = pexConfig.Field(
336 doc=
"Source flux field to use in source selection and to get fluxes from the catalog.",
343 msg =
"applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
344 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
346 msg = (
"Only doing astrometry, but Colorterms are not applied for astrometry;"
347 "applyColorTerms=True will be ignored.")
361 fluxField = f
"slot_{self.sourceFluxType}Flux_instFlux"
362 self.
sourceSelector[
'science'].signalToNoise.fluxField = fluxField
363 self.
sourceSelector[
'science'].signalToNoise.errField = fluxField +
"Err"
369 badFlags = [
'base_PixelFlags_flag_edge',
'base_PixelFlags_flag_saturated',
370 'base_PixelFlags_flag_interpolatedCenter',
'base_SdssCentroid_flag',
371 'base_PsfFlux_flag',
'base_PixelFlags_flag_suspectCenter']
383 """Write model to outfile."""
384 with open(filename,
"w")
as file:
385 file.write(repr(model))
386 log.info(
"Wrote %s to file: %s", model, filename)
390 """Jointly astrometrically and photometrically calibrate a group of images."""
392 ConfigClass = JointcalConfig
393 RunnerClass = JointcalRunner
394 _DefaultName =
"jointcal"
396 def __init__(self, butler=None, profile_jointcal=False, **kwargs):
398 Instantiate a JointcalTask.
402 butler : `lsst.daf.persistence.Butler`
403 The butler is passed to the refObjLoader constructor in case it is
404 needed. Ignored if the refObjLoader argument provides a loader directly.
405 Used to initialize the astrometry and photometry refObjLoaders.
406 profile_jointcal : `bool`
407 Set to True to profile different stages of this jointcal run.
409 pipeBase.CmdLineTask.__init__(self, **kwargs)
411 self.makeSubtask(
"sourceSelector")
412 if self.config.doAstrometry:
413 self.makeSubtask(
'astrometryRefObjLoader', butler=butler)
414 self.makeSubtask(
"astrometryReferenceSelector")
417 if self.config.doPhotometry:
418 self.makeSubtask(
'photometryRefObjLoader', butler=butler)
419 self.makeSubtask(
"photometryReferenceSelector")
424 self.
job = Job.load_metrics_package(subset=
'jointcal')
429 def _getMetadataName(self):
433 def _makeArgumentParser(cls):
434 """Create an argument parser"""
436 parser.add_argument(
"--profile_jointcal", default=
False, action=
"store_true",
437 help=
"Profile steps of jointcal separately.")
438 parser.add_id_argument(
"--id",
"calexp", help=
"data ID, e.g. --id visit=6789 ccd=0..9",
439 ContainerClass=PerTractCcdDataIdContainer)
442 def _build_ccdImage(self, dataRef, associations, jointcalControl):
444 Extract the necessary things from this dataRef to add a new ccdImage.
448 dataRef : `lsst.daf.persistence.ButlerDataRef`
449 DataRef to extract info from.
450 associations : `lsst.jointcal.Associations`
451 Object to add the info to, to construct a new CcdImage
452 jointcalControl : `jointcal.JointcalControl`
453 Control object for associations management
459 The TAN WCS of this image, read from the calexp
460 (`lsst.afw.geom.SkyWcs`).
462 A key to identify this dataRef by its visit and ccd ids
465 This calexp's filter (`str`).
467 if "visit" in dataRef.dataId.keys():
468 visit = dataRef.dataId[
"visit"]
470 visit = dataRef.getButler().queryMetadata(
"calexp", (
"visit"), dataRef.dataId)[0]
472 src = dataRef.get(
"src", flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS, immediate=
True)
474 visitInfo = dataRef.get(
'calexp_visitInfo')
475 detector = dataRef.get(
'calexp_detector')
476 ccdId = detector.getId()
477 photoCalib = dataRef.get(
'calexp_photoCalib')
478 tanWcs = dataRef.get(
'calexp_wcs')
479 bbox = dataRef.get(
'calexp_bbox')
480 filt = dataRef.get(
'calexp_filter')
481 filterName = filt.getName()
483 goodSrc = self.sourceSelector.
run(src)
485 if len(goodSrc.sourceCat) == 0:
486 self.log.
warn(
"No sources selected in visit %s ccd %s", visit, ccdId)
488 self.log.
info(
"%d sources selected in visit %d ccd %d", len(goodSrc.sourceCat), visit, ccdId)
489 associations.createCcdImage(goodSrc.sourceCat,
500 Result = collections.namedtuple(
'Result_from_build_CcdImage', (
'wcs',
'key',
'filter'))
501 Key = collections.namedtuple(
'Key', (
'visit',
'ccd'))
502 return Result(tanWcs, Key(visit, ccdId), filterName)
504 def _getDebugPath(self, filename):
505 """Constructs a path to filename using the configured debug path.
507 return os.path.join(self.config.debugOutputPath, filename)
512 Jointly calibrate the astrometry and photometry across a set of images.
516 dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef`
517 List of data references to the exposures to be fit.
518 profile_jointcal : `bool`
519 Profile the individual steps of jointcal.
523 result : `lsst.pipe.base.Struct`
524 Struct of metadata from the fit, containing:
527 The provided data references that were fit (with updated WCSs)
529 The original WCS from each dataRef
531 Dictionary of internally-computed metrics for testing/validation.
533 if len(dataRefs) == 0:
534 raise ValueError(
'Need a non-empty list of data references!')
538 sourceFluxField =
"slot_%sFlux" % (self.config.sourceFluxType,)
542 visit_ccd_to_dataRef = {}
545 load_cat_prof_file =
'jointcal_build_ccdImage.prof' if profile_jointcal
else ''
546 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
549 camera = dataRefs[0].get(
'camera', immediate=
True)
553 oldWcsList.append(result.wcs)
554 visit_ccd_to_dataRef[result.key] = ref
555 filters.append(result.filter)
556 filters = collections.Counter(filters)
558 associations.computeCommonTangentPoint()
560 boundingCircle = associations.computeBoundingCircle()
562 radius =
lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
564 self.log.
info(f
"Data has center={center} with radius={radius.asDegrees()} degrees.")
567 defaultFilter = filters.most_common(1)[0][0]
568 self.log.
debug(
"Using %s band for reference flux", defaultFilter)
571 tract = dataRefs[0].dataId[
'tract']
573 if self.config.doAstrometry:
577 referenceSelector=self.astrometryReferenceSelector,
579 profile_jointcal=profile_jointcal,
585 if self.config.doPhotometry:
589 referenceSelector=self.photometryReferenceSelector,
591 profile_jointcal=profile_jointcal,
594 reject_bad_fluxes=
True)
599 return pipeBase.Struct(dataRefs=dataRefs,
600 oldWcsList=oldWcsList,
604 defaultFilter=defaultFilter,
605 exitStatus=exitStatus)
607 def _get_refcat_coordinate_error_override(self, refCat, name):
608 """Check whether we should override the refcat coordinate errors, and
609 return the overridden error if necessary.
613 refCat : `lsst.afw.table.SimpleCatalog`
614 The reference catalog to check for a ``coord_raErr`` field.
616 Whether we are doing "astrometry" or "photometry".
620 refCoordErr : `float`
621 The refcat coordinate error to use, or NaN if we are not overriding
626 lsst.pex.config.FieldValidationError
627 Raised if the refcat does not contain coordinate errors and
628 ``config.astrometryReferenceErr`` is not set.
632 if name.lower() ==
"photometry":
633 if 'coord_raErr' not in refCat.schema:
638 if self.config.astrometryReferenceErr
is None and 'coord_raErr' not in refCat.schema:
639 msg = (
"Reference catalog does not contain coordinate errors, "
640 "and config.astrometryReferenceErr not supplied.")
641 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
645 if self.config.astrometryReferenceErr
is not None and 'coord_raErr' in refCat.schema:
646 self.log.
warn(
"Overriding reference catalog coordinate errors with %f/coordinate [mas]",
647 self.config.astrometryReferenceErr)
649 if self.config.astrometryReferenceErr
is None:
652 return self.config.astrometryReferenceErr
654 def _compute_proper_motion_epoch(self, ccdImageList):
655 """Return the proper motion correction epoch of the provided images.
659 ccdImageList : `list` [`lsst.jointcal.CcdImage`]
660 The images to compute the appropriate epoch for.
664 epoch : `astropy.time.Time`
665 The date to use for proper motion corrections.
667 mjds = [ccdImage.getMjd()
for ccdImage
in ccdImageList]
668 return astropy.time.Time(np.mean(mjds), format=
'mjd', scale=
"tai")
670 def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
672 tract="", profile_jointcal=False, match_cut=3.0,
673 reject_bad_fluxes=False, *,
674 name="", refObjLoader=None, referenceSelector=None,
676 """Load reference catalog, perform the fit, and return the result.
680 associations : `lsst.jointcal.Associations`
681 The star/reference star associations to fit.
682 defaultFilter : `str`
683 filter to load from reference catalog.
684 center : `lsst.geom.SpherePoint`
685 ICRS center of field to load from reference catalog.
686 radius : `lsst.geom.Angle`
687 On-sky radius to load from reference catalog.
689 Name of thing being fit: "astrometry" or "photometry".
690 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
691 Reference object loader to use to load a reference catalog.
692 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
693 Selector to use to pick objects from the loaded reference catalog.
694 fit_function : callable
695 Function to call to perform fit (takes Associations object).
696 filters : `list` [`str`], optional
697 List of filters to load from the reference catalog.
698 tract : `str`, optional
699 Name of tract currently being fit.
700 profile_jointcal : `bool`, optional
701 Separately profile the fitting step.
702 match_cut : `float`, optional
703 Radius in arcseconds to find cross-catalog matches to during
704 associations.associateCatalogs.
705 reject_bad_fluxes : `bool`, optional
706 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
710 result : `Photometry` or `Astrometry`
711 Result of `fit_function()`
713 self.log.
info(
"====== Now processing %s...", name)
716 associations.associateCatalogs(match_cut)
718 associations.fittedStarListSize())
720 applyColorterms =
False if name.lower() ==
"astrometry" else self.config.applyColorTerms
723 center, radius, defaultFilter,
724 applyColorterms=applyColorterms,
728 associations.collectRefStars(refCat,
729 self.config.matchCut*lsst.geom.arcseconds,
731 refCoordinateErr=refCoordErr,
732 rejectBadFluxes=reject_bad_fluxes)
734 associations.refStarListSize())
736 associations.prepareFittedStars(self.config.minMeasurements)
740 associations.nFittedStarsWithAssociatedRefStar())
742 associations.fittedStarListSize())
744 associations.nCcdImagesValidForFit())
746 load_cat_prof_file =
'jointcal_fit_%s.prof'%name
if profile_jointcal
else ''
747 dataName =
"{}_{}".
format(tract, defaultFilter)
748 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
749 result = fit_function(associations, dataName)
752 if self.config.writeChi2FilesInitialFinal:
753 baseName = self.
_getDebugPath(f
"{name}_final_chi2-{dataName}")
754 result.fit.saveChi2Contributions(baseName+
"{type}")
755 self.log.
info(
"Wrote chi2 contributions files: %s", baseName)
759 def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName,
760 applyColorterms=False, epoch=None):
761 """Load the necessary reference catalog sources, convert fluxes to
762 correct units, and apply color term corrections if requested.
766 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
767 The reference catalog loader to use to get the data.
768 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
769 Source selector to apply to loaded reference catalog.
770 center : `lsst.geom.SpherePoint`
771 The center around which to load sources.
772 radius : `lsst.geom.Angle`
773 The radius around ``center`` to load sources in.
775 The name of the camera filter to load fluxes for.
776 applyColorterms : `bool`
777 Apply colorterm corrections to the refcat for ``filterName``?
778 epoch : `astropy.time.Time`, optional
779 Epoch to which to correct refcat proper motion and parallax,
780 or `None` to not apply such corrections.
784 refCat : `lsst.afw.table.SimpleCatalog`
785 The loaded reference catalog.
787 The name of the reference catalog flux field appropriate for ``filterName``.
789 skyCircle = refObjLoader.loadSkyCircle(center,
794 selected = referenceSelector.run(skyCircle.refCat)
796 if not selected.sourceCat.isContiguous():
797 refCat = selected.sourceCat.copy(deep=
True)
799 refCat = selected.sourceCat
802 refCatName = refObjLoader.ref_dataset_name
803 self.log.
info(
"Applying color terms for filterName=%r reference catalog=%s",
804 filterName, refCatName)
805 colorterm = self.config.colorterms.getColorterm(
806 filterName=filterName, photoCatName=refCatName, doRaise=
True)
808 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat, filterName)
809 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
813 return refCat, skyCircle.fluxField
815 def _check_star_lists(self, associations, name):
817 if associations.nCcdImagesValidForFit() == 0:
818 raise RuntimeError(
'No images in the ccdImageList!')
819 if associations.fittedStarListSize() == 0:
820 raise RuntimeError(
'No stars in the {} fittedStarList!'.
format(name))
821 if associations.refStarListSize() == 0:
822 raise RuntimeError(
'No stars in the {} reference star list!'.
format(name))
824 def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
825 """Compute chi2, log it, validate the model, and return chi2.
829 associations : `lsst.jointcal.Associations`
830 The star/reference star associations to fit.
831 fit : `lsst.jointcal.FitterBase`
832 The fitter to use for minimization.
833 model : `lsst.jointcal.Model`
836 Label to describe the chi2 (e.g. "Initialized", "Final").
837 writeChi2Name : `str`, optional
838 Filename prefix to write the chi2 contributions to.
839 Do not supply an extension: an appropriate one will be added.
843 chi2: `lsst.jointcal.Chi2Accumulator`
844 The chi2 object for the current fitter and model.
849 Raised if chi2 is infinite or NaN.
851 Raised if the model is not valid.
853 if writeChi2Name
is not None:
855 fit.saveChi2Contributions(fullpath+
"{type}")
856 self.log.
info(
"Wrote chi2 contributions files: %s", fullpath)
858 chi2 = fit.computeChi2()
859 self.log.
info(
"%s %s", chi2Label, chi2)
861 if not np.isfinite(chi2.chi2):
862 raise FloatingPointError(f
'{chi2Label} chi2 is invalid: {chi2}')
863 if not model.validate(associations.getCcdImageList(), chi2.ndof):
864 raise ValueError(
"Model is not valid: check log messages for warnings.")
867 def _fit_photometry(self, associations, dataName=None):
869 Fit the photometric data.
873 associations : `lsst.jointcal.Associations`
874 The star/reference star associations to fit.
876 Name of the data being processed (e.g. "1234_HSC-Y"), for
877 identifying debugging files.
881 fit_result : `namedtuple`
882 fit : `lsst.jointcal.PhotometryFit`
883 The photometric fitter used to perform the fit.
884 model : `lsst.jointcal.PhotometryModel`
885 The photometric model that was fit.
887 self.log.
info(
"=== Starting photometric fitting...")
890 if self.config.photometryModel ==
"constrainedFlux":
893 visitOrder=self.config.photometryVisitOrder,
894 errorPedestal=self.config.photometryErrorPedestal)
896 doLineSearch = self.config.allowLineSearch
897 elif self.config.photometryModel ==
"constrainedMagnitude":
900 visitOrder=self.config.photometryVisitOrder,
901 errorPedestal=self.config.photometryErrorPedestal)
903 doLineSearch = self.config.allowLineSearch
904 elif self.config.photometryModel ==
"simpleFlux":
906 errorPedestal=self.config.photometryErrorPedestal)
908 elif self.config.photometryModel ==
"simpleMagnitude":
910 errorPedestal=self.config.photometryErrorPedestal)
916 if self.config.writeChi2FilesInitialFinal:
917 baseName = f
"photometry_initial_chi2-{dataName}"
920 if self.config.writeInitialModel:
925 def getChi2Name(whatToFit):
926 if self.config.writeChi2FilesOuterLoop:
927 return f
"photometry_init-%s_chi2-{dataName}" % whatToFit
933 dumpMatrixFile = self.
_getDebugPath(
"photometry_preinit")
if self.config.writeInitMatrix
else ""
934 if self.config.photometryModel.startswith(
"constrained"):
937 fit.minimize(
"ModelVisit", dumpMatrixFile=dumpMatrixFile)
939 writeChi2Name=getChi2Name(
"ModelVisit"))
942 fit.minimize(
"Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
944 writeChi2Name=getChi2Name(
"Model"))
946 fit.minimize(
"Fluxes")
948 writeChi2Name=getChi2Name(
"Fluxes"))
950 fit.minimize(
"Model Fluxes", doLineSearch=doLineSearch)
952 writeChi2Name=getChi2Name(
"ModelFluxes"))
954 model.freezeErrorTransform()
955 self.log.
debug(
"Photometry error scales are frozen.")
959 self.config.maxPhotometrySteps,
962 doRankUpdate=self.config.photometryDoRankUpdate,
963 doLineSearch=doLineSearch,
970 def _fit_astrometry(self, associations, dataName=None):
972 Fit the astrometric data.
976 associations : `lsst.jointcal.Associations`
977 The star/reference star associations to fit.
979 Name of the data being processed (e.g. "1234_HSC-Y"), for
980 identifying debugging files.
984 fit_result : `namedtuple`
985 fit : `lsst.jointcal.AstrometryFit`
986 The astrometric fitter used to perform the fit.
987 model : `lsst.jointcal.AstrometryModel`
988 The astrometric model that was fit.
989 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
990 The model for the sky to tangent plane projection that was used in the fit.
993 self.log.
info(
"=== Starting astrometric fitting...")
995 associations.deprojectFittedStars()
1002 if self.config.astrometryModel ==
"constrained":
1004 sky_to_tan_projection,
1005 chipOrder=self.config.astrometryChipOrder,
1006 visitOrder=self.config.astrometryVisitOrder)
1007 elif self.config.astrometryModel ==
"simple":
1009 sky_to_tan_projection,
1010 self.config.useInputWcs,
1012 order=self.config.astrometrySimpleOrder)
1017 if self.config.writeChi2FilesInitialFinal:
1018 baseName = f
"astrometry_initial_chi2-{dataName}"
1021 if self.config.writeInitialModel:
1026 def getChi2Name(whatToFit):
1027 if self.config.writeChi2FilesOuterLoop:
1028 return f
"astrometry_init-%s_chi2-{dataName}" % whatToFit
1032 dumpMatrixFile = self.
_getDebugPath(
"astrometry_preinit")
if self.config.writeInitMatrix
else ""
1035 if self.config.astrometryModel ==
"constrained":
1036 fit.minimize(
"DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1038 writeChi2Name=getChi2Name(
"DistortionsVisit"))
1041 fit.minimize(
"Distortions", dumpMatrixFile=dumpMatrixFile)
1043 writeChi2Name=getChi2Name(
"Distortions"))
1045 fit.minimize(
"Positions")
1047 writeChi2Name=getChi2Name(
"Positions"))
1049 fit.minimize(
"Distortions Positions")
1051 writeChi2Name=getChi2Name(
"DistortionsPositions"))
1055 self.config.maxAstrometrySteps,
1057 "Distortions Positions",
1058 doRankUpdate=self.config.astrometryDoRankUpdate,
1064 return Astrometry(fit, model, sky_to_tan_projection)
1066 def _check_stars(self, associations):
1067 """Count measured and reference stars per ccd and warn/log them."""
1068 for ccdImage
in associations.getCcdImageList():
1069 nMeasuredStars, nRefStars = ccdImage.countStars()
1070 self.log.
debug(
"ccdImage %s has %s measured and %s reference stars",
1071 ccdImage.getName(), nMeasuredStars, nRefStars)
1072 if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1073 self.log.
warn(
"ccdImage %s has only %s measuredStars (desired %s)",
1074 ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1075 if nRefStars < self.config.minRefStarsPerCcd:
1076 self.log.
warn(
"ccdImage %s has only %s RefStars (desired %s)",
1077 ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1079 def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1082 doLineSearch=False):
1083 """Run fitter.minimize up to max_steps times, returning the final chi2.
1087 associations : `lsst.jointcal.Associations`
1088 The star/reference star associations to fit.
1089 fitter : `lsst.jointcal.FitterBase`
1090 The fitter to use for minimization.
1092 Maximum number of steps to run outlier rejection before declaring
1093 convergence failure.
1094 name : {'photometry' or 'astrometry'}
1095 What type of data are we fitting (for logs and debugging files).
1097 Passed to ``fitter.minimize()`` to define the parameters to fit.
1098 dataName : `str`, optional
1099 Descriptive name for this dataset (e.g. tract and filter),
1101 doRankUpdate : `bool`, optional
1102 Do an Eigen rank update during minimization, or recompute the full
1103 matrix and gradient?
1104 doLineSearch : `bool`, optional
1105 Do a line search for the optimum step during minimization?
1109 chi2: `lsst.jointcal.Chi2Statistic`
1110 The final chi2 after the fit converges, or is forced to end.
1115 Raised if the fitter fails with a non-finite value.
1117 Raised if the fitter fails for some other reason;
1118 log messages will provide further details.
1120 dumpMatrixFile = self.
_getDebugPath(f
"{name}_postinit")
if self.config.writeInitMatrix
else ""
1122 oldChi2.chi2 = float(
"inf")
1123 for i
in range(max_steps):
1124 if self.config.writeChi2FilesOuterLoop:
1125 writeChi2Name = f
"{name}_iterate_{i}_chi2-{dataName}"
1127 writeChi2Name =
None
1128 result = fitter.minimize(whatToFit,
1129 self.config.outlierRejectSigma,
1130 doRankUpdate=doRankUpdate,
1131 doLineSearch=doLineSearch,
1132 dumpMatrixFile=dumpMatrixFile)
1135 f
"Fit iteration {i}", writeChi2Name=writeChi2Name)
1137 if result == MinimizeResult.Converged:
1139 self.log.
debug(
"fit has converged - no more outliers - redo minimization "
1140 "one more time in case we have lost accuracy in rank update.")
1142 result = fitter.minimize(whatToFit, self.config.outlierRejectSigma)
1146 if chi2.chi2/chi2.ndof >= 4.0:
1147 self.log.
error(
"Potentially bad fit: High chi-squared/ndof.")
1150 elif result == MinimizeResult.Chi2Increased:
1151 self.log.
warn(
"Still some outliers remaining but chi2 increased - retry")
1153 chi2Ratio = chi2.chi2 / oldChi2.chi2
1155 self.log.
warn(
'Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1156 chi2.chi2, oldChi2.chi2, chi2Ratio)
1163 msg = (
"Large chi2 increase between steps: fit likely cannot converge."
1164 " Try setting one or more of the `writeChi2*` config fields and looking"
1165 " at how individual star chi2-values evolve during the fit.")
1166 raise RuntimeError(msg)
1168 elif result == MinimizeResult.NonFinite:
1171 fitter.saveChi2Contributions(filename+
"{type}")
1172 msg =
"Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1173 raise FloatingPointError(msg.format(filename))
1174 elif result == MinimizeResult.Failed:
1175 raise RuntimeError(
"Chi2 minimization failure, cannot complete fit.")
1177 raise RuntimeError(
"Unxepected return code from minimize().")
1179 self.log.
error(
"%s failed to converge after %d steps"%(name, max_steps))
1183 def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1185 Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1189 associations : `lsst.jointcal.Associations`
1190 The star/reference star associations to fit.
1191 model : `lsst.jointcal.AstrometryModel`
1192 The astrometric model that was fit.
1193 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1194 Dict of ccdImage identifiers to dataRefs that were fit.
1197 ccdImageList = associations.getCcdImageList()
1198 for ccdImage
in ccdImageList:
1200 ccd = ccdImage.ccdId
1201 visit = ccdImage.visit
1202 dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1203 self.log.
info(
"Updating WCS for visit: %d, ccd: %d", visit, ccd)
1204 skyWcs = model.makeSkyWcs(ccdImage)
1206 dataRef.put(skyWcs,
'jointcal_wcs')
1208 self.log.
fatal(
'Failed to write updated Wcs: %s', str(e))
1211 def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1213 Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1217 associations : `lsst.jointcal.Associations`
1218 The star/reference star associations to fit.
1219 model : `lsst.jointcal.PhotometryModel`
1220 The photoometric model that was fit.
1221 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1222 Dict of ccdImage identifiers to dataRefs that were fit.
1225 ccdImageList = associations.getCcdImageList()
1226 for ccdImage
in ccdImageList:
1228 ccd = ccdImage.ccdId
1229 visit = ccdImage.visit
1230 dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1231 self.log.
info(
"Updating PhotoCalib for visit: %d, ccd: %d", visit, ccd)
1232 photoCalib = model.toPhotoCalib(ccdImage)
1234 dataRef.put(photoCalib,
'jointcal_photoCalib')
1236 self.log.
fatal(
'Failed to write updated PhotoCalib: %s', str(e))