24 import astropy.units
as u
34 from lsst.pipe.tasks.colorterms
import ColortermLibrary
35 from lsst.verify
import Job, Measurement
40 from .dataIds
import PerTractCcdDataIdContainer
45 __all__ = [
"JointcalConfig",
"JointcalRunner",
"JointcalTask"]
47 Photometry = collections.namedtuple(
'Photometry', (
'fit',
'model'))
48 Astrometry = collections.namedtuple(
'Astrometry', (
'fit',
'model',
'sky_to_tan_projection'))
53 meas = Measurement(job.metrics[name], value)
54 job.measurements.insert(meas)
58 """Subclass of TaskRunner for jointcalTask 60 jointcalTask.runDataRef() takes a number of arguments, one of which is a list of dataRefs 61 extracted from the command line (whereas most CmdLineTasks' runDataRef methods take 62 single dataRef, are are called repeatedly). This class transforms the processed 63 arguments generated by the ArgumentParser into the arguments expected by 64 Jointcal.runDataRef(). 66 See pipeBase.TaskRunner for more information. 72 Return a list of tuples per tract, each containing (dataRefs, kwargs). 74 Jointcal operates on lists of dataRefs simultaneously. 76 kwargs[
'profile_jointcal'] = parsedCmd.profile_jointcal
77 kwargs[
'butler'] = parsedCmd.butler
81 for ref
in parsedCmd.id.refList:
82 refListDict.setdefault(ref.dataId[
"tract"], []).append(ref)
84 result = [(refListDict[tract], kwargs)
for tract
in sorted(refListDict.keys())]
92 Arguments for Task.runDataRef() 97 if self.doReturnResults is False: 99 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise. 101 if self.doReturnResults is True: 103 - ``result``: the result of calling jointcal.runDataRef() 104 - ``exitStatus``: 0 if the task completed successfully, 1 otherwise. 109 dataRefList, kwargs = args
110 butler = kwargs.pop(
'butler')
111 task = self.TaskClass(config=self.config, log=self.log, butler=butler)
114 result = task.runDataRef(dataRefList, **kwargs)
115 exitStatus = result.exitStatus
116 job_path = butler.get(
'verify_job_filename')
117 result.job.write(job_path[0])
118 except Exception
as e:
123 eName = type(e).__name__
124 tract = dataRefList[0].dataId[
'tract']
125 task.log.fatal(
"Failed processing tract %s, %s: %s", tract, eName, e)
128 kwargs[
'butler'] = butler
129 if self.doReturnResults:
130 return pipeBase.Struct(result=result, exitStatus=exitStatus)
132 return pipeBase.Struct(exitStatus=exitStatus)
136 """Configuration for JointcalTask""" 138 doAstrometry = pexConfig.Field(
139 doc=
"Fit astrometry and write the fitted result.",
143 doPhotometry = pexConfig.Field(
144 doc=
"Fit photometry and write the fitted result.",
148 coaddName = pexConfig.Field(
149 doc=
"Type of coadd, typically deep or goodSeeing",
153 positionErrorPedestal = pexConfig.Field(
154 doc=
"Systematic term to apply to the measured position error (pixels)",
158 photometryErrorPedestal = pexConfig.Field(
159 doc=
"Systematic term to apply to the measured error on flux or magnitude as a " 160 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
165 matchCut = pexConfig.Field(
166 doc=
"Matching radius between fitted and reference stars (arcseconds)",
170 minMeasurements = pexConfig.Field(
171 doc=
"Minimum number of associated measured stars for a fitted star to be included in the fit",
175 minMeasuredStarsPerCcd = pexConfig.Field(
176 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
180 minRefStarsPerCcd = pexConfig.Field(
181 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
185 allowLineSearch = pexConfig.Field(
186 doc=
"Allow a line search during minimization, if it is reasonable for the model" 187 " (models with a significant non-linear component, e.g. constrainedPhotometry).",
191 astrometrySimpleOrder = pexConfig.Field(
192 doc=
"Polynomial order for fitting the simple astrometry model.",
196 astrometryChipOrder = pexConfig.Field(
197 doc=
"Order of the per-chip transform for the constrained astrometry model.",
201 astrometryVisitOrder = pexConfig.Field(
202 doc=
"Order of the per-visit transform for the constrained astrometry model.",
206 useInputWcs = pexConfig.Field(
207 doc=
"Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
211 astrometryModel = pexConfig.ChoiceField(
212 doc=
"Type of model to fit to astrometry",
214 default=
"constrained",
215 allowed={
"simple":
"One polynomial per ccd",
216 "constrained":
"One polynomial per ccd, and one polynomial per visit"}
218 photometryModel = pexConfig.ChoiceField(
219 doc=
"Type of model to fit to photometry",
221 default=
"constrainedMagnitude",
222 allowed={
"simpleFlux":
"One constant zeropoint per ccd and visit, fitting in flux space.",
223 "constrainedFlux":
"Constrained zeropoint per ccd, and one polynomial per visit," 224 " fitting in flux space.",
225 "simpleMagnitude":
"One constant zeropoint per ccd and visit," 226 " fitting in magnitude space.",
227 "constrainedMagnitude":
"Constrained zeropoint per ccd, and one polynomial per visit," 228 " fitting in magnitude space.",
231 applyColorTerms = pexConfig.Field(
232 doc=
"Apply photometric color terms to reference stars?" 233 "Requires that colorterms be set to a ColortermLibrary",
237 colorterms = pexConfig.ConfigField(
238 doc=
"Library of photometric reference catalog name to color term dict.",
239 dtype=ColortermLibrary,
241 photometryVisitOrder = pexConfig.Field(
242 doc=
"Order of the per-visit polynomial transform for the constrained photometry model.",
246 photometryDoRankUpdate = pexConfig.Field(
247 doc=
"Do the rank update step during minimization. " 248 "Skipping this can help deal with models that are too non-linear.",
252 astrometryDoRankUpdate = pexConfig.Field(
253 doc=
"Do the rank update step during minimization (should not change the astrometry fit). " 254 "Skipping this can help deal with models that are too non-linear.",
258 outlierRejectSigma = pexConfig.Field(
259 doc=
"How many sigma to reject outliers at during minimization.",
263 maxPhotometrySteps = pexConfig.Field(
264 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
268 maxAstrometrySteps = pexConfig.Field(
269 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
273 astrometryRefObjLoader = pexConfig.ConfigurableField(
274 target=LoadIndexedReferenceObjectsTask,
275 doc=
"Reference object loader for astrometric fit",
277 photometryRefObjLoader = pexConfig.ConfigurableField(
278 target=LoadIndexedReferenceObjectsTask,
279 doc=
"Reference object loader for photometric fit",
281 sourceSelector = sourceSelectorRegistry.makeField(
282 doc=
"How to select sources for cross-matching",
285 astrometryReferenceSelector = pexConfig.ConfigurableField(
286 target=ReferenceSourceSelectorTask,
287 doc=
"How to down-select the loaded astrometry reference catalog.",
289 photometryReferenceSelector = pexConfig.ConfigurableField(
290 target=ReferenceSourceSelectorTask,
291 doc=
"How to down-select the loaded photometry reference catalog.",
293 astrometryReferenceErr = pexConfig.Field(
294 doc=
"Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*_err` fields." 295 " If None, then raise an exception if the reference catalog is missing coordinate errors." 296 " If specified, overrides any existing `coord_*_err` values.",
301 writeInitMatrix = pexConfig.Field(
303 doc=
"Write the pre/post-initialization Hessian and gradient to text files, for debugging." 304 "The output files will be of the form 'astrometry_preinit-mat.txt', in the current directory." 305 "Note that these files are the dense versions of the matrix, and so may be very large.",
308 writeChi2FilesInitialFinal = pexConfig.Field(
310 doc=
"Write .csv files containing the contributions to chi2 for the initialization and final fit.",
313 writeChi2FilesOuterLoop = pexConfig.Field(
315 doc=
"Write .csv files containing the contributions to chi2 for the outer fit loop.",
318 sourceFluxType = pexConfig.Field(
320 doc=
"Source flux field to use in source selection and to get fluxes from the catalog.",
327 msg =
"applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary." 328 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
341 fluxField = f
"slot_{self.sourceFluxType}Flux_instFlux" 342 self.
sourceSelector[
'science'].signalToNoise.fluxField = fluxField
343 self.
sourceSelector[
'science'].signalToNoise.errField = fluxField +
"Err" 349 badFlags = [
'base_PixelFlags_flag_edge',
'base_PixelFlags_flag_saturated',
350 'base_PixelFlags_flag_interpolatedCenter',
'base_SdssCentroid_flag',
351 'base_PsfFlux_flag',
'base_PixelFlags_flag_suspectCenter']
356 """Jointly astrometrically and photometrically calibrate a group of images.""" 358 ConfigClass = JointcalConfig
359 RunnerClass = JointcalRunner
360 _DefaultName =
"jointcal" 362 def __init__(self, butler=None, profile_jointcal=False, **kwargs):
364 Instantiate a JointcalTask. 368 butler : `lsst.daf.persistence.Butler` 369 The butler is passed to the refObjLoader constructor in case it is 370 needed. Ignored if the refObjLoader argument provides a loader directly. 371 Used to initialize the astrometry and photometry refObjLoaders. 372 profile_jointcal : `bool` 373 Set to True to profile different stages of this jointcal run. 375 pipeBase.CmdLineTask.__init__(self, **kwargs)
377 self.makeSubtask(
"sourceSelector")
378 if self.config.doAstrometry:
379 self.makeSubtask(
'astrometryRefObjLoader', butler=butler)
380 self.makeSubtask(
"astrometryReferenceSelector")
383 if self.config.doPhotometry:
384 self.makeSubtask(
'photometryRefObjLoader', butler=butler)
385 self.makeSubtask(
"photometryReferenceSelector")
390 self.
job = Job.load_metrics_package(subset=
'jointcal')
395 def _getMetadataName(self):
399 def _makeArgumentParser(cls):
400 """Create an argument parser""" 402 parser.add_argument(
"--profile_jointcal", default=
False, action=
"store_true",
403 help=
"Profile steps of jointcal separately.")
404 parser.add_id_argument(
"--id",
"calexp", help=
"data ID, e.g. --id visit=6789 ccd=0..9",
405 ContainerClass=PerTractCcdDataIdContainer)
408 def _build_ccdImage(self, dataRef, associations, jointcalControl):
410 Extract the necessary things from this dataRef to add a new ccdImage. 414 dataRef : `lsst.daf.persistence.ButlerDataRef` 415 DataRef to extract info from. 416 associations : `lsst.jointcal.Associations` 417 Object to add the info to, to construct a new CcdImage 418 jointcalControl : `jointcal.JointcalControl` 419 Control object for associations management 425 The TAN WCS of this image, read from the calexp 426 (`lsst.afw.geom.SkyWcs`). 428 A key to identify this dataRef by its visit and ccd ids 431 This calexp's filter (`str`). 433 if "visit" in dataRef.dataId.keys():
434 visit = dataRef.dataId[
"visit"]
436 visit = dataRef.getButler().queryMetadata(
"calexp", (
"visit"), dataRef.dataId)[0]
438 src = dataRef.get(
"src", flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS, immediate=
True)
440 visitInfo = dataRef.get(
'calexp_visitInfo')
441 detector = dataRef.get(
'calexp_detector')
442 ccdId = detector.getId()
443 photoCalib = dataRef.get(
'calexp_photoCalib')
444 tanWcs = dataRef.get(
'calexp_wcs')
445 bbox = dataRef.get(
'calexp_bbox')
446 filt = dataRef.get(
'calexp_filter')
447 filterName = filt.getName()
449 goodSrc = self.sourceSelector.run(src)
451 if len(goodSrc.sourceCat) == 0:
452 self.log.warn(
"No sources selected in visit %s ccd %s", visit, ccdId)
454 self.log.info(
"%d sources selected in visit %d ccd %d", len(goodSrc.sourceCat), visit, ccdId)
455 associations.createCcdImage(goodSrc.sourceCat,
466 Result = collections.namedtuple(
'Result_from_build_CcdImage', (
'wcs',
'key',
'filter'))
467 Key = collections.namedtuple(
'Key', (
'visit',
'ccd'))
468 return Result(tanWcs, Key(visit, ccdId), filterName)
473 Jointly calibrate the astrometry and photometry across a set of images. 477 dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef` 478 List of data references to the exposures to be fit. 479 profile_jointcal : `bool` 480 Profile the individual steps of jointcal. 484 result : `lsst.pipe.base.Struct` 485 Struct of metadata from the fit, containing: 488 The provided data references that were fit (with updated WCSs) 490 The original WCS from each dataRef 492 Dictionary of internally-computed metrics for testing/validation. 494 if len(dataRefs) == 0:
495 raise ValueError(
'Need a non-empty list of data references!')
499 sourceFluxField =
"slot_%sFlux" % (self.config.sourceFluxType,)
503 visit_ccd_to_dataRef = {}
506 load_cat_prof_file =
'jointcal_build_ccdImage.prof' if profile_jointcal
else '' 507 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
510 camera = dataRefs[0].get(
'camera', immediate=
True)
514 oldWcsList.append(result.wcs)
515 visit_ccd_to_dataRef[result.key] = ref
516 filters.append(result.filter)
517 filters = collections.Counter(filters)
519 associations.computeCommonTangentPoint()
521 boundingCircle = associations.computeBoundingCircle()
523 radius =
lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
526 defaultFilter = filters.most_common(1)[0][0]
527 self.log.debug(
"Using %s band for reference flux", defaultFilter)
530 tract = dataRefs[0].dataId[
'tract']
532 if self.config.doAstrometry:
536 referenceSelector=self.astrometryReferenceSelector,
538 profile_jointcal=profile_jointcal,
544 if self.config.doPhotometry:
548 referenceSelector=self.photometryReferenceSelector,
550 profile_jointcal=profile_jointcal,
553 reject_bad_fluxes=
True)
558 return pipeBase.Struct(dataRefs=dataRefs,
559 oldWcsList=oldWcsList,
563 defaultFilter=defaultFilter,
564 exitStatus=exitStatus)
566 def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
567 name="", refObjLoader=None, referenceSelector=None,
568 filters=[], fit_function=None,
569 tract=None, profile_jointcal=False, match_cut=3.0,
570 reject_bad_fluxes=False):
571 """Load reference catalog, perform the fit, and return the result. 575 associations : `lsst.jointcal.Associations` 576 The star/reference star associations to fit. 577 defaultFilter : `str` 578 filter to load from reference catalog. 579 center : `lsst.afw.geom.SpherePoint` 580 ICRS center of field to load from reference catalog. 581 radius : `lsst.afw.geom.Angle` 582 On-sky radius to load from reference catalog. 584 Name of thing being fit: "Astrometry" or "Photometry". 585 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask` 586 Reference object loader to load from for fit. 587 filters : `list` of `str`, optional 588 List of filters to load from the reference catalog. 589 fit_function : callable 590 Function to call to perform fit (takes associations object). 592 Name of tract currently being fit. 593 profile_jointcal : `bool`, optional 594 Separately profile the fitting step. 595 match_cut : `float`, optional 596 Radius in arcseconds to find cross-catalog matches to during 597 associations.associateCatalogs. 598 reject_bad_fluxes : `bool`, optional 599 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr. 603 result : `Photometry` or `Astrometry` 604 Result of `fit_function()` 606 self.log.info(
"====== Now processing %s...", name)
609 associations.associateCatalogs(match_cut)
611 associations.fittedStarListSize())
613 applyColorterms =
False if name ==
"Astrometry" else self.config.applyColorTerms
614 if name ==
"Astrometry":
615 referenceSelector = self.config.astrometryReferenceSelector
616 elif name ==
"Photometry":
617 referenceSelector = self.config.photometryReferenceSelector
619 center, radius, defaultFilter,
620 applyColorterms=applyColorterms)
622 if self.config.astrometryReferenceErr
is None:
623 refCoordErr = float(
'nan')
625 refCoordErr = self.config.astrometryReferenceErr
627 associations.collectRefStars(refCat,
628 self.config.matchCut*lsst.geom.arcseconds,
630 refCoordinateErr=refCoordErr,
631 rejectBadFluxes=reject_bad_fluxes)
633 associations.refStarListSize())
635 associations.prepareFittedStars(self.config.minMeasurements)
639 associations.nFittedStarsWithAssociatedRefStar())
641 associations.fittedStarListSize())
643 associations.nCcdImagesValidForFit())
645 load_cat_prof_file =
'jointcal_fit_%s.prof'%name
if profile_jointcal
else '' 646 dataName =
"{}_{}".format(tract, defaultFilter)
647 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
648 result = fit_function(associations, dataName)
651 if self.config.writeChi2FilesInitialFinal:
652 baseName = f
"{name}_final_chi2-{dataName}" 653 result.fit.saveChi2Contributions(baseName+
"{type}")
657 def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName,
658 applyColorterms=False):
659 """Load the necessary reference catalog sources, convert fluxes to 660 correct units, and apply color term corrections if requested. 664 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask` 665 The reference catalog loader to use to get the data. 666 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask` 667 Source selector to apply to loaded reference catalog. 668 center : `lsst.geom.SpherePoint` 669 The center around which to load sources. 670 radius : `lsst.geom.Angle` 671 The radius around ``center`` to load sources in. 673 The name of the camera filter to load fluxes for. 674 applyColorterms : `bool` 675 Apply colorterm corrections to the refcat for ``filterName``? 679 refCat : `lsst.afw.table.SimpleCatalog` 680 The loaded reference catalog. 682 The name of the reference catalog flux field appropriate for ``filterName``. 684 skyCircle = refObjLoader.loadSkyCircle(center,
688 selected = referenceSelector.run(skyCircle.refCat)
690 if not selected.sourceCat.isContiguous():
691 refCat = selected.sourceCat.copy(deep=
True)
693 refCat = selected.sourceCat
695 if self.config.astrometryReferenceErr
is None and 'coord_ra_err' not in refCat.schema:
696 msg = (
"Reference catalog does not contain coordinate errors, " 697 "and config.astrometryReferenceErr not supplied.")
698 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
702 if self.config.astrometryReferenceErr
is not None and 'coord_ra_err' in refCat.schema:
703 self.log.warn(
"Overriding reference catalog coordinate errors with %f/coordinate [mas]",
704 self.config.astrometryReferenceErr)
708 refCatName = refObjLoader.ref_dataset_name
709 except AttributeError:
711 raise RuntimeError(
"Cannot perform colorterm corrections with a.net refcats.")
712 self.log.info(
"Applying color terms for filterName=%r reference catalog=%s",
713 filterName, refCatName)
714 colorterm = self.config.colorterms.getColorterm(
715 filterName=filterName, photoCatName=refCatName, doRaise=
True)
717 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat, filterName)
718 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
720 refCat[skyCircle.fluxField+
'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
722 return refCat, skyCircle.fluxField
724 def _check_star_lists(self, associations, name):
726 if associations.nCcdImagesValidForFit() == 0:
727 raise RuntimeError(
'No images in the ccdImageList!')
728 if associations.fittedStarListSize() == 0:
729 raise RuntimeError(
'No stars in the {} fittedStarList!'.format(name))
730 if associations.refStarListSize() == 0:
731 raise RuntimeError(
'No stars in the {} reference star list!'.format(name))
733 def _logChi2AndValidate(self, associations, fit, model, chi2Label="Model",
735 """Compute chi2, log it, validate the model, and return chi2. 739 associations : `lsst.jointcal.Associations` 740 The star/reference star associations to fit. 741 fit : `lsst.jointcal.FitterBase` 742 The fitter to use for minimization. 743 model : `lsst.jointcal.Model` 745 chi2Label : str, optional 746 Label to describe the chi2 (e.g. "Initialized", "Final"). 747 writeChi2Name : `str`, optional 748 Filename prefix to write the chi2 contributions to. 749 Do not supply an extension: an appropriate one will be added. 753 chi2: `lsst.jointcal.Chi2Accumulator` 754 The chi2 object for the current fitter and model. 759 Raised if chi2 is infinite or NaN. 761 Raised if the model is not valid. 763 if writeChi2Name
is not None:
764 fit.saveChi2Contributions(writeChi2Name+
"{type}")
765 self.log.info(
"Wrote chi2 contributions files: %s", writeChi2Name)
767 chi2 = fit.computeChi2()
768 self.log.info(
"%s %s", chi2Label, chi2)
770 if not np.isfinite(chi2.chi2):
771 raise FloatingPointError(f
'{chi2Label} chi2 is invalid: {chi2}')
772 if not model.validate(associations.getCcdImageList(), chi2.ndof):
773 raise ValueError(
"Model is not valid: check log messages for warnings.")
776 def _fit_photometry(self, associations, dataName=None):
778 Fit the photometric data. 782 associations : `lsst.jointcal.Associations` 783 The star/reference star associations to fit. 785 Name of the data being processed (e.g. "1234_HSC-Y"), for 786 identifying debugging files. 790 fit_result : `namedtuple` 791 fit : `lsst.jointcal.PhotometryFit` 792 The photometric fitter used to perform the fit. 793 model : `lsst.jointcal.PhotometryModel` 794 The photometric model that was fit. 796 self.log.info(
"=== Starting photometric fitting...")
799 if self.config.photometryModel ==
"constrainedFlux":
802 visitOrder=self.config.photometryVisitOrder,
803 errorPedestal=self.config.photometryErrorPedestal)
805 doLineSearch = self.config.allowLineSearch
806 elif self.config.photometryModel ==
"constrainedMagnitude":
809 visitOrder=self.config.photometryVisitOrder,
810 errorPedestal=self.config.photometryErrorPedestal)
812 doLineSearch = self.config.allowLineSearch
813 elif self.config.photometryModel ==
"simpleFlux":
815 errorPedestal=self.config.photometryErrorPedestal)
817 elif self.config.photometryModel ==
"simpleMagnitude":
819 errorPedestal=self.config.photometryErrorPedestal)
825 if self.config.writeChi2FilesInitialFinal:
826 baseName = f
"photometry_initial_chi2-{dataName}" 831 def getChi2Name(whatToFit):
832 if self.config.writeChi2FilesOuterLoop:
833 return f
"photometry_init-%s_chi2-{dataName}" % whatToFit
839 dumpMatrixFile =
"photometry_preinit" if self.config.writeInitMatrix
else "" 840 if self.config.photometryModel.startswith(
"constrained"):
843 fit.minimize(
"ModelVisit", dumpMatrixFile=dumpMatrixFile)
844 self.
_logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name(
"ModelVisit"))
847 fit.minimize(
"Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
850 fit.minimize(
"Fluxes")
853 fit.minimize(
"Model Fluxes", doLineSearch=doLineSearch)
855 writeChi2Name=getChi2Name(
"ModelFluxes"))
857 model.freezeErrorTransform()
858 self.log.debug(
"Photometry error scales are frozen.")
862 self.config.maxPhotometrySteps,
865 doRankUpdate=self.config.photometryDoRankUpdate,
866 doLineSearch=doLineSearch,
873 def _fit_astrometry(self, associations, dataName=None):
875 Fit the astrometric data. 879 associations : `lsst.jointcal.Associations` 880 The star/reference star associations to fit. 882 Name of the data being processed (e.g. "1234_HSC-Y"), for 883 identifying debugging files. 887 fit_result : `namedtuple` 888 fit : `lsst.jointcal.AstrometryFit` 889 The astrometric fitter used to perform the fit. 890 model : `lsst.jointcal.AstrometryModel` 891 The astrometric model that was fit. 892 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler` 893 The model for the sky to tangent plane projection that was used in the fit. 896 self.log.info(
"=== Starting astrometric fitting...")
898 associations.deprojectFittedStars()
905 if self.config.astrometryModel ==
"constrained":
907 sky_to_tan_projection,
908 chipOrder=self.config.astrometryChipOrder,
909 visitOrder=self.config.astrometryVisitOrder)
910 elif self.config.astrometryModel ==
"simple":
912 sky_to_tan_projection,
913 self.config.useInputWcs,
915 order=self.config.astrometrySimpleOrder)
920 if self.config.writeChi2FilesInitialFinal:
921 baseName = f
"astrometry_initial_chi2-{dataName}" 926 def getChi2Name(whatToFit):
927 if self.config.writeChi2FilesOuterLoop:
928 return f
"astrometry_init-%s_chi2-{dataName}" % whatToFit
932 dumpMatrixFile =
"astrometry_preinit" if self.config.writeInitMatrix
else "" 935 if self.config.astrometryModel ==
"constrained":
936 fit.minimize(
"DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
937 self.
_logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name(
"DistortionsVisit"))
940 fit.minimize(
"Distortions", dumpMatrixFile=dumpMatrixFile)
941 self.
_logChi2AndValidate(associations, fit, model, writeChi2Name=getChi2Name(
"Distortions"))
943 fit.minimize(
"Positions")
946 fit.minimize(
"Distortions Positions")
948 writeChi2Name=getChi2Name(
"DistortionsPositions"))
952 self.config.maxAstrometrySteps,
954 "Distortions Positions",
955 doRankUpdate=self.config.astrometryDoRankUpdate,
961 return Astrometry(fit, model, sky_to_tan_projection)
963 def _check_stars(self, associations):
964 """Count measured and reference stars per ccd and warn/log them.""" 965 for ccdImage
in associations.getCcdImageList():
966 nMeasuredStars, nRefStars = ccdImage.countStars()
967 self.log.debug(
"ccdImage %s has %s measured and %s reference stars",
968 ccdImage.getName(), nMeasuredStars, nRefStars)
969 if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
970 self.log.warn(
"ccdImage %s has only %s measuredStars (desired %s)",
971 ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
972 if nRefStars < self.config.minRefStarsPerCcd:
973 self.log.warn(
"ccdImage %s has only %s RefStars (desired %s)",
974 ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
976 def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
980 """Run fitter.minimize up to max_steps times, returning the final chi2. 984 associations : `lsst.jointcal.Associations` 985 The star/reference star associations to fit. 986 fitter : `lsst.jointcal.FitterBase` 987 The fitter to use for minimization. 989 Maximum number of steps to run outlier rejection before declaring 991 name : {'photometry' or 'astrometry'} 992 What type of data are we fitting (for logs and debugging files). 994 Passed to ``fitter.minimize()`` to define the parameters to fit. 995 dataName : `str`, optional 996 Descriptive name for this dataset (e.g. tract and filter), 998 doRankUpdate : `bool`, optional 999 Do an Eigen rank update during minimization, or recompute the full 1000 matrix and gradient? 1001 doLineSearch : `bool`, optional 1002 Do a line search for the optimum step during minimization? 1006 chi2: `lsst.jointcal.Chi2Statistic` 1007 The final chi2 after the fit converges, or is forced to end. 1012 Raised if the fitter fails with a non-finite value. 1014 Raised if the fitter fails for some other reason; 1015 log messages will provide further details. 1017 dumpMatrixFile =
"%s_postinit" % name
if self.config.writeInitMatrix
else "" 1018 for i
in range(max_steps):
1019 if self.config.writeChi2FilesOuterLoop:
1020 writeChi2Name = f
"{name}_iterate_{i}_chi2-{dataName}" 1022 writeChi2Name =
None 1023 result = fitter.minimize(whatToFit,
1024 self.config.outlierRejectSigma,
1025 doRankUpdate=doRankUpdate,
1026 doLineSearch=doLineSearch,
1027 dumpMatrixFile=dumpMatrixFile)
1030 writeChi2Name=writeChi2Name)
1032 if result == MinimizeResult.Converged:
1034 self.log.debug(
"fit has converged - no more outliers - redo minimization " 1035 "one more time in case we have lost accuracy in rank update.")
1037 result = fitter.minimize(whatToFit, self.config.outlierRejectSigma)
1041 if chi2.chi2/chi2.ndof >= 4.0:
1042 self.log.error(
"Potentially bad fit: High chi-squared/ndof.")
1045 elif result == MinimizeResult.Chi2Increased:
1046 self.log.warn(
"still some outliers but chi2 increases - retry")
1047 elif result == MinimizeResult.NonFinite:
1048 filename =
"{}_failure-nonfinite_chi2-{}.csv".format(name, dataName)
1050 fitter.saveChi2Contributions(filename)
1051 msg =
"Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}" 1052 raise FloatingPointError(msg.format(filename))
1053 elif result == MinimizeResult.Failed:
1054 raise RuntimeError(
"Chi2 minimization failure, cannot complete fit.")
1056 raise RuntimeError(
"Unxepected return code from minimize().")
1058 self.log.error(
"%s failed to converge after %d steps"%(name, max_steps))
1062 def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1064 Write the fitted astrometric results to a new 'jointcal_wcs' dataRef. 1068 associations : `lsst.jointcal.Associations` 1069 The star/reference star associations to fit. 1070 model : `lsst.jointcal.AstrometryModel` 1071 The astrometric model that was fit. 1072 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef` 1073 Dict of ccdImage identifiers to dataRefs that were fit. 1076 ccdImageList = associations.getCcdImageList()
1077 for ccdImage
in ccdImageList:
1079 ccd = ccdImage.ccdId
1080 visit = ccdImage.visit
1081 dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1082 self.log.info(
"Updating WCS for visit: %d, ccd: %d", visit, ccd)
1083 skyWcs = model.makeSkyWcs(ccdImage)
1085 dataRef.put(skyWcs,
'jointcal_wcs')
1086 except pexExceptions.Exception
as e:
1087 self.log.fatal(
'Failed to write updated Wcs: %s', str(e))
1090 def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1092 Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef. 1096 associations : `lsst.jointcal.Associations` 1097 The star/reference star associations to fit. 1098 model : `lsst.jointcal.PhotometryModel` 1099 The photoometric model that was fit. 1100 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef` 1101 Dict of ccdImage identifiers to dataRefs that were fit. 1104 ccdImageList = associations.getCcdImageList()
1105 for ccdImage
in ccdImageList:
1107 ccd = ccdImage.ccdId
1108 visit = ccdImage.visit
1109 dataRef = visit_ccd_to_dataRef[(visit, ccd)]
1110 self.log.info(
"Updating PhotoCalib for visit: %d, ccd: %d", visit, ccd)
1111 photoCalib = model.toPhotoCalib(ccdImage)
1113 dataRef.put(photoCalib,
'jointcal_photoCalib')
1114 except pexExceptions.Exception
as e:
1115 self.log.fatal(
'Failed to write updated PhotoCalib: %s', str(e))
def runDataRef(self, dataRefs, profile_jointcal=False)
def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterName, applyColorterms=False)
def _build_ccdImage(self, dataRef, associations, jointcalControl)
def _fit_photometry(self, associations, dataName=None)
def getTargetList(parsedCmd, kwargs)
def _fit_astrometry(self, associations, dataName=None)
def _check_star_lists(self, associations, name)
The class that implements the relations between MeasuredStar and FittedStar.
A projection handler in which all CCDs from the same visit have the same tangent point.
def _logChi2AndValidate(self, associations, fit, model, chi2Label="Model", writeChi2Name=None)
A model where there is one independent transform per CcdImage.
def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", doRankUpdate=True, doLineSearch=False)
def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef)
def _check_stars(self, associations)
Class that handles the photometric least squares problem.
Class that handles the astrometric least squares problem.
def add_measurement(job, name, value)
def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius, name="", refObjLoader=None, referenceSelector=None, filters=[], fit_function=None, tract=None, profile_jointcal=False, match_cut=3.0, reject_bad_fluxes=False)
A multi-component model, fitting mappings for sensors and visits simultaneously.
def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef)
def __init__(self, butler=None, profile_jointcal=False, kwargs)