23"""Make the final fgcmcal output products.
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
39from astropy import units
41import lsst.daf.base as dafBase
42import lsst.pex.config as pexConfig
43import lsst.pipe.base as pipeBase
44from lsst.pipe.base import connectionTypes
45from lsst.afw.image import TransmissionCurve
46from lsst.meas.algorithms import ReferenceObjectLoader, LoadReferenceObjectsConfig
47from lsst.pipe.tasks.photoCal import PhotoCalTask
49import lsst.afw.image as afwImage
50import lsst.afw.math as afwMath
51import lsst.afw.table as afwTable
53from .utilities import computeApproxPixelAreaFields
54from .utilities import lookupStaticCalibrations
55from .utilities import FGCM_ILLEGAL_VALUE
59__all__ = ['FgcmOutputProductsConfig', 'FgcmOutputProductsTask']
62class FgcmOutputProductsConnections(pipeBase.PipelineTaskConnections,
63 dimensions=("instrument",),
64 defaultTemplates={
"cycleNumber":
"0"}):
65 camera = connectionTypes.PrerequisiteInput(
66 doc=
"Camera instrument",
68 storageClass=
"Camera",
69 dimensions=(
"instrument",),
70 lookupFunction=lookupStaticCalibrations,
74 fgcmLookUpTable = connectionTypes.PrerequisiteInput(
75 doc=(
"Atmosphere + instrument look-up-table for FGCM throughput and "
76 "chromatic corrections."),
77 name=
"fgcmLookUpTable",
78 storageClass=
"Catalog",
79 dimensions=(
"instrument",),
83 fgcmVisitCatalog = connectionTypes.Input(
84 doc=
"Catalog of visit information for fgcm",
85 name=
"fgcmVisitCatalog",
86 storageClass=
"Catalog",
87 dimensions=(
"instrument",),
91 fgcmStandardStars = connectionTypes.Input(
92 doc=
"Catalog of standard star data from fgcm fit",
93 name=
"fgcmStandardStars{cycleNumber}",
94 storageClass=
"SimpleCatalog",
95 dimensions=(
"instrument",),
99 fgcmZeropoints = connectionTypes.Input(
100 doc=
"Catalog of zeropoints from fgcm fit",
101 name=
"fgcmZeropoints{cycleNumber}",
102 storageClass=
"Catalog",
103 dimensions=(
"instrument",),
107 fgcmAtmosphereParameters = connectionTypes.Input(
108 doc=
"Catalog of atmosphere parameters from fgcm fit",
109 name=
"fgcmAtmosphereParameters{cycleNumber}",
110 storageClass=
"Catalog",
111 dimensions=(
"instrument",),
115 refCat = connectionTypes.PrerequisiteInput(
116 doc=
"Reference catalog to use for photometric calibration",
118 storageClass=
"SimpleCatalog",
119 dimensions=(
"skypix",),
124 fgcmPhotoCalib = connectionTypes.Output(
125 doc=(
"Per-visit photometric calibrations derived from fgcm calibration. "
126 "These catalogs use detector id for the id and are sorted for "
127 "fast lookups of a detector."),
128 name=
"fgcmPhotoCalibCatalog",
129 storageClass=
"ExposureCatalog",
130 dimensions=(
"instrument",
"visit",),
134 fgcmTransmissionAtmosphere = connectionTypes.Output(
135 doc=
"Per-visit atmosphere transmission files produced from fgcm calibration",
136 name=
"transmission_atmosphere_fgcm",
137 storageClass=
"TransmissionCurve",
138 dimensions=(
"instrument",
143 fgcmOffsets = connectionTypes.Output(
144 doc=
"Per-band offsets computed from doReferenceCalibration",
145 name=
"fgcmReferenceCalibrationOffsets",
146 storageClass=
"Catalog",
147 dimensions=(
"instrument",),
151 def __init__(self, *, config=None):
152 super().__init__(config=config)
154 if str(int(config.connections.cycleNumber)) != config.connections.cycleNumber:
155 raise ValueError(
"cycleNumber must be of integer format")
157 if not config.doReferenceCalibration:
158 self.prerequisiteInputs.remove(
"refCat")
159 if not config.doAtmosphereOutput:
160 self.inputs.remove(
"fgcmAtmosphereParameters")
161 if not config.doZeropointOutput:
162 self.inputs.remove(
"fgcmZeropoints")
163 if not config.doReferenceCalibration:
164 self.outputs.remove(
"fgcmOffsets")
167class FgcmOutputProductsConfig(pipeBase.PipelineTaskConfig,
168 pipelineConnections=FgcmOutputProductsConnections):
169 """Config for FgcmOutputProductsTask"""
171 cycleNumber = pexConfig.Field(
172 doc=
"Final fit cycle from FGCM fit",
175 deprecated=(
"This config is no longer used, and will be removed after v25. "
176 "Please set config.connections.cycleNumber directly instead."),
178 physicalFilterMap = pexConfig.DictField(
179 doc=
"Mapping from 'physicalFilter' to band.",
186 doReferenceCalibration = pexConfig.Field(
187 doc=(
"Transfer 'absolute' calibration from reference catalog? "
188 "This afterburner step is unnecessary if reference stars "
189 "were used in the full fit in FgcmFitCycleTask."),
193 doAtmosphereOutput = pexConfig.Field(
194 doc=
"Output atmospheres in transmission_atmosphere_fgcm format",
198 doZeropointOutput = pexConfig.Field(
199 doc=
"Output zeropoints in fgcm_photoCalib format",
203 doComposeWcsJacobian = pexConfig.Field(
204 doc=
"Compose Jacobian of WCS with fgcm calibration for output photoCalib?",
208 doApplyMeanChromaticCorrection = pexConfig.Field(
209 doc=
"Apply the mean chromatic correction to the zeropoints?",
213 photoCal = pexConfig.ConfigurableField(
215 doc=
"task to perform 'absolute' calibration",
217 referencePixelizationNside = pexConfig.Field(
218 doc=
"Healpix nside to pixelize catalog to compare to reference catalog",
222 referencePixelizationMinStars = pexConfig.Field(
223 doc=(
"Minimum number of stars per healpix pixel to select for comparison"
224 "to the specified reference catalog"),
228 referenceMinMatch = pexConfig.Field(
229 doc=
"Minimum number of stars matched to reference catalog to be used in statistics",
233 referencePixelizationNPixels = pexConfig.Field(
234 doc=(
"Number of healpix pixels to sample to do comparison. "
235 "Doing too many will take a long time and not yield any more "
236 "precise results because the final number is the median offset "
237 "(per band) from the set of pixels."),
242 def setDefaults(self):
243 pexConfig.Config.setDefaults(self)
253 self.photoCal.applyColorTerms =
False
254 self.photoCal.fluxField =
'instFlux'
255 self.photoCal.magErrFloor = 0.003
256 self.photoCal.match.referenceSelection.doSignalToNoise =
True
257 self.photoCal.match.referenceSelection.signalToNoise.minimum = 10.0
258 self.photoCal.match.sourceSelection.doSignalToNoise =
True
259 self.photoCal.match.sourceSelection.signalToNoise.minimum = 10.0
260 self.photoCal.match.sourceSelection.signalToNoise.fluxField =
'instFlux'
261 self.photoCal.match.sourceSelection.signalToNoise.errField =
'instFluxErr'
262 self.photoCal.match.sourceSelection.doFlags =
True
263 self.photoCal.match.sourceSelection.flags.good = []
264 self.photoCal.match.sourceSelection.flags.bad = [
'flag_badStar']
265 self.photoCal.match.sourceSelection.doUnresolved =
False
268class FgcmOutputProductsTask(pipeBase.PipelineTask):
270 Output products from FGCM
global calibration.
273 ConfigClass = FgcmOutputProductsConfig
274 _DefaultName = "fgcmOutputProducts"
276 def __init__(self, **kwargs):
277 super().__init__(**kwargs)
279 def runQuantum(self, butlerQC, inputRefs, outputRefs):
281 handleDict[
'camera'] = butlerQC.get(inputRefs.camera)
282 handleDict[
'fgcmLookUpTable'] = butlerQC.get(inputRefs.fgcmLookUpTable)
283 handleDict[
'fgcmVisitCatalog'] = butlerQC.get(inputRefs.fgcmVisitCatalog)
284 handleDict[
'fgcmStandardStars'] = butlerQC.get(inputRefs.fgcmStandardStars)
286 if self.config.doZeropointOutput:
287 handleDict[
'fgcmZeropoints'] = butlerQC.get(inputRefs.fgcmZeropoints)
288 photoCalibRefDict = {photoCalibRef.dataId.byName()[
'visit']:
289 photoCalibRef
for photoCalibRef
in outputRefs.fgcmPhotoCalib}
291 if self.config.doAtmosphereOutput:
292 handleDict[
'fgcmAtmosphereParameters'] = butlerQC.get(inputRefs.fgcmAtmosphereParameters)
293 atmRefDict = {atmRef.dataId.byName()[
'visit']: atmRef
for
294 atmRef
in outputRefs.fgcmTransmissionAtmosphere}
296 if self.config.doReferenceCalibration:
297 refConfig = LoadReferenceObjectsConfig()
298 self.refObjLoader = ReferenceObjectLoader(dataIds=[ref.datasetRef.dataId
299 for ref
in inputRefs.refCat],
300 refCats=butlerQC.get(inputRefs.refCat),
301 name=self.config.connections.refCat,
305 self.refObjLoader =
None
307 struct = self.run(handleDict, self.config.physicalFilterMap)
310 if struct.photoCalibCatalogs
is not None:
311 self.log.info(
"Outputting photoCalib catalogs.")
312 for visit, expCatalog
in struct.photoCalibCatalogs:
313 butlerQC.put(expCatalog, photoCalibRefDict[visit])
314 self.log.info(
"Done outputting photoCalib catalogs.")
317 if struct.atmospheres
is not None:
318 self.log.info(
"Outputting atmosphere transmission files.")
319 for visit, atm
in struct.atmospheres:
320 butlerQC.put(atm, atmRefDict[visit])
321 self.log.info(
"Done outputting atmosphere files.")
323 if self.config.doReferenceCalibration:
325 schema = afwTable.Schema()
326 schema.addField(
'offset', type=np.float64,
327 doc=
"Post-process calibration offset (mag)")
328 offsetCat = afwTable.BaseCatalog(schema)
329 offsetCat.resize(len(struct.offsets))
330 offsetCat[
'offset'][:] = struct.offsets
332 butlerQC.put(offsetCat, outputRefs.fgcmOffsets)
336 def run(self, handleDict, physicalFilterMap):
337 """Run the output products task.
342 All handles are `lsst.daf.butler.DeferredDatasetHandle`
343 handle dictionary with keys:
346 Camera object (`lsst.afw.cameraGeom.Camera`)
347 ``
"fgcmLookUpTable"``
348 handle
for the FGCM look-up table.
349 ``
"fgcmVisitCatalog"``
350 handle
for visit summary catalog.
351 ``
"fgcmStandardStars"``
352 handle
for the output standard star catalog.
354 handle
for the zeropoint data catalog.
355 ``
"fgcmAtmosphereParameters"``
356 handle
for the atmosphere parameter catalog.
357 ``
"fgcmBuildStarsTableConfig"``
358 Config
for `lsst.fgcmcal.fgcmBuildStarsTableTask`.
359 physicalFilterMap : `dict`
360 Dictionary of mappings
from physical filter to FGCM band.
364 retStruct : `lsst.pipe.base.Struct`
365 Output structure
with keys:
367 offsets : `np.ndarray`
368 Final reference offsets, per band.
369 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
370 Generator that returns (visit, transmissionCurve) tuples.
371 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
372 Generator that returns (visit, exposureCatalog) tuples.
374 stdCat = handleDict['fgcmStandardStars'].get()
375 md = stdCat.getMetadata()
376 bands = md.getArray(
'BANDS')
378 if self.config.doReferenceCalibration:
379 lutCat = handleDict[
'fgcmLookUpTable'].get()
380 offsets = self._computeReferenceOffsets(stdCat, lutCat, physicalFilterMap, bands)
382 offsets = np.zeros(len(bands))
386 if self.config.doZeropointOutput:
387 zptCat = handleDict[
'fgcmZeropoints'].get()
388 visitCat = handleDict[
'fgcmVisitCatalog'].get()
390 pcgen = self._outputZeropoints(handleDict[
'camera'], zptCat, visitCat, offsets, bands,
395 if self.config.doAtmosphereOutput:
396 atmCat = handleDict[
'fgcmAtmosphereParameters'].get()
397 atmgen = self._outputAtmospheres(handleDict, atmCat)
401 retStruct = pipeBase.Struct(offsets=offsets,
403 retStruct.photoCalibCatalogs = pcgen
407 def generateTractOutputProducts(self, handleDict, tract,
408 visitCat, zptCat, atmCat, stdCat,
409 fgcmBuildStarsConfig):
411 Generate the output products for a given tract,
as specified
in the config.
413 This method
is here to have an alternate entry-point
for
419 All handles are `lsst.daf.butler.DeferredDatasetHandle`
420 handle dictionary
with keys:
423 Camera object (`lsst.afw.cameraGeom.Camera`)
424 ``
"fgcmLookUpTable"``
425 handle
for the FGCM look-up table.
428 visitCat : `lsst.afw.table.BaseCatalog`
429 FGCM visitCat
from `FgcmBuildStarsTask`
430 zptCat : `lsst.afw.table.BaseCatalog`
431 FGCM zeropoint catalog
from `FgcmFitCycleTask`
432 atmCat : `lsst.afw.table.BaseCatalog`
433 FGCM atmosphere parameter catalog
from `FgcmFitCycleTask`
434 stdCat : `lsst.afw.table.SimpleCatalog`
435 FGCM standard star catalog
from `FgcmFitCycleTask`
436 fgcmBuildStarsConfig : `lsst.fgcmcal.FgcmBuildStarsConfig`
437 Configuration object
from `FgcmBuildStarsTask`
441 retStruct : `lsst.pipe.base.Struct`
442 Output structure
with keys:
444 offsets : `np.ndarray`
445 Final reference offsets, per band.
446 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
447 Generator that returns (visit, transmissionCurve) tuples.
448 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
449 Generator that returns (visit, exposureCatalog) tuples.
451 physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap
453 md = stdCat.getMetadata()
454 bands = md.getArray('BANDS')
456 if self.config.doComposeWcsJacobian
and not fgcmBuildStarsConfig.doApplyWcsJacobian:
457 raise RuntimeError(
"Cannot compose the WCS jacobian if it hasn't been applied "
458 "in fgcmBuildStarsTask.")
460 if not self.config.doComposeWcsJacobian
and fgcmBuildStarsConfig.doApplyWcsJacobian:
461 self.log.warning(
"Jacobian was applied in build-stars but doComposeWcsJacobian is not set.")
463 if self.config.doReferenceCalibration:
464 lutCat = handleDict[
'fgcmLookUpTable'].get()
465 offsets = self._computeReferenceOffsets(stdCat, lutCat, bands, physicalFilterMap)
467 offsets = np.zeros(len(bands))
469 if self.config.doZeropointOutput:
470 pcgen = self._outputZeropoints(handleDict[
'camera'], zptCat, visitCat, offsets, bands,
475 if self.config.doAtmosphereOutput:
476 atmgen = self._outputAtmospheres(handleDict, atmCat)
480 retStruct = pipeBase.Struct(offsets=offsets,
482 retStruct.photoCalibCatalogs = pcgen
486 def _computeReferenceOffsets(self, stdCat, lutCat, physicalFilterMap, bands):
488 Compute offsets relative to a reference catalog.
490 This method splits the star catalog into healpix pixels
491 and computes the calibration transfer
for a sample of
492 these pixels to approximate the
'absolute' calibration
493 values (on
for each band) to apply to transfer the
498 stdCat : `lsst.afw.table.SimpleCatalog`
500 lutCat : `lsst.afw.table.SimpleCatalog`
502 physicalFilterMap : `dict`
503 Dictionary of mappings
from physical filter to FGCM band.
504 bands : `list` [`str`]
505 List of band names
from FGCM output
508 offsets : `numpy.array` of floats
509 Per band zeropoint offsets
515 minObs = stdCat[
'ngood'].min(axis=1)
517 goodStars = (minObs >= 1)
518 stdCat = stdCat[goodStars]
520 self.log.info(
"Found %d stars with at least 1 good observation in each band" %
527 lutPhysicalFilters = lutCat[0][
'physicalFilters'].split(
',')
528 lutStdPhysicalFilters = lutCat[0][
'stdPhysicalFilters'].split(
',')
529 physicalFilterMapBands = list(physicalFilterMap.values())
530 physicalFilterMapFilters = list(physicalFilterMap.keys())
534 physicalFilterMapIndex = physicalFilterMapBands.index(band)
535 physicalFilter = physicalFilterMapFilters[physicalFilterMapIndex]
537 lutPhysicalFilterIndex = lutPhysicalFilters.index(physicalFilter)
538 stdPhysicalFilter = lutStdPhysicalFilters[lutPhysicalFilterIndex]
539 filterLabels.append(afwImage.FilterLabel(band=band,
540 physical=stdPhysicalFilter))
549 sourceMapper = afwTable.SchemaMapper(stdCat.schema)
550 sourceMapper.addMinimalSchema(afwTable.SimpleTable.makeMinimalSchema())
551 sourceMapper.editOutputSchema().addField(
'instFlux', type=np.float64,
552 doc=
"instrumental flux (counts)")
553 sourceMapper.editOutputSchema().addField(
'instFluxErr', type=np.float64,
554 doc=
"instrumental flux error (counts)")
555 badStarKey = sourceMapper.editOutputSchema().addField(
'flag_badStar',
563 ipring = hpg.angle_to_pixel(
564 self.config.referencePixelizationNside,
569 h, rev = esutil.stat.histogram(ipring, rev=
True)
571 gdpix, = np.where(h >= self.config.referencePixelizationMinStars)
573 self.log.info(
"Found %d pixels (nside=%d) with at least %d good stars" %
575 self.config.referencePixelizationNside,
576 self.config.referencePixelizationMinStars))
578 if gdpix.size < self.config.referencePixelizationNPixels:
579 self.log.warning(
"Found fewer good pixels (%d) than preferred in configuration (%d)" %
580 (gdpix.size, self.config.referencePixelizationNPixels))
583 gdpix = np.random.choice(gdpix, size=self.config.referencePixelizationNPixels, replace=
False)
585 results = np.zeros(gdpix.size, dtype=[(
'hpix',
'i4'),
586 (
'nstar',
'i4', len(bands)),
587 (
'nmatch',
'i4', len(bands)),
588 (
'zp',
'f4', len(bands)),
589 (
'zpErr',
'f4', len(bands))])
590 results[
'hpix'] = ipring[rev[rev[gdpix]]]
593 selected = np.zeros(len(stdCat), dtype=bool)
595 refFluxFields = [
None]*len(bands)
597 for p_index, pix
in enumerate(gdpix):
598 i1a = rev[rev[pix]: rev[pix + 1]]
606 for b_index, filterLabel
in enumerate(filterLabels):
607 struct = self._computeOffsetOneBand(sourceMapper, badStarKey, b_index,
609 selected, refFluxFields)
610 results[
'nstar'][p_index, b_index] = len(i1a)
611 results[
'nmatch'][p_index, b_index] = len(struct.arrays.refMag)
612 results[
'zp'][p_index, b_index] = struct.zp
613 results[
'zpErr'][p_index, b_index] = struct.sigma
616 offsets = np.zeros(len(bands))
618 for b_index, band
in enumerate(bands):
620 ok, = np.where(results[
'nmatch'][:, b_index] >= self.config.referenceMinMatch)
621 offsets[b_index] = np.median(results[
'zp'][ok, b_index])
624 madSigma = 1.4826*np.median(np.abs(results[
'zp'][ok, b_index] - offsets[b_index]))
625 self.log.info(
"Reference catalog offset for %s band: %.12f +/- %.12f",
626 band, offsets[b_index], madSigma)
630 def _computeOffsetOneBand(self, sourceMapper, badStarKey,
631 b_index, filterLabel, stdCat, selected, refFluxFields):
633 Compute the zeropoint offset between the fgcm stdCat and the reference
634 stars
for one pixel
in one band
638 sourceMapper : `lsst.afw.table.SchemaMapper`
639 Mapper to go
from stdCat to calibratable catalog
640 badStarKey : `lsst.afw.table.Key`
641 Key
for the field
with bad stars
643 Index of the band
in the star catalog
644 filterLabel : `lsst.afw.image.FilterLabel`
645 filterLabel
with band
and physical filter
646 stdCat : `lsst.afw.table.SimpleCatalog`
648 selected : `numpy.array(dtype=bool)`
649 Boolean array of which stars are
in the pixel
650 refFluxFields : `list`
651 List of names of flux fields
for reference catalog
654 sourceCat = afwTable.SimpleCatalog(sourceMapper.getOutputSchema())
655 sourceCat.reserve(selected.sum())
656 sourceCat.extend(stdCat[selected], mapper=sourceMapper)
657 sourceCat['instFlux'] = 10.**(stdCat[
'mag_std_noabs'][selected, b_index]/(-2.5))
658 sourceCat[
'instFluxErr'] = (np.log(10.)/2.5)*(stdCat[
'magErr_std'][selected, b_index]
659 * sourceCat[
'instFlux'])
663 badStar = (stdCat[
'mag_std_noabs'][selected, b_index] > 90.0)
664 for rec
in sourceCat[badStar]:
665 rec.set(badStarKey,
True)
667 exposure = afwImage.ExposureF()
668 exposure.setFilter(filterLabel)
670 if refFluxFields[b_index]
is None:
673 ctr = stdCat[0].getCoord()
674 rad = 0.05*lsst.geom.degrees
675 refDataTest = self.refObjLoader.loadSkyCircle(ctr, rad, filterLabel.bandLabel)
676 refFluxFields[b_index] = refDataTest.fluxField
679 calConfig = copy.copy(self.config.photoCal.value)
680 calConfig.match.referenceSelection.signalToNoise.fluxField = refFluxFields[b_index]
681 calConfig.match.referenceSelection.signalToNoise.errField = refFluxFields[b_index] +
'Err'
682 calTask = self.config.photoCal.target(refObjLoader=self.refObjLoader,
684 schema=sourceCat.getSchema())
686 struct = calTask.run(exposure, sourceCat)
690 def _outputZeropoints(self, camera, zptCat, visitCat, offsets, bands,
691 physicalFilterMap, tract=None):
692 """Output the zeropoints in fgcm_photoCalib format.
696 camera : `lsst.afw.cameraGeom.Camera`
697 Camera from the butler.
698 zptCat : `lsst.afw.table.BaseCatalog`
699 FGCM zeropoint catalog
from `FgcmFitCycleTask`.
700 visitCat : `lsst.afw.table.BaseCatalog`
701 FGCM visitCat
from `FgcmBuildStarsTask`.
702 offsets : `numpy.array`
703 Float array of absolute calibration offsets, one
for each filter.
704 bands : `list` [`str`]
705 List of band names
from FGCM output.
706 physicalFilterMap : `dict`
707 Dictionary of mappings
from physical filter to FGCM band.
708 tract: `int`, optional
709 Tract number to output. Default
is None (
global calibration)
713 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
714 Generator that returns (visit, exposureCatalog) tuples.
719 cannot_compute = fgcm.fgcmUtilities.zpFlagDict[
'CANNOT_COMPUTE_ZEROPOINT']
720 selected = (((zptCat[
'fgcmFlag'] & cannot_compute) == 0)
721 & (zptCat[
'fgcmZptVar'] > 0.0)
722 & (zptCat[
'fgcmZpt'] > FGCM_ILLEGAL_VALUE))
725 badVisits = np.unique(zptCat[
'visit'][~selected])
726 goodVisits = np.unique(zptCat[
'visit'][selected])
727 allBadVisits = badVisits[~np.isin(badVisits, goodVisits)]
728 for allBadVisit
in allBadVisits:
729 self.log.warning(f
'No suitable photoCalib for visit {allBadVisit}')
733 for f
in physicalFilterMap:
735 if physicalFilterMap[f]
in bands:
736 offsetMapping[f] = offsets[bands.index(physicalFilterMap[f])]
740 for ccdIndex, detector
in enumerate(camera):
741 ccdMapping[detector.getId()] = ccdIndex
746 scalingMapping[rec[
'visit']] = rec[
'scaling']
748 if self.config.doComposeWcsJacobian:
749 approxPixelAreaFields = computeApproxPixelAreaFields(camera)
753 zptVisitCatalog =
None
755 metadata = dafBase.PropertyList()
756 metadata.add(
"COMMENT",
"Catalog id is detector id, sorted.")
757 metadata.add(
"COMMENT",
"Only detectors with data have entries.")
759 for rec
in zptCat[selected]:
761 scaling = scalingMapping[rec[
'visit']][ccdMapping[rec[
'detector']]]
768 postCalibrationOffset = offsetMapping[rec[
'filtername']]
769 if self.config.doApplyMeanChromaticCorrection:
770 postCalibrationOffset += rec[
'fgcmDeltaChrom']
772 fgcmSuperStarField = self._getChebyshevBoundedField(rec[
'fgcmfZptSstarCheb'],
773 rec[
'fgcmfZptChebXyMax'])
775 fgcmZptField = self._getChebyshevBoundedField((rec[
'fgcmfZptCheb']*units.AB).to_value(units.nJy),
776 rec[
'fgcmfZptChebXyMax'],
777 offset=postCalibrationOffset,
780 if self.config.doComposeWcsJacobian:
782 fgcmField = afwMath.ProductBoundedField([approxPixelAreaFields[rec[
'detector']],
788 fgcmField = afwMath.ProductBoundedField([fgcmSuperStarField, fgcmZptField])
791 calibCenter = fgcmField.evaluate(fgcmField.getBBox().getCenter())
792 calibErr = (np.log(10.0)/2.5)*calibCenter*np.sqrt(rec[
'fgcmZptVar'])
793 photoCalib = afwImage.PhotoCalib(calibrationMean=calibCenter,
794 calibrationErr=calibErr,
795 calibration=fgcmField,
799 if rec[
'visit'] != lastVisit:
804 zptVisitCatalog.sort()
805 yield (int(lastVisit), zptVisitCatalog)
808 zptExpCatSchema = afwTable.ExposureTable.makeMinimalSchema()
809 zptExpCatSchema.addField(
'visit', type=
'L', doc=
'Visit number')
812 zptVisitCatalog = afwTable.ExposureCatalog(zptExpCatSchema)
813 zptVisitCatalog.setMetadata(metadata)
815 lastVisit = int(rec[
'visit'])
817 catRecord = zptVisitCatalog.addNew()
818 catRecord[
'id'] = int(rec[
'detector'])
819 catRecord[
'visit'] = rec[
'visit']
820 catRecord.setPhotoCalib(photoCalib)
824 zptVisitCatalog.sort()
825 yield (int(lastVisit), zptVisitCatalog)
827 def _getChebyshevBoundedField(self, coefficients, xyMax, offset=0.0, scaling=1.0):
829 Make a ChebyshevBoundedField from fgcm coefficients,
with optional offset
834 coefficients: `numpy.array`
835 Flattened array of chebyshev coefficients
836 xyMax: `list` of length 2
837 Maximum x
and y of the chebyshev bounding box
838 offset: `float`, optional
839 Absolute calibration offset. Default
is 0.0
840 scaling: `float`, optional
841 Flat scaling value
from fgcmBuildStars. Default
is 1.0
845 boundedField: `lsst.afw.math.ChebyshevBoundedField`
848 orderPlus1 = int(np.sqrt(coefficients.size))
849 pars = np.zeros((orderPlus1, orderPlus1))
851 bbox = lsst.geom.Box2I(lsst.geom.Point2I(0.0, 0.0),
852 lsst.geom.Point2I(*xyMax))
854 pars[:, :] = (coefficients.reshape(orderPlus1, orderPlus1)
855 * (10.**(offset/-2.5))*scaling)
857 boundedField = afwMath.ChebyshevBoundedField(bbox, pars)
861 def _outputAtmospheres(self, handleDict, atmCat):
863 Output the atmospheres.
868 All data handles are `lsst.daf.butler.DeferredDatasetHandle`
869 The handleDict has the follownig keys:
871 ``"fgcmLookUpTable"``
872 handle
for the FGCM look-up table.
873 atmCat : `lsst.afw.table.BaseCatalog`
874 FGCM atmosphere parameter catalog
from fgcmFitCycleTask.
878 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
879 Generator that returns (visit, transmissionCurve) tuples.
882 lutCat = handleDict[
'fgcmLookUpTable'].get()
884 atmosphereTableName = lutCat[0][
'tablename']
885 elevation = lutCat[0][
'elevation']
886 atmLambda = lutCat[0][
'atmLambda']
891 atmTable = fgcm.FgcmAtmosphereTable.initWithTableName(atmosphereTableName)
899 modGen = fgcm.ModtranGenerator(elevation)
900 lambdaRange = np.array([atmLambda[0], atmLambda[-1]])/10.
901 lambdaStep = (atmLambda[1] - atmLambda[0])/10.
902 except (ValueError, IOError)
as e:
903 raise RuntimeError(
"FGCM look-up-table generated with modtran, "
904 "but modtran not configured to run.")
from e
906 zenith = np.degrees(np.arccos(1./atmCat[
'secZenith']))
908 for i, visit
in enumerate(atmCat[
'visit']):
909 if atmTable
is not None:
911 atmVals = atmTable.interpolateAtmosphere(pmb=atmCat[i][
'pmb'],
912 pwv=atmCat[i][
'pwv'],
914 tau=atmCat[i][
'tau'],
915 alpha=atmCat[i][
'alpha'],
917 ctranslamstd=[atmCat[i][
'cTrans'],
918 atmCat[i][
'lamStd']])
921 modAtm = modGen(pmb=atmCat[i][
'pmb'],
922 pwv=atmCat[i][
'pwv'],
924 tau=atmCat[i][
'tau'],
925 alpha=atmCat[i][
'alpha'],
927 lambdaRange=lambdaRange,
928 lambdaStep=lambdaStep,
929 ctranslamstd=[atmCat[i][
'cTrans'],
930 atmCat[i][
'lamStd']])
931 atmVals = modAtm[
'COMBINED']
934 curve = TransmissionCurve.makeSpatiallyConstant(throughput=atmVals,
935 wavelengths=atmLambda,
936 throughputAtMin=atmVals[0],
937 throughputAtMax=atmVals[-1])
939 yield (int(visit), curve)