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