lsst.fgcmcal g7d743f24c3+f67ed74f77
Loading...
Searching...
No Matches
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 copy
35
36import numpy as np
37import hpgeom as hpg
38import esutil
39from astropy import units
40
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
48import lsst.geom
49import lsst.afw.image as afwImage
50import lsst.afw.math as afwMath
51import lsst.afw.table as afwTable
52
53from .utilities import computeApproxPixelAreaFields
54from .utilities import FGCM_ILLEGAL_VALUE
55
56import fgcm
57
58__all__ = ['FgcmOutputProductsConfig', 'FgcmOutputProductsTask']
59
60
61class FgcmOutputProductsConnections(pipeBase.PipelineTaskConnections,
62 dimensions=("instrument",),
63 defaultTemplates={"cycleNumber": "0"}):
64 camera = connectionTypes.PrerequisiteInput(
65 doc="Camera instrument",
66 name="camera",
67 storageClass="Camera",
68 dimensions=("instrument",),
69 isCalibration=True,
70 )
71
72 fgcmLookUpTable = connectionTypes.PrerequisiteInput(
73 doc=("Atmosphere + instrument look-up-table for FGCM throughput and "
74 "chromatic corrections."),
75 name="fgcmLookUpTable",
76 storageClass="Catalog",
77 dimensions=("instrument",),
78 deferLoad=True,
79 )
80
81 fgcmVisitCatalog = connectionTypes.Input(
82 doc="Catalog of visit information for fgcm",
83 name="fgcmVisitCatalog",
84 storageClass="Catalog",
85 dimensions=("instrument",),
86 deferLoad=True,
87 )
88
89 fgcmStandardStars = connectionTypes.Input(
90 doc="Catalog of standard star data from fgcm fit",
91 name="fgcmStandardStars{cycleNumber}",
92 storageClass="SimpleCatalog",
93 dimensions=("instrument",),
94 deferLoad=True,
95 )
96
97 fgcmZeropoints = connectionTypes.Input(
98 doc="Catalog of zeropoints from fgcm fit",
99 name="fgcmZeropoints{cycleNumber}",
100 storageClass="Catalog",
101 dimensions=("instrument",),
102 deferLoad=True,
103 )
104
105 fgcmAtmosphereParameters = connectionTypes.Input(
106 doc="Catalog of atmosphere parameters from fgcm fit",
107 name="fgcmAtmosphereParameters{cycleNumber}",
108 storageClass="Catalog",
109 dimensions=("instrument",),
110 deferLoad=True,
111 )
112
113 refCat = connectionTypes.PrerequisiteInput(
114 doc="Reference catalog to use for photometric calibration",
115 name="cal_ref_cat",
116 storageClass="SimpleCatalog",
117 dimensions=("skypix",),
118 deferLoad=True,
119 multiple=True,
120 )
121
122 fgcmPhotoCalib = connectionTypes.Output(
123 doc=("Per-visit photometric calibrations derived from fgcm calibration. "
124 "These catalogs use detector id for the id and are sorted for "
125 "fast lookups of a detector."),
126 name="fgcmPhotoCalibCatalog",
127 storageClass="ExposureCatalog",
128 dimensions=("instrument", "visit",),
129 multiple=True,
130 )
131
132 fgcmTransmissionAtmosphere = connectionTypes.Output(
133 doc="Per-visit atmosphere transmission files produced from fgcm calibration",
134 name="transmission_atmosphere_fgcm",
135 storageClass="TransmissionCurve",
136 dimensions=("instrument",
137 "visit",),
138 multiple=True,
139 )
140
141 fgcmOffsets = connectionTypes.Output(
142 doc="Per-band offsets computed from doReferenceCalibration",
143 name="fgcmReferenceCalibrationOffsets",
144 storageClass="Catalog",
145 dimensions=("instrument",),
146 multiple=False,
147 )
148
149 def __init__(self, *, config=None):
150 super().__init__(config=config)
151
152 if str(int(config.connections.cycleNumber)) != config.connections.cycleNumber:
153 raise ValueError("cycleNumber must be of integer format")
154
155 if not config.doReferenceCalibration:
156 self.prerequisiteInputs.remove("refCat")
157 if not config.doAtmosphereOutput:
158 self.inputs.remove("fgcmAtmosphereParameters")
159 if not config.doZeropointOutput:
160 self.inputs.remove("fgcmZeropoints")
161 if not config.doReferenceCalibration:
162 self.outputs.remove("fgcmOffsets")
163
164 def getSpatialBoundsConnections(self):
165 return ("fgcmPhotoCalib",)
166
167
168class FgcmOutputProductsConfig(pipeBase.PipelineTaskConfig,
169 pipelineConnections=FgcmOutputProductsConnections):
170 """Config for FgcmOutputProductsTask"""
171
172 physicalFilterMap = pexConfig.DictField(
173 doc="Mapping from 'physicalFilter' to band.",
174 keytype=str,
175 itemtype=str,
176 default={},
177 )
178 # The following fields refer to calibrating from a reference
179 # catalog, but in the future this might need to be expanded
180 doReferenceCalibration = pexConfig.Field(
181 doc=("Transfer 'absolute' calibration from reference catalog? "
182 "This afterburner step is unnecessary if reference stars "
183 "were used in the full fit in FgcmFitCycleTask."),
184 dtype=bool,
185 default=False,
186 )
187 doAtmosphereOutput = pexConfig.Field(
188 doc="Output atmospheres in transmission_atmosphere_fgcm format",
189 dtype=bool,
190 default=True,
191 )
192 doZeropointOutput = pexConfig.Field(
193 doc="Output zeropoints in fgcm_photoCalib format",
194 dtype=bool,
195 default=True,
196 )
197 doComposeWcsJacobian = pexConfig.Field(
198 doc="Compose Jacobian of WCS with fgcm calibration for output photoCalib?",
199 dtype=bool,
200 default=True,
201 )
202 doApplyMeanChromaticCorrection = pexConfig.Field(
203 doc="Apply the mean chromatic correction to the zeropoints?",
204 dtype=bool,
205 default=True,
206 )
207 photoCal = pexConfig.ConfigurableField(
208 target=PhotoCalTask,
209 doc="task to perform 'absolute' calibration",
210 )
211 referencePixelizationNside = pexConfig.Field(
212 doc="Healpix nside to pixelize catalog to compare to reference catalog",
213 dtype=int,
214 default=64,
215 )
216 referencePixelizationMinStars = pexConfig.Field(
217 doc=("Minimum number of stars per healpix pixel to select for comparison"
218 "to the specified reference catalog"),
219 dtype=int,
220 default=200,
221 )
222 referenceMinMatch = pexConfig.Field(
223 doc="Minimum number of stars matched to reference catalog to be used in statistics",
224 dtype=int,
225 default=50,
226 )
227 referencePixelizationNPixels = pexConfig.Field(
228 doc=("Number of healpix pixels to sample to do comparison. "
229 "Doing too many will take a long time and not yield any more "
230 "precise results because the final number is the median offset "
231 "(per band) from the set of pixels."),
232 dtype=int,
233 default=100,
234 )
235
236 def setDefaults(self):
237 pexConfig.Config.setDefaults(self)
238
239 # In order to transfer the "absolute" calibration from a reference
240 # catalog to the relatively calibrated FGCM standard stars (one number
241 # per band), we use the PhotoCalTask to match stars in a sample of healpix
242 # pixels. These basic settings ensure that only well-measured, good stars
243 # from the source and reference catalogs are used for the matching.
244
245 # applyColorTerms needs to be False if doReferenceCalibration is False,
246 # as is the new default after DM-16702
247 self.photoCal.applyColorTerms = False
248 self.photoCal.fluxField = 'instFlux'
249 self.photoCal.magErrFloor = 0.003
250 self.photoCal.match.referenceSelection.doSignalToNoise = True
251 self.photoCal.match.referenceSelection.signalToNoise.minimum = 10.0
252 self.photoCal.match.sourceSelection.doSignalToNoise = True
253 self.photoCal.match.sourceSelection.signalToNoise.minimum = 10.0
254 self.photoCal.match.sourceSelection.signalToNoise.fluxField = 'instFlux'
255 self.photoCal.match.sourceSelection.signalToNoise.errField = 'instFluxErr'
256 self.photoCal.match.sourceSelection.doFlags = True
257 self.photoCal.match.sourceSelection.flags.good = []
258 self.photoCal.match.sourceSelection.flags.bad = ['flag_badStar']
259 self.photoCal.match.sourceSelection.doUnresolved = False
260 self.photoCal.match.sourceSelection.doRequirePrimary = False
261
262
263class FgcmOutputProductsTask(pipeBase.PipelineTask):
264 """
265 Output products from FGCM global calibration.
266 """
267
268 ConfigClass = FgcmOutputProductsConfig
269 _DefaultName = "fgcmOutputProducts"
270
271 def __init__(self, **kwargs):
272 super().__init__(**kwargs)
273
274 def runQuantum(self, butlerQC, inputRefs, outputRefs):
275 handleDict = {}
276 handleDict['camera'] = butlerQC.get(inputRefs.camera)
277 handleDict['fgcmLookUpTable'] = butlerQC.get(inputRefs.fgcmLookUpTable)
278 handleDict['fgcmVisitCatalog'] = butlerQC.get(inputRefs.fgcmVisitCatalog)
279 handleDict['fgcmStandardStars'] = butlerQC.get(inputRefs.fgcmStandardStars)
280
281 if self.config.doZeropointOutput:
282 handleDict['fgcmZeropoints'] = butlerQC.get(inputRefs.fgcmZeropoints)
283 photoCalibRefDict = {photoCalibRef.dataId['visit']:
284 photoCalibRef for photoCalibRef in outputRefs.fgcmPhotoCalib}
285
286 if self.config.doAtmosphereOutput:
287 handleDict['fgcmAtmosphereParameters'] = butlerQC.get(inputRefs.fgcmAtmosphereParameters)
288 atmRefDict = {atmRef.dataId['visit']: atmRef for
289 atmRef in outputRefs.fgcmTransmissionAtmosphere}
290
291 if self.config.doReferenceCalibration:
292 refConfig = LoadReferenceObjectsConfig()
293 self.refObjLoader = ReferenceObjectLoader(dataIds=[ref.datasetRef.dataId
294 for ref in inputRefs.refCat],
295 refCats=butlerQC.get(inputRefs.refCat),
296 name=self.config.connections.refCat,
297 log=self.log,
298 config=refConfig)
299 else:
300 self.refObjLoader = None
301
302 struct = self.run(handleDict, self.config.physicalFilterMap)
303
304 # Output the photoCalib exposure catalogs
305 if struct.photoCalibCatalogs is not None:
306 self.log.info("Outputting photoCalib catalogs.")
307 for visit, expCatalog in struct.photoCalibCatalogs:
308 butlerQC.put(expCatalog, photoCalibRefDict[visit])
309 self.log.info("Done outputting photoCalib catalogs.")
310
311 # Output the atmospheres
312 if struct.atmospheres is not None:
313 self.log.info("Outputting atmosphere transmission files.")
314 for visit, atm in struct.atmospheres:
315 butlerQC.put(atm, atmRefDict[visit])
316 self.log.info("Done outputting atmosphere files.")
317
318 if self.config.doReferenceCalibration:
319 # Turn offset into simple catalog for persistence if necessary
320 schema = afwTable.Schema()
321 schema.addField('offset', type=np.float64,
322 doc="Post-process calibration offset (mag)")
323 offsetCat = afwTable.BaseCatalog(schema)
324 offsetCat.resize(len(struct.offsets))
325 offsetCat['offset'][:] = struct.offsets
326
327 butlerQC.put(offsetCat, outputRefs.fgcmOffsets)
328
329 return
330
331 def run(self, handleDict, physicalFilterMap):
332 """Run the output products task.
333
334 Parameters
335 ----------
336 handleDict : `dict`
337 All handles are `lsst.daf.butler.DeferredDatasetHandle`
338 handle dictionary with keys:
339
340 ``"camera"``
341 Camera object (`lsst.afw.cameraGeom.Camera`)
342 ``"fgcmLookUpTable"``
343 handle for the FGCM look-up table.
344 ``"fgcmVisitCatalog"``
345 handle for visit summary catalog.
346 ``"fgcmStandardStars"``
347 handle for the output standard star catalog.
348 ``"fgcmZeropoints"``
349 handle for the zeropoint data catalog.
350 ``"fgcmAtmosphereParameters"``
351 handle for the atmosphere parameter catalog.
352 ``"fgcmBuildStarsTableConfig"``
353 Config for `lsst.fgcmcal.fgcmBuildStarsTableTask`.
354 physicalFilterMap : `dict`
355 Dictionary of mappings from physical filter to FGCM band.
356
357 Returns
358 -------
359 retStruct : `lsst.pipe.base.Struct`
360 Output structure with keys:
361
362 offsets : `np.ndarray`
363 Final reference offsets, per band.
364 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
365 Generator that returns (visit, transmissionCurve) tuples.
366 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
367 Generator that returns (visit, exposureCatalog) tuples.
368 """
369 stdCat = handleDict['fgcmStandardStars'].get()
370 md = stdCat.getMetadata()
371 bands = md.getArray('BANDS')
372
373 if self.config.doReferenceCalibration:
374 lutCat = handleDict['fgcmLookUpTable'].get()
375 offsets = self._computeReferenceOffsets(stdCat, lutCat, physicalFilterMap, bands)
376 else:
377 offsets = np.zeros(len(bands))
378
379 del stdCat
380
381 if self.config.doZeropointOutput:
382 zptCat = handleDict['fgcmZeropoints'].get()
383 visitCat = handleDict['fgcmVisitCatalog'].get()
384
385 pcgen = self._outputZeropoints(handleDict['camera'], zptCat, visitCat, offsets, bands,
386 physicalFilterMap)
387 else:
388 pcgen = None
389
390 if self.config.doAtmosphereOutput:
391 atmCat = handleDict['fgcmAtmosphereParameters'].get()
392 atmgen = self._outputAtmospheres(handleDict, atmCat)
393 else:
394 atmgen = None
395
396 retStruct = pipeBase.Struct(offsets=offsets,
397 atmospheres=atmgen)
398 retStruct.photoCalibCatalogs = pcgen
399
400 return retStruct
401
402 def generateTractOutputProducts(self, handleDict, tract,
403 visitCat, zptCat, atmCat, stdCat,
404 fgcmBuildStarsConfig):
405 """
406 Generate the output products for a given tract, as specified in the config.
407
408 This method is here to have an alternate entry-point for
409 FgcmCalibrateTract.
410
411 Parameters
412 ----------
413 handleDict : `dict`
414 All handles are `lsst.daf.butler.DeferredDatasetHandle`
415 handle dictionary with keys:
416
417 ``"camera"``
418 Camera object (`lsst.afw.cameraGeom.Camera`)
419 ``"fgcmLookUpTable"``
420 handle for the FGCM look-up table.
421 tract : `int`
422 Tract number
423 visitCat : `lsst.afw.table.BaseCatalog`
424 FGCM visitCat from `FgcmBuildStarsTask`
425 zptCat : `lsst.afw.table.BaseCatalog`
426 FGCM zeropoint catalog from `FgcmFitCycleTask`
427 atmCat : `lsst.afw.table.BaseCatalog`
428 FGCM atmosphere parameter catalog from `FgcmFitCycleTask`
429 stdCat : `lsst.afw.table.SimpleCatalog`
430 FGCM standard star catalog from `FgcmFitCycleTask`
431 fgcmBuildStarsConfig : `lsst.fgcmcal.FgcmBuildStarsConfig`
432 Configuration object from `FgcmBuildStarsTask`
433
434 Returns
435 -------
436 retStruct : `lsst.pipe.base.Struct`
437 Output structure with keys:
438
439 offsets : `np.ndarray`
440 Final reference offsets, per band.
441 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
442 Generator that returns (visit, transmissionCurve) tuples.
443 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
444 Generator that returns (visit, exposureCatalog) tuples.
445 """
446 physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap
447
448 md = stdCat.getMetadata()
449 bands = md.getArray('BANDS')
450
451 if self.config.doComposeWcsJacobian and not fgcmBuildStarsConfig.doApplyWcsJacobian:
452 raise RuntimeError("Cannot compose the WCS jacobian if it hasn't been applied "
453 "in fgcmBuildStarsTask.")
454
455 if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian:
456 self.log.warning("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.")
457
458 if self.config.doReferenceCalibration:
459 lutCat = handleDict['fgcmLookUpTable'].get()
460 offsets = self._computeReferenceOffsets(stdCat, lutCat, bands, physicalFilterMap)
461 else:
462 offsets = np.zeros(len(bands))
463
464 if self.config.doZeropointOutput:
465 pcgen = self._outputZeropoints(handleDict['camera'], zptCat, visitCat, offsets, bands,
466 physicalFilterMap)
467 else:
468 pcgen = None
469
470 if self.config.doAtmosphereOutput:
471 atmgen = self._outputAtmospheres(handleDict, atmCat)
472 else:
473 atmgen = None
474
475 retStruct = pipeBase.Struct(offsets=offsets,
476 atmospheres=atmgen)
477 retStruct.photoCalibCatalogs = pcgen
478
479 return retStruct
480
481 def _computeReferenceOffsets(self, stdCat, lutCat, physicalFilterMap, bands):
482 """
483 Compute offsets relative to a reference catalog.
484
485 This method splits the star catalog into healpix pixels
486 and computes the calibration transfer for a sample of
487 these pixels to approximate the 'absolute' calibration
488 values (on for each band) to apply to transfer the
489 absolute scale.
490
491 Parameters
492 ----------
493 stdCat : `lsst.afw.table.SimpleCatalog`
494 FGCM standard stars
495 lutCat : `lsst.afw.table.SimpleCatalog`
496 FGCM Look-up table
497 physicalFilterMap : `dict`
498 Dictionary of mappings from physical filter to FGCM band.
499 bands : `list` [`str`]
500 List of band names from FGCM output
501 Returns
502 -------
503 offsets : `numpy.array` of floats
504 Per band zeropoint offsets
505 """
506
507 # Only use stars that are observed in all the bands that were actually used
508 # This will ensure that we use the same healpix pixels for the absolute
509 # calibration of each band.
510 minObs = stdCat['ngood'].min(axis=1)
511
512 goodStars = (minObs >= 1)
513 stdCat = stdCat[goodStars]
514
515 self.log.info("Found %d stars with at least 1 good observation in each band" %
516 (len(stdCat)))
517
518 # Associate each band with the appropriate physicalFilter and make
519 # filterLabels
520 filterLabels = []
521
522 lutPhysicalFilters = lutCat[0]['physicalFilters'].split(',')
523 lutStdPhysicalFilters = lutCat[0]['stdPhysicalFilters'].split(',')
524 physicalFilterMapBands = list(physicalFilterMap.values())
525 physicalFilterMapFilters = list(physicalFilterMap.keys())
526 for band in bands:
527 # Find a physical filter associated from the band by doing
528 # a reverse lookup on the physicalFilterMap dict
529 physicalFilterMapIndex = physicalFilterMapBands.index(band)
530 physicalFilter = physicalFilterMapFilters[physicalFilterMapIndex]
531 # Find the appropriate fgcm standard physicalFilter
532 lutPhysicalFilterIndex = lutPhysicalFilters.index(physicalFilter)
533 stdPhysicalFilter = lutStdPhysicalFilters[lutPhysicalFilterIndex]
534 filterLabels.append(afwImage.FilterLabel(band=band,
535 physical=stdPhysicalFilter))
536
537 # We have to make a table for each pixel with flux/fluxErr
538 # This is a temporary table generated for input to the photoCal task.
539 # These fluxes are not instFlux (they are top-of-the-atmosphere approximate and
540 # have had chromatic corrections applied to get to the standard system
541 # specified by the atmosphere/instrumental parameters), nor are they
542 # in Jansky (since they don't have a proper absolute calibration: the overall
543 # zeropoint is estimated from the telescope size, etc.)
544 sourceMapper = afwTable.SchemaMapper(stdCat.schema)
545 sourceMapper.addMinimalSchema(afwTable.SimpleTable.makeMinimalSchema())
546 sourceMapper.editOutputSchema().addField('instFlux', type=np.float64,
547 doc="instrumental flux (counts)")
548 sourceMapper.editOutputSchema().addField('instFluxErr', type=np.float64,
549 doc="instrumental flux error (counts)")
550 badStarKey = sourceMapper.editOutputSchema().addField('flag_badStar',
551 type='Flag',
552 doc="bad flag")
553
554 # Split up the stars
555 # Note that there is an assumption here that the ra/dec coords stored
556 # on-disk are in radians, and therefore that starObs['coord_ra'] /
557 # starObs['coord_dec'] return radians when used as an array of numpy float64s.
558 ipring = hpg.angle_to_pixel(
559 self.config.referencePixelizationNside,
560 stdCat['coord_ra'],
561 stdCat['coord_dec'],
562 degrees=False,
563 )
564 h, rev = esutil.stat.histogram(ipring, rev=True)
565
566 gdpix, = np.where(h >= self.config.referencePixelizationMinStars)
567
568 self.log.info("Found %d pixels (nside=%d) with at least %d good stars" %
569 (gdpix.size,
570 self.config.referencePixelizationNside,
571 self.config.referencePixelizationMinStars))
572
573 if gdpix.size < self.config.referencePixelizationNPixels:
574 self.log.warning("Found fewer good pixels (%d) than preferred in configuration (%d)" %
575 (gdpix.size, self.config.referencePixelizationNPixels))
576 else:
577 # Sample out the pixels we want to use
578 gdpix = np.random.choice(gdpix, size=self.config.referencePixelizationNPixels, replace=False)
579
580 results = np.zeros(gdpix.size, dtype=[('hpix', 'i4'),
581 ('nstar', 'i4', len(bands)),
582 ('nmatch', 'i4', len(bands)),
583 ('zp', 'f4', len(bands)),
584 ('zpErr', 'f4', len(bands))])
585 results['hpix'] = ipring[rev[rev[gdpix]]]
586
587 # We need a boolean index to deal with catalogs...
588 selected = np.zeros(len(stdCat), dtype=bool)
589
590 refFluxFields = [None]*len(bands)
591
592 for p_index, pix in enumerate(gdpix):
593 i1a = rev[rev[pix]: rev[pix + 1]]
594
595 # the stdCat afwTable can only be indexed with boolean arrays,
596 # and not numpy index arrays (see DM-16497). This little trick
597 # converts the index array into a boolean array
598 selected[:] = False
599 selected[i1a] = True
600
601 for b_index, filterLabel in enumerate(filterLabels):
602 struct = self._computeOffsetOneBand(sourceMapper, badStarKey, b_index,
603 filterLabel, stdCat,
604 selected, refFluxFields)
605 results['nstar'][p_index, b_index] = len(i1a)
606 results['nmatch'][p_index, b_index] = len(struct.arrays.refMag)
607 results['zp'][p_index, b_index] = struct.zp
608 results['zpErr'][p_index, b_index] = struct.sigma
609
610 # And compute the summary statistics
611 offsets = np.zeros(len(bands))
612
613 for b_index, band in enumerate(bands):
614 # make configurable
615 ok, = np.where(results['nmatch'][:, b_index] >= self.config.referenceMinMatch)
616 offsets[b_index] = np.median(results['zp'][ok, b_index])
617 # use median absolute deviation to estimate Normal sigma
618 # see https://en.wikipedia.org/wiki/Median_absolute_deviation
619 madSigma = 1.4826*np.median(np.abs(results['zp'][ok, b_index] - offsets[b_index]))
620 self.log.info("Reference catalog offset for %s band: %.12f +/- %.12f",
621 band, offsets[b_index], madSigma)
622
623 return offsets
624
625 def _computeOffsetOneBand(self, sourceMapper, badStarKey,
626 b_index, filterLabel, stdCat, selected, refFluxFields):
627 """
628 Compute the zeropoint offset between the fgcm stdCat and the reference
629 stars for one pixel in one band
630
631 Parameters
632 ----------
633 sourceMapper : `lsst.afw.table.SchemaMapper`
634 Mapper to go from stdCat to calibratable catalog
635 badStarKey : `lsst.afw.table.Key`
636 Key for the field with bad stars
637 b_index : `int`
638 Index of the band in the star catalog
639 filterLabel : `lsst.afw.image.FilterLabel`
640 filterLabel with band and physical filter
641 stdCat : `lsst.afw.table.SimpleCatalog`
642 FGCM standard stars
643 selected : `numpy.array(dtype=bool)`
644 Boolean array of which stars are in the pixel
645 refFluxFields : `list`
646 List of names of flux fields for reference catalog
647 """
648
649 sourceCat = afwTable.SimpleCatalog(sourceMapper.getOutputSchema())
650 sourceCat.reserve(selected.sum())
651 sourceCat.extend(stdCat[selected], mapper=sourceMapper)
652 sourceCat['instFlux'] = 10.**(stdCat['mag_std_noabs'][selected, b_index]/(-2.5))
653 sourceCat['instFluxErr'] = (np.log(10.)/2.5)*(stdCat['magErr_std'][selected, b_index]
654 * sourceCat['instFlux'])
655 # Make sure we only use stars that have valid measurements
656 # (This is perhaps redundant with requirements above that the
657 # stars be observed in all bands, but it can't hurt)
658 badStar = (stdCat['mag_std_noabs'][selected, b_index] > 90.0)
659 for rec in sourceCat[badStar]:
660 rec.set(badStarKey, True)
661
662 exposure = afwImage.ExposureF()
663 exposure.setFilter(filterLabel)
664
665 if refFluxFields[b_index] is None:
666 # Need to find the flux field in the reference catalog
667 # to work around limitations of DirectMatch in PhotoCal
668 ctr = stdCat[0].getCoord()
669 rad = 0.05*lsst.geom.degrees
670 refDataTest = self.refObjLoader.loadSkyCircle(ctr, rad, filterLabel.bandLabel)
671 refFluxFields[b_index] = refDataTest.fluxField
672
673 # Make a copy of the config so that we can modify it
674 calConfig = copy.copy(self.config.photoCal.value)
675 calConfig.match.referenceSelection.signalToNoise.fluxField = refFluxFields[b_index]
676 calConfig.match.referenceSelection.signalToNoise.errField = refFluxFields[b_index] + 'Err'
677 calTask = self.config.photoCal.target(refObjLoader=self.refObjLoader,
678 config=calConfig,
679 schema=sourceCat.getSchema())
680
681 struct = calTask.run(exposure, sourceCat)
682
683 return struct
684
685 def _outputZeropoints(self, camera, zptCat, visitCat, offsets, bands,
686 physicalFilterMap, tract=None):
687 """Output the zeropoints in fgcm_photoCalib format.
688
689 Parameters
690 ----------
691 camera : `lsst.afw.cameraGeom.Camera`
692 Camera from the butler.
693 zptCat : `lsst.afw.table.BaseCatalog`
694 FGCM zeropoint catalog from `FgcmFitCycleTask`.
695 visitCat : `lsst.afw.table.BaseCatalog`
696 FGCM visitCat from `FgcmBuildStarsTask`.
697 offsets : `numpy.array`
698 Float array of absolute calibration offsets, one for each filter.
699 bands : `list` [`str`]
700 List of band names from FGCM output.
701 physicalFilterMap : `dict`
702 Dictionary of mappings from physical filter to FGCM band.
703 tract: `int`, optional
704 Tract number to output. Default is None (global calibration)
705
706 Returns
707 -------
708 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
709 Generator that returns (visit, exposureCatalog) tuples.
710 """
711 # Select visit/ccds where we have a calibration
712 # This includes ccds where we were able to interpolate from neighboring
713 # ccds.
714 cannot_compute = fgcm.fgcmUtilities.zpFlagDict['CANNOT_COMPUTE_ZEROPOINT']
715 selected = (((zptCat['fgcmFlag'] & cannot_compute) == 0)
716 & (zptCat['fgcmZptVar'] > 0.0)
717 & (zptCat['fgcmZpt'] > FGCM_ILLEGAL_VALUE))
718
719 # Log warnings for any visit which has no valid zeropoints
720 badVisits = np.unique(zptCat['visit'][~selected])
721 goodVisits = np.unique(zptCat['visit'][selected])
722 allBadVisits = badVisits[~np.isin(badVisits, goodVisits)]
723 for allBadVisit in allBadVisits:
724 self.log.warning(f'No suitable photoCalib for visit {allBadVisit}')
725
726 # Get a mapping from filtername to the offsets
727 offsetMapping = {}
728 for f in physicalFilterMap:
729 # Not every filter in the map will necesarily have a band.
730 if physicalFilterMap[f] in bands:
731 offsetMapping[f] = offsets[bands.index(physicalFilterMap[f])]
732
733 # Get a mapping from "ccd" to the ccd index used for the scaling
734 ccdMapping = {}
735 for ccdIndex, detector in enumerate(camera):
736 ccdMapping[detector.getId()] = ccdIndex
737
738 # And a mapping to get the flat-field scaling values
739 scalingMapping = {}
740 for rec in visitCat:
741 scalingMapping[rec['visit']] = rec['scaling']
742
743 if self.config.doComposeWcsJacobian:
744 approxPixelAreaFields = computeApproxPixelAreaFields(camera)
745
746 # The zptCat is sorted by visit, which is useful
747 lastVisit = -1
748 zptVisitCatalog = None
749
750 metadata = dafBase.PropertyList()
751 metadata.add("COMMENT", "Catalog id is detector id, sorted.")
752 metadata.add("COMMENT", "Only detectors with data have entries.")
753
754 for rec in zptCat[selected]:
755 # Retrieve overall scaling
756 scaling = scalingMapping[rec['visit']][ccdMapping[rec['detector']]]
757
758 # The postCalibrationOffset describe any zeropoint offsets
759 # to apply after the fgcm calibration. The first part comes
760 # from the reference catalog match (used in testing). The
761 # second part comes from the mean chromatic correction
762 # (if configured).
763 postCalibrationOffset = offsetMapping[rec['filtername']]
764 if self.config.doApplyMeanChromaticCorrection:
765 postCalibrationOffset += rec['fgcmDeltaChrom']
766
767 fgcmSuperStarField = self._getChebyshevBoundedField(rec['fgcmfZptSstarCheb'],
768 rec['fgcmfZptChebXyMax'])
769 # Convert from FGCM AB to nJy
770 fgcmZptField = self._getChebyshevBoundedField((rec['fgcmfZptCheb']*units.AB).to_value(units.nJy),
771 rec['fgcmfZptChebXyMax'],
772 offset=postCalibrationOffset,
773 scaling=scaling)
774
775 if self.config.doComposeWcsJacobian:
776
777 fgcmField = afwMath.ProductBoundedField([approxPixelAreaFields[rec['detector']],
778 fgcmSuperStarField,
779 fgcmZptField])
780 else:
781 # The photoCalib is just the product of the fgcmSuperStarField and the
782 # fgcmZptField
783 fgcmField = afwMath.ProductBoundedField([fgcmSuperStarField, fgcmZptField])
784
785 # The "mean" calibration will be set to the center of the ccd for reference
786 calibCenter = fgcmField.evaluate(fgcmField.getBBox().getCenter())
787 calibErr = (np.log(10.0)/2.5)*calibCenter*np.sqrt(rec['fgcmZptVar'])
788 photoCalib = afwImage.PhotoCalib(calibrationMean=calibCenter,
789 calibrationErr=calibErr,
790 calibration=fgcmField,
791 isConstant=False)
792
793 # Return full per-visit exposure catalogs
794 if rec['visit'] != lastVisit:
795 # This is a new visit. If the last visit was not -1, yield
796 # the ExposureCatalog
797 if lastVisit > -1:
798 # ensure that the detectors are in sorted order, for fast lookups
799 zptVisitCatalog.sort()
800 yield (int(lastVisit), zptVisitCatalog)
801 else:
802 # We need to create a new schema
803 zptExpCatSchema = afwTable.ExposureTable.makeMinimalSchema()
804 zptExpCatSchema.addField('visit', type='L', doc='Visit number')
805
806 # And start a new one
807 zptVisitCatalog = afwTable.ExposureCatalog(zptExpCatSchema)
808 zptVisitCatalog.setMetadata(metadata)
809
810 lastVisit = int(rec['visit'])
811
812 catRecord = zptVisitCatalog.addNew()
813 catRecord['id'] = int(rec['detector'])
814 catRecord['visit'] = rec['visit']
815 catRecord.setPhotoCalib(photoCalib)
816
817 # Final output of last exposure catalog
818 # ensure that the detectors are in sorted order, for fast lookups
819 zptVisitCatalog.sort()
820 yield (int(lastVisit), zptVisitCatalog)
821
822 def _getChebyshevBoundedField(self, coefficients, xyMax, offset=0.0, scaling=1.0):
823 """
824 Make a ChebyshevBoundedField from fgcm coefficients, with optional offset
825 and scaling.
826
827 Parameters
828 ----------
829 coefficients: `numpy.array`
830 Flattened array of chebyshev coefficients
831 xyMax: `list` of length 2
832 Maximum x and y of the chebyshev bounding box
833 offset: `float`, optional
834 Absolute calibration offset. Default is 0.0
835 scaling: `float`, optional
836 Flat scaling value from fgcmBuildStars. Default is 1.0
837
838 Returns
839 -------
840 boundedField: `lsst.afw.math.ChebyshevBoundedField`
841 """
842
843 orderPlus1 = int(np.sqrt(coefficients.size))
844 pars = np.zeros((orderPlus1, orderPlus1))
845
846 bbox = lsst.geom.Box2I(lsst.geom.Point2I(0.0, 0.0),
847 lsst.geom.Point2I(*xyMax))
848
849 pars[:, :] = (coefficients.reshape(orderPlus1, orderPlus1)
850 * (10.**(offset/-2.5))*scaling)
851
852 boundedField = afwMath.ChebyshevBoundedField(bbox, pars)
853
854 return boundedField
855
856 def _outputAtmospheres(self, handleDict, atmCat):
857 """
858 Output the atmospheres.
859
860 Parameters
861 ----------
862 handleDict : `dict`
863 All data handles are `lsst.daf.butler.DeferredDatasetHandle`
864 The handleDict has the follownig keys:
865
866 ``"fgcmLookUpTable"``
867 handle for the FGCM look-up table.
868 atmCat : `lsst.afw.table.BaseCatalog`
869 FGCM atmosphere parameter catalog from fgcmFitCycleTask.
870
871 Returns
872 -------
873 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
874 Generator that returns (visit, transmissionCurve) tuples.
875 """
876 # First, we need to grab the look-up table and key info
877 lutCat = handleDict['fgcmLookUpTable'].get()
878
879 atmosphereTableName = lutCat[0]['tablename']
880 elevation = lutCat[0]['elevation']
881 atmLambda = lutCat[0]['atmLambda']
882 lutCat = None
883
884 # Make the atmosphere table if possible
885 try:
886 atmTable = fgcm.FgcmAtmosphereTable.initWithTableName(atmosphereTableName)
887 atmTable.loadTable()
888 except IOError:
889 atmTable = None
890
891 if atmTable is None:
892 # Try to use MODTRAN instead
893 try:
894 modGen = fgcm.ModtranGenerator(elevation)
895 lambdaRange = np.array([atmLambda[0], atmLambda[-1]])/10.
896 lambdaStep = (atmLambda[1] - atmLambda[0])/10.
897 except (ValueError, IOError) as e:
898 raise RuntimeError("FGCM look-up-table generated with modtran, "
899 "but modtran not configured to run.") from e
900
901 zenith = np.degrees(np.arccos(1./atmCat['secZenith']))
902
903 for i, visit in enumerate(atmCat['visit']):
904 if atmTable is not None:
905 # Interpolate the atmosphere table
906 atmVals = atmTable.interpolateAtmosphere(pmb=atmCat[i]['pmb'],
907 pwv=atmCat[i]['pwv'],
908 o3=atmCat[i]['o3'],
909 tau=atmCat[i]['tau'],
910 alpha=atmCat[i]['alpha'],
911 zenith=zenith[i],
912 ctranslamstd=[atmCat[i]['cTrans'],
913 atmCat[i]['lamStd']])
914 else:
915 # Run modtran
916 modAtm = modGen(pmb=atmCat[i]['pmb'],
917 pwv=atmCat[i]['pwv'],
918 o3=atmCat[i]['o3'],
919 tau=atmCat[i]['tau'],
920 alpha=atmCat[i]['alpha'],
921 zenith=zenith[i],
922 lambdaRange=lambdaRange,
923 lambdaStep=lambdaStep,
924 ctranslamstd=[atmCat[i]['cTrans'],
925 atmCat[i]['lamStd']])
926 atmVals = modAtm['COMBINED']
927
928 # Now need to create something to persist...
929 curve = TransmissionCurve.makeSpatiallyConstant(throughput=atmVals,
930 wavelengths=atmLambda,
931 throughputAtMin=atmVals[0],
932 throughputAtMax=atmVals[-1])
933
934 yield (int(visit), curve)