Hide keyboard shortcuts

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# 

22 

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 

28 

29__all__ = ["DcrModel", "applyDcr", "calculateDcr", "calculateImageParallacticAngle"] 

30 

31 

32class DcrModel: 

33 """A model of the true sky after correcting chromatic effects. 

34 

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 

41 

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 """ 

51 

52 def __init__(self, modelImages, filterInfo=None, psf=None, mask=None, variance=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 

60 @classmethod 

61 def fromImage(cls, maskedImage, dcrNumSubfilters, filterInfo=None, psf=None): 

62 """Initialize a DcrModel by dividing a coadd between the subfilters. 

63 

64 Parameters 

65 ---------- 

66 maskedImage : `lsst.afw.image.MaskedImage` 

67 Input coadded image to divide equally between the subfilters. 

68 dcrNumSubfilters : `int` 

69 Number of sub-filters used to model chromatic effects within a band. 

70 filterInfo : `lsst.afw.image.Filter`, optional 

71 The filter definition, set in the current instruments' obs package. 

72 Required for any calculation of DCR, including making matched templates. 

73 psf : `lsst.afw.detection.Psf`, optional 

74 Point spread function (PSF) of the model. 

75 Required if the ``DcrModel`` will be persisted. 

76 

77 Returns 

78 ------- 

79 dcrModel : `lsst.pipe.tasks.DcrModel` 

80 Best fit model of the true sky after correcting chromatic effects. 

81 

82 Raises 

83 ------ 

84 ValueError 

85 If there are any unmasked NAN values in ``maskedImage``. 

86 """ 

87 # NANs will potentially contaminate the entire image, 

88 # depending on the shift or convolution type used. 

89 model = maskedImage.image.clone() 

90 mask = maskedImage.mask.clone() 

91 # We divide the variance by N and not N**2 because we will assume each 

92 # subfilter is independent. That means that the significance of 

93 # detected sources will be lower by a factor of sqrt(N) in the 

94 # subfilter images, but we will recover it when we combine the 

95 # subfilter images to construct matched templates. 

96 variance = maskedImage.variance.clone() 

97 variance /= dcrNumSubfilters 

98 model /= dcrNumSubfilters 

99 modelImages = [model, ] 

100 for subfilter in range(1, dcrNumSubfilters): 

101 modelImages.append(model.clone()) 

102 return cls(modelImages, filterInfo, psf, mask, variance) 

103 

104 @classmethod 

105 def fromDataRef(cls, dataRef, datasetType="dcrCoadd", numSubfilters=None, **kwargs): 

106 """Load an existing DcrModel from a repository. 

107 

108 Parameters 

109 ---------- 

110 dataRef : `lsst.daf.persistence.ButlerDataRef` 

111 Data reference defining the patch for coaddition and the 

112 reference Warp 

113 datasetType : `str`, optional 

114 Name of the DcrModel in the registry {"dcrCoadd", "dcrCoadd_sub"} 

115 numSubfilters : `int` 

116 Number of sub-filters used to model chromatic effects within a band. 

117 **kwargs 

118 Additional keyword arguments to pass to look up the model in the data registry. 

119 Common keywords and their types include: ``tract``:`str`, ``patch``:`str`, 

120 ``bbox``:`lsst.afw.geom.Box2I` 

121 

122 Returns 

123 ------- 

124 dcrModel : `lsst.pipe.tasks.DcrModel` 

125 Best fit model of the true sky after correcting chromatic effects. 

126 """ 

127 modelImages = [] 

128 filterInfo = None 

129 psf = None 

130 mask = None 

131 variance = None 

132 for subfilter in range(numSubfilters): 

133 dcrCoadd = dataRef.get(datasetType, subfilter=subfilter, 

134 numSubfilters=numSubfilters, **kwargs) 

135 if filterInfo is None: 

136 filterInfo = dcrCoadd.getFilter() 

137 if psf is None: 

138 psf = dcrCoadd.getPsf() 

139 if mask is None: 

140 mask = dcrCoadd.mask 

141 if variance is None: 

142 variance = dcrCoadd.variance 

143 modelImages.append(dcrCoadd.image) 

144 return cls(modelImages, filterInfo, psf, mask, variance) 

145 

146 def __len__(self): 

147 """Return the number of subfilters. 

148 

149 Returns 

150 ------- 

151 dcrNumSubfilters : `int` 

152 The number of DCR subfilters in the model. 

153 """ 

154 return self.dcrNumSubfilters 

155 

156 def __getitem__(self, subfilter): 

157 """Iterate over the subfilters of the DCR model. 

158 

159 Parameters 

160 ---------- 

161 subfilter : `int` 

162 Index of the current ``subfilter`` within the full band. 

163 Negative indices are allowed, and count in reverse order 

164 from the highest ``subfilter``. 

165 

166 Returns 

167 ------- 

168 modelImage : `lsst.afw.image.Image` 

169 The DCR model for the given ``subfilter``. 

170 

171 Raises 

172 ------ 

173 IndexError 

174 If the requested ``subfilter`` is greater or equal to the number 

175 of subfilters in the model. 

176 """ 

177 if np.abs(subfilter) >= len(self): 

178 raise IndexError("subfilter out of bounds.") 

179 return self.modelImages[subfilter] 

180 

181 def __setitem__(self, subfilter, maskedImage): 

182 """Update the model image for one subfilter. 

183 

184 Parameters 

185 ---------- 

186 subfilter : `int` 

187 Index of the current subfilter within the full band. 

188 maskedImage : `lsst.afw.image.Image` 

189 The DCR model to set for the given ``subfilter``. 

190 

191 Raises 

192 ------ 

193 IndexError 

194 If the requested ``subfilter`` is greater or equal to the number 

195 of subfilters in the model. 

196 ValueError 

197 If the bounding box of the new image does not match. 

198 """ 

199 if np.abs(subfilter) >= len(self): 

200 raise IndexError("subfilter out of bounds.") 

201 if maskedImage.getBBox() != self.bbox: 

202 raise ValueError("The bounding box of a subfilter must not change.") 

203 self.modelImages[subfilter] = maskedImage 

204 

205 @property 

206 def filter(self): 

207 """Return the filter of the model. 

208 

209 Returns 

210 ------- 

211 filter : `lsst.afw.image.Filter` 

212 The filter definition, set in the current instruments' obs package. 

213 """ 

214 return self._filter 

215 

216 @property 

217 def psf(self): 

218 """Return the psf of the model. 

219 

220 Returns 

221 ------- 

222 psf : `lsst.afw.detection.Psf` 

223 Point spread function (PSF) of the model. 

224 """ 

225 return self._psf 

226 

227 @property 

228 def bbox(self): 

229 """Return the common bounding box of each subfilter image. 

230 

231 Returns 

232 ------- 

233 bbox : `lsst.afw.geom.Box2I` 

234 Bounding box of the DCR model. 

235 """ 

236 return self[0].getBBox() 

237 

238 @property 

239 def mask(self): 

240 """Return the common mask of each subfilter image. 

241 

242 Returns 

243 ------- 

244 mask : `lsst.afw.image.Mask` 

245 Mask plane of the DCR model. 

246 """ 

247 return self._mask 

248 

249 @property 

250 def variance(self): 

251 """Return the common variance of each subfilter image. 

252 

253 Returns 

254 ------- 

255 variance : `lsst.afw.image.Image` 

256 Variance plane of the DCR model. 

257 """ 

258 return self._variance 

259 

260 def getReferenceImage(self, bbox=None): 

261 """Calculate a reference image from the average of the subfilter images. 

262 

263 Parameters 

264 ---------- 

265 bbox : `lsst.afw.geom.Box2I`, optional 

266 Sub-region of the coadd. Returns the entire image if `None`. 

267 

268 Returns 

269 ------- 

270 refImage : `numpy.ndarray` 

271 The reference image with no chromatic effects applied. 

272 """ 

273 bbox = bbox or self.bbox 

274 return np.mean([model[bbox].array for model in self], axis=0) 

275 

276 def assign(self, dcrSubModel, bbox=None): 

277 """Update a sub-region of the ``DcrModel`` with new values. 

278 

279 Parameters 

280 ---------- 

281 dcrSubModel : `lsst.pipe.tasks.DcrModel` 

282 New model of the true scene after correcting chromatic effects. 

283 bbox : `lsst.afw.geom.Box2I`, optional 

284 Sub-region of the coadd. 

285 Defaults to the bounding box of ``dcrSubModel``. 

286 

287 Raises 

288 ------ 

289 ValueError 

290 If the new model has a different number of subfilters. 

291 """ 

292 if len(dcrSubModel) != len(self): 

293 raise ValueError("The number of DCR subfilters must be the same " 

294 "between the old and new models.") 

295 bbox = bbox or self.bbox 

296 for model, subModel in zip(self, dcrSubModel): 

297 model.assign(subModel[bbox], bbox) 

298 

299 def buildMatchedTemplate(self, exposure=None, order=3, 

300 visitInfo=None, bbox=None, wcs=None, mask=None, 

301 splitSubfilters=True, splitThreshold=0., amplifyModel=1.): 

302 """Create a DCR-matched template image for an exposure. 

303 

304 Parameters 

305 ---------- 

306 exposure : `lsst.afw.image.Exposure`, optional 

307 The input exposure to build a matched template for. 

308 May be omitted if all of the metadata is supplied separately 

309 order : `int`, optional 

310 Interpolation order of the DCR shift. 

311 visitInfo : `lsst.afw.image.VisitInfo`, optional 

312 Metadata for the exposure. Ignored if ``exposure`` is set. 

313 bbox : `lsst.afw.geom.Box2I`, optional 

314 Sub-region of the coadd. Ignored if ``exposure`` is set. 

315 wcs : `lsst.afw.geom.SkyWcs`, optional 

316 Coordinate system definition (wcs) for the exposure. 

317 Ignored if ``exposure`` is set. 

318 mask : `lsst.afw.image.Mask`, optional 

319 reference mask to use for the template image. 

320 splitSubfilters : `bool`, optional 

321 Calculate DCR for two evenly-spaced wavelengths in each subfilter, 

322 instead of at the midpoint. Default: True 

323 splitThreshold : `float`, optional 

324 Minimum DCR difference within a subfilter required to use ``splitSubfilters`` 

325 amplifyModel : `float`, optional 

326 Multiplication factor to amplify differences between model planes. 

327 Used to speed convergence of iterative forward modeling. 

328 

329 Returns 

330 ------- 

331 templateImage : `lsst.afw.image.ImageF` 

332 The DCR-matched template 

333 

334 Raises 

335 ------ 

336 ValueError 

337 If neither ``exposure`` or all of ``visitInfo``, ``bbox``, and ``wcs`` are set. 

338 """ 

339 if self.filter is None: 

340 raise ValueError("'filterInfo' must be set for the DcrModel in order to calculate DCR.") 

341 if exposure is not None: 

342 visitInfo = exposure.getInfo().getVisitInfo() 

343 bbox = exposure.getBBox() 

344 wcs = exposure.getInfo().getWcs() 

345 elif visitInfo is None or bbox is None or wcs is None: 

346 raise ValueError("Either exposure or visitInfo, bbox, and wcs must be set.") 

347 dcrShift = calculateDcr(visitInfo, wcs, self.filter, len(self), splitSubfilters=splitSubfilters) 

348 templateImage = afwImage.ImageF(bbox) 

349 refModel = self.getReferenceImage(bbox) 

350 for subfilter, dcr in enumerate(dcrShift): 

351 if amplifyModel > 1: 

352 model = (self[subfilter][bbox].array - refModel)*amplifyModel + refModel 

353 else: 

354 model = self[subfilter][bbox].array 

355 templateImage.array += applyDcr(model, dcr, splitSubfilters=splitSubfilters, 

356 splitThreshold=splitThreshold, order=order) 

357 return templateImage 

358 

359 def buildMatchedExposure(self, exposure=None, 

360 visitInfo=None, bbox=None, wcs=None, mask=None): 

361 """Wrapper to create an exposure from a template image. 

362 

363 Parameters 

364 ---------- 

365 exposure : `lsst.afw.image.Exposure`, optional 

366 The input exposure to build a matched template for. 

367 May be omitted if all of the metadata is supplied separately 

368 visitInfo : `lsst.afw.image.VisitInfo`, optional 

369 Metadata for the exposure. Ignored if ``exposure`` is set. 

370 bbox : `lsst.afw.geom.Box2I`, optional 

371 Sub-region of the coadd. Ignored if ``exposure`` is set. 

372 wcs : `lsst.afw.geom.SkyWcs`, optional 

373 Coordinate system definition (wcs) for the exposure. 

374 Ignored if ``exposure`` is set. 

375 mask : `lsst.afw.image.Mask`, optional 

376 reference mask to use for the template image. 

377 

378 Returns 

379 ------- 

380 templateExposure : `lsst.afw.image.exposureF` 

381 The DCR-matched template 

382 """ 

383 if bbox is None: 

384 bbox = exposure.getBBox() 

385 templateImage = self.buildMatchedTemplate(exposure=exposure, visitInfo=visitInfo, 

386 bbox=bbox, wcs=wcs, mask=mask) 

387 maskedImage = afwImage.MaskedImageF(bbox) 

388 maskedImage.image = templateImage[bbox] 

389 maskedImage.mask = self.mask[bbox] 

390 maskedImage.variance = self.variance[bbox] 

391 templateExposure = afwImage.ExposureF(bbox, wcs) 

392 templateExposure.setMaskedImage(maskedImage[bbox]) 

393 templateExposure.setPsf(self.psf) 

394 templateExposure.setFilter(self.filter) 

395 return templateExposure 

396 

397 def conditionDcrModel(self, modelImages, bbox, gain=1.): 

398 """Average two iterations' solutions to reduce oscillations. 

399 

400 Parameters 

401 ---------- 

402 modelImages : `list` of `lsst.afw.image.Image` 

403 The new DCR model images from the current iteration. 

404 The values will be modified in place. 

405 bbox : `lsst.afw.geom.Box2I` 

406 Sub-region of the coadd 

407 gain : `float`, optional 

408 Relative weight to give the new solution when updating the model. 

409 Defaults to 1.0, which gives equal weight to both solutions. 

410 """ 

411 # Calculate weighted averages of the images. 

412 for model, newModel in zip(self, modelImages): 

413 newModel *= gain 

414 newModel += model[bbox] 

415 newModel /= 1. + gain 

416 

417 def regularizeModelIter(self, subfilter, newModel, bbox, regularizationFactor, 

418 regularizationWidth=2): 

419 """Restrict large variations in the model between iterations. 

420 

421 Parameters 

422 ---------- 

423 subfilter : `int` 

424 Index of the current subfilter within the full band. 

425 newModel : `lsst.afw.image.Image` 

426 The new DCR model for one subfilter from the current iteration. 

427 Values in ``newModel`` that are extreme compared with the last 

428 iteration are modified in place. 

429 bbox : `lsst.afw.geom.Box2I` 

430 Sub-region to coadd 

431 regularizationFactor : `float` 

432 Maximum relative change of the model allowed between iterations. 

433 regularizationWidth : int, optional 

434 Minimum radius of a region to include in regularization, in pixels. 

435 """ 

436 refImage = self[subfilter][bbox].array 

437 highThreshold = np.abs(refImage)*regularizationFactor 

438 lowThreshold = refImage/regularizationFactor 

439 newImage = newModel.array 

440 self.applyImageThresholds(newImage, highThreshold=highThreshold, lowThreshold=lowThreshold, 

441 regularizationWidth=regularizationWidth) 

442 

443 def regularizeModelFreq(self, modelImages, bbox, statsCtrl, regularizationFactor, 

444 regularizationWidth=2, mask=None, convergenceMaskPlanes="DETECTED"): 

445 """Restrict large variations in the model between subfilters. 

446 

447 Parameters 

448 ---------- 

449 modelImages : `list` of `lsst.afw.image.Image` 

450 The new DCR model images from the current iteration. 

451 The values will be modified in place. 

452 bbox : `lsst.afw.geom.Box2I` 

453 Sub-region to coadd 

454 statsCtrl : `lsst.afw.math.StatisticsControl` 

455 Statistics control object for coaddition. 

456 regularizationFactor : `float` 

457 Maximum relative change of the model allowed between subfilters. 

458 regularizationWidth : `int`, optional 

459 Minimum radius of a region to include in regularization, in pixels. 

460 mask : `lsst.afw.image.Mask`, optional 

461 Optional alternate mask 

462 convergenceMaskPlanes : `list` of `str`, or `str`, optional 

463 Mask planes to use to calculate convergence. 

464 

465 Notes 

466 ----- 

467 This implementation of frequency regularization restricts each subfilter 

468 image to be a smoothly-varying function times a reference image. 

469 """ 

470 # ``regularizationFactor`` is the maximum change between subfilter images, so the maximum difference 

471 # between one subfilter image and the average will be the square root of that. 

472 maxDiff = np.sqrt(regularizationFactor) 

473 noiseLevel = self.calculateNoiseCutoff(modelImages[0], statsCtrl, bufferSize=5, mask=mask, bbox=bbox) 

474 referenceImage = self.getReferenceImage(bbox) 

475 badPixels = np.isnan(referenceImage) | (referenceImage <= 0.) 

476 if np.sum(~badPixels) == 0: 

477 # Skip regularization if there are no valid pixels 

478 return 

479 referenceImage[badPixels] = 0. 

480 filterWidth = regularizationWidth 

481 fwhm = 2.*filterWidth 

482 # The noise should be lower in the smoothed image by sqrt(Nsmooth) ~ fwhm pixels 

483 noiseLevel /= fwhm 

484 smoothRef = ndimage.filters.gaussian_filter(referenceImage, filterWidth, mode='constant') 

485 # Add a three sigma offset to both the reference and model to prevent dividing by zero. 

486 # Note that this will also slightly suppress faint variations in color. 

487 smoothRef += 3.*noiseLevel 

488 

489 lowThreshold = smoothRef/maxDiff 

490 highThreshold = smoothRef*maxDiff 

491 for model in modelImages: 

492 self.applyImageThresholds(model.array, 

493 highThreshold=highThreshold, 

494 lowThreshold=lowThreshold, 

495 regularizationWidth=regularizationWidth) 

496 smoothModel = ndimage.filters.gaussian_filter(model.array, filterWidth, mode='constant') 

497 smoothModel += 3.*noiseLevel 

498 relativeModel = smoothModel/smoothRef 

499 # Now sharpen the smoothed relativeModel using an alpha of 3. 

500 alpha = 3. 

501 relativeModel2 = ndimage.filters.gaussian_filter(relativeModel, filterWidth/alpha) 

502 relativeModel += alpha*(relativeModel - relativeModel2) 

503 model.array = relativeModel*referenceImage 

504 

505 def calculateNoiseCutoff(self, image, statsCtrl, bufferSize, 

506 convergenceMaskPlanes="DETECTED", mask=None, bbox=None): 

507 """Helper function to calculate the background noise level of an image. 

508 

509 Parameters 

510 ---------- 

511 image : `lsst.afw.image.Image` 

512 The input image to evaluate the background noise properties. 

513 statsCtrl : `lsst.afw.math.StatisticsControl` 

514 Statistics control object for coaddition. 

515 bufferSize : `int` 

516 Number of additional pixels to exclude 

517 from the edges of the bounding box. 

518 convergenceMaskPlanes : `list` of `str`, or `str` 

519 Mask planes to use to calculate convergence. 

520 mask : `lsst.afw.image.Mask`, Optional 

521 Optional alternate mask 

522 bbox : `lsst.afw.geom.Box2I`, optional 

523 Sub-region of the masked image to calculate the noise level over. 

524 

525 Returns 

526 ------- 

527 noiseCutoff : `float` 

528 The threshold value to treat pixels as noise in an image.. 

529 """ 

530 if bbox is None: 

531 bbox = self.bbox 

532 if mask is None: 

533 mask = self.mask[bbox] 

534 bboxShrink = geom.Box2I(bbox) 

535 bboxShrink.grow(-bufferSize) 

536 convergeMask = mask.getPlaneBitMask(convergenceMaskPlanes) 

537 

538 backgroundPixels = mask[bboxShrink].array & (statsCtrl.getAndMask() | convergeMask) == 0 

539 noiseCutoff = np.std(image[bboxShrink].array[backgroundPixels]) 

540 return noiseCutoff 

541 

542 def applyImageThresholds(self, image, highThreshold=None, lowThreshold=None, regularizationWidth=2): 

543 """Restrict image values to be between upper and lower limits. 

544 

545 This method flags all pixels in an image that are outside of the given 

546 threshold values. The threshold values are taken from a reference image, 

547 so noisy pixels are likely to get flagged. In order to exclude those 

548 noisy pixels, the array of flags is eroded and dilated, which removes 

549 isolated pixels outside of the thresholds from the list of pixels to be 

550 modified. Pixels that remain flagged after this operation have their 

551 values set to the appropriate upper or lower threshold value. 

552 

553 Parameters 

554 ---------- 

555 image : `numpy.ndarray` 

556 The image to apply the thresholds to. 

557 The values will be modified in place. 

558 highThreshold : `numpy.ndarray`, optional 

559 Array of upper limit values for each pixel of ``image``. 

560 lowThreshold : `numpy.ndarray`, optional 

561 Array of lower limit values for each pixel of ``image``. 

562 regularizationWidth : `int`, optional 

563 Minimum radius of a region to include in regularization, in pixels. 

564 """ 

565 # Generate the structure for binary erosion and dilation, which is used to remove noise-like pixels. 

566 # Groups of pixels with a radius smaller than ``regularizationWidth`` 

567 # will be excluded from regularization. 

568 filterStructure = ndimage.iterate_structure(ndimage.generate_binary_structure(2, 1), 

569 regularizationWidth) 

570 if highThreshold is not None: 

571 highPixels = image > highThreshold 

572 if regularizationWidth > 0: 

573 # Erode and dilate ``highPixels`` to exclude noisy pixels. 

574 highPixels = ndimage.morphology.binary_opening(highPixels, structure=filterStructure) 

575 image[highPixels] = highThreshold[highPixels] 

576 if lowThreshold is not None: 

577 lowPixels = image < lowThreshold 

578 if regularizationWidth > 0: 

579 # Erode and dilate ``lowPixels`` to exclude noisy pixels. 

580 lowPixels = ndimage.morphology.binary_opening(lowPixels, structure=filterStructure) 

581 image[lowPixels] = lowThreshold[lowPixels] 

582 

583 

584def applyDcr(image, dcr, useInverse=False, splitSubfilters=False, splitThreshold=0., 

585 doPrefilter=True, order=3): 

586 """Shift an image along the X and Y directions. 

587 

588 Parameters 

589 ---------- 

590 image : `numpy.ndarray` 

591 The input image to shift. 

592 dcr : `tuple` 

593 Shift calculated with ``calculateDcr``. 

594 Uses numpy axes ordering (Y, X). 

595 If ``splitSubfilters`` is set, each element is itself a `tuple` 

596 of two `float`, corresponding to the DCR shift at the two wavelengths. 

597 Otherwise, each element is a `float` corresponding to the DCR shift at 

598 the effective wavelength of the subfilter. 

599 useInverse : `bool`, optional 

600 Apply the shift in the opposite direction. Default: False 

601 splitSubfilters : `bool`, optional 

602 Calculate DCR for two evenly-spaced wavelengths in each subfilter, 

603 instead of at the midpoint. Default: False 

604 splitThreshold : `float`, optional 

605 Minimum DCR difference within a subfilter required to use ``splitSubfilters`` 

606 doPrefilter : `bool`, optional 

607 Spline filter the image before shifting, if set. Filtering is required, 

608 so only set to False if the image is already filtered. 

609 Filtering takes ~20% of the time of shifting, so if `applyDcr` will be 

610 called repeatedly on the same image it is more efficient to precalculate 

611 the filter. 

612 order : `int`, optional 

613 The order of the spline interpolation, default is 3. 

614 

615 Returns 

616 ------- 

617 shiftedImage : `numpy.ndarray` 

618 A copy of the input image with the specified shift applied. 

619 """ 

620 if doPrefilter: 

621 prefilteredImage = ndimage.spline_filter(image, order=order) 

622 else: 

623 prefilteredImage = image 

624 if splitSubfilters: 

625 shiftAmp = np.max(np.abs([_dcr0 - _dcr1 for _dcr0, _dcr1 in zip(dcr[0], dcr[1])])) 

626 if shiftAmp >= splitThreshold: 

627 if useInverse: 

628 shift = [-1.*s for s in dcr[0]] 

629 shift1 = [-1.*s for s in dcr[1]] 

630 else: 

631 shift = dcr[0] 

632 shift1 = dcr[1] 

633 shiftedImage = ndimage.shift(prefilteredImage, shift, prefilter=False, order=order) 

634 shiftedImage += ndimage.shift(prefilteredImage, shift1, prefilter=False, order=order) 

635 shiftedImage /= 2. 

636 return shiftedImage 

637 else: 

638 # If the difference in the DCR shifts is less than the threshold, 

639 # then just use the average shift for efficiency. 

640 dcr = (np.mean(dcr[0]), np.mean(dcr[1])) 

641 if useInverse: 

642 shift = [-1.*s for s in dcr] 

643 else: 

644 shift = dcr 

645 shiftedImage = ndimage.shift(prefilteredImage, shift, prefilter=False, order=order) 

646 return shiftedImage 

647 

648 

649def calculateDcr(visitInfo, wcs, filterInfo, dcrNumSubfilters, splitSubfilters=False): 

650 """Calculate the shift in pixels of an exposure due to DCR. 

651 

652 Parameters 

653 ---------- 

654 visitInfo : `lsst.afw.image.VisitInfo` 

655 Metadata for the exposure. 

656 wcs : `lsst.afw.geom.SkyWcs` 

657 Coordinate system definition (wcs) for the exposure. 

658 filterInfo : `lsst.afw.image.Filter` 

659 The filter definition, set in the current instruments' obs package. 

660 dcrNumSubfilters : `int` 

661 Number of sub-filters used to model chromatic effects within a band. 

662 splitSubfilters : `bool`, optional 

663 Calculate DCR for two evenly-spaced wavelengths in each subfilter, 

664 instead of at the midpoint. Default: False 

665 

666 Returns 

667 ------- 

668 dcrShift : `tuple` of two `float` 

669 The 2D shift due to DCR, in pixels. 

670 Uses numpy axes ordering (Y, X). 

671 """ 

672 rotation = calculateImageParallacticAngle(visitInfo, wcs) 

673 dcrShift = [] 

674 weight = [0.75, 0.25] 

675 lambdaEff = filterInfo.getFilterProperty().getLambdaEff() 

676 for wl0, wl1 in wavelengthGenerator(filterInfo, dcrNumSubfilters): 

677 # Note that diffRefractAmp can be negative, since it's relative to the midpoint of the full band 

678 diffRefractAmp0 = differentialRefraction(wavelength=wl0, wavelengthRef=lambdaEff, 

679 elevation=visitInfo.getBoresightAzAlt().getLatitude(), 

680 observatory=visitInfo.getObservatory(), 

681 weather=visitInfo.getWeather()) 

682 diffRefractAmp1 = differentialRefraction(wavelength=wl1, wavelengthRef=lambdaEff, 

683 elevation=visitInfo.getBoresightAzAlt().getLatitude(), 

684 observatory=visitInfo.getObservatory(), 

685 weather=visitInfo.getWeather()) 

686 if splitSubfilters: 

687 diffRefractPix0 = diffRefractAmp0.asArcseconds()/wcs.getPixelScale().asArcseconds() 

688 diffRefractPix1 = diffRefractAmp1.asArcseconds()/wcs.getPixelScale().asArcseconds() 

689 diffRefractArr = [diffRefractPix0*weight[0] + diffRefractPix1*weight[1], 

690 diffRefractPix0*weight[1] + diffRefractPix1*weight[0]] 

691 shiftX = [diffRefractPix*np.sin(rotation.asRadians()) for diffRefractPix in diffRefractArr] 

692 shiftY = [diffRefractPix*np.cos(rotation.asRadians()) for diffRefractPix in diffRefractArr] 

693 dcrShift.append(((shiftY[0], shiftX[0]), (shiftY[1], shiftX[1]))) 

694 else: 

695 diffRefractAmp = (diffRefractAmp0 + diffRefractAmp1)/2. 

696 diffRefractPix = diffRefractAmp.asArcseconds()/wcs.getPixelScale().asArcseconds() 

697 shiftX = diffRefractPix*np.sin(rotation.asRadians()) 

698 shiftY = diffRefractPix*np.cos(rotation.asRadians()) 

699 dcrShift.append((shiftY, shiftX)) 

700 return dcrShift 

701 

702 

703def calculateImageParallacticAngle(visitInfo, wcs): 

704 """Calculate the total sky rotation angle of an exposure. 

705 

706 Parameters 

707 ---------- 

708 visitInfo : `lsst.afw.image.VisitInfo` 

709 Metadata for the exposure. 

710 wcs : `lsst.afw.geom.SkyWcs` 

711 Coordinate system definition (wcs) for the exposure. 

712 

713 Returns 

714 ------- 

715 `lsst.geom.Angle` 

716 The rotation of the image axis, East from North. 

717 Equal to the parallactic angle plus any additional rotation of the 

718 coordinate system. 

719 A rotation angle of 0 degrees is defined with 

720 North along the +y axis and East along the +x axis. 

721 A rotation angle of 90 degrees is defined with 

722 North along the +x axis and East along the -y axis. 

723 """ 

724 parAngle = visitInfo.getBoresightParAngle().asRadians() 

725 cd = wcs.getCdMatrix() 

726 if wcs.isFlipped: 

727 cdAngle = (np.arctan2(-cd[0, 1], cd[0, 0]) + np.arctan2(cd[1, 0], cd[1, 1]))/2. 

728 rotAngle = (cdAngle + parAngle)*geom.radians 

729 else: 

730 cdAngle = (np.arctan2(cd[0, 1], -cd[0, 0]) + np.arctan2(cd[1, 0], cd[1, 1]))/2. 

731 rotAngle = (cdAngle - parAngle)*geom.radians 

732 return rotAngle 

733 

734 

735def wavelengthGenerator(filterInfo, dcrNumSubfilters): 

736 """Iterate over the wavelength endpoints of subfilters. 

737 

738 Parameters 

739 ---------- 

740 filterInfo : `lsst.afw.image.Filter` 

741 The filter definition, set in the current instruments' obs package. 

742 dcrNumSubfilters : `int` 

743 Number of sub-filters used to model chromatic effects within a band. 

744 

745 Yields 

746 ------ 

747 `tuple` of two `float` 

748 The next set of wavelength endpoints for a subfilter, in nm. 

749 """ 

750 lambdaMin = filterInfo.getFilterProperty().getLambdaMin() 

751 lambdaMax = filterInfo.getFilterProperty().getLambdaMax() 

752 wlStep = (lambdaMax - lambdaMin)/dcrNumSubfilters 

753 for wl in np.linspace(lambdaMin, lambdaMax, dcrNumSubfilters, endpoint=False): 

754 yield (wl, wl + wlStep)