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