28 import astropy.units
as u
40 from lsst.pipe.tasks.colorterms
import ColortermLibrary
41 from lsst.verify
import Job, Measurement
46 from .dataIds
import PerTractCcdDataIdContainer
51 __all__ = [
"JointcalConfig",
"JointcalRunner",
"JointcalTask"]
53 Photometry = collections.namedtuple(
'Photometry', (
'fit',
'model'))
54 Astrometry = collections.namedtuple(
'Astrometry', (
'fit',
'model',
'sky_to_tan_projection'))
59 meas = Measurement(job.metrics[name], value)
60 job.measurements.insert(meas)
64 """Subclass of TaskRunner for jointcalTask
66 jointcalTask.runDataRef() takes a number of arguments, one of which is a list of dataRefs
67 extracted from the command line (whereas most CmdLineTasks' runDataRef methods take
68 single dataRef, are are called repeatedly). This class transforms the processed
69 arguments generated by the ArgumentParser into the arguments expected by
70 Jointcal.runDataRef().
72 See pipeBase.TaskRunner for more information.
78 Return a list of tuples per tract, each containing (dataRefs, kwargs).
80 Jointcal operates on lists of dataRefs simultaneously.
82 kwargs[
'profile_jointcal'] = parsedCmd.profile_jointcal
83 kwargs[
'butler'] = parsedCmd.butler
87 for ref
in parsedCmd.id.refList:
88 refListDict.setdefault(ref.dataId[
"tract"], []).append(ref)
90 result = [(refListDict[tract], kwargs)
for tract
in sorted(refListDict.keys())]
98 Arguments for Task.runDataRef()
103 if self.doReturnResults is False:
105 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
107 if self.doReturnResults is True:
109 - ``result``: the result of calling jointcal.runDataRef()
110 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise.
115 dataRefList, kwargs = args
116 butler = kwargs.pop(
'butler')
117 task = self.TaskClass(config=self.config, log=self.log, butler=butler)
120 result = task.runDataRef(dataRefList, **kwargs)
121 exitStatus = result.exitStatus
122 job_path = butler.get(
'verify_job_filename')
123 result.job.write(job_path[0])
124 except Exception
as e:
129 eName = type(e).__name__
130 tract = dataRefList[0].dataId[
'tract']
131 task.log.fatal(
"Failed processing tract %s, %s: %s", tract, eName, e)
134 kwargs[
'butler'] = butler
135 if self.doReturnResults:
136 return pipeBase.Struct(result=result, exitStatus=exitStatus)
138 return pipeBase.Struct(exitStatus=exitStatus)
142 """Configuration for JointcalTask"""
144 doAstrometry = pexConfig.Field(
145 doc=
"Fit astrometry and write the fitted result.",
149 doPhotometry = pexConfig.Field(
150 doc=
"Fit photometry and write the fitted result.",
154 coaddName = pexConfig.Field(
155 doc=
"Type of coadd, typically deep or goodSeeing",
159 positionErrorPedestal = pexConfig.Field(
160 doc=
"Systematic term to apply to the measured position error (pixels)",
164 photometryErrorPedestal = pexConfig.Field(
165 doc=
"Systematic term to apply to the measured error on flux or magnitude as a "
166 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
171 matchCut = pexConfig.Field(
172 doc=
"Matching radius between fitted and reference stars (arcseconds)",
176 minMeasurements = pexConfig.Field(
177 doc=
"Minimum number of associated measured stars for a fitted star to be included in the fit",
181 minMeasuredStarsPerCcd = pexConfig.Field(
182 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
186 minRefStarsPerCcd = pexConfig.Field(
187 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
191 allowLineSearch = pexConfig.Field(
192 doc=
"Allow a line search during minimization, if it is reasonable for the model"
193 " (models with a significant non-linear component, e.g. constrainedPhotometry).",
197 astrometrySimpleOrder = pexConfig.Field(
198 doc=
"Polynomial order for fitting the simple astrometry model.",
202 astrometryChipOrder = pexConfig.Field(
203 doc=
"Order of the per-chip transform for the constrained astrometry model.",
207 astrometryVisitOrder = pexConfig.Field(
208 doc=
"Order of the per-visit transform for the constrained astrometry model.",
212 useInputWcs = pexConfig.Field(
213 doc=
"Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
217 astrometryModel = pexConfig.ChoiceField(
218 doc=
"Type of model to fit to astrometry",
220 default=
"constrained",
221 allowed={
"simple":
"One polynomial per ccd",
222 "constrained":
"One polynomial per ccd, and one polynomial per visit"}
224 photometryModel = pexConfig.ChoiceField(
225 doc=
"Type of model to fit to photometry",
227 default=
"constrainedMagnitude",
228 allowed={
"simpleFlux":
"One constant zeropoint per ccd and visit, fitting in flux space.",
229 "constrainedFlux":
"Constrained zeropoint per ccd, and one polynomial per visit,"
230 " fitting in flux space.",
231 "simpleMagnitude":
"One constant zeropoint per ccd and visit,"
232 " fitting in magnitude space.",
233 "constrainedMagnitude":
"Constrained zeropoint per ccd, and one polynomial per visit,"
234 " fitting in magnitude space.",
237 applyColorTerms = pexConfig.Field(
238 doc=
"Apply photometric color terms to reference stars?"
239 "Requires that colorterms be set to a ColortermLibrary",
243 colorterms = pexConfig.ConfigField(
244 doc=
"Library of photometric reference catalog name to color term dict.",
245 dtype=ColortermLibrary,
247 photometryVisitOrder = pexConfig.Field(
248 doc=
"Order of the per-visit polynomial transform for the constrained photometry model.",
252 photometryDoRankUpdate = pexConfig.Field(
253 doc=(
"Do the rank update step during minimization. "
254 "Skipping this can help deal with models that are too non-linear."),
258 astrometryDoRankUpdate = pexConfig.Field(
259 doc=(
"Do the rank update step during minimization (should not change the astrometry fit). "
260 "Skipping this can help deal with models that are too non-linear."),
264 outlierRejectSigma = pexConfig.Field(
265 doc=
"How many sigma to reject outliers at during minimization.",
269 maxPhotometrySteps = pexConfig.Field(
270 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
274 maxAstrometrySteps = pexConfig.Field(
275 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
279 astrometryRefObjLoader = pexConfig.ConfigurableField(
280 target=LoadIndexedReferenceObjectsTask,
281 doc=
"Reference object loader for astrometric fit",
283 photometryRefObjLoader = pexConfig.ConfigurableField(
284 target=LoadIndexedReferenceObjectsTask,
285 doc=
"Reference object loader for photometric fit",
287 sourceSelector = sourceSelectorRegistry.makeField(
288 doc=
"How to select sources for cross-matching",
291 astrometryReferenceSelector = pexConfig.ConfigurableField(
292 target=ReferenceSourceSelectorTask,
293 doc=
"How to down-select the loaded astrometry reference catalog.",
295 photometryReferenceSelector = pexConfig.ConfigurableField(
296 target=ReferenceSourceSelectorTask,
297 doc=
"How to down-select the loaded photometry reference catalog.",
299 astrometryReferenceErr = pexConfig.Field(
300 doc=(
"Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
301 "If None, then raise an exception if the reference catalog is missing coordinate errors. "
302 "If specified, overrides any existing `coord_*Err` values."),
307 writeInitMatrix = pexConfig.Field(
309 doc=(
"Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
310 "The output files will be of the form 'astrometry_preinit-mat.txt', in the current directory. "
311 "Note that these files are the dense versions of the matrix, and so may be very large."),
314 writeChi2FilesInitialFinal = pexConfig.Field(
316 doc=
"Write .csv files containing the contributions to chi2 for the initialization and final fit.",
319 writeChi2FilesOuterLoop = pexConfig.Field(
321 doc=
"Write .csv files containing the contributions to chi2 for the outer fit loop.",
324 writeInitialModel = pexConfig.Field(
326 doc=(
"Write the pre-initialization model to text files, for debugging."
327 " Output is written to `initial[Astro|Photo]metryModel.txt` in the current working directory."),
330 debugOutputPath = pexConfig.Field(
333 doc=(
"Path to write debug output files to. Used by "
334 "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
336 sourceFluxType = pexConfig.Field(
338 doc=
"Source flux field to use in source selection and to get fluxes from the catalog.",
345 msg =
"applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
346 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
348 msg = (
"Only doing astrometry, but Colorterms are not applied for astrometry;"
349 "applyColorTerms=True will be ignored.")
358 self.
sourceSelectorsourceSelector[
'science'].doSignalToNoise =
True
361 self.
sourceSelectorsourceSelector[
'science'].signalToNoise.minimum = 10.
363 fluxField = f
"slot_{self.sourceFluxType}Flux_instFlux"
364 self.
sourceSelectorsourceSelector[
'science'].signalToNoise.fluxField = fluxField
365 self.
sourceSelectorsourceSelector[
'science'].signalToNoise.errField = fluxField +
"Err"
371 badFlags = [
'base_PixelFlags_flag_edge',
'base_PixelFlags_flag_saturated',
372 'base_PixelFlags_flag_interpolatedCenter',
'base_SdssCentroid_flag',
373 'base_PsfFlux_flag',
'base_PixelFlags_flag_suspectCenter']
385 """Write model to outfile."""
386 with open(filename,
"w")
as file:
387 file.write(repr(model))
388 log.info(
"Wrote %s to file: %s", model, filename)
391 @dataclasses.dataclass
393 """The input data jointcal needs for each detector/visit."""
395 """The visit identifier of this exposure."""
397 """The catalog derived from this exposure."""
399 """The VisitInfo of this exposure."""
401 """The detector of this exposure."""
403 """The photometric calibration of this exposure."""
405 """The WCS of this exposure."""
407 """The bounding box of this exposure."""
409 """The filter of this exposure."""
413 """Astrometricly and photometricly calibrate across multiple visits of the
418 butler : `lsst.daf.persistence.Butler`
419 The butler is passed to the refObjLoader constructor in case it is
420 needed. Ignored if the refObjLoader argument provides a loader directly.
421 Used to initialize the astrometry and photometry refObjLoaders.
422 profile_jointcal : `bool`
423 Set to True to profile different stages of this jointcal run.
426 ConfigClass = JointcalConfig
427 RunnerClass = JointcalRunner
428 _DefaultName =
"jointcal"
430 def __init__(self, butler=None, profile_jointcal=False, **kwargs):
433 self.makeSubtask(
"sourceSelector")
434 if self.config.doAstrometry:
435 self.makeSubtask(
'astrometryRefObjLoader', butler=butler)
436 self.makeSubtask(
"astrometryReferenceSelector")
439 if self.config.doPhotometry:
440 self.makeSubtask(
'photometryRefObjLoader', butler=butler)
441 self.makeSubtask(
"photometryReferenceSelector")
446 self.
jobjob = Job.load_metrics_package(subset=
'jointcal')
451 def _getMetadataName(self):
455 def _makeArgumentParser(cls):
456 """Create an argument parser"""
457 parser = pipeBase.ArgumentParser(name=cls.
_DefaultName_DefaultName)
458 parser.add_argument(
"--profile_jointcal", default=
False, action=
"store_true",
459 help=
"Profile steps of jointcal separately.")
460 parser.add_id_argument(
"--id",
"calexp", help=
"data ID, e.g. --id visit=6789 ccd=0..9",
461 ContainerClass=PerTractCcdDataIdContainer)
464 def _build_ccdImage(self, data, associations, jointcalControl):
466 Extract the necessary things from this dataRef to add a new ccdImage.
470 data : `JointcalInputData`
471 The loaded input data.
472 associations : `lsst.jointcal.Associations`
473 Object to add the info to, to construct a new CcdImage
474 jointcalControl : `jointcal.JointcalControl`
475 Control object for associations management
481 The TAN WCS of this image, read from the calexp
482 (`lsst.afw.geom.SkyWcs`).
484 A key to identify this dataRef by its visit and ccd ids
487 This calexp's filter band (`str`) (used to e.g. load refcats)
489 goodSrc = self.sourceSelector.run(data.catalog)
491 if len(goodSrc.sourceCat) == 0:
492 self.log.warn(
"No sources selected in visit %s ccd %s", data.visit, data.detector.getId())
494 self.log.info(
"%d sources selected in visit %d ccd %d", len(goodSrc.sourceCat),
496 data.detector.getId())
497 associations.createCcdImage(goodSrc.sourceCat,
501 data.filter.physicalLabel,
505 data.detector.getId(),
508 Result = collections.namedtuple(
'Result_from_build_CcdImage', (
'wcs',
'key',
'band'))
509 Key = collections.namedtuple(
'Key', (
'visit',
'ccd'))
510 return Result(data.wcs, Key(data.visit, data.detector.getId()), data.filter.bandLabel)
512 def _readDataId(self, butler, dataId):
513 """Read all of the data for one dataId from the butler. (gen2 version)"""
515 if "visit" in dataId.keys():
516 visit = dataId[
"visit"]
518 visit = butler.getButler().queryMetadata(
"calexp", (
"visit"), butler.dataId)[0]
520 catalog = butler.get(
'src',
521 flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS,
525 visitInfo=butler.get(
'calexp_visitInfo', dataId=dataId),
526 detector=butler.get(
'calexp_detector', dataId=dataId),
527 photoCalib=butler.get(
'calexp_photoCalib', dataId=dataId),
528 wcs=butler.get(
'calexp_wcs', dataId=dataId),
529 bbox=butler.get(
'calexp_bbox', dataId=dataId),
530 filter=butler.get(
'calexp_filterLabel', dataId=dataId))
532 def loadData(self, dataRefs, associations, jointcalControl, profile_jointcal=False):
533 """Read the data that jointcal needs to run. (Gen2 version)"""
534 visit_ccd_to_dataRef = {}
537 load_cat_prof_file =
'jointcal_loadData.prof' if profile_jointcal
else ''
538 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
540 camera = dataRefs[0].get(
'camera', immediate=
True)
542 for dataRef
in dataRefs:
543 data = self.
_readDataId_readDataId(dataRef.getButler(), dataRef.dataId)
544 result = self.
_build_ccdImage_build_ccdImage(data, associations, jointcalControl)
545 oldWcsList.append(result.wcs)
546 visit_ccd_to_dataRef[result.key] = dataRef
547 bands.append(result.band)
548 bands = collections.Counter(bands)
550 return oldWcsList, bands, visit_ccd_to_dataRef
552 def _getDebugPath(self, filename):
553 """Constructs a path to filename using the configured debug path.
555 return os.path.join(self.config.debugOutputPath, filename)
557 def _prep_sky(self, associations, bands):
558 """Prepare on-sky and other data that must be computed after data has
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 defaultBand = bands.most_common(1)[0][0]
571 self.log.debug(
"Using '%s' filter band for reference flux", defaultBand)
573 return boundingCircle, center, radius, defaultBand
578 Jointly calibrate the astrometry and photometry across a set of images.
580 NOTE: this is for gen2 middleware only.
584 dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef`
585 List of data references to the exposures to be fit.
586 profile_jointcal : `bool`
587 Profile the individual steps of jointcal.
591 result : `lsst.pipe.base.Struct`
592 Struct of metadata from the fit, containing:
595 The provided data references that were fit (with updated WCSs)
597 The original WCS from each dataRef
599 Dictionary of internally-computed metrics for testing/validation.
601 if len(dataRefs) == 0:
602 raise ValueError(
'Need a non-empty list of data references!')
606 sourceFluxField =
"slot_%sFlux" % (self.config.sourceFluxType,)
610 oldWcsList, bands, visit_ccd_to_dataRef = self.
loadDataloadData(dataRefs,
613 profile_jointcal=profile_jointcal)
615 boundingCircle, center, radius, defaultBand = self.
_prep_sky_prep_sky(associations, bands)
617 tract = dataRefs[0].dataId[
'tract']
619 if self.config.doAstrometry:
623 referenceSelector=self.astrometryReferenceSelector,
625 profile_jointcal=profile_jointcal,
631 if self.config.doPhotometry:
635 referenceSelector=self.photometryReferenceSelector,
637 profile_jointcal=profile_jointcal,
639 reject_bad_fluxes=
True)
644 return pipeBase.Struct(dataRefs=dataRefs,
645 oldWcsList=oldWcsList,
649 defaultBand=defaultBand,
650 exitStatus=exitStatus)
652 def _get_refcat_coordinate_error_override(self, refCat, name):
653 """Check whether we should override the refcat coordinate errors, and
654 return the overridden error if necessary.
658 refCat : `lsst.afw.table.SimpleCatalog`
659 The reference catalog to check for a ``coord_raErr`` field.
661 Whether we are doing "astrometry" or "photometry".
665 refCoordErr : `float`
666 The refcat coordinate error to use, or NaN if we are not overriding
671 lsst.pex.config.FieldValidationError
672 Raised if the refcat does not contain coordinate errors and
673 ``config.astrometryReferenceErr`` is not set.
677 if name.lower() ==
"photometry":
678 if 'coord_raErr' not in refCat.schema:
683 if self.config.astrometryReferenceErr
is None and 'coord_raErr' not in refCat.schema:
684 msg = (
"Reference catalog does not contain coordinate errors, "
685 "and config.astrometryReferenceErr not supplied.")
686 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
690 if self.config.astrometryReferenceErr
is not None and 'coord_raErr' in refCat.schema:
691 self.log.warn(
"Overriding reference catalog coordinate errors with %f/coordinate [mas]",
692 self.config.astrometryReferenceErr)
694 if self.config.astrometryReferenceErr
is None:
697 return self.config.astrometryReferenceErr
699 def _compute_proper_motion_epoch(self, ccdImageList):
700 """Return the proper motion correction epoch of the provided images.
704 ccdImageList : `list` [`lsst.jointcal.CcdImage`]
705 The images to compute the appropriate epoch for.
709 epoch : `astropy.time.Time`
710 The date to use for proper motion corrections.
712 mjds = [ccdImage.getMjd()
for ccdImage
in ccdImageList]
713 return astropy.time.Time(np.mean(mjds), format=
'mjd', scale=
"tai")
715 def _do_load_refcat_and_fit(self, associations, defaultBand, center, radius,
716 tract="", profile_jointcal=False, match_cut=3.0,
717 reject_bad_fluxes=False, *,
718 name="", refObjLoader=None, referenceSelector=None,
720 """Load reference catalog, perform the fit, and return the result.
724 associations : `lsst.jointcal.Associations`
725 The star/reference star associations to fit.
727 filter to load from reference catalog.
728 center : `lsst.geom.SpherePoint`
729 ICRS center of field to load from reference catalog.
730 radius : `lsst.geom.Angle`
731 On-sky radius to load from reference catalog.
733 Name of thing being fit: "astrometry" or "photometry".
734 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
735 Reference object loader to use to load a reference catalog.
736 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
737 Selector to use to pick objects from the loaded reference catalog.
738 fit_function : callable
739 Function to call to perform fit (takes Associations object).
740 tract : `str`, optional
741 Name of tract currently being fit.
742 profile_jointcal : `bool`, optional
743 Separately profile the fitting step.
744 match_cut : `float`, optional
745 Radius in arcseconds to find cross-catalog matches to during
746 associations.associateCatalogs.
747 reject_bad_fluxes : `bool`, optional
748 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
752 result : `Photometry` or `Astrometry`
753 Result of `fit_function()`
755 self.log.info(
"====== Now processing %s...", name)
758 associations.associateCatalogs(match_cut)
760 associations.fittedStarListSize())
762 applyColorterms =
False if name.lower() ==
"astrometry" else self.config.applyColorTerms
765 center, radius, defaultBand,
766 applyColorterms=applyColorterms,
770 associations.collectRefStars(refCat,
771 self.config.matchCut*lsst.geom.arcseconds,
773 refCoordinateErr=refCoordErr,
774 rejectBadFluxes=reject_bad_fluxes)
776 associations.refStarListSize())
778 associations.prepareFittedStars(self.config.minMeasurements)
782 associations.nFittedStarsWithAssociatedRefStar())
784 associations.fittedStarListSize())
786 associations.nCcdImagesValidForFit())
788 load_cat_prof_file =
'jointcal_fit_%s.prof'%name
if profile_jointcal
else ''
789 dataName =
"{}_{}".format(tract, defaultBand)
790 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
791 result = fit_function(associations, dataName)
794 if self.config.writeChi2FilesInitialFinal:
795 baseName = self.
_getDebugPath_getDebugPath(f
"{name}_final_chi2-{dataName}")
796 result.fit.saveChi2Contributions(baseName+
"{type}")
797 self.log.info(
"Wrote chi2 contributions files: %s", baseName)
801 def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName,
802 applyColorterms=False, epoch=None):
803 """Load the necessary reference catalog sources, convert fluxes to
804 correct units, and apply color term corrections if requested.
808 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
809 The reference catalog loader to use to get the data.
810 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
811 Source selector to apply to loaded reference catalog.
812 center : `lsst.geom.SpherePoint`
813 The center around which to load sources.
814 radius : `lsst.geom.Angle`
815 The radius around ``center`` to load sources in.
817 The name of the camera filter to load fluxes for.
818 applyColorterms : `bool`
819 Apply colorterm corrections to the refcat for ``filterName``?
820 epoch : `astropy.time.Time`, optional
821 Epoch to which to correct refcat proper motion and parallax,
822 or `None` to not apply such corrections.
826 refCat : `lsst.afw.table.SimpleCatalog`
827 The loaded reference catalog.
829 The name of the reference catalog flux field appropriate for ``filterName``.
831 skyCircle = refObjLoader.loadSkyCircle(center,
836 selected = referenceSelector.run(skyCircle.refCat)
838 if not selected.sourceCat.isContiguous():
839 refCat = selected.sourceCat.copy(deep=
True)
841 refCat = selected.sourceCat
844 refCatName = refObjLoader.ref_dataset_name
845 self.log.info(
"Applying color terms for filterName=%r reference catalog=%s",
846 filterName, refCatName)
847 colorterm = self.config.colorterms.getColorterm(
848 filterName=filterName, photoCatName=refCatName, doRaise=
True)
850 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat, filterName)
851 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
853 refCat[skyCircle.fluxField+
'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
855 return refCat, skyCircle.fluxField
857 def _check_star_lists(self, associations, name):
859 if associations.nCcdImagesValidForFit() == 0:
860 raise RuntimeError(
'No images in the ccdImageList!')
861 if associations.fittedStarListSize() == 0:
862 raise RuntimeError(
'No stars in the {} fittedStarList!'.format(name))
863 if associations.refStarListSize() == 0:
864 raise RuntimeError(
'No stars in the {} reference star list!'.format(name))
866 def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
867 """Compute chi2, log it, validate the model, and return chi2.
871 associations : `lsst.jointcal.Associations`
872 The star/reference star associations to fit.
873 fit : `lsst.jointcal.FitterBase`
874 The fitter to use for minimization.
875 model : `lsst.jointcal.Model`
878 Label to describe the chi2 (e.g. "Initialized", "Final").
879 writeChi2Name : `str`, optional
880 Filename prefix to write the chi2 contributions to.
881 Do not supply an extension: an appropriate one will be added.
885 chi2: `lsst.jointcal.Chi2Accumulator`
886 The chi2 object for the current fitter and model.
891 Raised if chi2 is infinite or NaN.
893 Raised if the model is not valid.
895 if writeChi2Name
is not None:
897 fit.saveChi2Contributions(fullpath+
"{type}")
898 self.log.info(
"Wrote chi2 contributions files: %s", fullpath)
900 chi2 = fit.computeChi2()
901 self.log.info(
"%s %s", chi2Label, chi2)
903 if not np.isfinite(chi2.chi2):
904 raise FloatingPointError(f
'{chi2Label} chi2 is invalid: {chi2}')
905 if not model.validate(associations.getCcdImageList(), chi2.ndof):
906 raise ValueError(
"Model is not valid: check log messages for warnings.")
909 def _fit_photometry(self, associations, dataName=None):
911 Fit the photometric data.
915 associations : `lsst.jointcal.Associations`
916 The star/reference star associations to fit.
918 Name of the data being processed (e.g. "1234_HSC-Y"), for
919 identifying debugging files.
923 fit_result : `namedtuple`
924 fit : `lsst.jointcal.PhotometryFit`
925 The photometric fitter used to perform the fit.
926 model : `lsst.jointcal.PhotometryModel`
927 The photometric model that was fit.
929 self.log.info(
"=== Starting photometric fitting...")
932 if self.config.photometryModel ==
"constrainedFlux":
935 visitOrder=self.config.photometryVisitOrder,
936 errorPedestal=self.config.photometryErrorPedestal)
938 doLineSearch = self.config.allowLineSearch
939 elif self.config.photometryModel ==
"constrainedMagnitude":
942 visitOrder=self.config.photometryVisitOrder,
943 errorPedestal=self.config.photometryErrorPedestal)
945 doLineSearch = self.config.allowLineSearch
946 elif self.config.photometryModel ==
"simpleFlux":
948 errorPedestal=self.config.photometryErrorPedestal)
950 elif self.config.photometryModel ==
"simpleMagnitude":
952 errorPedestal=self.config.photometryErrorPedestal)
958 if self.config.writeChi2FilesInitialFinal:
959 baseName = f
"photometry_initial_chi2-{dataName}"
962 if self.config.writeInitialModel:
963 fullpath = self.
_getDebugPath_getDebugPath(
"initialPhotometryModel.txt")
965 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialized", writeChi2Name=baseName)
967 def getChi2Name(whatToFit):
968 if self.config.writeChi2FilesOuterLoop:
969 return f
"photometry_init-%s_chi2-{dataName}" % whatToFit
975 dumpMatrixFile = self.
_getDebugPath_getDebugPath(
"photometry_preinit")
if self.config.writeInitMatrix
else ""
976 if self.config.photometryModel.startswith(
"constrained"):
979 fit.minimize(
"ModelVisit", dumpMatrixFile=dumpMatrixFile)
980 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize ModelVisit",
981 writeChi2Name=getChi2Name(
"ModelVisit"))
984 fit.minimize(
"Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
986 writeChi2Name=getChi2Name(
"Model"))
988 fit.minimize(
"Fluxes")
990 writeChi2Name=getChi2Name(
"Fluxes"))
992 fit.minimize(
"Model Fluxes", doLineSearch=doLineSearch)
993 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize ModelFluxes",
994 writeChi2Name=getChi2Name(
"ModelFluxes"))
996 model.freezeErrorTransform()
997 self.log.debug(
"Photometry error scales are frozen.")
1001 self.config.maxPhotometrySteps,
1004 doRankUpdate=self.config.photometryDoRankUpdate,
1005 doLineSearch=doLineSearch,
1012 def _fit_astrometry(self, associations, dataName=None):
1014 Fit the astrometric data.
1018 associations : `lsst.jointcal.Associations`
1019 The star/reference star associations to fit.
1021 Name of the data being processed (e.g. "1234_HSC-Y"), for
1022 identifying debugging files.
1026 fit_result : `namedtuple`
1027 fit : `lsst.jointcal.AstrometryFit`
1028 The astrometric fitter used to perform the fit.
1029 model : `lsst.jointcal.AstrometryModel`
1030 The astrometric model that was fit.
1031 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1032 The model for the sky to tangent plane projection that was used in the fit.
1035 self.log.info(
"=== Starting astrometric fitting...")
1037 associations.deprojectFittedStars()
1044 if self.config.astrometryModel ==
"constrained":
1046 sky_to_tan_projection,
1047 chipOrder=self.config.astrometryChipOrder,
1048 visitOrder=self.config.astrometryVisitOrder)
1049 elif self.config.astrometryModel ==
"simple":
1051 sky_to_tan_projection,
1052 self.config.useInputWcs,
1054 order=self.config.astrometrySimpleOrder)
1059 if self.config.writeChi2FilesInitialFinal:
1060 baseName = f
"astrometry_initial_chi2-{dataName}"
1063 if self.config.writeInitialModel:
1064 fullpath = self.
_getDebugPath_getDebugPath(
"initialAstrometryModel.txt")
1066 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initial", writeChi2Name=baseName)
1068 def getChi2Name(whatToFit):
1069 if self.config.writeChi2FilesOuterLoop:
1070 return f
"astrometry_init-%s_chi2-{dataName}" % whatToFit
1074 dumpMatrixFile = self.
_getDebugPath_getDebugPath(
"astrometry_preinit")
if self.config.writeInitMatrix
else ""
1077 if self.config.astrometryModel ==
"constrained":
1078 fit.minimize(
"DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1079 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize DistortionsVisit",
1080 writeChi2Name=getChi2Name(
"DistortionsVisit"))
1083 fit.minimize(
"Distortions", dumpMatrixFile=dumpMatrixFile)
1084 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize Distortions",
1085 writeChi2Name=getChi2Name(
"Distortions"))
1087 fit.minimize(
"Positions")
1088 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize Positions",
1089 writeChi2Name=getChi2Name(
"Positions"))
1091 fit.minimize(
"Distortions Positions")
1092 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize DistortionsPositions",
1093 writeChi2Name=getChi2Name(
"DistortionsPositions"))
1097 self.config.maxAstrometrySteps,
1099 "Distortions Positions",
1100 doRankUpdate=self.config.astrometryDoRankUpdate,
1106 return Astrometry(fit, model, sky_to_tan_projection)
1108 def _check_stars(self, associations):
1109 """Count measured and reference stars per ccd and warn/log them."""
1110 for ccdImage
in associations.getCcdImageList():
1111 nMeasuredStars, nRefStars = ccdImage.countStars()
1112 self.log.debug(
"ccdImage %s has %s measured and %s reference stars",
1113 ccdImage.getName(), nMeasuredStars, nRefStars)
1114 if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1115 self.log.warn(
"ccdImage %s has only %s measuredStars (desired %s)",
1116 ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1117 if nRefStars < self.config.minRefStarsPerCcd:
1118 self.log.warn(
"ccdImage %s has only %s RefStars (desired %s)",
1119 ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1121 def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1124 doLineSearch=False):
1125 """Run fitter.minimize up to max_steps times, returning the final chi2.
1129 associations : `lsst.jointcal.Associations`
1130 The star/reference star associations to fit.
1131 fitter : `lsst.jointcal.FitterBase`
1132 The fitter to use for minimization.
1134 Maximum number of steps to run outlier rejection before declaring
1135 convergence failure.
1136 name : {'photometry' or 'astrometry'}
1137 What type of data are we fitting (for logs and debugging files).
1139 Passed to ``fitter.minimize()`` to define the parameters to fit.
1140 dataName : `str`, optional
1141 Descriptive name for this dataset (e.g. tract and filter),
1143 doRankUpdate : `bool`, optional
1144 Do an Eigen rank update during minimization, or recompute the full
1145 matrix and gradient?
1146 doLineSearch : `bool`, optional
1147 Do a line search for the optimum step during minimization?
1151 chi2: `lsst.jointcal.Chi2Statistic`
1152 The final chi2 after the fit converges, or is forced to end.
1157 Raised if the fitter fails with a non-finite value.
1159 Raised if the fitter fails for some other reason;
1160 log messages will provide further details.
1162 dumpMatrixFile = self.
_getDebugPath_getDebugPath(f
"{name}_postinit")
if self.config.writeInitMatrix
else ""
1164 oldChi2.chi2 = float(
"inf")
1165 for i
in range(max_steps):
1166 if self.config.writeChi2FilesOuterLoop:
1167 writeChi2Name = f
"{name}_iterate_{i}_chi2-{dataName}"
1169 writeChi2Name =
None
1170 result = fitter.minimize(whatToFit,
1171 self.config.outlierRejectSigma,
1172 doRankUpdate=doRankUpdate,
1173 doLineSearch=doLineSearch,
1174 dumpMatrixFile=dumpMatrixFile)
1176 chi2 = self.
_logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(),
1177 f
"Fit iteration {i}", writeChi2Name=writeChi2Name)
1179 if result == MinimizeResult.Converged:
1181 self.log.debug(
"fit has converged - no more outliers - redo minimization "
1182 "one more time in case we have lost accuracy in rank update.")
1184 result = fitter.minimize(whatToFit, self.config.outlierRejectSigma)
1185 chi2 = self.
_logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(),
"Fit completed")
1188 if chi2.chi2/chi2.ndof >= 4.0:
1189 self.log.error(
"Potentially bad fit: High chi-squared/ndof.")
1192 elif result == MinimizeResult.Chi2Increased:
1193 self.log.warn(
"Still some outliers remaining but chi2 increased - retry")
1195 chi2Ratio = chi2.chi2 / oldChi2.chi2
1197 self.log.warn(
'Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1198 chi2.chi2, oldChi2.chi2, chi2Ratio)
1205 msg = (
"Large chi2 increase between steps: fit likely cannot converge."
1206 " Try setting one or more of the `writeChi2*` config fields and looking"
1207 " at how individual star chi2-values evolve during the fit.")
1208 raise RuntimeError(msg)
1210 elif result == MinimizeResult.NonFinite:
1211 filename = self.
_getDebugPath_getDebugPath(
"{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1213 fitter.saveChi2Contributions(filename+
"{type}")
1214 msg =
"Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1215 raise FloatingPointError(msg.format(filename))
1216 elif result == MinimizeResult.Failed:
1217 raise RuntimeError(
"Chi2 minimization failure, cannot complete fit.")
1219 raise RuntimeError(
"Unxepected return code from minimize().")
1221 self.log.error(
"%s failed to converge after %d steps"%(name, max_steps))
1225 def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1227 Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1231 associations : `lsst.jointcal.Associations`
1232 The star/reference star associations to fit.
1233 model : `lsst.jointcal.AstrometryModel`
1234 The astrometric model that was fit.
1235 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1236 Dict of ccdImage identifiers to dataRefs that were fit.
1239 ccdImageList = associations.getCcdImageList()
1240 for ccdImage
in ccdImageList:
1242 ccd = ccdImage.ccdId
1243 visit = ccdImage.visit
1244 dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1245 self.log.info(
"Updating WCS for visit: %d, ccd: %d", visit, ccd)
1246 skyWcs = model.makeSkyWcs(ccdImage)
1248 dataRef.put(skyWcs,
'jointcal_wcs')
1249 except pexExceptions.Exception
as e:
1250 self.log.fatal(
'Failed to write updated Wcs: %s', str(e))
1253 def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1255 Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1259 associations : `lsst.jointcal.Associations`
1260 The star/reference star associations to fit.
1261 model : `lsst.jointcal.PhotometryModel`
1262 The photoometric model that was fit.
1263 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1264 Dict of ccdImage identifiers to dataRefs that were fit.
1267 ccdImageList = associations.getCcdImageList()
1268 for ccdImage
in ccdImageList:
1270 ccd = ccdImage.ccdId
1271 visit = ccdImage.visit
1272 dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1273 self.log.info(
"Updating PhotoCalib for visit: %d, ccd: %d", visit, ccd)
1274 photoCalib = model.toPhotoCalib(ccdImage)
1276 dataRef.put(photoCalib,
'jointcal_photoCalib')
1277 except pexExceptions.Exception
as e:
1278 self.log.fatal(
'Failed to write updated PhotoCalib: %s', str(e))
The class that implements the relations between MeasuredStar and FittedStar.
Class that handles the astrometric least squares problem.
Simple structure to accumulate chi2 and ndof.
A multi-component model, fitting mappings for sensors and visits simultaneously.
A projection handler in which all CCDs from the same visit have the same tangent point.
Class that handles the photometric least squares problem.
A model where there is one independent transform per CcdImage.
def getTargetList(parsedCmd, **kwargs)
def _compute_proper_motion_epoch(self, ccdImageList)
def _get_refcat_coordinate_error_override(self, refCat, name)
def _getDebugPath(self, filename)
def _check_star_lists(self, associations, name)
def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef)
def _check_stars(self, associations)
def loadData(self, dataRefs, associations, jointcalControl, profile_jointcal=False)
def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName, applyColorterms=False, epoch=None)
def _readDataId(self, butler, dataId)
def __init__(self, butler=None, profile_jointcal=False, **kwargs)
def _build_ccdImage(self, data, associations, jointcalControl)
def _do_load_refcat_and_fit(self, associations, defaultBand, center, radius, tract="", profile_jointcal=False, match_cut=3.0, reject_bad_fluxes=False, *name="", refObjLoader=None, referenceSelector=None, fit_function=None)
def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", doRankUpdate=True, doLineSearch=False)
def _fit_photometry(self, associations, dataName=None)
def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef)
def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None)
def runDataRef(self, dataRefs, profile_jointcal=False)
def _fit_astrometry(self, associations, dataName=None)
def _prep_sky(self, associations, bands)
def add_measurement(job, name, value)
def writeModel(model, filename, log)