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

629 statements  

« prev     ^ index     » next       coverage.py v7.5.1, created at 2024-05-15 02:55 -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, MaskStreaksTask, 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.scaleZeroPoint import ScaleZeroPointTask 

53from lsst.skymap import BaseSkyMap 

54from lsst.utils.timer import timeMethod 

55 

56log = logging.getLogger(__name__) 

57 

58 

59class AssembleCoaddConnections( 

60 pipeBase.PipelineTaskConnections, 

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

62 defaultTemplates={ 

63 "inputCoaddName": "deep", 

64 "outputCoaddName": "deep", 

65 "warpType": "direct", 

66 "warpTypeSuffix": "", 

67 }, 

68): 

69 inputWarps = pipeBase.connectionTypes.Input( 

70 doc=( 

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

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

73 "warpType config parameter" 

74 ), 

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

76 storageClass="ExposureF", 

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

78 deferLoad=True, 

79 multiple=True, 

80 ) 

81 skyMap = pipeBase.connectionTypes.Input( 

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

83 name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME, 

84 storageClass="SkyMap", 

85 dimensions=("skymap",), 

86 ) 

87 selectedVisits = pipeBase.connectionTypes.Input( 

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

89 name="{outputCoaddName}Visits", 

90 storageClass="StructuredDataDict", 

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

92 ) 

93 brightObjectMask = pipeBase.connectionTypes.PrerequisiteInput( 

94 doc=( 

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

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

97 ), 

98 name="brightObjectMask", 

99 storageClass="ObjectMaskCatalog", 

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

101 minimum=0, 

102 ) 

103 coaddExposure = pipeBase.connectionTypes.Output( 

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

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

106 storageClass="ExposureF", 

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

108 ) 

109 nImage = pipeBase.connectionTypes.Output( 

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

111 name="{outputCoaddName}Coadd_nImage", 

112 storageClass="ImageU", 

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

114 ) 

115 inputMap = pipeBase.connectionTypes.Output( 

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

117 name="{outputCoaddName}Coadd_inputMap", 

118 storageClass="HealSparseMap", 

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

120 ) 

121 

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

123 super().__init__(config=config) 

124 

125 if not config.doMaskBrightObjects: 

126 self.prerequisiteInputs.remove("brightObjectMask") 

127 

128 if not config.doSelectVisits: 

129 self.inputs.remove("selectedVisits") 

130 

131 if not config.doNImage: 

132 self.outputs.remove("nImage") 

133 

134 if not self.config.doInputMap: 

135 self.outputs.remove("inputMap") 

136 

137 

138class AssembleCoaddConfig( 

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

140): 

141 warpType = pexConfig.Field( 

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

143 dtype=str, 

144 default="direct", 

145 ) 

146 subregionSize = pexConfig.ListField( 

147 dtype=int, 

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

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

150 " at once.", 

151 length=2, 

152 default=(2000, 2000), 

153 ) 

154 statistic = pexConfig.Field( 

155 dtype=str, 

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

157 default="MEANCLIP", 

158 ) 

159 doOnlineForMean = pexConfig.Field( 

160 dtype=bool, 

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

162 default=False, 

163 ) 

164 sigmaClip = pexConfig.Field( 

165 dtype=float, 

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

167 default=3.0, 

168 ) 

169 clipIter = pexConfig.Field( 

170 dtype=int, 

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

172 default=2, 

173 ) 

174 calcErrorFromInputVariance = pexConfig.Field( 

175 dtype=bool, 

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

177 "statistic. Passed to " 

178 "StatisticsControl.setCalcErrorFromInputVariance()", 

179 default=True, 

180 ) 

181 scaleZeroPoint = pexConfig.ConfigurableField( 

182 target=ScaleZeroPointTask, 

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

184 ) 

185 doInterp = pexConfig.Field( 

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

187 dtype=bool, 

188 default=True, 

189 ) 

190 interpImage = pexConfig.ConfigurableField( 

191 target=InterpImageTask, 

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

193 ) 

194 doWrite = pexConfig.Field( 

195 doc="Persist coadd?", 

196 dtype=bool, 

197 default=True, 

198 ) 

199 doNImage = pexConfig.Field( 

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

201 dtype=bool, 

202 default=False, 

203 ) 

204 doUsePsfMatchedPolygons = pexConfig.Field( 

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

206 "to True by CompareWarp only.", 

207 dtype=bool, 

208 default=False, 

209 ) 

210 maskPropagationThresholds = pexConfig.DictField( 

211 keytype=str, 

212 itemtype=float, 

213 doc=( 

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

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

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

217 "would have contributed exceeds this value." 

218 ), 

219 default={"SAT": 0.1}, 

220 ) 

221 removeMaskPlanes = pexConfig.ListField( 

222 dtype=str, 

223 default=["NOT_DEBLENDED"], 

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

225 ) 

226 doMaskBrightObjects = pexConfig.Field( 

227 dtype=bool, 

228 default=False, 

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

230 ) 

231 brightObjectMaskName = pexConfig.Field( 

232 dtype=str, 

233 default="BRIGHT_OBJECT", 

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

235 ) 

236 coaddPsf = pexConfig.ConfigField( 

237 doc="Configuration for CoaddPsf", 

238 dtype=measAlg.CoaddPsfConfig, 

239 ) 

240 doAttachTransmissionCurve = pexConfig.Field( 

241 dtype=bool, 

242 default=False, 

243 optional=False, 

244 doc=( 

245 "Attach a piecewise TransmissionCurve for the coadd? " 

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

247 ), 

248 ) 

249 hasFakes = pexConfig.Field( 

250 dtype=bool, 

251 default=False, 

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

253 ) 

254 doSelectVisits = pexConfig.Field( 

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

256 dtype=bool, 

257 default=False, 

258 ) 

259 doInputMap = pexConfig.Field( 

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

261 dtype=bool, 

262 default=False, 

263 ) 

264 inputMapper = pexConfig.ConfigurableField( 

265 doc="Input map creation subtask.", 

266 target=HealSparseInputMapTask, 

267 ) 

268 

269 def setDefaults(self): 

270 super().setDefaults() 

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

272 self.coaddPsf.cacheSize = 0 

273 

274 def validate(self): 

275 super().validate() 

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

277 raise ValueError( 

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

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

280 ) 

281 

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

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

284 stackableStats = [ 

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

286 ] 

287 raise ValueError( 

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

289 ) 

290 

291 

292class AssembleCoaddTask(CoaddBaseTask, pipeBase.PipelineTask): 

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

294 

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

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

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

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

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

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

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

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

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

304 

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

306 

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

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

309 zeropoint for each Warp 

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

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

312 

313 You can retarget these subtasks if you wish. 

314 

315 Raises 

316 ------ 

317 RuntimeError 

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

319 

320 Notes 

321 ----- 

322 Debugging: 

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

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

325 documentation for the subtasks for further information. 

326 """ 

327 

328 ConfigClass = AssembleCoaddConfig 

329 _DefaultName = "assembleCoadd" 

330 

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

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

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

334 if args: 

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

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

337 warnings.warn( 

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

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

340 FutureWarning, 

341 stacklevel=2, 

342 ) 

343 

344 super().__init__(**kwargs) 

345 self.makeSubtask("interpImage") 

346 self.makeSubtask("scaleZeroPoint") 

347 

348 if self.config.doMaskBrightObjects: 

349 mask = afwImage.Mask() 

350 try: 

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

352 except pexExceptions.LsstCppException: 

353 raise RuntimeError( 

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

355 % mask.getMaskPlaneDict().keys() 

356 ) 

357 del mask 

358 

359 if self.config.doInputMap: 

360 self.makeSubtask("inputMapper") 

361 

362 self.warpType = self.config.warpType 

363 

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

365 inputData = butlerQC.get(inputRefs) 

366 

367 # Construct skyInfo expected by run 

368 # Do not remove skyMap from inputData in case _makeSupplementaryData 

369 # needs it 

370 skyMap = inputData["skyMap"] 

371 outputDataId = butlerQC.quantum.dataId 

372 

373 inputData["skyInfo"] = makeSkyInfo( 

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

375 ) 

376 

377 if self.config.doSelectVisits: 

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

379 else: 

380 warpRefList = inputData["inputWarps"] 

381 

382 inputs = self.prepareInputs(warpRefList) 

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

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

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

386 

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

388 retStruct = self.run( 

389 inputData["skyInfo"], 

390 inputs.tempExpRefList, 

391 inputs.imageScalerList, 

392 inputs.weightList, 

393 supplementaryData=supplementaryData, 

394 ) 

395 

396 inputData.setdefault("brightObjectMask", None) 

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

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

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

400 

401 if self.config.doWrite: 

402 butlerQC.put(retStruct, outputRefs) 

403 return retStruct 

404 

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

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

407 

408 Parameters 

409 ---------- 

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

411 The coadded exposure to process. 

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

413 Table of bright objects to mask. 

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

415 Data identification. 

416 """ 

417 if self.config.doInterp: 

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

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

420 varArray = coaddExposure.variance.array 

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

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

423 

424 if self.config.doMaskBrightObjects: 

425 self.setBrightObjectMasks(coaddExposure, brightObjectMasks, dataId) 

426 

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

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

429 

430 Duplicates interface of `runQuantum` method. 

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

432 coadd dataRef for performing preliminary processing before 

433 assembling the coadd. 

434 

435 Parameters 

436 ---------- 

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

438 Gen3 Butler object for fetching additional data products before 

439 running the Task specialized for quantum being processed. 

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

441 Attributes are the names of the connections describing input 

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

443 corresponding dataset type. DataIds are guaranteed to match data 

444 objects in ``inputData``. 

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

446 Attributes are the names of the connections describing output 

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

448 for the corresponding dataset type. 

449 """ 

450 return pipeBase.Struct() 

451 

452 @deprecated( 

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

454 version="v25.0", 

455 category=FutureWarning, 

456 ) 

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

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

459 

460 def prepareInputs(self, refList): 

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

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

463 

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

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

466 normalize the photometric zeropoint and compute the weight for each 

467 Warp. 

468 

469 Parameters 

470 ---------- 

471 refList : `list` 

472 List of data references to tempExp. 

473 

474 Returns 

475 ------- 

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

477 Results as a struct with attributes: 

478 

479 ``tempExprefList`` 

480 `list` of data references to tempExp. 

481 ``weightList`` 

482 `list` of weightings. 

483 ``imageScalerList`` 

484 `list` of image scalers. 

485 """ 

486 statsCtrl = afwMath.StatisticsControl() 

487 statsCtrl.setNumSigmaClip(self.config.sigmaClip) 

488 statsCtrl.setNumIter(self.config.clipIter) 

489 statsCtrl.setAndMask(self.getBadPixelMask()) 

490 statsCtrl.setNanSafe(True) 

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

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

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

494 # tempExp. 

495 tempExpRefList = [] 

496 weightList = [] 

497 imageScalerList = [] 

498 tempExpName = self.getTempExpDatasetName(self.warpType) 

499 for tempExpRef in refList: 

500 tempExp = tempExpRef.get() 

501 # Ignore any input warp that is empty of data 

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

503 continue 

504 maskedImage = tempExp.getMaskedImage() 

505 imageScaler = self.scaleZeroPoint.computeImageScaler( 

506 exposure=tempExp, 

507 dataRef=tempExpRef, # FIXME 

508 ) 

509 try: 

510 imageScaler.scaleMaskedImage(maskedImage) 

511 except Exception as e: 

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

513 continue 

514 statObj = afwMath.makeStatistics( 

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

516 ) 

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

518 weight = 1.0 / float(meanVar) 

519 if not numpy.isfinite(weight): 

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

521 continue 

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

523 

524 del maskedImage 

525 del tempExp 

526 

527 tempExpRefList.append(tempExpRef) 

528 weightList.append(weight) 

529 imageScalerList.append(imageScaler) 

530 

531 return pipeBase.Struct( 

532 tempExpRefList=tempExpRefList, weightList=weightList, imageScalerList=imageScalerList 

533 ) 

534 

535 def prepareStats(self, mask=None): 

536 """Prepare the statistics for coadding images. 

537 

538 Parameters 

539 ---------- 

540 mask : `int`, optional 

541 Bit mask value to exclude from coaddition. 

542 

543 Returns 

544 ------- 

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

546 Statistics as a struct with attributes: 

547 

548 ``statsCtrl`` 

549 Statistics control object for coadd 

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

551 ``statsFlags`` 

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

553 """ 

554 if mask is None: 

555 mask = self.getBadPixelMask() 

556 statsCtrl = afwMath.StatisticsControl() 

557 statsCtrl.setNumSigmaClip(self.config.sigmaClip) 

558 statsCtrl.setNumIter(self.config.clipIter) 

559 statsCtrl.setAndMask(mask) 

560 statsCtrl.setNanSafe(True) 

561 statsCtrl.setWeighted(True) 

562 statsCtrl.setCalcErrorFromInputVariance(self.config.calcErrorFromInputVariance) 

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

564 bit = afwImage.Mask.getMaskPlane(plane) 

565 statsCtrl.setMaskPropagationThreshold(bit, threshold) 

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

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

568 

569 @timeMethod 

570 def run( 

571 self, 

572 skyInfo, 

573 tempExpRefList, 

574 imageScalerList, 

575 weightList, 

576 altMaskList=None, 

577 mask=None, 

578 supplementaryData=None, 

579 ): 

580 """Assemble a coadd from input warps. 

581 

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

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

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

585 conserve memory usage. Iterate over subregions within the outer 

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

587 subregions from the coaddTempExps with the statistic specified. 

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

589 

590 Parameters 

591 ---------- 

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

593 Struct with geometric information about the patch. 

594 tempExpRefList : `list` 

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

596 imageScalerList : `list` 

597 List of image scalers. 

598 weightList : `list` 

599 List of weights. 

600 altMaskList : `list`, optional 

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

602 tempExp. 

603 mask : `int`, optional 

604 Bit mask value to exclude from coaddition. 

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

606 Struct with additional data products needed to assemble coadd. 

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

608 and override `run`. 

609 

610 Returns 

611 ------- 

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

613 Results as a struct with attributes: 

614 

615 ``coaddExposure`` 

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

617 ``nImage`` 

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

619 ``inputMap`` 

620 Bit-wise map of inputs, if requested. 

621 ``warpRefList`` 

622 Input list of refs to the warps 

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

624 ``imageScalerList`` 

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

626 ``weightList`` 

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

628 

629 Raises 

630 ------ 

631 lsst.pipe.base.NoWorkFound 

632 Raised if no data references are provided. 

633 """ 

634 tempExpName = self.getTempExpDatasetName(self.warpType) 

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

636 if not tempExpRefList: 

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

638 

639 stats = self.prepareStats(mask=mask) 

640 

641 if altMaskList is None: 

642 altMaskList = [None] * len(tempExpRefList) 

643 

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

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

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

647 self.assembleMetadata(coaddExposure, tempExpRefList, weightList) 

648 coaddMaskedImage = coaddExposure.getMaskedImage() 

649 subregionSizeArr = self.config.subregionSize 

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

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

652 # assembleSubregion. 

653 if self.config.doNImage: 

654 nImage = afwImage.ImageU(skyInfo.bbox) 

655 else: 

656 nImage = None 

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

658 # masked in assembleSubregion. 

659 if self.config.doInputMap: 

660 self.inputMapper.build_ccd_input_map( 

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

662 ) 

663 

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

665 try: 

666 self.assembleOnlineMeanCoadd( 

667 coaddExposure, 

668 tempExpRefList, 

669 imageScalerList, 

670 weightList, 

671 altMaskList, 

672 stats.ctrl, 

673 nImage=nImage, 

674 ) 

675 except Exception as e: 

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

677 raise 

678 else: 

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

680 try: 

681 self.assembleSubregion( 

682 coaddExposure, 

683 subBBox, 

684 tempExpRefList, 

685 imageScalerList, 

686 weightList, 

687 altMaskList, 

688 stats.flags, 

689 stats.ctrl, 

690 nImage=nImage, 

691 ) 

692 except Exception as e: 

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

694 raise 

695 

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

697 # accumulation. 

698 if self.config.doInputMap: 

699 self.inputMapper.finalize_ccd_input_map_mask() 

700 inputMap = self.inputMapper.ccd_input_map 

701 else: 

702 inputMap = None 

703 

704 self.setInexactPsf(coaddMaskedImage.getMask()) 

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

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

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

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

709 return pipeBase.Struct( 

710 coaddExposure=coaddExposure, 

711 nImage=nImage, 

712 warpRefList=tempExpRefList, 

713 imageScalerList=imageScalerList, 

714 weightList=weightList, 

715 inputMap=inputMap, 

716 ) 

717 

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

719 """Set the metadata for the coadd. 

720 

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

722 

723 Parameters 

724 ---------- 

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

726 The target exposure for the coadd. 

727 tempExpRefList : `list` 

728 List of data references to tempExp. 

729 weightList : `list` 

730 List of weights. 

731 

732 Raises 

733 ------ 

734 AssertionError 

735 Raised if there is a length mismatch. 

736 """ 

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

738 

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

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

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

742 # (see #2777). 

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

744 

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

746 

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

748 

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

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

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

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

753 coaddInputs.ccds.reserve(numCcds) 

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

755 

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

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

758 

759 if self.config.doUsePsfMatchedPolygons: 

760 self.shrinkValidPolygons(coaddInputs) 

761 

762 coaddInputs.visits.sort() 

763 coaddInputs.ccds.sort() 

764 if self.warpType == "psfMatched": 

765 # The modelPsf BBox for a psfMatchedWarp/coaddTempExp was 

766 # dynamically defined by ModelPsfMatchTask as the square box 

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

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

769 # having the maximum width (sufficient because square) 

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

771 modelPsfWidthList = [ 

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

773 ] 

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

775 else: 

776 psf = measAlg.CoaddPsf( 

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

778 ) 

779 coaddExposure.setPsf(psf) 

780 apCorrMap = measAlg.makeCoaddApCorrMap( 

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

782 ) 

783 coaddExposure.getInfo().setApCorrMap(apCorrMap) 

784 if self.config.doAttachTransmissionCurve: 

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

786 coaddExposure.getInfo().setTransmissionCurve(transmissionCurve) 

787 

788 def assembleSubregion( 

789 self, 

790 coaddExposure, 

791 bbox, 

792 tempExpRefList, 

793 imageScalerList, 

794 weightList, 

795 altMaskList, 

796 statsFlags, 

797 statsCtrl, 

798 nImage=None, 

799 ): 

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

801 

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

803 if one is passed. Remove mask planes listed in 

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

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

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

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

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

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

810 

811 Parameters 

812 ---------- 

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

814 The target exposure for the coadd. 

815 bbox : `lsst.geom.Box` 

816 Sub-region to coadd. 

817 tempExpRefList : `list` 

818 List of data reference to tempExp. 

819 imageScalerList : `list` 

820 List of image scalers. 

821 weightList : `list` 

822 List of weights. 

823 altMaskList : `list` 

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

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

826 name to which to add the spans. 

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

828 Property object for statistic for coadd. 

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

830 Statistics control object for coadd. 

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

832 Keeps track of exposure count for each pixel. 

833 """ 

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

835 

836 coaddExposure.mask.addMaskPlane("REJECTED") 

837 coaddExposure.mask.addMaskPlane("CLIPPED") 

838 coaddExposure.mask.addMaskPlane("SENSOR_EDGE") 

839 maskMap = self.setRejectedMaskMapping(statsCtrl) 

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

841 maskedImageList = [] 

842 if nImage is not None: 

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

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

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

846 

847 maskedImage = exposure.getMaskedImage() 

848 mask = maskedImage.getMask() 

849 if altMask is not None: 

850 self.applyAltMaskPlanes(mask, altMask) 

851 imageScaler.scaleMaskedImage(maskedImage) 

852 

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

854 # In legacyCoadd, pixels may also be excluded by 

855 # afwMath.statisticsStack. 

856 if nImage is not None: 

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

858 if self.config.removeMaskPlanes: 

859 self.removeMaskPlanes(maskedImage) 

860 maskedImageList.append(maskedImage) 

861 

862 if self.config.doInputMap: 

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

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

865 

866 with self.timer("stack"): 

867 coaddSubregion = afwMath.statisticsStack( 

868 maskedImageList, 

869 statsFlags, 

870 statsCtrl, 

871 weightList, 

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

873 maskMap, 

874 ) 

875 coaddExposure.maskedImage.assign(coaddSubregion, bbox) 

876 if nImage is not None: 

877 nImage.assign(subNImage, bbox) 

878 

879 def assembleOnlineMeanCoadd( 

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

881 ): 

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

883 

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

885 It only works for MEAN statistics. 

886 

887 Parameters 

888 ---------- 

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

890 The target exposure for the coadd. 

891 tempExpRefList : `list` 

892 List of data reference to tempExp. 

893 imageScalerList : `list` 

894 List of image scalers. 

895 weightList : `list` 

896 List of weights. 

897 altMaskList : `list` 

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

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

900 name to which to add the spans. 

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

902 Statistics control object for coadd. 

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

904 Keeps track of exposure count for each pixel. 

905 """ 

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

907 

908 coaddExposure.mask.addMaskPlane("REJECTED") 

909 coaddExposure.mask.addMaskPlane("CLIPPED") 

910 coaddExposure.mask.addMaskPlane("SENSOR_EDGE") 

911 maskMap = self.setRejectedMaskMapping(statsCtrl) 

912 thresholdDict = AccumulatorMeanStack.stats_ctrl_to_threshold_dict(statsCtrl) 

913 

914 bbox = coaddExposure.maskedImage.getBBox() 

915 

916 stacker = AccumulatorMeanStack( 

917 coaddExposure.image.array.shape, 

918 statsCtrl.getAndMask(), 

919 mask_threshold_dict=thresholdDict, 

920 mask_map=maskMap, 

921 no_good_pixels_mask=statsCtrl.getNoGoodPixelsMask(), 

922 calc_error_from_input_variance=self.config.calcErrorFromInputVariance, 

923 compute_n_image=(nImage is not None), 

924 ) 

925 

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

927 tempExpRefList, imageScalerList, altMaskList, weightList 

928 ): 

929 exposure = tempExpRef.get() 

930 maskedImage = exposure.getMaskedImage() 

931 mask = maskedImage.getMask() 

932 if altMask is not None: 

933 self.applyAltMaskPlanes(mask, altMask) 

934 imageScaler.scaleMaskedImage(maskedImage) 

935 if self.config.removeMaskPlanes: 

936 self.removeMaskPlanes(maskedImage) 

937 

938 stacker.add_masked_image(maskedImage, weight=weight) 

939 

940 if self.config.doInputMap: 

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

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

943 

944 stacker.fill_stacked_masked_image(coaddExposure.maskedImage) 

945 

946 if nImage is not None: 

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

948 

949 def removeMaskPlanes(self, maskedImage): 

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

951 

952 Parameters 

953 ---------- 

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

955 The masked image to be modified. 

956 

957 Raises 

958 ------ 

959 InvalidParameterError 

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

961 """ 

962 mask = maskedImage.getMask() 

963 for maskPlane in self.config.removeMaskPlanes: 

964 try: 

965 mask &= ~mask.getPlaneBitMask(maskPlane) 

966 except pexExceptions.InvalidParameterError: 

967 self.log.debug( 

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

969 ) 

970 

971 @staticmethod 

972 def setRejectedMaskMapping(statsCtrl): 

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

974 

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

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

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

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

979 

980 Parameters 

981 ---------- 

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

983 Statistics control object for coadd. 

984 

985 Returns 

986 ------- 

987 maskMap : `list` of `tuple` of `int` 

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

989 mask planes of the coadd. 

990 """ 

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

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

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

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

995 maskMap = [ 

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

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

998 (clipped, clipped), 

999 ] 

1000 return maskMap 

1001 

1002 def applyAltMaskPlanes(self, mask, altMaskSpans): 

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

1004 

1005 Parameters 

1006 ---------- 

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

1008 Original mask. 

1009 altMaskSpans : `dict` 

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

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

1012 and list of SpanSets to apply to the mask. 

1013 

1014 Returns 

1015 ------- 

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

1017 Updated mask. 

1018 """ 

1019 if self.config.doUsePsfMatchedPolygons: 

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

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

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

1023 # still be rejected because of NO_DATA. 

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

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

1026 # *before* applying the new masks below. 

1027 for spanSet in altMaskSpans["NO_DATA"]: 

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

1029 

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

1031 maskClipValue = mask.addMaskPlane(plane) 

1032 for spanSet in spanSetList: 

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

1034 return mask 

1035 

1036 def shrinkValidPolygons(self, coaddInputs): 

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

1038 

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

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

1041 

1042 Parameters 

1043 ---------- 

1044 coaddInputs : `lsst.afw.image.coaddInputs` 

1045 Original mask. 

1046 """ 

1047 for ccd in coaddInputs.ccds: 

1048 polyOrig = ccd.getValidPolygon() 

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

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

1051 if polyOrig: 

1052 validPolygon = polyOrig.intersectionSingle(validPolyBBox) 

1053 else: 

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

1055 ccd.setValidPolygon(validPolygon) 

1056 

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

1058 """Set the bright object masks. 

1059 

1060 Parameters 

1061 ---------- 

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

1063 Exposure under consideration. 

1064 brightObjectMasks : `lsst.afw.table` 

1065 Table of bright objects to mask. 

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

1067 Data identifier dict for patch. 

1068 """ 

1069 if brightObjectMasks is None: 

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

1071 return 

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

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

1074 wcs = exposure.getWcs() 

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

1076 

1077 for rec in brightObjectMasks: 

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

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

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

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

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

1083 

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

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

1086 

1087 bbox = geom.BoxI( 

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

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

1090 ) 

1091 spans = afwGeom.SpanSet(bbox) 

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

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

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

1095 else: 

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

1097 continue 

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

1099 

1100 def setInexactPsf(self, mask): 

1101 """Set INEXACT_PSF mask plane. 

1102 

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

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

1105 these pixels. 

1106 

1107 Parameters 

1108 ---------- 

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

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

1111 """ 

1112 mask.addMaskPlane("INEXACT_PSF") 

1113 inexactPsf = mask.getPlaneBitMask("INEXACT_PSF") 

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

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

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

1117 array = mask.getArray() 

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

1119 array[selected] |= inexactPsf 

1120 

1121 def filterWarps(self, inputs, goodVisits): 

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

1123 goodVisit. 

1124 

1125 Parameters 

1126 ---------- 

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

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

1129 containing visit. 

1130 goodVisit : `dict` 

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

1132 

1133 Returns 

1134 ------- 

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

1136 Filtered and sorted list of inputRefs with visitId in goodVisits 

1137 ordered by goodVisit. 

1138 """ 

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

1140 filteredInputs = [] 

1141 for visit in goodVisits.keys(): 

1142 if visit in inputWarpDict: 

1143 filteredInputs.append(inputWarpDict[visit]) 

1144 return filteredInputs 

1145 

1146 

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

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

1149 footprint. 

1150 

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

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

1153 ignoreMask set. Return the count. 

1154 

1155 Parameters 

1156 ---------- 

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

1158 Mask to define intersection region by. 

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

1160 Footprint to define the intersection region by. 

1161 bitmask : `Unknown` 

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

1163 ignoreMask : `Unknown` 

1164 Pixels to not consider. 

1165 

1166 Returns 

1167 ------- 

1168 result : `int` 

1169 Number of pixels in footprint with specified mask. 

1170 """ 

1171 bbox = footprint.getBBox() 

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

1173 fp = afwImage.Mask(bbox) 

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

1175 footprint.spans.setMask(fp, bitmask) 

1176 return numpy.logical_and( 

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

1178 ).sum() 

1179 

1180 

1181class CompareWarpAssembleCoaddConnections(AssembleCoaddConnections): 

1182 psfMatchedWarps = pipeBase.connectionTypes.Input( 

1183 doc=( 

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

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

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

1187 ), 

1188 name="{inputCoaddName}Coadd_psfMatchedWarp", 

1189 storageClass="ExposureF", 

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

1191 deferLoad=True, 

1192 multiple=True, 

1193 ) 

1194 templateCoadd = pipeBase.connectionTypes.Output( 

1195 doc=( 

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

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

1198 ), 

1199 name="{outputCoaddName}CoaddPsfMatched", 

1200 storageClass="ExposureF", 

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

1202 ) 

1203 

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

1205 super().__init__(config=config) 

1206 if not config.assembleStaticSkyModel.doWrite: 

1207 self.outputs.remove("templateCoadd") 

1208 config.validate() 

1209 

1210 

1211class CompareWarpAssembleCoaddConfig( 

1212 AssembleCoaddConfig, pipelineConnections=CompareWarpAssembleCoaddConnections 

1213): 

1214 assembleStaticSkyModel = pexConfig.ConfigurableField( 

1215 target=AssembleCoaddTask, 

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

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

1218 ) 

1219 detect = pexConfig.ConfigurableField( 

1220 target=SourceDetectionTask, 

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

1222 ) 

1223 detectTemplate = pexConfig.ConfigurableField( 

1224 target=SourceDetectionTask, 

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

1226 ) 

1227 maskStreaks = pexConfig.ConfigurableField( 

1228 target=MaskStreaksTask, 

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

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

1231 "streakMaskName", 

1232 ) 

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

1234 maxNumEpochs = pexConfig.Field( 

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

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

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

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

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

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

1241 "than transient and not masked.", 

1242 dtype=int, 

1243 default=2, 

1244 ) 

1245 maxFractionEpochsLow = pexConfig.RangeField( 

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

1247 "Effective maxNumEpochs = " 

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

1249 dtype=float, 

1250 default=0.4, 

1251 min=0.0, 

1252 max=1.0, 

1253 ) 

1254 maxFractionEpochsHigh = pexConfig.RangeField( 

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

1256 "Effective maxNumEpochs = " 

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

1258 dtype=float, 

1259 default=0.03, 

1260 min=0.0, 

1261 max=1.0, 

1262 ) 

1263 spatialThreshold = pexConfig.RangeField( 

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

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

1266 dtype=float, 

1267 default=0.5, 

1268 min=0.0, 

1269 max=1.0, 

1270 inclusiveMin=True, 

1271 inclusiveMax=True, 

1272 ) 

1273 doScaleWarpVariance = pexConfig.Field( 

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

1275 dtype=bool, 

1276 default=True, 

1277 ) 

1278 scaleWarpVariance = pexConfig.ConfigurableField( 

1279 target=ScaleVarianceTask, 

1280 doc="Rescale variance on warps", 

1281 ) 

1282 doPreserveContainedBySource = pexConfig.Field( 

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

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

1285 dtype=bool, 

1286 default=True, 

1287 ) 

1288 doPrefilterArtifacts = pexConfig.Field( 

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

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

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

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

1293 dtype=bool, 

1294 default=True, 

1295 ) 

1296 prefilterArtifactsMaskPlanes = pexConfig.ListField( 

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

1298 dtype=str, 

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

1300 ) 

1301 prefilterArtifactsRatio = pexConfig.Field( 

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

1303 dtype=float, 

1304 default=0.05, 

1305 ) 

1306 doFilterMorphological = pexConfig.Field( 

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

1308 "be streaks.", 

1309 dtype=bool, 

1310 default=False, 

1311 ) 

1312 growStreakFp = pexConfig.Field( 

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

1314 ) 

1315 

1316 def setDefaults(self): 

1317 AssembleCoaddConfig.setDefaults(self) 

1318 self.statistic = "MEAN" 

1319 self.doUsePsfMatchedPolygons = True 

1320 

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

1322 # matching kernel. CompareWarp applies psfMatched EDGE pixels to 

1323 # directWarps before assembling. 

1324 if "EDGE" in self.badMaskPlanes: 

1325 self.badMaskPlanes.remove("EDGE") 

1326 self.removeMaskPlanes.append("EDGE") 

1327 self.assembleStaticSkyModel.badMaskPlanes = [ 

1328 "NO_DATA", 

1329 ] 

1330 self.assembleStaticSkyModel.warpType = "psfMatched" 

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

1332 self.assembleStaticSkyModel.statistic = "MEANCLIP" 

1333 self.assembleStaticSkyModel.sigmaClip = 2.5 

1334 self.assembleStaticSkyModel.clipIter = 3 

1335 self.assembleStaticSkyModel.calcErrorFromInputVariance = False 

1336 self.assembleStaticSkyModel.doWrite = False 

1337 self.detect.doTempLocalBackground = False 

1338 self.detect.reEstimateBackground = False 

1339 self.detect.returnOriginalFootprints = False 

1340 self.detect.thresholdPolarity = "both" 

1341 self.detect.thresholdValue = 5 

1342 self.detect.minPixels = 4 

1343 self.detect.isotropicGrow = True 

1344 self.detect.thresholdType = "pixel_stdev" 

1345 self.detect.nSigmaToGrow = 0.4 

1346 # The default nSigmaToGrow for SourceDetectionTask is already 2.4, 

1347 # Explicitly restating because ratio with detect.nSigmaToGrow matters 

1348 self.detectTemplate.nSigmaToGrow = 2.4 

1349 self.detectTemplate.doTempLocalBackground = False 

1350 self.detectTemplate.reEstimateBackground = False 

1351 self.detectTemplate.returnOriginalFootprints = False 

1352 

1353 def validate(self): 

1354 super().validate() 

1355 if self.assembleStaticSkyModel.doNImage: 

1356 raise ValueError( 

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

1358 "Please set assembleStaticSkyModel.doNImage=False" 

1359 ) 

1360 

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

1362 raise ValueError( 

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

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

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

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

1367 ) 

1368 

1369 

1370class CompareWarpAssembleCoaddTask(AssembleCoaddTask): 

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

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

1373 

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

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

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

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

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

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

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

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

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

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

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

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

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

1387 ``temporalThreshold`` and ``spatialThreshold``. The temporalThreshold is 

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

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

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

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

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

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

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

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

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

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

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

1399 surveys. 

1400 

1401 ``CompareWarpAssembleCoaddTask`` sub-classes 

1402 ``AssembleCoaddTask`` and instantiates ``AssembleCoaddTask`` 

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

1404 

1405 Notes 

1406 ----- 

1407 Debugging: 

1408 This task supports the following debug variables: 

1409 - ``saveCountIm`` 

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

1411 - ``figPath`` 

1412 Path to save the debug fits images and figures 

1413 """ 

1414 

1415 ConfigClass = CompareWarpAssembleCoaddConfig 

1416 _DefaultName = "compareWarpAssembleCoadd" 

1417 

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

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

1420 self.makeSubtask("assembleStaticSkyModel") 

1421 detectionSchema = afwTable.SourceTable.makeMinimalSchema() 

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

1423 if self.config.doPreserveContainedBySource: 

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

1425 if self.config.doScaleWarpVariance: 

1426 self.makeSubtask("scaleWarpVariance") 

1427 if self.config.doFilterMorphological: 

1428 self.makeSubtask("maskStreaks") 

1429 

1430 @utils.inheritDoc(AssembleCoaddTask) 

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

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

1433 subtract from PSF-Matched warps. 

1434 

1435 Returns 

1436 ------- 

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

1438 Results as a struct with attributes: 

1439 

1440 ``templateCoadd`` 

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

1442 ``nImage`` 

1443 Keeps track of exposure count for each pixel 

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

1445 

1446 Raises 

1447 ------ 

1448 RuntimeError 

1449 Raised if ``templateCoadd`` is `None`. 

1450 """ 

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

1452 # generation. 

1453 staticSkyModelInputRefs = copy.deepcopy(inputRefs) 

1454 staticSkyModelInputRefs.inputWarps = inputRefs.psfMatchedWarps 

1455 

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

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

1458 staticSkyModelOutputRefs = copy.deepcopy(outputRefs) 

1459 if self.config.assembleStaticSkyModel.doWrite: 

1460 staticSkyModelOutputRefs.coaddExposure = staticSkyModelOutputRefs.templateCoadd 

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

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

1463 del outputRefs.templateCoadd 

1464 del staticSkyModelOutputRefs.templateCoadd 

1465 

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

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

1468 del staticSkyModelOutputRefs.nImage 

1469 

1470 templateCoadd = self.assembleStaticSkyModel.runQuantum( 

1471 butlerQC, staticSkyModelInputRefs, staticSkyModelOutputRefs 

1472 ) 

1473 if templateCoadd is None: 

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

1475 

1476 return pipeBase.Struct( 

1477 templateCoadd=templateCoadd.coaddExposure, 

1478 nImage=templateCoadd.nImage, 

1479 warpRefList=templateCoadd.warpRefList, 

1480 imageScalerList=templateCoadd.imageScalerList, 

1481 weightList=templateCoadd.weightList, 

1482 ) 

1483 

1484 def _noTemplateMessage(self, warpType): 

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

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

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

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

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

1490 

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

1492 another algorithm like: 

1493 

1494 from lsst.pipe.tasks.assembleCoadd import SafeClipAssembleCoaddTask 

1495 config.assemble.retarget(SafeClipAssembleCoaddTask) 

1496 """ % { 

1497 "warpName": warpName 

1498 } 

1499 return message 

1500 

1501 @utils.inheritDoc(AssembleCoaddTask) 

1502 @timeMethod 

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

1504 """Notes 

1505 ----- 

1506 Assemble the coadd. 

1507 

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

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

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

1511 method. 

1512 """ 

1513 # Check and match the order of the supplementaryData 

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

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

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

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

1518 

1519 if dataIds != psfMatchedDataIds: 

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

1521 supplementaryData.warpRefList = reorderAndPadList( 

1522 supplementaryData.warpRefList, psfMatchedDataIds, dataIds 

1523 ) 

1524 supplementaryData.imageScalerList = reorderAndPadList( 

1525 supplementaryData.imageScalerList, psfMatchedDataIds, dataIds 

1526 ) 

1527 

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

1529 # artifacts. 

1530 spanSetMaskList = self.findArtifacts( 

1531 supplementaryData.templateCoadd, supplementaryData.warpRefList, supplementaryData.imageScalerList 

1532 ) 

1533 

1534 badMaskPlanes = self.config.badMaskPlanes[:] 

1535 badMaskPlanes.append("CLIPPED") 

1536 badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes) 

1537 

1538 result = AssembleCoaddTask.run( 

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

1540 ) 

1541 

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

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

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

1545 return result 

1546 

1547 def applyAltEdgeMask(self, mask, altMaskList): 

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

1549 

1550 Parameters 

1551 ---------- 

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

1553 Original mask. 

1554 altMaskList : `list` of `dict` 

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

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

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

1558 the mask. 

1559 """ 

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

1561 for visitMask in altMaskList: 

1562 if "EDGE" in visitMask: 

1563 for spanSet in visitMask["EDGE"]: 

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

1565 

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

1567 """Find artifacts. 

1568 

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

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

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

1572 difference image and filters the artifacts detected in each using 

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

1574 difficult to subtract cleanly. 

1575 

1576 Parameters 

1577 ---------- 

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

1579 Exposure to serve as model of static sky. 

1580 tempExpRefList : `list` 

1581 List of data references to warps. 

1582 imageScalerList : `list` 

1583 List of image scalers. 

1584 

1585 Returns 

1586 ------- 

1587 altMasks : `list` of `dict` 

1588 List of dicts containing information about CLIPPED 

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

1590 """ 

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

1592 coaddBBox = templateCoadd.getBBox() 

1593 slateIm = afwImage.ImageU(coaddBBox) 

1594 epochCountImage = afwImage.ImageU(coaddBBox) 

1595 nImage = afwImage.ImageU(coaddBBox) 

1596 spanSetArtifactList = [] 

1597 spanSetNoDataMaskList = [] 

1598 spanSetEdgeList = [] 

1599 spanSetBadMorphoList = [] 

1600 badPixelMask = self.getBadPixelMask() 

1601 

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

1603 templateCoadd.mask.clearAllMaskPlanes() 

1604 

1605 if self.config.doPreserveContainedBySource: 

1606 templateFootprints = self.detectTemplate.detectFootprints(templateCoadd) 

1607 else: 

1608 templateFootprints = None 

1609 

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

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

1612 if warpDiffExp is not None: 

1613 # This nImage only approximates the final nImage because it 

1614 # uses the PSF-matched mask. 

1615 nImage.array += ( 

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

1617 ).astype(numpy.uint16) 

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

1619 fpSet.positive.merge(fpSet.negative) 

1620 footprints = fpSet.positive 

1621 slateIm.set(0) 

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

1623 

1624 # Remove artifacts due to defects before they contribute to 

1625 # the epochCountImage. 

1626 if self.config.doPrefilterArtifacts: 

1627 spanSetList = self.prefilterArtifacts(spanSetList, warpDiffExp) 

1628 

1629 # Clear mask before adding prefiltered spanSets 

1630 self.detect.clearMask(warpDiffExp.mask) 

1631 for spans in spanSetList: 

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

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

1634 epochCountImage += slateIm 

1635 

1636 if self.config.doFilterMorphological: 

1637 maskName = self.config.streakMaskName 

1638 # clear single frame streak mask if it exists already 

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

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

1641 else: 

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

1643 

1644 _ = self.maskStreaks.run(warpDiffExp) 

1645 streakMask = warpDiffExp.mask 

1646 spanSetStreak = afwGeom.SpanSet.fromMask( 

1647 streakMask, streakMask.getPlaneBitMask(maskName) 

1648 ).split() 

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

1650 # wings. 

1651 psf = warpDiffExp.getPsf() 

1652 for s, sset in enumerate(spanSetStreak): 

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

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

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

1656 spanSetStreak[s] = sset_dilated 

1657 

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

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

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

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

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

1663 # the direct warp. 

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

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

1666 nansMask.setXY0(warpDiffExp.getXY0()) 

1667 edgeMask = warpDiffExp.mask 

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

1669 else: 

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

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

1672 nansMask = afwImage.MaskX(coaddBBox, 1) 

1673 spanSetList = [] 

1674 spanSetEdgeMask = [] 

1675 spanSetStreak = [] 

1676 

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

1678 

1679 spanSetNoDataMaskList.append(spanSetNoDataMask) 

1680 spanSetArtifactList.append(spanSetList) 

1681 spanSetEdgeList.append(spanSetEdgeMask) 

1682 if self.config.doFilterMorphological: 

1683 spanSetBadMorphoList.append(spanSetStreak) 

1684 

1685 if lsstDebug.Info(__name__).saveCountIm: 

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

1687 epochCountImage.writeFits(path) 

1688 

1689 for i, spanSetList in enumerate(spanSetArtifactList): 

1690 if spanSetList: 

1691 filteredSpanSetList = self.filterArtifacts( 

1692 spanSetList, epochCountImage, nImage, templateFootprints 

1693 ) 

1694 spanSetArtifactList[i] = filteredSpanSetList 

1695 if self.config.doFilterMorphological: 

1696 spanSetArtifactList[i] += spanSetBadMorphoList[i] 

1697 

1698 altMasks = [] 

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

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

1701 return altMasks 

1702 

1703 def prefilterArtifacts(self, spanSetList, exp): 

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

1705 

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

1707 temporal information should go in this method. 

1708 

1709 Parameters 

1710 ---------- 

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

1712 List of SpanSets representing artifact candidates. 

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

1714 Exposure containing mask planes used to prefilter. 

1715 

1716 Returns 

1717 ------- 

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

1719 List of SpanSets with artifacts. 

1720 """ 

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

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

1723 returnSpanSetList = [] 

1724 bbox = exp.getBBox() 

1725 x0, y0 = exp.getXY0() 

1726 for i, span in enumerate(spanSetList): 

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

1728 yIndexLocal = numpy.array(y) - y0 

1729 xIndexLocal = numpy.array(x) - x0 

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

1731 if goodRatio > self.config.prefilterArtifactsRatio: 

1732 returnSpanSetList.append(span) 

1733 return returnSpanSetList 

1734 

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

1736 """Filter artifact candidates. 

1737 

1738 Parameters 

1739 ---------- 

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

1741 List of SpanSets representing artifact candidates. 

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

1743 Image of accumulated number of warpDiff detections. 

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

1745 Image of the accumulated number of total epochs contributing. 

1746 

1747 Returns 

1748 ------- 

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

1750 List of SpanSets with artifacts. 

1751 """ 

1752 maskSpanSetList = [] 

1753 x0, y0 = epochCountImage.getXY0() 

1754 for i, span in enumerate(spanSetList): 

1755 y, x = span.indices() 

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

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

1758 outlierN = epochCountImage.array[yIdxLocal, xIdxLocal] 

1759 totalN = nImage.array[yIdxLocal, xIdxLocal] 

1760 

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

1762 # characteristic config.maxNumEpochs. 

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

1764 totalN 

1765 ) 

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

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

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

1769 percentBelowThreshold = nPixelsBelowThreshold / len(outlierN) 

1770 if percentBelowThreshold > self.config.spatialThreshold: 

1771 maskSpanSetList.append(span) 

1772 

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

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

1775 # do not clip. 

1776 filteredMaskSpanSetList = [] 

1777 for span in maskSpanSetList: 

1778 doKeep = True 

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

1780 if footprint.spans.contains(span): 

1781 doKeep = False 

1782 break 

1783 if doKeep: 

1784 filteredMaskSpanSetList.append(span) 

1785 maskSpanSetList = filteredMaskSpanSetList 

1786 

1787 return maskSpanSetList 

1788 

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

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

1791 

1792 Parameters 

1793 ---------- 

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

1795 Handle for the warp. 

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

1797 An image scaler object. 

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

1799 Exposure to be substracted from the scaled warp. 

1800 

1801 Returns 

1802 ------- 

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

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

1805 """ 

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

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

1808 if warpRef is None: 

1809 return None 

1810 

1811 warp = warpRef.get() 

1812 # direct image scaler OK for PSF-matched Warp 

1813 imageScaler.scaleMaskedImage(warp.getMaskedImage()) 

1814 mi = warp.getMaskedImage() 

1815 if self.config.doScaleWarpVariance: 

1816 try: 

1817 self.scaleWarpVariance.run(mi) 

1818 except Exception as exc: 

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

1820 mi -= templateCoadd.getMaskedImage() 

1821 return warp