LSSTApplications  18.0.0+106,18.0.0+50,19.0.0,19.0.0+1,19.0.0+10,19.0.0+11,19.0.0+13,19.0.0+17,19.0.0+2,19.0.0-1-g20d9b18+6,19.0.0-1-g425ff20,19.0.0-1-g5549ca4,19.0.0-1-g580fafe+6,19.0.0-1-g6fe20d0+1,19.0.0-1-g7011481+9,19.0.0-1-g8c57eb9+6,19.0.0-1-gb5175dc+11,19.0.0-1-gdc0e4a7+9,19.0.0-1-ge272bc4+6,19.0.0-1-ge3aa853,19.0.0-10-g448f008b,19.0.0-12-g6990b2c,19.0.0-2-g0d9f9cd+11,19.0.0-2-g3d9e4fb2+11,19.0.0-2-g5037de4,19.0.0-2-gb96a1c4+3,19.0.0-2-gd955cfd+15,19.0.0-3-g2d13df8,19.0.0-3-g6f3c7dc,19.0.0-4-g725f80e+11,19.0.0-4-ga671dab3b+1,19.0.0-4-gad373c5+3,19.0.0-5-ga2acb9c+2,19.0.0-5-gfe96e6c+2,w.2020.01
LSSTDataManagementBasePackage
simple_curve.py
Go to the documentation of this file.
1 #
2 # LSST Data Management System
3 #
4 # Copyright 2019 AURA/LSST.
5 #
6 # This product includes software developed by the
7 # LSST Project (http://www.lsst.org/).
8 #
9 # This program is free software: you can redistribute it and/or modify
10 # it under the terms of the GNU General Public License as published by
11 # the Free Software Foundation, either version 3 of the License, or
12 # (at your option) any later version.
13 #
14 # This program is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 # GNU General Public License for more details.
18 #
19 # You should have received a copy of the LSST License Statement and
20 # the GNU General Public License along with this program. If not,
21 # see <https://www.lsstcorp.org/LegalNotices/>.
22 #
23 
24 __all__ = ["Curve", "AmpCurve", "DetectorCurve", "ImageCurve"]
25 
26 from scipy.interpolate import interp1d
27 from astropy.table import QTable
28 import astropy.units as u
29 from abc import ABC, abstractmethod
30 import datetime
31 import os
32 import numpy
33 
34 import lsst.afw.cameraGeom.utils as cgUtils
35 from lsst.geom import Point2I
36 
37 
38 class Curve(ABC):
39  """ An abstract class to represent an arbitrary curve with
40  interpolation.
41  """
42  mode = ''
43  subclasses = dict()
44 
45  def __init__(self, wavelength, efficiency, metadata):
46  if not (isinstance(wavelength, u.Quantity) and wavelength.unit.physical_type == 'length'):
47  raise ValueError('The wavelength must be a quantity with a length sense.')
48  if not isinstance(efficiency, u.Quantity) or efficiency.unit != u.percent:
49  raise ValueError('The efficiency must be a quantity with units of percent.')
50  self.wavelength = wavelength
51  self.efficiency = efficiency
52  # make sure needed metadata is set if built directly from ctor.
53  metadata.update({'MODE': self.mode, 'TYPE': 'QE'})
54  self.metadata = metadata
55 
56  @classmethod
57  @abstractmethod
58  def fromTable(cls, table):
59  """Class method for constructing a `Curve` object.
60 
61  Parameters
62  ----------
63  table : `astropy.table.QTable`
64  Table containing metadata and columns necessary
65  for constructing a `Curve` object.
66 
67  Returns
68  -------
69  curve : `Curve`
70  A `Curve` subclass of the appropriate type according
71  to the table metadata
72  """
73  pass
74 
75  @abstractmethod
76  def toTable(self):
77  """Convert this `Curve` object to an `astropy.table.QTable`.
78 
79  Returns
80  -------
81  table : `astropy.table.QTable`
82  A table object containing the data from this `Curve`.
83  """
84  pass
85 
86  @abstractmethod
87  def evaluate(self, detector, position, wavelength, kind='linear'):
88  """Interpolate the curve at the specified position and wavelength.
89 
90  Parameters
91  ----------
92  detector : `lsst.afw.cameraGeom.Detector`
93  Is used to find the appropriate curve given the position for
94  curves that vary over the detector. Ignored in the case where
95  there is only a single curve per detector.
96  position : `lsst.geom.Point2D`
97  The position on the detector at which to evaluate the curve.
98  wavelength : `astropy.units.Quantity`
99  The wavelength(s) at which to make the interpolation.
100  kind : `str`
101  The type of interpolation to do. See documentation for
102  `scipy.interpolate.interp1d` for accepted values.
103 
104  Returns
105  -------
106  value : `astropy.units.Quantity`
107  Interpolated value(s). Number of values returned will match the
108  length of `wavelength`.
109  """
110  pass
111 
112  @classmethod
113  def __init_subclass__(cls, **kwargs):
114  """Register subclasses with the abstract base class"""
115  super().__init_subclass__(**kwargs)
116  if cls.mode in Curve.subclasses:
117  raise ValueError(f'Class for mode, {cls.mode}, already defined')
118  Curve.subclasses[cls.mode] = cls
119 
120  @abstractmethod
121  def __eq__(self, other):
122  """Define equality for this class"""
123  pass
124 
125  def compare_metadata(self, other,
126  keys_to_compare=['MODE', 'TYPE', 'CALIBDATE', 'INSTRUME', 'OBSTYPE', 'DETECTOR']):
127  """Compare metadata in this object to another.
128 
129  Parameters
130  ----------
131  other : `Curve`
132  The object with which to compare metadata.
133  keys_to_compare : `list`
134  List of metadata keys to compare.
135 
136  Returns
137  -------
138  same : `bool`
139  Are the metadata the same?
140  """
141  for k in keys_to_compare:
142  if self.metadata[k] != other.metadata[k]:
143  return False
144  return True
145 
146  def interpolate(self, wavelengths, values, wavelength, kind):
147  """Interplate the curve at the specified wavelength(s).
148 
149  Parameters
150  ----------
151  wavelengths : `astropy.units.Quantity`
152  The wavelength values for the curve.
153  values : `astropy.units.Quantity`
154  The y-values for the curve.
155  wavelength : `astropy.units.Quantity`
156  The wavelength(s) at which to make the interpolation.
157  kind : `str`
158  The type of interpolation to do. See documentation for
159  `scipy.interpolate.interp1d` for accepted values.
160 
161  Returns
162  -------
163  value : `astropy.units.Quantity`
164  Interpolated value(s)
165  """
166  if not isinstance(wavelength, u.Quantity):
167  raise ValueError("Wavelengths at which to interpolate must be astropy quantities")
168  if not (isinstance(wavelengths, u.Quantity) and isinstance(values, u.Quantity)):
169  raise ValueError("Model to be interpreted must be astropy quantities")
170  interp_wavelength = wavelength.to(wavelengths.unit)
171  f = interp1d(wavelengths, values, kind=kind)
172  return f(interp_wavelength.value)*values.unit
173 
174  def getMetadata(self):
175  """Return metadata
176 
177  Returns
178  -------
179  metadata : `dict`
180  Dictionary of metadata for this curve.
181  """
182  # Needed to duck type as an object that can be ingested
183  return self.metadata
184 
185  @classmethod
186  def readText(cls, filename):
187  """Class method for constructing a `Curve` object from
188  the standardized text format.
189 
190  Parameters
191  ----------
192  filename : `str`
193  Path to the text file to read.
194 
195  Returns
196  -------
197  curve : `Curve`
198  A `Curve` subclass of the appropriate type according
199  to the table metadata
200  """
201  table = QTable.read(filename, format='ascii.ecsv')
202  return cls.subclasses[table.meta['MODE']].fromTable(table)
203 
204  @classmethod
205  def readFits(cls, filename):
206  """Class method for constructing a `Curve` object from
207  the standardized FITS format.
208 
209  Parameters
210  ----------
211  filename : `str`
212  Path to the FITS file to read.
213 
214  Returns
215  -------
216  curve : `Curve`
217  A `Curve` subclass of the appropriate type according
218  to the table metadata
219  """
220  table = QTable.read(filename, format='fits')
221  return cls.subclasses[table.meta['MODE']].fromTable(table)
222 
223  @staticmethod
224  def _check_cols(cols, table):
225  """Check that the columns are in the table"""
226  for col in cols:
227  if col not in table.columns:
228  raise ValueError(f'The table must include a column named "{col}".')
229 
230  def _to_table_with_meta(self):
231  """Compute standard metadata before writing file out"""
232  now = datetime.datetime.utcnow()
233  table = self.toTable()
234  metadata = table.meta
235  metadata["DATE"] = now.isoformat()
236  metadata["CALIB_CREATION_DATE"] = now.strftime("%Y-%m-%d")
237  metadata["CALIB_CREATION_TIME"] = now.strftime("%T %Z").strip()
238  return table
239 
240  def writeText(self, filename):
241  """ Write the `Curve` out to a text file.
242 
243  Parameters
244  ----------
245  filename : `str`
246  Path to the text file to write.
247 
248  Returns
249  -------
250  filename : `str`
251  Because this method forces a particular extension return
252  the name of the file actually written.
253  """
254  table = self._to_table_with_meta()
255  # Force file extension to .ecsv
256  path, ext = os.path.splitext(filename)
257  filename = path + ".ecsv"
258  table.write(filename, format="ascii.ecsv")
259  return filename
260 
261  def writeFits(self, filename):
262  """ Write the `Curve` out to a FITS file.
263 
264  Parameters
265  ----------
266  filename : `str`
267  Path to the FITS file to write.
268 
269  Returns
270  -------
271  filename : `str`
272  Because this method forces a particular extension return
273  the name of the file actually written.
274  """
275  table = self._to_table_with_meta()
276  # Force file extension to .ecsv
277  path, ext = os.path.splitext(filename)
278  filename = path + ".fits"
279  table.write(filename, format="fits")
280  return filename
281 
282 
284  """Subclass of `Curve` that represents a single curve per detector.
285 
286  Parameters
287  ----------
288  wavelength : `astropy.units.Quantity`
289  Wavelength values for this curve
290  efficiency : `astropy.units.Quantity`
291  Quantum efficiency values for this curve
292  metadata : `dict`
293  Dictionary of metadata for this curve
294  """
295  mode = 'DETECTOR'
296 
297  def __eq__(self, other):
298  return (self.compare_metadata(other) and
299  numpy.array_equal(self.wavelength, other.wavelength) and
300  numpy.array_equal(self.wavelength, other.wavelength))
301 
302  @classmethod
303  def fromTable(cls, table):
304  # Docstring inherited from base classs
305  cls._check_cols(['wavelength', 'efficiency'], table)
306  return cls(table['wavelength'], table['efficiency'], table.meta)
307 
308  def toTable(self):
309  # Docstring inherited from base classs
310  return QTable({'wavelength': self.wavelength, 'efficiency': self.efficiency}, meta=self.metadata)
311 
312  def evaluate(self, detector, position, wavelength, kind='linear'):
313  # Docstring inherited from base classs
314  return self.interpolate(self.wavelength, self.efficiency, wavelength, kind=kind)
315 
316 
318  """Subclass of `Curve` that represents a curve per amp.
319 
320  Parameters
321  ----------
322  amp_name_list : iterable of `str`
323  The name of the amp for each entry
324  wavelength : `astropy.units.Quantity`
325  Wavelength values for this curve
326  efficiency : `astropy.units.Quantity`
327  Quantum efficiency values for this curve
328  metadata : `dict`
329  Dictionary of metadata for this curve
330  """
331  mode = 'AMP'
332 
333  def __init__(self, amp_name_list, wavelength, efficiency, metadata):
334  super().__init__(wavelength, efficiency, metadata)
335  amp_names = set(amp_name_list)
336  self.data = {}
337  for amp_name in amp_names:
338  idx = numpy.where(amp_name_list == amp_name)[0]
339  # Deal with the case where the keys are bytes from FITS
340  name = amp_name
341  if isinstance(name, bytes):
342  name = name.decode()
343  self.data[name] = (wavelength[idx], efficiency[idx])
344 
345  def __eq__(self, other):
346  if not self.compare_metadata(other):
347  return False
348  for k in self.data:
349  if not numpy.array_equal(self.data[k][0], other.data[k][0]):
350  return False
351  if not numpy.array_equal(self.data[k][1], other.data[k][1]):
352  return False
353  return True
354 
355  @classmethod
356  def fromTable(cls, table):
357  # Docstring inherited from base classs
358  cls._check_cols(['amp_name', 'wavelength', 'efficiency'], table)
359  return cls(table['amp_name'], table['wavelength'],
360  table['efficiency'], table.meta)
361 
362  def toTable(self):
363  # Docstring inherited from base classs
364  wavelength = None
365  efficiency = None
366  names = numpy.array([])
367  # Loop over the amps and concatenate into three same length columns to feed
368  # to the Table constructor.
369  for amp_name, val in self.data.items():
370  # This will preserve the quantity
371  if wavelength is None:
372  wunit = val[0].unit
373  wavelength = val[0].value
374  else:
375  wavelength = numpy.concatenate([wavelength, val[0].value])
376  if efficiency is None:
377  eunit = val[1].unit
378  efficiency = val[1].value
379  else:
380  efficiency = numpy.concatenate([efficiency, val[1].value])
381  names = numpy.concatenate([names, numpy.full(val[0].shape, amp_name)])
382  names = numpy.array(names)
383  # Note that in future, the astropy.unit should make it through concatenation
384  return QTable({'amp_name': names, 'wavelength': wavelength*wunit, 'efficiency': efficiency*eunit},
385  meta=self.metadata)
386 
387  def evaluate(self, detector, position, wavelength, kind='linear'):
388  # Docstring inherited from base classs
389  amp = cgUtils.findAmp(detector, Point2I(position)) # cast to Point2I if Point2D passed
390  w, e = self.data[amp.getName()]
391  return self.interpolate(w, e, wavelength, kind=kind)
392 
393 
395  mode = 'IMAGE'
396 
397  def fromTable(self, table):
398  # Docstring inherited from base classs
399  raise NotImplementedError()
400 
401  def toTable(self):
402  # Docstring inherited from base classs
403  raise NotImplementedError()
404 
405  def evaluate(self, detector, position, wavelength, kind='linear'):
406  # Docstring inherited from base classs
407  raise NotImplementedError()
def evaluate(self, detector, position, wavelength, kind='linear')
def evaluate(self, detector, position, wavelength, kind='linear')
Definition: simple_curve.py:87
def compare_metadata(self, other, keys_to_compare=['MODE', TYPE, CALIBDATE, INSTRUME, OBSTYPE, DETECTOR)
def evaluate(self, detector, position, wavelength, kind='linear')
std::vector< SchemaItem< Flag > > * items
daf::base::PropertySet * set
Definition: fits.cc:902
def __init__(self, wavelength, efficiency, metadata)
Definition: simple_curve.py:45
def evaluate(self, detector, position, wavelength, kind='linear')
Point< int, 2 > Point2I
Definition: Point.h:321
def __init__(self, amp_name_list, wavelength, efficiency, metadata)
def interpolate(self, wavelengths, values, wavelength, kind)
bool strip
Definition: fits.cc:901