26 import astropy.units
as u
30 import lsst.pex.config
as pexConfig
38 from lsst.verify
import Job, Measurement
43 from .dataIds
import PerTractCcdDataIdContainer
48 __all__ = [
"JointcalConfig",
"JointcalRunner",
"JointcalTask"]
50 Photometry = collections.namedtuple(
'Photometry', (
'fit',
'model'))
51 Astrometry = collections.namedtuple(
'Astrometry', (
'fit',
'model',
'sky_to_tan_projection'))
56 meas = Measurement(job.metrics[name], value)
57 job.measurements.insert(meas)
61 """Subclass of TaskRunner for jointcalTask
63 jointcalTask.runDataRef() takes a number of arguments, one of which is a list of dataRefs
64 extracted from the command line (whereas most CmdLineTasks' runDataRef methods take
65 single dataRef, are are called repeatedly). This class transforms the processed
66 arguments generated by the ArgumentParser into the arguments expected by
67 Jointcal.runDataRef().
69 See pipeBase.TaskRunner for more information.
75 Return a list of tuples per tract, each containing (dataRefs, kwargs).
77 Jointcal operates on lists of dataRefs simultaneously.
79 kwargs[
'profile_jointcal'] = parsedCmd.profile_jointcal
80 kwargs[
'butler'] = parsedCmd.butler
84 for ref
in parsedCmd.id.refList:
85 refListDict.setdefault(ref.dataId[
"tract"], []).
append(ref)
87 result = [(refListDict[tract], kwargs)
for tract
in sorted(refListDict.keys())]
95 Arguments for Task.runDataRef()
100 if self.doReturnResults is False:
102 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
104 if self.doReturnResults is True:
106 - ``result``: the result of calling jointcal.runDataRef()
107 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
112 dataRefList, kwargs = args
113 butler = kwargs.pop(
'butler')
114 task = self.TaskClass(config=self.config, log=self.log, butler=butler)
117 result = task.runDataRef(dataRefList, **kwargs)
118 exitStatus = result.exitStatus
119 job_path = butler.get(
'verify_job_filename')
120 result.job.write(job_path[0])
121 except Exception
as e:
126 eName =
type(e).__name__
127 tract = dataRefList[0].dataId[
'tract']
128 task.log.fatal(
"Failed processing tract %s, %s: %s", tract, eName, e)
131 kwargs[
'butler'] = butler
132 if self.doReturnResults:
133 return pipeBase.Struct(result=result, exitStatus=exitStatus)
135 return pipeBase.Struct(exitStatus=exitStatus)
139 """Configuration for JointcalTask"""
141 doAstrometry = pexConfig.Field(
142 doc=
"Fit astrometry and write the fitted result.",
146 doPhotometry = pexConfig.Field(
147 doc=
"Fit photometry and write the fitted result.",
151 coaddName = pexConfig.Field(
152 doc=
"Type of coadd, typically deep or goodSeeing",
156 positionErrorPedestal = pexConfig.Field(
157 doc=
"Systematic term to apply to the measured position error (pixels)",
161 photometryErrorPedestal = pexConfig.Field(
162 doc=
"Systematic term to apply to the measured error on flux or magnitude as a "
163 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
168 matchCut = pexConfig.Field(
169 doc=
"Matching radius between fitted and reference stars (arcseconds)",
173 minMeasurements = pexConfig.Field(
174 doc=
"Minimum number of associated measured stars for a fitted star to be included in the fit",
178 minMeasuredStarsPerCcd = pexConfig.Field(
179 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
183 minRefStarsPerCcd = pexConfig.Field(
184 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
188 allowLineSearch = pexConfig.Field(
189 doc=
"Allow a line search during minimization, if it is reasonable for the model"
190 " (models with a significant non-linear component, e.g. constrainedPhotometry).",
194 astrometrySimpleOrder = pexConfig.Field(
195 doc=
"Polynomial order for fitting the simple astrometry model.",
199 astrometryChipOrder = pexConfig.Field(
200 doc=
"Order of the per-chip transform for the constrained astrometry model.",
204 astrometryVisitOrder = pexConfig.Field(
205 doc=
"Order of the per-visit transform for the constrained astrometry model.",
209 useInputWcs = pexConfig.Field(
210 doc=
"Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
214 astrometryModel = pexConfig.ChoiceField(
215 doc=
"Type of model to fit to astrometry",
217 default=
"constrained",
218 allowed={
"simple":
"One polynomial per ccd",
219 "constrained":
"One polynomial per ccd, and one polynomial per visit"}
221 photometryModel = pexConfig.ChoiceField(
222 doc=
"Type of model to fit to photometry",
224 default=
"constrainedMagnitude",
225 allowed={
"simpleFlux":
"One constant zeropoint per ccd and visit, fitting in flux space.",
226 "constrainedFlux":
"Constrained zeropoint per ccd, and one polynomial per visit,"
227 " fitting in flux space.",
228 "simpleMagnitude":
"One constant zeropoint per ccd and visit,"
229 " fitting in magnitude space.",
230 "constrainedMagnitude":
"Constrained zeropoint per ccd, and one polynomial per visit,"
231 " fitting in magnitude space.",
234 applyColorTerms = pexConfig.Field(
235 doc=
"Apply photometric color terms to reference stars?"
236 "Requires that colorterms be set to a ColortermLibrary",
240 colorterms = pexConfig.ConfigField(
241 doc=
"Library of photometric reference catalog name to color term dict.",
242 dtype=ColortermLibrary,
244 photometryVisitOrder = pexConfig.Field(
245 doc=
"Order of the per-visit polynomial transform for the constrained photometry model.",
249 photometryDoRankUpdate = pexConfig.Field(
250 doc=(
"Do the rank update step during minimization. "
251 "Skipping this can help deal with models that are too non-linear."),
255 astrometryDoRankUpdate = pexConfig.Field(
256 doc=(
"Do the rank update step during minimization (should not change the astrometry fit). "
257 "Skipping this can help deal with models that are too non-linear."),
261 outlierRejectSigma = pexConfig.Field(
262 doc=
"How many sigma to reject outliers at during minimization.",
266 maxPhotometrySteps = pexConfig.Field(
267 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
271 maxAstrometrySteps = pexConfig.Field(
272 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
276 astrometryRefObjLoader = pexConfig.ConfigurableField(
277 target=LoadIndexedReferenceObjectsTask,
278 doc=
"Reference object loader for astrometric fit",
280 photometryRefObjLoader = pexConfig.ConfigurableField(
281 target=LoadIndexedReferenceObjectsTask,
282 doc=
"Reference object loader for photometric fit",
284 sourceSelector = sourceSelectorRegistry.makeField(
285 doc=
"How to select sources for cross-matching",
288 astrometryReferenceSelector = pexConfig.ConfigurableField(
289 target=ReferenceSourceSelectorTask,
290 doc=
"How to down-select the loaded astrometry reference catalog.",
292 photometryReferenceSelector = pexConfig.ConfigurableField(
293 target=ReferenceSourceSelectorTask,
294 doc=
"How to down-select the loaded photometry reference catalog.",
296 astrometryReferenceErr = pexConfig.Field(
297 doc=(
"Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
298 "If None, then raise an exception if the reference catalog is missing coordinate errors. "
299 "If specified, overrides any existing `coord_*Err` values."),
304 writeInitMatrix = pexConfig.Field(
306 doc=(
"Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
307 "The output files will be of the form 'astrometry_preinit-mat.txt', in the current directory. "
308 "Note that these files are the dense versions of the matrix, and so may be very large."),
311 writeChi2FilesInitialFinal = pexConfig.Field(
313 doc=
"Write .csv files containing the contributions to chi2 for the initialization and final fit.",
316 writeChi2FilesOuterLoop = pexConfig.Field(
318 doc=
"Write .csv files containing the contributions to chi2 for the outer fit loop.",
321 writeInitialModel = pexConfig.Field(
323 doc=(
"Write the pre-initialization model to text files, for debugging."
324 " Output is written to `initial[Astro|Photo]metryModel.txt` in the current working directory."),
327 debugOutputPath = pexConfig.Field(
330 doc=(
"Path to write debug output files to. Used by "
331 "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
333 sourceFluxType = pexConfig.Field(
335 doc=
"Source flux field to use in source selection and to get fluxes from the catalog.",
342 msg =
"applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
343 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
345 msg = (
"Only doing astrometry, but Colorterms are not applied for astrometry;"
346 "applyColorTerms=True will be ignored.")
360 fluxField = f
"slot_{self.sourceFluxType}Flux_instFlux"
361 self.
sourceSelector[
'science'].signalToNoise.fluxField = fluxField
362 self.
sourceSelector[
'science'].signalToNoise.errField = fluxField +
"Err"
368 badFlags = [
'base_PixelFlags_flag_edge',
'base_PixelFlags_flag_saturated',
369 'base_PixelFlags_flag_interpolatedCenter',
'base_SdssCentroid_flag',
370 'base_PsfFlux_flag',
'base_PixelFlags_flag_suspectCenter']
386 """Write model to outfile."""
387 with open(filename,
"w")
as file:
388 file.write(repr(model))
389 log.info(
"Wrote %s to file: %s", model, filename)
393 """Jointly astrometrically and photometrically calibrate a group of images."""
395 ConfigClass = JointcalConfig
396 RunnerClass = JointcalRunner
397 _DefaultName =
"jointcal"
399 def __init__(self, butler=None, profile_jointcal=False, **kwargs):
401 Instantiate a JointcalTask.
405 butler : `lsst.daf.persistence.Butler`
406 The butler is passed to the refObjLoader constructor in case it is
407 needed. Ignored if the refObjLoader argument provides a loader directly.
408 Used to initialize the astrometry and photometry refObjLoaders.
409 profile_jointcal : `bool`
410 Set to True to profile different stages of this jointcal run.
412 pipeBase.CmdLineTask.__init__(self, **kwargs)
414 self.makeSubtask(
"sourceSelector")
415 if self.config.doAstrometry:
416 self.makeSubtask(
'astrometryRefObjLoader', butler=butler)
417 self.makeSubtask(
"astrometryReferenceSelector")
420 if self.config.doPhotometry:
421 self.makeSubtask(
'photometryRefObjLoader', butler=butler)
422 self.makeSubtask(
"photometryReferenceSelector")
427 self.
job = Job.load_metrics_package(subset=
'jointcal')
432 def _getMetadataName(self):
436 def _makeArgumentParser(cls):
437 """Create an argument parser"""
439 parser.add_argument(
"--profile_jointcal", default=
False, action=
"store_true",
440 help=
"Profile steps of jointcal separately.")
441 parser.add_id_argument(
"--id",
"calexp", help=
"data ID, e.g. --id visit=6789 ccd=0..9",
442 ContainerClass=PerTractCcdDataIdContainer)
445 def _build_ccdImage(self, dataRef, associations, jointcalControl):
447 Extract the necessary things from this dataRef to add a new ccdImage.
451 dataRef : `lsst.daf.persistence.ButlerDataRef`
452 DataRef to extract info from.
453 associations : `lsst.jointcal.Associations`
454 Object to add the info to, to construct a new CcdImage
455 jointcalControl : `jointcal.JointcalControl`
456 Control object for associations management
462 The TAN WCS of this image, read from the calexp
463 (`lsst.afw.geom.SkyWcs`).
465 A key to identify this dataRef by its visit and ccd ids
468 This calexp's filter (`str`).
470 if "visit" in dataRef.dataId.keys():
471 visit = dataRef.dataId[
"visit"]
473 visit = dataRef.getButler().queryMetadata(
"calexp", (
"visit"), dataRef.dataId)[0]
475 src = dataRef.get(
"src", flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS, immediate=
True)
477 visitInfo = dataRef.get(
'calexp_visitInfo')
478 detector = dataRef.get(
'calexp_detector')
479 ccdId = detector.getId()
480 photoCalib = dataRef.get(
'calexp_photoCalib')
481 tanWcs = dataRef.get(
'calexp_wcs')
482 bbox = dataRef.get(
'calexp_bbox')
483 filt = dataRef.get(
'calexp_filter')
484 filterName = filt.getName()
486 goodSrc = self.sourceSelector.
run(src)
488 if len(goodSrc.sourceCat) == 0:
489 self.log.
warn(
"No sources selected in visit %s ccd %s", visit, ccdId)
491 self.log.
info(
"%d sources selected in visit %d ccd %d", len(goodSrc.sourceCat), visit, ccdId)
492 associations.createCcdImage(goodSrc.sourceCat,
503 Result = collections.namedtuple(
'Result_from_build_CcdImage', (
'wcs',
'key',
'filter'))
504 Key = collections.namedtuple(
'Key', (
'visit',
'ccd'))
505 return Result(tanWcs, Key(visit, ccdId), filterName)
507 def _getDebugPath(self, filename):
508 """Constructs a path to filename using the configured debug path.
510 return os.path.join(self.config.debugOutputPath, filename)
515 Jointly calibrate the astrometry and photometry across a set of images.
519 dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef`
520 List of data references to the exposures to be fit.
521 profile_jointcal : `bool`
522 Profile the individual steps of jointcal.
526 result : `lsst.pipe.base.Struct`
527 Struct of metadata from the fit, containing:
530 The provided data references that were fit (with updated WCSs)
532 The original WCS from each dataRef
534 Dictionary of internally-computed metrics for testing/validation.
536 if len(dataRefs) == 0:
537 raise ValueError(
'Need a non-empty list of data references!')
541 sourceFluxField =
"slot_%sFlux" % (self.config.sourceFluxType,)
545 visit_ccd_to_dataRef = {}
548 load_cat_prof_file =
'jointcal_build_ccdImage.prof' if profile_jointcal
else ''
549 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
552 camera = dataRefs[0].get(
'camera', immediate=
True)
556 oldWcsList.append(result.wcs)
557 visit_ccd_to_dataRef[result.key] = ref
558 filters.append(result.filter)
559 filters = collections.Counter(filters)
561 associations.computeCommonTangentPoint()
563 boundingCircle = associations.computeBoundingCircle()
565 radius =
lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
567 self.log.
info(f
"Data has center={center} with radius={radius.asDegrees()} degrees.")
570 defaultFilter = filters.most_common(1)[0][0]
571 self.log.
debug(
"Using %s band for reference flux", defaultFilter)
574 tract = dataRefs[0].dataId[
'tract']
576 if self.config.doAstrometry:
580 referenceSelector=self.astrometryReferenceSelector,
582 profile_jointcal=profile_jointcal,
588 if self.config.doPhotometry:
592 referenceSelector=self.photometryReferenceSelector,
594 profile_jointcal=profile_jointcal,
597 reject_bad_fluxes=
True)
602 return pipeBase.Struct(dataRefs=dataRefs,
603 oldWcsList=oldWcsList,
607 defaultFilter=defaultFilter,
608 exitStatus=exitStatus)
610 def _get_refcat_coordinate_error_override(self, refCat, name):
611 """Check whether we should override the refcat coordinate errors, and
612 return the overridden error if necessary.
616 refCat : `lsst.afw.table.SimpleCatalog`
617 The reference catalog to check for a ``coord_raErr`` field.
619 Whether we are doing "astrometry" or "photometry".
623 refCoordErr : `float`
624 The refcat coordinate error to use, or NaN if we are not overriding
629 lsst.pex.config.FieldValidationError
630 Raised if the refcat does not contain coordinate errors and
631 ``config.astrometryReferenceErr`` is not set.
635 if name.lower() ==
"photometry":
636 if 'coord_raErr' not in refCat.schema:
641 if self.config.astrometryReferenceErr
is None and 'coord_raErr' not in refCat.schema:
642 msg = (
"Reference catalog does not contain coordinate errors, "
643 "and config.astrometryReferenceErr not supplied.")
644 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
648 if self.config.astrometryReferenceErr
is not None and 'coord_raErr' in refCat.schema:
649 self.log.
warn(
"Overriding reference catalog coordinate errors with %f/coordinate [mas]",
650 self.config.astrometryReferenceErr)
652 if self.config.astrometryReferenceErr
is None:
655 return self.config.astrometryReferenceErr
657 def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
659 tract="", profile_jointcal=False, match_cut=3.0,
660 reject_bad_fluxes=False, *,
661 name="", refObjLoader=None, referenceSelector=None,
663 """Load reference catalog, perform the fit, and return the result.
667 associations : `lsst.jointcal.Associations`
668 The star/reference star associations to fit.
669 defaultFilter : `str`
670 filter to load from reference catalog.
671 center : `lsst.geom.SpherePoint`
672 ICRS center of field to load from reference catalog.
673 radius : `lsst.geom.Angle`
674 On-sky radius to load from reference catalog.
676 Name of thing being fit: "astrometry" or "photometry".
677 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
678 Reference object loader to use to load a reference catalog.
679 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
680 Selector to use to pick objects from the loaded reference catalog.
681 fit_function : callable
682 Function to call to perform fit (takes Associations object).
683 filters : `list` [`str`], optional
684 List of filters to load from the reference catalog.
685 tract : `str`, optional
686 Name of tract currently being fit.
687 profile_jointcal : `bool`, optional
688 Separately profile the fitting step.
689 match_cut : `float`, optional
690 Radius in arcseconds to find cross-catalog matches to during
691 associations.associateCatalogs.
692 reject_bad_fluxes : `bool`, optional
693 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
697 result : `Photometry` or `Astrometry`
698 Result of `fit_function()`
700 self.log.
info(
"====== Now processing %s...", name)
703 associations.associateCatalogs(match_cut)
705 associations.fittedStarListSize())
707 applyColorterms =
False if name.lower() ==
"astrometry" else self.config.applyColorTerms
709 center, radius, defaultFilter,
710 applyColorterms=applyColorterms)
713 associations.collectRefStars(refCat,
714 self.config.matchCut*lsst.geom.arcseconds,
716 refCoordinateErr=refCoordErr,
717 rejectBadFluxes=reject_bad_fluxes)
719 associations.refStarListSize())
721 associations.prepareFittedStars(self.config.minMeasurements)
725 associations.nFittedStarsWithAssociatedRefStar())
727 associations.fittedStarListSize())
729 associations.nCcdImagesValidForFit())
731 load_cat_prof_file =
'jointcal_fit_%s.prof'%name
if profile_jointcal
else ''
732 dataName =
"{}_{}".
format(tract, defaultFilter)
733 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
734 result = fit_function(associations, dataName)
737 if self.config.writeChi2FilesInitialFinal:
738 baseName = self.
_getDebugPath(f
"{name}_final_chi2-{dataName}")
739 result.fit.saveChi2Contributions(baseName+
"{type}")
740 self.log.
info(
"Wrote chi2 contributions files: %s", baseName)
744 def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName,
745 applyColorterms=False):
746 """Load the necessary reference catalog sources, convert fluxes to
747 correct units, and apply color term corrections if requested.
751 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
752 The reference catalog loader to use to get the data.
753 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
754 Source selector to apply to loaded reference catalog.
755 center : `lsst.geom.SpherePoint`
756 The center around which to load sources.
757 radius : `lsst.geom.Angle`
758 The radius around ``center`` to load sources in.
760 The name of the camera filter to load fluxes for.
761 applyColorterms : `bool`
762 Apply colorterm corrections to the refcat for ``filterName``?
766 refCat : `lsst.afw.table.SimpleCatalog`
767 The loaded reference catalog.
769 The name of the reference catalog flux field appropriate for ``filterName``.
771 skyCircle = refObjLoader.loadSkyCircle(center,
775 selected = referenceSelector.run(skyCircle.refCat)
777 if not selected.sourceCat.isContiguous():
778 refCat = selected.sourceCat.copy(deep=
True)
780 refCat = selected.sourceCat
784 refCatName = refObjLoader.ref_dataset_name
785 except AttributeError:
787 raise RuntimeError(
"Cannot perform colorterm corrections with a.net refcats.")
788 self.log.
info(
"Applying color terms for filterName=%r reference catalog=%s",
789 filterName, refCatName)
790 colorterm = self.config.colorterms.getColorterm(
791 filterName=filterName, photoCatName=refCatName, doRaise=
True)
793 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat, filterName)
794 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
798 return refCat, skyCircle.fluxField
800 def _check_star_lists(self, associations, name):
802 if associations.nCcdImagesValidForFit() == 0:
803 raise RuntimeError(
'No images in the ccdImageList!')
804 if associations.fittedStarListSize() == 0:
805 raise RuntimeError(
'No stars in the {} fittedStarList!'.
format(name))
806 if associations.refStarListSize() == 0:
807 raise RuntimeError(
'No stars in the {} reference star list!'.
format(name))
809 def _logChi2AndValidate(self, associations, fit, model, chi2Label="Model",
811 """Compute chi2, log it, validate the model, and return chi2.
815 associations : `lsst.jointcal.Associations`
816 The star/reference star associations to fit.
817 fit : `lsst.jointcal.FitterBase`
818 The fitter to use for minimization.
819 model : `lsst.jointcal.Model`
821 chi2Label : str, optional
822 Label to describe the chi2 (e.g. "Initialized", "Final").
823 writeChi2Name : `str`, optional
824 Filename prefix to write the chi2 contributions to.
825 Do not supply an extension: an appropriate one will be added.
829 chi2: `lsst.jointcal.Chi2Accumulator`
830 The chi2 object for the current fitter and model.
835 Raised if chi2 is infinite or NaN.
837 Raised if the model is not valid.
839 if writeChi2Name
is not None:
841 fit.saveChi2Contributions(fullpath+
"{type}")
842 self.log.
info(
"Wrote chi2 contributions files: %s", fullpath)
844 chi2 = fit.computeChi2()
845 self.log.
info(
"%s %s", chi2Label, chi2)
847 if not np.isfinite(chi2.chi2):
848 raise FloatingPointError(f
'{chi2Label} chi2 is invalid: {chi2}')
849 if not model.validate(associations.getCcdImageList(), chi2.ndof):
850 raise ValueError(
"Model is not valid: check log messages for warnings.")
853 def _fit_photometry(self, associations, dataName=None):
855 Fit the photometric data.
859 associations : `lsst.jointcal.Associations`
860 The star/reference star associations to fit.
862 Name of the data being processed (e.g. "1234_HSC-Y"), for
863 identifying debugging files.
867 fit_result : `namedtuple`
868 fit : `lsst.jointcal.PhotometryFit`
869 The photometric fitter used to perform the fit.
870 model : `lsst.jointcal.PhotometryModel`
871 The photometric model that was fit.
873 self.log.
info(
"=== Starting photometric fitting...")
876 if self.config.photometryModel ==
"constrainedFlux":
879 visitOrder=self.config.photometryVisitOrder,
880 errorPedestal=self.config.photometryErrorPedestal)
882 doLineSearch = self.config.allowLineSearch
883 elif self.config.photometryModel ==
"constrainedMagnitude":
886 visitOrder=self.config.photometryVisitOrder,
887 errorPedestal=self.config.photometryErrorPedestal)
889 doLineSearch = self.config.allowLineSearch
890 elif self.config.photometryModel ==
"simpleFlux":
892 errorPedestal=self.config.photometryErrorPedestal)
894 elif self.config.photometryModel ==
"simpleMagnitude":
896 errorPedestal=self.config.photometryErrorPedestal)
902 if self.config.writeChi2FilesInitialFinal:
903 baseName = f
"photometry_initial_chi2-{dataName}"
906 if self.config.writeInitialModel:
911 def getChi2Name(whatToFit):
912 if self.config.writeChi2FilesOuterLoop:
913 return f
"photometry_init-%s_chi2-{dataName}" % whatToFit
919 dumpMatrixFile = self.
_getDebugPath(
"photometry_preinit")
if self.config.writeInitMatrix
else ""
920 if self.config.photometryModel.startswith(
"constrained"):
923 fit.minimize(
"ModelVisit", dumpMatrixFile=dumpMatrixFile)
924 self.
_logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name(
"ModelVisit"))
927 fit.minimize(
"Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
930 fit.minimize(
"Fluxes")
933 fit.minimize(
"Model Fluxes", doLineSearch=doLineSearch)
935 writeChi2Name=getChi2Name(
"ModelFluxes"))
937 model.freezeErrorTransform()
938 self.log.
debug(
"Photometry error scales are frozen.")
942 self.config.maxPhotometrySteps,
945 doRankUpdate=self.config.photometryDoRankUpdate,
946 doLineSearch=doLineSearch,
953 def _fit_astrometry(self, associations, dataName=None):
955 Fit the astrometric data.
959 associations : `lsst.jointcal.Associations`
960 The star/reference star associations to fit.
962 Name of the data being processed (e.g. "1234_HSC-Y"), for
963 identifying debugging files.
967 fit_result : `namedtuple`
968 fit : `lsst.jointcal.AstrometryFit`
969 The astrometric fitter used to perform the fit.
970 model : `lsst.jointcal.AstrometryModel`
971 The astrometric model that was fit.
972 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
973 The model for the sky to tangent plane projection that was used in the fit.
976 self.log.
info(
"=== Starting astrometric fitting...")
978 associations.deprojectFittedStars()
985 if self.config.astrometryModel ==
"constrained":
987 sky_to_tan_projection,
988 chipOrder=self.config.astrometryChipOrder,
989 visitOrder=self.config.astrometryVisitOrder)
990 elif self.config.astrometryModel ==
"simple":
992 sky_to_tan_projection,
993 self.config.useInputWcs,
995 order=self.config.astrometrySimpleOrder)
1000 if self.config.writeChi2FilesInitialFinal:
1001 baseName = f
"astrometry_initial_chi2-{dataName}"
1004 if self.config.writeInitialModel:
1009 def getChi2Name(whatToFit):
1010 if self.config.writeChi2FilesOuterLoop:
1011 return f
"astrometry_init-%s_chi2-{dataName}" % whatToFit
1015 dumpMatrixFile = self.
_getDebugPath(
"astrometry_preinit")
if self.config.writeInitMatrix
else ""
1018 if self.config.astrometryModel ==
"constrained":
1019 fit.minimize(
"DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1020 self.
_logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name(
"DistortionsVisit"))
1023 fit.minimize(
"Distortions", dumpMatrixFile=dumpMatrixFile)
1024 self.
_logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name(
"Distortions"))
1026 fit.minimize(
"Positions")
1027 self.
_logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name(
"Positions"))
1029 fit.minimize(
"Distortions Positions")
1031 writeChi2Name=getChi2Name(
"DistortionsPositions"))
1035 self.config.maxAstrometrySteps,
1037 "Distortions Positions",
1038 doRankUpdate=self.config.astrometryDoRankUpdate,
1044 return Astrometry(fit, model, sky_to_tan_projection)
1046 def _check_stars(self, associations):
1047 """Count measured and reference stars per ccd and warn/log them."""
1048 for ccdImage
in associations.getCcdImageList():
1049 nMeasuredStars, nRefStars = ccdImage.countStars()
1050 self.log.
debug(
"ccdImage %s has %s measured and %s reference stars",
1051 ccdImage.getName(), nMeasuredStars, nRefStars)
1052 if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1053 self.log.
warn(
"ccdImage %s has only %s measuredStars (desired %s)",
1054 ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1055 if nRefStars < self.config.minRefStarsPerCcd:
1056 self.log.
warn(
"ccdImage %s has only %s RefStars (desired %s)",
1057 ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1059 def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1062 doLineSearch=False):
1063 """Run fitter.minimize up to max_steps times, returning the final chi2.
1067 associations : `lsst.jointcal.Associations`
1068 The star/reference star associations to fit.
1069 fitter : `lsst.jointcal.FitterBase`
1070 The fitter to use for minimization.
1072 Maximum number of steps to run outlier rejection before declaring
1073 convergence failure.
1074 name : {'photometry' or 'astrometry'}
1075 What type of data are we fitting (for logs and debugging files).
1077 Passed to ``fitter.minimize()`` to define the parameters to fit.
1078 dataName : `str`, optional
1079 Descriptive name for this dataset (e.g. tract and filter),
1081 doRankUpdate : `bool`, optional
1082 Do an Eigen rank update during minimization, or recompute the full
1083 matrix and gradient?
1084 doLineSearch : `bool`, optional
1085 Do a line search for the optimum step during minimization?
1089 chi2: `lsst.jointcal.Chi2Statistic`
1090 The final chi2 after the fit converges, or is forced to end.
1095 Raised if the fitter fails with a non-finite value.
1097 Raised if the fitter fails for some other reason;
1098 log messages will provide further details.
1100 dumpMatrixFile = self.
_getDebugPath(f
"{name}_postinit")
if self.config.writeInitMatrix
else ""
1101 for i
in range(max_steps):
1102 if self.config.writeChi2FilesOuterLoop:
1103 writeChi2Name = f
"{name}_iterate_{i}_chi2-{dataName}"
1105 writeChi2Name =
None
1106 result = fitter.minimize(whatToFit,
1107 self.config.outlierRejectSigma,
1108 doRankUpdate=doRankUpdate,
1109 doLineSearch=doLineSearch,
1110 dumpMatrixFile=dumpMatrixFile)
1113 writeChi2Name=writeChi2Name)
1115 if result == MinimizeResult.Converged:
1117 self.log.
debug(
"fit has converged - no more outliers - redo minimization "
1118 "one more time in case we have lost accuracy in rank update.")
1120 result = fitter.minimize(whatToFit, self.config.outlierRejectSigma)
1124 if chi2.chi2/chi2.ndof >= 4.0:
1125 self.log.
error(
"Potentially bad fit: High chi-squared/ndof.")
1128 elif result == MinimizeResult.Chi2Increased:
1129 self.log.
warn(
"still some outliers but chi2 increases - retry")
1130 elif result == MinimizeResult.NonFinite:
1133 fitter.saveChi2Contributions(filename+
"{type}")
1134 msg =
"Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1135 raise FloatingPointError(msg.format(filename))
1136 elif result == MinimizeResult.Failed:
1137 raise RuntimeError(
"Chi2 minimization failure, cannot complete fit.")
1139 raise RuntimeError(
"Unxepected return code from minimize().")
1141 self.log.
error(
"%s failed to converge after %d steps"%(name, max_steps))
1145 def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1147 Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1151 associations : `lsst.jointcal.Associations`
1152 The star/reference star associations to fit.
1153 model : `lsst.jointcal.AstrometryModel`
1154 The astrometric model that was fit.
1155 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1156 Dict of ccdImage identifiers to dataRefs that were fit.
1159 ccdImageList = associations.getCcdImageList()
1160 for ccdImage
in ccdImageList:
1162 ccd = ccdImage.ccdId
1163 visit = ccdImage.visit
1164 dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1165 self.log.
info(
"Updating WCS for visit: %d, ccd: %d", visit, ccd)
1166 skyWcs = model.makeSkyWcs(ccdImage)
1168 dataRef.put(skyWcs,
'jointcal_wcs')
1170 self.log.
fatal(
'Failed to write updated Wcs: %s', str(e))
1173 def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1175 Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1179 associations : `lsst.jointcal.Associations`
1180 The star/reference star associations to fit.
1181 model : `lsst.jointcal.PhotometryModel`
1182 The photoometric model that was fit.
1183 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1184 Dict of ccdImage identifiers to dataRefs that were fit.
1187 ccdImageList = associations.getCcdImageList()
1188 for ccdImage
in ccdImageList:
1190 ccd = ccdImage.ccdId
1191 visit = ccdImage.visit
1192 dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1193 self.log.
info(
"Updating PhotoCalib for visit: %d, ccd: %d", visit, ccd)
1194 photoCalib = model.toPhotoCalib(ccdImage)
1196 dataRef.put(photoCalib,
'jointcal_photoCalib')
1198 self.log.
fatal(
'Failed to write updated PhotoCalib: %s', str(e))