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
ptcDataset.py
Go to the documentation of this file.
2# LSST Data Management System
3# Copyright 2008-2017 AURA/LSST.
4#
5# This product includes software developed by the
6# LSST Project (http://www.lsst.org/).
7#
8# This program is free software: you can redistribute it and/or modify
9# it under the terms of the GNU General Public License as published by
10# the Free Software Foundation, either version 3 of the License, or
11# (at your option) any later version.
12#
13# This program is distributed in the hope that it will be useful,
14# but WITHOUT ANY WARRANTY; without even the implied warranty of
15# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16# GNU General Public License for more details.
17#
18# You should have received a copy of the LSST License Statement and
19# the GNU General Public License along with this program. If not,
20# see <https://www.lsstcorp.org/LegalNotices/>.
21#
22"""
23Define dataset class for MeasurePhotonTransferCurve task
24"""
25
26__all__ = ['PhotonTransferCurveDataset']
27
28import numpy as np
29from astropy.table import Table
30
31from lsst.ip.isr import IsrCalib
32
33
35 """A simple class to hold the output data from the PTC task.
36
37 The dataset is made up of a dictionary for each item, keyed by the
38 amplifiers' names, which much be supplied at construction time.
39 New items cannot be added to the class to save accidentally saving to the
40 wrong property, and the class can be frozen if desired.
41 inputExpIdPairs records the exposures used to produce the data.
42 When fitPtc() or fitCovariancesAstier() is run, a mask is built up, which
43 is by definition always the same length as inputExpIdPairs, rawExpTimes,
44 rawMeans and rawVars, and is a list of bools, which are incrementally set
45 to False as points are discarded from the fits.
46 PTC fit parameters for polynomials are stored in a list in ascending order
47 of polynomial term, i.e. par[0]*x^0 + par[1]*x + par[2]*x^2 etc
48 with the length of the list corresponding to the order of the polynomial
49 plus one.
50
51 Parameters
52 ----------
53 ampNames : `list`
54 List with the names of the amplifiers of the detector at hand.
55 ptcFitType : `str`
56 Type of model fitted to the PTC: "POLYNOMIAL", "EXPAPPROXIMATION",
57 or "FULLCOVARIANCE".
58 covMatrixSide : `int`
59 Maximum lag of covariances (size of square covariance matrices).
60 kwargs : `dict`, optional
61 Other keyword arguments to pass to the parent init.
62
63 Notes
64 -----
65 The stored attributes are:
66
67 badAmps : `list`
68 List with bad amplifiers names.
69 inputExpIdPairs : `dict`, [`str`, `list`]
70 Dictionary keyed by amp names containing the input exposures IDs.
71 expIdMask : `dict`, [`str`, `list`]
72 Dictionary keyed by amp names containing the mask produced after
73 outlier rejection. The mask produced by the "FULLCOVARIANCE"
74 option may differ from the one produced in the other two PTC
75 fit types.
76 rawExpTimes : `dict`, [`str`, `list`]
77 Dictionary keyed by amp names containing the unmasked exposure times.
78 rawMeans : `dict`, [`str`, `list`]
79 Dictionary keyed by amp namescontaining the unmasked average of the
80 means of the exposures in each flat pair.
81 rawVars : `dict`, [`str`, `list`]
82 Dictionary keyed by amp names containing the variance of the
83 difference image of the exposures in each flat pair.
84 gain : `dict`, [`str`, `list`]
85 Dictionary keyed by amp names containing the fitted gains.
86 gainErr : `dict`, [`str`, `list`]
87 Dictionary keyed by amp names containing the errors on the
88 fitted gains.
89 noise : `dict`, [`str`, `list`]
90 Dictionary keyed by amp names containing the fitted noise.
91 noiseErr : `dict`, [`str`, `list`]
92 Dictionary keyed by amp names containing the errors on the fitted
93 noise.
94 ptcFitPars : `dict`, [`str`, `list`]
95 Dictionary keyed by amp names containing the fitted parameters of the
96 PTC model for ptcFitTye in ["POLYNOMIAL", "EXPAPPROXIMATION"].
97 ptcFitParsError : `dict`, [`str`, `list`]
98 Dictionary keyed by amp names containing the errors on the fitted
99 parameters of the PTC model for ptcFitTye in
100 ["POLYNOMIAL", "EXPAPPROXIMATION"].
101 ptcFitChiSq : `dict`, [`str`, `list`]
102 Dictionary keyed by amp names containing the reduced chi squared
103 of the fit for ptcFitTye in ["POLYNOMIAL", "EXPAPPROXIMATION"].
104 ptcTurnoff : `float`
105 Flux value (in ADU) where the variance of the PTC curve starts
106 decreasing consistently.
107 covariances : `dict`, [`str`, `list`]
108 Dictionary keyed by amp names containing a list of measured
109 covariances per mean flux.
110 covariancesModel : `dict`, [`str`, `list`]
111 Dictionary keyed by amp names containinging covariances model
112 (Eq. 20 of Astier+19) per mean flux.
113 covariancesSqrtWeights : `dict`, [`str`, `list`]
114 Dictionary keyed by amp names containinging sqrt. of covariances
115 weights.
116 aMatrix : `dict`, [`str`, `list`]
117 Dictionary keyed by amp names containing the "a" parameters from
118 the model in Eq. 20 of Astier+19.
119 bMatrix : `dict`, [`str`, `list`]
120 Dictionary keyed by amp names containing the "b" parameters from
121 the model in Eq. 20 of Astier+19.
122 covariancesModelNoB : `dict`, [`str`, `list`]
123 Dictionary keyed by amp names containing covariances model
124 (with 'b'=0 in Eq. 20 of Astier+19)
125 per mean flux.
126 aMatrixNoB : `dict`, [`str`, `list`]
127 Dictionary keyed by amp names containing the "a" parameters from the
128 model in Eq. 20 of Astier+19
129 (and 'b' = 0).
130 finalVars : `dict`, [`str`, `list`]
131 Dictionary keyed by amp names containing the masked variance of the
132 difference image of each flat
133 pair. If needed, each array will be right-padded with
134 np.nan to match the length of rawExpTimes.
135 finalModelVars : `dict`, [`str`, `list`]
136 Dictionary keyed by amp names containing the masked modeled
137 variance of the difference image of each flat pair. If needed, each
138 array will be right-padded with np.nan to match the length of
139 rawExpTimes.
140 finalMeans : `dict`, [`str`, `list`]
141 Dictionary keyed by amp names containing the masked average of the
142 means of the exposures in each flat pair. If needed, each array
143 will be right-padded with np.nan to match the length of
144 rawExpTimes.
145 photoCharge : `dict`, [`str`, `list`]
146 Dictionary keyed by amp names containing the integrated photocharge
147 for linearity calibration.
148
149 Version 1.1 adds the `ptcTurnoff` attribute.
150 """
151
152 _OBSTYPE = 'PTC'
153 _SCHEMA = 'Gen3 Photon Transfer Curve'
154 _VERSION = 1.1
155
156 def __init__(self, ampNames=[], ptcFitType=None, covMatrixSide=1, **kwargs):
157 self.ptcFitType = ptcFitType
158 self.ampNames = ampNames
159 self.covMatrixSide = covMatrixSide
160
161 self.badAmps = [np.nan]
162
163 self.inputExpIdPairs = {ampName: [] for ampName in ampNames}
164 self.expIdMask = {ampName: [] for ampName in ampNames}
165 self.rawExpTimes = {ampName: [] for ampName in ampNames}
166 self.rawMeans = {ampName: [] for ampName in ampNames}
167 self.rawVars = {ampName: [] for ampName in ampNames}
168 self.photoCharge = {ampName: [] for ampName in ampNames}
169
170 self.gain = {ampName: np.nan for ampName in ampNames}
171 self.gainErr = {ampName: np.nan for ampName in ampNames}
172 self.noise = {ampName: np.nan for ampName in ampNames}
173 self.noiseErr = {ampName: np.nan for ampName in ampNames}
174
175 self.ptcFitPars = {ampName: [] for ampName in ampNames}
176 self.ptcFitParsError = {ampName: [] for ampName in ampNames}
177 self.ptcFitChiSq = {ampName: np.nan for ampName in ampNames}
178 self.ptcTurnoff = {ampName: np.nan for ampName in ampNames}
179
180 self.covariances = {ampName: [] for ampName in ampNames}
181 self.covariancesModel = {ampName: [] for ampName in ampNames}
182 self.covariancesSqrtWeights = {ampName: [] for ampName in ampNames}
183 self.aMatrix = {ampName: np.nan for ampName in ampNames}
184 self.bMatrix = {ampName: np.nan for ampName in ampNames}
185 self.covariancesModelNoB = {ampName: [] for ampName in ampNames}
186 self.aMatrixNoB = {ampName: np.nan for ampName in ampNames}
187
188 self.finalVars = {ampName: [] for ampName in ampNames}
189 self.finalModelVars = {ampName: [] for ampName in ampNames}
190 self.finalMeans = {ampName: [] for ampName in ampNames}
191
192 super().__init__(**kwargs)
193 self.requiredAttributesrequiredAttributesrequiredAttributes.update(['badAmps', 'inputExpIdPairs', 'expIdMask', 'rawExpTimes',
194 'rawMeans', 'rawVars', 'gain', 'gainErr', 'noise', 'noiseErr',
195 'ptcFitPars', 'ptcFitParsError', 'ptcFitChiSq', 'ptcTurnoff',
196 'aMatrixNoB', 'covariances', 'covariancesModel',
197 'covariancesSqrtWeights', 'covariancesModelNoB',
198 'aMatrix', 'bMatrix', 'finalVars', 'finalModelVars', 'finalMeans',
199 'photoCharge'])
200
201 self.updateMetadataupdateMetadata(setCalibInfo=True, setCalibId=True, **kwargs)
202
203 def setAmpValues(self, ampName, inputExpIdPair=[(np.nan, np.nan)], expIdMask=[np.nan],
204 rawExpTime=[np.nan], rawMean=[np.nan], rawVar=[np.nan], photoCharge=[np.nan],
205 gain=np.nan, gainErr=np.nan, noise=np.nan, noiseErr=np.nan, ptcFitPars=[np.nan],
206 ptcFitParsError=[np.nan], ptcFitChiSq=np.nan, ptcTurnoff=np.nan, covArray=[],
207 covArrayModel=[], covSqrtWeights=[], aMatrix=[], bMatrix=[], covArrayModelNoB=[],
208 aMatrixNoB=[], finalVar=[np.nan], finalModelVar=[np.nan], finalMean=[np.nan]):
209 """Function to initialize an amp of a PhotonTransferCurveDataset.
210
211 Notes
212 -----
213 The parameters are all documented in `init`.
214 """
215 nanMatrix = np.full((self.covMatrixSide, self.covMatrixSide), np.nan)
216 if len(covArray) == 0:
217 covArray = [nanMatrix]
218 if len(covArrayModel) == 0:
219 covArrayModel = [nanMatrix]
220 if len(covSqrtWeights) == 0:
221 covSqrtWeights = [nanMatrix]
222 if len(covArrayModelNoB) == 0:
223 covArrayModelNoB = [nanMatrix]
224 if len(aMatrix) == 0:
225 aMatrix = nanMatrix
226 if len(bMatrix) == 0:
227 bMatrix = nanMatrix
228 if len(aMatrixNoB) == 0:
229 aMatrixNoB = nanMatrix
230
231 self.inputExpIdPairs[ampName] = inputExpIdPair
232 self.expIdMask[ampName] = expIdMask
233 self.rawExpTimes[ampName] = rawExpTime
234 self.rawMeans[ampName] = rawMean
235 self.rawVars[ampName] = rawVar
236 self.photoCharge[ampName] = photoCharge
237 self.gain[ampName] = gain
238 self.gainErr[ampName] = gainErr
239 self.noise[ampName] = noise
240 self.noiseErr[ampName] = noiseErr
241 self.ptcFitPars[ampName] = ptcFitPars
242 self.ptcFitParsError[ampName] = ptcFitParsError
243 self.ptcFitChiSq[ampName] = ptcFitChiSq
244 self.ptcTurnoff[ampName] = ptcTurnoff
245 self.covariances[ampName] = covArray
246 self.covariancesSqrtWeights[ampName] = covSqrtWeights
247 self.covariancesModel[ampName] = covArrayModel
248 self.covariancesModelNoB[ampName] = covArrayModelNoB
249 self.aMatrix[ampName] = aMatrix
250 self.bMatrix[ampName] = bMatrix
251 self.aMatrixNoB[ampName] = aMatrixNoB
252 self.ptcFitPars[ampName] = ptcFitPars
253 self.ptcFitParsError[ampName] = ptcFitParsError
254 self.ptcFitChiSq[ampName] = ptcFitChiSq
255 self.finalVars[ampName] = finalVar
256 self.finalModelVars[ampName] = finalModelVar
257 self.finalMeans[ampName] = finalMean
258
259 def updateMetadata(self, **kwargs):
260 """Update calibration metadata.
261 This calls the base class's method after ensuring the required
262 calibration keywords will be saved.
263
264 Parameters
265 ----------
266 setDate : `bool`, optional
267 Update the CALIBDATE fields in the metadata to the current
268 time. Defaults to False.
269 kwargs :
270 Other keyword parameters to set in the metadata.
271 """
272 super().updateMetadata(PTC_FIT_TYPE=self.ptcFitType, **kwargs)
273
274 @classmethod
275 def fromDict(cls, dictionary):
276 """Construct a calibration from a dictionary of properties.
277 Must be implemented by the specific calibration subclasses.
278
279 Parameters
280 ----------
281 dictionary : `dict`
282 Dictionary of properties.
283
284 Returns
285 -------
287 Constructed calibration.
288
289 Raises
290 ------
291 RuntimeError
292 Raised if the supplied dictionary is for a different
293 calibration.
294 """
295 calib = cls()
296 if calib._OBSTYPE != dictionary['metadata']['OBSTYPE']:
297 raise RuntimeError(f"Incorrect Photon Transfer Curve dataset supplied. "
298 f"Expected {calib._OBSTYPE}, found {dictionary['metadata']['OBSTYPE']}")
299 calib.setMetadata(dictionary['metadata'])
300 calib.ptcFitType = dictionary['ptcFitType']
301 calib.covMatrixSide = dictionary['covMatrixSide']
302 calib.badAmps = np.array(dictionary['badAmps'], 'str').tolist()
303
304 # The cov matrices are square
305 covMatrixSide = calib.covMatrixSide
306 # Number of final signal levels
307 covDimensionsProduct = len(np.array(list(dictionary['covariances'].values())[0]).ravel())
308 nSignalPoints = int(covDimensionsProduct/(covMatrixSide*covMatrixSide))
309
310 for ampName in dictionary['ampNames']:
311 covsAmp = np.array(dictionary['covariances'][ampName]).reshape((nSignalPoints, covMatrixSide,
312 covMatrixSide))
313
314 # After cpPtcExtract runs in the PTC pipeline, the datasets
315 # created ('PARTIAL' and 'DUMMY') have a single measurement.
316 # Apply the maskign to the final ptcDataset, after running
317 # cpPtcSolve.
318 if len(covsAmp) > 1:
319 # Masks for covariances padding in `toTable`
320 maskCovsAmp = np.array([~np.isnan(entry).all() for entry in covsAmp])
321 maskAmp = ~np.isnan(np.array(dictionary['finalMeans'][ampName]))
322 else:
323 maskCovsAmp = np.array([True])
324 maskAmp = np.array([True])
325
326 calib.ampNames.append(ampName)
327 calib.inputExpIdPairs[ampName] = np.array(dictionary['inputExpIdPairs'][ampName]).tolist()
328 calib.expIdMask[ampName] = np.array(dictionary['expIdMask'][ampName]).tolist()
329 calib.rawExpTimes[ampName] = np.array(dictionary['rawExpTimes'][ampName]).tolist()
330 calib.rawMeans[ampName] = np.array(dictionary['rawMeans'][ampName]).tolist()
331 calib.rawVars[ampName] = np.array(dictionary['rawVars'][ampName]).tolist()
332 calib.gain[ampName] = np.array(dictionary['gain'][ampName]).tolist()
333 calib.gainErr[ampName] = np.array(dictionary['gainErr'][ampName]).tolist()
334 calib.noise[ampName] = np.array(dictionary['noise'][ampName]).tolist()
335 calib.noiseErr[ampName] = np.array(dictionary['noiseErr'][ampName]).tolist()
336 calib.ptcFitPars[ampName] = np.array(dictionary['ptcFitPars'][ampName]).tolist()
337 calib.ptcFitParsError[ampName] = np.array(dictionary['ptcFitParsError'][ampName]).tolist()
338 calib.ptcFitChiSq[ampName] = np.array(dictionary['ptcFitChiSq'][ampName]).tolist()
339 calib.ptcTurnoff[ampName] = np.array(dictionary['ptcTurnoff'][ampName]).tolist()
340 calib.covariances[ampName] = covsAmp[maskCovsAmp].tolist()
341 calib.covariancesModel[ampName] = np.array(
342 dictionary['covariancesModel'][ampName]).reshape(
343 (nSignalPoints, covMatrixSide, covMatrixSide))[maskCovsAmp].tolist()
344 calib.covariancesSqrtWeights[ampName] = np.array(
345 dictionary['covariancesSqrtWeights'][ampName]).reshape(
346 (nSignalPoints, covMatrixSide, covMatrixSide))[maskCovsAmp].tolist()
347 calib.aMatrix[ampName] = np.array(dictionary['aMatrix'][ampName]).reshape(
348 (covMatrixSide, covMatrixSide)).tolist()
349 calib.bMatrix[ampName] = np.array(dictionary['bMatrix'][ampName]).reshape(
350 (covMatrixSide, covMatrixSide)).tolist()
351 calib.covariancesModelNoB[ampName] = np.array(
352 dictionary['covariancesModelNoB'][ampName]).reshape(
353 (nSignalPoints, covMatrixSide, covMatrixSide))[maskCovsAmp].tolist()
354 calib.aMatrixNoB[ampName] = np.array(
355 dictionary['aMatrixNoB'][ampName]).reshape((covMatrixSide, covMatrixSide)).tolist()
356 calib.finalVars[ampName] = np.array(dictionary['finalVars'][ampName])[maskAmp].tolist()
357 calib.finalModelVars[ampName] = np.array(dictionary['finalModelVars'][ampName])[maskAmp].tolist()
358 calib.finalMeans[ampName] = np.array(dictionary['finalMeans'][ampName])[maskAmp].tolist()
359 calib.photoCharge[ampName] = np.array(dictionary['photoCharge'][ampName]).tolist()
360 calib.updateMetadata()
361 return calib
362
363 def toDict(self):
364 """Return a dictionary containing the calibration properties.
365 The dictionary should be able to be round-tripped through
366 `fromDict`.
367
368 Returns
369 -------
370 dictionary : `dict`
371 Dictionary of properties.
372 """
374
375 outDict = dict()
376 metadata = self.getMetadata()
377 outDict['metadata'] = metadata
378
379 outDict['ptcFitType'] = self.ptcFitType
380 outDict['covMatrixSide'] = self.covMatrixSide
381 outDict['ampNames'] = self.ampNames
382 outDict['badAmps'] = self.badAmps
383 outDict['inputExpIdPairs'] = self.inputExpIdPairs
384 outDict['expIdMask'] = self.expIdMask
385 outDict['rawExpTimes'] = self.rawExpTimes
386 outDict['rawMeans'] = self.rawMeans
387 outDict['rawVars'] = self.rawVars
388 outDict['gain'] = self.gain
389 outDict['gainErr'] = self.gainErr
390 outDict['noise'] = self.noise
391 outDict['noiseErr'] = self.noiseErr
392 outDict['ptcFitPars'] = self.ptcFitPars
393 outDict['ptcFitParsError'] = self.ptcFitParsError
394 outDict['ptcFitChiSq'] = self.ptcFitChiSq
395 outDict['ptcTurnoff'] = self.ptcTurnoff
396 outDict['covariances'] = self.covariances
397 outDict['covariancesModel'] = self.covariancesModel
398 outDict['covariancesSqrtWeights'] = self.covariancesSqrtWeights
399 outDict['aMatrix'] = self.aMatrix
400 outDict['bMatrix'] = self.bMatrix
401 outDict['covariancesModelNoB'] = self.covariancesModelNoB
402 outDict['aMatrixNoB'] = self.aMatrixNoB
403 outDict['finalVars'] = self.finalVars
404 outDict['finalModelVars'] = self.finalModelVars
405 outDict['finalMeans'] = self.finalMeans
406 outDict['photoCharge'] = self.photoCharge
407
408 return outDict
409
410 @classmethod
411 def fromTable(cls, tableList):
412 """Construct calibration from a list of tables.
413 This method uses the `fromDict` method to create the
414 calibration, after constructing an appropriate dictionary from
415 the input tables.
416
417 Parameters
418 ----------
419 tableList : `list` [`lsst.afw.table.Table`]
420 List of tables to use to construct the datasetPtc.
421
422 Returns
423 -------
425 The calibration defined in the tables.
426 """
427 ptcTable = tableList[0]
428
429 metadata = ptcTable.meta
430 inDict = dict()
431 inDict['metadata'] = metadata
432 inDict['ampNames'] = []
433 inDict['ptcFitType'] = []
434 inDict['covMatrixSide'] = []
435 inDict['inputExpIdPairs'] = dict()
436 inDict['expIdMask'] = dict()
437 inDict['rawExpTimes'] = dict()
438 inDict['rawMeans'] = dict()
439 inDict['rawVars'] = dict()
440 inDict['gain'] = dict()
441 inDict['gainErr'] = dict()
442 inDict['noise'] = dict()
443 inDict['noiseErr'] = dict()
444 inDict['ptcFitPars'] = dict()
445 inDict['ptcFitParsError'] = dict()
446 inDict['ptcFitChiSq'] = dict()
447 inDict['ptcTurnoff'] = dict()
448 inDict['covariances'] = dict()
449 inDict['covariancesModel'] = dict()
450 inDict['covariancesSqrtWeights'] = dict()
451 inDict['aMatrix'] = dict()
452 inDict['bMatrix'] = dict()
453 inDict['covariancesModelNoB'] = dict()
454 inDict['aMatrixNoB'] = dict()
455 inDict['finalVars'] = dict()
456 inDict['finalModelVars'] = dict()
457 inDict['finalMeans'] = dict()
458 inDict['badAmps'] = []
459 inDict['photoCharge'] = dict()
460
461 calibVersion = metadata['PTC_VERSION']
462 if calibVersion == 1.0:
463 cls().log.warning(f"Previous version found for PTC dataset: {calibVersion}. "
464 f"Setting 'ptcTurnoff' in all amps to last value in 'finalMeans'.")
465 for record in ptcTable:
466 ampName = record['AMPLIFIER_NAME']
467
468 inDict['ptcFitType'] = record['PTC_FIT_TYPE']
469 inDict['covMatrixSide'] = record['COV_MATRIX_SIDE']
470 inDict['ampNames'].append(ampName)
471 inDict['inputExpIdPairs'][ampName] = record['INPUT_EXP_ID_PAIRS']
472 inDict['expIdMask'][ampName] = record['EXP_ID_MASK']
473 inDict['rawExpTimes'][ampName] = record['RAW_EXP_TIMES']
474 inDict['rawMeans'][ampName] = record['RAW_MEANS']
475 inDict['rawVars'][ampName] = record['RAW_VARS']
476 inDict['gain'][ampName] = record['GAIN']
477 inDict['gainErr'][ampName] = record['GAIN_ERR']
478 inDict['noise'][ampName] = record['NOISE']
479 inDict['noiseErr'][ampName] = record['NOISE_ERR']
480 inDict['ptcFitPars'][ampName] = record['PTC_FIT_PARS']
481 inDict['ptcFitParsError'][ampName] = record['PTC_FIT_PARS_ERROR']
482 inDict['ptcFitChiSq'][ampName] = record['PTC_FIT_CHI_SQ']
483 inDict['covariances'][ampName] = record['COVARIANCES']
484 inDict['covariancesModel'][ampName] = record['COVARIANCES_MODEL']
485 inDict['covariancesSqrtWeights'][ampName] = record['COVARIANCES_SQRT_WEIGHTS']
486 inDict['aMatrix'][ampName] = record['A_MATRIX']
487 inDict['bMatrix'][ampName] = record['B_MATRIX']
488 inDict['covariancesModelNoB'][ampName] = record['COVARIANCES_MODEL_NO_B']
489 inDict['aMatrixNoB'][ampName] = record['A_MATRIX_NO_B']
490 inDict['finalVars'][ampName] = record['FINAL_VARS']
491 inDict['finalModelVars'][ampName] = record['FINAL_MODEL_VARS']
492 inDict['finalMeans'][ampName] = record['FINAL_MEANS']
493 inDict['badAmps'] = record['BAD_AMPS']
494 inDict['photoCharge'][ampName] = record['PHOTO_CHARGE']
495 if calibVersion == 1.0:
496 mask = record['FINAL_MEANS'].mask
497 array = record['FINAL_MEANS'][~mask]
498 if len(array) > 0:
499 inDict['ptcTurnoff'][ampName] = record['FINAL_MEANS'][~mask][-1]
500 else:
501 inDict['ptcTurnoff'][ampName] = np.nan
502 else:
503 inDict['ptcTurnoff'][ampName] = record['PTC_TURNOFF']
504 return cls().fromDict(inDict)
505
506 def toTable(self):
507 """Construct a list of tables containing the information in this
508 calibration.
509
510 The list of tables should create an identical calibration
511 after being passed to this class's fromTable method.
512
513 Returns
514 -------
515 tableList : `list` [`astropy.table.Table`]
516 List of tables containing the linearity calibration
517 information.
518 """
519 tableList = []
521 nPoints = []
522 for i, ampName in enumerate(self.ampNames):
523 nPoints.append(len(list(self.covariances.values())[i]))
524 nSignalPoints = max(nPoints)
525 nPadPoints = {}
526 for i, ampName in enumerate(self.ampNames):
527 nPadPoints[ampName] = nSignalPoints - len(list(self.covariances.values())[i])
528 covMatrixSide = self.covMatrixSide
529
530 catalog = Table([{'AMPLIFIER_NAME': ampName,
531 'PTC_FIT_TYPE': self.ptcFitType,
532 'COV_MATRIX_SIDE': self.covMatrixSide,
533 'INPUT_EXP_ID_PAIRS': self.inputExpIdPairs[ampName]
534 if len(self.expIdMask[ampName]) else np.nan,
535 'EXP_ID_MASK': self.expIdMask[ampName]
536 if len(self.expIdMask[ampName]) else np.nan,
537 'RAW_EXP_TIMES': np.array(self.rawExpTimes[ampName]).tolist()
538 if len(self.rawExpTimes[ampName]) else np.nan,
539 'RAW_MEANS': np.array(self.rawMeans[ampName]).tolist()
540 if len(self.rawMeans[ampName]) else np.nan,
541 'RAW_VARS': np.array(self.rawVars[ampName]).tolist()
542 if len(self.rawVars[ampName]) else np.nan,
543 'GAIN': self.gain[ampName],
544 'GAIN_ERR': self.gainErr[ampName],
545 'NOISE': self.noise[ampName],
546 'NOISE_ERR': self.noiseErr[ampName],
547 'PTC_FIT_PARS': np.array(self.ptcFitPars[ampName]).tolist(),
548 'PTC_FIT_PARS_ERROR': np.array(self.ptcFitParsError[ampName]).tolist(),
549 'PTC_FIT_CHI_SQ': self.ptcFitChiSq[ampName],
550 'PTC_TURNOFF': self.ptcTurnoff[ampName],
551 'COVARIANCES': np.pad(np.array(self.covariances[ampName]),
552 ((0, nPadPoints[ampName]), (0, 0), (0, 0)),
553 'constant', constant_values=np.nan).reshape(
554 nSignalPoints*covMatrixSide**2).tolist(),
555 'COVARIANCES_MODEL': np.pad(np.array(self.covariancesModel[ampName]),
556 ((0, nPadPoints[ampName]), (0, 0), (0, 0)),
557 'constant', constant_values=np.nan).reshape(
558 nSignalPoints*covMatrixSide**2).tolist(),
559 'COVARIANCES_SQRT_WEIGHTS': np.pad(np.array(self.covariancesSqrtWeights[ampName]),
560 ((0, nPadPoints[ampName]), (0, 0), (0, 0)),
561 'constant', constant_values=0.0).reshape(
562 nSignalPoints*covMatrixSide**2).tolist(),
563 'A_MATRIX': np.array(self.aMatrix[ampName]).reshape(covMatrixSide**2).tolist(),
564 'B_MATRIX': np.array(self.bMatrix[ampName]).reshape(covMatrixSide**2).tolist(),
565 'COVARIANCES_MODEL_NO_B':
566 np.pad(np.array(self.covariancesModelNoB[ampName]),
567 ((0, nPadPoints[ampName]), (0, 0), (0, 0)),
568 'constant', constant_values=np.nan).reshape(
569 nSignalPoints*covMatrixSide**2).tolist(),
570 'A_MATRIX_NO_B': np.array(self.aMatrixNoB[ampName]).reshape(
571 covMatrixSide**2).tolist(),
572 'FINAL_VARS': np.pad(np.array(self.finalVars[ampName]), (0, nPadPoints[ampName]),
573 'constant', constant_values=np.nan).tolist(),
574 'FINAL_MODEL_VARS': np.pad(np.array(self.finalModelVars[ampName]),
575 (0, nPadPoints[ampName]),
576 'constant', constant_values=np.nan).tolist(),
577 'FINAL_MEANS': np.pad(np.array(self.finalMeans[ampName]),
578 (0, nPadPoints[ampName]),
579 'constant', constant_values=np.nan).tolist(),
580 'BAD_AMPS': np.array(self.badAmps).tolist() if len(self.badAmps) else np.nan,
581 'PHOTO_CHARGE': np.array(self.photoCharge[ampName]).tolist(),
582 } for ampName in self.ampNames])
583 inMeta = self.getMetadata().toDict()
584 outMeta = {k: v for k, v in inMeta.items() if v is not None}
585 outMeta.update({k: "" for k, v in inMeta.items() if v is None})
586 catalog.meta = outMeta
587 tableList.append(catalog)
588
589 return(tableList)
590
591 def fromDetector(self, detector):
592 """Read metadata parameters from a detector.
593
594 Parameters
595 ----------
596 detector : `lsst.afw.cameraGeom.detector`
597 Input detector with parameters to use.
598
599 Returns
600 -------
602 The calibration constructed from the detector.
603 """
604
605 pass
606
607 def getExpIdsUsed(self, ampName):
608 """Get the exposures used, i.e. not discarded, for a given amp.
609 If no mask has been created yet, all exposures are returned.
610 """
611 if len(self.expIdMask[ampName]) == 0:
612 return self.inputExpIdPairs[ampName]
613
614 # if the mask exists it had better be the same length as the expIdPairs
615 assert len(self.expIdMask[ampName]) == len(self.inputExpIdPairs[ampName])
616
617 pairs = self.inputExpIdPairs[ampName]
618 mask = self.expIdMask[ampName]
619 # cast to bool required because numpy
620 return [(exp1, exp2) for ((exp1, exp2), m) in zip(pairs, mask) if bool(m) is True]
621
622 def getGoodAmps(self):
623 return [amp for amp in self.ampNames if amp not in self.badAmps]
int max
table::Key< int > to
def requiredAttributes(self, value)
Definition: calibType.py:149
def updateMetadata(self, camera=None, detector=None, filterName=None, setCalibId=False, setCalibInfo=False, setDate=False, **kwargs)
Definition: calibType.py:188
def __init__(self, ampNames=[], ptcFitType=None, covMatrixSide=1, **kwargs)
Definition: ptcDataset.py:156
def setAmpValues(self, ampName, inputExpIdPair=[(np.nan, np.nan)], expIdMask=[np.nan], rawExpTime=[np.nan], rawMean=[np.nan], rawVar=[np.nan], photoCharge=[np.nan], gain=np.nan, gainErr=np.nan, noise=np.nan, noiseErr=np.nan, ptcFitPars=[np.nan], ptcFitParsError=[np.nan], ptcFitChiSq=np.nan, ptcTurnoff=np.nan, covArray=[], covArrayModel=[], covSqrtWeights=[], aMatrix=[], bMatrix=[], covArrayModelNoB=[], aMatrixNoB=[], finalVar=[np.nan], finalModelVar=[np.nan], finalMean=[np.nan])
Definition: ptcDataset.py:208
daf::base::PropertyList * list
Definition: fits.cc:928