lsst.fgcmcal g9c3108eb99+601a40e866
fgcmOutputProducts.py
Go to the documentation of this file.
1# See COPYRIGHT file at the top of the source tree.
2#
3# This file is part of fgcmcal.
4#
5# Developed for the LSST Data Management System.
6# This product includes software developed by the LSST Project
7# (https://www.lsst.org).
8# See the COPYRIGHT file at the top-level directory of this distribution
9# for details of code ownership.
10#
11# This program is free software: you can redistribute it and/or modify
12# it under the terms of the GNU General Public License as published by
13# the Free Software Foundation, either version 3 of the License, or
14# (at your option) any later version.
15#
16# This program is distributed in the hope that it will be useful,
17# but WITHOUT ANY WARRANTY; without even the implied warranty of
18# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19# GNU General Public License for more details.
20#
21# You should have received a copy of the GNU General Public License
22# along with this program. If not, see <https://www.gnu.org/licenses/>.
23"""Make the final fgcmcal output products.
24
25This task takes the final output from fgcmFitCycle and produces the following
26outputs for use in the DM stack: the FGCM standard stars in a reference
27catalog format; the model atmospheres in "transmission_atmosphere_fgcm"
28format; and the zeropoints in "fgcm_photoCalib" format. Optionally, the
29task can transfer the 'absolute' calibration from a reference catalog
30to put the fgcm standard stars in units of Jansky. This is accomplished
31by matching stars in a sample of healpix pixels, and applying the median
32offset per band.
33"""
34import sys
35import traceback
36import copy
37
38import numpy as np
39import healpy as hp
40import esutil
41from astropy import units
42
43import lsst.daf.base as dafBase
44import lsst.pex.config as pexConfig
45import lsst.pipe.base as pipeBase
46from lsst.pipe.base import connectionTypes
47from lsst.afw.image import TransmissionCurve
48from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask
49from lsst.meas.algorithms import ReferenceObjectLoader
50from lsst.pipe.tasks.photoCal import PhotoCalTask
51import lsst.geom
52import lsst.afw.image as afwImage
53import lsst.afw.math as afwMath
54import lsst.afw.table as afwTable
55from lsst.meas.algorithms import IndexerRegistry
56from lsst.meas.algorithms import DatasetConfig
57from lsst.meas.algorithms.ingestIndexReferenceTask import addRefCatMetadata
58
59from .utilities import computeApproxPixelAreaFields
60from .utilities import lookupStaticCalibrations
61from .utilities import FGCM_ILLEGAL_VALUE
62
63import fgcm
64
65__all__ = ['FgcmOutputProductsConfig', 'FgcmOutputProductsTask', 'FgcmOutputProductsRunner']
66
67
68class FgcmOutputProductsConnections(pipeBase.PipelineTaskConnections,
69 dimensions=("instrument",),
70 defaultTemplates={"cycleNumber": "0"}):
71 camera = connectionTypes.PrerequisiteInput(
72 doc="Camera instrument",
73 name="camera",
74 storageClass="Camera",
75 dimensions=("instrument",),
76 lookupFunction=lookupStaticCalibrations,
77 isCalibration=True,
78 )
79
80 fgcmLookUpTable = connectionTypes.PrerequisiteInput(
81 doc=("Atmosphere + instrument look-up-table for FGCM throughput and "
82 "chromatic corrections."),
83 name="fgcmLookUpTable",
84 storageClass="Catalog",
85 dimensions=("instrument",),
86 deferLoad=True,
87 )
88
89 fgcmVisitCatalog = connectionTypes.Input(
90 doc="Catalog of visit information for fgcm",
91 name="fgcmVisitCatalog",
92 storageClass="Catalog",
93 dimensions=("instrument",),
94 deferLoad=True,
95 )
96
97 fgcmStandardStars = connectionTypes.Input(
98 doc="Catalog of standard star data from fgcm fit",
99 name="fgcmStandardStars{cycleNumber}",
100 storageClass="SimpleCatalog",
101 dimensions=("instrument",),
102 deferLoad=True,
103 )
104
105 fgcmZeropoints = connectionTypes.Input(
106 doc="Catalog of zeropoints from fgcm fit",
107 name="fgcmZeropoints{cycleNumber}",
108 storageClass="Catalog",
109 dimensions=("instrument",),
110 deferLoad=True,
111 )
112
113 fgcmAtmosphereParameters = connectionTypes.Input(
114 doc="Catalog of atmosphere parameters from fgcm fit",
115 name="fgcmAtmosphereParameters{cycleNumber}",
116 storageClass="Catalog",
117 dimensions=("instrument",),
118 deferLoad=True,
119 )
120
121 refCat = connectionTypes.PrerequisiteInput(
122 doc="Reference catalog to use for photometric calibration",
123 name="cal_ref_cat",
124 storageClass="SimpleCatalog",
125 dimensions=("skypix",),
126 deferLoad=True,
127 multiple=True,
128 )
129
130 fgcmPhotoCalib = connectionTypes.Output(
131 doc=("Per-visit photometric calibrations derived from fgcm calibration. "
132 "These catalogs use detector id for the id and are sorted for "
133 "fast lookups of a detector."),
134 name="fgcmPhotoCalibCatalog",
135 storageClass="ExposureCatalog",
136 dimensions=("instrument", "visit",),
137 multiple=True,
138 )
139
140 fgcmTransmissionAtmosphere = connectionTypes.Output(
141 doc="Per-visit atmosphere transmission files produced from fgcm calibration",
142 name="transmission_atmosphere_fgcm",
143 storageClass="TransmissionCurve",
144 dimensions=("instrument",
145 "visit",),
146 multiple=True,
147 )
148
149 fgcmOffsets = connectionTypes.Output(
150 doc="Per-band offsets computed from doReferenceCalibration",
151 name="fgcmReferenceCalibrationOffsets",
152 storageClass="Catalog",
153 dimensions=("instrument",),
154 multiple=False,
155 )
156
157 def __init__(self, *, config=None):
158 super().__init__(config=config)
159
160 if str(int(config.connections.cycleNumber)) != config.connections.cycleNumber:
161 raise ValueError("cycleNumber must be of integer format")
162 if config.connections.refCat != config.refObjLoader.ref_dataset_name:
163 raise ValueError("connections.refCat must be the same as refObjLoader.ref_dataset_name")
164
165 if config.doRefcatOutput:
166 raise ValueError("FgcmOutputProductsTask (Gen3) does not support doRefcatOutput")
167
168 if not config.doReferenceCalibration:
169 self.prerequisiteInputs.remove("refCat")
170 if not config.doAtmosphereOutput:
171 self.inputs.remove("fgcmAtmosphereParameters")
172 if not config.doZeropointOutput:
173 self.inputs.remove("fgcmZeropoints")
174 if not config.doReferenceCalibration:
175 self.outputs.remove("fgcmOffsets")
176
177
178class FgcmOutputProductsConfig(pipeBase.PipelineTaskConfig,
179 pipelineConnections=FgcmOutputProductsConnections):
180 """Config for FgcmOutputProductsTask"""
181
182 cycleNumber = pexConfig.Field(
183 doc="Final fit cycle from FGCM fit",
184 dtype=int,
185 default=None,
186 )
187 physicalFilterMap = pexConfig.DictField(
188 doc="Mapping from 'physicalFilter' to band.",
189 keytype=str,
190 itemtype=str,
191 default={},
192 )
193 # The following fields refer to calibrating from a reference
194 # catalog, but in the future this might need to be expanded
195 doReferenceCalibration = pexConfig.Field(
196 doc=("Transfer 'absolute' calibration from reference catalog? "
197 "This afterburner step is unnecessary if reference stars "
198 "were used in the full fit in FgcmFitCycleTask."),
199 dtype=bool,
200 default=False,
201 )
202 doRefcatOutput = pexConfig.Field(
203 doc="Output standard stars in reference catalog format",
204 dtype=bool,
205 default=True,
206 )
207 doAtmosphereOutput = pexConfig.Field(
208 doc="Output atmospheres in transmission_atmosphere_fgcm format",
209 dtype=bool,
210 default=True,
211 )
212 doZeropointOutput = pexConfig.Field(
213 doc="Output zeropoints in fgcm_photoCalib format",
214 dtype=bool,
215 default=True,
216 )
217 doComposeWcsJacobian = pexConfig.Field(
218 doc="Compose Jacobian of WCS with fgcm calibration for output photoCalib?",
219 dtype=bool,
220 default=True,
221 )
222 doApplyMeanChromaticCorrection = pexConfig.Field(
223 doc="Apply the mean chromatic correction to the zeropoints?",
224 dtype=bool,
225 default=True,
226 )
227 refObjLoader = pexConfig.ConfigurableField(
228 target=LoadIndexedReferenceObjectsTask,
229 doc="reference object loader for 'absolute' photometric calibration",
230 )
231 photoCal = pexConfig.ConfigurableField(
232 target=PhotoCalTask,
233 doc="task to perform 'absolute' calibration",
234 )
235 referencePixelizationNside = pexConfig.Field(
236 doc="Healpix nside to pixelize catalog to compare to reference catalog",
237 dtype=int,
238 default=64,
239 )
240 referencePixelizationMinStars = pexConfig.Field(
241 doc=("Minimum number of stars per healpix pixel to select for comparison"
242 "to the specified reference catalog"),
243 dtype=int,
244 default=200,
245 )
246 referenceMinMatch = pexConfig.Field(
247 doc="Minimum number of stars matched to reference catalog to be used in statistics",
248 dtype=int,
249 default=50,
250 )
251 referencePixelizationNPixels = pexConfig.Field(
252 doc=("Number of healpix pixels to sample to do comparison. "
253 "Doing too many will take a long time and not yield any more "
254 "precise results because the final number is the median offset "
255 "(per band) from the set of pixels."),
256 dtype=int,
257 default=100,
258 )
259 datasetConfig = pexConfig.ConfigField(
260 dtype=DatasetConfig,
261 doc="Configuration for writing/reading ingested catalog",
262 )
263
264 def setDefaults(self):
265 pexConfig.Config.setDefaults(self)
266
267 # In order to transfer the "absolute" calibration from a reference
268 # catalog to the relatively calibrated FGCM standard stars (one number
269 # per band), we use the PhotoCalTask to match stars in a sample of healpix
270 # pixels. These basic settings ensure that only well-measured, good stars
271 # from the source and reference catalogs are used for the matching.
272
273 # applyColorTerms needs to be False if doReferenceCalibration is False,
274 # as is the new default after DM-16702
275 self.photoCal.applyColorTerms = False
276 self.photoCal.fluxField = 'instFlux'
277 self.photoCal.magErrFloor = 0.003
278 self.photoCal.match.referenceSelection.doSignalToNoise = True
279 self.photoCal.match.referenceSelection.signalToNoise.minimum = 10.0
280 self.photoCal.match.sourceSelection.doSignalToNoise = True
281 self.photoCal.match.sourceSelection.signalToNoise.minimum = 10.0
282 self.photoCal.match.sourceSelection.signalToNoise.fluxField = 'instFlux'
283 self.photoCal.match.sourceSelection.signalToNoise.errField = 'instFluxErr'
284 self.photoCal.match.sourceSelection.doFlags = True
285 self.photoCal.match.sourceSelection.flags.good = []
286 self.photoCal.match.sourceSelection.flags.bad = ['flag_badStar']
287 self.photoCal.match.sourceSelection.doUnresolved = False
288 self.datasetConfig.ref_dataset_name = 'fgcm_stars'
289 self.datasetConfig.format_version = 1
290
291 def validate(self):
292 super().validate()
293
294 # Force the connections to conform with cycleNumber
295 self.connections.cycleNumber = str(self.cycleNumber)
296
297
298class FgcmOutputProductsRunner(pipeBase.ButlerInitializedTaskRunner):
299 """Subclass of TaskRunner for fgcmOutputProductsTask
300
301 fgcmOutputProductsTask.run() takes one argument, the butler, and
302 does not run on any data in the repository.
303 This runner does not use any parallelization.
304 """
305
306 @staticmethod
307 def getTargetList(parsedCmd):
308 """
309 Return a list with one element, the butler.
310 """
311 return [parsedCmd.butler]
312
313 def __call__(self, butler):
314 """
315 Parameters
316 ----------
317 butler: `lsst.daf.persistence.Butler`
318
319 Returns
320 -------
321 exitStatus: `list` with `pipeBase.Struct`
322 exitStatus (0: success; 1: failure)
323 if self.doReturnResults also
324 results (`np.array` with absolute zeropoint offsets)
325 """
326 task = self.TaskClass(butler=butler, config=self.config, log=self.log)
327
328 exitStatus = 0
329 if self.doRaise:
330 results = task.runDataRef(butler)
331 else:
332 try:
333 results = task.runDataRef(butler)
334 except Exception as e:
335 exitStatus = 1
336 task.log.fatal("Failed: %s" % e)
337 if not isinstance(e, pipeBase.TaskError):
338 traceback.print_exc(file=sys.stderr)
339
340 task.writeMetadata(butler)
341
342 if self.doReturnResults:
343 # The results here are the zeropoint offsets for each band
344 return [pipeBase.Struct(exitStatus=exitStatus,
345 results=results)]
346 else:
347 return [pipeBase.Struct(exitStatus=exitStatus)]
348
349 def run(self, parsedCmd):
350 """
351 Run the task, with no multiprocessing
352
353 Parameters
354 ----------
355 parsedCmd: `lsst.pipe.base.ArgumentParser` parsed command line
356 """
357
358 resultList = []
359
360 if self.precall(parsedCmd):
361 targetList = self.getTargetList(parsedCmd)
362 # make sure that we only get 1
363 resultList = self(targetList[0])
364
365 return resultList
366
367
368class FgcmOutputProductsTask(pipeBase.PipelineTask, pipeBase.CmdLineTask):
369 """
370 Output products from FGCM global calibration.
371 """
372
373 ConfigClass = FgcmOutputProductsConfig
374 RunnerClass = FgcmOutputProductsRunner
375 _DefaultName = "fgcmOutputProducts"
376
377 def __init__(self, butler=None, **kwargs):
378 super().__init__(**kwargs)
379
380 # no saving of metadata for now
381 def _getMetadataName(self):
382 return None
383
384 def runQuantum(self, butlerQC, inputRefs, outputRefs):
385 dataRefDict = {}
386 dataRefDict['camera'] = butlerQC.get(inputRefs.camera)
387 dataRefDict['fgcmLookUpTable'] = butlerQC.get(inputRefs.fgcmLookUpTable)
388 dataRefDict['fgcmVisitCatalog'] = butlerQC.get(inputRefs.fgcmVisitCatalog)
389 dataRefDict['fgcmStandardStars'] = butlerQC.get(inputRefs.fgcmStandardStars)
390
391 if self.config.doZeropointOutput:
392 dataRefDict['fgcmZeropoints'] = butlerQC.get(inputRefs.fgcmZeropoints)
393 photoCalibRefDict = {photoCalibRef.dataId.byName()['visit']:
394 photoCalibRef for photoCalibRef in outputRefs.fgcmPhotoCalib}
395
396 if self.config.doAtmosphereOutput:
397 dataRefDict['fgcmAtmosphereParameters'] = butlerQC.get(inputRefs.fgcmAtmosphereParameters)
398 atmRefDict = {atmRef.dataId.byName()['visit']: atmRef for
399 atmRef in outputRefs.fgcmTransmissionAtmosphere}
400
401 if self.config.doReferenceCalibration:
402 refConfig = self.config.refObjLoader
403 self.refObjLoader = ReferenceObjectLoader(dataIds=[ref.datasetRef.dataId
404 for ref in inputRefs.refCat],
405 refCats=butlerQC.get(inputRefs.refCat),
406 config=refConfig,
407 log=self.log)
408 else:
409 self.refObjLoader = None
410
411 struct = self.run(dataRefDict, self.config.physicalFilterMap, returnCatalogs=True)
412
413 # Output the photoCalib exposure catalogs
414 if struct.photoCalibCatalogs is not None:
415 self.log.info("Outputting photoCalib catalogs.")
416 for visit, expCatalog in struct.photoCalibCatalogs:
417 butlerQC.put(expCatalog, photoCalibRefDict[visit])
418 self.log.info("Done outputting photoCalib catalogs.")
419
420 # Output the atmospheres
421 if struct.atmospheres is not None:
422 self.log.info("Outputting atmosphere transmission files.")
423 for visit, atm in struct.atmospheres:
424 butlerQC.put(atm, atmRefDict[visit])
425 self.log.info("Done outputting atmosphere files.")
426
427 if self.config.doReferenceCalibration:
428 # Turn offset into simple catalog for persistence if necessary
429 schema = afwTable.Schema()
430 schema.addField('offset', type=np.float64,
431 doc="Post-process calibration offset (mag)")
432 offsetCat = afwTable.BaseCatalog(schema)
433 offsetCat.resize(len(struct.offsets))
434 offsetCat['offset'][:] = struct.offsets
435
436 butlerQC.put(offsetCat, outputRefs.fgcmOffsets)
437
438 return
439
440 @pipeBase.timeMethod
441 def runDataRef(self, butler):
442 """
443 Make FGCM output products for use in the stack
444
445 Parameters
446 ----------
447 butler: `lsst.daf.persistence.Butler`
448 cycleNumber: `int`
449 Final fit cycle number, override config.
450
451 Returns
452 -------
453 offsets: `lsst.pipe.base.Struct`
454 A structure with array of zeropoint offsets
455
456 Raises
457 ------
458 RuntimeError:
459 Raised if any one of the following is true:
460
461 - butler cannot find "fgcmBuildStars_config" or
462 "fgcmBuildStarsTable_config".
463 - butler cannot find "fgcmFitCycle_config".
464 - "fgcmFitCycle_config" does not refer to
465 `self.config.cycleNumber`.
466 - butler cannot find "fgcmAtmosphereParameters" and
467 `self.config.doAtmosphereOutput` is `True`.
468 - butler cannot find "fgcmStandardStars" and
469 `self.config.doReferenceCalibration` is `True` or
470 `self.config.doRefcatOutput` is `True`.
471 - butler cannot find "fgcmZeropoints" and
472 `self.config.doZeropointOutput` is `True`.
473 """
474 if self.config.doReferenceCalibration:
475 # We need the ref obj loader to get the flux field
476 self.makeSubtask("refObjLoader", butler=butler)
477
478 # Check to make sure that the fgcmBuildStars config exists, to retrieve
479 # the visit and ccd dataset tags
480 if not butler.datasetExists('fgcmBuildStarsTable_config') and \
481 not butler.datasetExists('fgcmBuildStars_config'):
482 raise RuntimeError("Cannot find fgcmBuildStarsTable_config or fgcmBuildStars_config, "
483 "which is prereq for fgcmOutputProducts")
484
485 if butler.datasetExists('fgcmBuildStarsTable_config'):
486 fgcmBuildStarsConfig = butler.get('fgcmBuildStarsTable_config')
487 else:
488 fgcmBuildStarsConfig = butler.get('fgcmBuildStars_config')
489 visitDataRefName = fgcmBuildStarsConfig.visitDataRefName
490 ccdDataRefName = fgcmBuildStarsConfig.ccdDataRefName
491 physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap
492
493 if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian:
494 raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied "
495 "in fgcmBuildStarsTask.")
496
497 if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian:
498 self.log.warning("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.")
499
500 # And make sure that the atmosphere was output properly
501 if (self.config.doAtmosphereOutput
502 and not butler.datasetExists('fgcmAtmosphereParameters', fgcmcycle=self.config.cycleNumber)):
503 raise RuntimeError(f"Atmosphere parameters are missing for cycle {self.config.cycleNumber}.")
504
505 if not butler.datasetExists('fgcmStandardStars',
506 fgcmcycle=self.config.cycleNumber):
507 raise RuntimeError("Standard stars are missing for cycle %d." %
508 (self.config.cycleNumber))
509
510 if (self.config.doZeropointOutput
511 and (not butler.datasetExists('fgcmZeropoints', fgcmcycle=self.config.cycleNumber))):
512 raise RuntimeError("Zeropoints are missing for cycle %d." %
513 (self.config.cycleNumber))
514
515 dataRefDict = {}
516 # This is the _actual_ camera
517 dataRefDict['camera'] = butler.get('camera')
518 dataRefDict['fgcmLookUpTable'] = butler.dataRef('fgcmLookUpTable')
519 dataRefDict['fgcmVisitCatalog'] = butler.dataRef('fgcmVisitCatalog')
520 dataRefDict['fgcmStandardStars'] = butler.dataRef('fgcmStandardStars',
521 fgcmcycle=self.config.cycleNumber)
522
523 if self.config.doZeropointOutput:
524 dataRefDict['fgcmZeropoints'] = butler.dataRef('fgcmZeropoints',
525 fgcmcycle=self.config.cycleNumber)
526 if self.config.doAtmosphereOutput:
527 dataRefDict['fgcmAtmosphereParameters'] = butler.dataRef('fgcmAtmosphereParameters',
528 fgcmcycle=self.config.cycleNumber)
529
530 struct = self.run(dataRefDict, physicalFilterMap, butler=butler, returnCatalogs=False)
531
532 if struct.photoCalibs is not None:
533 self.log.info("Outputting photoCalib files.")
534
535 for visit, detector, physicalFilter, photoCalib in struct.photoCalibs:
536 butler.put(photoCalib, 'fgcm_photoCalib',
537 dataId={visitDataRefName: visit,
538 ccdDataRefName: detector,
539 'filter': physicalFilter})
540
541 self.log.info("Done outputting photoCalib files.")
542
543 if struct.atmospheres is not None:
544 self.log.info("Outputting atmosphere transmission files.")
545 for visit, atm in struct.atmospheres:
546 butler.put(atm, "transmission_atmosphere_fgcm",
547 dataId={visitDataRefName: visit})
548 self.log.info("Done outputting atmosphere transmissions.")
549
550 return pipeBase.Struct(offsets=struct.offsets)
551
552 def run(self, dataRefDict, physicalFilterMap, returnCatalogs=True, butler=None):
553 """Run the output products task.
554
555 Parameters
556 ----------
557 dataRefDict : `dict`
558 All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or
559 `lsst.daf.butler.DeferredDatasetHandle` (gen3)
560 dataRef dictionary with keys:
561
562 ``"camera"``
563 Camera object (`lsst.afw.cameraGeom.Camera`)
564 ``"fgcmLookUpTable"``
565 dataRef for the FGCM look-up table.
566 ``"fgcmVisitCatalog"``
567 dataRef for visit summary catalog.
568 ``"fgcmStandardStars"``
569 dataRef for the output standard star catalog.
570 ``"fgcmZeropoints"``
571 dataRef for the zeropoint data catalog.
572 ``"fgcmAtmosphereParameters"``
573 dataRef for the atmosphere parameter catalog.
574 ``"fgcmBuildStarsTableConfig"``
575 Config for `lsst.fgcmcal.fgcmBuildStarsTableTask`.
576 physicalFilterMap : `dict`
577 Dictionary of mappings from physical filter to FGCM band.
578 returnCatalogs : `bool`, optional
579 Return photoCalibs as per-visit exposure catalogs.
580 butler : `lsst.daf.persistence.Butler`, optional
581 Gen2 butler used for reference star outputs
582
583 Returns
584 -------
585 retStruct : `lsst.pipe.base.Struct`
586 Output structure with keys:
587
588 offsets : `np.ndarray`
589 Final reference offsets, per band.
590 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
591 Generator that returns (visit, transmissionCurve) tuples.
592 photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)]
593 Generator that returns (visit, ccd, filtername, photoCalib) tuples.
594 (returned if returnCatalogs is False).
595 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
596 Generator that returns (visit, exposureCatalog) tuples.
597 (returned if returnCatalogs is True).
598 """
599 stdCat = dataRefDict['fgcmStandardStars'].get()
600 md = stdCat.getMetadata()
601 bands = md.getArray('BANDS')
602
603 if self.config.doReferenceCalibration:
604 lutCat = dataRefDict['fgcmLookUpTable'].get()
605 offsets = self._computeReferenceOffsets(stdCat, lutCat, physicalFilterMap, bands)
606 else:
607 offsets = np.zeros(len(bands))
608
609 # This is Gen2 only, and requires the butler.
610 if self.config.doRefcatOutput and butler is not None:
611 self._outputStandardStars(butler, stdCat, offsets, bands, self.config.datasetConfig)
612
613 del stdCat
614
615 if self.config.doZeropointOutput:
616 zptCat = dataRefDict['fgcmZeropoints'].get()
617 visitCat = dataRefDict['fgcmVisitCatalog'].get()
618
619 pcgen = self._outputZeropoints(dataRefDict['camera'], zptCat, visitCat, offsets, bands,
620 physicalFilterMap, returnCatalogs=returnCatalogs)
621 else:
622 pcgen = None
623
624 if self.config.doAtmosphereOutput:
625 atmCat = dataRefDict['fgcmAtmosphereParameters'].get()
626 atmgen = self._outputAtmospheres(dataRefDict, atmCat)
627 else:
628 atmgen = None
629
630 retStruct = pipeBase.Struct(offsets=offsets,
631 atmospheres=atmgen)
632 if returnCatalogs:
633 retStruct.photoCalibCatalogs = pcgen
634 else:
635 retStruct.photoCalibs = pcgen
636
637 return retStruct
638
639 def generateTractOutputProducts(self, dataRefDict, tract,
640 visitCat, zptCat, atmCat, stdCat,
641 fgcmBuildStarsConfig,
642 returnCatalogs=True,
643 butler=None):
644 """
645 Generate the output products for a given tract, as specified in the config.
646
647 This method is here to have an alternate entry-point for
648 FgcmCalibrateTract.
649
650 Parameters
651 ----------
652 dataRefDict : `dict`
653 All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or
654 `lsst.daf.butler.DeferredDatasetHandle` (gen3)
655 dataRef dictionary with keys:
656
657 ``"camera"``
658 Camera object (`lsst.afw.cameraGeom.Camera`)
659 ``"fgcmLookUpTable"``
660 dataRef for the FGCM look-up table.
661 tract : `int`
662 Tract number
663 visitCat : `lsst.afw.table.BaseCatalog`
664 FGCM visitCat from `FgcmBuildStarsTask`
665 zptCat : `lsst.afw.table.BaseCatalog`
666 FGCM zeropoint catalog from `FgcmFitCycleTask`
667 atmCat : `lsst.afw.table.BaseCatalog`
668 FGCM atmosphere parameter catalog from `FgcmFitCycleTask`
669 stdCat : `lsst.afw.table.SimpleCatalog`
670 FGCM standard star catalog from `FgcmFitCycleTask`
671 fgcmBuildStarsConfig : `lsst.fgcmcal.FgcmBuildStarsConfig`
672 Configuration object from `FgcmBuildStarsTask`
673 returnCatalogs : `bool`, optional
674 Return photoCalibs as per-visit exposure catalogs.
675 butler: `lsst.daf.persistence.Butler`, optional
676 Gen2 butler used for reference star outputs
677
678 Returns
679 -------
680 retStruct : `lsst.pipe.base.Struct`
681 Output structure with keys:
682
683 offsets : `np.ndarray`
684 Final reference offsets, per band.
685 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
686 Generator that returns (visit, transmissionCurve) tuples.
687 photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)]
688 Generator that returns (visit, ccd, filtername, photoCalib) tuples.
689 (returned if returnCatalogs is False).
690 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
691 Generator that returns (visit, exposureCatalog) tuples.
692 (returned if returnCatalogs is True).
693 """
694 physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap
695
696 md = stdCat.getMetadata()
697 bands = md.getArray('BANDS')
698
699 if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian:
700 raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied "
701 "in fgcmBuildStarsTask.")
702
703 if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian:
704 self.log.warning("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.")
705
706 if self.config.doReferenceCalibration:
707 lutCat = dataRefDict['fgcmLookUpTable'].get()
708 offsets = self._computeReferenceOffsets(stdCat, lutCat, bands, physicalFilterMap)
709 else:
710 offsets = np.zeros(len(bands))
711
712 if self.config.doRefcatOutput and butler is not None:
713 # Create a special config that has the tract number in it
714 datasetConfig = copy.copy(self.config.datasetConfig)
715 datasetConfig.ref_dataset_name = '%s_%d' % (self.config.datasetConfig.ref_dataset_name,
716 tract)
717 self._outputStandardStars(butler, stdCat, offsets, bands, datasetConfig)
718
719 if self.config.doZeropointOutput:
720 pcgen = self._outputZeropoints(dataRefDict['camera'], zptCat, visitCat, offsets, bands,
721 physicalFilterMap, returnCatalogs=returnCatalogs)
722 else:
723 pcgen = None
724
725 if self.config.doAtmosphereOutput:
726 atmgen = self._outputAtmospheres(dataRefDict, atmCat)
727 else:
728 atmgen = None
729
730 retStruct = pipeBase.Struct(offsets=offsets,
731 atmospheres=atmgen)
732 if returnCatalogs:
733 retStruct.photoCalibCatalogs = pcgen
734 else:
735 retStruct.photoCalibs = pcgen
736
737 return retStruct
738
739 def _computeReferenceOffsets(self, stdCat, lutCat, physicalFilterMap, bands):
740 """
741 Compute offsets relative to a reference catalog.
742
743 This method splits the star catalog into healpix pixels
744 and computes the calibration transfer for a sample of
745 these pixels to approximate the 'absolute' calibration
746 values (on for each band) to apply to transfer the
747 absolute scale.
748
749 Parameters
750 ----------
751 stdCat : `lsst.afw.table.SimpleCatalog`
752 FGCM standard stars
753 lutCat : `lsst.afw.table.SimpleCatalog`
754 FGCM Look-up table
755 physicalFilterMap : `dict`
756 Dictionary of mappings from physical filter to FGCM band.
757 bands : `list` [`str`]
758 List of band names from FGCM output
759 Returns
760 -------
761 offsets : `numpy.array` of floats
762 Per band zeropoint offsets
763 """
764
765 # Only use stars that are observed in all the bands that were actually used
766 # This will ensure that we use the same healpix pixels for the absolute
767 # calibration of each band.
768 minObs = stdCat['ngood'].min(axis=1)
769
770 goodStars = (minObs >= 1)
771 stdCat = stdCat[goodStars]
772
773 self.log.info("Found %d stars with at least 1 good observation in each band" %
774 (len(stdCat)))
775
776 # Associate each band with the appropriate physicalFilter and make
777 # filterLabels
778 filterLabels = []
779
780 lutPhysicalFilters = lutCat[0]['physicalFilters'].split(',')
781 lutStdPhysicalFilters = lutCat[0]['stdPhysicalFilters'].split(',')
782 physicalFilterMapBands = list(physicalFilterMap.values())
783 physicalFilterMapFilters = list(physicalFilterMap.keys())
784 for band in bands:
785 # Find a physical filter associated from the band by doing
786 # a reverse lookup on the physicalFilterMap dict
787 physicalFilterMapIndex = physicalFilterMapBands.index(band)
788 physicalFilter = physicalFilterMapFilters[physicalFilterMapIndex]
789 # Find the appropriate fgcm standard physicalFilter
790 lutPhysicalFilterIndex = lutPhysicalFilters.index(physicalFilter)
791 stdPhysicalFilter = lutStdPhysicalFilters[lutPhysicalFilterIndex]
792 filterLabels.append(afwImage.FilterLabel(band=band,
793 physical=stdPhysicalFilter))
794
795 # We have to make a table for each pixel with flux/fluxErr
796 # This is a temporary table generated for input to the photoCal task.
797 # These fluxes are not instFlux (they are top-of-the-atmosphere approximate and
798 # have had chromatic corrections applied to get to the standard system
799 # specified by the atmosphere/instrumental parameters), nor are they
800 # in Jansky (since they don't have a proper absolute calibration: the overall
801 # zeropoint is estimated from the telescope size, etc.)
802 sourceMapper = afwTable.SchemaMapper(stdCat.schema)
803 sourceMapper.addMinimalSchema(afwTable.SimpleTable.makeMinimalSchema())
804 sourceMapper.editOutputSchema().addField('instFlux', type=np.float64,
805 doc="instrumental flux (counts)")
806 sourceMapper.editOutputSchema().addField('instFluxErr', type=np.float64,
807 doc="instrumental flux error (counts)")
808 badStarKey = sourceMapper.editOutputSchema().addField('flag_badStar',
809 type='Flag',
810 doc="bad flag")
811
812 # Split up the stars
813 # Note that there is an assumption here that the ra/dec coords stored
814 # on-disk are in radians, and therefore that starObs['coord_ra'] /
815 # starObs['coord_dec'] return radians when used as an array of numpy float64s.
816 theta = np.pi/2. - stdCat['coord_dec']
817 phi = stdCat['coord_ra']
818
819 ipring = hp.ang2pix(self.config.referencePixelizationNside, theta, phi)
820 h, rev = esutil.stat.histogram(ipring, rev=True)
821
822 gdpix, = np.where(h >= self.config.referencePixelizationMinStars)
823
824 self.log.info("Found %d pixels (nside=%d) with at least %d good stars" %
825 (gdpix.size,
826 self.config.referencePixelizationNside,
827 self.config.referencePixelizationMinStars))
828
829 if gdpix.size < self.config.referencePixelizationNPixels:
830 self.log.warning("Found fewer good pixels (%d) than preferred in configuration (%d)" %
831 (gdpix.size, self.config.referencePixelizationNPixels))
832 else:
833 # Sample out the pixels we want to use
834 gdpix = np.random.choice(gdpix, size=self.config.referencePixelizationNPixels, replace=False)
835
836 results = np.zeros(gdpix.size, dtype=[('hpix', 'i4'),
837 ('nstar', 'i4', len(bands)),
838 ('nmatch', 'i4', len(bands)),
839 ('zp', 'f4', len(bands)),
840 ('zpErr', 'f4', len(bands))])
841 results['hpix'] = ipring[rev[rev[gdpix]]]
842
843 # We need a boolean index to deal with catalogs...
844 selected = np.zeros(len(stdCat), dtype=bool)
845
846 refFluxFields = [None]*len(bands)
847
848 for p_index, pix in enumerate(gdpix):
849 i1a = rev[rev[pix]: rev[pix + 1]]
850
851 # the stdCat afwTable can only be indexed with boolean arrays,
852 # and not numpy index arrays (see DM-16497). This little trick
853 # converts the index array into a boolean array
854 selected[:] = False
855 selected[i1a] = True
856
857 for b_index, filterLabel in enumerate(filterLabels):
858 struct = self._computeOffsetOneBand(sourceMapper, badStarKey, b_index,
859 filterLabel, stdCat,
860 selected, refFluxFields)
861 results['nstar'][p_index, b_index] = len(i1a)
862 results['nmatch'][p_index, b_index] = len(struct.arrays.refMag)
863 results['zp'][p_index, b_index] = struct.zp
864 results['zpErr'][p_index, b_index] = struct.sigma
865
866 # And compute the summary statistics
867 offsets = np.zeros(len(bands))
868
869 for b_index, band in enumerate(bands):
870 # make configurable
871 ok, = np.where(results['nmatch'][:, b_index] >= self.config.referenceMinMatch)
872 offsets[b_index] = np.median(results['zp'][ok, b_index])
873 # use median absolute deviation to estimate Normal sigma
874 # see https://en.wikipedia.org/wiki/Median_absolute_deviation
875 madSigma = 1.4826*np.median(np.abs(results['zp'][ok, b_index] - offsets[b_index]))
876 self.log.info("Reference catalog offset for %s band: %.12f +/- %.12f",
877 band, offsets[b_index], madSigma)
878
879 return offsets
880
881 def _computeOffsetOneBand(self, sourceMapper, badStarKey,
882 b_index, filterLabel, stdCat, selected, refFluxFields):
883 """
884 Compute the zeropoint offset between the fgcm stdCat and the reference
885 stars for one pixel in one band
886
887 Parameters
888 ----------
889 sourceMapper : `lsst.afw.table.SchemaMapper`
890 Mapper to go from stdCat to calibratable catalog
891 badStarKey : `lsst.afw.table.Key`
892 Key for the field with bad stars
893 b_index : `int`
894 Index of the band in the star catalog
895 filterLabel : `lsst.afw.image.FilterLabel`
896 filterLabel with band and physical filter
897 stdCat : `lsst.afw.table.SimpleCatalog`
898 FGCM standard stars
899 selected : `numpy.array(dtype=bool)`
900 Boolean array of which stars are in the pixel
901 refFluxFields : `list`
902 List of names of flux fields for reference catalog
903 """
904
905 sourceCat = afwTable.SimpleCatalog(sourceMapper.getOutputSchema())
906 sourceCat.reserve(selected.sum())
907 sourceCat.extend(stdCat[selected], mapper=sourceMapper)
908 sourceCat['instFlux'] = 10.**(stdCat['mag_std_noabs'][selected, b_index]/(-2.5))
909 sourceCat['instFluxErr'] = (np.log(10.)/2.5)*(stdCat['magErr_std'][selected, b_index]
910 * sourceCat['instFlux'])
911 # Make sure we only use stars that have valid measurements
912 # (This is perhaps redundant with requirements above that the
913 # stars be observed in all bands, but it can't hurt)
914 badStar = (stdCat['mag_std_noabs'][selected, b_index] > 90.0)
915 for rec in sourceCat[badStar]:
916 rec.set(badStarKey, True)
917
918 exposure = afwImage.ExposureF()
919 exposure.setFilterLabel(filterLabel)
920
921 if refFluxFields[b_index] is None:
922 # Need to find the flux field in the reference catalog
923 # to work around limitations of DirectMatch in PhotoCal
924 ctr = stdCat[0].getCoord()
925 rad = 0.05*lsst.geom.degrees
926 refDataTest = self.refObjLoader.loadSkyCircle(ctr, rad, filterLabel.bandLabel)
927 refFluxFields[b_index] = refDataTest.fluxField
928
929 # Make a copy of the config so that we can modify it
930 calConfig = copy.copy(self.config.photoCal.value)
931 calConfig.match.referenceSelection.signalToNoise.fluxField = refFluxFields[b_index]
932 calConfig.match.referenceSelection.signalToNoise.errField = refFluxFields[b_index] + 'Err'
933 calTask = self.config.photoCal.target(refObjLoader=self.refObjLoader,
934 config=calConfig,
935 schema=sourceCat.getSchema())
936
937 struct = calTask.run(exposure, sourceCat)
938
939 return struct
940
941 def _outputStandardStars(self, butler, stdCat, offsets, bands, datasetConfig):
942 """
943 Output standard stars in indexed reference catalog format.
944 This is not currently supported in Gen3.
945
946 Parameters
947 ----------
948 butler : `lsst.daf.persistence.Butler`
949 stdCat : `lsst.afw.table.SimpleCatalog`
950 FGCM standard star catalog from fgcmFitCycleTask
951 offsets : `numpy.array` of floats
952 Per band zeropoint offsets
953 bands : `list` [`str`]
954 List of band names from FGCM output
955 datasetConfig : `lsst.meas.algorithms.DatasetConfig`
956 Config for reference dataset
957 """
958
959 self.log.info("Outputting standard stars to %s" % (datasetConfig.ref_dataset_name))
960
961 indexer = IndexerRegistry[self.config.datasetConfig.indexer.name](
962 self.config.datasetConfig.indexer.active)
963
964 # We determine the conversion from the native units (typically radians) to
965 # degrees for the first star. This allows us to treat coord_ra/coord_dec as
966 # numpy arrays rather than Angles, which would we approximately 600x slower.
967 # TODO: Fix this after DM-16524 (HtmIndexer.indexPoints should take coords
968 # (as Angles) for input
969 conv = stdCat[0]['coord_ra'].asDegrees()/float(stdCat[0]['coord_ra'])
970 indices = np.array(indexer.indexPoints(stdCat['coord_ra']*conv,
971 stdCat['coord_dec']*conv))
972
973 formattedCat = self._formatCatalog(stdCat, offsets, bands)
974
975 # Write the master schema
976 dataId = indexer.makeDataId('master_schema',
977 datasetConfig.ref_dataset_name)
978 masterCat = afwTable.SimpleCatalog(formattedCat.schema)
979 addRefCatMetadata(masterCat)
980 butler.put(masterCat, 'ref_cat', dataId=dataId)
981
982 # Break up the pixels using a histogram
983 h, rev = esutil.stat.histogram(indices, rev=True)
984 gd, = np.where(h > 0)
985 selected = np.zeros(len(formattedCat), dtype=bool)
986 for i in gd:
987 i1a = rev[rev[i]: rev[i + 1]]
988
989 # the formattedCat afwTable can only be indexed with boolean arrays,
990 # and not numpy index arrays (see DM-16497). This little trick
991 # converts the index array into a boolean array
992 selected[:] = False
993 selected[i1a] = True
994
995 # Write the individual pixel
996 dataId = indexer.makeDataId(indices[i1a[0]],
997 datasetConfig.ref_dataset_name)
998 butler.put(formattedCat[selected], 'ref_cat', dataId=dataId)
999
1000 # And save the dataset configuration
1001 dataId = indexer.makeDataId(None, datasetConfig.ref_dataset_name)
1002 butler.put(datasetConfig, 'ref_cat_config', dataId=dataId)
1003
1004 self.log.info("Done outputting standard stars.")
1005
1006 def _formatCatalog(self, fgcmStarCat, offsets, bands):
1007 """
1008 Turn an FGCM-formatted star catalog, applying zeropoint offsets.
1009
1010 Parameters
1011 ----------
1012 fgcmStarCat : `lsst.afw.Table.SimpleCatalog`
1013 SimpleCatalog as output by fgcmcal
1014 offsets : `list` with len(self.bands) entries
1015 Zeropoint offsets to apply
1016 bands : `list` [`str`]
1017 List of band names from FGCM output
1018
1019 Returns
1020 -------
1021 formattedCat: `lsst.afw.table.SimpleCatalog`
1022 SimpleCatalog suitable for using as a reference catalog
1023 """
1024
1025 sourceMapper = afwTable.SchemaMapper(fgcmStarCat.schema)
1026 minSchema = LoadIndexedReferenceObjectsTask.makeMinimalSchema(bands,
1027 addCentroid=False,
1028 addIsResolved=True,
1029 coordErrDim=0)
1030 sourceMapper.addMinimalSchema(minSchema)
1031 for band in bands:
1032 sourceMapper.editOutputSchema().addField('%s_nGood' % (band), type=np.int32)
1033 sourceMapper.editOutputSchema().addField('%s_nTotal' % (band), type=np.int32)
1034 sourceMapper.editOutputSchema().addField('%s_nPsfCandidate' % (band), type=np.int32)
1035
1036 formattedCat = afwTable.SimpleCatalog(sourceMapper.getOutputSchema())
1037 formattedCat.reserve(len(fgcmStarCat))
1038 formattedCat.extend(fgcmStarCat, mapper=sourceMapper)
1039
1040 # Note that we don't have to set `resolved` because the default is False
1041
1042 for b, band in enumerate(bands):
1043 mag = fgcmStarCat['mag_std_noabs'][:, b].astype(np.float64) + offsets[b]
1044 # We want fluxes in nJy from calibrated AB magnitudes
1045 # (after applying offset). Updated after RFC-549 and RFC-575.
1046 flux = (mag*units.ABmag).to_value(units.nJy)
1047 fluxErr = (np.log(10.)/2.5)*flux*fgcmStarCat['magErr_std'][:, b].astype(np.float64)
1048
1049 formattedCat['%s_flux' % (band)][:] = flux
1050 formattedCat['%s_fluxErr' % (band)][:] = fluxErr
1051 formattedCat['%s_nGood' % (band)][:] = fgcmStarCat['ngood'][:, b]
1052 formattedCat['%s_nTotal' % (band)][:] = fgcmStarCat['ntotal'][:, b]
1053 formattedCat['%s_nPsfCandidate' % (band)][:] = fgcmStarCat['npsfcand'][:, b]
1054
1055 addRefCatMetadata(formattedCat)
1056
1057 return formattedCat
1058
1059 def _outputZeropoints(self, camera, zptCat, visitCat, offsets, bands,
1060 physicalFilterMap, returnCatalogs=True,
1061 tract=None):
1062 """Output the zeropoints in fgcm_photoCalib format.
1063
1064 Parameters
1065 ----------
1066 camera : `lsst.afw.cameraGeom.Camera`
1067 Camera from the butler.
1068 zptCat : `lsst.afw.table.BaseCatalog`
1069 FGCM zeropoint catalog from `FgcmFitCycleTask`.
1070 visitCat : `lsst.afw.table.BaseCatalog`
1071 FGCM visitCat from `FgcmBuildStarsTask`.
1072 offsets : `numpy.array`
1073 Float array of absolute calibration offsets, one for each filter.
1074 bands : `list` [`str`]
1075 List of band names from FGCM output.
1076 physicalFilterMap : `dict`
1077 Dictionary of mappings from physical filter to FGCM band.
1078 returnCatalogs : `bool`, optional
1079 Return photoCalibs as per-visit exposure catalogs.
1080 tract: `int`, optional
1081 Tract number to output. Default is None (global calibration)
1082
1083 Returns
1084 -------
1085 photoCalibs : `generator` [(`int`, `int`, `str`, `lsst.afw.image.PhotoCalib`)]
1086 Generator that returns (visit, ccd, filtername, photoCalib) tuples.
1087 (returned if returnCatalogs is False).
1088 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
1089 Generator that returns (visit, exposureCatalog) tuples.
1090 (returned if returnCatalogs is True).
1091 """
1092 # Select visit/ccds where we have a calibration
1093 # This includes ccds where we were able to interpolate from neighboring
1094 # ccds.
1095 cannot_compute = fgcm.fgcmUtilities.zpFlagDict['CANNOT_COMPUTE_ZEROPOINT']
1096 selected = (((zptCat['fgcmFlag'] & cannot_compute) == 0)
1097 & (zptCat['fgcmZptVar'] > 0.0)
1098 & (zptCat['fgcmZpt'] > FGCM_ILLEGAL_VALUE))
1099
1100 # Log warnings for any visit which has no valid zeropoints
1101 badVisits = np.unique(zptCat['visit'][~selected])
1102 goodVisits = np.unique(zptCat['visit'][selected])
1103 allBadVisits = badVisits[~np.isin(badVisits, goodVisits)]
1104 for allBadVisit in allBadVisits:
1105 self.log.warning(f'No suitable photoCalib for visit {allBadVisit}')
1106
1107 # Get a mapping from filtername to the offsets
1108 offsetMapping = {}
1109 for f in physicalFilterMap:
1110 # Not every filter in the map will necesarily have a band.
1111 if physicalFilterMap[f] in bands:
1112 offsetMapping[f] = offsets[bands.index(physicalFilterMap[f])]
1113
1114 # Get a mapping from "ccd" to the ccd index used for the scaling
1115 ccdMapping = {}
1116 for ccdIndex, detector in enumerate(camera):
1117 ccdMapping[detector.getId()] = ccdIndex
1118
1119 # And a mapping to get the flat-field scaling values
1120 scalingMapping = {}
1121 for rec in visitCat:
1122 scalingMapping[rec['visit']] = rec['scaling']
1123
1124 if self.config.doComposeWcsJacobian:
1125 approxPixelAreaFields = computeApproxPixelAreaFields(camera)
1126
1127 # The zptCat is sorted by visit, which is useful
1128 lastVisit = -1
1129 zptVisitCatalog = None
1130
1131 metadata = dafBase.PropertyList()
1132 metadata.add("COMMENT", "Catalog id is detector id, sorted.")
1133 metadata.add("COMMENT", "Only detectors with data have entries.")
1134
1135 for rec in zptCat[selected]:
1136 # Retrieve overall scaling
1137 scaling = scalingMapping[rec['visit']][ccdMapping[rec['detector']]]
1138
1139 # The postCalibrationOffset describe any zeropoint offsets
1140 # to apply after the fgcm calibration. The first part comes
1141 # from the reference catalog match (used in testing). The
1142 # second part comes from the mean chromatic correction
1143 # (if configured).
1144 postCalibrationOffset = offsetMapping[rec['filtername']]
1145 if self.config.doApplyMeanChromaticCorrection:
1146 postCalibrationOffset += rec['fgcmDeltaChrom']
1147
1148 fgcmSuperStarField = self._getChebyshevBoundedField(rec['fgcmfZptSstarCheb'],
1149 rec['fgcmfZptChebXyMax'])
1150 # Convert from FGCM AB to nJy
1151 fgcmZptField = self._getChebyshevBoundedField((rec['fgcmfZptCheb']*units.AB).to_value(units.nJy),
1152 rec['fgcmfZptChebXyMax'],
1153 offset=postCalibrationOffset,
1154 scaling=scaling)
1155
1156 if self.config.doComposeWcsJacobian:
1157
1158 fgcmField = afwMath.ProductBoundedField([approxPixelAreaFields[rec['detector']],
1159 fgcmSuperStarField,
1160 fgcmZptField])
1161 else:
1162 # The photoCalib is just the product of the fgcmSuperStarField and the
1163 # fgcmZptField
1164 fgcmField = afwMath.ProductBoundedField([fgcmSuperStarField, fgcmZptField])
1165
1166 # The "mean" calibration will be set to the center of the ccd for reference
1167 calibCenter = fgcmField.evaluate(fgcmField.getBBox().getCenter())
1168 calibErr = (np.log(10.0)/2.5)*calibCenter*np.sqrt(rec['fgcmZptVar'])
1169 photoCalib = afwImage.PhotoCalib(calibrationMean=calibCenter,
1170 calibrationErr=calibErr,
1171 calibration=fgcmField,
1172 isConstant=False)
1173
1174 if not returnCatalogs:
1175 # Return individual photoCalibs
1176 yield (int(rec['visit']), int(rec['detector']), rec['filtername'], photoCalib)
1177 else:
1178 # Return full per-visit exposure catalogs
1179 if rec['visit'] != lastVisit:
1180 # This is a new visit. If the last visit was not -1, yield
1181 # the ExposureCatalog
1182 if lastVisit > -1:
1183 # ensure that the detectors are in sorted order, for fast lookups
1184 zptVisitCatalog.sort()
1185 yield (int(lastVisit), zptVisitCatalog)
1186 else:
1187 # We need to create a new schema
1188 zptExpCatSchema = afwTable.ExposureTable.makeMinimalSchema()
1189 zptExpCatSchema.addField('visit', type='I', doc='Visit number')
1190
1191 # And start a new one
1192 zptVisitCatalog = afwTable.ExposureCatalog(zptExpCatSchema)
1193 zptVisitCatalog.setMetadata(metadata)
1194
1195 lastVisit = int(rec['visit'])
1196
1197 catRecord = zptVisitCatalog.addNew()
1198 catRecord['id'] = int(rec['detector'])
1199 catRecord['visit'] = rec['visit']
1200 catRecord.setPhotoCalib(photoCalib)
1201
1202 # Final output of last exposure catalog
1203 if returnCatalogs:
1204 # ensure that the detectors are in sorted order, for fast lookups
1205 zptVisitCatalog.sort()
1206 yield (int(lastVisit), zptVisitCatalog)
1207
1208 def _getChebyshevBoundedField(self, coefficients, xyMax, offset=0.0, scaling=1.0):
1209 """
1210 Make a ChebyshevBoundedField from fgcm coefficients, with optional offset
1211 and scaling.
1212
1213 Parameters
1214 ----------
1215 coefficients: `numpy.array`
1216 Flattened array of chebyshev coefficients
1217 xyMax: `list` of length 2
1218 Maximum x and y of the chebyshev bounding box
1219 offset: `float`, optional
1220 Absolute calibration offset. Default is 0.0
1221 scaling: `float`, optional
1222 Flat scaling value from fgcmBuildStars. Default is 1.0
1223
1224 Returns
1225 -------
1226 boundedField: `lsst.afw.math.ChebyshevBoundedField`
1227 """
1228
1229 orderPlus1 = int(np.sqrt(coefficients.size))
1230 pars = np.zeros((orderPlus1, orderPlus1))
1231
1232 bbox = lsst.geom.Box2I(lsst.geom.Point2I(0.0, 0.0),
1233 lsst.geom.Point2I(*xyMax))
1234
1235 pars[:, :] = (coefficients.reshape(orderPlus1, orderPlus1)
1236 * (10.**(offset/-2.5))*scaling)
1237
1238 boundedField = afwMath.ChebyshevBoundedField(bbox, pars)
1239
1240 return boundedField
1241
1242 def _outputAtmospheres(self, dataRefDict, atmCat):
1243 """
1244 Output the atmospheres.
1245
1246 Parameters
1247 ----------
1248 dataRefDict : `dict`
1249 All dataRefs are `lsst.daf.persistence.ButlerDataRef` (gen2) or
1250 `lsst.daf.butler.DeferredDatasetHandle` (gen3)
1251 dataRef dictionary with keys:
1252
1253 ``"fgcmLookUpTable"``
1254 dataRef for the FGCM look-up table.
1255 atmCat : `lsst.afw.table.BaseCatalog`
1256 FGCM atmosphere parameter catalog from fgcmFitCycleTask.
1257
1258 Returns
1259 -------
1260 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
1261 Generator that returns (visit, transmissionCurve) tuples.
1262 """
1263 # First, we need to grab the look-up table and key info
1264 lutCat = dataRefDict['fgcmLookUpTable'].get()
1265
1266 atmosphereTableName = lutCat[0]['tablename']
1267 elevation = lutCat[0]['elevation']
1268 atmLambda = lutCat[0]['atmLambda']
1269 lutCat = None
1270
1271 # Make the atmosphere table if possible
1272 try:
1273 atmTable = fgcm.FgcmAtmosphereTable.initWithTableName(atmosphereTableName)
1274 atmTable.loadTable()
1275 except IOError:
1276 atmTable = None
1277
1278 if atmTable is None:
1279 # Try to use MODTRAN instead
1280 try:
1281 modGen = fgcm.ModtranGenerator(elevation)
1282 lambdaRange = np.array([atmLambda[0], atmLambda[-1]])/10.
1283 lambdaStep = (atmLambda[1] - atmLambda[0])/10.
1284 except (ValueError, IOError) as e:
1285 raise RuntimeError("FGCM look-up-table generated with modtran, "
1286 "but modtran not configured to run.") from e
1287
1288 zenith = np.degrees(np.arccos(1./atmCat['secZenith']))
1289
1290 for i, visit in enumerate(atmCat['visit']):
1291 if atmTable is not None:
1292 # Interpolate the atmosphere table
1293 atmVals = atmTable.interpolateAtmosphere(pmb=atmCat[i]['pmb'],
1294 pwv=atmCat[i]['pwv'],
1295 o3=atmCat[i]['o3'],
1296 tau=atmCat[i]['tau'],
1297 alpha=atmCat[i]['alpha'],
1298 zenith=zenith[i],
1299 ctranslamstd=[atmCat[i]['cTrans'],
1300 atmCat[i]['lamStd']])
1301 else:
1302 # Run modtran
1303 modAtm = modGen(pmb=atmCat[i]['pmb'],
1304 pwv=atmCat[i]['pwv'],
1305 o3=atmCat[i]['o3'],
1306 tau=atmCat[i]['tau'],
1307 alpha=atmCat[i]['alpha'],
1308 zenith=zenith[i],
1309 lambdaRange=lambdaRange,
1310 lambdaStep=lambdaStep,
1311 ctranslamstd=[atmCat[i]['cTrans'],
1312 atmCat[i]['lamStd']])
1313 atmVals = modAtm['COMBINED']
1314
1315 # Now need to create something to persist...
1316 curve = TransmissionCurve.makeSpatiallyConstant(throughput=atmVals,
1317 wavelengths=atmLambda,
1318 throughputAtMin=atmVals[0],
1319 throughputAtMax=atmVals[-1])
1320
1321 yield (int(visit), curve)
def computeApproxPixelAreaFields(camera)
Definition: utilities.py:495