29import astropy.units
as u
39from lsst.verify
import Job, Measurement
42 LoadReferenceObjectsConfig)
48__all__ = [
"JointcalConfig",
"JointcalTask"]
50Photometry = collections.namedtuple(
'Photometry', (
'fit',
'model'))
51Astrometry = collections.namedtuple(
'Astrometry', (
'fit',
'model',
'sky_to_tan_projection'))
56 meas = Measurement(job.metrics[name], value)
57 job.measurements.insert(meas)
61 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter")):
62 """Middleware input/output connections for jointcal data."""
63 inputCamera = pipeBase.connectionTypes.PrerequisiteInput(
64 doc=
"The camera instrument that took these observations.",
66 storageClass=
"Camera",
67 dimensions=(
"instrument",),
70 inputSourceTableVisit = pipeBase.connectionTypes.Input(
71 doc=
"Source table in parquet format, per visit",
72 name=
"sourceTable_visit",
73 storageClass=
"DataFrame",
74 dimensions=(
"instrument",
"visit"),
78 inputVisitSummary = pipeBase.connectionTypes.Input(
79 doc=(
"Per-visit consolidated exposure metadata built from calexps. "
80 "These catalogs use detector id for the id and must be sorted for "
81 "fast lookups of a detector."),
83 storageClass=
"ExposureCatalog",
84 dimensions=(
"instrument",
"visit"),
88 astrometryRefCat = pipeBase.connectionTypes.PrerequisiteInput(
89 doc=
"The astrometry reference catalog to match to loaded input catalog sources.",
90 name=
"gaia_dr2_20200414",
91 storageClass=
"SimpleCatalog",
92 dimensions=(
"skypix",),
96 photometryRefCat = pipeBase.connectionTypes.PrerequisiteInput(
97 doc=
"The photometry reference catalog to match to loaded input catalog sources.",
98 name=
"ps1_pv3_3pi_20170110",
99 storageClass=
"SimpleCatalog",
100 dimensions=(
"skypix",),
105 outputWcs = pipeBase.connectionTypes.Output(
106 doc=(
"Per-tract, per-visit world coordinate systems derived from the fitted model."
107 " These catalogs only contain entries for detectors with an output, and use"
108 " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
109 name=
"jointcalSkyWcsCatalog",
110 storageClass=
"ExposureCatalog",
111 dimensions=(
"instrument",
"visit",
"skymap",
"tract"),
114 outputPhotoCalib = pipeBase.connectionTypes.Output(
115 doc=(
"Per-tract, per-visit photometric calibrations derived from the fitted model."
116 " These catalogs only contain entries for detectors with an output, and use"
117 " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
118 name=
"jointcalPhotoCalibCatalog",
119 storageClass=
"ExposureCatalog",
120 dimensions=(
"instrument",
"visit",
"skymap",
"tract"),
128 for name
in (
"astrometry",
"photometry"):
129 vars()[f
"{name}_matched_fittedStars"] = pipeBase.connectionTypes.Output(
130 doc=f
"The number of cross-matched fittedStars for {name}",
131 name=f
"metricvalue_jointcal_{name}_matched_fittedStars",
132 storageClass=
"MetricValue",
133 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter"),
135 vars()[f
"{name}_collected_refStars"] = pipeBase.connectionTypes.Output(
136 doc=f
"The number of {name} reference stars drawn from the reference catalog, before matching.",
137 name=f
"metricvalue_jointcal_{name}_collected_refStars",
138 storageClass=
"MetricValue",
139 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter"),
141 vars()[f
"{name}_prepared_refStars"] = pipeBase.connectionTypes.Output(
142 doc=f
"The number of {name} reference stars matched to fittedStars.",
143 name=f
"metricvalue_jointcal_{name}_prepared_refStars",
144 storageClass=
"MetricValue",
145 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter"),
147 vars()[f
"{name}_prepared_fittedStars"] = pipeBase.connectionTypes.Output(
148 doc=f
"The number of cross-matched fittedStars after cleanup, for {name}.",
149 name=f
"metricvalue_jointcal_{name}_prepared_fittedStars",
150 storageClass=
"MetricValue",
151 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter"),
153 vars()[f
"{name}_prepared_ccdImages"] = pipeBase.connectionTypes.Output(
154 doc=f
"The number of ccdImages that will be fit for {name}, after cleanup.",
155 name=f
"metricvalue_jointcal_{name}_prepared_ccdImages",
156 storageClass=
"MetricValue",
157 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter"),
159 vars()[f
"{name}_final_chi2"] = pipeBase.connectionTypes.Output(
160 doc=f
"The final chi2 of the {name} fit.",
161 name=f
"metricvalue_jointcal_{name}_final_chi2",
162 storageClass=
"MetricValue",
163 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter"),
165 vars()[f
"{name}_final_ndof"] = pipeBase.connectionTypes.Output(
166 doc=f
"The number of degrees of freedom of the fitted {name}.",
167 name=f
"metricvalue_jointcal_{name}_final_ndof",
168 storageClass=
"MetricValue",
169 dimensions=(
"skymap",
"tract",
"instrument",
"physical_filter"),
180 if not config.doAstrometry:
181 self.prerequisiteInputs.
remove(
"astrometryRefCat")
184 if "metricvalue_jointcal_astrometry" in key:
186 if not config.doPhotometry:
187 self.prerequisiteInputs.
remove(
"photometryRefCat")
190 if "metricvalue_jointcal_photometry" in key:
194 return (
"inputVisitSummary",)
198 pipelineConnections=JointcalTaskConnections):
199 """Configuration for JointcalTask"""
201 doAstrometry = pexConfig.Field(
202 doc=
"Fit astrometry and write the fitted result.",
206 doPhotometry = pexConfig.Field(
207 doc=
"Fit photometry and write the fitted result.",
211 sourceFluxType = pexConfig.Field(
213 doc=
"Source flux field to use in source selection and to get fluxes from the catalog.",
214 default=
'apFlux_12_0'
216 positionErrorPedestal = pexConfig.Field(
217 doc=
"Systematic term to apply to the measured position error (pixels)",
221 photometryErrorPedestal = pexConfig.Field(
222 doc=
"Systematic term to apply to the measured error on flux or magnitude as a "
223 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
228 matchCut = pexConfig.Field(
229 doc=
"Matching radius between fitted and reference stars (arcseconds)",
233 minMeasurements = pexConfig.Field(
234 doc=
"Minimum number of associated measured stars for a fitted star to be included in the fit",
238 minMeasuredStarsPerCcd = pexConfig.Field(
239 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
243 minRefStarsPerCcd = pexConfig.Field(
244 doc=
"Minimum number of measuredStars per ccdImage before printing warnings",
248 allowLineSearch = pexConfig.Field(
249 doc=
"Allow a line search during minimization, if it is reasonable for the model"
250 " (models with a significant non-linear component, e.g. constrainedPhotometry).",
254 astrometrySimpleOrder = pexConfig.Field(
255 doc=
"Polynomial order for fitting the simple astrometry model.",
259 astrometryChipOrder = pexConfig.Field(
260 doc=
"Order of the per-chip transform for the constrained astrometry model.",
264 astrometryVisitOrder = pexConfig.Field(
265 doc=
"Order of the per-visit transform for the constrained astrometry model.",
269 useInputWcs = pexConfig.Field(
270 doc=
"Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
274 astrometryModel = pexConfig.ChoiceField(
275 doc=
"Type of model to fit to astrometry",
277 default=
"constrained",
278 allowed={
"simple":
"One polynomial per ccd",
279 "constrained":
"One polynomial per ccd, and one polynomial per visit"}
281 photometryModel = pexConfig.ChoiceField(
282 doc=
"Type of model to fit to photometry",
284 default=
"constrainedMagnitude",
285 allowed={
"simpleFlux":
"One constant zeropoint per ccd and visit, fitting in flux space.",
286 "constrainedFlux":
"Constrained zeropoint per ccd, and one polynomial per visit,"
287 " fitting in flux space.",
288 "simpleMagnitude":
"One constant zeropoint per ccd and visit,"
289 " fitting in magnitude space.",
290 "constrainedMagnitude":
"Constrained zeropoint per ccd, and one polynomial per visit,"
291 " fitting in magnitude space.",
294 applyColorTerms = pexConfig.Field(
295 doc=
"Apply photometric color terms to reference stars?"
296 "Requires that colorterms be set to a ColortermLibrary",
300 colorterms = pexConfig.ConfigField(
301 doc=
"Library of photometric reference catalog name to color term dict.",
302 dtype=ColortermLibrary,
304 photometryVisitOrder = pexConfig.Field(
305 doc=
"Order of the per-visit polynomial transform for the constrained photometry model.",
309 photometryDoRankUpdate = pexConfig.Field(
310 doc=(
"Do the rank update step during minimization. "
311 "Skipping this can help deal with models that are too non-linear."),
315 astrometryDoRankUpdate = pexConfig.Field(
316 doc=(
"Do the rank update step during minimization (should not change the astrometry fit). "
317 "Skipping this can help deal with models that are too non-linear."),
321 outlierRejectSigma = pexConfig.Field(
322 doc=
"How many sigma to reject outliers at during minimization.",
326 astrometryOutlierRelativeTolerance = pexConfig.Field(
327 doc=(
"Convergence tolerance for outlier rejection threshold when fitting astrometry. Iterations will "
328 "stop when the fractional change in the chi2 cut level is below this value. If tolerance is set "
329 "to zero, iterations will continue until there are no more outliers. We suggest a value of 0.002"
330 "as a balance between a shorter minimization runtime and achieving a final fitted model that is"
331 "close to the solution found when removing all outliers."),
335 maxPhotometrySteps = pexConfig.Field(
336 doc=
"Maximum number of minimize iterations to take when fitting photometry.",
340 maxAstrometrySteps = pexConfig.Field(
341 doc=
"Maximum number of minimize iterations to take when fitting astrometry.",
345 astrometryRefObjLoader = pexConfig.ConfigField(
346 dtype=LoadReferenceObjectsConfig,
347 doc=
"Reference object loader for astrometric fit",
349 photometryRefObjLoader = pexConfig.ConfigField(
350 dtype=LoadReferenceObjectsConfig,
351 doc=
"Reference object loader for photometric fit",
353 sourceSelector = sourceSelectorRegistry.makeField(
354 doc=
"How to select sources for cross-matching",
357 astrometryReferenceSelector = pexConfig.ConfigurableField(
358 target=ReferenceSourceSelectorTask,
359 doc=
"How to down-select the loaded astrometry reference catalog.",
361 photometryReferenceSelector = pexConfig.ConfigurableField(
362 target=ReferenceSourceSelectorTask,
363 doc=
"How to down-select the loaded photometry reference catalog.",
365 astrometryReferenceErr = pexConfig.Field(
366 doc=(
"Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
367 "If None, then raise an exception if the reference catalog is missing coordinate errors. "
368 "If specified, overrides any existing `coord_*Err` values."),
375 writeInitMatrix = pexConfig.Field(
377 doc=(
"Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
378 "Output files will be written to `config.debugOutputPath` and will "
379 "be of the form 'astrometry_[pre|post]init-TRACT-FILTER-mat.txt'. "
380 "Note that these files are the dense versions of the matrix, and so may be very large."),
383 writeChi2FilesInitialFinal = pexConfig.Field(
385 doc=(
"Write .csv files containing the contributions to chi2 for the initialization and final fit. "
386 "Output files will be written to `config.debugOutputPath` and will "
387 "be of the form `astrometry_[initial|final]_chi2-TRACT-FILTER."),
390 writeChi2FilesOuterLoop = pexConfig.Field(
392 doc=(
"Write .csv files containing the contributions to chi2 for the outer fit loop. "
393 "Output files will be written to `config.debugOutputPath` and will "
394 "be of the form `astrometry_init-NN_chi2-TRACT-FILTER`."),
397 writeInitialModel = pexConfig.Field(
399 doc=(
"Write the pre-initialization model to text files, for debugging. "
400 "Output files will be written to `config.debugOutputPath` and will be "
401 "of the form `initial_astrometry_model-TRACT_FILTER.txt`."),
404 debugOutputPath = pexConfig.Field(
407 doc=(
"Path to write debug output files to. Used by "
408 "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
410 detailedProfile = pexConfig.Field(
413 doc=
"Output separate profiling information for different parts of jointcal, e.g. data read, fitting"
419 msg =
"applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
420 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
422 msg = (
"Only doing astrometry, but Colorterms are not applied for astrometry;"
423 "applyColorTerms=True will be ignored.")
424 logging.getLogger(
"lsst.jointcal").warning(msg)
431 self.
sourceSelector[
"science"].unresolved.name =
"sizeExtendedness"
438 self.
sourceSelector[
"science"].signalToNoise.fluxField = f
"{self.sourceFluxType}_instFlux"
439 self.
sourceSelector[
"science"].signalToNoise.errField = f
"{self.sourceFluxType}_instFluxErr"
442 self.
sourceSelector[
"science"].isolated.parentName =
"parentSourceId"
443 self.
sourceSelector[
"science"].isolated.nChildName =
"deblend_nChild"
447 badFlags = [
"pixelFlags_edge",
448 "pixelFlags_saturated",
449 "pixelFlags_interpolatedCenter",
450 "pixelFlags_interpolated",
451 "pixelFlags_crCenter",
453 "hsmPsfMoments_flag",
454 f
"{self.sourceFluxType}_flag",
458 self.
sourceSelector[
"science"].requireFiniteRaDec.raColName =
"ra"
459 self.
sourceSelector[
"science"].requireFiniteRaDec.decColName =
"dec"
468 """Write model to outfile."""
469 with open(filename,
"w")
as file:
470 file.write(repr(model))
471 log.info(
"Wrote %s to file: %s", model, filename)
474@dataclasses.dataclass
476 """The input data jointcal needs for each detector/visit."""
478 """The visit identifier of this exposure."""
480 """The catalog derived from this exposure."""
482 """The VisitInfo of this exposure."""
484 """The detector of this exposure."""
486 """The photometric calibration of this exposure."""
488 """The WCS of this exposure."""
490 """The bounding box of this exposure."""
492 """The filter of this exposure."""
496 """Astrometricly and photometricly calibrate across multiple visits of the
500 ConfigClass = JointcalConfig
501 _DefaultName =
"jointcal"
504 super().__init__(**kwargs)
505 self.makeSubtask(
"sourceSelector")
506 if self.
config.doAstrometry:
507 self.makeSubtask(
"astrometryReferenceSelector")
510 if self.
config.doPhotometry:
511 self.makeSubtask(
"photometryReferenceSelector")
516 self.
job = Job.load_metrics_package(subset=
'jointcal')
521 inputs = butlerQC.get(inputRefs)
523 tract = butlerQC.quantum.dataId[
'tract']
524 if self.
config.doAstrometry:
526 dataIds=[ref.datasetRef.dataId
527 for ref
in inputRefs.astrometryRefCat],
528 refCats=inputs.pop(
'astrometryRefCat'),
529 config=self.
config.astrometryRefObjLoader,
530 name=self.
config.connections.astrometryRefCat,
532 if self.
config.doPhotometry:
534 dataIds=[ref.datasetRef.dataId
535 for ref
in inputRefs.photometryRefCat],
536 refCats=inputs.pop(
'photometryRefCat'),
537 config=self.
config.photometryRefObjLoader,
538 name=self.
config.connections.photometryRefCat,
540 outputs = self.
run(**inputs, tract=tract)
542 if self.
config.doAstrometry:
543 self.
_put_output(butlerQC, outputs.outputWcs, outputRefs.outputWcs,
544 inputs[
'inputCamera'],
"setWcs")
545 if self.
config.doPhotometry:
546 self.
_put_output(butlerQC, outputs.outputPhotoCalib, outputRefs.outputPhotoCalib,
547 inputs[
'inputCamera'],
"setPhotoCalib")
550 """Persist all measured metrics stored in a job.
554 butlerQC : `lsst.pipe.base.QuantumContext`
555 A butler which is specialized to operate in the context of a
556 `lsst.daf.butler.Quantum`; This is the input to `runQuantum`.
557 job : `lsst.verify.job.Job`
558 Measurements of metrics to persist.
559 outputRefs : `list` [`lsst.pipe.base.OutputQuantizedConnection`]
560 The DatasetRefs to persist the data to.
562 for key
in job.measurements.keys():
563 butlerQC.put(job.measurements[key], getattr(outputRefs, key.fqn.replace(
'jointcal.',
'')))
565 def _put_output(self, butlerQC, outputs, outputRefs, camera, setter):
566 """Persist the output datasets to their appropriate datarefs.
570 butlerQC : `lsst.pipe.base.QuantumContext`
571 A butler which is specialized to operate in the context of a
572 `lsst.daf.butler.Quantum`; This is the input to `runQuantum`.
573 outputs : `dict` [`tuple`, `lsst.afw.geom.SkyWcs`] or
574 `dict` [`tuple, `lsst.afw.image.PhotoCalib`]
575 The fitted objects to persist.
576 outputRefs : `list` [`lsst.pipe.base.OutputQuantizedConnection`]
577 The DatasetRefs to persist the data to.
578 camera : `lsst.afw.cameraGeom.Camera`
579 The camera for this instrument, to get detector ids from.
581 The method to call on the ExposureCatalog to set each output.
584 schema.addField(
'visit', type=
'L', doc=
'Visit number')
586 def new_catalog(visit, size):
587 """Return an catalog ready to be filled with appropriate output."""
590 catalog[
'visit'] = visit
592 metadata.add(
"COMMENT",
"Catalog id is detector id, sorted.")
593 metadata.add(
"COMMENT",
"Only detectors with data have entries.")
597 detectors_per_visit = collections.defaultdict(int)
600 detectors_per_visit[key[0]] += 1
602 for ref
in outputRefs:
603 visit = ref.dataId[
'visit']
604 catalog = new_catalog(visit, detectors_per_visit[visit])
607 for detector
in camera:
608 detectorId = detector.getId()
609 key = (ref.dataId[
'visit'], detectorId)
610 if key
not in outputs:
612 self.
log.debug(
"No %s output for detector %s in visit %s",
613 setter[3:], detectorId, visit)
616 catalog[i].setId(detectorId)
617 getattr(catalog[i], setter)(outputs[key])
621 butlerQC.put(catalog, ref)
622 self.
log.info(
"Wrote %s detectors to %s", i, ref)
624 def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None):
630 sourceFluxField =
"flux"
634 oldWcsList, bands = self.
_load_data(inputSourceTableVisit,
640 boundingCircle, center, radius, defaultFilter, epoch = self.
_prep_sky(associations, bands)
642 if self.
config.doAstrometry:
646 referenceSelector=self.astrometryReferenceSelector,
650 astrometry_output = self.
_make_output(associations.getCcdImageList(),
654 astrometry_output =
None
656 if self.
config.doPhotometry:
660 referenceSelector=self.photometryReferenceSelector,
664 reject_bad_fluxes=
True)
665 photometry_output = self.
_make_output(associations.getCcdImageList(),
669 photometry_output =
None
671 return pipeBase.Struct(outputWcs=astrometry_output,
672 outputPhotoCalib=photometry_output,
677 def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations,
678 jointcalControl, camera):
679 """Read the data that jointcal needs to run.
681 Modifies ``associations`` in-place with the loaded data.
685 inputSourceTableVisit : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
686 References to visit-level DataFrames to load the catalog data from.
687 inputVisitSummary : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
688 Visit-level exposure summary catalog with metadata.
689 associations : `lsst.jointcal.Associations`
690 Object to add the loaded data to by constructing new CcdImages.
691 jointcalControl : `jointcal.JointcalControl`
692 Control object for C++ associations management.
693 camera : `lsst.afw.cameraGeom.Camera`
694 Camera object for detector geometry.
698 oldWcsList: `list` [`lsst.afw.geom.SkyWcs`]
699 The original WCS of the input data, to aid in writing tests.
700 bands : `list` [`str`]
701 The filter bands of each input dataset.
705 load_cat_profile_file =
'jointcal_load_data.prof' if self.
config.detailedProfile
else ''
706 with lsst.utils.timer.profile(load_cat_profile_file):
709 catalogMap = {ref.dataId[
'visit']: i
for i, ref
in enumerate(inputSourceTableVisit)}
710 detectorDict = {detector.getId(): detector
for detector
in camera}
714 for visitSummaryRef
in inputVisitSummary:
715 visitSummary = visitSummaryRef.get()
717 dataRef = inputSourceTableVisit[catalogMap[visitSummaryRef.dataId[
'visit']]]
719 inColumns = dataRef.get(component=
'columns')
723 visitCatalog = dataRef.get(parameters={
'columns': columns})
726 if len(selected.sourceCat) == 0:
727 self.
log.warning(
"No sources selected in visit %s. Skipping...",
728 visitSummary[
"visit"][0])
732 detectors = {id: index
for index, id
in enumerate(visitSummary[
'id'])}
733 for id, index
in detectors.items():
738 self.
config.sourceFluxType,
744 if result
is not None:
745 oldWcsList.append(result.wcs)
747 filters.append(data.filter)
748 if len(filters) == 0:
749 raise RuntimeError(
"No data to process: did source selector remove all sources?")
750 filters = collections.Counter(filters)
752 return oldWcsList, filters
755 """Return a data structure for this detector+visit."""
758 visitInfo=visitRecord.getVisitInfo(),
759 detector=detectorDict[visitRecord.getId()],
760 photoCalib=visitRecord.getPhotoCalib(),
761 wcs=visitRecord.getWcs(),
762 bbox=visitRecord.getBBox(),
766 physical=visitRecord[
'physical_filter']))
770 Extract the necessary things from this catalog+metadata to add a new
775 data : `JointcalInputData`
776 The loaded input data.
777 associations : `lsst.jointcal.Associations`
778 Object to add the info to, to construct a new CcdImage
779 jointcalControl : `jointcal.JointcalControl`
780 Control object for associations management
786 The TAN WCS of this image, read from the calexp
787 (`lsst.afw.geom.SkyWcs`).
789 A key to identify this dataRef by its visit and ccd ids
792 if there are no sources in the loaded catalog.
794 if len(data.catalog) == 0:
795 self.
log.warning(
"No sources selected in visit %s ccd %s", data.visit, data.detector.getId())
798 associations.createCcdImage(data.catalog,
802 data.filter.physicalLabel,
806 data.detector.getId(),
809 Result = collections.namedtuple(
'Result_from_build_CcdImage', (
'wcs',
'key'))
810 Key = collections.namedtuple(
'Key', (
'visit',
'ccd'))
811 return Result(data.wcs, Key(data.visit, data.detector.getId()))
814 """Constructs a path to filename using the configured debug path.
816 return os.path.join(self.
config.debugOutputPath, filename)
819 """Prepare on-sky and other data that must be computed after data has
822 associations.computeCommonTangentPoint()
824 boundingCircle = associations.computeBoundingCircle()
826 radius =
lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
828 self.
log.info(f
"Data has center={center} with radius={radius.asDegrees()} degrees.")
831 defaultFilter = filters.most_common(1)[0][0]
832 self.
log.debug(
"Using '%s' filter for reference flux", defaultFilter.physicalLabel)
836 associations.setEpoch(epoch.jyear)
838 return boundingCircle, center, radius, defaultFilter, epoch
841 """Check whether we should override the refcat coordinate errors, and
842 return the overridden error if necessary.
846 refCat : `lsst.afw.table.SimpleCatalog`
847 The reference catalog to check for a ``coord_raErr`` field.
849 Whether we are doing "astrometry" or "photometry".
853 refCoordErr : `float`
854 The refcat coordinate error to use, or NaN if we are not overriding
859 lsst.pex.config.FieldValidationError
860 Raised if the refcat does not contain coordinate errors and
861 ``config.astrometryReferenceErr`` is not set.
865 if name.lower() ==
"photometry":
866 if 'coord_raErr' not in refCat.schema:
871 if self.
config.astrometryReferenceErr
is None and 'coord_raErr' not in refCat.schema:
872 msg = (
"Reference catalog does not contain coordinate errors, "
873 "and config.astrometryReferenceErr not supplied.")
874 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
878 if self.
config.astrometryReferenceErr
is not None and 'coord_raErr' in refCat.schema:
879 self.
log.warning(
"Overriding reference catalog coordinate errors with %f/coordinate [mas]",
880 self.
config.astrometryReferenceErr)
882 if self.
config.astrometryReferenceErr
is None:
885 return self.
config.astrometryReferenceErr
888 """Return the proper motion correction epoch of the provided images.
892 ccdImageList : `list` [`lsst.jointcal.CcdImage`]
893 The images to compute the appropriate epoch for.
897 epoch : `astropy.time.Time`
898 The date to use for proper motion corrections.
900 return astropy.time.Time(np.mean([ccdImage.getEpoch()
for ccdImage
in ccdImageList]),
905 tract="", match_cut=3.0,
906 reject_bad_fluxes=False, *,
907 name="", refObjLoader=None, referenceSelector=None,
908 fit_function=None, epoch=None):
909 """Load reference catalog, perform the fit, and return the result.
913 associations : `lsst.jointcal.Associations`
914 The star/reference star associations to fit.
915 defaultFilter : `lsst.afw.image.FilterLabel`
916 filter to load from reference catalog.
917 center : `lsst.geom.SpherePoint`
918 ICRS center of field to load from reference catalog.
919 radius : `lsst.geom.Angle`
920 On-sky radius to load from reference catalog.
922 Name of thing being fit: "astrometry" or "photometry".
923 refObjLoader : `lsst.meas.algorithms.ReferenceObjectLoader`
924 Reference object loader to use to load a reference catalog.
925 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
926 Selector to use to pick objects from the loaded reference catalog.
927 fit_function : callable
928 Function to call to perform fit (takes Associations object).
929 tract : `str`, optional
930 Name of tract currently being fit.
931 match_cut : `float`, optional
932 Radius in arcseconds to find cross-catalog matches to during
933 associations.associateCatalogs.
934 reject_bad_fluxes : `bool`, optional
935 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
936 epoch : `astropy.time.Time`, optional
937 Epoch to which to correct refcat proper motion and parallax,
938 or `None` to not apply such corrections.
942 result : `Photometry` or `Astrometry`
943 Result of `fit_function()`
945 self.
log.info(
"====== Now processing %s...", name)
948 associations.associateCatalogs(match_cut)
950 associations.fittedStarListSize())
952 applyColorterms =
False if name.lower() ==
"astrometry" else self.
config.applyColorTerms
954 center, radius, defaultFilter,
955 applyColorterms=applyColorterms,
959 associations.collectRefStars(refCat,
960 self.
config.matchCut*lsst.geom.arcseconds,
962 refCoordinateErr=refCoordErr,
963 rejectBadFluxes=reject_bad_fluxes)
965 associations.refStarListSize())
967 associations.prepareFittedStars(self.
config.minMeasurements)
971 associations.nFittedStarsWithAssociatedRefStar())
973 associations.fittedStarListSize())
975 associations.nCcdImagesValidForFit())
977 fit_profile_file =
'jointcal_fit_%s.prof'%name
if self.
config.detailedProfile
else ''
978 dataName =
"{}_{}".format(tract, defaultFilter.physicalLabel)
979 with lsst.utils.timer.profile(fit_profile_file):
980 result = fit_function(associations, dataName)
983 if self.
config.writeChi2FilesInitialFinal:
984 baseName = self.
_getDebugPath(f
"{name}_final_chi2-{dataName}")
985 result.fit.saveChi2Contributions(baseName+
"{type}")
986 self.
log.info(
"Wrote chi2 contributions files: %s", baseName)
991 applyColorterms=False, epoch=None):
992 """Load the necessary reference catalog sources, convert fluxes to
993 correct units, and apply color term corrections if requested.
997 refObjLoader : `lsst.meas.algorithms.ReferenceObjectLoader`
998 The reference catalog loader to use to get the data.
999 referenceSelector : `lsst.meas.algorithms.ReferenceSourceSelectorTask`
1000 Source selector to apply to loaded reference catalog.
1001 center : `lsst.geom.SpherePoint`
1002 The center around which to load sources.
1003 radius : `lsst.geom.Angle`
1004 The radius around ``center`` to load sources in.
1005 filterLabel : `lsst.afw.image.FilterLabel`
1006 The camera filter to load fluxes for.
1007 applyColorterms : `bool`
1008 Apply colorterm corrections to the refcat for ``filterName``?
1009 epoch : `astropy.time.Time`, optional
1010 Epoch to which to correct refcat proper motion and parallax,
1011 or `None` to not apply such corrections.
1015 refCat : `lsst.afw.table.SimpleCatalog`
1016 The loaded reference catalog.
1018 The name of the reference catalog flux field appropriate for ``filterName``.
1020 skyCircle = refObjLoader.loadSkyCircle(center,
1022 filterLabel.bandLabel,
1025 selected = referenceSelector.run(skyCircle.refCat)
1027 if not selected.sourceCat.isContiguous():
1028 refCat = selected.sourceCat.copy(deep=
True)
1030 refCat = selected.sourceCat
1033 refCatName = refObjLoader.name
1034 self.
log.info(
"Applying color terms for physical filter=%r reference catalog=%s",
1035 filterLabel.physicalLabel, refCatName)
1036 colorterm = self.
config.colorterms.getColorterm(filterLabel.physicalLabel,
1040 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat)
1041 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
1043 refCat[skyCircle.fluxField+
'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
1045 return refCat, skyCircle.fluxField
1049 if associations.nCcdImagesValidForFit() == 0:
1050 raise RuntimeError(
'No images in the ccdImageList!')
1051 if associations.fittedStarListSize() == 0:
1052 raise RuntimeError(
'No stars in the {} fittedStarList!'.format(name))
1053 if associations.refStarListSize() == 0:
1054 raise RuntimeError(
'No stars in the {} reference star list!'.format(name))
1057 """Compute chi2, log it, validate the model, and return chi2.
1061 associations : `lsst.jointcal.Associations`
1062 The star/reference star associations to fit.
1063 fit : `lsst.jointcal.FitterBase`
1064 The fitter to use for minimization.
1065 model : `lsst.jointcal.Model`
1066 The model being fit.
1068 Label to describe the chi2 (e.g. "Initialized", "Final").
1069 writeChi2Name : `str`, optional
1070 Filename prefix to write the chi2 contributions to.
1071 Do not supply an extension: an appropriate one will be added.
1075 chi2: `lsst.jointcal.Chi2Accumulator`
1076 The chi2 object for the current fitter and model.
1081 Raised if chi2 is infinite or NaN.
1083 Raised if the model is not valid.
1085 if writeChi2Name
is not None:
1087 fit.saveChi2Contributions(fullpath+
"{type}")
1088 self.
log.info(
"Wrote chi2 contributions files: %s", fullpath)
1090 chi2 = fit.computeChi2()
1091 self.
log.info(
"%s %s", chi2Label, chi2)
1093 if not np.isfinite(chi2.chi2):
1094 raise FloatingPointError(f
'{chi2Label} chi2 is invalid: {chi2}')
1095 if not model.validate(associations.getCcdImageList(), chi2.ndof):
1096 raise ValueError(
"Model is not valid: check log messages for warnings.")
1101 Fit the photometric data.
1105 associations : `lsst.jointcal.Associations`
1106 The star/reference star associations to fit.
1108 Name of the data being processed (e.g. "1234_HSC-Y"), for
1109 identifying debugging files.
1113 fit_result : `namedtuple`
1114 fit : `lsst.jointcal.PhotometryFit`
1115 The photometric fitter used to perform the fit.
1116 model : `lsst.jointcal.PhotometryModel`
1117 The photometric model that was fit.
1119 self.
log.info(
"=== Starting photometric fitting...")
1122 if self.
config.photometryModel ==
"constrainedFlux":
1125 visitOrder=self.
config.photometryVisitOrder,
1126 errorPedestal=self.
config.photometryErrorPedestal)
1128 doLineSearch = self.
config.allowLineSearch
1129 elif self.
config.photometryModel ==
"constrainedMagnitude":
1132 visitOrder=self.
config.photometryVisitOrder,
1133 errorPedestal=self.
config.photometryErrorPedestal)
1135 doLineSearch = self.
config.allowLineSearch
1136 elif self.
config.photometryModel ==
"simpleFlux":
1138 errorPedestal=self.
config.photometryErrorPedestal)
1139 doLineSearch =
False
1140 elif self.
config.photometryModel ==
"simpleMagnitude":
1142 errorPedestal=self.
config.photometryErrorPedestal)
1143 doLineSearch =
False
1148 if self.
config.writeChi2FilesInitialFinal:
1149 baseName = f
"photometry_initial_chi2-{dataName}"
1152 if self.
config.writeInitialModel:
1153 fullpath = self.
_getDebugPath(f
"initial_photometry_model-{dataName}.txt")
1157 def getChi2Name(whatToFit):
1158 if self.
config.writeChi2FilesOuterLoop:
1159 return f
"photometry_init-%s_chi2-{dataName}" % whatToFit
1165 if self.
config.writeInitMatrix:
1166 dumpMatrixFile = self.
_getDebugPath(f
"photometry_preinit-{dataName}")
1169 if self.
config.photometryModel.startswith(
"constrained"):
1172 fit.minimize(
"ModelVisit", dumpMatrixFile=dumpMatrixFile)
1174 writeChi2Name=getChi2Name(
"ModelVisit"))
1177 fit.minimize(
"Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
1179 writeChi2Name=getChi2Name(
"Model"))
1181 fit.minimize(
"Fluxes")
1183 writeChi2Name=getChi2Name(
"Fluxes"))
1185 fit.minimize(
"Model Fluxes", doLineSearch=doLineSearch)
1187 writeChi2Name=getChi2Name(
"ModelFluxes"))
1189 model.freezeErrorTransform()
1190 self.
log.debug(
"Photometry error scales are frozen.")
1194 self.
config.maxPhotometrySteps,
1197 doRankUpdate=self.
config.photometryDoRankUpdate,
1198 doLineSearch=doLineSearch,
1207 Fit the astrometric data.
1211 associations : `lsst.jointcal.Associations`
1212 The star/reference star associations to fit.
1214 Name of the data being processed (e.g. "1234_HSC-Y"), for
1215 identifying debugging files.
1219 fit_result : `namedtuple`
1220 fit : `lsst.jointcal.AstrometryFit`
1221 The astrometric fitter used to perform the fit.
1222 model : `lsst.jointcal.AstrometryModel`
1223 The astrometric model that was fit.
1224 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1225 The model for the sky to tangent plane projection that was used in the fit.
1228 self.
log.info(
"=== Starting astrometric fitting...")
1230 associations.deprojectFittedStars()
1237 if self.
config.astrometryModel ==
"constrained":
1239 sky_to_tan_projection,
1240 chipOrder=self.
config.astrometryChipOrder,
1241 visitOrder=self.
config.astrometryVisitOrder)
1242 elif self.
config.astrometryModel ==
"simple":
1244 sky_to_tan_projection,
1247 order=self.
config.astrometrySimpleOrder)
1252 if self.
config.writeChi2FilesInitialFinal:
1253 baseName = f
"astrometry_initial_chi2-{dataName}"
1256 if self.
config.writeInitialModel:
1257 fullpath = self.
_getDebugPath(f
"initial_astrometry_model-{dataName}.txt")
1261 def getChi2Name(whatToFit):
1262 if self.
config.writeChi2FilesOuterLoop:
1263 return f
"astrometry_init-%s_chi2-{dataName}" % whatToFit
1267 if self.
config.writeInitMatrix:
1268 dumpMatrixFile = self.
_getDebugPath(f
"astrometry_preinit-{dataName}")
1273 if self.
config.astrometryModel ==
"constrained":
1274 fit.minimize(
"DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1276 writeChi2Name=getChi2Name(
"DistortionsVisit"))
1279 fit.minimize(
"Distortions", dumpMatrixFile=dumpMatrixFile)
1281 writeChi2Name=getChi2Name(
"Distortions"))
1283 fit.minimize(
"Positions")
1285 writeChi2Name=getChi2Name(
"Positions"))
1287 fit.minimize(
"Distortions Positions")
1289 writeChi2Name=getChi2Name(
"DistortionsPositions"))
1293 self.
config.maxAstrometrySteps,
1295 "Distortions Positions",
1296 sigmaRelativeTolerance=self.
config.astrometryOutlierRelativeTolerance,
1297 doRankUpdate=self.
config.astrometryDoRankUpdate,
1303 return Astrometry(fit, model, sky_to_tan_projection)
1306 """Count measured and reference stars per ccd and warn/log them."""
1307 for ccdImage
in associations.getCcdImageList():
1308 nMeasuredStars, nRefStars = ccdImage.countStars()
1309 self.
log.debug(
"ccdImage %s has %s measured and %s reference stars",
1310 ccdImage.getName(), nMeasuredStars, nRefStars)
1311 if nMeasuredStars < self.
config.minMeasuredStarsPerCcd:
1312 self.
log.warning(
"ccdImage %s has only %s measuredStars (desired %s)",
1313 ccdImage.getName(), nMeasuredStars, self.
config.minMeasuredStarsPerCcd)
1314 if nRefStars < self.
config.minRefStarsPerCcd:
1315 self.
log.warning(
"ccdImage %s has only %s RefStars (desired %s)",
1316 ccdImage.getName(), nRefStars, self.
config.minRefStarsPerCcd)
1320 sigmaRelativeTolerance=0,
1322 doLineSearch=False):
1323 """Run fitter.minimize up to max_steps times, returning the final chi2.
1327 associations : `lsst.jointcal.Associations`
1328 The star/reference star associations to fit.
1329 fitter : `lsst.jointcal.FitterBase`
1330 The fitter to use for minimization.
1332 Maximum number of steps to run outlier rejection before declaring
1333 convergence failure.
1334 name : {'photometry' or 'astrometry'}
1335 What type of data are we fitting (for logs and debugging files).
1337 Passed to ``fitter.minimize()`` to define the parameters to fit.
1338 dataName : `str`, optional
1339 Descriptive name for this dataset (e.g. tract and filter),
1341 sigmaRelativeTolerance : `float`, optional
1342 Convergence tolerance for the fractional change in the chi2 cut
1343 level for determining outliers. If set to zero, iterations will
1344 continue until there are no outliers.
1345 doRankUpdate : `bool`, optional
1346 Do an Eigen rank update during minimization, or recompute the full
1347 matrix and gradient?
1348 doLineSearch : `bool`, optional
1349 Do a line search for the optimum step during minimization?
1353 chi2: `lsst.jointcal.Chi2Statistic`
1354 The final chi2 after the fit converges, or is forced to end.
1359 Raised if the fitter fails with a non-finite value.
1361 Raised if the fitter fails for some other reason;
1362 log messages will provide further details.
1364 if self.
config.writeInitMatrix:
1365 dumpMatrixFile = self.
_getDebugPath(f
"{name}_postinit-{dataName}")
1369 oldChi2.chi2 = float(
"inf")
1370 for i
in range(max_steps):
1371 if self.
config.writeChi2FilesOuterLoop:
1372 writeChi2Name = f
"{name}_iterate_{i}_chi2-{dataName}"
1374 writeChi2Name =
None
1375 result = fitter.minimize(whatToFit,
1376 self.
config.outlierRejectSigma,
1377 sigmaRelativeTolerance=sigmaRelativeTolerance,
1378 doRankUpdate=doRankUpdate,
1379 doLineSearch=doLineSearch,
1380 dumpMatrixFile=dumpMatrixFile)
1383 f
"Fit iteration {i}", writeChi2Name=writeChi2Name)
1385 if result == MinimizeResult.Converged:
1387 self.
log.debug(
"fit has converged - no more outliers - redo minimization "
1388 "one more time in case we have lost accuracy in rank update.")
1390 result = fitter.minimize(whatToFit, self.
config.outlierRejectSigma,
1391 sigmaRelativeTolerance=sigmaRelativeTolerance)
1395 if chi2.chi2/chi2.ndof >= 4.0:
1396 self.
log.error(
"Potentially bad fit: High chi-squared/ndof.")
1399 elif result == MinimizeResult.Chi2Increased:
1400 self.
log.warning(
"Still some outliers remaining but chi2 increased - retry")
1402 chi2Ratio = chi2.chi2 / oldChi2.chi2
1404 self.
log.warning(
'Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1405 chi2.chi2, oldChi2.chi2, chi2Ratio)
1412 msg = (
"Large chi2 increase between steps: fit likely cannot converge."
1413 " Try setting one or more of the `writeChi2*` config fields and looking"
1414 " at how individual star chi2-values evolve during the fit.")
1415 raise RuntimeError(msg)
1417 elif result == MinimizeResult.NonFinite:
1418 filename = self.
_getDebugPath(
"{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1420 fitter.saveChi2Contributions(filename+
"{type}")
1421 msg =
"Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1422 raise FloatingPointError(msg.format(filename))
1423 elif result == MinimizeResult.Failed:
1424 raise RuntimeError(
"Chi2 minimization failure, cannot complete fit.")
1426 raise RuntimeError(
"Unxepected return code from minimize().")
1428 self.
log.error(
"%s failed to converge after %d steps"%(name, max_steps))
1433 """Return the internal jointcal models converted to the afw
1434 structures that will be saved to disk.
1438 ccdImageList : `lsst.jointcal.CcdImageList`
1439 The list of CcdImages to get the output for.
1440 model : `lsst.jointcal.AstrometryModel` or `lsst.jointcal.PhotometryModel`
1441 The internal jointcal model to convert for each `lsst.jointcal.CcdImage`.
1443 The name of the function to call on ``model`` to get the converted
1444 structure. Must accept an `lsst.jointcal.CcdImage`.
1448 output : `dict` [`tuple`, `lsst.jointcal.AstrometryModel`] or
1449 `dict` [`tuple`, `lsst.jointcal.PhotometryModel`]
1450 The data to be saved, keyed on (visit, detector).
1453 for ccdImage
in ccdImageList:
1454 ccd = ccdImage.ccdId
1455 visit = ccdImage.visit
1456 self.
log.debug(
"%s for visit: %d, ccd: %d", func, visit, ccd)
1457 output[(visit, ccd)] = getattr(model, func)(ccdImage)
1462 """Return an afw SourceTable to use as a base for creating the
1463 SourceCatalog to insert values from the dataFrame into.
1467 table : `lsst.afw.table.SourceTable`
1468 Table with schema and slots to use to make SourceCatalogs.
1471 schema.addField(
"centroid_x",
"D")
1472 schema.addField(
"centroid_y",
"D")
1473 schema.addField(
"centroid_xErr",
"F")
1474 schema.addField(
"centroid_yErr",
"F")
1475 schema.addField(
"shape_xx",
"D")
1476 schema.addField(
"shape_yy",
"D")
1477 schema.addField(
"shape_xy",
"D")
1478 schema.addField(
"flux_instFlux",
"D")
1479 schema.addField(
"flux_instFluxErr",
"D")
1481 table.defineCentroid(
"centroid")
1482 table.defineShape(
"shape")
1488 Get the sourceTable_visit columns to load from the catalogs.
1493 List of columns known to be available in the sourceTable_visit.
1494 config : `JointcalConfig`
1495 A filled-in config to to help define column names.
1496 sourceSelector : `lsst.meas.algorithms.BaseSourceSelectorTask`
1497 A configured source selector to define column names to load.
1502 List of columns to read from sourceTable_visit.
1504 Name of the ixx/iyy/ixy columns.
1506 columns = [
'visit',
'detector',
1507 'sourceId',
'x',
'xErr',
'y',
'yErr',
1508 config.sourceFluxType +
'_instFlux', config.sourceFluxType +
'_instFluxErr']
1510 if 'ixx' in inColumns:
1512 ixxColumns = [
'ixx',
'iyy',
'ixy']
1515 ixxColumns = [
'Ixx',
'Iyy',
'Ixy']
1516 columns.extend(ixxColumns)
1518 if sourceSelector.config.doFlags:
1519 columns.extend(sourceSelector.config.flags.bad)
1520 if sourceSelector.config.doUnresolved:
1521 columns.append(sourceSelector.config.unresolved.name)
1522 if sourceSelector.config.doIsolated:
1523 columns.append(sourceSelector.config.isolated.parentName)
1524 columns.append(sourceSelector.config.isolated.nChildName)
1525 if sourceSelector.config.doRequireFiniteRaDec:
1526 columns.append(sourceSelector.config.requireFiniteRaDec.raColName)
1527 columns.append(sourceSelector.config.requireFiniteRaDec.decColName)
1528 if sourceSelector.config.doRequirePrimary:
1529 columns.append(sourceSelector.config.requirePrimary.primaryColName)
1531 return columns, ixxColumns
1535 ixxColumns, sourceFluxType, log):
1536 """Return an afw SourceCatalog extracted from a visit-level dataframe,
1537 limited to just one detector.
1541 table : `lsst.afw.table.SourceTable`
1542 Table factory to use to make the SourceCatalog that will be
1543 populated with data from ``visitCatalog``.
1544 visitCatalog : `pandas.DataFrame`
1545 DataFrame to extract a detector catalog from.
1547 Numeric id of the detector to extract from ``visitCatalog``.
1548 ixxColumns : `list` [`str`]
1549 Names of the ixx/iyy/ixy columns in the catalog.
1550 sourceFluxType : `str`
1551 Name of the catalog field to load instFluxes from.
1552 log : `logging.Logger`
1553 Logging instance to log to.
1557 catalog : `lsst.afw.table.SourceCatalog`, or `None`
1558 Detector-level catalog extracted from ``visitCatalog``, or `None`
1559 if there was no data to load.
1562 mapping = {
'x':
'centroid_x',
1564 'xErr':
'centroid_xErr',
1565 'yErr':
'centroid_yErr',
1566 ixxColumns[0]:
'shape_xx',
1567 ixxColumns[1]:
'shape_yy',
1568 ixxColumns[2]:
'shape_xy',
1569 f
'{sourceFluxType}_instFlux':
'flux_instFlux',
1570 f
'{sourceFluxType}_instFluxErr':
'flux_instFluxErr',
1574 matched = visitCatalog[
'detector'] == detectorId
1578 catalog.resize(sum(matched))
1579 view = visitCatalog.loc[matched]
1580 catalog[
'id'] = view.index
1581 for dfCol, afwCol
in mapping.items():
1582 catalog[afwCol] = view[dfCol]
1584 log.debug(
"%d sources selected in visit %d detector %d",
1586 view[
'visit'].iloc[0],
A representation of a detector in a mosaic camera.
A group of labels for a filter in an exposure or coadd.
The photometric calibration of an exposure.
Information about a single exposure of an imaging camera.
Custom catalog class for ExposureRecord/Table.
static Schema makeMinimalSchema()
Return a minimal schema for Exposure tables and records.
static std::shared_ptr< SourceTable > make(Schema const &schema, std::shared_ptr< IdFactory > const &idFactory)
Construct a new table.
static Schema makeMinimalSchema()
Return a minimal schema for Source tables and records.
Class for storing ordered metadata with comments.
A class representing an angle.
An integer coordinate rectangle.
Point in an unspecified spherical coordinate system.
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.
__init__(self, *config=None)
getSpatialBoundsConnections(self)
_build_ccdImage(self, data, associations, jointcalControl)
_iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", sigmaRelativeTolerance=0, doRankUpdate=True, doLineSearch=False)
_put_metrics(self, butlerQC, job, outputRefs)
_getDebugPath(self, filename)
_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)
_check_stars(self, associations)
_prep_sky(self, associations, filters)
_compute_proper_motion_epoch(self, ccdImageList)
_make_one_input_data(self, visitRecord, catalog, detectorDict)
_check_star_lists(self, associations, name)
_fit_photometry(self, associations, dataName=None)
_put_output(self, butlerQC, outputs, outputRefs, camera, setter)
_load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel, applyColorterms=False, epoch=None)
_logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None)
_get_refcat_coordinate_error_override(self, refCat, name)
_load_data(self, inputSourceTableVisit, inputVisitSummary, associations, jointcalControl, camera)
_make_output(self, ccdImageList, model, func)
runQuantum(self, butlerQC, inputRefs, outputRefs)
_fit_astrometry(self, associations, dataName=None)
run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None)
extract_detector_catalog_from_visit_catalog(table, visitCatalog, detectorId, ixxColumns, sourceFluxType, log)
add_measurement(job, name, value)
writeModel(model, filename, log)
get_sourceTable_visit_columns(inColumns, config, sourceSelector)