28 import astropy.units
as u
33 import lsst.pipe.base
as pipeBase
39 from lsst.obs.base
import Instrument
40 from lsst.pipe.tasks.colorterms
import ColortermLibrary
41 from lsst.verify
import Job, Measurement
44 ReferenceObjectLoader)
47 from .dataIds
import PerTractCcdDataIdContainer
52 __all__ = [
"JointcalConfig",
"JointcalRunner",
"JointcalTask"]
54 Photometry = collections.namedtuple(
'Photometry', (
'fit',
'model'))
55 Astrometry = collections.namedtuple(
'Astrometry', (
'fit',
'model',
'sky_to_tan_projection'))
60 meas = Measurement(job.metrics[name], value)
61 job.measurements.insert(meas)
65 """Subclass of TaskRunner for jointcalTask (gen2)
67 jointcalTask.runDataRef() takes a number of arguments, one of which is a list of dataRefs
68 extracted from the command line (whereas most CmdLineTasks' runDataRef methods take
69 single dataRef, are are called repeatedly). This class transforms the processed
70 arguments generated by the ArgumentParser into the arguments expected by
71 Jointcal.runDataRef().
73 See pipeBase.TaskRunner for more information.
79 Return a list of tuples per tract, each containing (dataRefs, kwargs).
81 Jointcal operates on lists of dataRefs simultaneously.
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 """Lookup function that asserts/hopes that a static calibration dataset
143 exists in a particular collection, since this task can't provide a single
144 date/time to use to search for one properly.
146 This is mostly useful for the ``camera`` dataset, in cases where the task's
147 quantum dimensions do *not* include something temporal, like ``exposure``
152 datasetType : `lsst.daf.butler.DatasetType`
153 Type of dataset being searched for.
154 registry : `lsst.daf.butler.Registry`
155 Data repository registry to search.
156 quantumDataId : `lsst.daf.butler.DataCoordinate`
157 Data ID of the quantum this camera should match.
158 collections : `Iterable` [ `str` ]
159 Collections that should be searched - but this lookup function works
160 by ignoring this in favor of a more-or-less hard-coded value.
164 refs : `Iterator` [ `lsst.daf.butler.DatasetRef` ]
165 Iterator over dataset references; should have only one element.
169 This implementation duplicates one in fgcmcal, and is at least quite
170 similar to another in cp_pipe. This duplicate has the most documentation.
171 Fixing this is DM-29661.
173 instrument = Instrument.fromName(quantumDataId[
"instrument"], registry)
174 unboundedCollection = instrument.makeUnboundedCalibrationRunName()
175 return registry.queryDatasets(datasetType,
176 dataId=quantumDataId,
177 collections=[unboundedCollection],
182 """Lookup function that finds all refcats for all visits that overlap a
183 tract, rather than just the refcats that directly overlap the tract.
187 datasetType : `lsst.daf.butler.DatasetType`
188 Type of dataset being searched for.
189 registry : `lsst.daf.butler.Registry`
190 Data repository registry to search.
191 quantumDataId : `lsst.daf.butler.DataCoordinate`
192 Data ID of the quantum; expected to be something we can use as a
193 constraint to query for overlapping visits.
194 collections : `Iterable` [ `str` ]
195 Collections to search.
199 refs : `Iterator` [ `lsst.daf.butler.DatasetRef` ]
200 Iterator over refcat references.
209 for visit_data_id
in set(registry.queryDataIds(
"visit", dataId=quantumDataId).expanded()):
211 registry.queryDatasets(
213 collections=collections,
214 dataId=visit_data_id,
222 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter")):
223 """Middleware input/output connections for jointcal data."""
224 inputCamera = pipeBase.connectionTypes.PrerequisiteInput(
225 doc=
"The camera instrument that took these observations.",
227 storageClass=
"Camera",
228 dimensions=(
"instrument",),
230 lookupFunction=lookupStaticCalibrations,
232 inputSourceTableVisit = pipeBase.connectionTypes.Input(
233 doc=
"Source table in parquet format, per visit",
234 name=
"sourceTable_visit",
235 storageClass=
"DataFrame",
236 dimensions=(
"instrument",
"visit"),
240 inputVisitSummary = pipeBase.connectionTypes.Input(
241 doc=(
"Per-visit consolidated exposure metadata built from calexps. "
242 "These catalogs use detector id for the id and must be sorted for "
243 "fast lookups of a detector."),
245 storageClass=
"ExposureCatalog",
246 dimensions=(
"instrument",
"visit"),
250 astrometryRefCat = pipeBase.connectionTypes.PrerequisiteInput(
251 doc=
"The astrometry reference catalog to match to loaded input catalog sources.",
252 name=
"gaia_dr2_20200414",
253 storageClass=
"SimpleCatalog",
254 dimensions=(
"skypix",),
257 lookupFunction=lookupVisitRefCats,
259 photometryRefCat = pipeBase.connectionTypes.PrerequisiteInput(
260 doc=
"The photometry reference catalog to match to loaded input catalog sources.",
261 name=
"ps1_pv3_3pi_20170110",
262 storageClass=
"SimpleCatalog",
263 dimensions=(
"skypix",),
266 lookupFunction=lookupVisitRefCats,
269 outputWcs = pipeBase.connectionTypes.Output(
270 doc=(
"Per-tract, per-visit world coordinate systems derived from the fitted model."
271 " These catalogs only contain entries for detectors with an output, and use"
272 " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
273 name=
"jointcalSkyWcsCatalog",
274 storageClass=
"ExposureCatalog",
275 dimensions=(
"instrument",
"visit",
"skymap",
"tract"),
278 outputPhotoCalib = pipeBase.connectionTypes.Output(
279 doc=(
"Per-tract, per-visit photometric calibrations derived from the fitted model."
280 " These catalogs only contain entries for detectors with an output, and use"
281 " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
282 name=
"jointcalPhotoCalibCatalog",
283 storageClass=
"ExposureCatalog",
284 dimensions=(
"instrument",
"visit",
"skymap",
"tract"),
296 if not config.doAstrometry:
297 self.prerequisiteInputs.remove(
"astrometryRefCat")
298 self.outputs.remove(
"outputWcs")
299 if not config.doPhotometry:
300 self.prerequisiteInputs.remove(
"photometryRefCat")
301 self.outputs.remove(
"outputPhotoCalib")
305 pipelineConnections=JointcalTaskConnections):
306 """Configuration for JointcalTask"""
308 doAstrometry = pexConfig.Field(
309 doc=
"Fit astrometry and write the fitted result.",
313 doPhotometry = pexConfig.Field(
314 doc=
"Fit photometry and write the fitted result.",
318 coaddName = pexConfig.Field(
319 doc=
"Type of coadd, typically deep or goodSeeing",
324 sourceFluxType = pexConfig.Field(
326 doc=
"Source flux field to use in source selection and to get fluxes from the catalog.",
329 positionErrorPedestal = pexConfig.Field(
330 doc=
"Systematic term to apply to the measured position error (pixels)",
334 photometryErrorPedestal = pexConfig.Field(
335 doc=
"Systematic term to apply to the measured error on flux or magnitude as a "
336 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
341 matchCut = pexConfig.Field(
342 doc=
"Matching radius between fitted and reference stars (arcseconds)",
346 minMeasurements = pexConfig.Field(
347 doc=
"Minimum number of associated measured stars for a fitted star to be included in the fit",
351 minMeasuredStarsPerCcd = pexConfig.Field(
352 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
356 minRefStarsPerCcd = pexConfig.Field(
357 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
361 allowLineSearch = pexConfig.Field(
362 doc=
"Allow a line search during minimization, if it is reasonable for the model"
363 " (models with a significant non-linear component, e.g. constrainedPhotometry).",
367 astrometrySimpleOrder = pexConfig.Field(
368 doc=
"Polynomial order for fitting the simple astrometry model.",
372 astrometryChipOrder = pexConfig.Field(
373 doc=
"Order of the per-chip transform for the constrained astrometry model.",
377 astrometryVisitOrder = pexConfig.Field(
378 doc=
"Order of the per-visit transform for the constrained astrometry model.",
382 useInputWcs = pexConfig.Field(
383 doc=
"Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
387 astrometryModel = pexConfig.ChoiceField(
388 doc=
"Type of model to fit to astrometry",
390 default=
"constrained",
391 allowed={
"simple":
"One polynomial per ccd",
392 "constrained":
"One polynomial per ccd, and one polynomial per visit"}
394 photometryModel = pexConfig.ChoiceField(
395 doc=
"Type of model to fit to photometry",
397 default=
"constrainedMagnitude",
398 allowed={
"simpleFlux":
"One constant zeropoint per ccd and visit, fitting in flux space.",
399 "constrainedFlux":
"Constrained zeropoint per ccd, and one polynomial per visit,"
400 " fitting in flux space.",
401 "simpleMagnitude":
"One constant zeropoint per ccd and visit,"
402 " fitting in magnitude space.",
403 "constrainedMagnitude":
"Constrained zeropoint per ccd, and one polynomial per visit,"
404 " fitting in magnitude space.",
407 applyColorTerms = pexConfig.Field(
408 doc=
"Apply photometric color terms to reference stars?"
409 "Requires that colorterms be set to a ColortermLibrary",
413 colorterms = pexConfig.ConfigField(
414 doc=
"Library of photometric reference catalog name to color term dict.",
415 dtype=ColortermLibrary,
417 photometryVisitOrder = pexConfig.Field(
418 doc=
"Order of the per-visit polynomial transform for the constrained photometry model.",
422 photometryDoRankUpdate = pexConfig.Field(
423 doc=(
"Do the rank update step during minimization. "
424 "Skipping this can help deal with models that are too non-linear."),
428 astrometryDoRankUpdate = pexConfig.Field(
429 doc=(
"Do the rank update step during minimization (should not change the astrometry fit). "
430 "Skipping this can help deal with models that are too non-linear."),
434 outlierRejectSigma = pexConfig.Field(
435 doc=
"How many sigma to reject outliers at during minimization.",
439 astrometryOutlierRelativeTolerance = pexConfig.Field(
440 doc=(
"Convergence tolerance for outlier rejection threshold when fitting astrometry. Iterations will "
441 "stop when the fractional change in the chi2 cut level is below this value. If tolerance is set "
442 "to zero, iterations will continue until there are no more outliers. We suggest a value of 0.002"
443 "as a balance between a shorter minimization runtime and achieving a final fitted model that is"
444 "close to the solution found when removing all outliers."),
448 maxPhotometrySteps = pexConfig.Field(
449 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
453 maxAstrometrySteps = pexConfig.Field(
454 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
458 astrometryRefObjLoader = pexConfig.ConfigurableField(
459 target=LoadIndexedReferenceObjectsTask,
460 doc=
"Reference object loader for astrometric fit",
462 photometryRefObjLoader = pexConfig.ConfigurableField(
463 target=LoadIndexedReferenceObjectsTask,
464 doc=
"Reference object loader for photometric fit",
466 sourceSelector = sourceSelectorRegistry.makeField(
467 doc=
"How to select sources for cross-matching",
470 astrometryReferenceSelector = pexConfig.ConfigurableField(
471 target=ReferenceSourceSelectorTask,
472 doc=
"How to down-select the loaded astrometry reference catalog.",
474 photometryReferenceSelector = pexConfig.ConfigurableField(
475 target=ReferenceSourceSelectorTask,
476 doc=
"How to down-select the loaded photometry reference catalog.",
478 astrometryReferenceErr = pexConfig.Field(
479 doc=(
"Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
480 "If None, then raise an exception if the reference catalog is missing coordinate errors. "
481 "If specified, overrides any existing `coord_*Err` values."),
488 writeInitMatrix = pexConfig.Field(
490 doc=(
"Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
491 "Output files will be written to `config.debugOutputPath` and will "
492 "be of the form 'astrometry_[pre|post]init-TRACT-FILTER-mat.txt'. "
493 "Note that these files are the dense versions of the matrix, and so may be very large."),
496 writeChi2FilesInitialFinal = pexConfig.Field(
498 doc=(
"Write .csv files containing the contributions to chi2 for the initialization and final fit. "
499 "Output files will be written to `config.debugOutputPath` and will "
500 "be of the form `astrometry_[initial|final]_chi2-TRACT-FILTER1."),
503 writeChi2FilesOuterLoop = pexConfig.Field(
505 doc=(
"Write .csv files containing the contributions to chi2 for the outer fit loop. "
506 "Output files will be written to `config.debugOutputPath` and will "
507 "be of the form `astrometry_init-NN_chi2-TRACT-FILTER`."),
510 writeInitialModel = pexConfig.Field(
512 doc=(
"Write the pre-initialization model to text files, for debugging. "
513 "Output files will be written to `config.debugOutputPath` and will be "
514 "of the form `initial_astrometry_model-TRACT_FILTER.txt`."),
517 debugOutputPath = pexConfig.Field(
520 doc=(
"Path to write debug output files to. Used by "
521 "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
523 detailedProfile = pexConfig.Field(
526 doc=
"Output separate profiling information for different parts of jointcal, e.g. data read, fitting"
532 msg =
"applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
533 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
535 msg = (
"Only doing astrometry, but Colorterms are not applied for astrometry;"
536 "applyColorTerms=True will be ignored.")
545 self.
sourceSelectorsourceSelector[
'science'].doSignalToNoise =
True
548 self.
sourceSelectorsourceSelector[
'science'].signalToNoise.minimum = 10.
550 fluxField = f
"slot_{self.sourceFluxType}Flux_instFlux"
551 self.
sourceSelectorsourceSelector[
'science'].signalToNoise.fluxField = fluxField
552 self.
sourceSelectorsourceSelector[
'science'].signalToNoise.errField = fluxField +
"Err"
558 badFlags = [
'base_PixelFlags_flag_edge',
'base_PixelFlags_flag_saturated',
559 'base_PixelFlags_flag_interpolatedCenter',
'base_SdssCentroid_flag',
560 'base_PsfFlux_flag',
'base_PixelFlags_flag_suspectCenter']
572 """Write model to outfile."""
573 with open(filename,
"w")
as file:
574 file.write(repr(model))
575 log.info(
"Wrote %s to file: %s", model, filename)
578 @dataclasses.dataclass
580 """The input data jointcal needs for each detector/visit."""
582 """The visit identifier of this exposure."""
584 """The catalog derived from this exposure."""
586 """The VisitInfo of this exposure."""
588 """The detector of this exposure."""
590 """The photometric calibration of this exposure."""
592 """The WCS of this exposure."""
594 """The bounding box of this exposure."""
596 """The filter of this exposure."""
600 """Astrometricly and photometricly calibrate across multiple visits of the
605 butler : `lsst.daf.persistence.Butler`
606 The butler is passed to the refObjLoader constructor in case it is
607 needed. Ignored if the refObjLoader argument provides a loader directly.
608 Used to initialize the astrometry and photometry refObjLoaders.
609 initInputs : `dict`, optional
610 Dictionary used to initialize PipelineTasks (empty for jointcal).
613 ConfigClass = JointcalConfig
614 RunnerClass = JointcalRunner
615 _DefaultName =
"jointcal"
617 def __init__(self, butler=None, initInputs=None, **kwargs):
619 self.makeSubtask(
"sourceSelector")
620 if self.config.doAstrometry:
621 if initInputs
is None:
623 self.makeSubtask(
'astrometryRefObjLoader', butler=butler)
624 self.makeSubtask(
"astrometryReferenceSelector")
627 if self.config.doPhotometry:
628 if initInputs
is None:
630 self.makeSubtask(
'photometryRefObjLoader', butler=butler)
631 self.makeSubtask(
"photometryReferenceSelector")
636 self.
jobjob = Job.load_metrics_package(subset=
'jointcal')
641 inputs = butlerQC.get(inputRefs)
643 tract = butlerQC.quantum.dataId[
'tract']
644 if self.config.doAstrometry:
646 dataIds=[ref.datasetRef.dataId
647 for ref
in inputRefs.astrometryRefCat],
648 refCats=inputs.pop(
'astrometryRefCat'),
649 config=self.config.astrometryRefObjLoader,
651 if self.config.doPhotometry:
653 dataIds=[ref.datasetRef.dataId
654 for ref
in inputRefs.photometryRefCat],
655 refCats=inputs.pop(
'photometryRefCat'),
656 config=self.config.photometryRefObjLoader,
658 outputs = self.
runrun(**inputs, tract=tract)
659 if self.config.doAstrometry:
660 self.
_put_output_put_output(butlerQC, outputs.outputWcs, outputRefs.outputWcs,
661 inputs[
'inputCamera'],
"setWcs")
662 if self.config.doPhotometry:
663 self.
_put_output_put_output(butlerQC, outputs.outputPhotoCalib, outputRefs.outputPhotoCalib,
664 inputs[
'inputCamera'],
"setPhotoCalib")
666 def _put_output(self, butlerQC, outputs, outputRefs, camera, setter):
667 """Persist the output datasets to their appropriate datarefs.
671 butlerQC : `ButlerQuantumContext`
672 A butler which is specialized to operate in the context of a
673 `lsst.daf.butler.Quantum`; This is the input to `runQuantum`.
674 outputs : `dict` [`tuple`, `lsst.afw.geom.SkyWcs`] or
675 `dict` [`tuple, `lsst.afw.image.PhotoCalib`]
676 The fitted objects to persist.
677 outputRefs : `list` [`OutputQuantizedConnection`]
678 The DatasetRefs to persist the data to.
679 camera : `lsst.afw.cameraGeom.Camera`
680 The camera for this instrument, to get detector ids from.
682 The method to call on the ExposureCatalog to set each output.
685 schema.addField(
'visit', type=
'I', doc=
'Visit number')
687 def new_catalog(visit, size):
688 """Return an catalog ready to be filled with appropriate output."""
691 catalog[
'visit'] = visit
693 metadata.add(
"COMMENT",
"Catalog id is detector id, sorted.")
694 metadata.add(
"COMMENT",
"Only detectors with data have entries.")
698 detectors_per_visit = collections.defaultdict(int)
701 detectors_per_visit[key[0]] += 1
703 for ref
in outputRefs:
704 visit = ref.dataId[
'visit']
705 catalog = new_catalog(visit, detectors_per_visit[visit])
708 for detector
in camera:
709 detectorId = detector.getId()
710 key = (ref.dataId[
'visit'], detectorId)
711 if key
not in outputs:
713 self.log.debug(
"No %s output for detector %s in visit %s",
714 setter[3:], detectorId, visit)
717 catalog[i].setId(detectorId)
718 getattr(catalog[i], setter)(outputs[key])
722 butlerQC.put(catalog, ref)
723 self.log.info(
"Wrote %s detectors to %s", i, ref)
725 def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None):
731 sourceFluxField =
"flux"
735 oldWcsList, bands = self.
_load_data_load_data(inputSourceTableVisit,
741 boundingCircle, center, radius, defaultFilter = self.
_prep_sky_prep_sky(associations, bands)
744 if self.config.doAstrometry:
748 referenceSelector=self.astrometryReferenceSelector,
752 astrometry_output = self.
_make_output_make_output(associations.getCcdImageList(),
756 astrometry_output =
None
758 if self.config.doPhotometry:
762 referenceSelector=self.photometryReferenceSelector,
766 reject_bad_fluxes=
True)
767 photometry_output = self.
_make_output_make_output(associations.getCcdImageList(),
771 photometry_output =
None
773 return pipeBase.Struct(outputWcs=astrometry_output,
774 outputPhotoCalib=photometry_output,
779 def _make_schema_table(self):
780 """Return an afw SourceTable to use as a base for creating the
781 SourceCatalog to insert values from the dataFrame into.
785 table : `lsst.afw.table.SourceTable`
786 Table with schema and slots to use to make SourceCatalogs.
789 schema.addField(
"centroid_x",
"D")
790 schema.addField(
"centroid_y",
"D")
791 schema.addField(
"centroid_xErr",
"F")
792 schema.addField(
"centroid_yErr",
"F")
793 schema.addField(
"shape_xx",
"D")
794 schema.addField(
"shape_yy",
"D")
795 schema.addField(
"shape_xy",
"D")
796 schema.addField(
"flux_instFlux",
"D")
797 schema.addField(
"flux_instFluxErr",
"D")
799 table.defineCentroid(
"centroid")
800 table.defineShape(
"shape")
803 def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId):
804 """Return an afw SourceCatalog extracted from a visit-level dataframe,
805 limited to just one detector.
809 table : `lsst.afw.table.SourceTable`
810 Table factory to use to make the SourceCatalog that will be
811 populated with data from ``visitCatalog``.
812 visitCatalog : `pandas.DataFrame`
813 DataFrame to extract a detector catalog from.
815 Numeric id of the detector to extract from ``visitCatalog``.
819 catalog : `lsst.afw.table.SourceCatalog`
820 Detector-level catalog extracted from ``visitCatalog``.
823 mapping = {
'sourceId':
'id',
826 'xErr':
'centroid_xErr',
827 'yErr':
'centroid_yErr',
831 f
'{self.config.sourceFluxType}_instFlux':
'flux_instFlux',
832 f
'{self.config.sourceFluxType}_instFluxErr':
'flux_instFluxErr',
840 detector_column =
"detector" if "detector" in visitCatalog.columns
else "ccd"
842 matched = visitCatalog[detector_column] == detectorId
843 catalog.resize(sum(matched))
844 view = visitCatalog.loc[matched]
845 for dfCol, afwCol
in mapping.items():
846 catalog[afwCol] = view[dfCol]
848 self.log.debug(
"%d sources selected in visit %d detector %d",
850 view[
'visit'].iloc[0],
854 def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations,
855 jointcalControl, camera):
856 """Read the data that jointcal needs to run. (Gen3 version)
858 Modifies ``associations`` in-place with the loaded data.
862 inputSourceTableVisit : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
863 References to visit-level DataFrames to load the catalog data from.
864 inputVisitSummary : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
865 Visit-level exposure summary catalog with metadata.
866 associations : `lsst.jointcal.Associations`
867 Object to add the loaded data to by constructing new CcdImages.
868 jointcalControl : `jointcal.JointcalControl`
869 Control object for C++ associations management.
870 camera : `lsst.afw.cameraGeom.Camera`
871 Camera object for detector geometry.
875 oldWcsList: `list` [`lsst.afw.geom.SkyWcs`]
876 The original WCS of the input data, to aid in writing tests.
877 bands : `list` [`str`]
878 The filter bands of each input dataset.
882 load_cat_prof_file =
'jointcal_load_data.prof' if self.config.detailedProfile
else ''
883 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
886 catalogMap = {ref.dataId[
'visit']: i
for i, ref
in enumerate(inputSourceTableVisit)}
887 detectorDict = {detector.getId(): detector
for detector
in camera}
889 for visitSummaryRef
in inputVisitSummary:
890 visitSummary = visitSummaryRef.get()
891 visitCatalog = inputSourceTableVisit[catalogMap[visitSummaryRef.dataId[
'visit']]].get()
892 selected = self.sourceSelector.
run(visitCatalog)
895 detectors = {id: index
for index, id
in enumerate(visitSummary[
'id'])}
896 for id, index
in detectors.items():
898 data = self.
_make_one_input_data_make_one_input_data(visitSummary[index], catalog, detectorDict)
899 result = self.
_build_ccdImage_build_ccdImage(data, associations, jointcalControl)
900 if result
is not None:
901 oldWcsList.append(result.wcs)
903 filters.append(data.filter)
904 if len(filters) == 0:
905 raise RuntimeError(
"No data to process: did source selector remove all sources?")
906 filters = collections.Counter(filters)
908 return oldWcsList, filters
910 def _make_one_input_data(self, visitRecord, catalog, detectorDict):
911 """Return a data structure for this detector+visit."""
914 visitInfo=visitRecord.getVisitInfo(),
915 detector=detectorDict[visitRecord.getId()],
916 photoCalib=visitRecord.getPhotoCalib(),
917 wcs=visitRecord.getWcs(),
918 bbox=visitRecord.getBBox(),
922 physical=visitRecord[
'physical_filter']))
927 def _getMetadataName(self):
931 def _makeArgumentParser(cls):
932 """Create an argument parser"""
933 parser = pipeBase.ArgumentParser(name=cls.
_DefaultName_DefaultName)
934 parser.add_id_argument(
"--id",
"calexp", help=
"data ID, e.g. --id visit=6789 ccd=0..9",
935 ContainerClass=PerTractCcdDataIdContainer)
938 def _build_ccdImage(self, data, associations, jointcalControl):
940 Extract the necessary things from this catalog+metadata to add a new
945 data : `JointcalInputData`
946 The loaded input data.
947 associations : `lsst.jointcal.Associations`
948 Object to add the info to, to construct a new CcdImage
949 jointcalControl : `jointcal.JointcalControl`
950 Control object for associations management
956 The TAN WCS of this image, read from the calexp
957 (`lsst.afw.geom.SkyWcs`).
959 A key to identify this dataRef by its visit and ccd ids
962 if there are no sources in the loaded catalog.
964 if len(data.catalog) == 0:
965 self.log.warn(
"No sources selected in visit %s ccd %s", data.visit, data.detector.getId())
968 associations.createCcdImage(data.catalog,
972 data.filter.physicalLabel,
976 data.detector.getId(),
979 Result = collections.namedtuple(
'Result_from_build_CcdImage', (
'wcs',
'key'))
980 Key = collections.namedtuple(
'Key', (
'visit',
'ccd'))
981 return Result(data.wcs, Key(data.visit, data.detector.getId()))
983 def _readDataId(self, butler, dataId):
984 """Read all of the data for one dataId from the butler. (gen2 version)"""
986 if "visit" in dataId.keys():
987 visit = dataId[
"visit"]
989 visit = butler.getButler().queryMetadata(
"calexp", (
"visit"), butler.dataId)[0]
990 detector = butler.get(
'calexp_detector', dataId=dataId)
992 catalog = butler.get(
'src',
993 flags=lsst.afw.table.SOURCE_IO_NO_FOOTPRINTS,
995 goodSrc = self.sourceSelector.
run(catalog)
996 self.log.debug(
"%d sources selected in visit %d detector %d",
997 len(goodSrc.sourceCat),
1001 catalog=goodSrc.sourceCat,
1002 visitInfo=butler.get(
'calexp_visitInfo', dataId=dataId),
1004 photoCalib=butler.get(
'calexp_photoCalib', dataId=dataId),
1005 wcs=butler.get(
'calexp_wcs', dataId=dataId),
1006 bbox=butler.get(
'calexp_bbox', dataId=dataId),
1007 filter=butler.get(
'calexp_filterLabel', dataId=dataId))
1009 def loadData(self, dataRefs, associations, jointcalControl):
1010 """Read the data that jointcal needs to run. (Gen2 version)"""
1011 visit_ccd_to_dataRef = {}
1014 load_cat_prof_file =
'jointcal_loadData.prof' if self.config.detailedProfile
else ''
1015 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
1017 camera = dataRefs[0].get(
'camera', immediate=
True)
1019 for dataRef
in dataRefs:
1020 data = self.
_readDataId_readDataId(dataRef.getButler(), dataRef.dataId)
1021 result = self.
_build_ccdImage_build_ccdImage(data, associations, jointcalControl)
1024 oldWcsList.append(result.wcs)
1025 visit_ccd_to_dataRef[result.key] = dataRef
1026 filters.append(data.filter)
1027 if len(filters) == 0:
1028 raise RuntimeError(
"No data to process: did source selector remove all sources?")
1029 filters = collections.Counter(filters)
1031 return oldWcsList, filters, visit_ccd_to_dataRef
1033 def _getDebugPath(self, filename):
1034 """Constructs a path to filename using the configured debug path.
1036 return os.path.join(self.config.debugOutputPath, filename)
1038 def _prep_sky(self, associations, filters):
1039 """Prepare on-sky and other data that must be computed after data has
1042 associations.computeCommonTangentPoint()
1044 boundingCircle = associations.computeBoundingCircle()
1046 radius =
lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
1048 self.log.info(f
"Data has center={center} with radius={radius.asDegrees()} degrees.")
1051 defaultFilter = filters.most_common(1)[0][0]
1052 self.log.debug(
"Using '%s' filter for reference flux", defaultFilter.physicalLabel)
1054 return boundingCircle, center, radius, defaultFilter
1056 @pipeBase.timeMethod
1059 Jointly calibrate the astrometry and photometry across a set of images.
1061 NOTE: this is for gen2 middleware only.
1065 dataRefs : `list` of `lsst.daf.persistence.ButlerDataRef`
1066 List of data references to the exposures to be fit.
1070 result : `lsst.pipe.base.Struct`
1071 Struct of metadata from the fit, containing:
1074 The provided data references that were fit (with updated WCSs)
1076 The original WCS from each dataRef
1078 Dictionary of internally-computed metrics for testing/validation.
1080 if len(dataRefs) == 0:
1081 raise ValueError(
'Need a non-empty list of data references!')
1085 sourceFluxField =
"slot_%sFlux" % (self.config.sourceFluxType,)
1089 oldWcsList, filters, visit_ccd_to_dataRef = self.
loadDataloadData(dataRefs,
1093 boundingCircle, center, radius, defaultFilter = self.
_prep_sky_prep_sky(associations, filters)
1096 tract = dataRefs[0].dataId[
'tract']
1098 if self.config.doAstrometry:
1102 referenceSelector=self.astrometryReferenceSelector,
1110 if self.config.doPhotometry:
1114 referenceSelector=self.photometryReferenceSelector,
1118 reject_bad_fluxes=
True)
1123 return pipeBase.Struct(dataRefs=dataRefs,
1124 oldWcsList=oldWcsList,
1128 defaultFilter=defaultFilter,
1130 exitStatus=exitStatus)
1132 def _get_refcat_coordinate_error_override(self, refCat, name):
1133 """Check whether we should override the refcat coordinate errors, and
1134 return the overridden error if necessary.
1138 refCat : `lsst.afw.table.SimpleCatalog`
1139 The reference catalog to check for a ``coord_raErr`` field.
1141 Whether we are doing "astrometry" or "photometry".
1145 refCoordErr : `float`
1146 The refcat coordinate error to use, or NaN if we are not overriding
1151 lsst.pex.config.FieldValidationError
1152 Raised if the refcat does not contain coordinate errors and
1153 ``config.astrometryReferenceErr`` is not set.
1157 if name.lower() ==
"photometry":
1158 if 'coord_raErr' not in refCat.schema:
1163 if self.config.astrometryReferenceErr
is None and 'coord_raErr' not in refCat.schema:
1164 msg = (
"Reference catalog does not contain coordinate errors, "
1165 "and config.astrometryReferenceErr not supplied.")
1166 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
1170 if self.config.astrometryReferenceErr
is not None and 'coord_raErr' in refCat.schema:
1171 self.log.warn(
"Overriding reference catalog coordinate errors with %f/coordinate [mas]",
1172 self.config.astrometryReferenceErr)
1174 if self.config.astrometryReferenceErr
is None:
1177 return self.config.astrometryReferenceErr
1179 def _compute_proper_motion_epoch(self, ccdImageList):
1180 """Return the proper motion correction epoch of the provided images.
1184 ccdImageList : `list` [`lsst.jointcal.CcdImage`]
1185 The images to compute the appropriate epoch for.
1189 epoch : `astropy.time.Time`
1190 The date to use for proper motion corrections.
1192 mjds = [ccdImage.getMjd()
for ccdImage
in ccdImageList]
1193 return astropy.time.Time(np.mean(mjds), format=
'mjd', scale=
"tai")
1195 def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
1196 tract="", match_cut=3.0,
1197 reject_bad_fluxes=False, *,
1198 name="", refObjLoader=None, referenceSelector=None,
1199 fit_function=None, epoch=None):
1200 """Load reference catalog, perform the fit, and return the result.
1204 associations : `lsst.jointcal.Associations`
1205 The star/reference star associations to fit.
1206 defaultFilter : `lsst.afw.image.FilterLabel`
1207 filter to load from reference catalog.
1208 center : `lsst.geom.SpherePoint`
1209 ICRS center of field to load from reference catalog.
1210 radius : `lsst.geom.Angle`
1211 On-sky radius to load from reference catalog.
1213 Name of thing being fit: "astrometry" or "photometry".
1214 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
1215 Reference object loader to use to load a reference catalog.
1216 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1217 Selector to use to pick objects from the loaded reference catalog.
1218 fit_function : callable
1219 Function to call to perform fit (takes Associations object).
1220 tract : `str`, optional
1221 Name of tract currently being fit.
1222 match_cut : `float`, optional
1223 Radius in arcseconds to find cross-catalog matches to during
1224 associations.associateCatalogs.
1225 reject_bad_fluxes : `bool`, optional
1226 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
1227 epoch : `astropy.time.Time`, optional
1228 Epoch to which to correct refcat proper motion and parallax,
1229 or `None` to not apply such corrections.
1233 result : `Photometry` or `Astrometry`
1234 Result of `fit_function()`
1236 self.log.info(
"====== Now processing %s...", name)
1239 associations.associateCatalogs(match_cut)
1241 associations.fittedStarListSize())
1243 applyColorterms =
False if name.lower() ==
"astrometry" else self.config.applyColorTerms
1245 center, radius, defaultFilter,
1246 applyColorterms=applyColorterms,
1250 associations.collectRefStars(refCat,
1251 self.config.matchCut*lsst.geom.arcseconds,
1253 refCoordinateErr=refCoordErr,
1254 rejectBadFluxes=reject_bad_fluxes)
1256 associations.refStarListSize())
1258 associations.prepareFittedStars(self.config.minMeasurements)
1262 associations.nFittedStarsWithAssociatedRefStar())
1264 associations.fittedStarListSize())
1266 associations.nCcdImagesValidForFit())
1268 load_cat_prof_file =
'jointcal_fit_%s.prof'%name
if self.config.detailedProfile
else ''
1269 dataName =
"{}_{}".format(tract, defaultFilter.physicalLabel)
1270 with pipeBase.cmdLineTask.profile(load_cat_prof_file):
1271 result = fit_function(associations, dataName)
1274 if self.config.writeChi2FilesInitialFinal:
1275 baseName = self.
_getDebugPath_getDebugPath(f
"{name}_final_chi2-{dataName}")
1276 result.fit.saveChi2Contributions(baseName+
"{type}")
1277 self.log.info(
"Wrote chi2 contributions files: %s", baseName)
1281 def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel,
1282 applyColorterms=False, epoch=None):
1283 """Load the necessary reference catalog sources, convert fluxes to
1284 correct units, and apply color term corrections if requested.
1288 refObjLoader : `lsst.meas.algorithms.LoadReferenceObjectsTask`
1289 The reference catalog loader to use to get the data.
1290 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1291 Source selector to apply to loaded reference catalog.
1292 center : `lsst.geom.SpherePoint`
1293 The center around which to load sources.
1294 radius : `lsst.geom.Angle`
1295 The radius around ``center`` to load sources in.
1296 filterLabel : `lsst.afw.image.FilterLabel`
1297 The camera filter to load fluxes for.
1298 applyColorterms : `bool`
1299 Apply colorterm corrections to the refcat for ``filterName``?
1300 epoch : `astropy.time.Time`, optional
1301 Epoch to which to correct refcat proper motion and parallax,
1302 or `None` to not apply such corrections.
1306 refCat : `lsst.afw.table.SimpleCatalog`
1307 The loaded reference catalog.
1309 The name of the reference catalog flux field appropriate for ``filterName``.
1311 skyCircle = refObjLoader.loadSkyCircle(center,
1313 filterLabel.bandLabel,
1316 selected = referenceSelector.run(skyCircle.refCat)
1318 if not selected.sourceCat.isContiguous():
1319 refCat = selected.sourceCat.copy(deep=
True)
1321 refCat = selected.sourceCat
1324 refCatName = refObjLoader.ref_dataset_name
1325 self.log.info(
"Applying color terms for physical filter=%r reference catalog=%s",
1326 filterLabel.physicalLabel, refCatName)
1327 colorterm = self.config.colorterms.getColorterm(filterLabel.physicalLabel,
1331 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat)
1332 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
1334 refCat[skyCircle.fluxField+
'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
1336 return refCat, skyCircle.fluxField
1338 def _check_star_lists(self, associations, name):
1340 if associations.nCcdImagesValidForFit() == 0:
1341 raise RuntimeError(
'No images in the ccdImageList!')
1342 if associations.fittedStarListSize() == 0:
1343 raise RuntimeError(
'No stars in the {} fittedStarList!'.format(name))
1344 if associations.refStarListSize() == 0:
1345 raise RuntimeError(
'No stars in the {} reference star list!'.format(name))
1347 def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
1348 """Compute chi2, log it, validate the model, and return chi2.
1352 associations : `lsst.jointcal.Associations`
1353 The star/reference star associations to fit.
1354 fit : `lsst.jointcal.FitterBase`
1355 The fitter to use for minimization.
1356 model : `lsst.jointcal.Model`
1357 The model being fit.
1359 Label to describe the chi2 (e.g. "Initialized", "Final").
1360 writeChi2Name : `str`, optional
1361 Filename prefix to write the chi2 contributions to.
1362 Do not supply an extension: an appropriate one will be added.
1366 chi2: `lsst.jointcal.Chi2Accumulator`
1367 The chi2 object for the current fitter and model.
1372 Raised if chi2 is infinite or NaN.
1374 Raised if the model is not valid.
1376 if writeChi2Name
is not None:
1378 fit.saveChi2Contributions(fullpath+
"{type}")
1379 self.log.info(
"Wrote chi2 contributions files: %s", fullpath)
1381 chi2 = fit.computeChi2()
1382 self.log.info(
"%s %s", chi2Label, chi2)
1384 if not np.isfinite(chi2.chi2):
1385 raise FloatingPointError(f
'{chi2Label} chi2 is invalid: {chi2}')
1386 if not model.validate(associations.getCcdImageList(), chi2.ndof):
1387 raise ValueError(
"Model is not valid: check log messages for warnings.")
1390 def _fit_photometry(self, associations, dataName=None):
1392 Fit the photometric data.
1396 associations : `lsst.jointcal.Associations`
1397 The star/reference star associations to fit.
1399 Name of the data being processed (e.g. "1234_HSC-Y"), for
1400 identifying debugging files.
1404 fit_result : `namedtuple`
1405 fit : `lsst.jointcal.PhotometryFit`
1406 The photometric fitter used to perform the fit.
1407 model : `lsst.jointcal.PhotometryModel`
1408 The photometric model that was fit.
1410 self.log.info(
"=== Starting photometric fitting...")
1413 if self.config.photometryModel ==
"constrainedFlux":
1416 visitOrder=self.config.photometryVisitOrder,
1417 errorPedestal=self.config.photometryErrorPedestal)
1419 doLineSearch = self.config.allowLineSearch
1420 elif self.config.photometryModel ==
"constrainedMagnitude":
1423 visitOrder=self.config.photometryVisitOrder,
1424 errorPedestal=self.config.photometryErrorPedestal)
1426 doLineSearch = self.config.allowLineSearch
1427 elif self.config.photometryModel ==
"simpleFlux":
1429 errorPedestal=self.config.photometryErrorPedestal)
1430 doLineSearch =
False
1431 elif self.config.photometryModel ==
"simpleMagnitude":
1433 errorPedestal=self.config.photometryErrorPedestal)
1434 doLineSearch =
False
1439 if self.config.writeChi2FilesInitialFinal:
1440 baseName = f
"photometry_initial_chi2-{dataName}"
1443 if self.config.writeInitialModel:
1444 fullpath = self.
_getDebugPath_getDebugPath(f
"initial_photometry_model-{dataName}.txt")
1446 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialized", writeChi2Name=baseName)
1448 def getChi2Name(whatToFit):
1449 if self.config.writeChi2FilesOuterLoop:
1450 return f
"photometry_init-%s_chi2-{dataName}" % whatToFit
1456 if self.config.writeInitMatrix:
1457 dumpMatrixFile = self.
_getDebugPath_getDebugPath(f
"photometry_preinit-{dataName}")
1460 if self.config.photometryModel.startswith(
"constrained"):
1463 fit.minimize(
"ModelVisit", dumpMatrixFile=dumpMatrixFile)
1464 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize ModelVisit",
1465 writeChi2Name=getChi2Name(
"ModelVisit"))
1468 fit.minimize(
"Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
1470 writeChi2Name=getChi2Name(
"Model"))
1472 fit.minimize(
"Fluxes")
1473 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize Fluxes",
1474 writeChi2Name=getChi2Name(
"Fluxes"))
1476 fit.minimize(
"Model Fluxes", doLineSearch=doLineSearch)
1477 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize ModelFluxes",
1478 writeChi2Name=getChi2Name(
"ModelFluxes"))
1480 model.freezeErrorTransform()
1481 self.log.debug(
"Photometry error scales are frozen.")
1485 self.config.maxPhotometrySteps,
1488 doRankUpdate=self.config.photometryDoRankUpdate,
1489 doLineSearch=doLineSearch,
1496 def _fit_astrometry(self, associations, dataName=None):
1498 Fit the astrometric data.
1502 associations : `lsst.jointcal.Associations`
1503 The star/reference star associations to fit.
1505 Name of the data being processed (e.g. "1234_HSC-Y"), for
1506 identifying debugging files.
1510 fit_result : `namedtuple`
1511 fit : `lsst.jointcal.AstrometryFit`
1512 The astrometric fitter used to perform the fit.
1513 model : `lsst.jointcal.AstrometryModel`
1514 The astrometric model that was fit.
1515 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1516 The model for the sky to tangent plane projection that was used in the fit.
1519 self.log.info(
"=== Starting astrometric fitting...")
1521 associations.deprojectFittedStars()
1528 if self.config.astrometryModel ==
"constrained":
1530 sky_to_tan_projection,
1531 chipOrder=self.config.astrometryChipOrder,
1532 visitOrder=self.config.astrometryVisitOrder)
1533 elif self.config.astrometryModel ==
"simple":
1535 sky_to_tan_projection,
1536 self.config.useInputWcs,
1538 order=self.config.astrometrySimpleOrder)
1543 if self.config.writeChi2FilesInitialFinal:
1544 baseName = f
"astrometry_initial_chi2-{dataName}"
1547 if self.config.writeInitialModel:
1548 fullpath = self.
_getDebugPath_getDebugPath(f
"initial_astrometry_model-{dataName}.txt")
1550 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initial", writeChi2Name=baseName)
1552 def getChi2Name(whatToFit):
1553 if self.config.writeChi2FilesOuterLoop:
1554 return f
"astrometry_init-%s_chi2-{dataName}" % whatToFit
1558 if self.config.writeInitMatrix:
1559 dumpMatrixFile = self.
_getDebugPath_getDebugPath(f
"astrometry_preinit-{dataName}")
1564 if self.config.astrometryModel ==
"constrained":
1565 fit.minimize(
"DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1566 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize DistortionsVisit",
1567 writeChi2Name=getChi2Name(
"DistortionsVisit"))
1570 fit.minimize(
"Distortions", dumpMatrixFile=dumpMatrixFile)
1571 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize Distortions",
1572 writeChi2Name=getChi2Name(
"Distortions"))
1574 fit.minimize(
"Positions")
1575 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize Positions",
1576 writeChi2Name=getChi2Name(
"Positions"))
1578 fit.minimize(
"Distortions Positions")
1579 self.
_logChi2AndValidate_logChi2AndValidate(associations, fit, model,
"Initialize DistortionsPositions",
1580 writeChi2Name=getChi2Name(
"DistortionsPositions"))
1584 self.config.maxAstrometrySteps,
1586 "Distortions Positions",
1587 sigmaRelativeTolerance=self.config.astrometryOutlierRelativeTolerance,
1588 doRankUpdate=self.config.astrometryDoRankUpdate,
1594 return Astrometry(fit, model, sky_to_tan_projection)
1596 def _check_stars(self, associations):
1597 """Count measured and reference stars per ccd and warn/log them."""
1598 for ccdImage
in associations.getCcdImageList():
1599 nMeasuredStars, nRefStars = ccdImage.countStars()
1600 self.log.debug(
"ccdImage %s has %s measured and %s reference stars",
1601 ccdImage.getName(), nMeasuredStars, nRefStars)
1602 if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1603 self.log.warn(
"ccdImage %s has only %s measuredStars (desired %s)",
1604 ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1605 if nRefStars < self.config.minRefStarsPerCcd:
1606 self.log.warn(
"ccdImage %s has only %s RefStars (desired %s)",
1607 ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1609 def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1611 sigmaRelativeTolerance=0,
1613 doLineSearch=False):
1614 """Run fitter.minimize up to max_steps times, returning the final chi2.
1618 associations : `lsst.jointcal.Associations`
1619 The star/reference star associations to fit.
1620 fitter : `lsst.jointcal.FitterBase`
1621 The fitter to use for minimization.
1623 Maximum number of steps to run outlier rejection before declaring
1624 convergence failure.
1625 name : {'photometry' or 'astrometry'}
1626 What type of data are we fitting (for logs and debugging files).
1628 Passed to ``fitter.minimize()`` to define the parameters to fit.
1629 dataName : `str`, optional
1630 Descriptive name for this dataset (e.g. tract and filter),
1632 sigmaRelativeTolerance : `float`, optional
1633 Convergence tolerance for the fractional change in the chi2 cut
1634 level for determining outliers. If set to zero, iterations will
1635 continue until there are no outliers.
1636 doRankUpdate : `bool`, optional
1637 Do an Eigen rank update during minimization, or recompute the full
1638 matrix and gradient?
1639 doLineSearch : `bool`, optional
1640 Do a line search for the optimum step during minimization?
1644 chi2: `lsst.jointcal.Chi2Statistic`
1645 The final chi2 after the fit converges, or is forced to end.
1650 Raised if the fitter fails with a non-finite value.
1652 Raised if the fitter fails for some other reason;
1653 log messages will provide further details.
1655 if self.config.writeInitMatrix:
1656 dumpMatrixFile = self.
_getDebugPath_getDebugPath(f
"{name}_postinit-{dataName}")
1660 oldChi2.chi2 = float(
"inf")
1661 for i
in range(max_steps):
1662 if self.config.writeChi2FilesOuterLoop:
1663 writeChi2Name = f
"{name}_iterate_{i}_chi2-{dataName}"
1665 writeChi2Name =
None
1666 result = fitter.minimize(whatToFit,
1667 self.config.outlierRejectSigma,
1668 sigmaRelativeTolerance=sigmaRelativeTolerance,
1669 doRankUpdate=doRankUpdate,
1670 doLineSearch=doLineSearch,
1671 dumpMatrixFile=dumpMatrixFile)
1673 chi2 = self.
_logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(),
1674 f
"Fit iteration {i}", writeChi2Name=writeChi2Name)
1676 if result == MinimizeResult.Converged:
1678 self.log.debug(
"fit has converged - no more outliers - redo minimization "
1679 "one more time in case we have lost accuracy in rank update.")
1681 result = fitter.minimize(whatToFit, self.config.outlierRejectSigma,
1682 sigmaRelativeTolerance=sigmaRelativeTolerance)
1683 chi2 = self.
_logChi2AndValidate_logChi2AndValidate(associations, fitter, fitter.getModel(),
"Fit completed")
1686 if chi2.chi2/chi2.ndof >= 4.0:
1687 self.log.error(
"Potentially bad fit: High chi-squared/ndof.")
1690 elif result == MinimizeResult.Chi2Increased:
1691 self.log.warn(
"Still some outliers remaining but chi2 increased - retry")
1693 chi2Ratio = chi2.chi2 / oldChi2.chi2
1695 self.log.warn(
'Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1696 chi2.chi2, oldChi2.chi2, chi2Ratio)
1703 msg = (
"Large chi2 increase between steps: fit likely cannot converge."
1704 " Try setting one or more of the `writeChi2*` config fields and looking"
1705 " at how individual star chi2-values evolve during the fit.")
1706 raise RuntimeError(msg)
1708 elif result == MinimizeResult.NonFinite:
1709 filename = self.
_getDebugPath_getDebugPath(
"{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1711 fitter.saveChi2Contributions(filename+
"{type}")
1712 msg =
"Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1713 raise FloatingPointError(msg.format(filename))
1714 elif result == MinimizeResult.Failed:
1715 raise RuntimeError(
"Chi2 minimization failure, cannot complete fit.")
1717 raise RuntimeError(
"Unxepected return code from minimize().")
1719 self.log.error(
"%s failed to converge after %d steps"%(name, max_steps))
1723 def _make_output(self, ccdImageList, model, func):
1724 """Return the internal jointcal models converted to the afw
1725 structures that will be saved to disk.
1729 ccdImageList : `lsst.jointcal.CcdImageList`
1730 The list of CcdImages to get the output for.
1731 model : `lsst.jointcal.AstrometryModel` or `lsst.jointcal.PhotometryModel`
1732 The internal jointcal model to convert for each `lsst.jointcal.CcdImage`.
1734 The name of the function to call on ``model`` to get the converted
1735 structure. Must accept an `lsst.jointcal.CcdImage`.
1739 output : `dict` [`tuple`, `lsst.jointcal.AstrometryModel`] or
1740 `dict` [`tuple`, `lsst.jointcal.PhotometryModel`]
1741 The data to be saved, keyed on (visit, detector).
1744 for ccdImage
in ccdImageList:
1745 ccd = ccdImage.ccdId
1746 visit = ccdImage.visit
1747 self.log.debug(
"%s for visit: %d, ccd: %d", func, visit, ccd)
1748 output[(visit, ccd)] = getattr(model, func)(ccdImage)
1751 def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef):
1753 Write the fitted astrometric results to a new 'jointcal_wcs' dataRef.
1757 associations : `lsst.jointcal.Associations`
1758 The star/reference star associations to fit.
1759 model : `lsst.jointcal.AstrometryModel`
1760 The astrometric model that was fit.
1761 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1762 Dict of ccdImage identifiers to dataRefs that were fit.
1764 ccdImageList = associations.getCcdImageList()
1765 output = self.
_make_output_make_output(ccdImageList, model,
"makeSkyWcs")
1766 for key, skyWcs
in output.items():
1767 dataRef = visit_ccd_to_dataRef[key]
1769 dataRef.put(skyWcs,
'jointcal_wcs')
1770 except pexExceptions.Exception
as e:
1771 self.log.fatal(
'Failed to write updated Wcs: %s', str(e))
1774 def _write_photometry_results(self, associations, model, visit_ccd_to_dataRef):
1776 Write the fitted photometric results to a new 'jointcal_photoCalib' dataRef.
1780 associations : `lsst.jointcal.Associations`
1781 The star/reference star associations to fit.
1782 model : `lsst.jointcal.PhotometryModel`
1783 The photoometric model that was fit.
1784 visit_ccd_to_dataRef : `dict` of Key: `lsst.daf.persistence.ButlerDataRef`
1785 Dict of ccdImage identifiers to dataRefs that were fit.
1788 ccdImageList = associations.getCcdImageList()
1789 output = self.
_make_output_make_output(ccdImageList, model,
"toPhotoCalib")
1790 for key, photoCalib
in output.items():
1791 dataRef = visit_ccd_to_dataRef[key]
1793 dataRef.put(photoCalib,
'jointcal_photoCalib')
1794 except pexExceptions.Exception
as e:
1795 self.log.fatal(
'Failed to write updated PhotoCalib: %s', str(e))
static Schema makeMinimalSchema()
static std::shared_ptr< SourceTable > make(Schema const &schema, std::shared_ptr< IdFactory > const &idFactory)
static Schema makeMinimalSchema()
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 __init__(self, *config=None)
def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations, jointcalControl, camera)
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 _make_one_input_data(self, visitRecord, catalog, detectorDict)
def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", sigmaRelativeTolerance=0, doRankUpdate=True, doLineSearch=False)
def _write_astrometry_results(self, associations, model, visit_ccd_to_dataRef)
def _check_stars(self, associations)
def _prep_sky(self, associations, filters)
def _readDataId(self, butler, dataId)
def runQuantum(self, butlerQC, inputRefs, outputRefs)
def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius, tract="", match_cut=3.0, reject_bad_fluxes=False, *name="", refObjLoader=None, referenceSelector=None, fit_function=None, epoch=None)
def _build_ccdImage(self, data, associations, jointcalControl)
def runDataRef(self, dataRefs)
def __init__(self, butler=None, initInputs=None, **kwargs)
def _extract_detector_catalog_from_visit_catalog(self, table, visitCatalog, detectorId)
def _make_schema_table(self)
def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None)
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 _make_output(self, ccdImageList, model, func)
def _put_output(self, butlerQC, outputs, outputRefs, camera, setter)
def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel, applyColorterms=False, epoch=None)
def _fit_astrometry(self, associations, dataName=None)
def loadData(self, dataRefs, associations, jointcalControl)
def lookupStaticCalibrations(datasetType, registry, quantumDataId, collections)
def add_measurement(job, name, value)
def lookupVisitRefCats(datasetType, registry, quantumDataId, collections)
def writeModel(model, filename, log)