lsst.ip.diffim  18.0.0-4-gd60ad06
dcrModel.py
Go to the documentation of this file.
1 # This file is part of ip_diffim.
2 #
3 # LSST Data Management System
4 # This product includes software developed by the
5 # LSST Project (http://www.lsst.org/).
6 # See COPYRIGHT file at the top of the source tree.
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 
23 import numpy as np
24 from scipy import ndimage
25 from lsst.afw.coord.refraction import differentialRefraction
26 import lsst.afw.geom as afwGeom
27 import lsst.afw.image as afwImage
28 from lsst.geom import radians
29 
30 __all__ = ["DcrModel", "applyDcr", "calculateDcr", "calculateImageParallacticAngle"]
31 
32 
33 class DcrModel:
34  """A model of the true sky after correcting chromatic effects.
35 
36  Attributes
37  ----------
38  dcrNumSubfilters : `int`
39  Number of sub-filters used to model chromatic effects within a band.
40  modelImages : `list` of `lsst.afw.image.Image`
41  A list of masked images, each containing the model for one subfilter
42 
43  Notes
44  -----
45  The ``DcrModel`` contains an estimate of the true sky, at a higher
46  wavelength resolution than the input observations. It can be forward-
47  modeled to produce Differential Chromatic Refraction (DCR) matched
48  templates for a given ``Exposure``, and provides utilities for conditioning
49  the model in ``dcrAssembleCoadd`` to avoid oscillating solutions between
50  iterations of forward modeling or between the subfilters of the model.
51  """
52 
53  def __init__(self, modelImages, filterInfo=None, psf=None, mask=None, variance=None):
54  self.dcrNumSubfilters = len(modelImages)
55  self.modelImages = modelImages
56  self._filter = filterInfo
57  self._psf = psf
58  self._mask = mask
59  self._variance = variance
60 
61  @classmethod
62  def fromImage(cls, maskedImage, dcrNumSubfilters, filterInfo=None, psf=None):
63  """Initialize a DcrModel by dividing a coadd between the subfilters.
64 
65  Parameters
66  ----------
67  maskedImage : `lsst.afw.image.MaskedImage`
68  Input coadded image to divide equally between the subfilters.
69  dcrNumSubfilters : `int`
70  Number of sub-filters used to model chromatic effects within a band.
71  filterInfo : `lsst.afw.image.Filter`, optional
72  The filter definition, set in the current instruments' obs package.
73  Required for any calculation of DCR, including making matched templates.
74  psf : `lsst.afw.detection.Psf`, optional
75  Point spread function (PSF) of the model.
76  Required if the ``DcrModel`` will be persisted.
77 
78  Returns
79  -------
80  dcrModel : `lsst.pipe.tasks.DcrModel`
81  Best fit model of the true sky after correcting chromatic effects.
82 
83  Raises
84  ------
85  ValueError
86  If there are any unmasked NAN values in ``maskedImage``.
87  """
88  # NANs will potentially contaminate the entire image,
89  # depending on the shift or convolution type used.
90  model = maskedImage.image.clone()
91  mask = maskedImage.mask.clone()
92  # We divide the variance by N and not N**2 because we will assume each
93  # subfilter is independent. That means that the significance of
94  # detected sources will be lower by a factor of sqrt(N) in the
95  # subfilter images, but we will recover it when we combine the
96  # subfilter images to construct matched templates.
97  variance = maskedImage.variance.clone()
98  variance /= dcrNumSubfilters
99  model /= dcrNumSubfilters
100  modelImages = [model, ]
101  for subfilter in range(1, dcrNumSubfilters):
102  modelImages.append(model.clone())
103  return cls(modelImages, filterInfo, psf, mask, variance)
104 
105  @classmethod
106  def fromDataRef(cls, dataRef, datasetType="dcrCoadd", numSubfilters=None, **kwargs):
107  """Load an existing DcrModel from a repository.
108 
109  Parameters
110  ----------
111  dataRef : `lsst.daf.persistence.ButlerDataRef`
112  Data reference defining the patch for coaddition and the
113  reference Warp
114  datasetType : `str`, optional
115  Name of the DcrModel in the registry {"dcrCoadd", "dcrCoadd_sub"}
116  numSubfilters : `int`
117  Number of sub-filters used to model chromatic effects within a band.
118  **kwargs
119  Additional keyword arguments to pass to look up the model in the data registry.
120  Common keywords and their types include: ``tract``:`str`, ``patch``:`str`,
121  ``bbox``:`lsst.afw.geom.Box2I`
122 
123  Returns
124  -------
125  dcrModel : `lsst.pipe.tasks.DcrModel`
126  Best fit model of the true sky after correcting chromatic effects.
127  """
128  modelImages = []
129  filterInfo = None
130  psf = None
131  mask = None
132  variance = None
133  for subfilter in range(numSubfilters):
134  dcrCoadd = dataRef.get(datasetType, subfilter=subfilter,
135  numSubfilters=numSubfilters, **kwargs)
136  if filterInfo is None:
137  filterInfo = dcrCoadd.getFilter()
138  if psf is None:
139  psf = dcrCoadd.getPsf()
140  if mask is None:
141  mask = dcrCoadd.mask
142  if variance is None:
143  variance = dcrCoadd.variance
144  modelImages.append(dcrCoadd.image)
145  return cls(modelImages, filterInfo, psf, mask, variance)
146 
147  def __len__(self):
148  """Return the number of subfilters.
149 
150  Returns
151  -------
152  dcrNumSubfilters : `int`
153  The number of DCR subfilters in the model.
154  """
155  return self.dcrNumSubfilters
156 
157  def __getitem__(self, subfilter):
158  """Iterate over the subfilters of the DCR model.
159 
160  Parameters
161  ----------
162  subfilter : `int`
163  Index of the current ``subfilter`` within the full band.
164  Negative indices are allowed, and count in reverse order
165  from the highest ``subfilter``.
166 
167  Returns
168  -------
169  modelImage : `lsst.afw.image.Image`
170  The DCR model for the given ``subfilter``.
171 
172  Raises
173  ------
174  IndexError
175  If the requested ``subfilter`` is greater or equal to the number
176  of subfilters in the model.
177  """
178  if np.abs(subfilter) >= len(self):
179  raise IndexError("subfilter out of bounds.")
180  return self.modelImages[subfilter]
181 
182  def __setitem__(self, subfilter, maskedImage):
183  """Update the model image for one subfilter.
184 
185  Parameters
186  ----------
187  subfilter : `int`
188  Index of the current subfilter within the full band.
189  maskedImage : `lsst.afw.image.Image`
190  The DCR model to set for the given ``subfilter``.
191 
192  Raises
193  ------
194  IndexError
195  If the requested ``subfilter`` is greater or equal to the number
196  of subfilters in the model.
197  ValueError
198  If the bounding box of the new image does not match.
199  """
200  if np.abs(subfilter) >= len(self):
201  raise IndexError("subfilter out of bounds.")
202  if maskedImage.getBBox() != self.bbox:
203  raise ValueError("The bounding box of a subfilter must not change.")
204  self.modelImages[subfilter] = maskedImage
205 
206  @property
207  def filter(self):
208  """Return the filter of the model.
209 
210  Returns
211  -------
212  filter : `lsst.afw.image.Filter`
213  The filter definition, set in the current instruments' obs package.
214  """
215  return self._filter
216 
217  @property
218  def psf(self):
219  """Return the psf of the model.
220 
221  Returns
222  -------
223  psf : `lsst.afw.detection.Psf`
224  Point spread function (PSF) of the model.
225  """
226  return self._psf
227 
228  @property
229  def bbox(self):
230  """Return the common bounding box of each subfilter image.
231 
232  Returns
233  -------
234  bbox : `lsst.afw.geom.Box2I`
235  Bounding box of the DCR model.
236  """
237  return self[0].getBBox()
238 
239  @property
240  def mask(self):
241  """Return the common mask of each subfilter image.
242 
243  Returns
244  -------
245  mask : `lsst.afw.image.Mask`
246  Mask plane of the DCR model.
247  """
248  return self._mask
249 
250  @property
251  def variance(self):
252  """Return the common variance of each subfilter image.
253 
254  Returns
255  -------
256  variance : `lsst.afw.image.Image`
257  Variance plane of the DCR model.
258  """
259  return self._variance
260 
261  def getReferenceImage(self, bbox=None):
262  """Calculate a reference image from the average of the subfilter images.
263 
264  Parameters
265  ----------
266  bbox : `lsst.afw.geom.Box2I`, optional
267  Sub-region of the coadd. Returns the entire image if `None`.
268 
269  Returns
270  -------
271  refImage : `numpy.ndarray`
272  The reference image with no chromatic effects applied.
273  """
274  bbox = bbox or self.bbox
275  return np.mean([model[bbox].array for model in self], axis=0)
276 
277  def assign(self, dcrSubModel, bbox=None):
278  """Update a sub-region of the ``DcrModel`` with new values.
279 
280  Parameters
281  ----------
282  dcrSubModel : `lsst.pipe.tasks.DcrModel`
283  New model of the true scene after correcting chromatic effects.
284  bbox : `lsst.afw.geom.Box2I`, optional
285  Sub-region of the coadd.
286  Defaults to the bounding box of ``dcrSubModel``.
287 
288  Raises
289  ------
290  ValueError
291  If the new model has a different number of subfilters.
292  """
293  if len(dcrSubModel) != len(self):
294  raise ValueError("The number of DCR subfilters must be the same "
295  "between the old and new models.")
296  bbox = bbox or self.bbox
297  for model, subModel in zip(self, dcrSubModel):
298  model.assign(subModel[bbox], bbox)
299 
300  def buildMatchedTemplate(self, exposure=None, order=3,
301  visitInfo=None, bbox=None, wcs=None, mask=None,
302  splitSubfilters=True, splitThreshold=0., amplifyModel=1.):
303  """Create a DCR-matched template image for an exposure.
304 
305  Parameters
306  ----------
307  exposure : `lsst.afw.image.Exposure`, optional
308  The input exposure to build a matched template for.
309  May be omitted if all of the metadata is supplied separately
310  order : `int`, optional
311  Interpolation order of the DCR shift.
312  visitInfo : `lsst.afw.image.VisitInfo`, optional
313  Metadata for the exposure. Ignored if ``exposure`` is set.
314  bbox : `lsst.afw.geom.Box2I`, optional
315  Sub-region of the coadd. Ignored if ``exposure`` is set.
316  wcs : `lsst.afw.geom.SkyWcs`, optional
317  Coordinate system definition (wcs) for the exposure.
318  Ignored if ``exposure`` is set.
319  mask : `lsst.afw.image.Mask`, optional
320  reference mask to use for the template image.
321  splitSubfilters : `bool`, optional
322  Calculate DCR for two evenly-spaced wavelengths in each subfilter,
323  instead of at the midpoint. Default: True
324  splitThreshold : `float`, optional
325  Minimum DCR difference within a subfilter required to use ``splitSubfilters``
326  amplifyModel : `float`, optional
327  Multiplication factor to amplify differences between model planes.
328  Used to speed convergence of iterative forward modeling.
329 
330  Returns
331  -------
332  templateImage : `lsst.afw.image.ImageF`
333  The DCR-matched template
334 
335  Raises
336  ------
337  ValueError
338  If neither ``exposure`` or all of ``visitInfo``, ``bbox``, and ``wcs`` are set.
339  """
340  if self.filter is None:
341  raise ValueError("'filterInfo' must be set for the DcrModel in order to calculate DCR.")
342  if exposure is not None:
343  visitInfo = exposure.getInfo().getVisitInfo()
344  bbox = exposure.getBBox()
345  wcs = exposure.getInfo().getWcs()
346  elif visitInfo is None or bbox is None or wcs is None:
347  raise ValueError("Either exposure or visitInfo, bbox, and wcs must be set.")
348  dcrShift = calculateDcr(visitInfo, wcs, self.filter, len(self), splitSubfilters=splitSubfilters)
349  templateImage = afwImage.ImageF(bbox)
350  refModel = self.getReferenceImage(bbox)
351  for subfilter, dcr in enumerate(dcrShift):
352  if amplifyModel > 1:
353  model = (self[subfilter][bbox].array - refModel)*amplifyModel + refModel
354  else:
355  model = self[subfilter][bbox].array
356  templateImage.array += applyDcr(model, dcr, splitSubfilters=splitSubfilters,
357  splitThreshold=splitThreshold, order=order)
358  return templateImage
359 
360  def buildMatchedExposure(self, exposure=None,
361  visitInfo=None, bbox=None, wcs=None, mask=None):
362  """Wrapper to create an exposure from a template image.
363 
364  Parameters
365  ----------
366  exposure : `lsst.afw.image.Exposure`, optional
367  The input exposure to build a matched template for.
368  May be omitted if all of the metadata is supplied separately
369  visitInfo : `lsst.afw.image.VisitInfo`, optional
370  Metadata for the exposure. Ignored if ``exposure`` is set.
371  bbox : `lsst.afw.geom.Box2I`, optional
372  Sub-region of the coadd. Ignored if ``exposure`` is set.
373  wcs : `lsst.afw.geom.SkyWcs`, optional
374  Coordinate system definition (wcs) for the exposure.
375  Ignored if ``exposure`` is set.
376  mask : `lsst.afw.image.Mask`, optional
377  reference mask to use for the template image.
378 
379  Returns
380  -------
381  templateExposure : `lsst.afw.image.exposureF`
382  The DCR-matched template
383  """
384  if bbox is None:
385  bbox = exposure.getBBox()
386  templateImage = self.buildMatchedTemplate(exposure=exposure, visitInfo=visitInfo,
387  bbox=bbox, wcs=wcs, mask=mask)
388  maskedImage = afwImage.MaskedImageF(bbox)
389  maskedImage.image = templateImage[bbox]
390  maskedImage.mask = self.mask[bbox]
391  maskedImage.variance = self.variance[bbox]
392  templateExposure = afwImage.ExposureF(bbox, wcs)
393  templateExposure.setMaskedImage(maskedImage[bbox])
394  templateExposure.setPsf(self.psf)
395  templateExposure.setFilter(self.filter)
396  return templateExposure
397 
398  def conditionDcrModel(self, modelImages, bbox, gain=1.):
399  """Average two iterations' solutions to reduce oscillations.
400 
401  Parameters
402  ----------
403  modelImages : `list` of `lsst.afw.image.Image`
404  The new DCR model images from the current iteration.
405  The values will be modified in place.
406  bbox : `lsst.afw.geom.Box2I`
407  Sub-region of the coadd
408  gain : `float`, optional
409  Relative weight to give the new solution when updating the model.
410  Defaults to 1.0, which gives equal weight to both solutions.
411  """
412  # Calculate weighted averages of the images.
413  for model, newModel in zip(self, modelImages):
414  newModel *= gain
415  newModel += model[bbox]
416  newModel /= 1. + gain
417 
418  def regularizeModelIter(self, subfilter, newModel, bbox, regularizationFactor,
419  regularizationWidth=2):
420  """Restrict large variations in the model between iterations.
421 
422  Parameters
423  ----------
424  subfilter : `int`
425  Index of the current subfilter within the full band.
426  newModel : `lsst.afw.image.Image`
427  The new DCR model for one subfilter from the current iteration.
428  Values in ``newModel`` that are extreme compared with the last
429  iteration are modified in place.
430  bbox : `lsst.afw.geom.Box2I`
431  Sub-region to coadd
432  regularizationFactor : `float`
433  Maximum relative change of the model allowed between iterations.
434  regularizationWidth : int, optional
435  Minimum radius of a region to include in regularization, in pixels.
436  """
437  refImage = self[subfilter][bbox].array
438  highThreshold = np.abs(refImage)*regularizationFactor
439  lowThreshold = refImage/regularizationFactor
440  newImage = newModel.array
441  self.applyImageThresholds(newImage, highThreshold=highThreshold, lowThreshold=lowThreshold,
442  regularizationWidth=regularizationWidth)
443 
444  def regularizeModelFreq(self, modelImages, bbox, statsCtrl, regularizationFactor,
445  regularizationWidth=2, mask=None, convergenceMaskPlanes="DETECTED"):
446  """Restrict large variations in the model between subfilters.
447 
448  Parameters
449  ----------
450  modelImages : `list` of `lsst.afw.image.Image`
451  The new DCR model images from the current iteration.
452  The values will be modified in place.
453  bbox : `lsst.afw.geom.Box2I`
454  Sub-region to coadd
455  statsCtrl : `lsst.afw.math.StatisticsControl`
456  Statistics control object for coaddition.
457  regularizationFactor : `float`
458  Maximum relative change of the model allowed between subfilters.
459  regularizationWidth : `int`, optional
460  Minimum radius of a region to include in regularization, in pixels.
461  mask : `lsst.afw.image.Mask`, optional
462  Optional alternate mask
463  convergenceMaskPlanes : `list` of `str`, or `str`, optional
464  Mask planes to use to calculate convergence.
465 
466  Notes
467  -----
468  This implementation of frequency regularization restricts each subfilter
469  image to be a smoothly-varying function times a reference image.
470  """
471  # ``regularizationFactor`` is the maximum change between subfilter images, so the maximum difference
472  # between one subfilter image and the average will be the square root of that.
473  maxDiff = np.sqrt(regularizationFactor)
474  noiseLevel = self.calculateNoiseCutoff(modelImages[0], statsCtrl, bufferSize=5, mask=mask, bbox=bbox)
475  referenceImage = self.getReferenceImage(bbox)
476  badPixels = np.isnan(referenceImage) | (referenceImage <= 0.)
477  if np.sum(~badPixels) == 0:
478  # Skip regularization if there are no valid pixels
479  return
480  referenceImage[badPixels] = 0.
481  filterWidth = regularizationWidth
482  fwhm = 2.*filterWidth
483  # The noise should be lower in the smoothed image by sqrt(Nsmooth) ~ fwhm pixels
484  noiseLevel /= fwhm
485  smoothRef = ndimage.filters.gaussian_filter(referenceImage, filterWidth, mode='constant')
486  # Add a three sigma offset to both the reference and model to prevent dividing by zero.
487  # Note that this will also slightly suppress faint variations in color.
488  smoothRef += 3.*noiseLevel
489 
490  lowThreshold = smoothRef/maxDiff
491  highThreshold = smoothRef*maxDiff
492  for model in modelImages:
493  self.applyImageThresholds(model.array,
494  highThreshold=highThreshold,
495  lowThreshold=lowThreshold,
496  regularizationWidth=regularizationWidth)
497  smoothModel = ndimage.filters.gaussian_filter(model.array, filterWidth, mode='constant')
498  smoothModel += 3.*noiseLevel
499  relativeModel = smoothModel/smoothRef
500  # Now sharpen the smoothed relativeModel using an alpha of 3.
501  alpha = 3.
502  relativeModel2 = ndimage.filters.gaussian_filter(relativeModel, filterWidth/alpha)
503  relativeModel += alpha*(relativeModel - relativeModel2)
504  model.array = relativeModel*referenceImage
505 
506  def calculateNoiseCutoff(self, image, statsCtrl, bufferSize,
507  convergenceMaskPlanes="DETECTED", mask=None, bbox=None):
508  """Helper function to calculate the background noise level of an image.
509 
510  Parameters
511  ----------
512  image : `lsst.afw.image.Image`
513  The input image to evaluate the background noise properties.
514  statsCtrl : `lsst.afw.math.StatisticsControl`
515  Statistics control object for coaddition.
516  bufferSize : `int`
517  Number of additional pixels to exclude
518  from the edges of the bounding box.
519  convergenceMaskPlanes : `list` of `str`, or `str`
520  Mask planes to use to calculate convergence.
521  mask : `lsst.afw.image.Mask`, Optional
522  Optional alternate mask
523  bbox : `lsst.afw.geom.Box2I`, optional
524  Sub-region of the masked image to calculate the noise level over.
525 
526  Returns
527  -------
528  noiseCutoff : `float`
529  The threshold value to treat pixels as noise in an image..
530  """
531  if bbox is None:
532  bbox = self.bbox
533  if mask is None:
534  mask = self.mask[bbox]
535  bboxShrink = afwGeom.Box2I(bbox)
536  bboxShrink.grow(-bufferSize)
537  convergeMask = mask.getPlaneBitMask(convergenceMaskPlanes)
538 
539  backgroundPixels = mask[bboxShrink].array & (statsCtrl.getAndMask() | convergeMask) == 0
540  noiseCutoff = np.std(image[bboxShrink].array[backgroundPixels])
541  return noiseCutoff
542 
543  def applyImageThresholds(self, image, highThreshold=None, lowThreshold=None, regularizationWidth=2):
544  """Restrict image values to be between upper and lower limits.
545 
546  This method flags all pixels in an image that are outside of the given
547  threshold values. The threshold values are taken from a reference image,
548  so noisy pixels are likely to get flagged. In order to exclude those
549  noisy pixels, the array of flags is eroded and dilated, which removes
550  isolated pixels outside of the thresholds from the list of pixels to be
551  modified. Pixels that remain flagged after this operation have their
552  values set to the appropriate upper or lower threshold value.
553 
554  Parameters
555  ----------
556  image : `numpy.ndarray`
557  The image to apply the thresholds to.
558  The values will be modified in place.
559  highThreshold : `numpy.ndarray`, optional
560  Array of upper limit values for each pixel of ``image``.
561  lowThreshold : `numpy.ndarray`, optional
562  Array of lower limit values for each pixel of ``image``.
563  regularizationWidth : `int`, optional
564  Minimum radius of a region to include in regularization, in pixels.
565  """
566  # Generate the structure for binary erosion and dilation, which is used to remove noise-like pixels.
567  # Groups of pixels with a radius smaller than ``regularizationWidth``
568  # will be excluded from regularization.
569  filterStructure = ndimage.iterate_structure(ndimage.generate_binary_structure(2, 1),
570  regularizationWidth)
571  if highThreshold is not None:
572  highPixels = image > highThreshold
573  if regularizationWidth > 0:
574  # Erode and dilate ``highPixels`` to exclude noisy pixels.
575  highPixels = ndimage.morphology.binary_opening(highPixels, structure=filterStructure)
576  image[highPixels] = highThreshold[highPixels]
577  if lowThreshold is not None:
578  lowPixels = image < lowThreshold
579  if regularizationWidth > 0:
580  # Erode and dilate ``lowPixels`` to exclude noisy pixels.
581  lowPixels = ndimage.morphology.binary_opening(lowPixels, structure=filterStructure)
582  image[lowPixels] = lowThreshold[lowPixels]
583 
584 
585 def applyDcr(image, dcr, useInverse=False, splitSubfilters=False, splitThreshold=0.,
586  doPrefilter=True, order=3):
587  """Shift an image along the X and Y directions.
588 
589  Parameters
590  ----------
591  image : `numpy.ndarray`
592  The input image to shift.
593  dcr : `tuple`
594  Shift calculated with ``calculateDcr``.
595  Uses numpy axes ordering (Y, X).
596  If ``splitSubfilters`` is set, each element is itself a `tuple`
597  of two `float`, corresponding to the DCR shift at the two wavelengths.
598  Otherwise, each element is a `float` corresponding to the DCR shift at
599  the effective wavelength of the subfilter.
600  useInverse : `bool`, optional
601  Apply the shift in the opposite direction. Default: False
602  splitSubfilters : `bool`, optional
603  Calculate DCR for two evenly-spaced wavelengths in each subfilter,
604  instead of at the midpoint. Default: False
605  splitThreshold : `float`, optional
606  Minimum DCR difference within a subfilter required to use ``splitSubfilters``
607  doPrefilter : `bool`, optional
608  Spline filter the image before shifting, if set. Filtering is required,
609  so only set to False if the image is already filtered.
610  Filtering takes ~20% of the time of shifting, so if `applyDcr` will be
611  called repeatedly on the same image it is more efficient to precalculate
612  the filter.
613  order : `int`, optional
614  The order of the spline interpolation, default is 3.
615 
616  Returns
617  -------
618  shiftedImage : `numpy.ndarray`
619  A copy of the input image with the specified shift applied.
620  """
621  if doPrefilter:
622  prefilteredImage = ndimage.spline_filter(image, order=order)
623  else:
624  prefilteredImage = image
625  if splitSubfilters:
626  shiftAmp = np.max(np.abs([_dcr0 - _dcr1 for _dcr0, _dcr1 in zip(dcr[0], dcr[1])]))
627  if shiftAmp >= splitThreshold:
628  if useInverse:
629  shift = [-1.*s for s in dcr[0]]
630  shift1 = [-1.*s for s in dcr[1]]
631  else:
632  shift = dcr[0]
633  shift1 = dcr[1]
634  shiftedImage = ndimage.shift(prefilteredImage, shift, prefilter=False, order=order)
635  shiftedImage += ndimage.shift(prefilteredImage, shift1, prefilter=False, order=order)
636  shiftedImage /= 2.
637  return shiftedImage
638  else:
639  # If the difference in the DCR shifts is less than the threshold,
640  # then just use the average shift for efficiency.
641  dcr = (np.mean(dcr[0]), np.mean(dcr[1]))
642  if useInverse:
643  shift = [-1.*s for s in dcr]
644  else:
645  shift = dcr
646  shiftedImage = ndimage.shift(prefilteredImage, shift, prefilter=False, order=order)
647  return shiftedImage
648 
649 
650 def calculateDcr(visitInfo, wcs, filterInfo, dcrNumSubfilters, splitSubfilters=False):
651  """Calculate the shift in pixels of an exposure due to DCR.
652 
653  Parameters
654  ----------
655  visitInfo : `lsst.afw.image.VisitInfo`
656  Metadata for the exposure.
657  wcs : `lsst.afw.geom.SkyWcs`
658  Coordinate system definition (wcs) for the exposure.
659  filterInfo : `lsst.afw.image.Filter`
660  The filter definition, set in the current instruments' obs package.
661  dcrNumSubfilters : `int`
662  Number of sub-filters used to model chromatic effects within a band.
663  splitSubfilters : `bool`, optional
664  Calculate DCR for two evenly-spaced wavelengths in each subfilter,
665  instead of at the midpoint. Default: False
666 
667  Returns
668  -------
669  dcrShift : `tuple` of two `float`
670  The 2D shift due to DCR, in pixels.
671  Uses numpy axes ordering (Y, X).
672  """
673  rotation = calculateImageParallacticAngle(visitInfo, wcs)
674  dcrShift = []
675  weight = [0.75, 0.25]
676  lambdaEff = filterInfo.getFilterProperty().getLambdaEff()
677  for wl0, wl1 in wavelengthGenerator(filterInfo, dcrNumSubfilters):
678  # Note that diffRefractAmp can be negative, since it's relative to the midpoint of the full band
679  diffRefractAmp0 = differentialRefraction(wavelength=wl0, wavelengthRef=lambdaEff,
680  elevation=visitInfo.getBoresightAzAlt().getLatitude(),
681  observatory=visitInfo.getObservatory(),
682  weather=visitInfo.getWeather())
683  diffRefractAmp1 = differentialRefraction(wavelength=wl1, wavelengthRef=lambdaEff,
684  elevation=visitInfo.getBoresightAzAlt().getLatitude(),
685  observatory=visitInfo.getObservatory(),
686  weather=visitInfo.getWeather())
687  if splitSubfilters:
688  diffRefractPix0 = diffRefractAmp0.asArcseconds()/wcs.getPixelScale().asArcseconds()
689  diffRefractPix1 = diffRefractAmp1.asArcseconds()/wcs.getPixelScale().asArcseconds()
690  diffRefractArr = [diffRefractPix0*weight[0] + diffRefractPix1*weight[1],
691  diffRefractPix0*weight[1] + diffRefractPix1*weight[0]]
692  shiftX = [diffRefractPix*np.sin(rotation.asRadians()) for diffRefractPix in diffRefractArr]
693  shiftY = [diffRefractPix*np.cos(rotation.asRadians()) for diffRefractPix in diffRefractArr]
694  dcrShift.append(((shiftY[0], shiftX[0]), (shiftY[1], shiftX[1])))
695  else:
696  diffRefractAmp = (diffRefractAmp0 + diffRefractAmp1)/2.
697  diffRefractPix = diffRefractAmp.asArcseconds()/wcs.getPixelScale().asArcseconds()
698  shiftX = diffRefractPix*np.sin(rotation.asRadians())
699  shiftY = diffRefractPix*np.cos(rotation.asRadians())
700  dcrShift.append((shiftY, shiftX))
701  return dcrShift
702 
703 
704 def calculateImageParallacticAngle(visitInfo, wcs):
705  """Calculate the total sky rotation angle of an exposure.
706 
707  Parameters
708  ----------
709  visitInfo : `lsst.afw.image.VisitInfo`
710  Metadata for the exposure.
711  wcs : `lsst.afw.geom.SkyWcs`
712  Coordinate system definition (wcs) for the exposure.
713 
714  Returns
715  -------
716  `lsst.geom.Angle`
717  The rotation of the image axis, East from North.
718  Equal to the parallactic angle plus any additional rotation of the
719  coordinate system.
720  A rotation angle of 0 degrees is defined with
721  North along the +y axis and East along the +x axis.
722  A rotation angle of 90 degrees is defined with
723  North along the +x axis and East along the -y axis.
724  """
725  parAngle = visitInfo.getBoresightParAngle().asRadians()
726  cd = wcs.getCdMatrix()
727  if wcs.isFlipped:
728  cdAngle = (np.arctan2(-cd[0, 1], cd[0, 0]) + np.arctan2(cd[1, 0], cd[1, 1]))/2.
729  rotAngle = (cdAngle + parAngle)*radians
730  else:
731  cdAngle = (np.arctan2(cd[0, 1], -cd[0, 0]) + np.arctan2(cd[1, 0], cd[1, 1]))/2.
732  rotAngle = (cdAngle - parAngle)*radians
733  return rotAngle
734 
735 
736 def wavelengthGenerator(filterInfo, dcrNumSubfilters):
737  """Iterate over the wavelength endpoints of subfilters.
738 
739  Parameters
740  ----------
741  filterInfo : `lsst.afw.image.Filter`
742  The filter definition, set in the current instruments' obs package.
743  dcrNumSubfilters : `int`
744  Number of sub-filters used to model chromatic effects within a band.
745 
746  Yields
747  ------
748  `tuple` of two `float`
749  The next set of wavelength endpoints for a subfilter, in nm.
750  """
751  lambdaMin = filterInfo.getFilterProperty().getLambdaMin()
752  lambdaMax = filterInfo.getFilterProperty().getLambdaMax()
753  wlStep = (lambdaMax - lambdaMin)/dcrNumSubfilters
754  for wl in np.linspace(lambdaMin, lambdaMax, dcrNumSubfilters, endpoint=False):
755  yield (wl, wl + wlStep)
def fromImage(cls, maskedImage, dcrNumSubfilters, filterInfo=None, psf=None)
Definition: dcrModel.py:62
def __setitem__(self, subfilter, maskedImage)
Definition: dcrModel.py:182
def calculateImageParallacticAngle(visitInfo, wcs)
Definition: dcrModel.py:704
def applyDcr(image, dcr, useInverse=False, splitSubfilters=False, splitThreshold=0., doPrefilter=True, order=3)
Definition: dcrModel.py:586
def regularizeModelIter(self, subfilter, newModel, bbox, regularizationFactor, regularizationWidth=2)
Definition: dcrModel.py:419
def __getitem__(self, subfilter)
Definition: dcrModel.py:157
def calculateNoiseCutoff(self, image, statsCtrl, bufferSize, convergenceMaskPlanes="DETECTED", mask=None, bbox=None)
Definition: dcrModel.py:507
def __init__(self, modelImages, filterInfo=None, psf=None, mask=None, variance=None)
Definition: dcrModel.py:53
def wavelengthGenerator(filterInfo, dcrNumSubfilters)
Definition: dcrModel.py:736
def applyImageThresholds(self, image, highThreshold=None, lowThreshold=None, regularizationWidth=2)
Definition: dcrModel.py:543
def buildMatchedExposure(self, exposure=None, visitInfo=None, bbox=None, wcs=None, mask=None)
Definition: dcrModel.py:361
def fromDataRef(cls, dataRef, datasetType="dcrCoadd", numSubfilters=None, kwargs)
Definition: dcrModel.py:106
def getReferenceImage(self, bbox=None)
Definition: dcrModel.py:261
def buildMatchedTemplate(self, exposure=None, order=3, visitInfo=None, bbox=None, wcs=None, mask=None, splitSubfilters=True, splitThreshold=0., amplifyModel=1.)
Definition: dcrModel.py:302
def assign(self, dcrSubModel, bbox=None)
Definition: dcrModel.py:277
def calculateDcr(visitInfo, wcs, filterInfo, dcrNumSubfilters, splitSubfilters=False)
Definition: dcrModel.py:650
def regularizeModelFreq(self, modelImages, bbox, statsCtrl, regularizationFactor, regularizationWidth=2, mask=None, convergenceMaskPlanes="DETECTED")
Definition: dcrModel.py:445
def conditionDcrModel(self, modelImages, bbox, gain=1.)
Definition: dcrModel.py:398