lsst.fgcmcal g2ffcdf413f+c922cbd9a8
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 cycleNumber = pexConfig.Field(
172 doc="Final fit cycle from FGCM fit",
173 dtype=int,
174 default=0,
175 deprecated=("This config is no longer used, and will be removed after v25. "
176 "Please set config.connections.cycleNumber directly instead."),
177 )
178 physicalFilterMap = pexConfig.DictField(
179 doc="Mapping from 'physicalFilter' to band.",
180 keytype=str,
181 itemtype=str,
182 default={},
183 )
184 # The following fields refer to calibrating from a reference
185 # catalog, but in the future this might need to be expanded
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."),
190 dtype=bool,
191 default=False,
192 )
193 doAtmosphereOutput = pexConfig.Field(
194 doc="Output atmospheres in transmission_atmosphere_fgcm format",
195 dtype=bool,
196 default=True,
197 )
198 doZeropointOutput = pexConfig.Field(
199 doc="Output zeropoints in fgcm_photoCalib format",
200 dtype=bool,
201 default=True,
202 )
203 doComposeWcsJacobian = pexConfig.Field(
204 doc="Compose Jacobian of WCS with fgcm calibration for output photoCalib?",
205 dtype=bool,
206 default=True,
207 )
208 doApplyMeanChromaticCorrection = pexConfig.Field(
209 doc="Apply the mean chromatic correction to the zeropoints?",
210 dtype=bool,
211 default=True,
212 )
213 photoCal = pexConfig.ConfigurableField(
214 target=PhotoCalTask,
215 doc="task to perform 'absolute' calibration",
216 )
217 referencePixelizationNside = pexConfig.Field(
218 doc="Healpix nside to pixelize catalog to compare to reference catalog",
219 dtype=int,
220 default=64,
221 )
222 referencePixelizationMinStars = pexConfig.Field(
223 doc=("Minimum number of stars per healpix pixel to select for comparison"
224 "to the specified reference catalog"),
225 dtype=int,
226 default=200,
227 )
228 referenceMinMatch = pexConfig.Field(
229 doc="Minimum number of stars matched to reference catalog to be used in statistics",
230 dtype=int,
231 default=50,
232 )
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."),
238 dtype=int,
239 default=100,
240 )
241
242 def setDefaults(self):
243 pexConfig.Config.setDefaults(self)
244
245 # In order to transfer the "absolute" calibration from a reference
246 # catalog to the relatively calibrated FGCM standard stars (one number
247 # per band), we use the PhotoCalTask to match stars in a sample of healpix
248 # pixels. These basic settings ensure that only well-measured, good stars
249 # from the source and reference catalogs are used for the matching.
250
251 # applyColorTerms needs to be False if doReferenceCalibration is False,
252 # as is the new default after DM-16702
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
266
267
268class FgcmOutputProductsTask(pipeBase.PipelineTask):
269 """
270 Output products from FGCM global calibration.
271 """
272
273 ConfigClass = FgcmOutputProductsConfig
274 _DefaultName = "fgcmOutputProducts"
275
276 def __init__(self, **kwargs):
277 super().__init__(**kwargs)
278
279 def runQuantum(self, butlerQC, inputRefs, outputRefs):
280 handleDict = {}
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)
285
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}
290
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}
295
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,
302 log=self.log,
303 config=refConfig)
304 else:
305 self.refObjLoader = None
306
307 struct = self.run(handleDict, self.config.physicalFilterMap)
308
309 # Output the photoCalib exposure catalogs
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.")
315
316 # Output the atmospheres
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.")
322
323 if self.config.doReferenceCalibration:
324 # Turn offset into simple catalog for persistence if necessary
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
331
332 butlerQC.put(offsetCat, outputRefs.fgcmOffsets)
333
334 return
335
336 def run(self, handleDict, physicalFilterMap):
337 """Run the output products task.
338
339 Parameters
340 ----------
341 handleDict : `dict`
342 All handles are `lsst.daf.butler.DeferredDatasetHandle`
343 handle dictionary with keys:
344
345 ``"camera"``
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.
353 ``"fgcmZeropoints"``
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.
361
362 Returns
363 -------
364 retStruct : `lsst.pipe.base.Struct`
365 Output structure with keys:
366
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.
373 """
374 stdCat = handleDict['fgcmStandardStars'].get()
375 md = stdCat.getMetadata()
376 bands = md.getArray('BANDS')
377
378 if self.config.doReferenceCalibration:
379 lutCat = handleDict['fgcmLookUpTable'].get()
380 offsets = self._computeReferenceOffsets(stdCat, lutCat, physicalFilterMap, bands)
381 else:
382 offsets = np.zeros(len(bands))
383
384 del stdCat
385
386 if self.config.doZeropointOutput:
387 zptCat = handleDict['fgcmZeropoints'].get()
388 visitCat = handleDict['fgcmVisitCatalog'].get()
389
390 pcgen = self._outputZeropoints(handleDict['camera'], zptCat, visitCat, offsets, bands,
391 physicalFilterMap)
392 else:
393 pcgen = None
394
395 if self.config.doAtmosphereOutput:
396 atmCat = handleDict['fgcmAtmosphereParameters'].get()
397 atmgen = self._outputAtmospheres(handleDict, atmCat)
398 else:
399 atmgen = None
400
401 retStruct = pipeBase.Struct(offsets=offsets,
402 atmospheres=atmgen)
403 retStruct.photoCalibCatalogs = pcgen
404
405 return retStruct
406
407 def generateTractOutputProducts(self, handleDict, tract,
408 visitCat, zptCat, atmCat, stdCat,
409 fgcmBuildStarsConfig):
410 """
411 Generate the output products for a given tract, as specified in the config.
412
413 This method is here to have an alternate entry-point for
414 FgcmCalibrateTract.
415
416 Parameters
417 ----------
418 handleDict : `dict`
419 All handles are `lsst.daf.butler.DeferredDatasetHandle`
420 handle dictionary with keys:
421
422 ``"camera"``
423 Camera object (`lsst.afw.cameraGeom.Camera`)
424 ``"fgcmLookUpTable"``
425 handle for the FGCM look-up table.
426 tract : `int`
427 Tract number
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`
438
439 Returns
440 -------
441 retStruct : `lsst.pipe.base.Struct`
442 Output structure with keys:
443
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.
450 """
451 physicalFilterMap = fgcmBuildStarsConfig.physicalFilterMap
452
453 md = stdCat.getMetadata()
454 bands = md.getArray('BANDS')
455
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.")
459
460 if not self.config.doComposeWcsJacobian and fgcmBuildStarsConfig.doApplyWcsJacobian:
461 self.log.warning("Jacobian was applied in build-stars but doComposeWcsJacobian is not set.")
462
463 if self.config.doReferenceCalibration:
464 lutCat = handleDict['fgcmLookUpTable'].get()
465 offsets = self._computeReferenceOffsets(stdCat, lutCat, bands, physicalFilterMap)
466 else:
467 offsets = np.zeros(len(bands))
468
469 if self.config.doZeropointOutput:
470 pcgen = self._outputZeropoints(handleDict['camera'], zptCat, visitCat, offsets, bands,
471 physicalFilterMap)
472 else:
473 pcgen = None
474
475 if self.config.doAtmosphereOutput:
476 atmgen = self._outputAtmospheres(handleDict, atmCat)
477 else:
478 atmgen = None
479
480 retStruct = pipeBase.Struct(offsets=offsets,
481 atmospheres=atmgen)
482 retStruct.photoCalibCatalogs = pcgen
483
484 return retStruct
485
486 def _computeReferenceOffsets(self, stdCat, lutCat, physicalFilterMap, bands):
487 """
488 Compute offsets relative to a reference catalog.
489
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
494 absolute scale.
495
496 Parameters
497 ----------
498 stdCat : `lsst.afw.table.SimpleCatalog`
499 FGCM standard stars
500 lutCat : `lsst.afw.table.SimpleCatalog`
501 FGCM Look-up table
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
506 Returns
507 -------
508 offsets : `numpy.array` of floats
509 Per band zeropoint offsets
510 """
511
512 # Only use stars that are observed in all the bands that were actually used
513 # This will ensure that we use the same healpix pixels for the absolute
514 # calibration of each band.
515 minObs = stdCat['ngood'].min(axis=1)
516
517 goodStars = (minObs >= 1)
518 stdCat = stdCat[goodStars]
519
520 self.log.info("Found %d stars with at least 1 good observation in each band" %
521 (len(stdCat)))
522
523 # Associate each band with the appropriate physicalFilter and make
524 # filterLabels
525 filterLabels = []
526
527 lutPhysicalFilters = lutCat[0]['physicalFilters'].split(',')
528 lutStdPhysicalFilters = lutCat[0]['stdPhysicalFilters'].split(',')
529 physicalFilterMapBands = list(physicalFilterMap.values())
530 physicalFilterMapFilters = list(physicalFilterMap.keys())
531 for band in bands:
532 # Find a physical filter associated from the band by doing
533 # a reverse lookup on the physicalFilterMap dict
534 physicalFilterMapIndex = physicalFilterMapBands.index(band)
535 physicalFilter = physicalFilterMapFilters[physicalFilterMapIndex]
536 # Find the appropriate fgcm standard physicalFilter
537 lutPhysicalFilterIndex = lutPhysicalFilters.index(physicalFilter)
538 stdPhysicalFilter = lutStdPhysicalFilters[lutPhysicalFilterIndex]
539 filterLabels.append(afwImage.FilterLabel(band=band,
540 physical=stdPhysicalFilter))
541
542 # We have to make a table for each pixel with flux/fluxErr
543 # This is a temporary table generated for input to the photoCal task.
544 # These fluxes are not instFlux (they are top-of-the-atmosphere approximate and
545 # have had chromatic corrections applied to get to the standard system
546 # specified by the atmosphere/instrumental parameters), nor are they
547 # in Jansky (since they don't have a proper absolute calibration: the overall
548 # zeropoint is estimated from the telescope size, etc.)
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',
556 type='Flag',
557 doc="bad flag")
558
559 # Split up the stars
560 # Note that there is an assumption here that the ra/dec coords stored
561 # on-disk are in radians, and therefore that starObs['coord_ra'] /
562 # starObs['coord_dec'] return radians when used as an array of numpy float64s.
563 ipring = hpg.angle_to_pixel(
564 self.config.referencePixelizationNside,
565 stdCat['coord_ra'],
566 stdCat['coord_dec'],
567 degrees=False,
568 )
569 h, rev = esutil.stat.histogram(ipring, rev=True)
570
571 gdpix, = np.where(h >= self.config.referencePixelizationMinStars)
572
573 self.log.info("Found %d pixels (nside=%d) with at least %d good stars" %
574 (gdpix.size,
575 self.config.referencePixelizationNside,
576 self.config.referencePixelizationMinStars))
577
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))
581 else:
582 # Sample out the pixels we want to use
583 gdpix = np.random.choice(gdpix, size=self.config.referencePixelizationNPixels, replace=False)
584
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]]]
591
592 # We need a boolean index to deal with catalogs...
593 selected = np.zeros(len(stdCat), dtype=bool)
594
595 refFluxFields = [None]*len(bands)
596
597 for p_index, pix in enumerate(gdpix):
598 i1a = rev[rev[pix]: rev[pix + 1]]
599
600 # the stdCat afwTable can only be indexed with boolean arrays,
601 # and not numpy index arrays (see DM-16497). This little trick
602 # converts the index array into a boolean array
603 selected[:] = False
604 selected[i1a] = True
605
606 for b_index, filterLabel in enumerate(filterLabels):
607 struct = self._computeOffsetOneBand(sourceMapper, badStarKey, b_index,
608 filterLabel, stdCat,
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
614
615 # And compute the summary statistics
616 offsets = np.zeros(len(bands))
617
618 for b_index, band in enumerate(bands):
619 # make configurable
620 ok, = np.where(results['nmatch'][:, b_index] >= self.config.referenceMinMatch)
621 offsets[b_index] = np.median(results['zp'][ok, b_index])
622 # use median absolute deviation to estimate Normal sigma
623 # see https://en.wikipedia.org/wiki/Median_absolute_deviation
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)
627
628 return offsets
629
630 def _computeOffsetOneBand(self, sourceMapper, badStarKey,
631 b_index, filterLabel, stdCat, selected, refFluxFields):
632 """
633 Compute the zeropoint offset between the fgcm stdCat and the reference
634 stars for one pixel in one band
635
636 Parameters
637 ----------
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
642 b_index : `int`
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`
647 FGCM standard stars
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
652 """
653
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'])
660 # Make sure we only use stars that have valid measurements
661 # (This is perhaps redundant with requirements above that the
662 # stars be observed in all bands, but it can't hurt)
663 badStar = (stdCat['mag_std_noabs'][selected, b_index] > 90.0)
664 for rec in sourceCat[badStar]:
665 rec.set(badStarKey, True)
666
667 exposure = afwImage.ExposureF()
668 exposure.setFilter(filterLabel)
669
670 if refFluxFields[b_index] is None:
671 # Need to find the flux field in the reference catalog
672 # to work around limitations of DirectMatch in PhotoCal
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
677
678 # Make a copy of the config so that we can modify it
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,
683 config=calConfig,
684 schema=sourceCat.getSchema())
685
686 struct = calTask.run(exposure, sourceCat)
687
688 return struct
689
690 def _outputZeropoints(self, camera, zptCat, visitCat, offsets, bands,
691 physicalFilterMap, tract=None):
692 """Output the zeropoints in fgcm_photoCalib format.
693
694 Parameters
695 ----------
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)
710
711 Returns
712 -------
713 photoCalibCatalogs : `generator` [(`int`, `lsst.afw.table.ExposureCatalog`)]
714 Generator that returns (visit, exposureCatalog) tuples.
715 """
716 # Select visit/ccds where we have a calibration
717 # This includes ccds where we were able to interpolate from neighboring
718 # ccds.
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))
723
724 # Log warnings for any visit which has no valid zeropoints
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}')
730
731 # Get a mapping from filtername to the offsets
732 offsetMapping = {}
733 for f in physicalFilterMap:
734 # Not every filter in the map will necesarily have a band.
735 if physicalFilterMap[f] in bands:
736 offsetMapping[f] = offsets[bands.index(physicalFilterMap[f])]
737
738 # Get a mapping from "ccd" to the ccd index used for the scaling
739 ccdMapping = {}
740 for ccdIndex, detector in enumerate(camera):
741 ccdMapping[detector.getId()] = ccdIndex
742
743 # And a mapping to get the flat-field scaling values
744 scalingMapping = {}
745 for rec in visitCat:
746 scalingMapping[rec['visit']] = rec['scaling']
747
748 if self.config.doComposeWcsJacobian:
749 approxPixelAreaFields = computeApproxPixelAreaFields(camera)
750
751 # The zptCat is sorted by visit, which is useful
752 lastVisit = -1
753 zptVisitCatalog = None
754
755 metadata = dafBase.PropertyList()
756 metadata.add("COMMENT", "Catalog id is detector id, sorted.")
757 metadata.add("COMMENT", "Only detectors with data have entries.")
758
759 for rec in zptCat[selected]:
760 # Retrieve overall scaling
761 scaling = scalingMapping[rec['visit']][ccdMapping[rec['detector']]]
762
763 # The postCalibrationOffset describe any zeropoint offsets
764 # to apply after the fgcm calibration. The first part comes
765 # from the reference catalog match (used in testing). The
766 # second part comes from the mean chromatic correction
767 # (if configured).
768 postCalibrationOffset = offsetMapping[rec['filtername']]
769 if self.config.doApplyMeanChromaticCorrection:
770 postCalibrationOffset += rec['fgcmDeltaChrom']
771
772 fgcmSuperStarField = self._getChebyshevBoundedField(rec['fgcmfZptSstarCheb'],
773 rec['fgcmfZptChebXyMax'])
774 # Convert from FGCM AB to nJy
775 fgcmZptField = self._getChebyshevBoundedField((rec['fgcmfZptCheb']*units.AB).to_value(units.nJy),
776 rec['fgcmfZptChebXyMax'],
777 offset=postCalibrationOffset,
778 scaling=scaling)
779
780 if self.config.doComposeWcsJacobian:
781
782 fgcmField = afwMath.ProductBoundedField([approxPixelAreaFields[rec['detector']],
783 fgcmSuperStarField,
784 fgcmZptField])
785 else:
786 # The photoCalib is just the product of the fgcmSuperStarField and the
787 # fgcmZptField
788 fgcmField = afwMath.ProductBoundedField([fgcmSuperStarField, fgcmZptField])
789
790 # The "mean" calibration will be set to the center of the ccd for reference
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,
796 isConstant=False)
797
798 # Return full per-visit exposure catalogs
799 if rec['visit'] != lastVisit:
800 # This is a new visit. If the last visit was not -1, yield
801 # the ExposureCatalog
802 if lastVisit > -1:
803 # ensure that the detectors are in sorted order, for fast lookups
804 zptVisitCatalog.sort()
805 yield (int(lastVisit), zptVisitCatalog)
806 else:
807 # We need to create a new schema
808 zptExpCatSchema = afwTable.ExposureTable.makeMinimalSchema()
809 zptExpCatSchema.addField('visit', type='L', doc='Visit number')
810
811 # And start a new one
812 zptVisitCatalog = afwTable.ExposureCatalog(zptExpCatSchema)
813 zptVisitCatalog.setMetadata(metadata)
814
815 lastVisit = int(rec['visit'])
816
817 catRecord = zptVisitCatalog.addNew()
818 catRecord['id'] = int(rec['detector'])
819 catRecord['visit'] = rec['visit']
820 catRecord.setPhotoCalib(photoCalib)
821
822 # Final output of last exposure catalog
823 # ensure that the detectors are in sorted order, for fast lookups
824 zptVisitCatalog.sort()
825 yield (int(lastVisit), zptVisitCatalog)
826
827 def _getChebyshevBoundedField(self, coefficients, xyMax, offset=0.0, scaling=1.0):
828 """
829 Make a ChebyshevBoundedField from fgcm coefficients, with optional offset
830 and scaling.
831
832 Parameters
833 ----------
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
842
843 Returns
844 -------
845 boundedField: `lsst.afw.math.ChebyshevBoundedField`
846 """
847
848 orderPlus1 = int(np.sqrt(coefficients.size))
849 pars = np.zeros((orderPlus1, orderPlus1))
850
851 bbox = lsst.geom.Box2I(lsst.geom.Point2I(0.0, 0.0),
852 lsst.geom.Point2I(*xyMax))
853
854 pars[:, :] = (coefficients.reshape(orderPlus1, orderPlus1)
855 * (10.**(offset/-2.5))*scaling)
856
857 boundedField = afwMath.ChebyshevBoundedField(bbox, pars)
858
859 return boundedField
860
861 def _outputAtmospheres(self, handleDict, atmCat):
862 """
863 Output the atmospheres.
864
865 Parameters
866 ----------
867 handleDict : `dict`
868 All data handles are `lsst.daf.butler.DeferredDatasetHandle`
869 The handleDict has the follownig keys:
870
871 ``"fgcmLookUpTable"``
872 handle for the FGCM look-up table.
873 atmCat : `lsst.afw.table.BaseCatalog`
874 FGCM atmosphere parameter catalog from fgcmFitCycleTask.
875
876 Returns
877 -------
878 atmospheres : `generator` [(`int`, `lsst.afw.image.TransmissionCurve`)]
879 Generator that returns (visit, transmissionCurve) tuples.
880 """
881 # First, we need to grab the look-up table and key info
882 lutCat = handleDict['fgcmLookUpTable'].get()
883
884 atmosphereTableName = lutCat[0]['tablename']
885 elevation = lutCat[0]['elevation']
886 atmLambda = lutCat[0]['atmLambda']
887 lutCat = None
888
889 # Make the atmosphere table if possible
890 try:
891 atmTable = fgcm.FgcmAtmosphereTable.initWithTableName(atmosphereTableName)
892 atmTable.loadTable()
893 except IOError:
894 atmTable = None
895
896 if atmTable is None:
897 # Try to use MODTRAN instead
898 try:
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
905
906 zenith = np.degrees(np.arccos(1./atmCat['secZenith']))
907
908 for i, visit in enumerate(atmCat['visit']):
909 if atmTable is not None:
910 # Interpolate the atmosphere table
911 atmVals = atmTable.interpolateAtmosphere(pmb=atmCat[i]['pmb'],
912 pwv=atmCat[i]['pwv'],
913 o3=atmCat[i]['o3'],
914 tau=atmCat[i]['tau'],
915 alpha=atmCat[i]['alpha'],
916 zenith=zenith[i],
917 ctranslamstd=[atmCat[i]['cTrans'],
918 atmCat[i]['lamStd']])
919 else:
920 # Run modtran
921 modAtm = modGen(pmb=atmCat[i]['pmb'],
922 pwv=atmCat[i]['pwv'],
923 o3=atmCat[i]['o3'],
924 tau=atmCat[i]['tau'],
925 alpha=atmCat[i]['alpha'],
926 zenith=zenith[i],
927 lambdaRange=lambdaRange,
928 lambdaStep=lambdaStep,
929 ctranslamstd=[atmCat[i]['cTrans'],
930 atmCat[i]['lamStd']])
931 atmVals = modAtm['COMBINED']
932
933 # Now need to create something to persist...
934 curve = TransmissionCurve.makeSpatiallyConstant(throughput=atmVals,
935 wavelengths=atmLambda,
936 throughputAtMin=atmVals[0],
937 throughputAtMax=atmVals[-1])
938
939 yield (int(visit), curve)