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