lsst.pipe.tasks gc8e401de96+07adebc1ff
Loading...
Searching...
No Matches
assembleCoadd.py
Go to the documentation of this file.
1# This file is part of pipe_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__ = ["AssembleCoaddTask", "AssembleCoaddConnections", "AssembleCoaddConfig",
23 "CompareWarpAssembleCoaddTask", "CompareWarpAssembleCoaddConfig"]
24
25import copy
26import numpy
27import warnings
28import logging
29import lsst.pex.config as pexConfig
30import lsst.pex.exceptions as pexExceptions
31import lsst.geom as geom
32import lsst.afw.geom as afwGeom
33import lsst.afw.image as afwImage
34import lsst.afw.math as afwMath
35import lsst.afw.table as afwTable
36import lsst.coadd.utils as coaddUtils
37import lsst.pipe.base as pipeBase
38import lsst.meas.algorithms as measAlg
39import lsstDebug
40import lsst.utils as utils
41from lsst.skymap import BaseSkyMap
42from .coaddBase import CoaddBaseTask, makeSkyInfo, reorderAndPadList
43from .interpImage import InterpImageTask
44from .scaleZeroPoint import ScaleZeroPointTask
45from .maskStreaks import MaskStreaksTask
46from .healSparseMapping import HealSparseInputMapTask
47from lsst.meas.algorithms import SourceDetectionTask, AccumulatorMeanStack, ScaleVarianceTask
48from lsst.utils.timer import timeMethod
49from deprecated.sphinx import deprecated
50
51log = logging.getLogger(__name__)
52
53
54class AssembleCoaddConnections(pipeBase.PipelineTaskConnections,
55 dimensions=("tract", "patch", "band", "skymap"),
56 defaultTemplates={"inputCoaddName": "deep",
57 "outputCoaddName": "deep",
58 "warpType": "direct",
59 "warpTypeSuffix": ""}):
60
61 inputWarps = pipeBase.connectionTypes.Input(
62 doc=("Input list of warps to be assemebled i.e. stacked."
63 "WarpType (e.g. direct, psfMatched) is controlled by the warpType config parameter"),
64 name="{inputCoaddName}Coadd_{warpType}Warp",
65 storageClass="ExposureF",
66 dimensions=("tract", "patch", "skymap", "visit", "instrument"),
67 deferLoad=True,
68 multiple=True
69 )
70 skyMap = pipeBase.connectionTypes.Input(
71 doc="Input definition of geometry/bbox and projection/wcs for coadded exposures",
72 name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME,
73 storageClass="SkyMap",
74 dimensions=("skymap", ),
75 )
76 selectedVisits = pipeBase.connectionTypes.Input(
77 doc="Selected visits to be coadded.",
78 name="{outputCoaddName}Visits",
79 storageClass="StructuredDataDict",
80 dimensions=("instrument", "tract", "patch", "skymap", "band")
81 )
82 brightObjectMask = pipeBase.connectionTypes.PrerequisiteInput(
83 doc=("Input Bright Object Mask mask produced with external catalogs to be applied to the mask plane"
84 " BRIGHT_OBJECT."),
85 name="brightObjectMask",
86 storageClass="ObjectMaskCatalog",
87 dimensions=("tract", "patch", "skymap", "band"),
88 minimum=0,
89 )
90 coaddExposure = pipeBase.connectionTypes.Output(
91 doc="Output coadded exposure, produced by stacking input warps",
92 name="{outputCoaddName}Coadd{warpTypeSuffix}",
93 storageClass="ExposureF",
94 dimensions=("tract", "patch", "skymap", "band"),
95 )
96 nImage = pipeBase.connectionTypes.Output(
97 doc="Output image of number of input images per pixel",
98 name="{outputCoaddName}Coadd_nImage",
99 storageClass="ImageU",
100 dimensions=("tract", "patch", "skymap", "band"),
101 )
102 inputMap = pipeBase.connectionTypes.Output(
103 doc="Output healsparse map of input images",
104 name="{outputCoaddName}Coadd_inputMap",
105 storageClass="HealSparseMap",
106 dimensions=("tract", "patch", "skymap", "band"),
107 )
108
109 def __init__(self, *, config=None):
110 super().__init__(config=config)
111
112 if not config.doMaskBrightObjects:
113 self.prerequisiteInputs.remove("brightObjectMask")
114
115 if not config.doSelectVisits:
116 self.inputs.remove("selectedVisits")
117
118 if not config.doNImage:
119 self.outputs.remove("nImage")
120
121 if not self.config.doInputMap:
122 self.outputs.remove("inputMap")
123
124
125class AssembleCoaddConfig(CoaddBaseTask.ConfigClass, pipeBase.PipelineTaskConfig,
126 pipelineConnections=AssembleCoaddConnections):
127 warpType = pexConfig.Field(
128 doc="Warp name: one of 'direct' or 'psfMatched'",
129 dtype=str,
130 default="direct",
131 )
132 subregionSize = pexConfig.ListField(
133 dtype=int,
134 doc="Width, height of stack subregion size; "
135 "make small enough that a full stack of images will fit into memory at once.",
136 length=2,
137 default=(2000, 2000),
138 )
139 statistic = pexConfig.Field(
140 dtype=str,
141 doc="Main stacking statistic for aggregating over the epochs.",
142 default="MEANCLIP",
143 )
144 doOnlineForMean = pexConfig.Field(
145 dtype=bool,
146 doc="Perform online coaddition when statistic=\"MEAN\" to save memory?",
147 default=False,
148 )
149 doSigmaClip = pexConfig.Field(
150 dtype=bool,
151 doc="Perform sigma clipped outlier rejection with MEANCLIP statistic? (DEPRECATED)",
152 default=False,
153 )
154 sigmaClip = pexConfig.Field(
155 dtype=float,
156 doc="Sigma for outlier rejection; ignored if non-clipping statistic selected.",
157 default=3.0,
158 )
159 clipIter = pexConfig.Field(
160 dtype=int,
161 doc="Number of iterations of outlier rejection; ignored if non-clipping statistic selected.",
162 default=2,
163 )
164 calcErrorFromInputVariance = pexConfig.Field(
165 dtype=bool,
166 doc="Calculate coadd variance from input variance by stacking statistic."
167 "Passed to StatisticsControl.setCalcErrorFromInputVariance()",
168 default=True,
169 )
170 scaleZeroPoint = pexConfig.ConfigurableField(
171 target=ScaleZeroPointTask,
172 doc="Task to adjust the photometric zero point of the coadd temp exposures",
173 )
174 doInterp = pexConfig.Field(
175 doc="Interpolate over NaN pixels? Also extrapolate, if necessary, but the results are ugly.",
176 dtype=bool,
177 default=True,
178 )
179 interpImage = pexConfig.ConfigurableField(
180 target=InterpImageTask,
181 doc="Task to interpolate (and extrapolate) over NaN pixels",
182 )
183 doWrite = pexConfig.Field(
184 doc="Persist coadd?",
185 dtype=bool,
186 default=True,
187 )
188 doNImage = pexConfig.Field(
189 doc="Create image of number of contributing exposures for each pixel",
190 dtype=bool,
191 default=False,
192 )
193 doUsePsfMatchedPolygons = pexConfig.Field(
194 doc="Use ValidPolygons from shrunk Psf-Matched Calexps? Should be set to True by CompareWarp only.",
195 dtype=bool,
196 default=False,
197 )
198 maskPropagationThresholds = pexConfig.DictField(
199 keytype=str,
200 itemtype=float,
201 doc=("Threshold (in fractional weight) of rejection at which we propagate a mask plane to "
202 "the coadd; that is, we set the mask bit on the coadd if the fraction the rejected frames "
203 "would have contributed exceeds this value."),
204 default={"SAT": 0.1},
205 )
206 removeMaskPlanes = pexConfig.ListField(dtype=str, default=["NOT_DEBLENDED"],
207 doc="Mask planes to remove before coadding")
208 doMaskBrightObjects = pexConfig.Field(dtype=bool, default=False,
209 doc="Set mask and flag bits for bright objects?")
210 brightObjectMaskName = pexConfig.Field(dtype=str, default="BRIGHT_OBJECT",
211 doc="Name of mask bit used for bright objects")
212 coaddPsf = pexConfig.ConfigField(
213 doc="Configuration for CoaddPsf",
214 dtype=measAlg.CoaddPsfConfig,
215 )
216 doAttachTransmissionCurve = pexConfig.Field(
217 dtype=bool, default=False, optional=False,
218 doc=("Attach a piecewise TransmissionCurve for the coadd? "
219 "(requires all input Exposures to have TransmissionCurves).")
220 )
221 hasFakes = pexConfig.Field(
222 dtype=bool,
223 default=False,
224 doc="Should be set to True if fake sources have been inserted into the input data."
225 )
226 doSelectVisits = pexConfig.Field(
227 doc="Coadd only visits selected by a SelectVisitsTask",
228 dtype=bool,
229 default=False,
230 )
231 doInputMap = pexConfig.Field(
232 doc="Create a bitwise map of coadd inputs",
233 dtype=bool,
234 default=False,
235 )
236 inputMapper = pexConfig.ConfigurableField(
237 doc="Input map creation subtask.",
238 target=HealSparseInputMapTask,
239 )
240
241 def setDefaults(self):
242 super().setDefaults()
243 self.badMaskPlanes = ["NO_DATA", "BAD", "SAT", "EDGE"]
244
245 def validate(self):
246 super().validate()
247 if self.doPsfMatch:
248 # Backwards compatibility.
249 # Configs do not have loggers
250 log.warning("Config doPsfMatch deprecated. Setting warpType='psfMatched'")
251 self.warpType = 'psfMatched'
252 if self.doSigmaClip and self.statistic != "MEANCLIP":
253 log.warning('doSigmaClip deprecated. To replicate behavior, setting statistic to "MEANCLIP"')
254 self.statistic = "MEANCLIP"
255 if self.doInterp and self.statistic not in ['MEAN', 'MEDIAN', 'MEANCLIP', 'VARIANCE', 'VARIANCECLIP']:
256 raise ValueError("Must set doInterp=False for statistic=%s, which does not "
257 "compute and set a non-zero coadd variance estimate." % (self.statistic))
258
259 unstackableStats = ['NOTHING', 'ERROR', 'ORMASK']
260 if not hasattr(afwMath.Property, self.statistic) or self.statistic in unstackableStats:
261 stackableStats = [str(k) for k in afwMath.Property.__members__.keys()
262 if str(k) not in unstackableStats]
263 raise ValueError("statistic %s is not allowed. Please choose one of %s."
264 % (self.statistic, stackableStats))
265
266
267class AssembleCoaddTask(CoaddBaseTask, pipeBase.PipelineTask):
268 """Assemble a coadded image from a set of warps.
269
270 Each Warp that goes into a coadd will typically have an independent
271 photometric zero-point. Therefore, we must scale each Warp to set it to
272 a common photometric zeropoint. WarpType may be one of 'direct' or
273 'psfMatched', and the boolean configs `config.makeDirect` and
274 `config.makePsfMatched` set which of the warp types will be coadded.
275 The coadd is computed as a mean with optional outlier rejection.
276 Criteria for outlier rejection are set in `AssembleCoaddConfig`.
277 Finally, Warps can have bad 'NaN' pixels which received no input from the
278 source calExps. We interpolate over these bad (NaN) pixels.
279
280 `AssembleCoaddTask` uses several sub-tasks. These are
281
282 - `~lsst.pipe.tasks.ScaleZeroPointTask`
283 - create and use an ``imageScaler`` object to scale the photometric zeropoint for each Warp
284 - `~lsst.pipe.tasks.InterpImageTask`
285 - interpolate across bad pixels (NaN) in the final coadd
286
287 You can retarget these subtasks if you wish.
288
289 Raises
290 ------
291 RuntimeError
292 Raised if unable to define mask plane for bright objects.
293
294 Notes
295 -----
296 Debugging:
297 `AssembleCoaddTask` has no debug variables of its own. Some of the
298 subtasks may support `~lsst.base.lsstDebug` variables. See the
299 documentation for the subtasks for further information.
300
301 Examples
302 --------
303 `AssembleCoaddTask` assembles a set of warped images into a coadded image.
304 The `AssembleCoaddTask` can be invoked by running ``assembleCoadd.py``
305 with the flag '--legacyCoadd'. Usage of assembleCoadd.py expects two
306 inputs: a data reference to the tract patch and filter to be coadded, and
307 a list of Warps to attempt to coadd. These are specified using ``--id`` and
308 ``--selectId``, respectively:
309
310 .. code-block:: none
311
312 --id = [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]
313 --selectId [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]
314
315 Only the Warps that cover the specified tract and patch will be coadded.
316 A list of the available optional arguments can be obtained by calling
317 ``assembleCoadd.py`` with the ``--help`` command line argument:
318
319 .. code-block:: none
320
321 assembleCoadd.py --help
322
323 To demonstrate usage of the `AssembleCoaddTask` in the larger context of
324 multi-band processing, we will generate the HSC-I & -R band coadds from
325 HSC engineering test data provided in the ``ci_hsc`` package. To begin,
326 assuming that the lsst stack has been already set up, we must set up the
327 obs_subaru and ``ci_hsc`` packages. This defines the environment variable
328 ``$CI_HSC_DIR`` and points at the location of the package. The raw HSC
329 data live in the ``$CI_HSC_DIR/raw directory``. To begin assembling the
330 coadds, we must first run:
331
332 - processCcd
333 - process the individual ccds in $CI_HSC_RAW to produce calibrated exposures
334 - makeSkyMap
335 - create a skymap that covers the area of the sky present in the raw exposures
336 - makeCoaddTempExp
337 - warp the individual calibrated exposures to the tangent plane of the coadd
338
339 We can perform all of these steps by running
340
341 .. code-block:: none
342
343 $CI_HSC_DIR scons warp-903986 warp-904014 warp-903990 warp-904010 warp-903988
344
345 This will produce warped exposures for each visit. To coadd the warped
346 data, we call assembleCoadd.py as follows:
347
348 .. code-block:: none
349
350 assembleCoadd.py --legacyCoadd $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-I \
351 --selectId visit=903986 ccd=16 --selectId visit=903986 ccd=22 --selectId visit=903986 ccd=23 \
352 --selectId visit=903986 ccd=100 --selectId visit=904014 ccd=1 --selectId visit=904014 ccd=6 \
353 --selectId visit=904014 ccd=12 --selectId visit=903990 ccd=18 --selectId visit=903990 ccd=25 \
354 --selectId visit=904010 ccd=4 --selectId visit=904010 ccd=10 --selectId visit=904010 ccd=100 \
355 --selectId visit=903988 ccd=16 --selectId visit=903988 ccd=17 --selectId visit=903988 ccd=23 \
356 --selectId visit=903988 ccd=24
357
358 that will process the HSC-I band data. The results are written in
359 ``$CI_HSC_DIR/DATA/deepCoadd-results/HSC-I``.
360
361 You may also choose to run:
362
363 .. code-block:: none
364
365 scons warp-903334 warp-903336 warp-903338 warp-903342 warp-903344 warp-903346
366 assembleCoadd.py --legacyCoadd $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-R \
367 --selectId visit=903334 ccd=16 --selectId visit=903334 ccd=22 --selectId visit=903334 ccd=23 \
368 --selectId visit=903334 ccd=100 --selectId visit=903336 ccd=17 --selectId visit=903336 ccd=24 \
369 --selectId visit=903338 ccd=18 --selectId visit=903338 ccd=25 --selectId visit=903342 ccd=4 \
370 --selectId visit=903342 ccd=10 --selectId visit=903342 ccd=100 --selectId visit=903344 ccd=0 \
371 --selectId visit=903344 ccd=5 --selectId visit=903344 ccd=11 --selectId visit=903346 ccd=1 \
372 --selectId visit=903346 ccd=6 --selectId visit=903346 ccd=12
373
374 to generate the coadd for the HSC-R band if you are interested in
375 following multiBand Coadd processing as discussed in `pipeTasks_multiBand`
376 (but note that normally, one would use the `SafeClipAssembleCoaddTask`
377 rather than `AssembleCoaddTask` to make the coadd.
378 """
379
380 ConfigClass = AssembleCoaddConfig
381 _DefaultName = "assembleCoadd"
382
383 def __init__(self, *args, **kwargs):
384 # TODO: DM-17415 better way to handle previously allowed passed args e.g.`AssembleCoaddTask(config)`
385 if args:
386 argNames = ["config", "name", "parentTask", "log"]
387 kwargs.update({k: v for k, v in zip(argNames, args)})
388 warnings.warn("AssembleCoadd received positional args, and casting them as kwargs: %s. "
389 "PipelineTask will not take positional args" % argNames, FutureWarning,
390 stacklevel=2)
391
392 super().__init__(**kwargs)
393 self.makeSubtask("interpImage")
394 self.makeSubtask("scaleZeroPoint")
395
396 if self.config.doMaskBrightObjects:
397 mask = afwImage.Mask()
398 try:
399 self.brightObjectBitmask = 1 << mask.addMaskPlane(self.config.brightObjectMaskName)
400 except pexExceptions.LsstCppException:
401 raise RuntimeError("Unable to define mask plane for bright objects; planes used are %s" %
402 mask.getMaskPlaneDict().keys())
403 del mask
404
405 if self.config.doInputMap:
406 self.makeSubtask("inputMapper")
407
408 self.warpType = self.config.warpType
409
410 @utils.inheritDoc(pipeBase.PipelineTask)
411 def runQuantum(self, butlerQC, inputRefs, outputRefs):
412 inputData = butlerQC.get(inputRefs)
413
414 # Construct skyInfo expected by run
415 # Do not remove skyMap from inputData in case _makeSupplementaryData needs it
416 skyMap = inputData["skyMap"]
417 outputDataId = butlerQC.quantum.dataId
418
419 inputData['skyInfo'] = makeSkyInfo(skyMap,
420 tractId=outputDataId['tract'],
421 patchId=outputDataId['patch'])
422
423 if self.config.doSelectVisits:
424 warpRefList = self.filterWarps(inputData['inputWarps'], inputData['selectedVisits'])
425 else:
426 warpRefList = inputData['inputWarps']
427
428 inputs = self.prepareInputs(warpRefList)
429 self.log.info("Found %d %s", len(inputs.tempExpRefList),
431 if len(inputs.tempExpRefList) == 0:
432 raise pipeBase.NoWorkFound("No coadd temporary exposures found")
433
434 supplementaryData = self._makeSupplementaryData(butlerQC, inputRefs, outputRefs)
435 retStruct = self.run(inputData['skyInfo'], inputs.tempExpRefList, inputs.imageScalerList,
436 inputs.weightList, supplementaryData=supplementaryData)
437
438 inputData.setdefault('brightObjectMask', None)
439 if self.config.doMaskBrightObjects and inputData["brightObjectMask"] is None:
440 log.warning("doMaskBrightObjects is set to True, but brightObjectMask not loaded")
441 self.processResults(retStruct.coaddExposure, inputData['brightObjectMask'], outputDataId)
442
443 if self.config.doWrite:
444 butlerQC.put(retStruct, outputRefs)
445 return retStruct
446
447 def processResults(self, coaddExposure, brightObjectMasks=None, dataId=None):
448 """Interpolate over missing data and mask bright stars.
449
450 Parameters
451 ----------
452 coaddExposure : `lsst.afw.image.Exposure`
453 The coadded exposure to process.
454 brightObjectMasks : `lsst.afw.table` or `None`, optional
455 Table of bright objects to mask.
456 dataId : `lsst.daf.butler.DataId` or `None`, optional
457 Data identification.
458 """
459 if self.config.doInterp:
460 self.interpImage.run(coaddExposure.getMaskedImage(), planeName="NO_DATA")
461 # The variance must be positive; work around for DM-3201.
462 varArray = coaddExposure.variance.array
463 with numpy.errstate(invalid="ignore"):
464 varArray[:] = numpy.where(varArray > 0, varArray, numpy.inf)
465
466 if self.config.doMaskBrightObjects:
467 self.setBrightObjectMasks(coaddExposure, brightObjectMasks, dataId)
468
469 def _makeSupplementaryData(self, butlerQC, inputRefs, outputRefs):
470 """Make additional inputs to run() specific to subclasses (Gen3).
471
472 Duplicates interface of `runQuantum` method.
473 Available to be implemented by subclasses only if they need the
474 coadd dataRef for performing preliminary processing before
475 assembling the coadd.
476
477 Parameters
478 ----------
479 butlerQC : `~lsst.pipe.base.ButlerQuantumContext`
480 Gen3 Butler object for fetching additional data products before
481 running the Task specialized for quantum being processed.
482 inputRefs : `~lsst.pipe.base.InputQuantizedConnection`
483 Attributes are the names of the connections describing input dataset types.
484 Values are DatasetRefs that task consumes for corresponding dataset type.
485 DataIds are guaranteed to match data objects in ``inputData``.
486 outputRefs : `~lsst.pipe.base.OutputQuantizedConnection`
487 Attributes are the names of the connections describing output dataset types.
488 Values are DatasetRefs that task is to produce
489 for corresponding dataset type.
490 """
491 return pipeBase.Struct()
492
493 @deprecated(
494 reason="makeSupplementaryDataGen3 is deprecated in favor of _makeSupplementaryData",
495 version="v25.0",
496 category=FutureWarning
497 )
498 def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs):
499 return self._makeSupplementaryData(butlerQC, inputRefs, outputRefs)
500
501 def prepareInputs(self, refList):
502 """Prepare the input warps for coaddition by measuring the weight for
503 each warp and the scaling for the photometric zero point.
504
505 Each Warp has its own photometric zeropoint and background variance.
506 Before coadding these Warps together, compute a scale factor to
507 normalize the photometric zeropoint and compute the weight for each Warp.
508
509 Parameters
510 ----------
511 refList : `list`
512 List of data references to tempExp.
513
514 Returns
515 -------
516 result : `~lsst.pipe.base.Struct`
517 Results as a struct with attributes:
518
519 ``tempExprefList``
520 `list` of data references to tempExp.
521 ``weightList``
522 `list` of weightings.
523 ``imageScalerList``
524 `list` of image scalers.
525 """
526 statsCtrl = afwMath.StatisticsControl()
527 statsCtrl.setNumSigmaClip(self.config.sigmaClip)
528 statsCtrl.setNumIter(self.config.clipIter)
529 statsCtrl.setAndMask(self.getBadPixelMask())
530 statsCtrl.setNanSafe(True)
531 # compute tempExpRefList: a list of tempExpRef that actually exist
532 # and weightList: a list of the weight of the associated coadd tempExp
533 # and imageScalerList: a list of scale factors for the associated coadd tempExp
534 tempExpRefList = []
535 weightList = []
536 imageScalerList = []
537 tempExpName = self.getTempExpDatasetName(self.warpType)
538 for tempExpRef in refList:
539 tempExp = tempExpRef.get()
540 # Ignore any input warp that is empty of data
541 if numpy.isnan(tempExp.image.array).all():
542 continue
543 maskedImage = tempExp.getMaskedImage()
544 imageScaler = self.scaleZeroPoint.computeImageScaler(
545 exposure=tempExp,
546 dataRef=tempExpRef, # FIXME
547 )
548 try:
549 imageScaler.scaleMaskedImage(maskedImage)
550 except Exception as e:
551 self.log.warning("Scaling failed for %s (skipping it): %s", tempExpRef.dataId, e)
552 continue
553 statObj = afwMath.makeStatistics(maskedImage.getVariance(), maskedImage.getMask(),
554 afwMath.MEANCLIP, statsCtrl)
555 meanVar, meanVarErr = statObj.getResult(afwMath.MEANCLIP)
556 weight = 1.0 / float(meanVar)
557 if not numpy.isfinite(weight):
558 self.log.warning("Non-finite weight for %s: skipping", tempExpRef.dataId)
559 continue
560 self.log.info("Weight of %s %s = %0.3f", tempExpName, tempExpRef.dataId, weight)
561
562 del maskedImage
563 del tempExp
564
565 tempExpRefList.append(tempExpRef)
566 weightList.append(weight)
567 imageScalerList.append(imageScaler)
568
569 return pipeBase.Struct(tempExpRefList=tempExpRefList, weightList=weightList,
570 imageScalerList=imageScalerList)
571
572 def prepareStats(self, mask=None):
573 """Prepare the statistics for coadding images.
574
575 Parameters
576 ----------
577 mask : `int`, optional
578 Bit mask value to exclude from coaddition.
579
580 Returns
581 -------
582 stats : `~lsst.pipe.base.Struct`
583 Statistics as a struct with attributes:
584
585 ``statsCtrl``
586 Statistics control object for coadd (`~lsst.afw.math.StatisticsControl`).
587 ``statsFlags``
588 Statistic for coadd (`~lsst.afw.math.Property`).
589 """
590 if mask is None:
591 mask = self.getBadPixelMask()
592 statsCtrl = afwMath.StatisticsControl()
593 statsCtrl.setNumSigmaClip(self.config.sigmaClip)
594 statsCtrl.setNumIter(self.config.clipIter)
595 statsCtrl.setAndMask(mask)
596 statsCtrl.setNanSafe(True)
597 statsCtrl.setWeighted(True)
598 statsCtrl.setCalcErrorFromInputVariance(self.config.calcErrorFromInputVariance)
599 for plane, threshold in self.config.maskPropagationThresholds.items():
600 bit = afwImage.Mask.getMaskPlane(plane)
601 statsCtrl.setMaskPropagationThreshold(bit, threshold)
602 statsFlags = afwMath.stringToStatisticsProperty(self.config.statistic)
603 return pipeBase.Struct(ctrl=statsCtrl, flags=statsFlags)
604
605 @timeMethod
606 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
607 altMaskList=None, mask=None, supplementaryData=None):
608 """Assemble a coadd from input warps.
609
610 Assemble the coadd using the provided list of coaddTempExps. Since
611 the full coadd covers a patch (a large area), the assembly is
612 performed over small areas on the image at a time in order to
613 conserve memory usage. Iterate over subregions within the outer
614 bbox of the patch using `assembleSubregion` to stack the corresponding
615 subregions from the coaddTempExps with the statistic specified.
616 Set the edge bits the coadd mask based on the weight map.
617
618 Parameters
619 ----------
620 skyInfo : `~lsst.pipe.base.Struct`
621 Struct with geometric information about the patch.
622 tempExpRefList : `list`
623 List of data references to Warps (previously called CoaddTempExps).
624 imageScalerList : `list`
625 List of image scalers.
626 weightList : `list`
627 List of weights.
628 altMaskList : `list`, optional
629 List of alternate masks to use rather than those stored with
630 tempExp.
631 mask : `int`, optional
632 Bit mask value to exclude from coaddition.
633 supplementaryData : `~lsst.pipe.base.Struct`, optional
634 Struct with additional data products needed to assemble coadd.
635 Only used by subclasses that implement ``_makeSupplementaryData``
636 and override `run`.
637
638 Returns
639 -------
640 result : `~lsst.pipe.base.Struct`
641 Results as a struct with attributes:
642
643 ``coaddExposure``
644 Coadded exposure (``lsst.afw.image.Exposure``).
645 ``nImage``
646 Exposure count image (``lsst.afw.image.Image``), if requested.
647 ``inputMap``
648 Bit-wise map of inputs, if requested.
649 ``warpRefList``
650 Input list of refs to the warps (``lsst.daf.butler.DeferredDatasetHandle``)
651 (unmodified).
652 ``imageScalerList``
653 Input list of image scalers (`list`) (unmodified).
654 ``weightList``
655 Input list of weights (`list`) (unmodified).
656 """
657 tempExpName = self.getTempExpDatasetName(self.warpType)
658 self.log.info("Assembling %s %s", len(tempExpRefList), tempExpName)
659 stats = self.prepareStats(mask=mask)
660
661 if altMaskList is None:
662 altMaskList = [None]*len(tempExpRefList)
663
664 coaddExposure = afwImage.ExposureF(skyInfo.bbox, skyInfo.wcs)
665 coaddExposure.setPhotoCalib(self.scaleZeroPoint.getPhotoCalib())
666 coaddExposure.getInfo().setCoaddInputs(self.inputRecorder.makeCoaddInputs())
667 self.assembleMetadata(coaddExposure, tempExpRefList, weightList)
668 coaddMaskedImage = coaddExposure.getMaskedImage()
669 subregionSizeArr = self.config.subregionSize
670 subregionSize = geom.Extent2I(subregionSizeArr[0], subregionSizeArr[1])
671 # if nImage is requested, create a zero one which can be passed to assembleSubregion
672 if self.config.doNImage:
673 nImage = afwImage.ImageU(skyInfo.bbox)
674 else:
675 nImage = None
676 # If inputMap is requested, create the initial version that can be masked in
677 # assembleSubregion.
678 if self.config.doInputMap:
679 self.inputMapper.build_ccd_input_map(skyInfo.bbox,
680 skyInfo.wcs,
681 coaddExposure.getInfo().getCoaddInputs().ccds)
682
683 if self.config.doOnlineForMean and self.config.statistic == "MEAN":
684 try:
685 self.assembleOnlineMeanCoadd(coaddExposure, tempExpRefList, imageScalerList,
686 weightList, altMaskList, stats.ctrl,
687 nImage=nImage)
688 except Exception as e:
689 self.log.exception("Cannot compute online coadd %s", e)
690 raise
691 else:
692 for subBBox in self._subBBoxIter(skyInfo.bbox, subregionSize):
693 try:
694 self.assembleSubregion(coaddExposure, subBBox, tempExpRefList, imageScalerList,
695 weightList, altMaskList, stats.flags, stats.ctrl,
696 nImage=nImage)
697 except Exception as e:
698 self.log.exception("Cannot compute coadd %s: %s", subBBox, e)
699 raise
700
701 # If inputMap is requested, we must finalize the map after the accumulation.
702 if self.config.doInputMap:
703 self.inputMapper.finalize_ccd_input_map_mask()
704 inputMap = self.inputMapper.ccd_input_map
705 else:
706 inputMap = None
707
708 self.setInexactPsf(coaddMaskedImage.getMask())
709 # Despite the name, the following doesn't really deal with "EDGE" pixels: it identifies
710 # pixels that didn't receive any unmasked inputs (as occurs around the edge of the field).
711 coaddUtils.setCoaddEdgeBits(coaddMaskedImage.getMask(), coaddMaskedImage.getVariance())
712 return pipeBase.Struct(coaddExposure=coaddExposure, nImage=nImage,
713 warpRefList=tempExpRefList, imageScalerList=imageScalerList,
714 weightList=weightList, inputMap=inputMap)
715
716 def assembleMetadata(self, coaddExposure, tempExpRefList, weightList):
717 """Set the metadata for the coadd.
718
719 This basic implementation sets the filter from the first input.
720
721 Parameters
722 ----------
723 coaddExposure : `lsst.afw.image.Exposure`
724 The target exposure for the coadd.
725 tempExpRefList : `list`
726 List of data references to tempExp.
727 weightList : `list`
728 List of weights.
729
730 Raises
731 ------
732 AssertionError
733 Raised if there is a length mismatch.
734 """
735 assert len(tempExpRefList) == len(weightList), "Length mismatch"
736
737 # We load a single pixel of each coaddTempExp, because we just want to get at the metadata
738 # (and we need more than just the PropertySet that contains the header), which is not possible
739 # with the current butler (see #2777).
740 bbox = geom.Box2I(coaddExposure.getBBox().getMin(), geom.Extent2I(1, 1))
741
742 tempExpList = [tempExpRef.get(parameters={'bbox': bbox}) for tempExpRef in tempExpRefList]
743
744 numCcds = sum(len(tempExp.getInfo().getCoaddInputs().ccds) for tempExp in tempExpList)
745
746 # Set the coadd FilterLabel to the band of the first input exposure:
747 # Coadds are calibrated, so the physical label is now meaningless.
748 coaddExposure.setFilter(afwImage.FilterLabel(tempExpList[0].getFilter().bandLabel))
749 coaddInputs = coaddExposure.getInfo().getCoaddInputs()
750 coaddInputs.ccds.reserve(numCcds)
751 coaddInputs.visits.reserve(len(tempExpList))
752
753 for tempExp, weight in zip(tempExpList, weightList):
754 self.inputRecorder.addVisitToCoadd(coaddInputs, tempExp, weight)
755
756 if self.config.doUsePsfMatchedPolygons:
757 self.shrinkValidPolygons(coaddInputs)
758
759 coaddInputs.visits.sort()
760 coaddInputs.ccds.sort()
761 if self.warpType == "psfMatched":
762 # The modelPsf BBox for a psfMatchedWarp/coaddTempExp was dynamically defined by
763 # ModelPsfMatchTask as the square box bounding its spatially-variable, pre-matched WarpedPsf.
764 # Likewise, set the PSF of a PSF-Matched Coadd to the modelPsf
765 # having the maximum width (sufficient because square)
766 modelPsfList = [tempExp.getPsf() for tempExp in tempExpList]
767 modelPsfWidthList = [modelPsf.computeBBox(modelPsf.getAveragePosition()).getWidth()
768 for modelPsf in modelPsfList]
769 psf = modelPsfList[modelPsfWidthList.index(max(modelPsfWidthList))]
770 else:
771 psf = measAlg.CoaddPsf(coaddInputs.ccds, coaddExposure.getWcs(),
772 self.config.coaddPsf.makeControl())
773 coaddExposure.setPsf(psf)
774 apCorrMap = measAlg.makeCoaddApCorrMap(coaddInputs.ccds, coaddExposure.getBBox(afwImage.PARENT),
775 coaddExposure.getWcs())
776 coaddExposure.getInfo().setApCorrMap(apCorrMap)
777 if self.config.doAttachTransmissionCurve:
778 transmissionCurve = measAlg.makeCoaddTransmissionCurve(coaddExposure.getWcs(), coaddInputs.ccds)
779 coaddExposure.getInfo().setTransmissionCurve(transmissionCurve)
780
781 def assembleSubregion(self, coaddExposure, bbox, tempExpRefList, imageScalerList, weightList,
782 altMaskList, statsFlags, statsCtrl, nImage=None):
783 """Assemble the coadd for a sub-region.
784
785 For each coaddTempExp, check for (and swap in) an alternative mask
786 if one is passed. Remove mask planes listed in
787 `config.removeMaskPlanes`. Finally, stack the actual exposures using
788 `lsst.afw.math.statisticsStack` with the statistic specified by
789 statsFlags. Typically, the statsFlag will be one of lsst.afw.math.MEAN for
790 a mean-stack or `lsst.afw.math.MEANCLIP` for outlier rejection using
791 an N-sigma clipped mean where N and iterations are specified by
792 statsCtrl. Assign the stacked subregion back to the coadd.
793
794 Parameters
795 ----------
796 coaddExposure : `lsst.afw.image.Exposure`
797 The target exposure for the coadd.
798 bbox : `lsst.geom.Box`
799 Sub-region to coadd.
800 tempExpRefList : `list`
801 List of data reference to tempExp.
802 imageScalerList : `list`
803 List of image scalers.
804 weightList : `list`
805 List of weights.
806 altMaskList : `list`
807 List of alternate masks to use rather than those stored with
808 tempExp, or None. Each element is dict with keys = mask plane
809 name to which to add the spans.
810 statsFlags : `lsst.afw.math.Property`
811 Property object for statistic for coadd.
813 Statistics control object for coadd.
814 nImage : `lsst.afw.image.ImageU`, optional
815 Keeps track of exposure count for each pixel.
816 """
817 self.log.debug("Computing coadd over %s", bbox)
818
819 coaddExposure.mask.addMaskPlane("REJECTED")
820 coaddExposure.mask.addMaskPlane("CLIPPED")
821 coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
822 maskMap = self.setRejectedMaskMapping(statsCtrl)
823 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
824 maskedImageList = []
825 if nImage is not None:
826 subNImage = afwImage.ImageU(bbox.getWidth(), bbox.getHeight())
827 for tempExpRef, imageScaler, altMask in zip(tempExpRefList, imageScalerList, altMaskList):
828
829 exposure = tempExpRef.get(parameters={'bbox': bbox})
830
831 maskedImage = exposure.getMaskedImage()
832 mask = maskedImage.getMask()
833 if altMask is not None:
834 self.applyAltMaskPlanes(mask, altMask)
835 imageScaler.scaleMaskedImage(maskedImage)
836
837 # Add 1 for each pixel which is not excluded by the exclude mask.
838 # In legacyCoadd, pixels may also be excluded by afwMath.statisticsStack.
839 if nImage is not None:
840 subNImage.getArray()[maskedImage.getMask().getArray() & statsCtrl.getAndMask() == 0] += 1
841 if self.config.removeMaskPlanes:
842 self.removeMaskPlanes(maskedImage)
843 maskedImageList.append(maskedImage)
844
845 if self.config.doInputMap:
846 visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
847 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
848
849 with self.timer("stack"):
850 coaddSubregion = afwMath.statisticsStack(maskedImageList, statsFlags, statsCtrl, weightList,
851 clipped, # also set output to CLIPPED if sigma-clipped
852 maskMap)
853 coaddExposure.maskedImage.assign(coaddSubregion, bbox)
854 if nImage is not None:
855 nImage.assign(subNImage, bbox)
856
857 def assembleOnlineMeanCoadd(self, coaddExposure, tempExpRefList, imageScalerList, weightList,
858 altMaskList, statsCtrl, nImage=None):
859 """Assemble the coadd using the "online" method.
860
861 This method takes a running sum of images and weights to save memory.
862 It only works for MEAN statistics.
863
864 Parameters
865 ----------
866 coaddExposure : `lsst.afw.image.Exposure`
867 The target exposure for the coadd.
868 tempExpRefList : `list`
869 List of data reference to tempExp.
870 imageScalerList : `list`
871 List of image scalers.
872 weightList : `list`
873 List of weights.
874 altMaskList : `list`
875 List of alternate masks to use rather than those stored with
876 tempExp, or None. Each element is dict with keys = mask plane
877 name to which to add the spans.
879 Statistics control object for coadd.
880 nImage : `lsst.afw.image.ImageU`, optional
881 Keeps track of exposure count for each pixel.
882 """
883 self.log.debug("Computing online coadd.")
884
885 coaddExposure.mask.addMaskPlane("REJECTED")
886 coaddExposure.mask.addMaskPlane("CLIPPED")
887 coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
888 maskMap = self.setRejectedMaskMapping(statsCtrl)
889 thresholdDict = AccumulatorMeanStack.stats_ctrl_to_threshold_dict(statsCtrl)
890
891 bbox = coaddExposure.maskedImage.getBBox()
892
893 stacker = AccumulatorMeanStack(
894 coaddExposure.image.array.shape,
895 statsCtrl.getAndMask(),
896 mask_threshold_dict=thresholdDict,
897 mask_map=maskMap,
898 no_good_pixels_mask=statsCtrl.getNoGoodPixelsMask(),
899 calc_error_from_input_variance=self.config.calcErrorFromInputVariance,
900 compute_n_image=(nImage is not None)
901 )
902
903 for tempExpRef, imageScaler, altMask, weight in zip(tempExpRefList,
904 imageScalerList,
905 altMaskList,
906 weightList):
907 exposure = tempExpRef.get()
908 maskedImage = exposure.getMaskedImage()
909 mask = maskedImage.getMask()
910 if altMask is not None:
911 self.applyAltMaskPlanes(mask, altMask)
912 imageScaler.scaleMaskedImage(maskedImage)
913 if self.config.removeMaskPlanes:
914 self.removeMaskPlanes(maskedImage)
915
916 stacker.add_masked_image(maskedImage, weight=weight)
917
918 if self.config.doInputMap:
919 visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
920 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
921
922 stacker.fill_stacked_masked_image(coaddExposure.maskedImage)
923
924 if nImage is not None:
925 nImage.array[:, :] = stacker.n_image
926
927 def removeMaskPlanes(self, maskedImage):
928 """Unset the mask of an image for mask planes specified in the config.
929
930 Parameters
931 ----------
932 maskedImage : `lsst.afw.image.MaskedImage`
933 The masked image to be modified.
934
935 Raises
936 ------
937 InvalidParameterError
938 Raised if no mask plane with that name was found.
939 """
940 mask = maskedImage.getMask()
941 for maskPlane in self.config.removeMaskPlanes:
942 try:
943 mask &= ~mask.getPlaneBitMask(maskPlane)
944 except pexExceptions.InvalidParameterError:
945 self.log.debug("Unable to remove mask plane %s: no mask plane with that name was found.",
946 maskPlane)
947
948 @staticmethod
950 """Map certain mask planes of the warps to new planes for the coadd.
951
952 If a pixel is rejected due to a mask value other than EDGE, NO_DATA,
953 or CLIPPED, set it to REJECTED on the coadd.
954 If a pixel is rejected due to EDGE, set the coadd pixel to SENSOR_EDGE.
955 If a pixel is rejected due to CLIPPED, set the coadd pixel to CLIPPED.
956
957 Parameters
958 ----------
960 Statistics control object for coadd.
961
962 Returns
963 -------
964 maskMap : `list` of `tuple` of `int`
965 A list of mappings of mask planes of the warped exposures to
966 mask planes of the coadd.
967 """
968 edge = afwImage.Mask.getPlaneBitMask("EDGE")
969 noData = afwImage.Mask.getPlaneBitMask("NO_DATA")
970 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
971 toReject = statsCtrl.getAndMask() & (~noData) & (~edge) & (~clipped)
972 maskMap = [(toReject, afwImage.Mask.getPlaneBitMask("REJECTED")),
973 (edge, afwImage.Mask.getPlaneBitMask("SENSOR_EDGE")),
974 (clipped, clipped)]
975 return maskMap
976
977 def applyAltMaskPlanes(self, mask, altMaskSpans):
978 """Apply in place alt mask formatted as SpanSets to a mask.
979
980 Parameters
981 ----------
982 mask : `lsst.afw.image.Mask`
983 Original mask.
984 altMaskSpans : `dict`
985 SpanSet lists to apply. Each element contains the new mask
986 plane name (e.g. "CLIPPED and/or "NO_DATA") as the key,
987 and list of SpanSets to apply to the mask.
988
989 Returns
990 -------
991 mask : `lsst.afw.image.Mask`
992 Updated mask.
993 """
994 if self.config.doUsePsfMatchedPolygons:
995 if ("NO_DATA" in altMaskSpans) and ("NO_DATA" in self.config.badMaskPlanes):
996 # Clear away any other masks outside the validPolygons. These pixels are no longer
997 # contributing to inexact PSFs, and will still be rejected because of NO_DATA
998 # self.config.doUsePsfMatchedPolygons should be True only in CompareWarpAssemble
999 # This mask-clearing step must only occur *before* applying the new masks below
1000 for spanSet in altMaskSpans['NO_DATA']:
1001 spanSet.clippedTo(mask.getBBox()).clearMask(mask, self.getBadPixelMask())
1002
1003 for plane, spanSetList in altMaskSpans.items():
1004 maskClipValue = mask.addMaskPlane(plane)
1005 for spanSet in spanSetList:
1006 spanSet.clippedTo(mask.getBBox()).setMask(mask, 2**maskClipValue)
1007 return mask
1008
1009 def shrinkValidPolygons(self, coaddInputs):
1010 """Shrink coaddInputs' ccds' ValidPolygons in place.
1011
1012 Either modify each ccd's validPolygon in place, or if CoaddInputs
1013 does not have a validPolygon, create one from its bbox.
1014
1015 Parameters
1016 ----------
1017 coaddInputs : `lsst.afw.image.coaddInputs`
1018 Original mask.
1019 """
1020 for ccd in coaddInputs.ccds:
1021 polyOrig = ccd.getValidPolygon()
1022 validPolyBBox = polyOrig.getBBox() if polyOrig else ccd.getBBox()
1023 validPolyBBox.grow(-self.config.matchingKernelSize//2)
1024 if polyOrig:
1025 validPolygon = polyOrig.intersectionSingle(validPolyBBox)
1026 else:
1027 validPolygon = afwGeom.polygon.Polygon(geom.Box2D(validPolyBBox))
1028 ccd.setValidPolygon(validPolygon)
1029
1030 def setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None):
1031 """Set the bright object masks.
1032
1033 Parameters
1034 ----------
1035 exposure : `lsst.afw.image.Exposure`
1036 Exposure under consideration.
1037 brightObjectMasks : `lsst.afw.table`
1038 Table of bright objects to mask.
1039 dataId : `lsst.daf.butler.DataId`, optional
1040 Data identifier dict for patch.
1041 """
1042 if brightObjectMasks is None:
1043 self.log.warning("Unable to apply bright object mask: none supplied")
1044 return
1045 self.log.info("Applying %d bright object masks to %s", len(brightObjectMasks), dataId)
1046 mask = exposure.getMaskedImage().getMask()
1047 wcs = exposure.getWcs()
1048 plateScale = wcs.getPixelScale().asArcseconds()
1049
1050 for rec in brightObjectMasks:
1051 center = geom.PointI(wcs.skyToPixel(rec.getCoord()))
1052 if rec["type"] == "box":
1053 assert rec["angle"] == 0.0, ("Angle != 0 for mask object %s" % rec["id"])
1054 width = rec["width"].asArcseconds()/plateScale # convert to pixels
1055 height = rec["height"].asArcseconds()/plateScale # convert to pixels
1056
1057 halfSize = geom.ExtentI(0.5*width, 0.5*height)
1058 bbox = geom.Box2I(center - halfSize, center + halfSize)
1059
1060 bbox = geom.BoxI(geom.PointI(int(center[0] - 0.5*width), int(center[1] - 0.5*height)),
1061 geom.PointI(int(center[0] + 0.5*width), int(center[1] + 0.5*height)))
1062 spans = afwGeom.SpanSet(bbox)
1063 elif rec["type"] == "circle":
1064 radius = int(rec["radius"].asArcseconds()/plateScale) # convert to pixels
1065 spans = afwGeom.SpanSet.fromShape(radius, offset=center)
1066 else:
1067 self.log.warning("Unexpected region type %s at %s", rec["type"], center)
1068 continue
1069 spans.clippedTo(mask.getBBox()).setMask(mask, self.brightObjectBitmask)
1070
1071 def setInexactPsf(self, mask):
1072 """Set INEXACT_PSF mask plane.
1073
1074 If any of the input images isn't represented in the coadd (due to
1075 clipped pixels or chip gaps), the `CoaddPsf` will be inexact. Flag
1076 these pixels.
1077
1078 Parameters
1079 ----------
1080 mask : `lsst.afw.image.Mask`
1081 Coadded exposure's mask, modified in-place.
1082 """
1083 mask.addMaskPlane("INEXACT_PSF")
1084 inexactPsf = mask.getPlaneBitMask("INEXACT_PSF")
1085 sensorEdge = mask.getPlaneBitMask("SENSOR_EDGE") # chip edges (so PSF is discontinuous)
1086 clipped = mask.getPlaneBitMask("CLIPPED") # pixels clipped from coadd
1087 rejected = mask.getPlaneBitMask("REJECTED") # pixels rejected from coadd due to masks
1088 array = mask.getArray()
1089 selected = array & (sensorEdge | clipped | rejected) > 0
1090 array[selected] |= inexactPsf
1091
1092 @staticmethod
1093 def _subBBoxIter(bbox, subregionSize):
1094 """Iterate over subregions of a bbox.
1095
1096 Parameters
1097 ----------
1098 bbox : `lsst.geom.Box2I`
1099 Bounding box over which to iterate.
1100 subregionSize : `lsst.geom.Extent2I`
1101 Size of sub-bboxes.
1102
1103 Yields
1104 ------
1105 subBBox : `lsst.geom.Box2I`
1106 Next sub-bounding box of size ``subregionSize`` or smaller; each ``subBBox``
1107 is contained within ``bbox``, so it may be smaller than ``subregionSize`` at
1108 the edges of ``bbox``, but it will never be empty.
1109
1110 Raises
1111 ------
1112 RuntimeError
1113 Raised if any of the following occur:
1114 - The given bbox is empty.
1115 - The subregionSize is 0.
1116 """
1117 if bbox.isEmpty():
1118 raise RuntimeError("bbox %s is empty" % (bbox,))
1119 if subregionSize[0] < 1 or subregionSize[1] < 1:
1120 raise RuntimeError("subregionSize %s must be nonzero" % (subregionSize,))
1121
1122 for rowShift in range(0, bbox.getHeight(), subregionSize[1]):
1123 for colShift in range(0, bbox.getWidth(), subregionSize[0]):
1124 subBBox = geom.Box2I(bbox.getMin() + geom.Extent2I(colShift, rowShift), subregionSize)
1125 subBBox.clip(bbox)
1126 if subBBox.isEmpty():
1127 raise RuntimeError("Bug: empty bbox! bbox=%s, subregionSize=%s, "
1128 "colShift=%s, rowShift=%s" %
1129 (bbox, subregionSize, colShift, rowShift))
1130 yield subBBox
1131
1132 def filterWarps(self, inputs, goodVisits):
1133 """Return list of only inputRefs with visitId in goodVisits ordered by goodVisit.
1134
1135 Parameters
1136 ----------
1137 inputs : `list` of `~lsst.pipe.base.connections.DeferredDatasetRef`
1138 List of `lsst.pipe.base.connections.DeferredDatasetRef` with dataId containing visit.
1139 goodVisit : `dict`
1140 Dictionary with good visitIds as the keys. Value ignored.
1141
1142 Returns
1143 -------
1144 filteredInputs : `list` of `~lsst.pipe.base.connections.DeferredDatasetRef`
1145 Filtered and sorted list of inputRefs with visitId in goodVisits ordered by goodVisit.
1146 """
1147 inputWarpDict = {inputRef.ref.dataId['visit']: inputRef for inputRef in inputs}
1148 filteredInputs = []
1149 for visit in goodVisits.keys():
1150 if visit in inputWarpDict:
1151 filteredInputs.append(inputWarpDict[visit])
1152 return filteredInputs
1153
1154
1155def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask):
1156 """Function to count the number of pixels with a specific mask in a
1157 footprint.
1158
1159 Find the intersection of mask & footprint. Count all pixels in the mask
1160 that are in the intersection that have bitmask set but do not have
1161 ignoreMask set. Return the count.
1162
1163 Parameters
1164 ----------
1165 mask : `lsst.afw.image.Mask`
1166 Mask to define intersection region by.
1167 footprint : `lsst.afw.detection.Footprint`
1168 Footprint to define the intersection region by.
1169 bitmask : `Unknown`
1170 Specific mask that we wish to count the number of occurances of.
1171 ignoreMask : `Unknown`
1172 Pixels to not consider.
1173
1174 Returns
1175 -------
1176 result : `int`
1177 Number of pixels in footprint with specified mask.
1178 """
1179 bbox = footprint.getBBox()
1180 bbox.clip(mask.getBBox(afwImage.PARENT))
1181 fp = afwImage.Mask(bbox)
1182 subMask = mask.Factory(mask, bbox, afwImage.PARENT)
1183 footprint.spans.setMask(fp, bitmask)
1184 return numpy.logical_and((subMask.getArray() & fp.getArray()) > 0,
1185 (subMask.getArray() & ignoreMask) == 0).sum()
1186
1187
1189 psfMatchedWarps = pipeBase.connectionTypes.Input(
1190 doc=("PSF-Matched Warps are required by CompareWarp regardless of the coadd type requested. "
1191 "Only PSF-Matched Warps make sense for image subtraction. "
1192 "Therefore, they must be an additional declared input."),
1193 name="{inputCoaddName}Coadd_psfMatchedWarp",
1194 storageClass="ExposureF",
1195 dimensions=("tract", "patch", "skymap", "visit"),
1196 deferLoad=True,
1197 multiple=True
1198 )
1199 templateCoadd = pipeBase.connectionTypes.Output(
1200 doc=("Model of the static sky, used to find temporal artifacts. Typically a PSF-Matched, "
1201 "sigma-clipped coadd. Written if and only if assembleStaticSkyModel.doWrite=True"),
1202 name="{outputCoaddName}CoaddPsfMatched",
1203 storageClass="ExposureF",
1204 dimensions=("tract", "patch", "skymap", "band"),
1205 )
1206
1207 def __init__(self, *, config=None):
1208 super().__init__(config=config)
1209 if not config.assembleStaticSkyModel.doWrite:
1210 self.outputs.remove("templateCoadd")
1211 config.validate()
1212
1213
1214class CompareWarpAssembleCoaddConfig(AssembleCoaddConfig,
1215 pipelineConnections=CompareWarpAssembleCoaddConnections):
1216 assembleStaticSkyModel = pexConfig.ConfigurableField(
1217 target=AssembleCoaddTask,
1218 doc="Task to assemble an artifact-free, PSF-matched Coadd to serve as a"
1219 " naive/first-iteration model of the static sky.",
1220 )
1221 detect = pexConfig.ConfigurableField(
1222 target=SourceDetectionTask,
1223 doc="Detect outlier sources on difference between each psfMatched warp and static sky model"
1224 )
1225 detectTemplate = pexConfig.ConfigurableField(
1226 target=SourceDetectionTask,
1227 doc="Detect sources on static sky model. Only used if doPreserveContainedBySource is True"
1228 )
1229 maskStreaks = pexConfig.ConfigurableField(
1230 target=MaskStreaksTask,
1231 doc="Detect streaks on difference between each psfMatched warp and static sky model. Only used if "
1232 "doFilterMorphological is True. Adds a mask plane to an exposure, with the mask plane name set by"
1233 "streakMaskName"
1234 )
1235 streakMaskName = pexConfig.Field(
1236 dtype=str,
1237 default="STREAK",
1238 doc="Name of mask bit used for streaks"
1239 )
1240 maxNumEpochs = pexConfig.Field(
1241 doc="Charactistic maximum local number of epochs/visits in which an artifact candidate can appear "
1242 "and still be masked. The effective maxNumEpochs is a broken linear function of local "
1243 "number of epochs (N): min(maxFractionEpochsLow*N, maxNumEpochs + maxFractionEpochsHigh*N). "
1244 "For each footprint detected on the image difference between the psfMatched warp and static sky "
1245 "model, if a significant fraction of pixels (defined by spatialThreshold) are residuals in more "
1246 "than the computed effective maxNumEpochs, the artifact candidate is deemed persistant rather "
1247 "than transient and not masked.",
1248 dtype=int,
1249 default=2
1250 )
1251 maxFractionEpochsLow = pexConfig.RangeField(
1252 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for low N. "
1253 "Effective maxNumEpochs = "
1254 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1255 dtype=float,
1256 default=0.4,
1257 min=0., max=1.,
1258 )
1259 maxFractionEpochsHigh = pexConfig.RangeField(
1260 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for high N. "
1261 "Effective maxNumEpochs = "
1262 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1263 dtype=float,
1264 default=0.03,
1265 min=0., max=1.,
1266 )
1267 spatialThreshold = pexConfig.RangeField(
1268 doc="Unitless fraction of pixels defining how much of the outlier region has to meet the "
1269 "temporal criteria. If 0, clip all. If 1, clip none.",
1270 dtype=float,
1271 default=0.5,
1272 min=0., max=1.,
1273 inclusiveMin=True, inclusiveMax=True
1274 )
1275 doScaleWarpVariance = pexConfig.Field(
1276 doc="Rescale Warp variance plane using empirical noise?",
1277 dtype=bool,
1278 default=True,
1279 )
1280 scaleWarpVariance = pexConfig.ConfigurableField(
1281 target=ScaleVarianceTask,
1282 doc="Rescale variance on warps",
1283 )
1284 doPreserveContainedBySource = pexConfig.Field(
1285 doc="Rescue artifacts from clipping that completely lie within a footprint detected"
1286 "on the PsfMatched Template Coadd. Replicates a behavior of SafeClip.",
1287 dtype=bool,
1288 default=True,
1289 )
1290 doPrefilterArtifacts = pexConfig.Field(
1291 doc="Ignore artifact candidates that are mostly covered by the bad pixel mask, "
1292 "because they will be excluded anyway. This prevents them from contributing "
1293 "to the outlier epoch count image and potentially being labeled as persistant."
1294 "'Mostly' is defined by the config 'prefilterArtifactsRatio'.",
1295 dtype=bool,
1296 default=True
1297 )
1298 prefilterArtifactsMaskPlanes = pexConfig.ListField(
1299 doc="Prefilter artifact candidates that are mostly covered by these bad mask planes.",
1300 dtype=str,
1301 default=('NO_DATA', 'BAD', 'SAT', 'SUSPECT'),
1302 )
1303 prefilterArtifactsRatio = pexConfig.Field(
1304 doc="Prefilter artifact candidates with less than this fraction overlapping good pixels",
1305 dtype=float,
1306 default=0.05
1307 )
1308 doFilterMorphological = pexConfig.Field(
1309 doc="Filter artifact candidates based on morphological criteria, i.g. those that appear to "
1310 "be streaks.",
1311 dtype=bool,
1312 default=False
1313 )
1314 growStreakFp = pexConfig.Field(
1315 doc="Grow streak footprints by this number multiplied by the PSF width",
1316 dtype=float,
1317 default=5
1318 )
1319
1320 def setDefaults(self):
1321 AssembleCoaddConfig.setDefaults(self)
1322 self.statistic = 'MEAN'
1324
1325 # Real EDGE removed by psfMatched NO_DATA border half the width of the matching kernel
1326 # CompareWarp applies psfMatched EDGE pixels to directWarps before assembling
1327 if "EDGE" in self.badMaskPlanes:
1328 self.badMaskPlanes.remove('EDGE')
1329 self.removeMaskPlanes.append('EDGE')
1330 self.assembleStaticSkyModel.badMaskPlanes = ["NO_DATA", ]
1331 self.assembleStaticSkyModel.warpType = 'psfMatched'
1332 self.assembleStaticSkyModel.connections.warpType = 'psfMatched'
1333 self.assembleStaticSkyModel.statistic = 'MEANCLIP'
1334 self.assembleStaticSkyModel.sigmaClip = 2.5
1335 self.assembleStaticSkyModel.clipIter = 3
1336 self.assembleStaticSkyModel.calcErrorFromInputVariance = False
1337 self.assembleStaticSkyModel.doWrite = False
1338 self.detect.doTempLocalBackground = False
1339 self.detect.reEstimateBackground = False
1340 self.detect.returnOriginalFootprints = False
1341 self.detect.thresholdPolarity = "both"
1342 self.detect.thresholdValue = 5
1343 self.detect.minPixels = 4
1344 self.detect.isotropicGrow = True
1345 self.detect.thresholdType = "pixel_stdev"
1346 self.detect.nSigmaToGrow = 0.4
1347 # The default nSigmaToGrow for SourceDetectionTask is already 2.4,
1348 # Explicitly restating because ratio with detect.nSigmaToGrow matters
1349 self.detectTemplate.nSigmaToGrow = 2.4
1350 self.detectTemplate.doTempLocalBackground = False
1351 self.detectTemplate.reEstimateBackground = False
1352 self.detectTemplate.returnOriginalFootprints = False
1353
1354 def validate(self):
1355 super().validate()
1356 if self.assembleStaticSkyModel.doNImage:
1357 raise ValueError("No dataset type exists for a PSF-Matched Template N Image."
1358 "Please set assembleStaticSkyModel.doNImage=False")
1359
1360 if self.assembleStaticSkyModel.doWrite and (self.warpType == self.assembleStaticSkyModel.warpType):
1361 raise ValueError("warpType (%s) == assembleStaticSkyModel.warpType (%s) and will compete for "
1362 "the same dataset name. Please set assembleStaticSkyModel.doWrite to False "
1363 "or warpType to 'direct'. assembleStaticSkyModel.warpType should ways be "
1364 "'PsfMatched'" % (self.warpType, self.assembleStaticSkyModel.warpType))
1365
1366
1368 """Assemble a compareWarp coadded image from a set of warps
1369 by masking artifacts detected by comparing PSF-matched warps.
1370
1371 In ``AssembleCoaddTask``, we compute the coadd as an clipped mean (i.e.,
1372 we clip outliers). The problem with doing this is that when computing the
1373 coadd PSF at a given location, individual visit PSFs from visits with
1374 outlier pixels contribute to the coadd PSF and cannot be treated correctly.
1375 In this task, we correct for this behavior by creating a new badMaskPlane
1376 'CLIPPED' which marks pixels in the individual warps suspected to contain
1377 an artifact. We populate this plane on the input warps by comparing
1378 PSF-matched warps with a PSF-matched median coadd which serves as a
1379 model of the static sky. Any group of pixels that deviates from the
1380 PSF-matched template coadd by more than config.detect.threshold sigma,
1381 is an artifact candidate. The candidates are then filtered to remove
1382 variable sources and sources that are difficult to subtract such as
1383 bright stars. This filter is configured using the config parameters
1384 ``temporalThreshold`` and ``spatialThreshold``. The temporalThreshold is
1385 the maximum fraction of epochs that the deviation can appear in and still
1386 be considered an artifact. The spatialThreshold is the maximum fraction of
1387 pixels in the footprint of the deviation that appear in other epochs
1388 (where other epochs is defined by the temporalThreshold). If the deviant
1389 region meets this criteria of having a significant percentage of pixels
1390 that deviate in only a few epochs, these pixels have the 'CLIPPED' bit
1391 set in the mask. These regions will not contribute to the final coadd.
1392 Furthermore, any routine to determine the coadd PSF can now be cognizant
1393 of clipped regions. Note that the algorithm implemented by this task is
1394 preliminary and works correctly for HSC data. Parameter modifications and
1395 or considerable redesigning of the algorithm is likley required for other
1396 surveys.
1397
1398 ``CompareWarpAssembleCoaddTask`` sub-classes
1399 ``AssembleCoaddTask`` and instantiates ``AssembleCoaddTask``
1400 as a subtask to generate the TemplateCoadd (the model of the static sky).
1401
1402 Notes
1403 -----
1404 Debugging:
1405 This task supports the following debug variables:
1406 - ``saveCountIm``
1407 If True then save the Epoch Count Image as a fits file in the `figPath`
1408 - ``figPath``
1409 Path to save the debug fits images and figures
1410 """
1411
1412 ConfigClass = CompareWarpAssembleCoaddConfig
1413 _DefaultName = "compareWarpAssembleCoadd"
1414
1415 def __init__(self, *args, **kwargs):
1416 AssembleCoaddTask.__init__(self, *args, **kwargs)
1417 self.makeSubtask("assembleStaticSkyModel")
1418 detectionSchema = afwTable.SourceTable.makeMinimalSchema()
1419 self.makeSubtask("detect", schema=detectionSchema)
1420 if self.config.doPreserveContainedBySource:
1421 self.makeSubtask("detectTemplate", schema=afwTable.SourceTable.makeMinimalSchema())
1422 if self.config.doScaleWarpVariance:
1423 self.makeSubtask("scaleWarpVariance")
1424 if self.config.doFilterMorphological:
1425 self.makeSubtask("maskStreaks")
1426
1427 @utils.inheritDoc(AssembleCoaddTask)
1428 def _makeSupplementaryData(self, butlerQC, inputRefs, outputRefs):
1429 """Generate a templateCoadd to use as a naive model of static sky to
1430 subtract from PSF-Matched warps.
1431
1432 Returns
1433 -------
1434 result : `~lsst.pipe.base.Struct`
1435 Results as a struct with attributes:
1436
1437 ``templateCoadd``
1438 Coadded exposure (`lsst.afw.image.Exposure`).
1439 ``nImage``
1440 Keeps track of exposure count for each pixel (`lsst.afw.image.ImageU`).
1441
1442 Raises
1443 ------
1444 RuntimeError
1445 Raised if ``templateCoadd`` is `None`.
1446 """
1447 # Ensure that psfMatchedWarps are used as input warps for template generation
1448 staticSkyModelInputRefs = copy.deepcopy(inputRefs)
1449 staticSkyModelInputRefs.inputWarps = inputRefs.psfMatchedWarps
1450
1451 # Because subtasks don't have connections we have to make one.
1452 # The main task's `templateCoadd` is the subtask's `coaddExposure`
1453 staticSkyModelOutputRefs = copy.deepcopy(outputRefs)
1454 if self.config.assembleStaticSkyModel.doWrite:
1455 staticSkyModelOutputRefs.coaddExposure = staticSkyModelOutputRefs.templateCoadd
1456 # Remove template coadd from both subtask's and main tasks outputs,
1457 # because it is handled by the subtask as `coaddExposure`
1458 del outputRefs.templateCoadd
1459 del staticSkyModelOutputRefs.templateCoadd
1460
1461 # A PSF-Matched nImage does not exist as a dataset type
1462 if 'nImage' in staticSkyModelOutputRefs.keys():
1463 del staticSkyModelOutputRefs.nImage
1464
1465 templateCoadd = self.assembleStaticSkyModel.runQuantum(butlerQC, staticSkyModelInputRefs,
1466 staticSkyModelOutputRefs)
1467 if templateCoadd is None:
1468 raise RuntimeError(self._noTemplateMessage(self.assembleStaticSkyModel.warpType))
1469
1470 return pipeBase.Struct(templateCoadd=templateCoadd.coaddExposure,
1471 nImage=templateCoadd.nImage,
1472 warpRefList=templateCoadd.warpRefList,
1473 imageScalerList=templateCoadd.imageScalerList,
1474 weightList=templateCoadd.weightList)
1475
1476 def _noTemplateMessage(self, warpType):
1477 warpName = (warpType[0].upper() + warpType[1:])
1478 message = """No %(warpName)s warps were found to build the template coadd which is
1479 required to run CompareWarpAssembleCoaddTask. To continue assembling this type of coadd,
1480 first either rerun makeCoaddTempExp with config.make%(warpName)s=True or
1481 coaddDriver with config.makeCoadTempExp.make%(warpName)s=True, before assembleCoadd.
1482
1483 Alternatively, to use another algorithm with existing warps, retarget the CoaddDriverConfig to
1484 another algorithm like:
1485
1486 from lsst.pipe.tasks.assembleCoadd import SafeClipAssembleCoaddTask
1487 config.assemble.retarget(SafeClipAssembleCoaddTask)
1488 """ % {"warpName": warpName}
1489 return message
1490
1491 @utils.inheritDoc(AssembleCoaddTask)
1492 @timeMethod
1493 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
1494 supplementaryData):
1495 """Assemble the coadd.
1496
1497 Find artifacts and apply them to the warps' masks creating a list of
1498 alternative masks with a new "CLIPPED" plane and updated "NO_DATA"
1499 plane. Then pass these alternative masks to the base class's ``run``
1500 method.
1501 """
1502 # Check and match the order of the supplementaryData
1503 # (PSF-matched) inputs to the order of the direct inputs,
1504 # so that the artifact mask is applied to the right warp
1505 dataIds = [ref.dataId for ref in tempExpRefList]
1506 psfMatchedDataIds = [ref.dataId for ref in supplementaryData.warpRefList]
1507
1508 if dataIds != psfMatchedDataIds:
1509 self.log.info("Reordering and or/padding PSF-matched visit input list")
1510 supplementaryData.warpRefList = reorderAndPadList(supplementaryData.warpRefList,
1511 psfMatchedDataIds, dataIds)
1512 supplementaryData.imageScalerList = reorderAndPadList(supplementaryData.imageScalerList,
1513 psfMatchedDataIds, dataIds)
1514
1515 # Use PSF-Matched Warps (and corresponding scalers) and coadd to find artifacts
1516 spanSetMaskList = self.findArtifacts(supplementaryData.templateCoadd,
1517 supplementaryData.warpRefList,
1518 supplementaryData.imageScalerList)
1519
1520 badMaskPlanes = self.config.badMaskPlanes[:]
1521 badMaskPlanes.append("CLIPPED")
1522 badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes)
1523
1524 result = AssembleCoaddTask.run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
1525 spanSetMaskList, mask=badPixelMask)
1526
1527 # Propagate PSF-matched EDGE pixels to coadd SENSOR_EDGE and INEXACT_PSF
1528 # Psf-Matching moves the real edge inwards
1529 self.applyAltEdgeMask(result.coaddExposure.maskedImage.mask, spanSetMaskList)
1530 return result
1531
1532 def applyAltEdgeMask(self, mask, altMaskList):
1533 """Propagate alt EDGE mask to SENSOR_EDGE AND INEXACT_PSF planes.
1534
1535 Parameters
1536 ----------
1537 mask : `lsst.afw.image.Mask`
1538 Original mask.
1539 altMaskList : `list` of `dict`
1540 List of Dicts containing ``spanSet`` lists.
1541 Each element contains the new mask plane name (e.g. "CLIPPED
1542 and/or "NO_DATA") as the key, and list of ``SpanSets`` to apply to
1543 the mask.
1544 """
1545 maskValue = mask.getPlaneBitMask(["SENSOR_EDGE", "INEXACT_PSF"])
1546 for visitMask in altMaskList:
1547 if "EDGE" in visitMask:
1548 for spanSet in visitMask['EDGE']:
1549 spanSet.clippedTo(mask.getBBox()).setMask(mask, maskValue)
1550
1551 def findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList):
1552 """Find artifacts.
1553
1554 Loop through warps twice. The first loop builds a map with the count
1555 of how many epochs each pixel deviates from the templateCoadd by more
1556 than ``config.chiThreshold`` sigma. The second loop takes each
1557 difference image and filters the artifacts detected in each using
1558 count map to filter out variable sources and sources that are
1559 difficult to subtract cleanly.
1560
1561 Parameters
1562 ----------
1563 templateCoadd : `lsst.afw.image.Exposure`
1564 Exposure to serve as model of static sky.
1565 tempExpRefList : `list`
1566 List of data references to warps.
1567 imageScalerList : `list`
1568 List of image scalers.
1569
1570 Returns
1571 -------
1572 altMasks : `list` of `dict`
1573 List of dicts containing information about CLIPPED
1574 (i.e., artifacts), NO_DATA, and EDGE pixels.
1575 """
1576 self.log.debug("Generating Count Image, and mask lists.")
1577 coaddBBox = templateCoadd.getBBox()
1578 slateIm = afwImage.ImageU(coaddBBox)
1579 epochCountImage = afwImage.ImageU(coaddBBox)
1580 nImage = afwImage.ImageU(coaddBBox)
1581 spanSetArtifactList = []
1582 spanSetNoDataMaskList = []
1583 spanSetEdgeList = []
1584 spanSetBadMorphoList = []
1585 badPixelMask = self.getBadPixelMask()
1586
1587 # mask of the warp diffs should = that of only the warp
1588 templateCoadd.mask.clearAllMaskPlanes()
1589
1590 if self.config.doPreserveContainedBySource:
1591 templateFootprints = self.detectTemplate.detectFootprints(templateCoadd)
1592 else:
1593 templateFootprints = None
1594
1595 for warpRef, imageScaler in zip(tempExpRefList, imageScalerList):
1596 warpDiffExp = self._readAndComputeWarpDiff(warpRef, imageScaler, templateCoadd)
1597 if warpDiffExp is not None:
1598 # This nImage only approximates the final nImage because it uses the PSF-matched mask
1599 nImage.array += (numpy.isfinite(warpDiffExp.image.array)
1600 * ((warpDiffExp.mask.array & badPixelMask) == 0)).astype(numpy.uint16)
1601 fpSet = self.detect.detectFootprints(warpDiffExp, doSmooth=False, clearMask=True)
1602 fpSet.positive.merge(fpSet.negative)
1603 footprints = fpSet.positive
1604 slateIm.set(0)
1605 spanSetList = [footprint.spans for footprint in footprints.getFootprints()]
1606
1607 # Remove artifacts due to defects before they contribute to the epochCountImage
1608 if self.config.doPrefilterArtifacts:
1609 spanSetList = self.prefilterArtifacts(spanSetList, warpDiffExp)
1610
1611 # Clear mask before adding prefiltered spanSets
1612 self.detect.clearMask(warpDiffExp.mask)
1613 for spans in spanSetList:
1614 spans.setImage(slateIm, 1, doClip=True)
1615 spans.setMask(warpDiffExp.mask, warpDiffExp.mask.getPlaneBitMask("DETECTED"))
1616 epochCountImage += slateIm
1617
1618 if self.config.doFilterMorphological:
1619 maskName = self.config.streakMaskName
1620 _ = self.maskStreaks.run(warpDiffExp)
1621 streakMask = warpDiffExp.mask
1622 spanSetStreak = afwGeom.SpanSet.fromMask(streakMask,
1623 streakMask.getPlaneBitMask(maskName)).split()
1624 # Pad the streaks to account for low-surface brightness wings
1625 psf = warpDiffExp.getPsf()
1626 for s, sset in enumerate(spanSetStreak):
1627 psfShape = psf.computeShape(sset.computeCentroid())
1628 dilation = self.config.growStreakFp * psfShape.getDeterminantRadius()
1629 sset_dilated = sset.dilated(int(dilation))
1630 spanSetStreak[s] = sset_dilated
1631
1632 # PSF-Matched warps have less available area (~the matching kernel) because the calexps
1633 # undergo a second convolution. Pixels with data in the direct warp
1634 # but not in the PSF-matched warp will not have their artifacts detected.
1635 # NaNs from the PSF-matched warp therefore must be masked in the direct warp
1636 nans = numpy.where(numpy.isnan(warpDiffExp.maskedImage.image.array), 1, 0)
1637 nansMask = afwImage.makeMaskFromArray(nans.astype(afwImage.MaskPixel))
1638 nansMask.setXY0(warpDiffExp.getXY0())
1639 edgeMask = warpDiffExp.mask
1640 spanSetEdgeMask = afwGeom.SpanSet.fromMask(edgeMask,
1641 edgeMask.getPlaneBitMask("EDGE")).split()
1642 else:
1643 # If the directWarp has <1% coverage, the psfMatchedWarp can have 0% and not exist
1644 # In this case, mask the whole epoch
1645 nansMask = afwImage.MaskX(coaddBBox, 1)
1646 spanSetList = []
1647 spanSetEdgeMask = []
1648 spanSetStreak = []
1649
1650 spanSetNoDataMask = afwGeom.SpanSet.fromMask(nansMask).split()
1651
1652 spanSetNoDataMaskList.append(spanSetNoDataMask)
1653 spanSetArtifactList.append(spanSetList)
1654 spanSetEdgeList.append(spanSetEdgeMask)
1655 if self.config.doFilterMorphological:
1656 spanSetBadMorphoList.append(spanSetStreak)
1657
1658 if lsstDebug.Info(__name__).saveCountIm:
1659 path = self._dataRef2DebugPath("epochCountIm", tempExpRefList[0], coaddLevel=True)
1660 epochCountImage.writeFits(path)
1661
1662 for i, spanSetList in enumerate(spanSetArtifactList):
1663 if spanSetList:
1664 filteredSpanSetList = self.filterArtifacts(spanSetList, epochCountImage, nImage,
1665 templateFootprints)
1666 spanSetArtifactList[i] = filteredSpanSetList
1667 if self.config.doFilterMorphological:
1668 spanSetArtifactList[i] += spanSetBadMorphoList[i]
1669
1670 altMasks = []
1671 for artifacts, noData, edge in zip(spanSetArtifactList, spanSetNoDataMaskList, spanSetEdgeList):
1672 altMasks.append({'CLIPPED': artifacts,
1673 'NO_DATA': noData,
1674 'EDGE': edge})
1675 return altMasks
1676
1677 def prefilterArtifacts(self, spanSetList, exp):
1678 """Remove artifact candidates covered by bad mask plane.
1679
1680 Any future editing of the candidate list that does not depend on
1681 temporal information should go in this method.
1682
1683 Parameters
1684 ----------
1685 spanSetList : `list` of `lsst.afw.geom.SpanSet`
1686 List of SpanSets representing artifact candidates.
1688 Exposure containing mask planes used to prefilter.
1689
1690 Returns
1691 -------
1692 returnSpanSetList : `list` of `lsst.afw.geom.SpanSet`
1693 List of SpanSets with artifacts.
1694 """
1695 badPixelMask = exp.mask.getPlaneBitMask(self.config.prefilterArtifactsMaskPlanes)
1696 goodArr = (exp.mask.array & badPixelMask) == 0
1697 returnSpanSetList = []
1698 bbox = exp.getBBox()
1699 x0, y0 = exp.getXY0()
1700 for i, span in enumerate(spanSetList):
1701 y, x = span.clippedTo(bbox).indices()
1702 yIndexLocal = numpy.array(y) - y0
1703 xIndexLocal = numpy.array(x) - x0
1704 goodRatio = numpy.count_nonzero(goodArr[yIndexLocal, xIndexLocal])/span.getArea()
1705 if goodRatio > self.config.prefilterArtifactsRatio:
1706 returnSpanSetList.append(span)
1707 return returnSpanSetList
1708
1709 def filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None):
1710 """Filter artifact candidates.
1711
1712 Parameters
1713 ----------
1714 spanSetList : `list` of `lsst.afw.geom.SpanSet`
1715 List of SpanSets representing artifact candidates.
1716 epochCountImage : `lsst.afw.image.Image`
1717 Image of accumulated number of warpDiff detections.
1718 nImage : `lsst.afw.image.ImageU`
1719 Image of the accumulated number of total epochs contributing.
1720
1721 Returns
1722 -------
1723 maskSpanSetList : `list`
1724 List of SpanSets with artifacts.
1725 """
1726 maskSpanSetList = []
1727 x0, y0 = epochCountImage.getXY0()
1728 for i, span in enumerate(spanSetList):
1729 y, x = span.indices()
1730 yIdxLocal = [y1 - y0 for y1 in y]
1731 xIdxLocal = [x1 - x0 for x1 in x]
1732 outlierN = epochCountImage.array[yIdxLocal, xIdxLocal]
1733 totalN = nImage.array[yIdxLocal, xIdxLocal]
1734
1735 # effectiveMaxNumEpochs is broken line (fraction of N) with characteristic config.maxNumEpochs
1736 effMaxNumEpochsHighN = (self.config.maxNumEpochs
1737 + self.config.maxFractionEpochsHigh*numpy.mean(totalN))
1738 effMaxNumEpochsLowN = self.config.maxFractionEpochsLow * numpy.mean(totalN)
1739 effectiveMaxNumEpochs = int(min(effMaxNumEpochsLowN, effMaxNumEpochsHighN))
1740 nPixelsBelowThreshold = numpy.count_nonzero((outlierN > 0)
1741 & (outlierN <= effectiveMaxNumEpochs))
1742 percentBelowThreshold = nPixelsBelowThreshold / len(outlierN)
1743 if percentBelowThreshold > self.config.spatialThreshold:
1744 maskSpanSetList.append(span)
1745
1746 if self.config.doPreserveContainedBySource and footprintsToExclude is not None:
1747 # If a candidate is contained by a footprint on the template coadd, do not clip
1748 filteredMaskSpanSetList = []
1749 for span in maskSpanSetList:
1750 doKeep = True
1751 for footprint in footprintsToExclude.positive.getFootprints():
1752 if footprint.spans.contains(span):
1753 doKeep = False
1754 break
1755 if doKeep:
1756 filteredMaskSpanSetList.append(span)
1757 maskSpanSetList = filteredMaskSpanSetList
1758
1759 return maskSpanSetList
1760
1761 def _readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd):
1762 """Fetch a warp from the butler and return a warpDiff.
1763
1764 Parameters
1765 ----------
1766 warpRef : `lsst.daf.butler.DeferredDatasetHandle`
1767 Handle for the warp.
1769 An image scaler object.
1770 templateCoadd : `lsst.afw.image.Exposure`
1771 Exposure to be substracted from the scaled warp.
1772
1773 Returns
1774 -------
1776 Exposure of the image difference between the warp and template.
1777 """
1778 # If the PSF-Matched warp did not exist for this direct warp
1779 # None is holding its place to maintain order in Gen 3
1780 if warpRef is None:
1781 return None
1782
1783 warp = warpRef.get()
1784 # direct image scaler OK for PSF-matched Warp
1785 imageScaler.scaleMaskedImage(warp.getMaskedImage())
1786 mi = warp.getMaskedImage()
1787 if self.config.doScaleWarpVariance:
1788 try:
1789 self.scaleWarpVariance.run(mi)
1790 except Exception as exc:
1791 self.log.warning("Unable to rescale variance of warp (%s); leaving it as-is", exc)
1792 mi -= templateCoadd.getMaskedImage()
1793 return warp
def assembleMetadata(self, coaddExposure, tempExpRefList, weightList)
def assembleOnlineMeanCoadd(self, coaddExposure, tempExpRefList, imageScalerList, weightList, altMaskList, statsCtrl, nImage=None)
def runQuantum(self, butlerQC, inputRefs, outputRefs)
def _makeSupplementaryData(self, butlerQC, inputRefs, outputRefs)
def processResults(self, coaddExposure, brightObjectMasks=None, dataId=None)
def setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None)
def assembleSubregion(self, coaddExposure, bbox, tempExpRefList, imageScalerList, weightList, altMaskList, statsFlags, statsCtrl, nImage=None)
def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, altMaskList=None, mask=None, supplementaryData=None)
def applyAltMaskPlanes(self, mask, altMaskSpans)
def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs)
def findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList)
def _readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd)
def filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None)
def getTempExpDatasetName(self, warpType="direct")
Definition: coaddBase.py:143
def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask)