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 cp_pipe. 

2# 

3# Developed for the LSST Data Management System. 

4# This product includes software developed by the LSST Project 

5# (https://www.lsst.org). 

6# See the COPYRIGHT file at the top-level directory of this distribution 

7# for details of code ownership. 

8# 

9# This program is free software: you can redistribute it and/or modify 

10# it under the terms of the GNU General Public License as published by 

11# the Free Software Foundation, either version 3 of the License, or 

12# (at your option) any later version. 

13# 

14# This program is distributed in the hope that it will be useful, 

15# but WITHOUT ANY WARRANTY; without even the implied warranty of 

16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 

17# GNU General Public License for more details. 

18# 

19# You should have received a copy of the GNU General Public License 

20# along with this program. If not, see <https://www.gnu.org/licenses/>. 

21# 

22"""Calculation of brighter-fatter effect correlations and kernels.""" 

23 

24__all__ = ['BrighterFatterKernelSolveTask', 

25 'BrighterFatterKernelSolveConfig'] 

26 

27import numpy as np 

28 

29import lsst.afw.math as afwMath 

30import lsst.pex.config as pexConfig 

31import lsst.pipe.base as pipeBase 

32import lsst.pipe.base.connectionTypes as cT 

33 

34from lsst.ip.isr import (BrighterFatterKernel) 

35from .utils import (funcPolynomial, irlsFit) 

36from ._lookupStaticCalibration import lookupStaticCalibration 

37 

38 

39class BrighterFatterKernelSolveConnections(pipeBase.PipelineTaskConnections, 

40 dimensions=("instrument", "exposure", "detector")): 

41 dummy = cT.Input( 

42 name="raw", 

43 doc="Dummy exposure.", 

44 storageClass='Exposure', 

45 dimensions=("instrument", "exposure", "detector"), 

46 multiple=True, 

47 deferLoad=True, 

48 ) 

49 camera = cT.PrerequisiteInput( 

50 name="camera", 

51 doc="Camera associated with this data.", 

52 storageClass="Camera", 

53 dimensions=("instrument", ), 

54 isCalibration=True, 

55 lookupFunction=lookupStaticCalibration, 

56 ) 

57 inputPtc = cT.PrerequisiteInput( 

58 name="ptc", 

59 doc="Photon transfer curve dataset.", 

60 storageClass="PhotonTransferCurveDataset", 

61 dimensions=("instrument", "detector"), 

62 isCalibration=True, 

63 ) 

64 

65 outputBFK = cT.Output( 

66 name="brighterFatterKernel", 

67 doc="Output measured brighter-fatter kernel.", 

68 storageClass="BrighterFatterKernel", 

69 dimensions=("instrument", "detector"), 

70 isCalibration=True, 

71 ) 

72 

73 

74class BrighterFatterKernelSolveConfig(pipeBase.PipelineTaskConfig, 

75 pipelineConnections=BrighterFatterKernelSolveConnections): 

76 level = pexConfig.ChoiceField( 

77 doc="The level at which to calculate the brighter-fatter kernels", 

78 dtype=str, 

79 default="AMP", 

80 allowed={ 

81 "AMP": "Every amplifier treated separately", 

82 "DETECTOR": "One kernel per detector", 

83 } 

84 ) 

85 ignoreAmpsForAveraging = pexConfig.ListField( 

86 dtype=str, 

87 doc="List of amp names to ignore when averaging the amplifier kernels into the detector" 

88 " kernel. Only relevant for level = DETECTOR", 

89 default=[] 

90 ) 

91 xcorrCheckRejectLevel = pexConfig.Field( 

92 dtype=float, 

93 doc="Rejection level for the sum of the input cross-correlations. Arrays which " 

94 "sum to greater than this are discarded before the clipped mean is calculated.", 

95 default=2.0 

96 ) 

97 nSigmaClip = pexConfig.Field( 

98 dtype=float, 

99 doc="Number of sigma to clip when calculating means for the cross-correlation", 

100 default=5 

101 ) 

102 forceZeroSum = pexConfig.Field( 

103 dtype=bool, 

104 doc="Force the correlation matrix to have zero sum by adjusting the (0,0) value?", 

105 default=False, 

106 ) 

107 useAmatrix = pexConfig.Field( 

108 dtype=bool, 

109 doc="Use the PTC 'a' matrix (Astier et al. 2019 equation 20) " 

110 "instead of the average of measured covariances?", 

111 default=False, 

112 ) 

113 

114 maxIterSuccessiveOverRelaxation = pexConfig.Field( 

115 dtype=int, 

116 doc="The maximum number of iterations allowed for the successive over-relaxation method", 

117 default=10000 

118 ) 

119 eLevelSuccessiveOverRelaxation = pexConfig.Field( 

120 dtype=float, 

121 doc="The target residual error for the successive over-relaxation method", 

122 default=5.0e-14 

123 ) 

124 

125 correlationQuadraticFit = pexConfig.Field( 

126 dtype=bool, 

127 doc="Use a quadratic fit to find the correlations instead of simple averaging?", 

128 default=False, 

129 ) 

130 correlationModelRadius = pexConfig.Field( 

131 dtype=int, 

132 doc="Build a model of the correlation coefficients for radii larger than this value in pixels?", 

133 default=100, 

134 ) 

135 correlationModelSlope = pexConfig.Field( 

136 dtype=float, 

137 doc="Slope of the correlation model for radii larger than correlationModelRadius", 

138 default=-1.35, 

139 ) 

140 

141 

142class BrighterFatterKernelSolveTask(pipeBase.PipelineTask, pipeBase.CmdLineTask): 

143 """Measure appropriate Brighter-Fatter Kernel from the PTC dataset. 

144 

145 """ 

146 ConfigClass = BrighterFatterKernelSolveConfig 

147 _DefaultName = 'cpBfkMeasure' 

148 

149 def runQuantum(self, butlerQC, inputRefs, outputRefs): 

150 """Ensure that the input and output dimensions are passed along. 

151 

152 Parameters 

153 ---------- 

154 butlerQC : `lsst.daf.butler.butlerQuantumContext.ButlerQuantumContext` 

155 Butler to operate on. 

156 inputRefs : `lsst.pipe.base.connections.InputQuantizedConnection` 

157 Input data refs to load. 

158 ouptutRefs : `lsst.pipe.base.connections.OutputQuantizedConnection` 

159 Output data refs to persist. 

160 """ 

161 inputs = butlerQC.get(inputRefs) 

162 

163 # Use the dimensions to set calib/provenance information. 

164 inputs['inputDims'] = inputRefs.inputPtc.dataId.byName() 

165 

166 outputs = self.run(**inputs) 

167 butlerQC.put(outputs, outputRefs) 

168 

169 def run(self, inputPtc, dummy, camera, inputDims): 

170 """Combine covariance information from PTC into brighter-fatter kernels. 

171 

172 Parameters 

173 ---------- 

174 inputPtc : `lsst.ip.isr.PhotonTransferCurveDataset` 

175 PTC data containing per-amplifier covariance measurements. 

176 dummy : `lsst.afw.image.Exposure 

177 The exposure used to select the appropriate PTC dataset. 

178 camera : `lsst.afw.cameraGeom.Camera` 

179 Camera to use for camera geometry information. 

180 inputDims : `lsst.daf.butler.DataCoordinate` or `dict` 

181 DataIds to use to populate the output calibration. 

182 

183 Returns 

184 ------- 

185 results : `lsst.pipe.base.Struct` 

186 The resulst struct containing: 

187 

188 ``outputBfk`` : `lsst.ip.isr.BrighterFatterKernel` 

189 Resulting Brighter-Fatter Kernel. 

190 """ 

191 if len(dummy) == 0: 

192 self.log.warn("No dummy exposure found.") 

193 

194 detector = camera[inputDims['detector']] 

195 detName = detector.getName() 

196 

197 if self.config.level == 'DETECTOR': 

198 detectorCorrList = list() 

199 

200 bfk = BrighterFatterKernel(camera=camera, detectorId=detector.getId(), level=self.config.level) 

201 bfk.means = inputPtc.finalMeans # ADU 

202 bfk.variances = inputPtc.finalVars # ADU^2 

203 # Use the PTC covariances as the cross-correlations. These 

204 # are scaled before the kernel is generated, which performs 

205 # the conversion. 

206 bfk.rawXcorrs = inputPtc.covariances # ADU^2 

207 bfk.badAmps = inputPtc.badAmps 

208 bfk.shape = (inputPtc.covMatrixSide*2 + 1, inputPtc.covMatrixSide*2 + 1) 

209 bfk.gain = inputPtc.gain 

210 bfk.noise = inputPtc.noise 

211 bfk.meanXcorrs = dict() 

212 bfk.valid = dict() 

213 

214 for amp in detector: 

215 ampName = amp.getName() 

216 mask = np.array(inputPtc.expIdMask[ampName], dtype=bool) 

217 

218 gain = bfk.gain[ampName] 

219 fluxes = np.array(bfk.means[ampName])[mask] 

220 variances = np.array(bfk.variances[ampName])[mask] 

221 xCorrList = [np.array(xcorr) for xcorr in bfk.rawXcorrs[ampName]] 

222 xCorrList = np.array(xCorrList)[mask] 

223 

224 if gain <= 0: 

225 # We've received very bad data. 

226 self.log.warn("Impossible gain recieved from PTC for %s: %f. Skipping amplifier.", 

227 ampName, gain) 

228 bfk.meanXcorrs[ampName] = np.zeros(bfk.shape) 

229 bfk.ampKernels[ampName] = np.zeros(bfk.shape) 

230 bfk.valid[ampName] = False 

231 continue 

232 

233 fluxes = np.array([flux*gain for flux in fluxes]) # Now in e^- 

234 variances = np.array([variance*gain*gain for variance in variances]) # Now in e^2- 

235 

236 # This should duplicate Coulton et al. 2017 Equation 22-29 (arxiv:1711.06273) 

237 scaledCorrList = list() 

238 for xcorrNum, (xcorr, flux, var) in enumerate(zip(xCorrList, fluxes, variances), 1): 

239 q = np.array(xcorr) * gain * gain # xcorr now in e^- 

240 q *= 2.0 # Remove factor of 1/2 applied in PTC. 

241 self.log.info("Amp: %s %d/%d Flux: %f Var: %f Q(0,0): %g Q(1,0): %g Q(0,1): %g", 

242 ampName, xcorrNum, len(xCorrList), flux, var, q[0][0], q[1][0], q[0][1]) 

243 

244 # Normalize by the flux, which removes the (0,0) 

245 # component attributable to Poisson noise. This 

246 # contains the two "t I delta(x - x')" terms in Coulton et al. 2017 equation 29 

247 q[0][0] -= 2.0*(flux) 

248 

249 if q[0][0] > 0.0: 

250 self.log.warn("Amp: %s %d skipped due to value of (variance-mean)=%f", 

251 ampName, xcorrNum, q[0][0]) 

252 continue 

253 

254 # This removes the "t (I_a^2 + I_b^2)" factor in Coulton et al. 2017 equation 29. 

255 q /= -2.0*(flux**2) 

256 scaled = self._tileArray(q) 

257 

258 xcorrCheck = np.abs(np.sum(scaled))/np.sum(np.abs(scaled)) 

259 if (xcorrCheck > self.config.xcorrCheckRejectLevel) or not (np.isfinite(xcorrCheck)): 

260 self.log.warn("Amp: %s %d skipped due to value of triangle-inequality sum %f", 

261 ampName, xcorrNum, xcorrCheck) 

262 continue 

263 

264 scaledCorrList.append(scaled) 

265 self.log.info("Amp: %s %d/%d Final: %g XcorrCheck: %f", 

266 ampName, xcorrNum, len(xCorrList), q[0][0], xcorrCheck) 

267 

268 if len(scaledCorrList) == 0: 

269 self.log.warn("Amp: %s All inputs rejected for amp!", ampName) 

270 bfk.meanXcorrs[ampName] = np.zeros(bfk.shape) 

271 bfk.ampKernels[ampName] = np.zeros(bfk.shape) 

272 bfk.valid[ampName] = False 

273 continue 

274 

275 if self.config.useAmatrix: 

276 # Use the aMatrix, ignoring the meanXcorr generated above. 

277 preKernel = np.pad(self._tileArray(np.array(inputPtc.aMatrix[ampName])), ((1, 1))) 

278 elif self.config.correlationQuadraticFit: 

279 # Use a quadratic fit to the correlations as a function of flux. 

280 preKernel = self.quadraticCorrelations(scaledCorrList, fluxes, f"Amp: {ampName}") 

281 else: 

282 # Use a simple average of the measured correlations. 

283 preKernel = self.averageCorrelations(scaledCorrList, f"Amp: {ampName}") 

284 

285 center = int((bfk.shape[0] - 1) / 2) 

286 

287 if self.config.forceZeroSum: 

288 totalSum = np.sum(preKernel) 

289 

290 if self.config.correlationModelRadius < (preKernel.shape[0] - 1) / 2: 

291 # Assume a correlation model of Corr(r) = -preFactor * r^(2 * slope) 

292 preFactor = np.sqrt(preKernel[center, center + 1] * preKernel[center + 1, center]) 

293 slopeFactor = 2.0 * np.abs(self.config.correlationModelSlope) 

294 totalSum += 2.0*np.pi*(preFactor / (slopeFactor*(center + 0.5))**slopeFactor) 

295 

296 preKernel[center, center] -= totalSum 

297 self.log.info("%s Zero-Sum Scale: %g", ampName, totalSum) 

298 

299 finalSum = np.sum(preKernel) 

300 bfk.meanXcorrs[ampName] = preKernel 

301 

302 postKernel = self.successiveOverRelax(preKernel) 

303 bfk.ampKernels[ampName] = postKernel 

304 if self.config.level == 'DETECTOR': 

305 detectorCorrList.extend(scaledCorrList) 

306 bfk.valid[ampName] = True 

307 self.log.info("Amp: %s Sum: %g Center Info Pre: %g Post: %g", 

308 ampName, finalSum, preKernel[center, center], postKernel[center, center]) 

309 

310 # Assemble a detector kernel? 

311 if self.config.level == 'DETECTOR': 

312 preKernel = self.averageCorrelations(detectorCorrList, f"Det: {detName}") 

313 finalSum = np.sum(preKernel) 

314 center = int((bfk.shape[0] - 1) / 2) 

315 

316 postKernel = self.successiveOverRelax(preKernel) 

317 bfk.detKernels[detName] = postKernel 

318 self.log.info("Det: %s Sum: %g Center Info Pre: %g Post: %g", 

319 detName, finalSum, preKernel[center, center], postKernel[center, center]) 

320 

321 return pipeBase.Struct( 

322 outputBFK=bfk, 

323 ) 

324 

325 def averageCorrelations(self, xCorrList, name): 

326 """Average input correlations. 

327 

328 Parameters 

329 ---------- 

330 xCorrList : `list` [`numpy.array`] 

331 List of cross-correlations. 

332 name : `str` 

333 Name for log messages. 

334 

335 Returns 

336 ------- 

337 meanXcorr : `numpy.array` 

338 The averaged cross-correlation. 

339 """ 

340 meanXcorr = np.zeros_like(xCorrList[0]) 

341 xCorrList = np.transpose(xCorrList) 

342 sctrl = afwMath.StatisticsControl() 

343 sctrl.setNumSigmaClip(self.config.nSigmaClip) 

344 for i in range(np.shape(meanXcorr)[0]): 

345 for j in range(np.shape(meanXcorr)[1]): 

346 meanXcorr[i, j] = afwMath.makeStatistics(xCorrList[i, j], 

347 afwMath.MEANCLIP, sctrl).getValue() 

348 

349 # To match previous definitions, pad by one element. 

350 meanXcorr = np.pad(meanXcorr, ((1, 1))) 

351 

352 return meanXcorr 

353 

354 def quadraticCorrelations(self, xCorrList, fluxList, name): 

355 """Measure a quadratic correlation model. 

356 

357 Parameters 

358 ---------- 

359 xCorrList : `list` [`numpy.array`] 

360 List of cross-correlations. 

361 fluxList : `numpy.array` 

362 Associated list of fluxes. 

363 name : `str` 

364 Name for log messages. 

365 

366 Returns 

367 ------- 

368 meanXcorr : `numpy.array` 

369 The averaged cross-correlation. 

370 """ 

371 meanXcorr = np.zeros_like(xCorrList[0]) 

372 fluxList = np.square(fluxList) 

373 xCorrList = np.array(xCorrList) 

374 

375 for i in range(np.shape(meanXcorr)[0]): 

376 for j in range(np.shape(meanXcorr)[1]): 

377 # Fit corrlation_i(x, y) = a0 + a1 * (flux_i)^2 The 

378 # i,j indices are inverted to apply the transposition, 

379 # as is done in the averaging case. 

380 linearFit, linearFitErr, chiSq, weights = irlsFit([0.0, 1e-4], fluxList, 

381 xCorrList[:, j, i], funcPolynomial) 

382 meanXcorr[i, j] = linearFit[1] # Discard the intercept. 

383 self.log.debug("Quad fit meanXcorr[%d,%d] = %g", i, j, linearFit[1]) 

384 

385 # To match previous definitions, pad by one element. 

386 meanXcorr = np.pad(meanXcorr, ((1, 1))) 

387 

388 return meanXcorr 

389 

390 @staticmethod 

391 def _tileArray(in_array): 

392 """Given an input quarter-image, tile/mirror it and return full image. 

393 

394 Given a square input of side-length n, of the form 

395 

396 input = array([[1, 2, 3], 

397 [4, 5, 6], 

398 [7, 8, 9]]) 

399 

400 return an array of size 2n-1 as 

401 

402 output = array([[ 9, 8, 7, 8, 9], 

403 [ 6, 5, 4, 5, 6], 

404 [ 3, 2, 1, 2, 3], 

405 [ 6, 5, 4, 5, 6], 

406 [ 9, 8, 7, 8, 9]]) 

407 

408 Parameters: 

409 ----------- 

410 input : `np.array` 

411 The square input quarter-array 

412 

413 Returns: 

414 -------- 

415 output : `np.array` 

416 The full, tiled array 

417 """ 

418 assert(in_array.shape[0] == in_array.shape[1]) 

419 length = in_array.shape[0] - 1 

420 output = np.zeros((2*length + 1, 2*length + 1)) 

421 

422 for i in range(length + 1): 

423 for j in range(length + 1): 

424 output[i + length, j + length] = in_array[i, j] 

425 output[-i + length, j + length] = in_array[i, j] 

426 output[i + length, -j + length] = in_array[i, j] 

427 output[-i + length, -j + length] = in_array[i, j] 

428 return output 

429 

430 def successiveOverRelax(self, source, maxIter=None, eLevel=None): 

431 """An implementation of the successive over relaxation (SOR) method. 

432 

433 A numerical method for solving a system of linear equations 

434 with faster convergence than the Gauss-Seidel method. 

435 

436 Parameters: 

437 ----------- 

438 source : `numpy.ndarray` 

439 The input array. 

440 maxIter : `int`, optional 

441 Maximum number of iterations to attempt before aborting. 

442 eLevel : `float`, optional 

443 The target error level at which we deem convergence to have 

444 occurred. 

445 

446 Returns: 

447 -------- 

448 output : `numpy.ndarray` 

449 The solution. 

450 """ 

451 if not maxIter: 

452 maxIter = self.config.maxIterSuccessiveOverRelaxation 

453 if not eLevel: 

454 eLevel = self.config.eLevelSuccessiveOverRelaxation 

455 

456 assert source.shape[0] == source.shape[1], "Input array must be square" 

457 # initialize, and set boundary conditions 

458 func = np.zeros([source.shape[0] + 2, source.shape[1] + 2]) 

459 resid = np.zeros([source.shape[0] + 2, source.shape[1] + 2]) 

460 rhoSpe = np.cos(np.pi/source.shape[0]) # Here a square grid is assumed 

461 

462 # Calculate the initial error 

463 for i in range(1, func.shape[0] - 1): 

464 for j in range(1, func.shape[1] - 1): 

465 resid[i, j] = (func[i, j - 1] + func[i, j + 1] + func[i - 1, j] 

466 + func[i + 1, j] - 4*func[i, j] - source[i - 1, j - 1]) 

467 inError = np.sum(np.abs(resid)) 

468 

469 # Iterate until convergence 

470 # We perform two sweeps per cycle, 

471 # updating 'odd' and 'even' points separately 

472 nIter = 0 

473 omega = 1.0 

474 dx = 1.0 

475 while nIter < maxIter*2: 

476 outError = 0 

477 if nIter%2 == 0: 

478 for i in range(1, func.shape[0] - 1, 2): 

479 for j in range(1, func.shape[1] - 1, 2): 

480 resid[i, j] = float(func[i, j-1] + func[i, j + 1] + func[i - 1, j] 

481 + func[i + 1, j] - 4.0*func[i, j] - dx*dx*source[i - 1, j - 1]) 

482 func[i, j] += omega*resid[i, j]*.25 

483 for i in range(2, func.shape[0] - 1, 2): 

484 for j in range(2, func.shape[1] - 1, 2): 

485 resid[i, j] = float(func[i, j - 1] + func[i, j + 1] + func[i - 1, j] 

486 + func[i + 1, j] - 4.0*func[i, j] - dx*dx*source[i - 1, j - 1]) 

487 func[i, j] += omega*resid[i, j]*.25 

488 else: 

489 for i in range(1, func.shape[0] - 1, 2): 

490 for j in range(2, func.shape[1] - 1, 2): 

491 resid[i, j] = float(func[i, j - 1] + func[i, j + 1] + func[i - 1, j] 

492 + func[i + 1, j] - 4.0*func[i, j] - dx*dx*source[i - 1, j - 1]) 

493 func[i, j] += omega*resid[i, j]*.25 

494 for i in range(2, func.shape[0] - 1, 2): 

495 for j in range(1, func.shape[1] - 1, 2): 

496 resid[i, j] = float(func[i, j - 1] + func[i, j + 1] + func[i - 1, j] 

497 + func[i + 1, j] - 4.0*func[i, j] - dx*dx*source[i - 1, j - 1]) 

498 func[i, j] += omega*resid[i, j]*.25 

499 outError = np.sum(np.abs(resid)) 

500 if outError < inError*eLevel: 

501 break 

502 if nIter == 0: 

503 omega = 1.0/(1 - rhoSpe*rhoSpe/2.0) 

504 else: 

505 omega = 1.0/(1 - rhoSpe*rhoSpe*omega/4.0) 

506 nIter += 1 

507 

508 if nIter >= maxIter*2: 

509 self.log.warn("Failure: SuccessiveOverRelaxation did not converge in %s iterations." 

510 "\noutError: %s, inError: %s," % (nIter//2, outError, inError*eLevel)) 

511 else: 

512 self.log.info("Success: SuccessiveOverRelaxation converged in %s iterations." 

513 "\noutError: %s, inError: %s", nIter//2, outError, inError*eLevel) 

514 return func[1: -1, 1: -1]