Coverage for python/lsst/ip/diffim/dcrModel.py : 10%

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