LSST Applications g1635faa6d4+215bc75b8c,g1653933729+a8ce1bb630,g22ce9dc20b+d972d8df89,g28da252d5a+0fcf840c6d,g29321ee8c0+e558be0e74,g2bbee38e9b+9634bc57db,g2bc492864f+9634bc57db,g2cdde0e794+c2c89b37c4,g3156d2b45e+41e33cbcdc,g347aa1857d+9634bc57db,g35bb328faa+a8ce1bb630,g3a166c0a6a+9634bc57db,g3e281a1b8c+9f2c4e2fc3,g414038480c+077ccc18e7,g41af890bb2+2a6f257a1d,g5fbc88fb19+17cd334064,g781aacb6e4+a8ce1bb630,g7ab3e175f3+59ce30aec6,g80478fca09+f8b2ab54e1,g82479be7b0+ba9d578ff8,g858d7b2824+d972d8df89,g9125e01d80+a8ce1bb630,g9726552aa6+10f999ec6a,ga5288a1d22+630363936d,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gb9c6c11c1e+9553554aa7,gbd46683f8f+0c4209622a,gc28159a63d+9634bc57db,gcf0d15dbbd+2db122af0a,gda3e153d99+d972d8df89,gda6a2b7d83+2db122af0a,gdaeeff99f8+1711a396fd,ge2409df99d+d1dc2f3b25,ge33fd446bb+d972d8df89,ge79ae78c31+9634bc57db,gf0baf85859+147a0692ba,gf3967379c6+02b11634a5,w.2024.45
LSST Data Management Base Package
Loading...
Searching...
No Matches
fgcmMakeLut.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 a look-up-table (LUT) for FGCM calibration.
24
25This task computes a look-up-table for the range in expected atmosphere
26variation and variation in instrumental throughput (as tracked by the
27transmission_filter products). By pre-computing linearized integrals,
28the FGCM fit is orders of magnitude faster for stars with a broad range
29of colors and observing bands, yielding precision at the 1-2 mmag level.
30
31Computing a LUT requires running MODTRAN or with a pre-generated
32atmosphere table packaged with fgcm.
33"""
34
35import numpy as np
36
37import lsst.pex.config as pexConfig
38import lsst.pipe.base as pipeBase
39from lsst.pipe.base import connectionTypes
40import lsst.afw.table as afwTable
41import lsst.afw.cameraGeom as afwCameraGeom
42from lsst.afw.image import TransmissionCurve
43from .utilities import lookupStaticCalibrations
44
45import fgcm
46
47__all__ = ['FgcmMakeLutParametersConfig', 'FgcmMakeLutConfig', 'FgcmMakeLutTask', 'SensorCorrectionTerms']
48
49
50class SensorCorrectionTerms(pexConfig.Config):
51 refLambda = pexConfig.Field(
52 doc="Reference wavelength for first-order correction terms.",
53 dtype=float,
54 optional=False,
55 )
56 correctionTermDict = pexConfig.DictField(
57 doc="Mapping of detector number to first-order correction term.",
58 keytype=int,
59 itemtype=float,
60 default={},
61 )
62
63
64class FgcmMakeLutConnections(pipeBase.PipelineTaskConnections,
65 dimensions=('instrument',),
66 defaultTemplates={}):
67 camera = connectionTypes.PrerequisiteInput(
68 doc="Camera instrument",
69 name="camera",
70 storageClass="Camera",
71 dimensions=("instrument",),
72 isCalibration=True,
73 )
74
75 transmission_optics = connectionTypes.PrerequisiteInput(
76 doc="Optics transmission curve information",
77 name="transmission_optics",
78 storageClass="TransmissionCurve",
79 dimensions=("instrument",),
80 isCalibration=True,
81 deferLoad=True,
82 )
83
84 transmission_sensor = connectionTypes.PrerequisiteInput(
85 doc="Sensor transmission curve information",
86 name="transmission_sensor",
87 storageClass="TransmissionCurve",
88 dimensions=("instrument", "detector",),
89 lookupFunction=lookupStaticCalibrations,
90 isCalibration=True,
91 deferLoad=True,
92 multiple=True,
93 )
94
95 transmission_filter = connectionTypes.PrerequisiteInput(
96 doc="Filter transmission curve information",
97 name="transmission_filter",
98 storageClass="TransmissionCurve",
99 dimensions=("band", "instrument", "physical_filter",),
100 lookupFunction=lookupStaticCalibrations,
101 isCalibration=True,
102 deferLoad=True,
103 multiple=True,
104 )
105
106 fgcmLookUpTable = connectionTypes.Output(
107 doc=("Atmosphere + instrument look-up-table for FGCM throughput and "
108 "chromatic corrections."),
109 name="fgcmLookUpTable",
110 storageClass="Catalog",
111 dimensions=("instrument",),
112 )
113
114 fgcmStandardAtmosphere = connectionTypes.Output(
115 doc="Standard atmosphere used for FGCM calibration.",
116 name="fgcm_standard_atmosphere",
117 storageClass="TransmissionCurve",
118 dimensions=("instrument",),
119 )
120
121 fgcmStandardPassbands = connectionTypes.Output(
122 doc="Standard passbands used for FGCM calibration.",
123 name="fgcm_standard_passband",
124 storageClass="TransmissionCurve",
125 dimensions=("instrument", "physical_filter"),
126 multiple=True,
127 )
128
129 def __init__(self, *, config=None):
130 if not config.doOpticsTransmission:
131 del self.transmission_optics
132 if not config.doSensorTransmission:
133 del self.transmission_sensor
134
135
136class FgcmMakeLutParametersConfig(pexConfig.Config):
137 """Config for parameters if atmosphereTableName not available"""
138 # TODO: When DM-16511 is done, it will be possible to get the
139 # telescope elevation directly from the camera.
140 elevation = pexConfig.Field(
141 doc="Telescope elevation (m)",
142 dtype=float,
143 default=None,
144 )
145 pmbRange = pexConfig.ListField(
146 doc=("Barometric Pressure range (millibar) "
147 "Recommended range depends on the site."),
148 dtype=float,
149 default=None,
150 )
151 pmbSteps = pexConfig.Field(
152 doc="Barometric Pressure number of steps",
153 dtype=int,
154 default=5,
155 )
156 pwvRange = pexConfig.ListField(
157 doc=("Precipitable Water Vapor range (mm) "
158 "Recommended range depends on the site."),
159 dtype=float,
160 default=None,
161 )
162 pwvSteps = pexConfig.Field(
163 doc="Precipitable Water Vapor number of steps",
164 dtype=int,
165 default=15,
166 )
167 o3Range = pexConfig.ListField(
168 doc="Ozone range (dob)",
169 dtype=float,
170 default=[220.0, 310.0],
171 )
172 o3Steps = pexConfig.Field(
173 doc="Ozone number of steps",
174 dtype=int,
175 default=3,
176 )
177 tauRange = pexConfig.ListField(
178 doc="Aerosol Optical Depth range (unitless)",
179 dtype=float,
180 default=[0.002, 0.35],
181 )
182 tauSteps = pexConfig.Field(
183 doc="Aerosol Optical Depth number of steps",
184 dtype=int,
185 default=11,
186 )
187 alphaRange = pexConfig.ListField(
188 doc="Aerosol alpha range (unitless)",
189 dtype=float,
190 default=[0.0, 2.0],
191 )
192 alphaSteps = pexConfig.Field(
193 doc="Aerosol alpha number of steps",
194 dtype=int,
195 default=9,
196 )
197 zenithRange = pexConfig.ListField(
198 doc="Zenith angle range (degree)",
199 dtype=float,
200 default=[0.0, 70.0],
201 )
202 zenithSteps = pexConfig.Field(
203 doc="Zenith angle number of steps",
204 dtype=int,
205 default=21,
206 )
207 # Note that the standard atmosphere parameters depend on the observatory
208 # and elevation, and so these should be set on a per-camera basis.
209 pmbStd = pexConfig.Field(
210 doc=("Standard Atmosphere pressure (millibar); "
211 "Recommended default depends on the site."),
212 dtype=float,
213 default=None,
214 )
215 pwvStd = pexConfig.Field(
216 doc=("Standard Atmosphere PWV (mm); "
217 "Recommended default depends on the site."),
218 dtype=float,
219 default=None,
220 )
221 o3Std = pexConfig.Field(
222 doc="Standard Atmosphere O3 (dob)",
223 dtype=float,
224 default=263.0,
225 )
226 tauStd = pexConfig.Field(
227 doc="Standard Atmosphere aerosol optical depth",
228 dtype=float,
229 default=0.03,
230 )
231 alphaStd = pexConfig.Field(
232 doc="Standard Atmosphere aerosol alpha",
233 dtype=float,
234 default=1.0,
235 )
236 airmassStd = pexConfig.Field(
237 doc=("Standard Atmosphere airmass; "
238 "Recommended default depends on the survey strategy."),
239 dtype=float,
240 default=None,
241 )
242 lambdaNorm = pexConfig.Field(
243 doc="Aerosol Optical Depth normalization wavelength (Angstrom)",
244 dtype=float,
245 default=7750.0,
246 )
247 lambdaStep = pexConfig.Field(
248 doc="Wavelength step for generating atmospheres (nm)",
249 dtype=float,
250 default=0.5,
251 )
252 lambdaRange = pexConfig.ListField(
253 doc="Wavelength range for LUT (Angstrom)",
254 dtype=float,
255 default=[3000.0, 11000.0],
256 )
257
258
259class FgcmMakeLutConfig(pipeBase.PipelineTaskConfig,
260 pipelineConnections=FgcmMakeLutConnections):
261 """Config for FgcmMakeLutTask"""
262 physicalFilters = pexConfig.ListField(
263 doc="List of physicalFilter labels to generate look-up table.",
264 dtype=str,
265 default=[],
266 )
267 stdPhysicalFilterOverrideMap = pexConfig.DictField(
268 doc=("Override mapping from physical filter labels to 'standard' physical "
269 "filter labels. The 'standard' physical filter defines the transmission "
270 "curve that the FGCM standard bandpass will be based on. "
271 "Any filter not listed here will be mapped to "
272 "itself (e.g. g->g or HSC-G->HSC-G). Use this override for cross-"
273 "filter calibration such as HSC-R->HSC-R2 and HSC-I->HSC-I2."),
274 keytype=str,
275 itemtype=str,
276 default={},
277 )
278 atmosphereTableName = pexConfig.Field(
279 doc="FGCM name or filename of precomputed atmospheres",
280 dtype=str,
281 default=None,
282 optional=True,
283 )
284 doOpticsTransmission = pexConfig.Field(
285 doc="Include optics transmission?",
286 dtype=bool,
287 default=True,
288 )
289 doSensorTransmission = pexConfig.Field(
290 doc="Include sensor transmission?",
291 dtype=bool,
292 default=True,
293 )
294 sensorCorrectionTermDict = pexConfig.ConfigDictField(
295 doc="Mapping of filter name to sensor correction terms.",
296 keytype=str,
297 itemtype=SensorCorrectionTerms,
298 default={},
299 )
300 parameters = pexConfig.ConfigField(
301 doc="Atmosphere parameters (required if no atmosphereTableName)",
302 dtype=FgcmMakeLutParametersConfig,
303 default=None,
304 check=None)
305
306 def validate(self):
307 """
308 Validate the config parameters.
309
310 This method behaves differently from the parent validate in the case
311 that atmosphereTableName is set. In this case, the config values
312 for standard values, step sizes, and ranges are loaded
313 directly from the specified atmosphereTableName.
314 """
315 # check that filterNames and stdFilterNames are okay
316 self._fields['physicalFilters'].validate(self)
317 self._fields['stdPhysicalFilterOverrideMap'].validate(self)
318
319 if self.atmosphereTableName is None:
320 # Validate the parameters
321 self._fields['parameters'].validate(self)
322
323
324class FgcmMakeLutTask(pipeBase.PipelineTask):
325 """
326 Make Look-Up Table for FGCM.
327
328 This task computes a look-up-table for the range in expected atmosphere
329 variation and variation in instrumental throughput (as tracked by the
330 transmission_filter products). By pre-computing linearized integrals,
331 the FGCM fit is orders of magnitude faster for stars with a broad range
332 of colors and observing bands, yielding precision at the 1-2 mmag level.
333
334 Computing a LUT requires running MODTRAN or with a pre-generated
335 atmosphere table packaged with fgcm.
336 """
337
338 ConfigClass = FgcmMakeLutConfig
339 _DefaultName = "fgcmMakeLut"
340
341 def __init__(self, initInputs=None, **kwargs):
342 super().__init__(**kwargs)
343
344 def runQuantum(self, butlerQC, inputRefs, outputRefs):
345 camera = butlerQC.get(inputRefs.camera)
346
347 if self.config.doOpticsTransmission:
348 opticsHandle = butlerQC.get(inputRefs.transmission_optics)
349 else:
350 opticsHandle = None
351
352 if self.config.doSensorTransmission:
353 sensorHandles = butlerQC.get(inputRefs.transmission_sensor)
354 sensorHandleDict = {sensorHandle.dataId['detector']: sensorHandle for
355 sensorHandle in sensorHandles}
356 else:
357 sensorHandleDict = {}
358
359 filterHandles = butlerQC.get(inputRefs.transmission_filter)
360 filterHandleDict = {filterHandle.dataId['physical_filter']: filterHandle for
361 filterHandle in filterHandles}
362
363 struct = self._fgcmMakeLut(camera,
364 opticsHandle,
365 sensorHandleDict,
366 filterHandleDict)
367
368 butlerQC.put(struct.fgcmLookUpTable, outputRefs.fgcmLookUpTable)
369 butlerQC.put(struct.fgcmStandardAtmosphere, outputRefs.fgcmStandardAtmosphere)
370
371 refDict = {passbandRef.dataId['physical_filter']: passbandRef for
372 passbandRef in outputRefs.fgcmStandardPassbands}
373 for physical_filter, passband in struct.fgcmStandardPassbands.items():
374 butlerQC.put(passband, refDict[physical_filter])
375
376 def _fgcmMakeLut(self, camera, opticsHandle, sensorHandleDict,
377 filterHandleDict):
378 """
379 Make a FGCM Look-up Table
380
381 Parameters
382 ----------
383 camera : `lsst.afw.cameraGeom.Camera`
384 Camera from the butler.
385 opticsHandle : `lsst.daf.butler.DeferredDatasetHandle`
386 Reference to optics transmission curve.
387 sensorHandleDict : `dict` of [`int`, `lsst.daf.butler.DeferredDatasetHandle`]
388 Dictionary of references to sensor transmission curves. Key will
389 be detector id.
390 filterHandleDict : `dict` of [`str`, `lsst.daf.butler.DeferredDatasetHandle`]
391 Dictionary of references to filter transmission curves. Key will
392 be physical filter label.
393
394 Returns
395 -------
396 retStruct : `lsst.pipe.base.Struct`
397 Output structure with keys:
398
399 fgcmLookUpTable : `BaseCatalog`
400 The FGCM look-up table.
401 fgcmStandardAtmosphere : `lsst.afw.image.TransmissionCurve`
402 Transmission curve for the FGCM standard atmosphere.
403 fgcmStandardPassbands : `dict` [`str`, `lsst.afw.image.TransmissionCurve`]
404 Dictionary of standard passbands, with the key as the
405 physical filter name.
406 """
407 # number of ccds from the length of the camera iterator
408 nCcd = len(camera)
409 self.log.info("Found %d ccds for look-up table" % (nCcd))
410
411 # Load in optics, etc.
412 self._loadThroughputs(camera,
413 opticsHandle,
414 sensorHandleDict,
415 filterHandleDict)
416
417 lutConfig = self._createLutConfig(nCcd)
418
419 # make the lut object
420 self.log.info("Making the LUT maker object")
421 self.fgcmLutMaker = fgcm.FgcmLUTMaker(lutConfig)
422
423 # generate the throughput dictionary.
424
425 # these will be in Angstroms
426 # note that lambdaStep is currently in nm, because of historical
427 # reasons in the code. Convert to Angstroms here.
428 throughputLambda = np.arange(self.fgcmLutMaker.lambdaRange[0],
429 self.fgcmLutMaker.lambdaRange[1]+self.fgcmLutMaker.lambdaStep*10,
430 self.fgcmLutMaker.lambdaStep*10.)
431
432 self.log.info("Built throughput lambda, %.1f-%.1f, step %.2f" %
433 (throughputLambda[0], throughputLambda[-1],
434 throughputLambda[1] - throughputLambda[0]))
435
436 throughputDict = {}
437 for i, physicalFilter in enumerate(self.config.physicalFilters):
438 tDict = {}
439 tDict['LAMBDA'] = throughputLambda
440 for ccdIndex, detector in enumerate(camera):
441 tDict[ccdIndex] = self._getThroughputDetector(detector, physicalFilter, throughputLambda)
442 throughputDict[physicalFilter] = tDict
443
444 # set the throughputs
445 self.fgcmLutMaker.setThroughputs(throughputDict)
446
447 # make the LUT
448 self.log.info("Making LUT")
449 self.fgcmLutMaker.makeLUT()
450
451 # and save the LUT
452
453 # build the index values
454 comma = ','
455 physicalFilterString = comma.join(self.config.physicalFilters)
456 stdPhysicalFilterString = comma.join(self._getStdPhysicalFilterList())
457
458 atmosphereTableName = 'NoTableWasUsed'
459 if self.config.atmosphereTableName is not None:
460 atmosphereTableName = self.config.atmosphereTableName
461
462 lutSchema = self._makeLutSchema(physicalFilterString, stdPhysicalFilterString,
463 atmosphereTableName)
464
465 lutCat = self._makeLutCat(lutSchema, physicalFilterString,
466 stdPhysicalFilterString, atmosphereTableName)
467
468 atmStd = TransmissionCurve.makeSpatiallyConstant(
469 throughput=self.fgcmLutMaker.atmStdTrans.astype(np.float64),
470 wavelengths=self.fgcmLutMaker.atmLambda.astype(np.float64),
471 throughputAtMin=self.fgcmLutMaker.atmStdTrans[0],
472 throughputAtMax=self.fgcmLutMaker.atmStdTrans[1],
473 )
474
475 fgcmStandardPassbands = {}
476 for i, physical_filter in enumerate(self.fgcmLutMaker.filterNames):
477 passband = self.fgcmLutMaker.throughputs[i]['THROUGHPUT_AVG']*self.fgcmLutMaker.atmStdTrans
478 fgcmStandardPassbands[physical_filter] = TransmissionCurve.makeSpatiallyConstant(
479 throughput=passband.astype(np.float64),
480 wavelengths=self.fgcmLutMaker.atmLambda.astype(np.float64),
481 throughputAtMin=passband[0],
482 throughputAtMax=passband[-1],
483 )
484
485 retStruct = pipeBase.Struct(
486 fgcmLookUpTable=lutCat,
487 fgcmStandardAtmosphere=atmStd,
488 fgcmStandardPassbands=fgcmStandardPassbands,
489 )
490
491 return retStruct
492
494 """Get the standard physical filter lists from config.physicalFilters
495 and config.stdPhysicalFilterOverrideMap
496
497 Returns
498 -------
499 stdPhysicalFilters : `list`
500 """
501 override = self.config.stdPhysicalFilterOverrideMap
502 return [override.get(physicalFilter, physicalFilter) for
503 physicalFilter in self.config.physicalFilters]
504
505 def _createLutConfig(self, nCcd):
506 """
507 Create the fgcmLut config dictionary
508
509 Parameters
510 ----------
511 nCcd: `int`
512 Number of CCDs in the camera
513 """
514
515 # create the common stub of the lutConfig
516 lutConfig = {}
517 lutConfig['logger'] = self.log
518 lutConfig['filterNames'] = self.config.physicalFilters
519 lutConfig['stdFilterNames'] = self._getStdPhysicalFilterList()
520 lutConfig['nCCD'] = nCcd
521
522 # atmosphereTable already validated if available
523 if self.config.atmosphereTableName is not None:
524 lutConfig['atmosphereTableName'] = self.config.atmosphereTableName
525 else:
526 # use the regular paramters (also validated if needed)
527 lutConfig['elevation'] = self.config.parameters.elevation
528 lutConfig['pmbRange'] = self.config.parameters.pmbRange
529 lutConfig['pmbSteps'] = self.config.parameters.pmbSteps
530 lutConfig['pwvRange'] = self.config.parameters.pwvRange
531 lutConfig['pwvSteps'] = self.config.parameters.pwvSteps
532 lutConfig['o3Range'] = self.config.parameters.o3Range
533 lutConfig['o3Steps'] = self.config.parameters.o3Steps
534 lutConfig['tauRange'] = self.config.parameters.tauRange
535 lutConfig['tauSteps'] = self.config.parameters.tauSteps
536 lutConfig['alphaRange'] = self.config.parameters.alphaRange
537 lutConfig['alphaSteps'] = self.config.parameters.alphaSteps
538 lutConfig['zenithRange'] = self.config.parameters.zenithRange
539 lutConfig['zenithSteps'] = self.config.parameters.zenithSteps
540 lutConfig['pmbStd'] = self.config.parameters.pmbStd
541 lutConfig['pwvStd'] = self.config.parameters.pwvStd
542 lutConfig['o3Std'] = self.config.parameters.o3Std
543 lutConfig['tauStd'] = self.config.parameters.tauStd
544 lutConfig['alphaStd'] = self.config.parameters.alphaStd
545 lutConfig['airmassStd'] = self.config.parameters.airmassStd
546 lutConfig['lambdaRange'] = self.config.parameters.lambdaRange
547 lutConfig['lambdaStep'] = self.config.parameters.lambdaStep
548 lutConfig['lambdaNorm'] = self.config.parameters.lambdaNorm
549
550 # Add any per-filter correction term updates if necessary.
551 # Note that sensorCTerms is the name of the config field in fgcm.
552 if self.config.sensorCorrectionTermDict:
553 lutConfig['sensorCTerms'] = {}
554 for key, value in self.config.sensorCorrectionTermDict.items():
555 lutConfig['sensorCTerms'][key] = (
556 value.refLambda,
557 dict(value.correctionTermDict),
558 )
559
560 return lutConfig
561
562 def _loadThroughputs(self, camera, opticsHandle, sensorHandleDict, filterHandleDict):
563 """Internal method to load throughput data for filters
564
565 Parameters
566 ----------
567 camera: `lsst.afw.cameraGeom.Camera`
568 Camera from the butler
569 opticsHandle : `lsst.daf.butler.DeferredDatasetHandle`
570 Reference to optics transmission curve.
571 sensorHandleDict : `dict` of [`int`, `lsst.daf.butler.DeferredDatasetHandle`]
572 Dictionary of references to sensor transmission curves. Key will
573 be detector id.
574 filterHandleDict : `dict` of [`str`, `lsst.daf.butler.DeferredDatasetHandle`]
575 Dictionary of references to filter transmission curves. Key will
576 be physical filter label.
577
578 Raises
579 ------
580 ValueError : Raised if configured filter name does not match any of the
581 available filter transmission curves.
582 """
583 if self.config.doOpticsTransmission:
584 self._opticsTransmission = opticsHandle.get()
585 else:
586 self._opticsTransmission = TransmissionCurve.makeSpatiallyConstant(
587 throughput=np.ones(100),
588 wavelengths=np.linspace(
589 self.config.parameters.lambdaRange[0],
590 self.config.parameters.lambdaRange[1],
591 100,
592 ),
593 throughputAtMin=1.0,
594 throughputAtMax=1.0,
595 )
596
598 for detector in camera:
599 if self.config.doSensorTransmission:
600 self._sensorsTransmission[detector.getId()] = sensorHandleDict[detector.getId()].get()
601 else:
602 self._sensorsTransmission[detector.getId()] = TransmissionCurve.makeSpatiallyConstant(
603 throughput=np.ones(100),
604 wavelengths=np.linspace(
605 self.config.parameters.lambdaRange[0],
606 self.config.parameters.lambdaRange[1],
607 100,
608 ),
609 throughputAtMin=1.0,
610 throughputAtMax=1.0,
611 )
612
614 for physicalFilter in self.config.physicalFilters:
615 self._filtersTransmission[physicalFilter] = filterHandleDict[physicalFilter].get()
616
617 def _getThroughputDetector(self, detector, physicalFilter, throughputLambda):
618 """Internal method to get throughput for a detector.
619
620 Returns the throughput at the center of the detector for a given filter.
621
622 Parameters
623 ----------
624 detector: `lsst.afw.cameraGeom._detector.Detector`
625 Detector on camera
626 physicalFilter: `str`
627 Physical filter label
628 throughputLambda: `np.array(dtype=np.float64)`
629 Wavelength steps (Angstrom)
630
631 Returns
632 -------
633 throughput: `np.array(dtype=np.float64)`
634 Throughput (max 1.0) at throughputLambda
635 """
636
637 c = detector.getCenter(afwCameraGeom.FOCAL_PLANE)
638 c.scale(1.0/detector.getPixelSize()[0]) # Assumes x and y pixel sizes in arcsec are the same
639
640 throughput = self._opticsTransmission.sampleAt(position=c,
641 wavelengths=throughputLambda)
642
643 throughput *= self._sensorsTransmission[detector.getId()].sampleAt(position=c,
644 wavelengths=throughputLambda)
645
646 throughput *= self._filtersTransmission[physicalFilter].sampleAt(position=c,
647 wavelengths=throughputLambda)
648
649 # Clip the throughput from 0 to 1
650 throughput = np.clip(throughput, 0.0, 1.0)
651
652 return throughput
653
654 def _makeLutSchema(self, physicalFilterString, stdPhysicalFilterString,
655 atmosphereTableName):
656 """
657 Make the LUT schema
658
659 Parameters
660 ----------
661 physicalFilterString: `str`
662 Combined string of all the physicalFilters
663 stdPhysicalFilterString: `str`
664 Combined string of all the standard physicalFilters
665 atmosphereTableName: `str`
666 Name of the atmosphere table used to generate LUT
667
668 Returns
669 -------
670 lutSchema: `afwTable.schema`
671 """
672
673 lutSchema = afwTable.Schema()
674
675 lutSchema.addField('tablename', type=str, doc='Atmosphere table name',
676 size=len(atmosphereTableName))
677 lutSchema.addField('elevation', type=float, doc="Telescope elevation used for LUT")
678 lutSchema.addField('physicalFilters', type=str, doc='physicalFilters in LUT',
679 size=len(physicalFilterString))
680 lutSchema.addField('stdPhysicalFilters', type=str, doc='Standard physicalFilters in LUT',
681 size=len(stdPhysicalFilterString))
682 lutSchema.addField('pmb', type='ArrayD', doc='Barometric Pressure',
683 size=self.fgcmLutMaker.pmb.size)
684 lutSchema.addField('pmbFactor', type='ArrayD', doc='PMB scaling factor',
685 size=self.fgcmLutMaker.pmb.size)
686 lutSchema.addField('pmbElevation', type=np.float64, doc='PMB Scaling at elevation')
687 lutSchema.addField('pwv', type='ArrayD', doc='Preciptable Water Vapor',
688 size=self.fgcmLutMaker.pwv.size)
689 lutSchema.addField('o3', type='ArrayD', doc='Ozone',
690 size=self.fgcmLutMaker.o3.size)
691 lutSchema.addField('tau', type='ArrayD', doc='Aerosol optical depth',
692 size=self.fgcmLutMaker.tau.size)
693 lutSchema.addField('lambdaNorm', type=np.float64, doc='AOD wavelength')
694 lutSchema.addField('alpha', type='ArrayD', doc='Aerosol alpha',
695 size=self.fgcmLutMaker.alpha.size)
696 lutSchema.addField('zenith', type='ArrayD', doc='Zenith angle',
697 size=self.fgcmLutMaker.zenith.size)
698 lutSchema.addField('nCcd', type=np.int32, doc='Number of CCDs')
699
700 # and the standard values
701 lutSchema.addField('pmbStd', type=np.float64, doc='PMB Standard')
702 lutSchema.addField('pwvStd', type=np.float64, doc='PWV Standard')
703 lutSchema.addField('o3Std', type=np.float64, doc='O3 Standard')
704 lutSchema.addField('tauStd', type=np.float64, doc='Tau Standard')
705 lutSchema.addField('alphaStd', type=np.float64, doc='Alpha Standard')
706 lutSchema.addField('zenithStd', type=np.float64, doc='Zenith angle Standard')
707 lutSchema.addField('lambdaRange', type='ArrayD', doc='Wavelength range',
708 size=2)
709 lutSchema.addField('lambdaStep', type=np.float64, doc='Wavelength step')
710 lutSchema.addField('lambdaStd', type='ArrayD', doc='Standard Wavelength',
711 size=len(self.fgcmLutMaker.filterNames))
712 lutSchema.addField('lambdaStdFilter', type='ArrayD', doc='Standard Wavelength (raw)',
713 size=len(self.fgcmLutMaker.filterNames))
714 lutSchema.addField('i0Std', type='ArrayD', doc='I0 Standard',
715 size=len(self.fgcmLutMaker.filterNames))
716 lutSchema.addField('i1Std', type='ArrayD', doc='I1 Standard',
717 size=len(self.fgcmLutMaker.filterNames))
718 lutSchema.addField('i10Std', type='ArrayD', doc='I10 Standard',
719 size=len(self.fgcmLutMaker.filterNames))
720 lutSchema.addField('i2Std', type='ArrayD', doc='I2 Standard',
721 size=len(self.fgcmLutMaker.filterNames))
722 lutSchema.addField('lambdaB', type='ArrayD', doc='Wavelength for passband (no atm)',
723 size=len(self.fgcmLutMaker.filterNames))
724 lutSchema.addField('atmLambda', type='ArrayD', doc='Atmosphere wavelengths (Angstrom)',
725 size=self.fgcmLutMaker.atmLambda.size)
726 lutSchema.addField('atmStdTrans', type='ArrayD', doc='Standard Atmosphere Throughput',
727 size=self.fgcmLutMaker.atmStdTrans.size)
728
729 # and the look-up-tables
730 lutSchema.addField('luttype', type=str, size=20, doc='Look-up table type')
731 lutSchema.addField('lut', type='ArrayF', doc='Look-up table for luttype',
732 size=self.fgcmLutMaker.lut['I0'].size)
733
734 return lutSchema
735
736 def _makeLutCat(self, lutSchema, physicalFilterString, stdPhysicalFilterString,
737 atmosphereTableName):
738 """
739 Make the LUT schema
740
741 Parameters
742 ----------
743 lutSchema: `afwTable.schema`
744 Lut catalog schema
745 physicalFilterString: `str`
746 Combined string of all the physicalFilters
747 stdPhysicalFilterString: `str`
748 Combined string of all the standard physicalFilters
749 atmosphereTableName: `str`
750 Name of the atmosphere table used to generate LUT
751
752 Returns
753 -------
754 lutCat: `afwTable.BaseCatalog`
755 Look-up table catalog for persistence.
756 """
757
758 # The somewhat strange format is to make sure that
759 # the rows of the afwTable do not get too large
760 # (see DM-11419)
761
762 lutCat = afwTable.BaseCatalog(lutSchema)
763 lutCat.table.preallocate(14)
764
765 # first fill the first index
766 rec = lutCat.addNew()
767
768 rec['tablename'] = atmosphereTableName
769 rec['elevation'] = self.fgcmLutMaker.atmosphereTable.elevation
770 rec['physicalFilters'] = physicalFilterString
771 rec['stdPhysicalFilters'] = stdPhysicalFilterString
772 rec['pmb'][:] = self.fgcmLutMaker.pmb
773 rec['pmbFactor'][:] = self.fgcmLutMaker.pmbFactor
774 rec['pmbElevation'] = self.fgcmLutMaker.pmbElevation
775 rec['pwv'][:] = self.fgcmLutMaker.pwv
776 rec['o3'][:] = self.fgcmLutMaker.o3
777 rec['tau'][:] = self.fgcmLutMaker.tau
778 rec['lambdaNorm'] = self.fgcmLutMaker.lambdaNorm
779 rec['alpha'][:] = self.fgcmLutMaker.alpha
780 rec['zenith'][:] = self.fgcmLutMaker.zenith
781 rec['nCcd'] = self.fgcmLutMaker.nCCD
782
783 rec['pmbStd'] = self.fgcmLutMaker.pmbStd
784 rec['pwvStd'] = self.fgcmLutMaker.pwvStd
785 rec['o3Std'] = self.fgcmLutMaker.o3Std
786 rec['tauStd'] = self.fgcmLutMaker.tauStd
787 rec['alphaStd'] = self.fgcmLutMaker.alphaStd
788 rec['zenithStd'] = self.fgcmLutMaker.zenithStd
789 rec['lambdaRange'][:] = self.fgcmLutMaker.lambdaRange
790 rec['lambdaStep'] = self.fgcmLutMaker.lambdaStep
791 rec['lambdaStd'][:] = self.fgcmLutMaker.lambdaStd
792 rec['lambdaStdFilter'][:] = self.fgcmLutMaker.lambdaStdFilter
793 rec['i0Std'][:] = self.fgcmLutMaker.I0Std
794 rec['i1Std'][:] = self.fgcmLutMaker.I1Std
795 rec['i10Std'][:] = self.fgcmLutMaker.I10Std
796 rec['i2Std'][:] = self.fgcmLutMaker.I2Std
797 rec['lambdaB'][:] = self.fgcmLutMaker.lambdaB
798 rec['atmLambda'][:] = self.fgcmLutMaker.atmLambda
799 rec['atmStdTrans'][:] = self.fgcmLutMaker.atmStdTrans
800
801 rec['luttype'] = 'I0'
802 rec['lut'][:] = self.fgcmLutMaker.lut['I0'].flatten()
803
804 # and add the rest
805 rec = lutCat.addNew()
806 rec['luttype'] = 'I1'
807 rec['lut'][:] = self.fgcmLutMaker.lut['I1'].flatten()
808
809 derivTypes = ['D_PMB', 'D_LNPWV', 'D_O3', 'D_LNTAU', 'D_ALPHA', 'D_SECZENITH',
810 'D_PMB_I1', 'D_LNPWV_I1', 'D_O3_I1', 'D_LNTAU_I1', 'D_ALPHA_I1',
811 'D_SECZENITH_I1']
812 for derivType in derivTypes:
813 rec = lutCat.addNew()
814 rec['luttype'] = derivType
815 rec['lut'][:] = self.fgcmLutMaker.lutDeriv[derivType].flatten()
816
817 return lutCat
Defines the fields and offsets for a table.
Definition Schema.h:51
_makeLutSchema(self, physicalFilterString, stdPhysicalFilterString, atmosphereTableName)
_fgcmMakeLut(self, camera, opticsHandle, sensorHandleDict, filterHandleDict)
_getThroughputDetector(self, detector, physicalFilter, throughputLambda)
_loadThroughputs(self, camera, opticsHandle, sensorHandleDict, filterHandleDict)
runQuantum(self, butlerQC, inputRefs, outputRefs)
__init__(self, initInputs=None, **kwargs)
_makeLutCat(self, lutSchema, physicalFilterString, stdPhysicalFilterString, atmosphereTableName)