LSST Applications 26.0.0,g0265f82a02+6660c170cc,g07994bdeae+30b05a742e,g0a0026dc87+17526d298f,g0a60f58ba1+17526d298f,g0e4bf8285c+96dd2c2ea9,g0ecae5effc+c266a536c8,g1e7d6db67d+6f7cb1f4bb,g26482f50c6+6346c0633c,g2bbee38e9b+6660c170cc,g2cc88a2952+0a4e78cd49,g3273194fdb+f6908454ef,g337abbeb29+6660c170cc,g337c41fc51+9a8f8f0815,g37c6e7c3d5+7bbafe9d37,g44018dc512+6660c170cc,g4a941329ef+4f7594a38e,g4c90b7bd52+5145c320d2,g58be5f913a+bea990ba40,g635b316a6c+8d6b3a3e56,g67924a670a+bfead8c487,g6ae5381d9b+81bc2a20b4,g93c4d6e787+26b17396bd,g98cecbdb62+ed2cb6d659,g98ffbb4407+81bc2a20b4,g9ddcbc5298+7f7571301f,ga1e77700b3+99e9273977,gae46bcf261+6660c170cc,gb2715bf1a1+17526d298f,gc86a011abf+17526d298f,gcf0d15dbbd+96dd2c2ea9,gdaeeff99f8+0d8dbea60f,gdb4ec4c597+6660c170cc,ge23793e450+96dd2c2ea9,gf041782ebf+171108ac67
LSST Data Management Base Package
Loading...
Searching...
No Matches
jointcal.py
Go to the documentation of this file.
1# This file is part of jointcal.
2#
3# Developed for the LSST Data Management System.
4# This product includes software developed by the LSST Project
5# (https://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <https://www.gnu.org/licenses/>.
21
22import dataclasses
23import collections
24import os
25import logging
26
27import astropy.time
28import numpy as np
29import astropy.units as u
30
31import lsst.geom
32import lsst.utils
33import lsst.pex.config as pexConfig
34import lsst.pipe.base as pipeBase
35from lsst.afw.image import fluxErrFromABMagErr
37import lsst.afw.table
38from lsst.pipe.base import Instrument
39from lsst.pipe.tasks.colorterms import ColortermLibrary
40from lsst.verify import Job, Measurement
41
42from lsst.meas.algorithms import (ReferenceObjectLoader, ReferenceSourceSelectorTask,
43 LoadReferenceObjectsConfig)
44from lsst.meas.algorithms.sourceSelector import sourceSelectorRegistry
45
46import lsst.jointcal
47from lsst.jointcal import MinimizeResult
48
49__all__ = ["JointcalConfig", "JointcalTask"]
50
51Photometry = collections.namedtuple('Photometry', ('fit', 'model'))
52Astrometry = collections.namedtuple('Astrometry', ('fit', 'model', 'sky_to_tan_projection'))
53
54
55# TODO: move this to MeasurementSet in lsst.verify per DM-12655.
56def add_measurement(job, name, value):
57 meas = Measurement(job.metrics[name], value)
58 job.measurements.insert(meas)
59
60
61def lookupStaticCalibrations(datasetType, registry, quantumDataId, collections):
62 """Lookup function that asserts/hopes that a static calibration dataset
63 exists in a particular collection, since this task can't provide a single
64 date/time to use to search for one properly.
65
66 This is mostly useful for the ``camera`` dataset, in cases where the task's
67 quantum dimensions do *not* include something temporal, like ``exposure``
68 or ``visit``.
69
70 Parameters
71 ----------
72 datasetType : `lsst.daf.butler.DatasetType`
73 Type of dataset being searched for.
74 registry : `lsst.daf.butler.Registry`
75 Data repository registry to search.
76 quantumDataId : `lsst.daf.butler.DataCoordinate`
77 Data ID of the quantum this camera should match.
78 collections : `Iterable` [ `str` ]
79 Collections that should be searched - but this lookup function works
80 by ignoring this in favor of a more-or-less hard-coded value.
81
82 Returns
83 -------
84 refs : `Iterator` [ `lsst.daf.butler.DatasetRef` ]
85 Iterator over dataset references; should have only one element.
86
87 Notes
88 -----
89 This implementation duplicates one in fgcmcal, and is at least quite
90 similar to another in cp_pipe. This duplicate has the most documentation.
91 Fixing this is DM-29661.
92 """
93 instrument = Instrument.fromName(quantumDataId["instrument"], registry)
94 unboundedCollection = instrument.makeUnboundedCalibrationRunName()
95 return registry.queryDatasets(datasetType,
96 dataId=quantumDataId,
97 collections=[unboundedCollection],
98 findFirst=True)
99
100
101def lookupVisitRefCats(datasetType, registry, quantumDataId, collections):
102 """Lookup function that finds all refcats for all visits that overlap a
103 tract, rather than just the refcats that directly overlap the tract.
104
105 Parameters
106 ----------
107 datasetType : `lsst.daf.butler.DatasetType`
108 Type of dataset being searched for.
109 registry : `lsst.daf.butler.Registry`
110 Data repository registry to search.
111 quantumDataId : `lsst.daf.butler.DataCoordinate`
112 Data ID of the quantum; expected to be something we can use as a
113 constraint to query for overlapping visits.
114 collections : `Iterable` [ `str` ]
115 Collections to search.
116
117 Returns
118 -------
119 refs : `Iterator` [ `lsst.daf.butler.DatasetRef` ]
120 Iterator over refcat references.
121 """
122 refs = set()
123 # Use .expanded() on the query methods below because we need data IDs with
124 # regions, both in the outer loop over visits (queryDatasets will expand
125 # any data ID we give it, but doing it up-front in bulk is much more
126 # efficient) and in the data IDs of the DatasetRefs this function yields
127 # (because the RefCatLoader relies on them to do some of its own
128 # filtering).
129 for visit_data_id in set(registry.queryDataIds("visit", dataId=quantumDataId).expanded()):
130 refs.update(
131 registry.queryDatasets(
132 datasetType,
133 collections=collections,
134 dataId=visit_data_id,
135 findFirst=True,
136 ).expanded()
137 )
138 yield from refs
139
140
141class JointcalTaskConnections(pipeBase.PipelineTaskConnections,
142 dimensions=("skymap", "tract", "instrument", "physical_filter")):
143 """Middleware input/output connections for jointcal data."""
144 inputCamera = pipeBase.connectionTypes.PrerequisiteInput(
145 doc="The camera instrument that took these observations.",
146 name="camera",
147 storageClass="Camera",
148 dimensions=("instrument",),
149 isCalibration=True,
150 lookupFunction=lookupStaticCalibrations,
151 )
152 inputSourceTableVisit = pipeBase.connectionTypes.Input(
153 doc="Source table in parquet format, per visit",
154 name="sourceTable_visit",
155 storageClass="DataFrame",
156 dimensions=("instrument", "visit"),
157 deferLoad=True,
158 multiple=True,
159 )
160 inputVisitSummary = pipeBase.connectionTypes.Input(
161 doc=("Per-visit consolidated exposure metadata built from calexps. "
162 "These catalogs use detector id for the id and must be sorted for "
163 "fast lookups of a detector."),
164 name="visitSummary",
165 storageClass="ExposureCatalog",
166 dimensions=("instrument", "visit"),
167 deferLoad=True,
168 multiple=True,
169 )
170 astrometryRefCat = pipeBase.connectionTypes.PrerequisiteInput(
171 doc="The astrometry reference catalog to match to loaded input catalog sources.",
172 name="gaia_dr2_20200414",
173 storageClass="SimpleCatalog",
174 dimensions=("skypix",),
175 deferLoad=True,
176 multiple=True,
177 lookupFunction=lookupVisitRefCats,
178 )
179 photometryRefCat = pipeBase.connectionTypes.PrerequisiteInput(
180 doc="The photometry reference catalog to match to loaded input catalog sources.",
181 name="ps1_pv3_3pi_20170110",
182 storageClass="SimpleCatalog",
183 dimensions=("skypix",),
184 deferLoad=True,
185 multiple=True,
186 lookupFunction=lookupVisitRefCats,
187 )
188
189 outputWcs = pipeBase.connectionTypes.Output(
190 doc=("Per-tract, per-visit world coordinate systems derived from the fitted model."
191 " These catalogs only contain entries for detectors with an output, and use"
192 " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
193 name="jointcalSkyWcsCatalog",
194 storageClass="ExposureCatalog",
195 dimensions=("instrument", "visit", "skymap", "tract"),
196 multiple=True
197 )
198 outputPhotoCalib = pipeBase.connectionTypes.Output(
199 doc=("Per-tract, per-visit photometric calibrations derived from the fitted model."
200 " These catalogs only contain entries for detectors with an output, and use"
201 " the detector id for the catalog id, sorted on id for fast lookups of a detector."),
202 name="jointcalPhotoCalibCatalog",
203 storageClass="ExposureCatalog",
204 dimensions=("instrument", "visit", "skymap", "tract"),
205 multiple=True
206 )
207
208 # measurements of metrics
209 # The vars() trick used here allows us to set class attributes
210 # programatically. Taken from:
211 # https://stackoverflow.com/questions/2519807/setting-a-class-attribute-with-a-given-name-in-python-while-defining-the-class
212 for name in ("astrometry", "photometry"):
213 vars()[f"{name}_matched_fittedStars"] = pipeBase.connectionTypes.Output(
214 doc=f"The number of cross-matched fittedStars for {name}",
215 name=f"metricvalue_jointcal_{name}_matched_fittedStars",
216 storageClass="MetricValue",
217 dimensions=("skymap", "tract", "instrument", "physical_filter"),
218 )
219 vars()[f"{name}_collected_refStars"] = pipeBase.connectionTypes.Output(
220 doc=f"The number of {name} reference stars drawn from the reference catalog, before matching.",
221 name=f"metricvalue_jointcal_{name}_collected_refStars",
222 storageClass="MetricValue",
223 dimensions=("skymap", "tract", "instrument", "physical_filter"),
224 )
225 vars()[f"{name}_prepared_refStars"] = pipeBase.connectionTypes.Output(
226 doc=f"The number of {name} reference stars matched to fittedStars.",
227 name=f"metricvalue_jointcal_{name}_prepared_refStars",
228 storageClass="MetricValue",
229 dimensions=("skymap", "tract", "instrument", "physical_filter"),
230 )
231 vars()[f"{name}_prepared_fittedStars"] = pipeBase.connectionTypes.Output(
232 doc=f"The number of cross-matched fittedStars after cleanup, for {name}.",
233 name=f"metricvalue_jointcal_{name}_prepared_fittedStars",
234 storageClass="MetricValue",
235 dimensions=("skymap", "tract", "instrument", "physical_filter"),
236 )
237 vars()[f"{name}_prepared_ccdImages"] = pipeBase.connectionTypes.Output(
238 doc=f"The number of ccdImages that will be fit for {name}, after cleanup.",
239 name=f"metricvalue_jointcal_{name}_prepared_ccdImages",
240 storageClass="MetricValue",
241 dimensions=("skymap", "tract", "instrument", "physical_filter"),
242 )
243 vars()[f"{name}_final_chi2"] = pipeBase.connectionTypes.Output(
244 doc=f"The final chi2 of the {name} fit.",
245 name=f"metricvalue_jointcal_{name}_final_chi2",
246 storageClass="MetricValue",
247 dimensions=("skymap", "tract", "instrument", "physical_filter"),
248 )
249 vars()[f"{name}_final_ndof"] = pipeBase.connectionTypes.Output(
250 doc=f"The number of degrees of freedom of the fitted {name}.",
251 name=f"metricvalue_jointcal_{name}_final_ndof",
252 storageClass="MetricValue",
253 dimensions=("skymap", "tract", "instrument", "physical_filter"),
254 )
255
256 def __init__(self, *, config=None):
257 super().__init__(config=config)
258 # When we are only doing one of astrometry or photometry, we don't
259 # need the reference catalog or produce the outputs for the other.
260 # This informs the middleware of that when the QuantumGraph is
261 # generated, so we don't block on getting something we won't need or
262 # create an expectation that downstream tasks will be able to consume
263 # something we won't produce.
264 if not config.doAstrometry:
265 self.prerequisiteInputs.remove("astrometryRefCat")
266 self.outputs.remove("outputWcs")
267 for key in list(self.outputs):
268 if "metricvalue_jointcal_astrometry" in key:
269 self.outputs.remove(key)
270 if not config.doPhotometry:
271 self.prerequisiteInputs.remove("photometryRefCat")
272 self.outputs.remove("outputPhotoCalib")
273 for key in list(self.outputs):
274 if "metricvalue_jointcal_photometry" in key:
275 self.outputs.remove(key)
276
278 return ("inputVisitSummary",)
279
280
281class JointcalConfig(pipeBase.PipelineTaskConfig,
282 pipelineConnections=JointcalTaskConnections):
283 """Configuration for JointcalTask"""
284
285 doAstrometry = pexConfig.Field(
286 doc="Fit astrometry and write the fitted result.",
287 dtype=bool,
288 default=True
289 )
290 doPhotometry = pexConfig.Field(
291 doc="Fit photometry and write the fitted result.",
292 dtype=bool,
293 default=True
294 )
295 sourceFluxType = pexConfig.Field(
296 dtype=str,
297 doc="Source flux field to use in source selection and to get fluxes from the catalog.",
298 default='apFlux_12_0'
299 )
300 positionErrorPedestal = pexConfig.Field(
301 doc="Systematic term to apply to the measured position error (pixels)",
302 dtype=float,
303 default=0.02,
304 )
305 photometryErrorPedestal = pexConfig.Field(
306 doc="Systematic term to apply to the measured error on flux or magnitude as a "
307 "fraction of source flux or magnitude delta (e.g. 0.05 is 5% of flux or +50 millimag).",
308 dtype=float,
309 default=0.0,
310 )
311 # TODO: DM-6885 matchCut should be an geom.Angle
312 matchCut = pexConfig.Field(
313 doc="Matching radius between fitted and reference stars (arcseconds)",
314 dtype=float,
315 default=3.0,
316 )
317 minMeasurements = pexConfig.Field(
318 doc="Minimum number of associated measured stars for a fitted star to be included in the fit",
319 dtype=int,
320 default=2,
321 )
322 minMeasuredStarsPerCcd = pexConfig.Field(
323 doc="Minimum number of measuredStars per ccdImage before printing warnings",
324 dtype=int,
325 default=100,
326 )
327 minRefStarsPerCcd = pexConfig.Field(
328 doc="Minimum number of measuredStars per ccdImage before printing warnings",
329 dtype=int,
330 default=30,
331 )
332 allowLineSearch = pexConfig.Field(
333 doc="Allow a line search during minimization, if it is reasonable for the model"
334 " (models with a significant non-linear component, e.g. constrainedPhotometry).",
335 dtype=bool,
336 default=False
337 )
338 astrometrySimpleOrder = pexConfig.Field(
339 doc="Polynomial order for fitting the simple astrometry model.",
340 dtype=int,
341 default=3,
342 )
343 astrometryChipOrder = pexConfig.Field(
344 doc="Order of the per-chip transform for the constrained astrometry model.",
345 dtype=int,
346 default=1,
347 )
348 astrometryVisitOrder = pexConfig.Field(
349 doc="Order of the per-visit transform for the constrained astrometry model.",
350 dtype=int,
351 default=5,
352 )
353 useInputWcs = pexConfig.Field(
354 doc="Use the input calexp WCSs to initialize a SimpleAstrometryModel.",
355 dtype=bool,
356 default=True,
357 )
358 astrometryModel = pexConfig.ChoiceField(
359 doc="Type of model to fit to astrometry",
360 dtype=str,
361 default="constrained",
362 allowed={"simple": "One polynomial per ccd",
363 "constrained": "One polynomial per ccd, and one polynomial per visit"}
364 )
365 photometryModel = pexConfig.ChoiceField(
366 doc="Type of model to fit to photometry",
367 dtype=str,
368 default="constrainedMagnitude",
369 allowed={"simpleFlux": "One constant zeropoint per ccd and visit, fitting in flux space.",
370 "constrainedFlux": "Constrained zeropoint per ccd, and one polynomial per visit,"
371 " fitting in flux space.",
372 "simpleMagnitude": "One constant zeropoint per ccd and visit,"
373 " fitting in magnitude space.",
374 "constrainedMagnitude": "Constrained zeropoint per ccd, and one polynomial per visit,"
375 " fitting in magnitude space.",
376 }
377 )
378 applyColorTerms = pexConfig.Field(
379 doc="Apply photometric color terms to reference stars?"
380 "Requires that colorterms be set to a ColortermLibrary",
381 dtype=bool,
382 default=False
383 )
384 colorterms = pexConfig.ConfigField(
385 doc="Library of photometric reference catalog name to color term dict.",
386 dtype=ColortermLibrary,
387 )
388 photometryVisitOrder = pexConfig.Field(
389 doc="Order of the per-visit polynomial transform for the constrained photometry model.",
390 dtype=int,
391 default=7,
392 )
393 photometryDoRankUpdate = pexConfig.Field(
394 doc=("Do the rank update step during minimization. "
395 "Skipping this can help deal with models that are too non-linear."),
396 dtype=bool,
397 default=True,
398 )
399 astrometryDoRankUpdate = pexConfig.Field(
400 doc=("Do the rank update step during minimization (should not change the astrometry fit). "
401 "Skipping this can help deal with models that are too non-linear."),
402 dtype=bool,
403 default=True,
404 )
405 outlierRejectSigma = pexConfig.Field(
406 doc="How many sigma to reject outliers at during minimization.",
407 dtype=float,
408 default=5.0,
409 )
410 astrometryOutlierRelativeTolerance = pexConfig.Field(
411 doc=("Convergence tolerance for outlier rejection threshold when fitting astrometry. Iterations will "
412 "stop when the fractional change in the chi2 cut level is below this value. If tolerance is set "
413 "to zero, iterations will continue until there are no more outliers. We suggest a value of 0.002"
414 "as a balance between a shorter minimization runtime and achieving a final fitted model that is"
415 "close to the solution found when removing all outliers."),
416 dtype=float,
417 default=0,
418 )
419 maxPhotometrySteps = pexConfig.Field(
420 doc="Maximum number of minimize iterations to take when fitting photometry.",
421 dtype=int,
422 default=20,
423 )
424 maxAstrometrySteps = pexConfig.Field(
425 doc="Maximum number of minimize iterations to take when fitting astrometry.",
426 dtype=int,
427 default=20,
428 )
429 astrometryRefObjLoader = pexConfig.ConfigField(
430 dtype=LoadReferenceObjectsConfig,
431 doc="Reference object loader for astrometric fit",
432 )
433 photometryRefObjLoader = pexConfig.ConfigField(
434 dtype=LoadReferenceObjectsConfig,
435 doc="Reference object loader for photometric fit",
436 )
437 sourceSelector = sourceSelectorRegistry.makeField(
438 doc="How to select sources for cross-matching",
439 default="science"
440 )
441 astrometryReferenceSelector = pexConfig.ConfigurableField(
442 target=ReferenceSourceSelectorTask,
443 doc="How to down-select the loaded astrometry reference catalog.",
444 )
445 photometryReferenceSelector = pexConfig.ConfigurableField(
446 target=ReferenceSourceSelectorTask,
447 doc="How to down-select the loaded photometry reference catalog.",
448 )
449 astrometryReferenceErr = pexConfig.Field(
450 doc=("Uncertainty on reference catalog coordinates [mas] to use in place of the `coord_*Err` fields. "
451 "If None, then raise an exception if the reference catalog is missing coordinate errors. "
452 "If specified, overrides any existing `coord_*Err` values."),
453 dtype=float,
454 default=None,
455 optional=True
456 )
457
458 # configs for outputting debug information
459 writeInitMatrix = pexConfig.Field(
460 dtype=bool,
461 doc=("Write the pre/post-initialization Hessian and gradient to text files, for debugging. "
462 "Output files will be written to `config.debugOutputPath` and will "
463 "be of the form 'astrometry_[pre|post]init-TRACT-FILTER-mat.txt'. "
464 "Note that these files are the dense versions of the matrix, and so may be very large."),
465 default=False
466 )
467 writeChi2FilesInitialFinal = pexConfig.Field(
468 dtype=bool,
469 doc=("Write .csv files containing the contributions to chi2 for the initialization and final fit. "
470 "Output files will be written to `config.debugOutputPath` and will "
471 "be of the form `astrometry_[initial|final]_chi2-TRACT-FILTER."),
472 default=False
473 )
474 writeChi2FilesOuterLoop = pexConfig.Field(
475 dtype=bool,
476 doc=("Write .csv files containing the contributions to chi2 for the outer fit loop. "
477 "Output files will be written to `config.debugOutputPath` and will "
478 "be of the form `astrometry_init-NN_chi2-TRACT-FILTER`."),
479 default=False
480 )
481 writeInitialModel = pexConfig.Field(
482 dtype=bool,
483 doc=("Write the pre-initialization model to text files, for debugging. "
484 "Output files will be written to `config.debugOutputPath` and will be "
485 "of the form `initial_astrometry_model-TRACT_FILTER.txt`."),
486 default=False
487 )
488 debugOutputPath = pexConfig.Field(
489 dtype=str,
490 default=".",
491 doc=("Path to write debug output files to. Used by "
492 "`writeInitialModel`, `writeChi2Files*`, `writeInitMatrix`.")
493 )
494 detailedProfile = pexConfig.Field(
495 dtype=bool,
496 default=False,
497 doc="Output separate profiling information for different parts of jointcal, e.g. data read, fitting"
498 )
499
500 def validate(self):
501 super().validate()
502 if self.doPhotometry and self.applyColorTerms and len(self.colorterms.data) == 0:
503 msg = "applyColorTerms=True requires the `colorterms` field be set to a ColortermLibrary."
504 raise pexConfig.FieldValidationError(JointcalConfig.colorterms, self, msg)
505 if self.doAstrometry and not self.doPhotometry and self.applyColorTerms:
506 msg = ("Only doing astrometry, but Colorterms are not applied for astrometry;"
507 "applyColorTerms=True will be ignored.")
508 logging.getLogger("lsst.jointcal").warning(msg)
509
510 def setDefaults(self):
511 # Use only primary stars.
512 self.sourceSelector["science"].doRequirePrimary = True
513 # Use only stars because aperture fluxes of galaxies are biased and depend on seeing.
514 self.sourceSelector["science"].doUnresolved = True
515 self.sourceSelector["science"].unresolved.name = "extendedness"
516 # with dependable signal to noise ratio.
517 self.sourceSelector["science"].doSignalToNoise = True
518 # Min SNR must be > 0 because jointcal cannot handle negative fluxes,
519 # and S/N > 10 to use sources that are not too faint, and thus better measured.
520 self.sourceSelector["science"].signalToNoise.minimum = 10.
521 # Base SNR selection on `sourceFluxType` because that is the flux that jointcal fits.
522 self.sourceSelector["science"].signalToNoise.fluxField = f"{self.sourceFluxType}_instFlux"
523 self.sourceSelector["science"].signalToNoise.errField = f"{self.sourceFluxType}_instFluxErr"
524 # Do not trust blended sources" aperture fluxes which also depend on seeing.
525 self.sourceSelector["science"].doIsolated = True
526 self.sourceSelector["science"].isolated.parentName = "parentSourceId"
527 self.sourceSelector["science"].isolated.nChildName = "deblend_nChild"
528 # Do not trust either flux or centroid measurements with flags,
529 # chosen from the usual QA flags for stars)
530 self.sourceSelector["science"].doFlags = True
531 badFlags = ["pixelFlags_edge",
532 "pixelFlags_saturated",
533 "pixelFlags_interpolatedCenter",
534 "pixelFlags_interpolated",
535 "pixelFlags_crCenter",
536 "pixelFlags_bad",
537 "hsmPsfMoments_flag",
538 f"{self.sourceFluxType}_flag",
539 ]
540 self.sourceSelector["science"].flags.bad = badFlags
541 self.sourceSelector["science"].doRequireFiniteRaDec = True
542 self.sourceSelector["science"].requireFiniteRaDec.raColName = "ra"
543 self.sourceSelector["science"].requireFiniteRaDec.decColName = "dec"
544
545 # Use Gaia-DR2 with proper motions for astrometry; phot_g_mean is the
546 # primary Gaia band, but is not like any normal photometric band.
547 self.astrometryRefObjLoader.requireProperMotion = True
548 self.astrometryRefObjLoader.anyFilterMapsToThis = "phot_g_mean"
549
550
551def writeModel(model, filename, log):
552 """Write model to outfile."""
553 with open(filename, "w") as file:
554 file.write(repr(model))
555 log.info("Wrote %s to file: %s", model, filename)
556
557
558@dataclasses.dataclass
560 """The input data jointcal needs for each detector/visit."""
561 visit: int
562 """The visit identifier of this exposure."""
564 """The catalog derived from this exposure."""
566 """The VisitInfo of this exposure."""
568 """The detector of this exposure."""
570 """The photometric calibration of this exposure."""
572 """The WCS of this exposure."""
574 """The bounding box of this exposure."""
576 """The filter of this exposure."""
577
578
579class JointcalTask(pipeBase.PipelineTask):
580 """Astrometricly and photometricly calibrate across multiple visits of the
581 same field.
582 """
583
584 ConfigClass = JointcalConfig
585 _DefaultName = "jointcal"
586
587 def __init__(self, **kwargs):
588 super().__init__(**kwargs)
589 self.makeSubtask("sourceSelector")
590 if self.config.doAstrometry:
591 self.makeSubtask("astrometryReferenceSelector")
592 else:
594 if self.config.doPhotometry:
595 self.makeSubtask("photometryReferenceSelector")
596 else:
598
599 # To hold various computed metrics for use by tests
600 self.job = Job.load_metrics_package(subset='jointcal')
601
602 def runQuantum(self, butlerQC, inputRefs, outputRefs):
603 # We override runQuantum to set up the refObjLoaders and write the
604 # outputs to the correct refs.
605 inputs = butlerQC.get(inputRefs)
606 # We want the tract number for writing debug files
607 tract = butlerQC.quantum.dataId['tract']
608 if self.config.doAstrometry:
610 dataIds=[ref.datasetRef.dataId
611 for ref in inputRefs.astrometryRefCat],
612 refCats=inputs.pop('astrometryRefCat'),
613 config=self.config.astrometryRefObjLoader,
614 name=self.config.connections.astrometryRefCat,
615 log=self.log)
616 if self.config.doPhotometry:
618 dataIds=[ref.datasetRef.dataId
619 for ref in inputRefs.photometryRefCat],
620 refCats=inputs.pop('photometryRefCat'),
621 config=self.config.photometryRefObjLoader,
622 name=self.config.connections.photometryRefCat,
623 log=self.log)
624 outputs = self.run(**inputs, tract=tract)
625 self._put_metrics(butlerQC, outputs.job, outputRefs)
626 if self.config.doAstrometry:
627 self._put_output(butlerQC, outputs.outputWcs, outputRefs.outputWcs,
628 inputs['inputCamera'], "setWcs")
629 if self.config.doPhotometry:
630 self._put_output(butlerQC, outputs.outputPhotoCalib, outputRefs.outputPhotoCalib,
631 inputs['inputCamera'], "setPhotoCalib")
632
633 def _put_metrics(self, butlerQC, job, outputRefs):
634 """Persist all measured metrics stored in a job.
635
636 Parameters
637 ----------
638 butlerQC : `lsst.pipe.base.ButlerQuantumContext`
639 A butler which is specialized to operate in the context of a
640 `lsst.daf.butler.Quantum`; This is the input to `runQuantum`.
641 job : `lsst.verify.job.Job`
642 Measurements of metrics to persist.
643 outputRefs : `list` [`lsst.pipe.base.connectionTypes.OutputQuantizedConnection`]
644 The DatasetRefs to persist the data to.
645 """
646 for key in job.measurements.keys():
647 butlerQC.put(job.measurements[key], getattr(outputRefs, key.fqn.replace('jointcal.', '')))
648
649 def _put_output(self, butlerQC, outputs, outputRefs, camera, setter):
650 """Persist the output datasets to their appropriate datarefs.
651
652 Parameters
653 ----------
654 butlerQC : `lsst.pipe.base.ButlerQuantumContext`
655 A butler which is specialized to operate in the context of a
656 `lsst.daf.butler.Quantum`; This is the input to `runQuantum`.
657 outputs : `dict` [`tuple`, `lsst.afw.geom.SkyWcs`] or
658 `dict` [`tuple, `lsst.afw.image.PhotoCalib`]
659 The fitted objects to persist.
660 outputRefs : `list` [`lsst.pipe.base.connectionTypes.OutputQuantizedConnection`]
661 The DatasetRefs to persist the data to.
663 The camera for this instrument, to get detector ids from.
664 setter : `str`
665 The method to call on the ExposureCatalog to set each output.
666 """
668 schema.addField('visit', type='L', doc='Visit number')
669
670 def new_catalog(visit, size):
671 """Return an catalog ready to be filled with appropriate output."""
672 catalog = lsst.afw.table.ExposureCatalog(schema)
673 catalog.resize(size)
674 catalog['visit'] = visit
675 metadata = lsst.daf.base.PropertyList()
676 metadata.add("COMMENT", "Catalog id is detector id, sorted.")
677 metadata.add("COMMENT", "Only detectors with data have entries.")
678 return catalog
679
680 # count how many detectors have output for each visit
681 detectors_per_visit = collections.defaultdict(int)
682 for key in outputs:
683 # key is (visit, detector_id), and we only need visit here
684 detectors_per_visit[key[0]] += 1
685
686 for ref in outputRefs:
687 visit = ref.dataId['visit']
688 catalog = new_catalog(visit, detectors_per_visit[visit])
689 # Iterate over every detector and skip the ones we don't have output for.
690 i = 0
691 for detector in camera:
692 detectorId = detector.getId()
693 key = (ref.dataId['visit'], detectorId)
694 if key not in outputs:
695 # skip detectors we don't have output for
696 self.log.debug("No %s output for detector %s in visit %s",
697 setter[3:], detectorId, visit)
698 continue
699
700 catalog[i].setId(detectorId)
701 getattr(catalog[i], setter)(outputs[key])
702 i += 1
703
704 catalog.sort() # ensure that the detectors are in sorted order, for fast lookups
705 butlerQC.put(catalog, ref)
706 self.log.info("Wrote %s detectors to %s", i, ref)
707
708 def run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None):
709 # Docstring inherited.
710
711 # We take values out of the Parquet table, and put them in "flux_",
712 # and the config.sourceFluxType field is used during that extraction,
713 # so just use "flux" here.
714 sourceFluxField = "flux"
715 jointcalControl = lsst.jointcal.JointcalControl(sourceFluxField)
716 associations = lsst.jointcal.Associations()
717 self.focalPlaneBBox = inputCamera.getFpBBox()
718 oldWcsList, bands = self._load_data(inputSourceTableVisit,
719 inputVisitSummary,
720 associations,
721 jointcalControl,
722 inputCamera)
723
724 boundingCircle, center, radius, defaultFilter, epoch = self._prep_sky(associations, bands)
725
726 if self.config.doAstrometry:
727 astrometry = self._do_load_refcat_and_fit(associations, defaultFilter, center, radius,
728 name="astrometry",
729 refObjLoader=self.astrometryRefObjLoader,
730 referenceSelector=self.astrometryReferenceSelector,
731 fit_function=self._fit_astrometry,
732 tract=tract,
733 epoch=epoch)
734 astrometry_output = self._make_output(associations.getCcdImageList(),
735 astrometry.model,
736 "makeSkyWcs")
737 else:
738 astrometry_output = None
739
740 if self.config.doPhotometry:
741 photometry = self._do_load_refcat_and_fit(associations, defaultFilter, center, radius,
742 name="photometry",
743 refObjLoader=self.photometryRefObjLoader,
744 referenceSelector=self.photometryReferenceSelector,
745 fit_function=self._fit_photometry,
746 tract=tract,
747 epoch=epoch,
748 reject_bad_fluxes=True)
749 photometry_output = self._make_output(associations.getCcdImageList(),
750 photometry.model,
751 "toPhotoCalib")
752 else:
753 photometry_output = None
754
755 return pipeBase.Struct(outputWcs=astrometry_output,
756 outputPhotoCalib=photometry_output,
757 job=self.job,
758 astrometryRefObjLoader=self.astrometryRefObjLoader,
759 photometryRefObjLoader=self.photometryRefObjLoader)
760
761 def _load_data(self, inputSourceTableVisit, inputVisitSummary, associations,
762 jointcalControl, camera):
763 """Read the data that jointcal needs to run.
764
765 Modifies ``associations`` in-place with the loaded data.
766
767 Parameters
768 ----------
769 inputSourceTableVisit : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
770 References to visit-level DataFrames to load the catalog data from.
771 inputVisitSummary : `list` [`lsst.daf.butler.DeferredDatasetHandle`]
772 Visit-level exposure summary catalog with metadata.
773 associations : `lsst.jointcal.Associations`
774 Object to add the loaded data to by constructing new CcdImages.
775 jointcalControl : `jointcal.JointcalControl`
776 Control object for C++ associations management.
778 Camera object for detector geometry.
779
780 Returns
781 -------
782 oldWcsList: `list` [`lsst.afw.geom.SkyWcs`]
783 The original WCS of the input data, to aid in writing tests.
784 bands : `list` [`str`]
785 The filter bands of each input dataset.
786 """
787 oldWcsList = []
788 filters = []
789 load_cat_profile_file = 'jointcal_load_data.prof' if self.config.detailedProfile else ''
790 with lsst.utils.timer.profile(load_cat_profile_file):
791 table = make_schema_table() # every detector catalog has the same layout
792 # No guarantee that the input is in the same order of visits, so we have to map one of them.
793 catalogMap = {ref.dataId['visit']: i for i, ref in enumerate(inputSourceTableVisit)}
794 detectorDict = {detector.getId(): detector for detector in camera}
795
796 columns = None
797
798 for visitSummaryRef in inputVisitSummary:
799 visitSummary = visitSummaryRef.get()
800
801 dataRef = inputSourceTableVisit[catalogMap[visitSummaryRef.dataId['visit']]]
802 if columns is None:
803 inColumns = dataRef.get(component='columns')
804 columns, ixxColumns = get_sourceTable_visit_columns(inColumns,
805 self.config,
807 visitCatalog = dataRef.get(parameters={'columns': columns})
808
809 selected = self.sourceSelector.run(visitCatalog)
810 if len(selected.sourceCat) == 0:
811 self.log.warning("No sources selected in visit %s. Skipping...",
812 visitSummary["visit"][0])
813 continue
814
815 # Build a CcdImage for each detector in this visit.
816 detectors = {id: index for index, id in enumerate(visitSummary['id'])}
817 for id, index in detectors.items():
819 selected.sourceCat,
820 id,
821 ixxColumns,
822 self.config.sourceFluxType,
823 self.log)
824 if catalog is None:
825 continue
826 data = self._make_one_input_data(visitSummary[index], catalog, detectorDict)
827 result = self._build_ccdImage(data, associations, jointcalControl)
828 if result is not None:
829 oldWcsList.append(result.wcs)
830 # A visit has only one band, so we can just use the first.
831 filters.append(data.filter)
832 if len(filters) == 0:
833 raise RuntimeError("No data to process: did source selector remove all sources?")
834 filters = collections.Counter(filters)
835
836 return oldWcsList, filters
837
838 def _make_one_input_data(self, visitRecord, catalog, detectorDict):
839 """Return a data structure for this detector+visit."""
840 return JointcalInputData(visit=visitRecord['visit'],
841 catalog=catalog,
842 visitInfo=visitRecord.getVisitInfo(),
843 detector=detectorDict[visitRecord.getId()],
844 photoCalib=visitRecord.getPhotoCalib(),
845 wcs=visitRecord.getWcs(),
846 bbox=visitRecord.getBBox(),
847 # ExposureRecord doesn't have a FilterLabel yet,
848 # so we have to make one.
849 filter=lsst.afw.image.FilterLabel(band=visitRecord['band'],
850 physical=visitRecord['physical_filter']))
851
852 def _build_ccdImage(self, data, associations, jointcalControl):
853 """
854 Extract the necessary things from this catalog+metadata to add a new
855 ccdImage.
856
857 Parameters
858 ----------
859 data : `JointcalInputData`
860 The loaded input data.
861 associations : `lsst.jointcal.Associations`
862 Object to add the info to, to construct a new CcdImage
863 jointcalControl : `jointcal.JointcalControl`
864 Control object for associations management
865
866 Returns
867 ------
868 namedtuple or `None`
869 ``wcs``
870 The TAN WCS of this image, read from the calexp
872 ``key``
873 A key to identify this dataRef by its visit and ccd ids
874 (`namedtuple`).
875 `None`
876 if there are no sources in the loaded catalog.
877 """
878 if len(data.catalog) == 0:
879 self.log.warning("No sources selected in visit %s ccd %s", data.visit, data.detector.getId())
880 return None
881
882 associations.createCcdImage(data.catalog,
883 data.wcs,
884 data.visitInfo,
885 data.bbox,
886 data.filter.physicalLabel,
887 data.photoCalib,
888 data.detector,
889 data.visit,
890 data.detector.getId(),
891 jointcalControl)
892
893 Result = collections.namedtuple('Result_from_build_CcdImage', ('wcs', 'key'))
894 Key = collections.namedtuple('Key', ('visit', 'ccd'))
895 return Result(data.wcs, Key(data.visit, data.detector.getId()))
896
897 def _getDebugPath(self, filename):
898 """Constructs a path to filename using the configured debug path.
899 """
900 return os.path.join(self.config.debugOutputPath, filename)
901
902 def _prep_sky(self, associations, filters):
903 """Prepare on-sky and other data that must be computed after data has
904 been read.
905 """
906 associations.computeCommonTangentPoint()
907
908 boundingCircle = associations.computeBoundingCircle()
909 center = lsst.geom.SpherePoint(boundingCircle.getCenter())
910 radius = lsst.geom.Angle(boundingCircle.getOpeningAngle().asRadians(), lsst.geom.radians)
911
912 self.log.info(f"Data has center={center} with radius={radius.asDegrees()} degrees.")
913
914 # Determine a default filter band associated with the catalog. See DM-9093
915 defaultFilter = filters.most_common(1)[0][0]
916 self.log.debug("Using '%s' filter for reference flux", defaultFilter.physicalLabel)
917
918 # compute and set the reference epoch of the observations, for proper motion corrections
919 epoch = self._compute_proper_motion_epoch(associations.getCcdImageList())
920 associations.setEpoch(epoch.jyear)
921
922 return boundingCircle, center, radius, defaultFilter, epoch
923
925 """Check whether we should override the refcat coordinate errors, and
926 return the overridden error if necessary.
927
928 Parameters
929 ----------
931 The reference catalog to check for a ``coord_raErr`` field.
932 name : `str`
933 Whether we are doing "astrometry" or "photometry".
934
935 Returns
936 -------
937 refCoordErr : `float`
938 The refcat coordinate error to use, or NaN if we are not overriding
939 those fields.
940
941 Raises
942 ------
944 Raised if the refcat does not contain coordinate errors and
945 ``config.astrometryReferenceErr`` is not set.
946 """
947 # This value doesn't matter for photometry, so just set something to
948 # keep old refcats from causing problems.
949 if name.lower() == "photometry":
950 if 'coord_raErr' not in refCat.schema:
951 return 100
952 else:
953 return float('nan')
954
955 if self.config.astrometryReferenceErr is None and 'coord_raErr' not in refCat.schema:
956 msg = ("Reference catalog does not contain coordinate errors, "
957 "and config.astrometryReferenceErr not supplied.")
958 raise pexConfig.FieldValidationError(JointcalConfig.astrometryReferenceErr,
959 self.config,
960 msg)
961
962 if self.config.astrometryReferenceErr is not None and 'coord_raErr' in refCat.schema:
963 self.log.warning("Overriding reference catalog coordinate errors with %f/coordinate [mas]",
964 self.config.astrometryReferenceErr)
965
966 if self.config.astrometryReferenceErr is None:
967 return float('nan')
968 else:
969 return self.config.astrometryReferenceErr
970
971 def _compute_proper_motion_epoch(self, ccdImageList):
972 """Return the proper motion correction epoch of the provided images.
973
974 Parameters
975 ----------
976 ccdImageList : `list` [`lsst.jointcal.CcdImage`]
977 The images to compute the appropriate epoch for.
978
979 Returns
980 -------
981 epoch : `astropy.time.Time`
982 The date to use for proper motion corrections.
983 """
984 return astropy.time.Time(np.mean([ccdImage.getEpoch() for ccdImage in ccdImageList]),
985 format="jyear",
986 scale="tai")
987
988 def _do_load_refcat_and_fit(self, associations, defaultFilter, center, radius,
989 tract="", match_cut=3.0,
990 reject_bad_fluxes=False, *,
991 name="", refObjLoader=None, referenceSelector=None,
992 fit_function=None, epoch=None):
993 """Load reference catalog, perform the fit, and return the result.
994
995 Parameters
996 ----------
997 associations : `lsst.jointcal.Associations`
998 The star/reference star associations to fit.
999 defaultFilter : `lsst.afw.image.FilterLabel`
1000 filter to load from reference catalog.
1001 center : `lsst.geom.SpherePoint`
1002 ICRS center of field to load from reference catalog.
1003 radius : `lsst.geom.Angle`
1004 On-sky radius to load from reference catalog.
1005 name : `str`
1006 Name of thing being fit: "astrometry" or "photometry".
1008 Reference object loader to use to load a reference catalog.
1010 Selector to use to pick objects from the loaded reference catalog.
1011 fit_function : callable
1012 Function to call to perform fit (takes Associations object).
1013 tract : `str`, optional
1014 Name of tract currently being fit.
1015 match_cut : `float`, optional
1016 Radius in arcseconds to find cross-catalog matches to during
1017 associations.associateCatalogs.
1018 reject_bad_fluxes : `bool`, optional
1019 Reject refCat sources with NaN/inf flux or NaN/0 fluxErr.
1020 epoch : `astropy.time.Time`, optional
1021 Epoch to which to correct refcat proper motion and parallax,
1022 or `None` to not apply such corrections.
1023
1024 Returns
1025 -------
1026 result : `Photometry` or `Astrometry`
1027 Result of `fit_function()`
1028 """
1029 self.log.info("====== Now processing %s...", name)
1030 # TODO: this should not print "trying to invert a singular transformation:"
1031 # if it does that, something's not right about the WCS...
1032 associations.associateCatalogs(match_cut)
1033 add_measurement(self.job, 'jointcal.%s_matched_fittedStars' % name,
1034 associations.fittedStarListSize())
1035
1036 applyColorterms = False if name.lower() == "astrometry" else self.config.applyColorTerms
1037 refCat, fluxField = self._load_reference_catalog(refObjLoader, referenceSelector,
1038 center, radius, defaultFilter,
1039 applyColorterms=applyColorterms,
1040 epoch=epoch)
1041 refCoordErr = self._get_refcat_coordinate_error_override(refCat, name)
1042
1043 associations.collectRefStars(refCat,
1044 self.config.matchCut*lsst.geom.arcseconds,
1045 fluxField,
1046 refCoordinateErr=refCoordErr,
1047 rejectBadFluxes=reject_bad_fluxes)
1048 add_measurement(self.job, 'jointcal.%s_collected_refStars' % name,
1049 associations.refStarListSize())
1050
1051 associations.prepareFittedStars(self.config.minMeasurements)
1052
1053 self._check_star_lists(associations, name)
1054 add_measurement(self.job, 'jointcal.%s_prepared_refStars' % name,
1055 associations.nFittedStarsWithAssociatedRefStar())
1056 add_measurement(self.job, 'jointcal.%s_prepared_fittedStars' % name,
1057 associations.fittedStarListSize())
1058 add_measurement(self.job, 'jointcal.%s_prepared_ccdImages' % name,
1059 associations.nCcdImagesValidForFit())
1060
1061 fit_profile_file = 'jointcal_fit_%s.prof'%name if self.config.detailedProfile else ''
1062 dataName = "{}_{}".format(tract, defaultFilter.physicalLabel)
1063 with lsst.utils.timer.profile(fit_profile_file):
1064 result = fit_function(associations, dataName)
1065 # TODO DM-12446: turn this into a "butler save" somehow.
1066 # Save reference and measurement chi2 contributions for this data
1067 if self.config.writeChi2FilesInitialFinal:
1068 baseName = self._getDebugPath(f"{name}_final_chi2-{dataName}")
1069 result.fit.saveChi2Contributions(baseName+"{type}")
1070 self.log.info("Wrote chi2 contributions files: %s", baseName)
1071
1072 return result
1073
1074 def _load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel,
1075 applyColorterms=False, epoch=None):
1076 """Load the necessary reference catalog sources, convert fluxes to
1077 correct units, and apply color term corrections if requested.
1078
1079 Parameters
1080 ----------
1082 The reference catalog loader to use to get the data.
1084 Source selector to apply to loaded reference catalog.
1085 center : `lsst.geom.SpherePoint`
1086 The center around which to load sources.
1087 radius : `lsst.geom.Angle`
1088 The radius around ``center`` to load sources in.
1089 filterLabel : `lsst.afw.image.FilterLabel`
1090 The camera filter to load fluxes for.
1091 applyColorterms : `bool`
1092 Apply colorterm corrections to the refcat for ``filterName``?
1093 epoch : `astropy.time.Time`, optional
1094 Epoch to which to correct refcat proper motion and parallax,
1095 or `None` to not apply such corrections.
1096
1097 Returns
1098 -------
1100 The loaded reference catalog.
1101 fluxField : `str`
1102 The name of the reference catalog flux field appropriate for ``filterName``.
1103 """
1104 skyCircle = refObjLoader.loadSkyCircle(center,
1105 radius,
1106 filterLabel.bandLabel,
1107 epoch=epoch)
1108
1109 selected = referenceSelector.run(skyCircle.refCat)
1110 # Need memory contiguity to get reference filters as a vector.
1111 if not selected.sourceCat.isContiguous():
1112 refCat = selected.sourceCat.copy(deep=True)
1113 else:
1114 refCat = selected.sourceCat
1115
1116 if applyColorterms:
1117 refCatName = refObjLoader.name
1118 self.log.info("Applying color terms for physical filter=%r reference catalog=%s",
1119 filterLabel.physicalLabel, refCatName)
1120 colorterm = self.config.colorterms.getColorterm(filterLabel.physicalLabel,
1121 refCatName,
1122 doRaise=True)
1123
1124 refMag, refMagErr = colorterm.getCorrectedMagnitudes(refCat)
1125 refCat[skyCircle.fluxField] = u.Magnitude(refMag, u.ABmag).to_value(u.nJy)
1126 # TODO: I didn't want to use this, but I'll deal with it in DM-16903
1127 refCat[skyCircle.fluxField+'Err'] = fluxErrFromABMagErr(refMagErr, refMag) * 1e9
1128
1129 return refCat, skyCircle.fluxField
1130
1131 def _check_star_lists(self, associations, name):
1132 # TODO: these should be len(blah), but we need this properly wrapped first.
1133 if associations.nCcdImagesValidForFit() == 0:
1134 raise RuntimeError('No images in the ccdImageList!')
1135 if associations.fittedStarListSize() == 0:
1136 raise RuntimeError('No stars in the {} fittedStarList!'.format(name))
1137 if associations.refStarListSize() == 0:
1138 raise RuntimeError('No stars in the {} reference star list!'.format(name))
1139
1140 def _logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None):
1141 """Compute chi2, log it, validate the model, and return chi2.
1142
1143 Parameters
1144 ----------
1145 associations : `lsst.jointcal.Associations`
1146 The star/reference star associations to fit.
1148 The fitter to use for minimization.
1149 model : `lsst.jointcal.Model`
1150 The model being fit.
1151 chi2Label : `str`
1152 Label to describe the chi2 (e.g. "Initialized", "Final").
1153 writeChi2Name : `str`, optional
1154 Filename prefix to write the chi2 contributions to.
1155 Do not supply an extension: an appropriate one will be added.
1156
1157 Returns
1158 -------
1160 The chi2 object for the current fitter and model.
1161
1162 Raises
1163 ------
1164 FloatingPointError
1165 Raised if chi2 is infinite or NaN.
1166 ValueError
1167 Raised if the model is not valid.
1168 """
1169 if writeChi2Name is not None:
1170 fullpath = self._getDebugPath(writeChi2Name)
1171 fit.saveChi2Contributions(fullpath+"{type}")
1172 self.log.info("Wrote chi2 contributions files: %s", fullpath)
1173
1174 chi2 = fit.computeChi2()
1175 self.log.info("%s %s", chi2Label, chi2)
1176 self._check_stars(associations)
1177 if not np.isfinite(chi2.chi2):
1178 raise FloatingPointError(f'{chi2Label} chi2 is invalid: {chi2}')
1179 if not model.validate(associations.getCcdImageList(), chi2.ndof):
1180 raise ValueError("Model is not valid: check log messages for warnings.")
1181 return chi2
1182
1183 def _fit_photometry(self, associations, dataName=None):
1184 """
1185 Fit the photometric data.
1186
1187 Parameters
1188 ----------
1189 associations : `lsst.jointcal.Associations`
1190 The star/reference star associations to fit.
1191 dataName : `str`
1192 Name of the data being processed (e.g. "1234_HSC-Y"), for
1193 identifying debugging files.
1194
1195 Returns
1196 -------
1197 fit_result : `namedtuple`
1199 The photometric fitter used to perform the fit.
1201 The photometric model that was fit.
1202 """
1203 self.log.info("=== Starting photometric fitting...")
1204
1205 # TODO: should use pex.config.RegistryField here (see DM-9195)
1206 if self.config.photometryModel == "constrainedFlux":
1207 model = lsst.jointcal.ConstrainedFluxModel(associations.getCcdImageList(),
1208 self.focalPlaneBBox,
1209 visitOrder=self.config.photometryVisitOrder,
1210 errorPedestal=self.config.photometryErrorPedestal)
1211 # potentially nonlinear problem, so we may need a line search to converge.
1212 doLineSearch = self.config.allowLineSearch
1213 elif self.config.photometryModel == "constrainedMagnitude":
1214 model = lsst.jointcal.ConstrainedMagnitudeModel(associations.getCcdImageList(),
1215 self.focalPlaneBBox,
1216 visitOrder=self.config.photometryVisitOrder,
1217 errorPedestal=self.config.photometryErrorPedestal)
1218 # potentially nonlinear problem, so we may need a line search to converge.
1219 doLineSearch = self.config.allowLineSearch
1220 elif self.config.photometryModel == "simpleFlux":
1221 model = lsst.jointcal.SimpleFluxModel(associations.getCcdImageList(),
1222 errorPedestal=self.config.photometryErrorPedestal)
1223 doLineSearch = False # purely linear in model parameters, so no line search needed
1224 elif self.config.photometryModel == "simpleMagnitude":
1225 model = lsst.jointcal.SimpleMagnitudeModel(associations.getCcdImageList(),
1226 errorPedestal=self.config.photometryErrorPedestal)
1227 doLineSearch = False # purely linear in model parameters, so no line search needed
1228
1229 fit = lsst.jointcal.PhotometryFit(associations, model)
1230 # TODO DM-12446: turn this into a "butler save" somehow.
1231 # Save reference and measurement chi2 contributions for this data
1232 if self.config.writeChi2FilesInitialFinal:
1233 baseName = f"photometry_initial_chi2-{dataName}"
1234 else:
1235 baseName = None
1236 if self.config.writeInitialModel:
1237 fullpath = self._getDebugPath(f"initial_photometry_model-{dataName}.txt")
1238 writeModel(model, fullpath, self.log)
1239 self._logChi2AndValidate(associations, fit, model, "Initialized", writeChi2Name=baseName)
1240
1241 def getChi2Name(whatToFit):
1242 if self.config.writeChi2FilesOuterLoop:
1243 return f"photometry_init-%s_chi2-{dataName}" % whatToFit
1244 else:
1245 return None
1246
1247 # The constrained model needs the visit transform fit first; the chip
1248 # transform is initialized from the singleFrame PhotoCalib, so it's close.
1249 if self.config.writeInitMatrix:
1250 dumpMatrixFile = self._getDebugPath(f"photometry_preinit-{dataName}")
1251 else:
1252 dumpMatrixFile = ""
1253 if self.config.photometryModel.startswith("constrained"):
1254 # no line search: should be purely (or nearly) linear,
1255 # and we want a large step size to initialize with.
1256 fit.minimize("ModelVisit", dumpMatrixFile=dumpMatrixFile)
1257 self._logChi2AndValidate(associations, fit, model, "Initialize ModelVisit",
1258 writeChi2Name=getChi2Name("ModelVisit"))
1259 dumpMatrixFile = "" # so we don't redo the output on the next step
1260
1261 fit.minimize("Model", doLineSearch=doLineSearch, dumpMatrixFile=dumpMatrixFile)
1262 self._logChi2AndValidate(associations, fit, model, "Initialize Model",
1263 writeChi2Name=getChi2Name("Model"))
1264
1265 fit.minimize("Fluxes") # no line search: always purely linear.
1266 self._logChi2AndValidate(associations, fit, model, "Initialize Fluxes",
1267 writeChi2Name=getChi2Name("Fluxes"))
1268
1269 fit.minimize("Model Fluxes", doLineSearch=doLineSearch)
1270 self._logChi2AndValidate(associations, fit, model, "Initialize ModelFluxes",
1271 writeChi2Name=getChi2Name("ModelFluxes"))
1272
1273 model.freezeErrorTransform()
1274 self.log.debug("Photometry error scales are frozen.")
1275
1276 chi2 = self._iterate_fit(associations,
1277 fit,
1278 self.config.maxPhotometrySteps,
1279 "photometry",
1280 "Model Fluxes",
1281 doRankUpdate=self.config.photometryDoRankUpdate,
1282 doLineSearch=doLineSearch,
1283 dataName=dataName)
1284
1285 add_measurement(self.job, 'jointcal.photometry_final_chi2', chi2.chi2)
1286 add_measurement(self.job, 'jointcal.photometry_final_ndof', chi2.ndof)
1287 return Photometry(fit, model)
1288
1289 def _fit_astrometry(self, associations, dataName=None):
1290 """
1291 Fit the astrometric data.
1292
1293 Parameters
1294 ----------
1295 associations : `lsst.jointcal.Associations`
1296 The star/reference star associations to fit.
1297 dataName : `str`
1298 Name of the data being processed (e.g. "1234_HSC-Y"), for
1299 identifying debugging files.
1300
1301 Returns
1302 -------
1303 fit_result : `namedtuple`
1305 The astrometric fitter used to perform the fit.
1307 The astrometric model that was fit.
1308 sky_to_tan_projection : `lsst.jointcal.ProjectionHandler`
1309 The model for the sky to tangent plane projection that was used in the fit.
1310 """
1311
1312 self.log.info("=== Starting astrometric fitting...")
1313
1314 associations.deprojectFittedStars()
1315
1316 # NOTE: need to return sky_to_tan_projection so that it doesn't get garbage collected.
1317 # TODO: could we package sky_to_tan_projection and model together so we don't have to manage
1318 # them so carefully?
1319 sky_to_tan_projection = lsst.jointcal.OneTPPerVisitHandler(associations.getCcdImageList())
1320
1321 if self.config.astrometryModel == "constrained":
1322 model = lsst.jointcal.ConstrainedAstrometryModel(associations.getCcdImageList(),
1323 sky_to_tan_projection,
1324 chipOrder=self.config.astrometryChipOrder,
1325 visitOrder=self.config.astrometryVisitOrder)
1326 elif self.config.astrometryModel == "simple":
1327 model = lsst.jointcal.SimpleAstrometryModel(associations.getCcdImageList(),
1328 sky_to_tan_projection,
1329 self.config.useInputWcs,
1330 nNotFit=0,
1331 order=self.config.astrometrySimpleOrder)
1332
1333 fit = lsst.jointcal.AstrometryFit(associations, model, self.config.positionErrorPedestal)
1334 # TODO DM-12446: turn this into a "butler save" somehow.
1335 # Save reference and measurement chi2 contributions for this data
1336 if self.config.writeChi2FilesInitialFinal:
1337 baseName = f"astrometry_initial_chi2-{dataName}"
1338 else:
1339 baseName = None
1340 if self.config.writeInitialModel:
1341 fullpath = self._getDebugPath(f"initial_astrometry_model-{dataName}.txt")
1342 writeModel(model, fullpath, self.log)
1343 self._logChi2AndValidate(associations, fit, model, "Initial", writeChi2Name=baseName)
1344
1345 def getChi2Name(whatToFit):
1346 if self.config.writeChi2FilesOuterLoop:
1347 return f"astrometry_init-%s_chi2-{dataName}" % whatToFit
1348 else:
1349 return None
1350
1351 if self.config.writeInitMatrix:
1352 dumpMatrixFile = self._getDebugPath(f"astrometry_preinit-{dataName}")
1353 else:
1354 dumpMatrixFile = ""
1355 # The constrained model needs the visit transform fit first; the chip
1356 # transform is initialized from the detector's cameraGeom, so it's close.
1357 if self.config.astrometryModel == "constrained":
1358 fit.minimize("DistortionsVisit", dumpMatrixFile=dumpMatrixFile)
1359 self._logChi2AndValidate(associations, fit, model, "Initialize DistortionsVisit",
1360 writeChi2Name=getChi2Name("DistortionsVisit"))
1361 dumpMatrixFile = "" # so we don't redo the output on the next step
1362
1363 fit.minimize("Distortions", dumpMatrixFile=dumpMatrixFile)
1364 self._logChi2AndValidate(associations, fit, model, "Initialize Distortions",
1365 writeChi2Name=getChi2Name("Distortions"))
1366
1367 fit.minimize("Positions")
1368 self._logChi2AndValidate(associations, fit, model, "Initialize Positions",
1369 writeChi2Name=getChi2Name("Positions"))
1370
1371 fit.minimize("Distortions Positions")
1372 self._logChi2AndValidate(associations, fit, model, "Initialize DistortionsPositions",
1373 writeChi2Name=getChi2Name("DistortionsPositions"))
1374
1375 chi2 = self._iterate_fit(associations,
1376 fit,
1377 self.config.maxAstrometrySteps,
1378 "astrometry",
1379 "Distortions Positions",
1380 sigmaRelativeTolerance=self.config.astrometryOutlierRelativeTolerance,
1381 doRankUpdate=self.config.astrometryDoRankUpdate,
1382 dataName=dataName)
1383
1384 add_measurement(self.job, 'jointcal.astrometry_final_chi2', chi2.chi2)
1385 add_measurement(self.job, 'jointcal.astrometry_final_ndof', chi2.ndof)
1386
1387 return Astrometry(fit, model, sky_to_tan_projection)
1388
1389 def _check_stars(self, associations):
1390 """Count measured and reference stars per ccd and warn/log them."""
1391 for ccdImage in associations.getCcdImageList():
1392 nMeasuredStars, nRefStars = ccdImage.countStars()
1393 self.log.debug("ccdImage %s has %s measured and %s reference stars",
1394 ccdImage.getName(), nMeasuredStars, nRefStars)
1395 if nMeasuredStars < self.config.minMeasuredStarsPerCcd:
1396 self.log.warning("ccdImage %s has only %s measuredStars (desired %s)",
1397 ccdImage.getName(), nMeasuredStars, self.config.minMeasuredStarsPerCcd)
1398 if nRefStars < self.config.minRefStarsPerCcd:
1399 self.log.warning("ccdImage %s has only %s RefStars (desired %s)",
1400 ccdImage.getName(), nRefStars, self.config.minRefStarsPerCcd)
1401
1402 def _iterate_fit(self, associations, fitter, max_steps, name, whatToFit,
1403 dataName="",
1404 sigmaRelativeTolerance=0,
1405 doRankUpdate=True,
1406 doLineSearch=False):
1407 """Run fitter.minimize up to max_steps times, returning the final chi2.
1408
1409 Parameters
1410 ----------
1411 associations : `lsst.jointcal.Associations`
1412 The star/reference star associations to fit.
1413 fitter : `lsst.jointcal.FitterBase`
1414 The fitter to use for minimization.
1415 max_steps : `int`
1416 Maximum number of steps to run outlier rejection before declaring
1417 convergence failure.
1418 name : {'photometry' or 'astrometry'}
1419 What type of data are we fitting (for logs and debugging files).
1420 whatToFit : `str`
1421 Passed to ``fitter.minimize()`` to define the parameters to fit.
1422 dataName : `str`, optional
1423 Descriptive name for this dataset (e.g. tract and filter),
1424 for debugging.
1425 sigmaRelativeTolerance : `float`, optional
1426 Convergence tolerance for the fractional change in the chi2 cut
1427 level for determining outliers. If set to zero, iterations will
1428 continue until there are no outliers.
1429 doRankUpdate : `bool`, optional
1430 Do an Eigen rank update during minimization, or recompute the full
1431 matrix and gradient?
1432 doLineSearch : `bool`, optional
1433 Do a line search for the optimum step during minimization?
1434
1435 Returns
1436 -------
1438 The final chi2 after the fit converges, or is forced to end.
1439
1440 Raises
1441 ------
1442 FloatingPointError
1443 Raised if the fitter fails with a non-finite value.
1444 RuntimeError
1445 Raised if the fitter fails for some other reason;
1446 log messages will provide further details.
1447 """
1448 if self.config.writeInitMatrix:
1449 dumpMatrixFile = self._getDebugPath(f"{name}_postinit-{dataName}")
1450 else:
1451 dumpMatrixFile = ""
1452 oldChi2 = lsst.jointcal.Chi2Statistic()
1453 oldChi2.chi2 = float("inf")
1454 for i in range(max_steps):
1455 if self.config.writeChi2FilesOuterLoop:
1456 writeChi2Name = f"{name}_iterate_{i}_chi2-{dataName}"
1457 else:
1458 writeChi2Name = None
1459 result = fitter.minimize(whatToFit,
1460 self.config.outlierRejectSigma,
1461 sigmaRelativeTolerance=sigmaRelativeTolerance,
1462 doRankUpdate=doRankUpdate,
1463 doLineSearch=doLineSearch,
1464 dumpMatrixFile=dumpMatrixFile)
1465 dumpMatrixFile = "" # clear it so we don't write the matrix again.
1466 chi2 = self._logChi2AndValidate(associations, fitter, fitter.getModel(),
1467 f"Fit iteration {i}", writeChi2Name=writeChi2Name)
1468
1469 if result == MinimizeResult.Converged:
1470 if doRankUpdate:
1471 self.log.debug("fit has converged - no more outliers - redo minimization "
1472 "one more time in case we have lost accuracy in rank update.")
1473 # Redo minimization one more time in case we have lost accuracy in rank update
1474 result = fitter.minimize(whatToFit, self.config.outlierRejectSigma,
1475 sigmaRelativeTolerance=sigmaRelativeTolerance)
1476 chi2 = self._logChi2AndValidate(associations, fitter, fitter.getModel(), "Fit completed")
1477
1478 # log a message for a large final chi2, TODO: DM-15247 for something better
1479 if chi2.chi2/chi2.ndof >= 4.0:
1480 self.log.error("Potentially bad fit: High chi-squared/ndof.")
1481
1482 break
1483 elif result == MinimizeResult.Chi2Increased:
1484 self.log.warning("Still some outliers remaining but chi2 increased - retry")
1485 # Check whether the increase was large enough to cause trouble.
1486 chi2Ratio = chi2.chi2 / oldChi2.chi2
1487 if chi2Ratio > 1.5:
1488 self.log.warning('Significant chi2 increase by a factor of %.4g / %.4g = %.4g',
1489 chi2.chi2, oldChi2.chi2, chi2Ratio)
1490 # Based on a variety of HSC jointcal logs (see DM-25779), it
1491 # appears that chi2 increases more than a factor of ~2 always
1492 # result in the fit diverging rapidly and ending at chi2 > 1e10.
1493 # Using 10 as the "failure" threshold gives some room between
1494 # leaving a warning and bailing early.
1495 if chi2Ratio > 10:
1496 msg = ("Large chi2 increase between steps: fit likely cannot converge."
1497 " Try setting one or more of the `writeChi2*` config fields and looking"
1498 " at how individual star chi2-values evolve during the fit.")
1499 raise RuntimeError(msg)
1500 oldChi2 = chi2
1501 elif result == MinimizeResult.NonFinite:
1502 filename = self._getDebugPath("{}_failure-nonfinite_chi2-{}.csv".format(name, dataName))
1503 # TODO DM-12446: turn this into a "butler save" somehow.
1504 fitter.saveChi2Contributions(filename+"{type}")
1505 msg = "Nonfinite value in chi2 minimization, cannot complete fit. Dumped star tables to: {}"
1506 raise FloatingPointError(msg.format(filename))
1507 elif result == MinimizeResult.Failed:
1508 raise RuntimeError("Chi2 minimization failure, cannot complete fit.")
1509 else:
1510 raise RuntimeError("Unxepected return code from minimize().")
1511 else:
1512 self.log.error("%s failed to converge after %d steps"%(name, max_steps))
1513
1514 return chi2
1515
1516 def _make_output(self, ccdImageList, model, func):
1517 """Return the internal jointcal models converted to the afw
1518 structures that will be saved to disk.
1519
1520 Parameters
1521 ----------
1522 ccdImageList : `lsst.jointcal.CcdImageList`
1523 The list of CcdImages to get the output for.
1525 The internal jointcal model to convert for each `lsst.jointcal.CcdImage`.
1526 func : `str`
1527 The name of the function to call on ``model`` to get the converted
1528 structure. Must accept an `lsst.jointcal.CcdImage`.
1529
1530 Returns
1531 -------
1532 output : `dict` [`tuple`, `lsst.jointcal.AstrometryModel`] or
1533 `dict` [`tuple`, `lsst.jointcal.PhotometryModel`]
1534 The data to be saved, keyed on (visit, detector).
1535 """
1536 output = {}
1537 for ccdImage in ccdImageList:
1538 ccd = ccdImage.ccdId
1539 visit = ccdImage.visit
1540 self.log.debug("%s for visit: %d, ccd: %d", func, visit, ccd)
1541 output[(visit, ccd)] = getattr(model, func)(ccdImage)
1542 return output
1543
1544
1546 """Return an afw SourceTable to use as a base for creating the
1547 SourceCatalog to insert values from the dataFrame into.
1548
1549 Returns
1550 -------
1552 Table with schema and slots to use to make SourceCatalogs.
1553 """
1555 schema.addField("centroid_x", "D")
1556 schema.addField("centroid_y", "D")
1557 schema.addField("centroid_xErr", "F")
1558 schema.addField("centroid_yErr", "F")
1559 schema.addField("shape_xx", "D")
1560 schema.addField("shape_yy", "D")
1561 schema.addField("shape_xy", "D")
1562 schema.addField("flux_instFlux", "D")
1563 schema.addField("flux_instFluxErr", "D")
1564 table = lsst.afw.table.SourceTable.make(schema)
1565 table.defineCentroid("centroid")
1566 table.defineShape("shape")
1567 return table
1568
1569
1570def get_sourceTable_visit_columns(inColumns, config, sourceSelector):
1571 """
1572 Get the sourceTable_visit columns to load from the catalogs.
1573
1574 Parameters
1575 ----------
1576 inColumns : `list`
1577 List of columns known to be available in the sourceTable_visit.
1578 config : `JointcalConfig`
1579 A filled-in config to to help define column names.
1581 A configured source selector to define column names to load.
1582
1583 Returns
1584 -------
1585 columns : `list`
1586 List of columns to read from sourceTable_visit.
1587 ixxColumns : `list`
1588 Name of the ixx/iyy/ixy columns.
1589 """
1590 columns = ['visit', 'detector',
1591 'sourceId', 'x', 'xErr', 'y', 'yErr',
1592 config.sourceFluxType + '_instFlux', config.sourceFluxType + '_instFluxErr']
1593
1594 if 'ixx' in inColumns:
1595 # New columns post-DM-31825
1596 ixxColumns = ['ixx', 'iyy', 'ixy']
1597 else:
1598 # Old columns pre-DM-31825
1599 ixxColumns = ['Ixx', 'Iyy', 'Ixy']
1600 columns.extend(ixxColumns)
1601
1602 if sourceSelector.config.doFlags:
1603 columns.extend(sourceSelector.config.flags.bad)
1604 if sourceSelector.config.doUnresolved:
1605 columns.append(sourceSelector.config.unresolved.name)
1606 if sourceSelector.config.doIsolated:
1607 columns.append(sourceSelector.config.isolated.parentName)
1608 columns.append(sourceSelector.config.isolated.nChildName)
1609 if sourceSelector.config.doRequireFiniteRaDec:
1610 columns.append(sourceSelector.config.requireFiniteRaDec.raColName)
1611 columns.append(sourceSelector.config.requireFiniteRaDec.decColName)
1612 if sourceSelector.config.doRequirePrimary:
1613 columns.append(sourceSelector.config.requirePrimary.primaryColName)
1614
1615 return columns, ixxColumns
1616
1617
1618def extract_detector_catalog_from_visit_catalog(table, visitCatalog, detectorId,
1619 ixxColumns, sourceFluxType, log):
1620 """Return an afw SourceCatalog extracted from a visit-level dataframe,
1621 limited to just one detector.
1622
1623 Parameters
1624 ----------
1626 Table factory to use to make the SourceCatalog that will be
1627 populated with data from ``visitCatalog``.
1628 visitCatalog : `pandas.DataFrame`
1629 DataFrame to extract a detector catalog from.
1630 detectorId : `int`
1631 Numeric id of the detector to extract from ``visitCatalog``.
1632 ixxColumns : `list` [`str`]
1633 Names of the ixx/iyy/ixy columns in the catalog.
1634 sourceFluxType : `str`
1635 Name of the catalog field to load instFluxes from.
1636 log : `logging.Logger`
1637 Logging instance to log to.
1638
1639 Returns
1640 -------
1641 catalog : `lsst.afw.table.SourceCatalog`, or `None`
1642 Detector-level catalog extracted from ``visitCatalog``, or `None`
1643 if there was no data to load.
1644 """
1645 # map from dataFrame column to afw table column
1646 mapping = {'x': 'centroid_x',
1647 'y': 'centroid_y',
1648 'xErr': 'centroid_xErr',
1649 'yErr': 'centroid_yErr',
1650 ixxColumns[0]: 'shape_xx',
1651 ixxColumns[1]: 'shape_yy',
1652 ixxColumns[2]: 'shape_xy',
1653 f'{sourceFluxType}_instFlux': 'flux_instFlux',
1654 f'{sourceFluxType}_instFluxErr': 'flux_instFluxErr',
1655 }
1656
1657 catalog = lsst.afw.table.SourceCatalog(table)
1658 matched = visitCatalog['detector'] == detectorId
1659 n = sum(matched)
1660 if n == 0:
1661 return None
1662 catalog.resize(sum(matched))
1663 view = visitCatalog.loc[matched]
1664 catalog['id'] = view.index
1665 for dfCol, afwCol in mapping.items():
1666 catalog[afwCol] = view[dfCol]
1667
1668 log.debug("%d sources selected in visit %d detector %d",
1669 len(catalog),
1670 view['visit'].iloc[0], # all visits in this catalog are the same, so take the first
1671 detectorId)
1672 return catalog
An immutable representation of a camera.
Definition Camera.h:43
A representation of a detector in a mosaic camera.
Definition Detector.h:185
A 2-dimensional celestial WCS that transform pixels to ICRS RA/Dec, using the LSST standard for pixel...
Definition SkyWcs.h:117
A group of labels for a filter in an exposure or coadd.
Definition FilterLabel.h:58
The photometric calibration of an exposure.
Definition PhotoCalib.h:114
Information about a single exposure of an imaging camera.
Definition VisitInfo.h:68
Custom catalog class for ExposureRecord/Table.
Definition Exposure.h:311
static Schema makeMinimalSchema()
Return a minimal schema for Exposure tables and records.
Definition Exposure.h:216
Custom catalog class for record/table subclasses that are guaranteed to have an ID,...
Table class that contains measurements made on a single exposure.
Definition Source.h:217
static std::shared_ptr< SourceTable > make(Schema const &schema, std::shared_ptr< IdFactory > const &idFactory)
Construct a new table.
Definition Source.cc:382
static Schema makeMinimalSchema()
Return a minimal schema for Source tables and records.
Definition Source.h:258
Class for storing ordered metadata with comments.
A class representing an angle.
Definition Angle.h:128
An integer coordinate rectangle.
Definition Box.h:55
Point in an unspecified spherical coordinate system.
Definition SpherePoint.h:57
The class that implements the relations between MeasuredStar and FittedStar.
Class that handles the astrometric least squares problem.
Interface between AstrometryFit and the combinations of Mappings from pixels to some tangent plane (a...
Handler of an actual image from a single CCD.
Definition CcdImage.h:64
Base class for Chi2Statistic and Chi2List, to allow addEntry inside Fitter for either class.
Definition Chi2.h:44
Simple structure to accumulate chi2 and ndof.
Definition Chi2.h:52
A multi-component model, fitting mappings for sensors and visits simultaneously.
Base class for fitters.
Definition FitterBase.h:55
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.
_build_ccdImage(self, data, associations, jointcalControl)
Definition jointcal.py:852
_iterate_fit(self, associations, fitter, max_steps, name, whatToFit, dataName="", sigmaRelativeTolerance=0, doRankUpdate=True, doLineSearch=False)
Definition jointcal.py:1406
_put_metrics(self, butlerQC, job, outputRefs)
Definition jointcal.py:633
_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)
Definition jointcal.py:992
_prep_sky(self, associations, filters)
Definition jointcal.py:902
_compute_proper_motion_epoch(self, ccdImageList)
Definition jointcal.py:971
_make_one_input_data(self, visitRecord, catalog, detectorDict)
Definition jointcal.py:838
_check_star_lists(self, associations, name)
Definition jointcal.py:1131
_fit_photometry(self, associations, dataName=None)
Definition jointcal.py:1183
_put_output(self, butlerQC, outputs, outputRefs, camera, setter)
Definition jointcal.py:649
_load_reference_catalog(self, refObjLoader, referenceSelector, center, radius, filterLabel, applyColorterms=False, epoch=None)
Definition jointcal.py:1075
_logChi2AndValidate(self, associations, fit, model, chi2Label, writeChi2Name=None)
Definition jointcal.py:1140
_get_refcat_coordinate_error_override(self, refCat, name)
Definition jointcal.py:924
_load_data(self, inputSourceTableVisit, inputVisitSummary, associations, jointcalControl, camera)
Definition jointcal.py:762
_make_output(self, ccdImageList, model, func)
Definition jointcal.py:1516
runQuantum(self, butlerQC, inputRefs, outputRefs)
Definition jointcal.py:602
_fit_astrometry(self, associations, dataName=None)
Definition jointcal.py:1289
run(self, inputSourceTableVisit, inputVisitSummary, inputCamera, tract=None)
Definition jointcal.py:708
extract_detector_catalog_from_visit_catalog(table, visitCatalog, detectorId, ixxColumns, sourceFluxType, log)
Definition jointcal.py:1619
add_measurement(job, name, value)
Definition jointcal.py:56
writeModel(model, filename, log)
Definition jointcal.py:551
get_sourceTable_visit_columns(inColumns, config, sourceSelector)
Definition jointcal.py:1570
lookupVisitRefCats(datasetType, registry, quantumDataId, collections)
Definition jointcal.py:101
T remove(T... args)
This is a virtual class that allows a lot of freedom in the choice of the projection from "Sky" (wher...