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 

23__all__ = ['PairedVisitListTaskRunner', 'SingleVisitListTaskRunner', 

24 'NonexistentDatasetTaskDataIdContainer', 'parseCmdlineNumberString', 

25 'countMaskedPixels', 'checkExpLengthEqual', 'ddict2dict'] 

26 

27import re 

28import numpy as np 

29from scipy.optimize import leastsq 

30import numpy.polynomial.polynomial as poly 

31from scipy.stats import norm 

32 

33import lsst.pipe.base as pipeBase 

34import lsst.ip.isr as ipIsr 

35from lsst.ip.isr import isrMock 

36import lsst.log 

37import lsst.afw.image 

38 

39import galsim 

40 

41 

42def sigmaClipCorrection(nSigClip): 

43 """Correct measured sigma to account for clipping. 

44 

45 If we clip our input data and then measure sigma, then the 

46 measured sigma is smaller than the true value because real 

47 points beyond the clip threshold have been removed. This is a 

48 small (1.5% at nSigClip=3) effect when nSigClip >~ 3, but the 

49 default parameters for measure crosstalk use nSigClip=2.0. 

50 This causes the measured sigma to be about 15% smaller than 

51 real. This formula corrects the issue, for the symmetric case 

52 (upper clip threshold equal to lower clip threshold). 

53 

54 Parameters 

55 ---------- 

56 nSigClip : `float` 

57 Number of sigma the measurement was clipped by. 

58 

59 Returns 

60 ------- 

61 scaleFactor : `float` 

62 Scale factor to increase the measured sigma by. 

63 

64 """ 

65 varFactor = 1.0 - (2 * nSigClip * norm.pdf(nSigClip)) / (norm.cdf(nSigClip) - norm.cdf(-nSigClip)) 

66 return 1.0 / np.sqrt(varFactor) 

67 

68 

69def calculateWeightedReducedChi2(measured, model, weightsMeasured, nData, nParsModel): 

70 """Calculate weighted reduced chi2. 

71 

72 Parameters 

73 ---------- 

74 

75 measured : `list` 

76 List with measured data. 

77 

78 model : `list` 

79 List with modeled data. 

80 

81 weightsMeasured : `list` 

82 List with weights for the measured data. 

83 

84 nData : `int` 

85 Number of data points. 

86 

87 nParsModel : `int` 

88 Number of parameters in the model. 

89 

90 Returns 

91 ------- 

92 

93 redWeightedChi2 : `float` 

94 Reduced weighted chi2. 

95 """ 

96 

97 wRes = (measured - model)*weightsMeasured 

98 return ((wRes*wRes).sum())/(nData-nParsModel) 

99 

100 

101def makeMockFlats(expTime, gain=1.0, readNoiseElectrons=5, fluxElectrons=1000, 

102 randomSeedFlat1=1984, randomSeedFlat2=666, powerLawBfParams=[], 

103 expId1=0, expId2=1): 

104 """Create a pair or mock flats with isrMock. 

105 

106 Parameters 

107 ---------- 

108 expTime : `float` 

109 Exposure time of the flats. 

110 

111 gain : `float`, optional 

112 Gain, in e/ADU. 

113 

114 readNoiseElectrons : `float`, optional 

115 Read noise rms, in electrons. 

116 

117 fluxElectrons : `float`, optional 

118 Flux of flats, in electrons per second. 

119 

120 randomSeedFlat1 : `int`, optional 

121 Random seed for the normal distrubutions for the mean signal and noise (flat1). 

122 

123 randomSeedFlat2 : `int`, optional 

124 Random seed for the normal distrubutions for the mean signal and noise (flat2). 

125 

126 powerLawBfParams : `list`, optional 

127 Parameters for `galsim.cdmodel.PowerLawCD` to simulate the brightter-fatter effect. 

128 

129 expId1 : `int`, optional 

130 Exposure ID for first flat. 

131 

132 expId2 : `int`, optional 

133 Exposure ID for second flat. 

134 

135 Returns 

136 ------- 

137 

138 flatExp1 : `lsst.afw.image.exposure.exposure.ExposureF` 

139 First exposure of flat field pair. 

140 

141 flatExp2 : `lsst.afw.image.exposure.exposure.ExposureF` 

142 Second exposure of flat field pair. 

143 

144 Notes 

145 ----- 

146 The parameters of `galsim.cdmodel.PowerLawCD` are `n, r0, t0, rx, tx, r, t, alpha`. For more 

147 information about their meaning, see the Galsim documentation 

148 https://galsim-developers.github.io/GalSim/_build/html/_modules/galsim/cdmodel.html 

149 and Gruen+15 (1501.02802). 

150 

151 Example: galsim.cdmodel.PowerLawCD(8, 1.1e-7, 1.1e-7, 1.0e-8, 1.0e-8, 1.0e-9, 1.0e-9, 2.0) 

152 """ 

153 flatFlux = fluxElectrons # e/s 

154 flatMean = flatFlux*expTime # e 

155 readNoise = readNoiseElectrons # e 

156 

157 mockImageConfig = isrMock.IsrMock.ConfigClass() 

158 

159 mockImageConfig.flatDrop = 0.99999 

160 mockImageConfig.isTrimmed = True 

161 

162 flatExp1 = isrMock.FlatMock(config=mockImageConfig).run() 

163 flatExp2 = flatExp1.clone() 

164 (shapeY, shapeX) = flatExp1.getDimensions() 

165 flatWidth = np.sqrt(flatMean) 

166 

167 rng1 = np.random.RandomState(randomSeedFlat1) 

168 flatData1 = rng1.normal(flatMean, flatWidth, (shapeX, shapeY)) + rng1.normal(0.0, readNoise, 

169 (shapeX, shapeY)) 

170 rng2 = np.random.RandomState(randomSeedFlat2) 

171 flatData2 = rng2.normal(flatMean, flatWidth, (shapeX, shapeY)) + rng2.normal(0.0, readNoise, 

172 (shapeX, shapeY)) 

173 # Simulate BF with power law model in galsim 

174 if len(powerLawBfParams): 

175 if not len(powerLawBfParams) == 8: 

176 raise RuntimeError("Wrong number of parameters for `galsim.cdmodel.PowerLawCD`. " 

177 f"Expected 8; passed {len(powerLawBfParams)}.") 

178 cd = galsim.cdmodel.PowerLawCD(*powerLawBfParams) 

179 tempFlatData1 = galsim.Image(flatData1) 

180 temp2FlatData1 = cd.applyForward(tempFlatData1) 

181 

182 tempFlatData2 = galsim.Image(flatData2) 

183 temp2FlatData2 = cd.applyForward(tempFlatData2) 

184 

185 flatExp1.image.array[:] = temp2FlatData1.array/gain # ADU 

186 flatExp2.image.array[:] = temp2FlatData2.array/gain # ADU 

187 else: 

188 flatExp1.image.array[:] = flatData1/gain # ADU 

189 flatExp2.image.array[:] = flatData2/gain # ADU 

190 

191 visitInfoExp1 = lsst.afw.image.VisitInfo(exposureId=expId1, exposureTime=expTime) 

192 visitInfoExp2 = lsst.afw.image.VisitInfo(exposureId=expId2, exposureTime=expTime) 

193 

194 flatExp1.getInfo().setVisitInfo(visitInfoExp1) 

195 flatExp2.getInfo().setVisitInfo(visitInfoExp2) 

196 

197 return flatExp1, flatExp2 

198 

199 

200def countMaskedPixels(maskedIm, maskPlane): 

201 """Count the number of pixels in a given mask plane.""" 

202 maskBit = maskedIm.mask.getPlaneBitMask(maskPlane) 

203 nPix = np.where(np.bitwise_and(maskedIm.mask.array, maskBit))[0].flatten().size 

204 return nPix 

205 

206 

207class PairedVisitListTaskRunner(pipeBase.TaskRunner): 

208 """Subclass of TaskRunner for handling intrinsically paired visits. 

209 

210 This transforms the processed arguments generated by the ArgumentParser 

211 into the arguments expected by tasks which take visit pairs for their 

212 run() methods. 

213 

214 Such tasks' run() methods tend to take two arguments, 

215 one of which is the dataRef (as usual), and the other is the list 

216 of visit-pairs, in the form of a list of tuples. 

217 This list is supplied on the command line as documented, 

218 and this class parses that, and passes the parsed version 

219 to the run() method. 

220 

221 See pipeBase.TaskRunner for more information. 

222 """ 

223 

224 @staticmethod 

225 def getTargetList(parsedCmd, **kwargs): 

226 """Parse the visit list and pass through explicitly.""" 

227 visitPairs = [] 

228 for visitStringPair in parsedCmd.visitPairs: 

229 visitStrings = visitStringPair.split(",") 

230 if len(visitStrings) != 2: 

231 raise RuntimeError("Found {} visits in {} instead of 2".format(len(visitStrings), 

232 visitStringPair)) 

233 try: 

234 visits = [int(visit) for visit in visitStrings] 

235 except Exception: 

236 raise RuntimeError("Could not parse {} as two integer visit numbers".format(visitStringPair)) 

237 visitPairs.append(visits) 

238 

239 return pipeBase.TaskRunner.getTargetList(parsedCmd, visitPairs=visitPairs, **kwargs) 

240 

241 

242def parseCmdlineNumberString(inputString): 

243 """Parse command line numerical expression sytax and return as list of int 

244 

245 Take an input of the form "'1..5:2^123..126'" as a string, and return 

246 a list of ints as [1, 3, 5, 123, 124, 125, 126] 

247 """ 

248 outList = [] 

249 for subString in inputString.split("^"): 

250 mat = re.search(r"^(\d+)\.\.(\d+)(?::(\d+))?$", subString) 

251 if mat: 

252 v1 = int(mat.group(1)) 

253 v2 = int(mat.group(2)) 

254 v3 = mat.group(3) 

255 v3 = int(v3) if v3 else 1 

256 for v in range(v1, v2 + 1, v3): 

257 outList.append(int(v)) 

258 else: 

259 outList.append(int(subString)) 

260 return outList 

261 

262 

263class SingleVisitListTaskRunner(pipeBase.TaskRunner): 

264 """Subclass of TaskRunner for tasks requiring a list of visits per dataRef. 

265 

266 This transforms the processed arguments generated by the ArgumentParser 

267 into the arguments expected by tasks which require a list of visits 

268 to be supplied for each dataRef, as is common in `lsst.cp.pipe` code. 

269 

270 Such tasks' run() methods tend to take two arguments, 

271 one of which is the dataRef (as usual), and the other is the list 

272 of visits. 

273 This list is supplied on the command line as documented, 

274 and this class parses that, and passes the parsed version 

275 to the run() method. 

276 

277 See `lsst.pipe.base.TaskRunner` for more information. 

278 """ 

279 

280 @staticmethod 

281 def getTargetList(parsedCmd, **kwargs): 

282 """Parse the visit list and pass through explicitly.""" 

283 # if this has been pre-parsed and therefore doesn't have length of one 

284 # then something has gone wrong, so execution should stop here. 

285 assert len(parsedCmd.visitList) == 1, 'visitList parsing assumptions violated' 

286 visits = parseCmdlineNumberString(parsedCmd.visitList[0]) 

287 

288 return pipeBase.TaskRunner.getTargetList(parsedCmd, visitList=visits, **kwargs) 

289 

290 

291class NonexistentDatasetTaskDataIdContainer(pipeBase.DataIdContainer): 

292 """A DataIdContainer for the tasks for which the output does 

293 not yet exist.""" 

294 

295 def makeDataRefList(self, namespace): 

296 """Compute refList based on idList. 

297 

298 This method must be defined as the dataset does not exist before this 

299 task is run. 

300 

301 Parameters 

302 ---------- 

303 namespace 

304 Results of parsing the command-line. 

305 

306 Notes 

307 ----- 

308 Not called if ``add_id_argument`` called 

309 with ``doMakeDataRefList=False``. 

310 Note that this is almost a copy-and-paste of the vanilla 

311 implementation, but without checking if the datasets already exist, 

312 as this task exists to make them. 

313 """ 

314 if self.datasetType is None: 

315 raise RuntimeError("Must call setDatasetType first") 

316 butler = namespace.butler 

317 for dataId in self.idList: 

318 refList = list(butler.subset(datasetType=self.datasetType, level=self.level, dataId=dataId)) 

319 # exclude nonexistent data 

320 # this is a recursive test, e.g. for the sake of "raw" data 

321 if not refList: 

322 namespace.log.warn("No data found for dataId=%s", dataId) 

323 continue 

324 self.refList += refList 

325 

326 

327def irlsFit(initialParams, dataX, dataY, function, weightsY=None): 

328 """Iteratively reweighted least squares fit. 

329 

330 This uses the `lsst.cp.pipe.utils.fitLeastSq`, but applies 

331 weights based on the Cauchy distribution to the fitter. See 

332 e.g. Holland and Welsch, 1977, doi:10.1080/03610927708827533 

333 

334 Parameters 

335 ---------- 

336 initialParams : `list` [`float`] 

337 Starting parameters. 

338 dataX : `numpy.array` [`float`] 

339 Abscissa data. 

340 dataY : `numpy.array` [`float`] 

341 Ordinate data. 

342 function : callable 

343 Function to fit. 

344 weightsY : `numpy.array` [`float`] 

345 Weights to apply to the data. 

346 

347 Returns 

348 ------- 

349 polyFit : `list` [`float`] 

350 Final best fit parameters. 

351 polyFitErr : `list` [`float`] 

352 Final errors on fit parameters. 

353 chiSq : `float` 

354 Reduced chi squared. 

355 weightsY : `list` [`float`] 

356 Final weights used for each point. 

357 

358 """ 

359 if not weightsY: 

360 weightsY = np.ones_like(dataX) 

361 

362 polyFit, polyFitErr, chiSq = fitLeastSq(initialParams, dataX, dataY, function, weightsY=weightsY) 

363 for iteration in range(10): 

364 # Use Cauchy weights 

365 resid = np.abs(dataY - function(polyFit, dataX)) / np.sqrt(dataY) 

366 weightsY = 1.0 / (1.0 + np.sqrt(resid / 2.385)) 

367 polyFit, polyFitErr, chiSq = fitLeastSq(initialParams, dataX, dataY, function, weightsY=weightsY) 

368 

369 return polyFit, polyFitErr, chiSq, weightsY 

370 

371 

372def fitLeastSq(initialParams, dataX, dataY, function, weightsY=None): 

373 """Do a fit and estimate the parameter errors using using scipy.optimize.leastq. 

374 

375 optimize.leastsq returns the fractional covariance matrix. To estimate the 

376 standard deviation of the fit parameters, multiply the entries of this matrix 

377 by the unweighted reduced chi squared and take the square root of the diagonal elements. 

378 

379 Parameters 

380 ---------- 

381 initialParams : `list` of `float` 

382 initial values for fit parameters. For ptcFitType=POLYNOMIAL, its length 

383 determines the degree of the polynomial. 

384 

385 dataX : `numpy.array` of `float` 

386 Data in the abscissa axis. 

387 

388 dataY : `numpy.array` of `float` 

389 Data in the ordinate axis. 

390 

391 function : callable object (function) 

392 Function to fit the data with. 

393 

394 weightsY : `numpy.array` of `float` 

395 Weights of the data in the ordinate axis. 

396 

397 Return 

398 ------ 

399 pFitSingleLeastSquares : `list` of `float` 

400 List with fitted parameters. 

401 

402 pErrSingleLeastSquares : `list` of `float` 

403 List with errors for fitted parameters. 

404 

405 reducedChiSqSingleLeastSquares : `float` 

406 Reduced chi squared, unweighted if weightsY is not provided. 

407 """ 

408 if weightsY is None: 

409 weightsY = np.ones(len(dataX)) 

410 

411 def errFunc(p, x, y, weightsY=None): 

412 if weightsY is None: 

413 weightsY = np.ones(len(x)) 

414 return (function(p, x) - y)*weightsY 

415 

416 pFit, pCov, infoDict, errMessage, success = leastsq(errFunc, initialParams, 

417 args=(dataX, dataY, weightsY), full_output=1, 

418 epsfcn=0.0001) 

419 

420 if (len(dataY) > len(initialParams)) and pCov is not None: 

421 reducedChiSq = calculateWeightedReducedChi2(dataY, function(pFit, dataX), weightsY, len(dataY), 

422 len(initialParams)) 

423 pCov *= reducedChiSq 

424 else: 

425 pCov = np.zeros((len(initialParams), len(initialParams))) 

426 pCov[:, :] = np.nan 

427 reducedChiSq = np.nan 

428 

429 errorVec = [] 

430 for i in range(len(pFit)): 

431 errorVec.append(np.fabs(pCov[i][i])**0.5) 

432 

433 pFitSingleLeastSquares = pFit 

434 pErrSingleLeastSquares = np.array(errorVec) 

435 

436 return pFitSingleLeastSquares, pErrSingleLeastSquares, reducedChiSq 

437 

438 

439def fitBootstrap(initialParams, dataX, dataY, function, weightsY=None, confidenceSigma=1.): 

440 """Do a fit using least squares and bootstrap to estimate parameter errors. 

441 

442 The bootstrap error bars are calculated by fitting 100 random data sets. 

443 

444 Parameters 

445 ---------- 

446 initialParams : `list` of `float` 

447 initial values for fit parameters. For ptcFitType=POLYNOMIAL, its length 

448 determines the degree of the polynomial. 

449 

450 dataX : `numpy.array` of `float` 

451 Data in the abscissa axis. 

452 

453 dataY : `numpy.array` of `float` 

454 Data in the ordinate axis. 

455 

456 function : callable object (function) 

457 Function to fit the data with. 

458 

459 weightsY : `numpy.array` of `float`, optional. 

460 Weights of the data in the ordinate axis. 

461 

462 confidenceSigma : `float`, optional. 

463 Number of sigmas that determine confidence interval for the bootstrap errors. 

464 

465 Return 

466 ------ 

467 pFitBootstrap : `list` of `float` 

468 List with fitted parameters. 

469 

470 pErrBootstrap : `list` of `float` 

471 List with errors for fitted parameters. 

472 

473 reducedChiSqBootstrap : `float` 

474 Reduced chi squared, unweighted if weightsY is not provided. 

475 """ 

476 if weightsY is None: 

477 weightsY = np.ones(len(dataX)) 

478 

479 def errFunc(p, x, y, weightsY): 

480 if weightsY is None: 

481 weightsY = np.ones(len(x)) 

482 return (function(p, x) - y)*weightsY 

483 

484 # Fit first time 

485 pFit, _ = leastsq(errFunc, initialParams, args=(dataX, dataY, weightsY), full_output=0) 

486 

487 # Get the stdev of the residuals 

488 residuals = errFunc(pFit, dataX, dataY, weightsY) 

489 # 100 random data sets are generated and fitted 

490 pars = [] 

491 for i in range(100): 

492 randomDelta = np.random.normal(0., np.fabs(residuals), len(dataY)) 

493 randomDataY = dataY + randomDelta 

494 randomFit, _ = leastsq(errFunc, initialParams, 

495 args=(dataX, randomDataY, weightsY), full_output=0) 

496 pars.append(randomFit) 

497 pars = np.array(pars) 

498 meanPfit = np.mean(pars, 0) 

499 

500 # confidence interval for parameter estimates 

501 errPfit = confidenceSigma*np.std(pars, 0) 

502 pFitBootstrap = meanPfit 

503 pErrBootstrap = errPfit 

504 

505 reducedChiSq = calculateWeightedReducedChi2(dataY, function(pFitBootstrap, dataX), weightsY, len(dataY), 

506 len(initialParams)) 

507 return pFitBootstrap, pErrBootstrap, reducedChiSq 

508 

509 

510def funcPolynomial(pars, x): 

511 """Polynomial function definition 

512 Parameters 

513 ---------- 

514 params : `list` 

515 Polynomial coefficients. Its length determines the polynomial order. 

516 

517 x : `numpy.array` 

518 Abscisa array. 

519 

520 Returns 

521 ------- 

522 Ordinate array after evaluating polynomial of order len(pars)-1 at `x`. 

523 """ 

524 return poly.polyval(x, [*pars]) 

525 

526 

527def funcAstier(pars, x): 

528 """Single brighter-fatter parameter model for PTC; Equation 16 of Astier+19. 

529 

530 Parameters 

531 ---------- 

532 params : `list` 

533 Parameters of the model: a00 (brightter-fatter), gain (e/ADU), and noise (e^2). 

534 

535 x : `numpy.array` 

536 Signal mu (ADU). 

537 

538 Returns 

539 ------- 

540 C_00 (variance) in ADU^2. 

541 """ 

542 a00, gain, noise = pars 

543 return 0.5/(a00*gain*gain)*(np.exp(2*a00*x*gain)-1) + noise/(gain*gain) # C_00 

544 

545 

546def arrangeFlatsByExpTime(exposureList): 

547 """Arrange exposures by exposure time. 

548 

549 Parameters 

550 ---------- 

551 exposureList : `list`[`lsst.afw.image.exposure.exposure.ExposureF`] 

552 Input list of exposures. 

553 

554 Returns 

555 ------ 

556 flatsAtExpTime : `dict` [`float`, 

557 `list`[`lsst.afw.image.exposure.exposure.ExposureF`]] 

558 Dictionary that groups flat-field exposures that have the same 

559 exposure time (seconds). 

560 """ 

561 flatsAtExpTime = {} 

562 for exp in exposureList: 

563 tempFlat = exp 

564 expTime = tempFlat.getInfo().getVisitInfo().getExposureTime() 

565 listAtExpTime = flatsAtExpTime.setdefault(expTime, []) 

566 listAtExpTime.append(tempFlat) 

567 

568 return flatsAtExpTime 

569 

570 

571def arrangeFlatsByExpId(exposureList): 

572 """Arrange exposures by exposure ID. 

573 

574 There is no guarantee that this will properly group exposures, but 

575 allows a sequence of flats that have different illumination 

576 (despite having the same exposure time) to be processed. 

577 

578 Parameters 

579 ---------- 

580 exposureList : `list`[`lsst.afw.image.exposure.exposure.ExposureF`] 

581 Input list of exposures. 

582 

583 Returns 

584 ------ 

585 flatsAtExpId : `dict` [`float`, 

586 `list`[`lsst.afw.image.exposure.exposure.ExposureF`]] 

587 Dictionary that groups flat-field exposures sequentially by 

588 their exposure id. 

589 

590 Notes 

591 ----- 

592 

593 This algorithm sorts the input exposures by their exposure id, and 

594 then assigns each pair of exposures (exp_j, exp_{j+1}) to pair k, 

595 such that 2*k = j, where j is the python index of one of the 

596 exposures (starting from zero). By checking for the IndexError 

597 while appending, we can ensure that there will only ever be fully 

598 populated pairs. 

599 """ 

600 flatsAtExpId = {} 

601 sortedExposures = sorted(exposureList, key=lambda exp: exp.getInfo().getVisitInfo().getExposureId()) 

602 

603 for jPair, exp in enumerate(sortedExposures): 

604 if (jPair + 1) % 2: 

605 kPair = jPair // 2 

606 listAtExpId = flatsAtExpId.setdefault(kPair, []) 

607 try: 

608 listAtExpId.append(exp) 

609 listAtExpId.append(sortedExposures[jPair + 1]) 

610 except IndexError: 

611 pass 

612 

613 return flatsAtExpId 

614 

615 

616def checkExpLengthEqual(exp1, exp2, v1=None, v2=None, raiseWithMessage=False): 

617 """Check the exposure lengths of two exposures are equal. 

618 

619 Parameters: 

620 ----------- 

621 exp1 : `lsst.afw.image.exposure.ExposureF` 

622 First exposure to check 

623 exp2 : `lsst.afw.image.exposure.ExposureF` 

624 Second exposure to check 

625 v1 : `int` or `str`, optional 

626 First visit of the visit pair 

627 v2 : `int` or `str`, optional 

628 Second visit of the visit pair 

629 raiseWithMessage : `bool` 

630 If True, instead of returning a bool, raise a RuntimeError if exposure 

631 times are not equal, with a message about which visits mismatch if the 

632 information is available. 

633 

634 Raises: 

635 ------- 

636 RuntimeError 

637 Raised if the exposure lengths of the two exposures are not equal 

638 """ 

639 expTime1 = exp1.getInfo().getVisitInfo().getExposureTime() 

640 expTime2 = exp2.getInfo().getVisitInfo().getExposureTime() 

641 if expTime1 != expTime2: 

642 if raiseWithMessage: 

643 msg = "Exposure lengths for visit pairs must be equal. " + \ 

644 "Found %s and %s" % (expTime1, expTime2) 

645 if v1 and v2: 

646 msg += " for visit pair %s, %s" % (v1, v2) 

647 raise RuntimeError(msg) 

648 else: 

649 return False 

650 return True 

651 

652 

653def validateIsrConfig(isrTask, mandatory=None, forbidden=None, desirable=None, undesirable=None, 

654 checkTrim=True, logName=None): 

655 """Check that appropriate ISR settings have been selected for the task. 

656 

657 Note that this checks that the task itself is configured correctly rather 

658 than checking a config. 

659 

660 Parameters 

661 ---------- 

662 isrTask : `lsst.ip.isr.IsrTask` 

663 The task whose config is to be validated 

664 

665 mandatory : `iterable` of `str` 

666 isr steps that must be set to True. Raises if False or missing 

667 

668 forbidden : `iterable` of `str` 

669 isr steps that must be set to False. Raises if True, warns if missing 

670 

671 desirable : `iterable` of `str` 

672 isr steps that should probably be set to True. Warns is False, info if 

673 missing 

674 

675 undesirable : `iterable` of `str` 

676 isr steps that should probably be set to False. Warns is True, info if 

677 missing 

678 

679 checkTrim : `bool` 

680 Check to ensure the isrTask's assembly subtask is trimming the images. 

681 This is a separate config as it is very ugly to do this within the 

682 normal configuration lists as it is an option of a sub task. 

683 

684 Raises 

685 ------ 

686 RuntimeError 

687 Raised if ``mandatory`` config parameters are False, 

688 or if ``forbidden`` parameters are True. 

689 

690 TypeError 

691 Raised if parameter ``isrTask`` is an invalid type. 

692 

693 Notes 

694 ----- 

695 Logs warnings using an isrValidation logger for desirable/undesirable 

696 options that are of the wrong polarity or if keys are missing. 

697 """ 

698 if not isinstance(isrTask, ipIsr.IsrTask): 

699 raise TypeError(f'Must supply an instance of lsst.ip.isr.IsrTask not {type(isrTask)}') 

700 

701 configDict = isrTask.config.toDict() 

702 

703 if logName and isinstance(logName, str): 

704 log = lsst.log.getLogger(logName) 

705 else: 

706 log = lsst.log.getLogger("isrValidation") 

707 

708 if mandatory: 

709 for configParam in mandatory: 

710 if configParam not in configDict: 

711 raise RuntimeError(f"Mandatory parameter {configParam} not found in the isr configuration.") 

712 if configDict[configParam] is False: 

713 raise RuntimeError(f"Must set config.isr.{configParam} to True for this task.") 

714 

715 if forbidden: 

716 for configParam in forbidden: 

717 if configParam not in configDict: 

718 log.warn(f"Failed to find forbidden key {configParam} in the isr config. The keys in the" 

719 " forbidden list should each have an associated Field in IsrConfig:" 

720 " check that there is not a typo in this case.") 

721 continue 

722 if configDict[configParam] is True: 

723 raise RuntimeError(f"Must set config.isr.{configParam} to False for this task.") 

724 

725 if desirable: 

726 for configParam in desirable: 

727 if configParam not in configDict: 

728 log.info(f"Failed to find key {configParam} in the isr config. You probably want" 

729 " to set the equivalent for your obs_package to True.") 

730 continue 

731 if configDict[configParam] is False: 

732 log.warn(f"Found config.isr.{configParam} set to False for this task." 

733 " The cp_pipe Config recommends setting this to True.") 

734 if undesirable: 

735 for configParam in undesirable: 

736 if configParam not in configDict: 

737 log.info(f"Failed to find key {configParam} in the isr config. You probably want" 

738 " to set the equivalent for your obs_package to False.") 

739 continue 

740 if configDict[configParam] is True: 

741 log.warn(f"Found config.isr.{configParam} set to True for this task." 

742 " The cp_pipe Config recommends setting this to False.") 

743 

744 if checkTrim: # subtask setting, seems non-trivial to combine with above lists 

745 if not isrTask.assembleCcd.config.doTrim: 

746 raise RuntimeError("Must trim when assembling CCDs. Set config.isr.assembleCcd.doTrim to True") 

747 

748 

749def ddict2dict(d): 

750 """Convert nested default dictionaries to regular dictionaries. 

751 

752 This is needed to prevent yaml persistence issues. 

753 

754 Parameters 

755 ---------- 

756 d : `defaultdict` 

757 A possibly nested set of `defaultdict`. 

758 

759 Returns 

760 ------- 

761 dict : `dict` 

762 A possibly nested set of `dict`. 

763 """ 

764 for k, v in d.items(): 

765 if isinstance(v, dict): 

766 d[k] = ddict2dict(v) 

767 return dict(d)