Coverage for python/lsst/drp/tasks/assemble_coadd.py: 16%

633 statements  

« prev     ^ index     » next       coverage.py v7.4.4, created at 2024-04-15 04:01 -0700

1# This file is part of drp_tasks. 

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__all__ = [ 

23 "AssembleCoaddTask", 

24 "AssembleCoaddConnections", 

25 "AssembleCoaddConfig", 

26 "CompareWarpAssembleCoaddTask", 

27 "CompareWarpAssembleCoaddConfig", 

28] 

29 

30import copy 

31import logging 

32import warnings 

33 

34import lsst.afw.geom as afwGeom 

35import lsst.afw.image as afwImage 

36import lsst.afw.math as afwMath 

37import lsst.afw.table as afwTable 

38import lsst.coadd.utils as coaddUtils 

39import lsst.geom as geom 

40import lsst.meas.algorithms as measAlg 

41import lsst.pex.config as pexConfig 

42import lsst.pex.exceptions as pexExceptions 

43import lsst.pipe.base as pipeBase 

44import lsst.utils as utils 

45import lsstDebug 

46import numpy 

47from deprecated.sphinx import deprecated 

48from lsst.meas.algorithms import AccumulatorMeanStack, ScaleVarianceTask, SourceDetectionTask 

49from lsst.pipe.tasks.coaddBase import CoaddBaseTask, makeSkyInfo, reorderAndPadList, subBBoxIter 

50from lsst.pipe.tasks.healSparseMapping import HealSparseInputMapTask 

51from lsst.pipe.tasks.interpImage import InterpImageTask 

52from lsst.pipe.tasks.maskStreaks import MaskStreaksTask 

53from lsst.pipe.tasks.scaleZeroPoint import ScaleZeroPointTask 

54from lsst.skymap import BaseSkyMap 

55from lsst.utils.timer import timeMethod 

56 

57log = logging.getLogger(__name__) 

58 

59 

60class AssembleCoaddConnections( 

61 pipeBase.PipelineTaskConnections, 

62 dimensions=("tract", "patch", "band", "skymap"), 

63 defaultTemplates={ 

64 "inputCoaddName": "deep", 

65 "outputCoaddName": "deep", 

66 "warpType": "direct", 

67 "warpTypeSuffix": "", 

68 }, 

69): 

70 inputWarps = pipeBase.connectionTypes.Input( 

71 doc=( 

72 "Input list of warps to be assemebled i.e. stacked." 

73 "WarpType (e.g. direct, psfMatched) is controlled by the " 

74 "warpType config parameter" 

75 ), 

76 name="{inputCoaddName}Coadd_{warpType}Warp", 

77 storageClass="ExposureF", 

78 dimensions=("tract", "patch", "skymap", "visit", "instrument"), 

79 deferLoad=True, 

80 multiple=True, 

81 ) 

82 skyMap = pipeBase.connectionTypes.Input( 

83 doc="Input definition of geometry/bbox and projection/wcs for coadded " "exposures", 

84 name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME, 

85 storageClass="SkyMap", 

86 dimensions=("skymap",), 

87 ) 

88 selectedVisits = pipeBase.connectionTypes.Input( 

89 doc="Selected visits to be coadded.", 

90 name="{outputCoaddName}Visits", 

91 storageClass="StructuredDataDict", 

92 dimensions=("instrument", "tract", "patch", "skymap", "band"), 

93 ) 

94 brightObjectMask = pipeBase.connectionTypes.PrerequisiteInput( 

95 doc=( 

96 "Input Bright Object Mask mask produced with external catalogs " 

97 "to be applied to the mask plane BRIGHT_OBJECT." 

98 ), 

99 name="brightObjectMask", 

100 storageClass="ObjectMaskCatalog", 

101 dimensions=("tract", "patch", "skymap", "band"), 

102 minimum=0, 

103 ) 

104 coaddExposure = pipeBase.connectionTypes.Output( 

105 doc="Output coadded exposure, produced by stacking input warps", 

106 name="{outputCoaddName}Coadd{warpTypeSuffix}", 

107 storageClass="ExposureF", 

108 dimensions=("tract", "patch", "skymap", "band"), 

109 ) 

110 nImage = pipeBase.connectionTypes.Output( 

111 doc="Output image of number of input images per pixel", 

112 name="{outputCoaddName}Coadd_nImage", 

113 storageClass="ImageU", 

114 dimensions=("tract", "patch", "skymap", "band"), 

115 ) 

116 inputMap = pipeBase.connectionTypes.Output( 

117 doc="Output healsparse map of input images", 

118 name="{outputCoaddName}Coadd_inputMap", 

119 storageClass="HealSparseMap", 

120 dimensions=("tract", "patch", "skymap", "band"), 

121 ) 

122 

123 def __init__(self, *, config=None): 

124 super().__init__(config=config) 

125 

126 if not config.doMaskBrightObjects: 

127 self.prerequisiteInputs.remove("brightObjectMask") 

128 

129 if not config.doSelectVisits: 

130 self.inputs.remove("selectedVisits") 

131 

132 if not config.doNImage: 

133 self.outputs.remove("nImage") 

134 

135 if not self.config.doInputMap: 

136 self.outputs.remove("inputMap") 

137 

138 

139class AssembleCoaddConfig( 

140 CoaddBaseTask.ConfigClass, pipeBase.PipelineTaskConfig, pipelineConnections=AssembleCoaddConnections 

141): 

142 warpType = pexConfig.Field( 

143 doc="Warp name: one of 'direct' or 'psfMatched'", 

144 dtype=str, 

145 default="direct", 

146 ) 

147 subregionSize = pexConfig.ListField( 

148 dtype=int, 

149 doc="Width, height of stack subregion size; " 

150 "make small enough that a full stack of images will fit into memory " 

151 " at once.", 

152 length=2, 

153 default=(2000, 2000), 

154 ) 

155 statistic = pexConfig.Field( 

156 dtype=str, 

157 doc="Main stacking statistic for aggregating over the epochs.", 

158 default="MEANCLIP", 

159 ) 

160 doOnlineForMean = pexConfig.Field( 

161 dtype=bool, 

162 doc='Perform online coaddition when statistic="MEAN" to save memory?', 

163 default=False, 

164 ) 

165 doSigmaClip = pexConfig.Field( 

166 dtype=bool, 

167 doc="Perform sigma clipped outlier rejection with MEANCLIP statistic?", 

168 deprecated=True, 

169 default=False, 

170 ) 

171 sigmaClip = pexConfig.Field( 

172 dtype=float, 

173 doc="Sigma for outlier rejection; ignored if non-clipping statistic " "selected.", 

174 default=3.0, 

175 ) 

176 clipIter = pexConfig.Field( 

177 dtype=int, 

178 doc="Number of iterations of outlier rejection; ignored if " "non-clipping statistic selected.", 

179 default=2, 

180 ) 

181 calcErrorFromInputVariance = pexConfig.Field( 

182 dtype=bool, 

183 doc="Calculate coadd variance from input variance by stacking " 

184 "statistic. Passed to " 

185 "StatisticsControl.setCalcErrorFromInputVariance()", 

186 default=True, 

187 ) 

188 scaleZeroPoint = pexConfig.ConfigurableField( 

189 target=ScaleZeroPointTask, 

190 doc="Task to adjust the photometric zero point of the coadd temp " "exposures", 

191 ) 

192 doInterp = pexConfig.Field( 

193 doc="Interpolate over NaN pixels? Also extrapolate, if necessary, but " "the results are ugly.", 

194 dtype=bool, 

195 default=True, 

196 ) 

197 interpImage = pexConfig.ConfigurableField( 

198 target=InterpImageTask, 

199 doc="Task to interpolate (and extrapolate) over NaN pixels", 

200 ) 

201 doWrite = pexConfig.Field( 

202 doc="Persist coadd?", 

203 dtype=bool, 

204 default=True, 

205 ) 

206 doNImage = pexConfig.Field( 

207 doc="Create image of number of contributing exposures for each pixel", 

208 dtype=bool, 

209 default=False, 

210 ) 

211 doUsePsfMatchedPolygons = pexConfig.Field( 

212 doc="Use ValidPolygons from shrunk Psf-Matched Calexps? Should be set " 

213 "to True by CompareWarp only.", 

214 dtype=bool, 

215 default=False, 

216 ) 

217 maskPropagationThresholds = pexConfig.DictField( 

218 keytype=str, 

219 itemtype=float, 

220 doc=( 

221 "Threshold (in fractional weight) of rejection at which we " 

222 "propagate a mask plane to the coadd; that is, we set the mask " 

223 "bit on the coadd if the fraction the rejected frames " 

224 "would have contributed exceeds this value." 

225 ), 

226 default={"SAT": 0.1}, 

227 ) 

228 removeMaskPlanes = pexConfig.ListField( 

229 dtype=str, 

230 default=["NOT_DEBLENDED"], 

231 doc="Mask planes to remove before coadding", 

232 ) 

233 doMaskBrightObjects = pexConfig.Field( 

234 dtype=bool, 

235 default=False, 

236 doc="Set mask and flag bits for bright objects?", 

237 ) 

238 brightObjectMaskName = pexConfig.Field( 

239 dtype=str, 

240 default="BRIGHT_OBJECT", 

241 doc="Name of mask bit used for bright objects", 

242 ) 

243 coaddPsf = pexConfig.ConfigField( 

244 doc="Configuration for CoaddPsf", 

245 dtype=measAlg.CoaddPsfConfig, 

246 ) 

247 doAttachTransmissionCurve = pexConfig.Field( 

248 dtype=bool, 

249 default=False, 

250 optional=False, 

251 doc=( 

252 "Attach a piecewise TransmissionCurve for the coadd? " 

253 "(requires all input Exposures to have TransmissionCurves)." 

254 ), 

255 ) 

256 hasFakes = pexConfig.Field( 

257 dtype=bool, 

258 default=False, 

259 doc="Should be set to True if fake sources have been inserted into the input data.", 

260 ) 

261 doSelectVisits = pexConfig.Field( 

262 doc="Coadd only visits selected by a SelectVisitsTask", 

263 dtype=bool, 

264 default=False, 

265 ) 

266 doInputMap = pexConfig.Field( 

267 doc="Create a bitwise map of coadd inputs", 

268 dtype=bool, 

269 default=False, 

270 ) 

271 inputMapper = pexConfig.ConfigurableField( 

272 doc="Input map creation subtask.", 

273 target=HealSparseInputMapTask, 

274 ) 

275 

276 def setDefaults(self): 

277 super().setDefaults() 

278 self.badMaskPlanes = ["NO_DATA", "BAD", "SAT", "EDGE"] 

279 

280 def validate(self): 

281 super().validate() 

282 if self.doSigmaClip and self.statistic != "MEANCLIP": 

283 log.warning('doSigmaClip deprecated. To replicate behavior, setting statistic to "MEANCLIP"') 

284 self.statistic = "MEANCLIP" 

285 if self.doInterp and self.statistic not in ["MEAN", "MEDIAN", "MEANCLIP", "VARIANCE", "VARIANCECLIP"]: 

286 raise ValueError( 

287 "Must set doInterp=False for statistic=%s, which does not " 

288 "compute and set a non-zero coadd variance estimate." % (self.statistic) 

289 ) 

290 

291 unstackableStats = ["NOTHING", "ERROR", "ORMASK"] 

292 if not hasattr(afwMath.Property, self.statistic) or self.statistic in unstackableStats: 

293 stackableStats = [ 

294 str(k) for k in afwMath.Property.__members__.keys() if str(k) not in unstackableStats 

295 ] 

296 raise ValueError( 

297 "statistic %s is not allowed. Please choose one of %s." % (self.statistic, stackableStats) 

298 ) 

299 

300 

301class AssembleCoaddTask(CoaddBaseTask, pipeBase.PipelineTask): 

302 """Assemble a coadded image from a set of warps. 

303 

304 Each Warp that goes into a coadd will typically have an independent 

305 photometric zero-point. Therefore, we must scale each Warp to set it to 

306 a common photometric zeropoint. WarpType may be one of 'direct' or 

307 'psfMatched', and the boolean configs `config.makeDirect` and 

308 `config.makePsfMatched` set which of the warp types will be coadded. 

309 The coadd is computed as a mean with optional outlier rejection. 

310 Criteria for outlier rejection are set in `AssembleCoaddConfig`. 

311 Finally, Warps can have bad 'NaN' pixels which received no input from the 

312 source calExps. We interpolate over these bad (NaN) pixels. 

313 

314 `AssembleCoaddTask` uses several sub-tasks. These are 

315 

316 - `~lsst.pipe.tasks.ScaleZeroPointTask` 

317 - create and use an ``imageScaler`` object to scale the photometric 

318 zeropoint for each Warp 

319 - `~lsst.pipe.tasks.InterpImageTask` 

320 - interpolate across bad pixels (NaN) in the final coadd 

321 

322 You can retarget these subtasks if you wish. 

323 

324 Raises 

325 ------ 

326 RuntimeError 

327 Raised if unable to define mask plane for bright objects. 

328 

329 Notes 

330 ----- 

331 Debugging: 

332 `AssembleCoaddTask` has no debug variables of its own. Some of the 

333 subtasks may support `~lsst.base.lsstDebug` variables. See the 

334 documentation for the subtasks for further information. 

335 """ 

336 

337 ConfigClass = AssembleCoaddConfig 

338 _DefaultName = "assembleCoadd" 

339 

340 def __init__(self, *args, **kwargs): 

341 # TODO: DM-17415 better way to handle previously allowed passed args 

342 # e.g.`AssembleCoaddTask(config)` 

343 if args: 

344 argNames = ["config", "name", "parentTask", "log"] 

345 kwargs.update({k: v for k, v in zip(argNames, args)}) 

346 warnings.warn( 

347 "AssembleCoadd received positional args, and casting them as kwargs: %s. " 

348 "PipelineTask will not take positional args" % argNames, 

349 FutureWarning, 

350 stacklevel=2, 

351 ) 

352 

353 super().__init__(**kwargs) 

354 self.makeSubtask("interpImage") 

355 self.makeSubtask("scaleZeroPoint") 

356 

357 if self.config.doMaskBrightObjects: 

358 mask = afwImage.Mask() 

359 try: 

360 self.brightObjectBitmask = 1 << mask.addMaskPlane(self.config.brightObjectMaskName) 

361 except pexExceptions.LsstCppException: 

362 raise RuntimeError( 

363 "Unable to define mask plane for bright objects; planes used are %s" 

364 % mask.getMaskPlaneDict().keys() 

365 ) 

366 del mask 

367 

368 if self.config.doInputMap: 

369 self.makeSubtask("inputMapper") 

370 

371 self.warpType = self.config.warpType 

372 

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

374 inputData = butlerQC.get(inputRefs) 

375 

376 # Construct skyInfo expected by run 

377 # Do not remove skyMap from inputData in case _makeSupplementaryData 

378 # needs it 

379 skyMap = inputData["skyMap"] 

380 outputDataId = butlerQC.quantum.dataId 

381 

382 inputData["skyInfo"] = makeSkyInfo( 

383 skyMap, tractId=outputDataId["tract"], patchId=outputDataId["patch"] 

384 ) 

385 

386 if self.config.doSelectVisits: 

387 warpRefList = self.filterWarps(inputData["inputWarps"], inputData["selectedVisits"]) 

388 else: 

389 warpRefList = inputData["inputWarps"] 

390 

391 inputs = self.prepareInputs(warpRefList) 

392 self.log.info("Found %d %s", len(inputs.tempExpRefList), self.getTempExpDatasetName(self.warpType)) 

393 if len(inputs.tempExpRefList) == 0: 

394 raise pipeBase.NoWorkFound("No coadd temporary exposures found") 

395 

396 supplementaryData = self._makeSupplementaryData(butlerQC, inputRefs, outputRefs) 

397 retStruct = self.run( 

398 inputData["skyInfo"], 

399 inputs.tempExpRefList, 

400 inputs.imageScalerList, 

401 inputs.weightList, 

402 supplementaryData=supplementaryData, 

403 ) 

404 

405 inputData.setdefault("brightObjectMask", None) 

406 if self.config.doMaskBrightObjects and inputData["brightObjectMask"] is None: 

407 log.warning("doMaskBrightObjects is set to True, but brightObjectMask not loaded") 

408 self.processResults(retStruct.coaddExposure, inputData["brightObjectMask"], outputDataId) 

409 

410 if self.config.doWrite: 

411 butlerQC.put(retStruct, outputRefs) 

412 return retStruct 

413 

414 def processResults(self, coaddExposure, brightObjectMasks=None, dataId=None): 

415 """Interpolate over missing data and mask bright stars. 

416 

417 Parameters 

418 ---------- 

419 coaddExposure : `lsst.afw.image.Exposure` 

420 The coadded exposure to process. 

421 brightObjectMasks : `lsst.afw.table` or `None`, optional 

422 Table of bright objects to mask. 

423 dataId : `lsst.daf.butler.DataId` or `None`, optional 

424 Data identification. 

425 """ 

426 if self.config.doInterp: 

427 self.interpImage.run(coaddExposure.getMaskedImage(), planeName="NO_DATA") 

428 # The variance must be positive; work around for DM-3201. 

429 varArray = coaddExposure.variance.array 

430 with numpy.errstate(invalid="ignore"): 

431 varArray[:] = numpy.where(varArray > 0, varArray, numpy.inf) 

432 

433 if self.config.doMaskBrightObjects: 

434 self.setBrightObjectMasks(coaddExposure, brightObjectMasks, dataId) 

435 

436 def _makeSupplementaryData(self, butlerQC, inputRefs, outputRefs): 

437 """Make additional inputs to run() specific to subclasses (Gen3). 

438 

439 Duplicates interface of `runQuantum` method. 

440 Available to be implemented by subclasses only if they need the 

441 coadd dataRef for performing preliminary processing before 

442 assembling the coadd. 

443 

444 Parameters 

445 ---------- 

446 butlerQC : `~lsst.pipe.base.QuantumContext` 

447 Gen3 Butler object for fetching additional data products before 

448 running the Task specialized for quantum being processed. 

449 inputRefs : `~lsst.pipe.base.InputQuantizedConnection` 

450 Attributes are the names of the connections describing input 

451 dataset types. Values are DatasetRefs that task consumes for the 

452 corresponding dataset type. DataIds are guaranteed to match data 

453 objects in ``inputData``. 

454 outputRefs : `~lsst.pipe.base.OutputQuantizedConnection` 

455 Attributes are the names of the connections describing output 

456 dataset types. Values are DatasetRefs that task is to produce 

457 for the corresponding dataset type. 

458 """ 

459 return pipeBase.Struct() 

460 

461 @deprecated( 

462 reason="makeSupplementaryDataGen3 is deprecated in favor of _makeSupplementaryData", 

463 version="v25.0", 

464 category=FutureWarning, 

465 ) 

466 def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs): 

467 return self._makeSupplementaryData(butlerQC, inputRefs, outputRefs) 

468 

469 def prepareInputs(self, refList): 

470 """Prepare the input warps for coaddition by measuring the weight for 

471 each warp and the scaling for the photometric zero point. 

472 

473 Each Warp has its own photometric zeropoint and background variance. 

474 Before coadding these Warps together, compute a scale factor to 

475 normalize the photometric zeropoint and compute the weight for each 

476 Warp. 

477 

478 Parameters 

479 ---------- 

480 refList : `list` 

481 List of data references to tempExp. 

482 

483 Returns 

484 ------- 

485 result : `~lsst.pipe.base.Struct` 

486 Results as a struct with attributes: 

487 

488 ``tempExprefList`` 

489 `list` of data references to tempExp. 

490 ``weightList`` 

491 `list` of weightings. 

492 ``imageScalerList`` 

493 `list` of image scalers. 

494 """ 

495 statsCtrl = afwMath.StatisticsControl() 

496 statsCtrl.setNumSigmaClip(self.config.sigmaClip) 

497 statsCtrl.setNumIter(self.config.clipIter) 

498 statsCtrl.setAndMask(self.getBadPixelMask()) 

499 statsCtrl.setNanSafe(True) 

500 # compute tempExpRefList: a list of tempExpRef that actually exist 

501 # and weightList: a list of the weight of the associated coadd tempExp 

502 # and imageScalerList: a list of scale factors for the associated coadd 

503 # tempExp. 

504 tempExpRefList = [] 

505 weightList = [] 

506 imageScalerList = [] 

507 tempExpName = self.getTempExpDatasetName(self.warpType) 

508 for tempExpRef in refList: 

509 tempExp = tempExpRef.get() 

510 # Ignore any input warp that is empty of data 

511 if numpy.isnan(tempExp.image.array).all(): 

512 continue 

513 maskedImage = tempExp.getMaskedImage() 

514 imageScaler = self.scaleZeroPoint.computeImageScaler( 

515 exposure=tempExp, 

516 dataRef=tempExpRef, # FIXME 

517 ) 

518 try: 

519 imageScaler.scaleMaskedImage(maskedImage) 

520 except Exception as e: 

521 self.log.warning("Scaling failed for %s (skipping it): %s", tempExpRef.dataId, e) 

522 continue 

523 statObj = afwMath.makeStatistics( 

524 maskedImage.getVariance(), maskedImage.getMask(), afwMath.MEANCLIP, statsCtrl 

525 ) 

526 meanVar, meanVarErr = statObj.getResult(afwMath.MEANCLIP) 

527 weight = 1.0 / float(meanVar) 

528 if not numpy.isfinite(weight): 

529 self.log.warning("Non-finite weight for %s: skipping", tempExpRef.dataId) 

530 continue 

531 self.log.info("Weight of %s %s = %0.3f", tempExpName, tempExpRef.dataId, weight) 

532 

533 del maskedImage 

534 del tempExp 

535 

536 tempExpRefList.append(tempExpRef) 

537 weightList.append(weight) 

538 imageScalerList.append(imageScaler) 

539 

540 return pipeBase.Struct( 

541 tempExpRefList=tempExpRefList, weightList=weightList, imageScalerList=imageScalerList 

542 ) 

543 

544 def prepareStats(self, mask=None): 

545 """Prepare the statistics for coadding images. 

546 

547 Parameters 

548 ---------- 

549 mask : `int`, optional 

550 Bit mask value to exclude from coaddition. 

551 

552 Returns 

553 ------- 

554 stats : `~lsst.pipe.base.Struct` 

555 Statistics as a struct with attributes: 

556 

557 ``statsCtrl`` 

558 Statistics control object for coadd 

559 (`~lsst.afw.math.StatisticsControl`). 

560 ``statsFlags`` 

561 Statistic for coadd (`~lsst.afw.math.Property`). 

562 """ 

563 if mask is None: 

564 mask = self.getBadPixelMask() 

565 statsCtrl = afwMath.StatisticsControl() 

566 statsCtrl.setNumSigmaClip(self.config.sigmaClip) 

567 statsCtrl.setNumIter(self.config.clipIter) 

568 statsCtrl.setAndMask(mask) 

569 statsCtrl.setNanSafe(True) 

570 statsCtrl.setWeighted(True) 

571 statsCtrl.setCalcErrorFromInputVariance(self.config.calcErrorFromInputVariance) 

572 for plane, threshold in self.config.maskPropagationThresholds.items(): 

573 bit = afwImage.Mask.getMaskPlane(plane) 

574 statsCtrl.setMaskPropagationThreshold(bit, threshold) 

575 statsFlags = afwMath.stringToStatisticsProperty(self.config.statistic) 

576 return pipeBase.Struct(ctrl=statsCtrl, flags=statsFlags) 

577 

578 @timeMethod 

579 def run( 

580 self, 

581 skyInfo, 

582 tempExpRefList, 

583 imageScalerList, 

584 weightList, 

585 altMaskList=None, 

586 mask=None, 

587 supplementaryData=None, 

588 ): 

589 """Assemble a coadd from input warps. 

590 

591 Assemble the coadd using the provided list of coaddTempExps. Since 

592 the full coadd covers a patch (a large area), the assembly is 

593 performed over small areas on the image at a time in order to 

594 conserve memory usage. Iterate over subregions within the outer 

595 bbox of the patch using `assembleSubregion` to stack the corresponding 

596 subregions from the coaddTempExps with the statistic specified. 

597 Set the edge bits the coadd mask based on the weight map. 

598 

599 Parameters 

600 ---------- 

601 skyInfo : `~lsst.pipe.base.Struct` 

602 Struct with geometric information about the patch. 

603 tempExpRefList : `list` 

604 List of data references to Warps (previously called CoaddTempExps). 

605 imageScalerList : `list` 

606 List of image scalers. 

607 weightList : `list` 

608 List of weights. 

609 altMaskList : `list`, optional 

610 List of alternate masks to use rather than those stored with 

611 tempExp. 

612 mask : `int`, optional 

613 Bit mask value to exclude from coaddition. 

614 supplementaryData : `~lsst.pipe.base.Struct`, optional 

615 Struct with additional data products needed to assemble coadd. 

616 Only used by subclasses that implement ``_makeSupplementaryData`` 

617 and override `run`. 

618 

619 Returns 

620 ------- 

621 result : `~lsst.pipe.base.Struct` 

622 Results as a struct with attributes: 

623 

624 ``coaddExposure`` 

625 Coadded exposure (`~lsst.afw.image.Exposure`). 

626 ``nImage`` 

627 Exposure count image (`~lsst.afw.image.Image`), if requested. 

628 ``inputMap`` 

629 Bit-wise map of inputs, if requested. 

630 ``warpRefList`` 

631 Input list of refs to the warps 

632 (`~lsst.daf.butler.DeferredDatasetHandle`) (unmodified). 

633 ``imageScalerList`` 

634 Input list of image scalers (`list`) (unmodified). 

635 ``weightList`` 

636 Input list of weights (`list`) (unmodified). 

637 

638 Raises 

639 ------ 

640 lsst.pipe.base.NoWorkFound 

641 Raised if no data references are provided. 

642 """ 

643 tempExpName = self.getTempExpDatasetName(self.warpType) 

644 self.log.info("Assembling %s %s", len(tempExpRefList), tempExpName) 

645 if not tempExpRefList: 

646 raise pipeBase.NoWorkFound("No exposures provided for co-addition.") 

647 

648 stats = self.prepareStats(mask=mask) 

649 

650 if altMaskList is None: 

651 altMaskList = [None] * len(tempExpRefList) 

652 

653 coaddExposure = afwImage.ExposureF(skyInfo.bbox, skyInfo.wcs) 

654 coaddExposure.setPhotoCalib(self.scaleZeroPoint.getPhotoCalib()) 

655 coaddExposure.getInfo().setCoaddInputs(self.inputRecorder.makeCoaddInputs()) 

656 self.assembleMetadata(coaddExposure, tempExpRefList, weightList) 

657 coaddMaskedImage = coaddExposure.getMaskedImage() 

658 subregionSizeArr = self.config.subregionSize 

659 subregionSize = geom.Extent2I(subregionSizeArr[0], subregionSizeArr[1]) 

660 # if nImage is requested, create a zero one which can be passed to 

661 # assembleSubregion. 

662 if self.config.doNImage: 

663 nImage = afwImage.ImageU(skyInfo.bbox) 

664 else: 

665 nImage = None 

666 # If inputMap is requested, create the initial version that can be 

667 # masked in assembleSubregion. 

668 if self.config.doInputMap: 

669 self.inputMapper.build_ccd_input_map( 

670 skyInfo.bbox, skyInfo.wcs, coaddExposure.getInfo().getCoaddInputs().ccds 

671 ) 

672 

673 if self.config.doOnlineForMean and self.config.statistic == "MEAN": 

674 try: 

675 self.assembleOnlineMeanCoadd( 

676 coaddExposure, 

677 tempExpRefList, 

678 imageScalerList, 

679 weightList, 

680 altMaskList, 

681 stats.ctrl, 

682 nImage=nImage, 

683 ) 

684 except Exception as e: 

685 self.log.exception("Cannot compute online coadd %s", e) 

686 raise 

687 else: 

688 for subBBox in subBBoxIter(skyInfo.bbox, subregionSize): 

689 try: 

690 self.assembleSubregion( 

691 coaddExposure, 

692 subBBox, 

693 tempExpRefList, 

694 imageScalerList, 

695 weightList, 

696 altMaskList, 

697 stats.flags, 

698 stats.ctrl, 

699 nImage=nImage, 

700 ) 

701 except Exception as e: 

702 self.log.exception("Cannot compute coadd %s: %s", subBBox, e) 

703 raise 

704 

705 # If inputMap is requested, we must finalize the map after the 

706 # accumulation. 

707 if self.config.doInputMap: 

708 self.inputMapper.finalize_ccd_input_map_mask() 

709 inputMap = self.inputMapper.ccd_input_map 

710 else: 

711 inputMap = None 

712 

713 self.setInexactPsf(coaddMaskedImage.getMask()) 

714 # Despite the name, the following doesn't really deal with "EDGE" 

715 # pixels: it identifies pixels that didn't receive any unmasked inputs 

716 # (as occurs around the edge of the field). 

717 coaddUtils.setCoaddEdgeBits(coaddMaskedImage.getMask(), coaddMaskedImage.getVariance()) 

718 return pipeBase.Struct( 

719 coaddExposure=coaddExposure, 

720 nImage=nImage, 

721 warpRefList=tempExpRefList, 

722 imageScalerList=imageScalerList, 

723 weightList=weightList, 

724 inputMap=inputMap, 

725 ) 

726 

727 def assembleMetadata(self, coaddExposure, tempExpRefList, weightList): 

728 """Set the metadata for the coadd. 

729 

730 This basic implementation sets the filter from the first input. 

731 

732 Parameters 

733 ---------- 

734 coaddExposure : `lsst.afw.image.Exposure` 

735 The target exposure for the coadd. 

736 tempExpRefList : `list` 

737 List of data references to tempExp. 

738 weightList : `list` 

739 List of weights. 

740 

741 Raises 

742 ------ 

743 AssertionError 

744 Raised if there is a length mismatch. 

745 """ 

746 assert len(tempExpRefList) == len(weightList), "Length mismatch" 

747 

748 # We load a single pixel of each coaddTempExp, because we just want to 

749 # get at the metadata (and we need more than just the PropertySet that 

750 # contains the header), which is not possible with the current butler 

751 # (see #2777). 

752 bbox = geom.Box2I(coaddExposure.getBBox().getMin(), geom.Extent2I(1, 1)) 

753 

754 tempExpList = [tempExpRef.get(parameters={"bbox": bbox}) for tempExpRef in tempExpRefList] 

755 

756 numCcds = sum(len(tempExp.getInfo().getCoaddInputs().ccds) for tempExp in tempExpList) 

757 

758 # Set the coadd FilterLabel to the band of the first input exposure: 

759 # Coadds are calibrated, so the physical label is now meaningless. 

760 coaddExposure.setFilter(afwImage.FilterLabel(tempExpList[0].getFilter().bandLabel)) 

761 coaddInputs = coaddExposure.getInfo().getCoaddInputs() 

762 coaddInputs.ccds.reserve(numCcds) 

763 coaddInputs.visits.reserve(len(tempExpList)) 

764 

765 for tempExp, weight in zip(tempExpList, weightList): 

766 self.inputRecorder.addVisitToCoadd(coaddInputs, tempExp, weight) 

767 

768 if self.config.doUsePsfMatchedPolygons: 

769 self.shrinkValidPolygons(coaddInputs) 

770 

771 coaddInputs.visits.sort() 

772 coaddInputs.ccds.sort() 

773 if self.warpType == "psfMatched": 

774 # The modelPsf BBox for a psfMatchedWarp/coaddTempExp was 

775 # dynamically defined by ModelPsfMatchTask as the square box 

776 # bounding its spatially-variable, pre-matched WarpedPsf. 

777 # Likewise, set the PSF of a PSF-Matched Coadd to the modelPsf 

778 # having the maximum width (sufficient because square) 

779 modelPsfList = [tempExp.getPsf() for tempExp in tempExpList] 

780 modelPsfWidthList = [ 

781 modelPsf.computeBBox(modelPsf.getAveragePosition()).getWidth() for modelPsf in modelPsfList 

782 ] 

783 psf = modelPsfList[modelPsfWidthList.index(max(modelPsfWidthList))] 

784 else: 

785 psf = measAlg.CoaddPsf( 

786 coaddInputs.ccds, coaddExposure.getWcs(), self.config.coaddPsf.makeControl() 

787 ) 

788 coaddExposure.setPsf(psf) 

789 apCorrMap = measAlg.makeCoaddApCorrMap( 

790 coaddInputs.ccds, coaddExposure.getBBox(afwImage.PARENT), coaddExposure.getWcs() 

791 ) 

792 coaddExposure.getInfo().setApCorrMap(apCorrMap) 

793 if self.config.doAttachTransmissionCurve: 

794 transmissionCurve = measAlg.makeCoaddTransmissionCurve(coaddExposure.getWcs(), coaddInputs.ccds) 

795 coaddExposure.getInfo().setTransmissionCurve(transmissionCurve) 

796 

797 def assembleSubregion( 

798 self, 

799 coaddExposure, 

800 bbox, 

801 tempExpRefList, 

802 imageScalerList, 

803 weightList, 

804 altMaskList, 

805 statsFlags, 

806 statsCtrl, 

807 nImage=None, 

808 ): 

809 """Assemble the coadd for a sub-region. 

810 

811 For each coaddTempExp, check for (and swap in) an alternative mask 

812 if one is passed. Remove mask planes listed in 

813 `config.removeMaskPlanes`. Finally, stack the actual exposures using 

814 `lsst.afw.math.statisticsStack` with the statistic specified by 

815 statsFlags. Typically, the statsFlag will be one of lsst.afw.math.MEAN 

816 for a mean-stack or `lsst.afw.math.MEANCLIP` for outlier rejection 

817 using an N-sigma clipped mean where N and iterations are specified by 

818 statsCtrl. Assign the stacked subregion back to the coadd. 

819 

820 Parameters 

821 ---------- 

822 coaddExposure : `lsst.afw.image.Exposure` 

823 The target exposure for the coadd. 

824 bbox : `lsst.geom.Box` 

825 Sub-region to coadd. 

826 tempExpRefList : `list` 

827 List of data reference to tempExp. 

828 imageScalerList : `list` 

829 List of image scalers. 

830 weightList : `list` 

831 List of weights. 

832 altMaskList : `list` 

833 List of alternate masks to use rather than those stored with 

834 tempExp, or None. Each element is dict with keys = mask plane 

835 name to which to add the spans. 

836 statsFlags : `lsst.afw.math.Property` 

837 Property object for statistic for coadd. 

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

839 Statistics control object for coadd. 

840 nImage : `lsst.afw.image.ImageU`, optional 

841 Keeps track of exposure count for each pixel. 

842 """ 

843 self.log.debug("Computing coadd over %s", bbox) 

844 

845 coaddExposure.mask.addMaskPlane("REJECTED") 

846 coaddExposure.mask.addMaskPlane("CLIPPED") 

847 coaddExposure.mask.addMaskPlane("SENSOR_EDGE") 

848 maskMap = self.setRejectedMaskMapping(statsCtrl) 

849 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED") 

850 maskedImageList = [] 

851 if nImage is not None: 

852 subNImage = afwImage.ImageU(bbox.getWidth(), bbox.getHeight()) 

853 for tempExpRef, imageScaler, altMask in zip(tempExpRefList, imageScalerList, altMaskList): 

854 exposure = tempExpRef.get(parameters={"bbox": bbox}) 

855 

856 maskedImage = exposure.getMaskedImage() 

857 mask = maskedImage.getMask() 

858 if altMask is not None: 

859 self.applyAltMaskPlanes(mask, altMask) 

860 imageScaler.scaleMaskedImage(maskedImage) 

861 

862 # Add 1 for each pixel which is not excluded by the exclude mask. 

863 # In legacyCoadd, pixels may also be excluded by 

864 # afwMath.statisticsStack. 

865 if nImage is not None: 

866 subNImage.getArray()[maskedImage.getMask().getArray() & statsCtrl.getAndMask() == 0] += 1 

867 if self.config.removeMaskPlanes: 

868 self.removeMaskPlanes(maskedImage) 

869 maskedImageList.append(maskedImage) 

870 

871 if self.config.doInputMap: 

872 visit = exposure.getInfo().getCoaddInputs().visits[0].getId() 

873 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask()) 

874 

875 with self.timer("stack"): 

876 coaddSubregion = afwMath.statisticsStack( 

877 maskedImageList, 

878 statsFlags, 

879 statsCtrl, 

880 weightList, 

881 clipped, # also set output to CLIPPED if sigma-clipped 

882 maskMap, 

883 ) 

884 coaddExposure.maskedImage.assign(coaddSubregion, bbox) 

885 if nImage is not None: 

886 nImage.assign(subNImage, bbox) 

887 

888 def assembleOnlineMeanCoadd( 

889 self, coaddExposure, tempExpRefList, imageScalerList, weightList, altMaskList, statsCtrl, nImage=None 

890 ): 

891 """Assemble the coadd using the "online" method. 

892 

893 This method takes a running sum of images and weights to save memory. 

894 It only works for MEAN statistics. 

895 

896 Parameters 

897 ---------- 

898 coaddExposure : `lsst.afw.image.Exposure` 

899 The target exposure for the coadd. 

900 tempExpRefList : `list` 

901 List of data reference to tempExp. 

902 imageScalerList : `list` 

903 List of image scalers. 

904 weightList : `list` 

905 List of weights. 

906 altMaskList : `list` 

907 List of alternate masks to use rather than those stored with 

908 tempExp, or None. Each element is dict with keys = mask plane 

909 name to which to add the spans. 

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

911 Statistics control object for coadd. 

912 nImage : `lsst.afw.image.ImageU`, optional 

913 Keeps track of exposure count for each pixel. 

914 """ 

915 self.log.debug("Computing online coadd.") 

916 

917 coaddExposure.mask.addMaskPlane("REJECTED") 

918 coaddExposure.mask.addMaskPlane("CLIPPED") 

919 coaddExposure.mask.addMaskPlane("SENSOR_EDGE") 

920 maskMap = self.setRejectedMaskMapping(statsCtrl) 

921 thresholdDict = AccumulatorMeanStack.stats_ctrl_to_threshold_dict(statsCtrl) 

922 

923 bbox = coaddExposure.maskedImage.getBBox() 

924 

925 stacker = AccumulatorMeanStack( 

926 coaddExposure.image.array.shape, 

927 statsCtrl.getAndMask(), 

928 mask_threshold_dict=thresholdDict, 

929 mask_map=maskMap, 

930 no_good_pixels_mask=statsCtrl.getNoGoodPixelsMask(), 

931 calc_error_from_input_variance=self.config.calcErrorFromInputVariance, 

932 compute_n_image=(nImage is not None), 

933 ) 

934 

935 for tempExpRef, imageScaler, altMask, weight in zip( 

936 tempExpRefList, imageScalerList, altMaskList, weightList 

937 ): 

938 exposure = tempExpRef.get() 

939 maskedImage = exposure.getMaskedImage() 

940 mask = maskedImage.getMask() 

941 if altMask is not None: 

942 self.applyAltMaskPlanes(mask, altMask) 

943 imageScaler.scaleMaskedImage(maskedImage) 

944 if self.config.removeMaskPlanes: 

945 self.removeMaskPlanes(maskedImage) 

946 

947 stacker.add_masked_image(maskedImage, weight=weight) 

948 

949 if self.config.doInputMap: 

950 visit = exposure.getInfo().getCoaddInputs().visits[0].getId() 

951 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask()) 

952 

953 stacker.fill_stacked_masked_image(coaddExposure.maskedImage) 

954 

955 if nImage is not None: 

956 nImage.array[:, :] = stacker.n_image 

957 

958 def removeMaskPlanes(self, maskedImage): 

959 """Unset the mask of an image for mask planes specified in the config. 

960 

961 Parameters 

962 ---------- 

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

964 The masked image to be modified. 

965 

966 Raises 

967 ------ 

968 InvalidParameterError 

969 Raised if no mask plane with that name was found. 

970 """ 

971 mask = maskedImage.getMask() 

972 for maskPlane in self.config.removeMaskPlanes: 

973 try: 

974 mask &= ~mask.getPlaneBitMask(maskPlane) 

975 except pexExceptions.InvalidParameterError: 

976 self.log.debug( 

977 "Unable to remove mask plane %s: no mask plane with that name was found.", maskPlane 

978 ) 

979 

980 @staticmethod 

981 def setRejectedMaskMapping(statsCtrl): 

982 """Map certain mask planes of the warps to new planes for the coadd. 

983 

984 If a pixel is rejected due to a mask value other than EDGE, NO_DATA, 

985 or CLIPPED, set it to REJECTED on the coadd. 

986 If a pixel is rejected due to EDGE, set the coadd pixel to SENSOR_EDGE. 

987 If a pixel is rejected due to CLIPPED, set the coadd pixel to CLIPPED. 

988 

989 Parameters 

990 ---------- 

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

992 Statistics control object for coadd. 

993 

994 Returns 

995 ------- 

996 maskMap : `list` of `tuple` of `int` 

997 A list of mappings of mask planes of the warped exposures to 

998 mask planes of the coadd. 

999 """ 

1000 edge = afwImage.Mask.getPlaneBitMask("EDGE") 

1001 noData = afwImage.Mask.getPlaneBitMask("NO_DATA") 

1002 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED") 

1003 toReject = statsCtrl.getAndMask() & (~noData) & (~edge) & (~clipped) 

1004 maskMap = [ 

1005 (toReject, afwImage.Mask.getPlaneBitMask("REJECTED")), 

1006 (edge, afwImage.Mask.getPlaneBitMask("SENSOR_EDGE")), 

1007 (clipped, clipped), 

1008 ] 

1009 return maskMap 

1010 

1011 def applyAltMaskPlanes(self, mask, altMaskSpans): 

1012 """Apply in place alt mask formatted as SpanSets to a mask. 

1013 

1014 Parameters 

1015 ---------- 

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

1017 Original mask. 

1018 altMaskSpans : `dict` 

1019 SpanSet lists to apply. Each element contains the new mask 

1020 plane name (e.g. "CLIPPED and/or "NO_DATA") as the key, 

1021 and list of SpanSets to apply to the mask. 

1022 

1023 Returns 

1024 ------- 

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

1026 Updated mask. 

1027 """ 

1028 if self.config.doUsePsfMatchedPolygons: 

1029 if ("NO_DATA" in altMaskSpans) and ("NO_DATA" in self.config.badMaskPlanes): 

1030 # Clear away any other masks outside the validPolygons. These 

1031 # pixels are no longer contributing to inexact PSFs, and will 

1032 # still be rejected because of NO_DATA. 

1033 # self.config.doUsePsfMatchedPolygons should be True only in 

1034 # CompareWarpAssemble. This mask-clearing step must only occur 

1035 # *before* applying the new masks below. 

1036 for spanSet in altMaskSpans["NO_DATA"]: 

1037 spanSet.clippedTo(mask.getBBox()).clearMask(mask, self.getBadPixelMask()) 

1038 

1039 for plane, spanSetList in altMaskSpans.items(): 

1040 maskClipValue = mask.addMaskPlane(plane) 

1041 for spanSet in spanSetList: 

1042 spanSet.clippedTo(mask.getBBox()).setMask(mask, 2**maskClipValue) 

1043 return mask 

1044 

1045 def shrinkValidPolygons(self, coaddInputs): 

1046 """Shrink coaddInputs' ccds' ValidPolygons in place. 

1047 

1048 Either modify each ccd's validPolygon in place, or if CoaddInputs 

1049 does not have a validPolygon, create one from its bbox. 

1050 

1051 Parameters 

1052 ---------- 

1053 coaddInputs : `lsst.afw.image.coaddInputs` 

1054 Original mask. 

1055 """ 

1056 for ccd in coaddInputs.ccds: 

1057 polyOrig = ccd.getValidPolygon() 

1058 validPolyBBox = polyOrig.getBBox() if polyOrig else ccd.getBBox() 

1059 validPolyBBox.grow(-self.config.matchingKernelSize // 2) 

1060 if polyOrig: 

1061 validPolygon = polyOrig.intersectionSingle(validPolyBBox) 

1062 else: 

1063 validPolygon = afwGeom.polygon.Polygon(geom.Box2D(validPolyBBox)) 

1064 ccd.setValidPolygon(validPolygon) 

1065 

1066 def setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None): 

1067 """Set the bright object masks. 

1068 

1069 Parameters 

1070 ---------- 

1071 exposure : `lsst.afw.image.Exposure` 

1072 Exposure under consideration. 

1073 brightObjectMasks : `lsst.afw.table` 

1074 Table of bright objects to mask. 

1075 dataId : `lsst.daf.butler.DataId`, optional 

1076 Data identifier dict for patch. 

1077 """ 

1078 if brightObjectMasks is None: 

1079 self.log.warning("Unable to apply bright object mask: none supplied") 

1080 return 

1081 self.log.info("Applying %d bright object masks to %s", len(brightObjectMasks), dataId) 

1082 mask = exposure.getMaskedImage().getMask() 

1083 wcs = exposure.getWcs() 

1084 plateScale = wcs.getPixelScale().asArcseconds() 

1085 

1086 for rec in brightObjectMasks: 

1087 center = geom.PointI(wcs.skyToPixel(rec.getCoord())) 

1088 if rec["type"] == "box": 

1089 assert rec["angle"] == 0.0, "Angle != 0 for mask object %s" % rec["id"] 

1090 width = rec["width"].asArcseconds() / plateScale # convert to pixels 

1091 height = rec["height"].asArcseconds() / plateScale # convert to pixels 

1092 

1093 halfSize = geom.ExtentI(0.5 * width, 0.5 * height) 

1094 bbox = geom.Box2I(center - halfSize, center + halfSize) 

1095 

1096 bbox = geom.BoxI( 

1097 geom.PointI(int(center[0] - 0.5 * width), int(center[1] - 0.5 * height)), 

1098 geom.PointI(int(center[0] + 0.5 * width), int(center[1] + 0.5 * height)), 

1099 ) 

1100 spans = afwGeom.SpanSet(bbox) 

1101 elif rec["type"] == "circle": 

1102 radius = int(rec["radius"].asArcseconds() / plateScale) # convert to pixels 

1103 spans = afwGeom.SpanSet.fromShape(radius, offset=center) 

1104 else: 

1105 self.log.warning("Unexpected region type %s at %s", rec["type"], center) 

1106 continue 

1107 spans.clippedTo(mask.getBBox()).setMask(mask, self.brightObjectBitmask) 

1108 

1109 def setInexactPsf(self, mask): 

1110 """Set INEXACT_PSF mask plane. 

1111 

1112 If any of the input images isn't represented in the coadd (due to 

1113 clipped pixels or chip gaps), the `CoaddPsf` will be inexact. Flag 

1114 these pixels. 

1115 

1116 Parameters 

1117 ---------- 

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

1119 Coadded exposure's mask, modified in-place. 

1120 """ 

1121 mask.addMaskPlane("INEXACT_PSF") 

1122 inexactPsf = mask.getPlaneBitMask("INEXACT_PSF") 

1123 sensorEdge = mask.getPlaneBitMask("SENSOR_EDGE") # chip edges (so PSF is discontinuous) 

1124 clipped = mask.getPlaneBitMask("CLIPPED") # pixels clipped from coadd 

1125 rejected = mask.getPlaneBitMask("REJECTED") # pixels rejected from coadd due to masks 

1126 array = mask.getArray() 

1127 selected = array & (sensorEdge | clipped | rejected) > 0 

1128 array[selected] |= inexactPsf 

1129 

1130 def filterWarps(self, inputs, goodVisits): 

1131 """Return list of only inputRefs with visitId in goodVisits ordered by 

1132 goodVisit. 

1133 

1134 Parameters 

1135 ---------- 

1136 inputs : `list` of `~lsst.pipe.base.connections.DeferredDatasetRef` 

1137 List of `lsst.pipe.base.connections.DeferredDatasetRef` with dataId 

1138 containing visit. 

1139 goodVisit : `dict` 

1140 Dictionary with good visitIds as the keys. Value ignored. 

1141 

1142 Returns 

1143 ------- 

1144 filterInputs : `list` [`lsst.pipe.base.connections.DeferredDatasetRef`] 

1145 Filtered and sorted list of inputRefs with visitId in goodVisits 

1146 ordered by goodVisit. 

1147 """ 

1148 inputWarpDict = {inputRef.ref.dataId["visit"]: inputRef for inputRef in inputs} 

1149 filteredInputs = [] 

1150 for visit in goodVisits.keys(): 

1151 if visit in inputWarpDict: 

1152 filteredInputs.append(inputWarpDict[visit]) 

1153 return filteredInputs 

1154 

1155 

1156def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask): 

1157 """Function to count the number of pixels with a specific mask in a 

1158 footprint. 

1159 

1160 Find the intersection of mask & footprint. Count all pixels in the mask 

1161 that are in the intersection that have bitmask set but do not have 

1162 ignoreMask set. Return the count. 

1163 

1164 Parameters 

1165 ---------- 

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

1167 Mask to define intersection region by. 

1168 footprint : `lsst.afw.detection.Footprint` 

1169 Footprint to define the intersection region by. 

1170 bitmask : `Unknown` 

1171 Specific mask that we wish to count the number of occurances of. 

1172 ignoreMask : `Unknown` 

1173 Pixels to not consider. 

1174 

1175 Returns 

1176 ------- 

1177 result : `int` 

1178 Number of pixels in footprint with specified mask. 

1179 """ 

1180 bbox = footprint.getBBox() 

1181 bbox.clip(mask.getBBox(afwImage.PARENT)) 

1182 fp = afwImage.Mask(bbox) 

1183 subMask = mask.Factory(mask, bbox, afwImage.PARENT) 

1184 footprint.spans.setMask(fp, bitmask) 

1185 return numpy.logical_and( 

1186 (subMask.getArray() & fp.getArray()) > 0, (subMask.getArray() & ignoreMask) == 0 

1187 ).sum() 

1188 

1189 

1190class CompareWarpAssembleCoaddConnections(AssembleCoaddConnections): 

1191 psfMatchedWarps = pipeBase.connectionTypes.Input( 

1192 doc=( 

1193 "PSF-Matched Warps are required by CompareWarp regardless of the coadd type requested. " 

1194 "Only PSF-Matched Warps make sense for image subtraction. " 

1195 "Therefore, they must be an additional declared input." 

1196 ), 

1197 name="{inputCoaddName}Coadd_psfMatchedWarp", 

1198 storageClass="ExposureF", 

1199 dimensions=("tract", "patch", "skymap", "visit"), 

1200 deferLoad=True, 

1201 multiple=True, 

1202 ) 

1203 templateCoadd = pipeBase.connectionTypes.Output( 

1204 doc=( 

1205 "Model of the static sky, used to find temporal artifacts. Typically a PSF-Matched, " 

1206 "sigma-clipped coadd. Written if and only if assembleStaticSkyModel.doWrite=True" 

1207 ), 

1208 name="{outputCoaddName}CoaddPsfMatched", 

1209 storageClass="ExposureF", 

1210 dimensions=("tract", "patch", "skymap", "band"), 

1211 ) 

1212 

1213 def __init__(self, *, config=None): 

1214 super().__init__(config=config) 

1215 if not config.assembleStaticSkyModel.doWrite: 

1216 self.outputs.remove("templateCoadd") 

1217 config.validate() 

1218 

1219 

1220class CompareWarpAssembleCoaddConfig( 

1221 AssembleCoaddConfig, pipelineConnections=CompareWarpAssembleCoaddConnections 

1222): 

1223 assembleStaticSkyModel = pexConfig.ConfigurableField( 

1224 target=AssembleCoaddTask, 

1225 doc="Task to assemble an artifact-free, PSF-matched Coadd to serve as " 

1226 "a naive/first-iteration model of the static sky.", 

1227 ) 

1228 detect = pexConfig.ConfigurableField( 

1229 target=SourceDetectionTask, 

1230 doc="Detect outlier sources on difference between each psfMatched warp and static sky model", 

1231 ) 

1232 detectTemplate = pexConfig.ConfigurableField( 

1233 target=SourceDetectionTask, 

1234 doc="Detect sources on static sky model. Only used if doPreserveContainedBySource is True", 

1235 ) 

1236 maskStreaks = pexConfig.ConfigurableField( 

1237 target=MaskStreaksTask, 

1238 doc="Detect streaks on difference between each psfMatched warp and static sky model. Only used if " 

1239 "doFilterMorphological is True. Adds a mask plane to an exposure, with the mask plane name set by" 

1240 "streakMaskName", 

1241 ) 

1242 streakMaskName = pexConfig.Field(dtype=str, default="STREAK", doc="Name of mask bit used for streaks") 

1243 maxNumEpochs = pexConfig.Field( 

1244 doc="Charactistic maximum local number of epochs/visits in which an artifact candidate can appear " 

1245 "and still be masked. The effective maxNumEpochs is a broken linear function of local " 

1246 "number of epochs (N): min(maxFractionEpochsLow*N, maxNumEpochs + maxFractionEpochsHigh*N). " 

1247 "For each footprint detected on the image difference between the psfMatched warp and static sky " 

1248 "model, if a significant fraction of pixels (defined by spatialThreshold) are residuals in more " 

1249 "than the computed effective maxNumEpochs, the artifact candidate is deemed persistant rather " 

1250 "than transient and not masked.", 

1251 dtype=int, 

1252 default=2, 

1253 ) 

1254 maxFractionEpochsLow = pexConfig.RangeField( 

1255 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for low N. " 

1256 "Effective maxNumEpochs = " 

1257 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)", 

1258 dtype=float, 

1259 default=0.4, 

1260 min=0.0, 

1261 max=1.0, 

1262 ) 

1263 maxFractionEpochsHigh = pexConfig.RangeField( 

1264 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for high N. " 

1265 "Effective maxNumEpochs = " 

1266 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)", 

1267 dtype=float, 

1268 default=0.03, 

1269 min=0.0, 

1270 max=1.0, 

1271 ) 

1272 spatialThreshold = pexConfig.RangeField( 

1273 doc="Unitless fraction of pixels defining how much of the outlier region has to meet the " 

1274 "temporal criteria. If 0, clip all. If 1, clip none.", 

1275 dtype=float, 

1276 default=0.5, 

1277 min=0.0, 

1278 max=1.0, 

1279 inclusiveMin=True, 

1280 inclusiveMax=True, 

1281 ) 

1282 doScaleWarpVariance = pexConfig.Field( 

1283 doc="Rescale Warp variance plane using empirical noise?", 

1284 dtype=bool, 

1285 default=True, 

1286 ) 

1287 scaleWarpVariance = pexConfig.ConfigurableField( 

1288 target=ScaleVarianceTask, 

1289 doc="Rescale variance on warps", 

1290 ) 

1291 doPreserveContainedBySource = pexConfig.Field( 

1292 doc="Rescue artifacts from clipping that completely lie within a footprint detected" 

1293 "on the PsfMatched Template Coadd. Replicates a behavior of SafeClip.", 

1294 dtype=bool, 

1295 default=True, 

1296 ) 

1297 doPrefilterArtifacts = pexConfig.Field( 

1298 doc="Ignore artifact candidates that are mostly covered by the bad pixel mask, " 

1299 "because they will be excluded anyway. This prevents them from contributing " 

1300 "to the outlier epoch count image and potentially being labeled as persistant." 

1301 "'Mostly' is defined by the config 'prefilterArtifactsRatio'.", 

1302 dtype=bool, 

1303 default=True, 

1304 ) 

1305 prefilterArtifactsMaskPlanes = pexConfig.ListField( 

1306 doc="Prefilter artifact candidates that are mostly covered by these bad mask planes.", 

1307 dtype=str, 

1308 default=("NO_DATA", "BAD", "SAT", "SUSPECT"), 

1309 ) 

1310 prefilterArtifactsRatio = pexConfig.Field( 

1311 doc="Prefilter artifact candidates with less than this fraction overlapping good pixels", 

1312 dtype=float, 

1313 default=0.05, 

1314 ) 

1315 doFilterMorphological = pexConfig.Field( 

1316 doc="Filter artifact candidates based on morphological criteria, i.g. those that appear to " 

1317 "be streaks.", 

1318 dtype=bool, 

1319 default=False, 

1320 ) 

1321 growStreakFp = pexConfig.Field( 

1322 doc="Grow streak footprints by this number multiplied by the PSF width", dtype=float, default=5 

1323 ) 

1324 

1325 def setDefaults(self): 

1326 AssembleCoaddConfig.setDefaults(self) 

1327 self.statistic = "MEAN" 

1328 self.doUsePsfMatchedPolygons = True 

1329 

1330 # Real EDGE removed by psfMatched NO_DATA border half the width of the 

1331 # matching kernel. CompareWarp applies psfMatched EDGE pixels to 

1332 # directWarps before assembling. 

1333 if "EDGE" in self.badMaskPlanes: 

1334 self.badMaskPlanes.remove("EDGE") 

1335 self.removeMaskPlanes.append("EDGE") 

1336 self.assembleStaticSkyModel.badMaskPlanes = [ 

1337 "NO_DATA", 

1338 ] 

1339 self.assembleStaticSkyModel.warpType = "psfMatched" 

1340 self.assembleStaticSkyModel.connections.warpType = "psfMatched" 

1341 self.assembleStaticSkyModel.statistic = "MEANCLIP" 

1342 self.assembleStaticSkyModel.sigmaClip = 2.5 

1343 self.assembleStaticSkyModel.clipIter = 3 

1344 self.assembleStaticSkyModel.calcErrorFromInputVariance = False 

1345 self.assembleStaticSkyModel.doWrite = False 

1346 self.detect.doTempLocalBackground = False 

1347 self.detect.reEstimateBackground = False 

1348 self.detect.returnOriginalFootprints = False 

1349 self.detect.thresholdPolarity = "both" 

1350 self.detect.thresholdValue = 5 

1351 self.detect.minPixels = 4 

1352 self.detect.isotropicGrow = True 

1353 self.detect.thresholdType = "pixel_stdev" 

1354 self.detect.nSigmaToGrow = 0.4 

1355 # The default nSigmaToGrow for SourceDetectionTask is already 2.4, 

1356 # Explicitly restating because ratio with detect.nSigmaToGrow matters 

1357 self.detectTemplate.nSigmaToGrow = 2.4 

1358 self.detectTemplate.doTempLocalBackground = False 

1359 self.detectTemplate.reEstimateBackground = False 

1360 self.detectTemplate.returnOriginalFootprints = False 

1361 

1362 def validate(self): 

1363 super().validate() 

1364 if self.assembleStaticSkyModel.doNImage: 

1365 raise ValueError( 

1366 "No dataset type exists for a PSF-Matched Template N Image." 

1367 "Please set assembleStaticSkyModel.doNImage=False" 

1368 ) 

1369 

1370 if self.assembleStaticSkyModel.doWrite and (self.warpType == self.assembleStaticSkyModel.warpType): 

1371 raise ValueError( 

1372 "warpType (%s) == assembleStaticSkyModel.warpType (%s) and will compete for " 

1373 "the same dataset name. Please set assembleStaticSkyModel.doWrite to False " 

1374 "or warpType to 'direct'. assembleStaticSkyModel.warpType should ways be " 

1375 "'PsfMatched'" % (self.warpType, self.assembleStaticSkyModel.warpType) 

1376 ) 

1377 

1378 

1379class CompareWarpAssembleCoaddTask(AssembleCoaddTask): 

1380 """Assemble a compareWarp coadded image from a set of warps 

1381 by masking artifacts detected by comparing PSF-matched warps. 

1382 

1383 In ``AssembleCoaddTask``, we compute the coadd as an clipped mean (i.e., 

1384 we clip outliers). The problem with doing this is that when computing the 

1385 coadd PSF at a given location, individual visit PSFs from visits with 

1386 outlier pixels contribute to the coadd PSF and cannot be treated correctly. 

1387 In this task, we correct for this behavior by creating a new badMaskPlane 

1388 'CLIPPED' which marks pixels in the individual warps suspected to contain 

1389 an artifact. We populate this plane on the input warps by comparing 

1390 PSF-matched warps with a PSF-matched median coadd which serves as a 

1391 model of the static sky. Any group of pixels that deviates from the 

1392 PSF-matched template coadd by more than config.detect.threshold sigma, 

1393 is an artifact candidate. The candidates are then filtered to remove 

1394 variable sources and sources that are difficult to subtract such as 

1395 bright stars. This filter is configured using the config parameters 

1396 ``temporalThreshold`` and ``spatialThreshold``. The temporalThreshold is 

1397 the maximum fraction of epochs that the deviation can appear in and still 

1398 be considered an artifact. The spatialThreshold is the maximum fraction of 

1399 pixels in the footprint of the deviation that appear in other epochs 

1400 (where other epochs is defined by the temporalThreshold). If the deviant 

1401 region meets this criteria of having a significant percentage of pixels 

1402 that deviate in only a few epochs, these pixels have the 'CLIPPED' bit 

1403 set in the mask. These regions will not contribute to the final coadd. 

1404 Furthermore, any routine to determine the coadd PSF can now be cognizant 

1405 of clipped regions. Note that the algorithm implemented by this task is 

1406 preliminary and works correctly for HSC data. Parameter modifications and 

1407 or considerable redesigning of the algorithm is likley required for other 

1408 surveys. 

1409 

1410 ``CompareWarpAssembleCoaddTask`` sub-classes 

1411 ``AssembleCoaddTask`` and instantiates ``AssembleCoaddTask`` 

1412 as a subtask to generate the TemplateCoadd (the model of the static sky). 

1413 

1414 Notes 

1415 ----- 

1416 Debugging: 

1417 This task supports the following debug variables: 

1418 - ``saveCountIm`` 

1419 If True then save the Epoch Count Image as a fits file in the `figPath` 

1420 - ``figPath`` 

1421 Path to save the debug fits images and figures 

1422 """ 

1423 

1424 ConfigClass = CompareWarpAssembleCoaddConfig 

1425 _DefaultName = "compareWarpAssembleCoadd" 

1426 

1427 def __init__(self, *args, **kwargs): 

1428 AssembleCoaddTask.__init__(self, *args, **kwargs) 

1429 self.makeSubtask("assembleStaticSkyModel") 

1430 detectionSchema = afwTable.SourceTable.makeMinimalSchema() 

1431 self.makeSubtask("detect", schema=detectionSchema) 

1432 if self.config.doPreserveContainedBySource: 

1433 self.makeSubtask("detectTemplate", schema=afwTable.SourceTable.makeMinimalSchema()) 

1434 if self.config.doScaleWarpVariance: 

1435 self.makeSubtask("scaleWarpVariance") 

1436 if self.config.doFilterMorphological: 

1437 self.makeSubtask("maskStreaks") 

1438 

1439 @utils.inheritDoc(AssembleCoaddTask) 

1440 def _makeSupplementaryData(self, butlerQC, inputRefs, outputRefs): 

1441 """Generate a templateCoadd to use as a naive model of static sky to 

1442 subtract from PSF-Matched warps. 

1443 

1444 Returns 

1445 ------- 

1446 result : `~lsst.pipe.base.Struct` 

1447 Results as a struct with attributes: 

1448 

1449 ``templateCoadd`` 

1450 Coadded exposure (`lsst.afw.image.Exposure`). 

1451 ``nImage`` 

1452 Keeps track of exposure count for each pixel 

1453 (`lsst.afw.image.ImageU`). 

1454 

1455 Raises 

1456 ------ 

1457 RuntimeError 

1458 Raised if ``templateCoadd`` is `None`. 

1459 """ 

1460 # Ensure that psfMatchedWarps are used as input warps for template 

1461 # generation. 

1462 staticSkyModelInputRefs = copy.deepcopy(inputRefs) 

1463 staticSkyModelInputRefs.inputWarps = inputRefs.psfMatchedWarps 

1464 

1465 # Because subtasks don't have connections we have to make one. 

1466 # The main task's `templateCoadd` is the subtask's `coaddExposure` 

1467 staticSkyModelOutputRefs = copy.deepcopy(outputRefs) 

1468 if self.config.assembleStaticSkyModel.doWrite: 

1469 staticSkyModelOutputRefs.coaddExposure = staticSkyModelOutputRefs.templateCoadd 

1470 # Remove template coadd from both subtask's and main tasks outputs, 

1471 # because it is handled by the subtask as `coaddExposure` 

1472 del outputRefs.templateCoadd 

1473 del staticSkyModelOutputRefs.templateCoadd 

1474 

1475 # A PSF-Matched nImage does not exist as a dataset type 

1476 if "nImage" in staticSkyModelOutputRefs.keys(): 

1477 del staticSkyModelOutputRefs.nImage 

1478 

1479 templateCoadd = self.assembleStaticSkyModel.runQuantum( 

1480 butlerQC, staticSkyModelInputRefs, staticSkyModelOutputRefs 

1481 ) 

1482 if templateCoadd is None: 

1483 raise RuntimeError(self._noTemplateMessage(self.assembleStaticSkyModel.warpType)) 

1484 

1485 return pipeBase.Struct( 

1486 templateCoadd=templateCoadd.coaddExposure, 

1487 nImage=templateCoadd.nImage, 

1488 warpRefList=templateCoadd.warpRefList, 

1489 imageScalerList=templateCoadd.imageScalerList, 

1490 weightList=templateCoadd.weightList, 

1491 ) 

1492 

1493 def _noTemplateMessage(self, warpType): 

1494 warpName = warpType[0].upper() + warpType[1:] 

1495 message = """No %(warpName)s warps were found to build the template coadd which is 

1496 required to run CompareWarpAssembleCoaddTask. To continue assembling this type of coadd, 

1497 first either rerun makeCoaddTempExp with config.make%(warpName)s=True or 

1498 coaddDriver with config.makeCoadTempExp.make%(warpName)s=True, before assembleCoadd. 

1499 

1500 Alternatively, to use another algorithm with existing warps, retarget the CoaddDriverConfig to 

1501 another algorithm like: 

1502 

1503 from lsst.pipe.tasks.assembleCoadd import SafeClipAssembleCoaddTask 

1504 config.assemble.retarget(SafeClipAssembleCoaddTask) 

1505 """ % { 

1506 "warpName": warpName 

1507 } 

1508 return message 

1509 

1510 @utils.inheritDoc(AssembleCoaddTask) 

1511 @timeMethod 

1512 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, supplementaryData): 

1513 """Notes 

1514 ----- 

1515 Assemble the coadd. 

1516 

1517 Find artifacts and apply them to the warps' masks creating a list of 

1518 alternative masks with a new "CLIPPED" plane and updated "NO_DATA" 

1519 plane. Then pass these alternative masks to the base class's ``run`` 

1520 method. 

1521 """ 

1522 # Check and match the order of the supplementaryData 

1523 # (PSF-matched) inputs to the order of the direct inputs, 

1524 # so that the artifact mask is applied to the right warp 

1525 dataIds = [ref.dataId for ref in tempExpRefList] 

1526 psfMatchedDataIds = [ref.dataId for ref in supplementaryData.warpRefList] 

1527 

1528 if dataIds != psfMatchedDataIds: 

1529 self.log.info("Reordering and or/padding PSF-matched visit input list") 

1530 supplementaryData.warpRefList = reorderAndPadList( 

1531 supplementaryData.warpRefList, psfMatchedDataIds, dataIds 

1532 ) 

1533 supplementaryData.imageScalerList = reorderAndPadList( 

1534 supplementaryData.imageScalerList, psfMatchedDataIds, dataIds 

1535 ) 

1536 

1537 # Use PSF-Matched Warps (and corresponding scalers) and coadd to find 

1538 # artifacts. 

1539 spanSetMaskList = self.findArtifacts( 

1540 supplementaryData.templateCoadd, supplementaryData.warpRefList, supplementaryData.imageScalerList 

1541 ) 

1542 

1543 badMaskPlanes = self.config.badMaskPlanes[:] 

1544 badMaskPlanes.append("CLIPPED") 

1545 badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes) 

1546 

1547 result = AssembleCoaddTask.run( 

1548 self, skyInfo, tempExpRefList, imageScalerList, weightList, spanSetMaskList, mask=badPixelMask 

1549 ) 

1550 

1551 # Propagate PSF-matched EDGE pixels to coadd SENSOR_EDGE and 

1552 # INEXACT_PSF. Psf-Matching moves the real edge inwards. 

1553 self.applyAltEdgeMask(result.coaddExposure.maskedImage.mask, spanSetMaskList) 

1554 return result 

1555 

1556 def applyAltEdgeMask(self, mask, altMaskList): 

1557 """Propagate alt EDGE mask to SENSOR_EDGE AND INEXACT_PSF planes. 

1558 

1559 Parameters 

1560 ---------- 

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

1562 Original mask. 

1563 altMaskList : `list` of `dict` 

1564 List of Dicts containing ``spanSet`` lists. 

1565 Each element contains the new mask plane name (e.g. "CLIPPED 

1566 and/or "NO_DATA") as the key, and list of ``SpanSets`` to apply to 

1567 the mask. 

1568 """ 

1569 maskValue = mask.getPlaneBitMask(["SENSOR_EDGE", "INEXACT_PSF"]) 

1570 for visitMask in altMaskList: 

1571 if "EDGE" in visitMask: 

1572 for spanSet in visitMask["EDGE"]: 

1573 spanSet.clippedTo(mask.getBBox()).setMask(mask, maskValue) 

1574 

1575 def findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList): 

1576 """Find artifacts. 

1577 

1578 Loop through warps twice. The first loop builds a map with the count 

1579 of how many epochs each pixel deviates from the templateCoadd by more 

1580 than ``config.chiThreshold`` sigma. The second loop takes each 

1581 difference image and filters the artifacts detected in each using 

1582 count map to filter out variable sources and sources that are 

1583 difficult to subtract cleanly. 

1584 

1585 Parameters 

1586 ---------- 

1587 templateCoadd : `lsst.afw.image.Exposure` 

1588 Exposure to serve as model of static sky. 

1589 tempExpRefList : `list` 

1590 List of data references to warps. 

1591 imageScalerList : `list` 

1592 List of image scalers. 

1593 

1594 Returns 

1595 ------- 

1596 altMasks : `list` of `dict` 

1597 List of dicts containing information about CLIPPED 

1598 (i.e., artifacts), NO_DATA, and EDGE pixels. 

1599 """ 

1600 self.log.debug("Generating Count Image, and mask lists.") 

1601 coaddBBox = templateCoadd.getBBox() 

1602 slateIm = afwImage.ImageU(coaddBBox) 

1603 epochCountImage = afwImage.ImageU(coaddBBox) 

1604 nImage = afwImage.ImageU(coaddBBox) 

1605 spanSetArtifactList = [] 

1606 spanSetNoDataMaskList = [] 

1607 spanSetEdgeList = [] 

1608 spanSetBadMorphoList = [] 

1609 badPixelMask = self.getBadPixelMask() 

1610 

1611 # mask of the warp diffs should = that of only the warp 

1612 templateCoadd.mask.clearAllMaskPlanes() 

1613 

1614 if self.config.doPreserveContainedBySource: 

1615 templateFootprints = self.detectTemplate.detectFootprints(templateCoadd) 

1616 else: 

1617 templateFootprints = None 

1618 

1619 for warpRef, imageScaler in zip(tempExpRefList, imageScalerList): 

1620 warpDiffExp = self._readAndComputeWarpDiff(warpRef, imageScaler, templateCoadd) 

1621 if warpDiffExp is not None: 

1622 # This nImage only approximates the final nImage because it 

1623 # uses the PSF-matched mask. 

1624 nImage.array += ( 

1625 numpy.isfinite(warpDiffExp.image.array) * ((warpDiffExp.mask.array & badPixelMask) == 0) 

1626 ).astype(numpy.uint16) 

1627 fpSet = self.detect.detectFootprints(warpDiffExp, doSmooth=False, clearMask=True) 

1628 fpSet.positive.merge(fpSet.negative) 

1629 footprints = fpSet.positive 

1630 slateIm.set(0) 

1631 spanSetList = [footprint.spans for footprint in footprints.getFootprints()] 

1632 

1633 # Remove artifacts due to defects before they contribute to 

1634 # the epochCountImage. 

1635 if self.config.doPrefilterArtifacts: 

1636 spanSetList = self.prefilterArtifacts(spanSetList, warpDiffExp) 

1637 

1638 # Clear mask before adding prefiltered spanSets 

1639 self.detect.clearMask(warpDiffExp.mask) 

1640 for spans in spanSetList: 

1641 spans.setImage(slateIm, 1, doClip=True) 

1642 spans.setMask(warpDiffExp.mask, warpDiffExp.mask.getPlaneBitMask("DETECTED")) 

1643 epochCountImage += slateIm 

1644 

1645 if self.config.doFilterMorphological: 

1646 maskName = self.config.streakMaskName 

1647 # clear single frame streak mask if it exists already 

1648 if maskName in warpDiffExp.mask.getMaskPlaneDict(): 

1649 warpDiffExp.mask.clearMaskPlane(warpDiffExp.mask.getMaskPlane(maskName)) 

1650 else: 

1651 self.log.debug(f"Did not (need to) clear {maskName} mask because it didn't exist") 

1652 

1653 _ = self.maskStreaks.run(warpDiffExp) 

1654 streakMask = warpDiffExp.mask 

1655 spanSetStreak = afwGeom.SpanSet.fromMask( 

1656 streakMask, streakMask.getPlaneBitMask(maskName) 

1657 ).split() 

1658 # Pad the streaks to account for low-surface brightness 

1659 # wings. 

1660 psf = warpDiffExp.getPsf() 

1661 for s, sset in enumerate(spanSetStreak): 

1662 psfShape = psf.computeShape(sset.computeCentroid()) 

1663 dilation = self.config.growStreakFp * psfShape.getDeterminantRadius() 

1664 sset_dilated = sset.dilated(int(dilation)) 

1665 spanSetStreak[s] = sset_dilated 

1666 

1667 # PSF-Matched warps have less available area (~the matching 

1668 # kernel) because the calexps undergo a second convolution. 

1669 # Pixels with data in the direct warp but not in the 

1670 # PSF-matched warp will not have their artifacts detected. 

1671 # NaNs from the PSF-matched warp therefore must be masked in 

1672 # the direct warp. 

1673 nans = numpy.where(numpy.isnan(warpDiffExp.maskedImage.image.array), 1, 0) 

1674 nansMask = afwImage.makeMaskFromArray(nans.astype(afwImage.MaskPixel)) 

1675 nansMask.setXY0(warpDiffExp.getXY0()) 

1676 edgeMask = warpDiffExp.mask 

1677 spanSetEdgeMask = afwGeom.SpanSet.fromMask(edgeMask, edgeMask.getPlaneBitMask("EDGE")).split() 

1678 else: 

1679 # If the directWarp has <1% coverage, the psfMatchedWarp can 

1680 # have 0% and not exist. In this case, mask the whole epoch. 

1681 nansMask = afwImage.MaskX(coaddBBox, 1) 

1682 spanSetList = [] 

1683 spanSetEdgeMask = [] 

1684 spanSetStreak = [] 

1685 

1686 spanSetNoDataMask = afwGeom.SpanSet.fromMask(nansMask).split() 

1687 

1688 spanSetNoDataMaskList.append(spanSetNoDataMask) 

1689 spanSetArtifactList.append(spanSetList) 

1690 spanSetEdgeList.append(spanSetEdgeMask) 

1691 if self.config.doFilterMorphological: 

1692 spanSetBadMorphoList.append(spanSetStreak) 

1693 

1694 if lsstDebug.Info(__name__).saveCountIm: 

1695 path = self._dataRef2DebugPath("epochCountIm", tempExpRefList[0], coaddLevel=True) 

1696 epochCountImage.writeFits(path) 

1697 

1698 for i, spanSetList in enumerate(spanSetArtifactList): 

1699 if spanSetList: 

1700 filteredSpanSetList = self.filterArtifacts( 

1701 spanSetList, epochCountImage, nImage, templateFootprints 

1702 ) 

1703 spanSetArtifactList[i] = filteredSpanSetList 

1704 if self.config.doFilterMorphological: 

1705 spanSetArtifactList[i] += spanSetBadMorphoList[i] 

1706 

1707 altMasks = [] 

1708 for artifacts, noData, edge in zip(spanSetArtifactList, spanSetNoDataMaskList, spanSetEdgeList): 

1709 altMasks.append({"CLIPPED": artifacts, "NO_DATA": noData, "EDGE": edge}) 

1710 return altMasks 

1711 

1712 def prefilterArtifacts(self, spanSetList, exp): 

1713 """Remove artifact candidates covered by bad mask plane. 

1714 

1715 Any future editing of the candidate list that does not depend on 

1716 temporal information should go in this method. 

1717 

1718 Parameters 

1719 ---------- 

1720 spanSetList : `list` [`lsst.afw.geom.SpanSet`] 

1721 List of SpanSets representing artifact candidates. 

1722 exp : `lsst.afw.image.Exposure` 

1723 Exposure containing mask planes used to prefilter. 

1724 

1725 Returns 

1726 ------- 

1727 returnSpanSetList : `list` [`lsst.afw.geom.SpanSet`] 

1728 List of SpanSets with artifacts. 

1729 """ 

1730 badPixelMask = exp.mask.getPlaneBitMask(self.config.prefilterArtifactsMaskPlanes) 

1731 goodArr = (exp.mask.array & badPixelMask) == 0 

1732 returnSpanSetList = [] 

1733 bbox = exp.getBBox() 

1734 x0, y0 = exp.getXY0() 

1735 for i, span in enumerate(spanSetList): 

1736 y, x = span.clippedTo(bbox).indices() 

1737 yIndexLocal = numpy.array(y) - y0 

1738 xIndexLocal = numpy.array(x) - x0 

1739 goodRatio = numpy.count_nonzero(goodArr[yIndexLocal, xIndexLocal]) / span.getArea() 

1740 if goodRatio > self.config.prefilterArtifactsRatio: 

1741 returnSpanSetList.append(span) 

1742 return returnSpanSetList 

1743 

1744 def filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None): 

1745 """Filter artifact candidates. 

1746 

1747 Parameters 

1748 ---------- 

1749 spanSetList : `list` [`lsst.afw.geom.SpanSet`] 

1750 List of SpanSets representing artifact candidates. 

1751 epochCountImage : `lsst.afw.image.Image` 

1752 Image of accumulated number of warpDiff detections. 

1753 nImage : `lsst.afw.image.ImageU` 

1754 Image of the accumulated number of total epochs contributing. 

1755 

1756 Returns 

1757 ------- 

1758 maskSpanSetList : `list` [`lsst.afw.geom.SpanSet`] 

1759 List of SpanSets with artifacts. 

1760 """ 

1761 maskSpanSetList = [] 

1762 x0, y0 = epochCountImage.getXY0() 

1763 for i, span in enumerate(spanSetList): 

1764 y, x = span.indices() 

1765 yIdxLocal = [y1 - y0 for y1 in y] 

1766 xIdxLocal = [x1 - x0 for x1 in x] 

1767 outlierN = epochCountImage.array[yIdxLocal, xIdxLocal] 

1768 totalN = nImage.array[yIdxLocal, xIdxLocal] 

1769 

1770 # effectiveMaxNumEpochs is broken line (fraction of N) with 

1771 # characteristic config.maxNumEpochs. 

1772 effMaxNumEpochsHighN = self.config.maxNumEpochs + self.config.maxFractionEpochsHigh * numpy.mean( 

1773 totalN 

1774 ) 

1775 effMaxNumEpochsLowN = self.config.maxFractionEpochsLow * numpy.mean(totalN) 

1776 effectiveMaxNumEpochs = int(min(effMaxNumEpochsLowN, effMaxNumEpochsHighN)) 

1777 nPixelsBelowThreshold = numpy.count_nonzero((outlierN > 0) & (outlierN <= effectiveMaxNumEpochs)) 

1778 percentBelowThreshold = nPixelsBelowThreshold / len(outlierN) 

1779 if percentBelowThreshold > self.config.spatialThreshold: 

1780 maskSpanSetList.append(span) 

1781 

1782 if self.config.doPreserveContainedBySource and footprintsToExclude is not None: 

1783 # If a candidate is contained by a footprint on the template coadd, 

1784 # do not clip. 

1785 filteredMaskSpanSetList = [] 

1786 for span in maskSpanSetList: 

1787 doKeep = True 

1788 for footprint in footprintsToExclude.positive.getFootprints(): 

1789 if footprint.spans.contains(span): 

1790 doKeep = False 

1791 break 

1792 if doKeep: 

1793 filteredMaskSpanSetList.append(span) 

1794 maskSpanSetList = filteredMaskSpanSetList 

1795 

1796 return maskSpanSetList 

1797 

1798 def _readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd): 

1799 """Fetch a warp from the butler and return a warpDiff. 

1800 

1801 Parameters 

1802 ---------- 

1803 warpRef : `lsst.daf.butler.DeferredDatasetHandle` 

1804 Handle for the warp. 

1805 imageScaler : `lsst.pipe.tasks.scaleZeroPoint.ImageScaler` 

1806 An image scaler object. 

1807 templateCoadd : `lsst.afw.image.Exposure` 

1808 Exposure to be substracted from the scaled warp. 

1809 

1810 Returns 

1811 ------- 

1812 warp : `lsst.afw.image.Exposure` 

1813 Exposure of the image difference between the warp and template. 

1814 """ 

1815 # If the PSF-Matched warp did not exist for this direct warp 

1816 # None is holding its place to maintain order in Gen 3 

1817 if warpRef is None: 

1818 return None 

1819 

1820 warp = warpRef.get() 

1821 # direct image scaler OK for PSF-matched Warp 

1822 imageScaler.scaleMaskedImage(warp.getMaskedImage()) 

1823 mi = warp.getMaskedImage() 

1824 if self.config.doScaleWarpVariance: 

1825 try: 

1826 self.scaleWarpVariance.run(mi) 

1827 except Exception as exc: 

1828 self.log.warning("Unable to rescale variance of warp (%s); leaving it as-is", exc) 

1829 mi -= templateCoadd.getMaskedImage() 

1830 return warp