lsst.pipe.tasks g7ab1c79000+374aa26984
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
391 super().__init__(**kwargs)
392 self.makeSubtask("interpImage")
393 self.makeSubtask("scaleZeroPoint")
394
395 if self.config.doMaskBrightObjects:
396 mask = afwImage.Mask()
397 try:
398 self.brightObjectBitmask = 1 << mask.addMaskPlane(self.config.brightObjectMaskName)
399 except pexExceptions.LsstCppException:
400 raise RuntimeError("Unable to define mask plane for bright objects; planes used are %s" %
401 mask.getMaskPlaneDict().keys())
402 del mask
403
404 if self.config.doInputMap:
405 self.makeSubtask("inputMapper")
406
407 self.warpType = self.config.warpType
408
409 @utils.inheritDoc(pipeBase.PipelineTask)
410 def runQuantum(self, butlerQC, inputRefs, outputRefs):
411 inputData = butlerQC.get(inputRefs)
412
413 # Construct skyInfo expected by run
414 # Do not remove skyMap from inputData in case _makeSupplementaryData needs it
415 skyMap = inputData["skyMap"]
416 outputDataId = butlerQC.quantum.dataId
417
418 inputData['skyInfo'] = makeSkyInfo(skyMap,
419 tractId=outputDataId['tract'],
420 patchId=outputDataId['patch'])
421
422 if self.config.doSelectVisits:
423 warpRefList = self.filterWarps(inputData['inputWarps'], inputData['selectedVisits'])
424 else:
425 warpRefList = inputData['inputWarps']
426
427 inputs = self.prepareInputs(warpRefList)
428 self.log.info("Found %d %s", len(inputs.tempExpRefList),
430 if len(inputs.tempExpRefList) == 0:
431 raise pipeBase.NoWorkFound("No coadd temporary exposures found")
432
433 supplementaryData = self._makeSupplementaryData(butlerQC, inputRefs, outputRefs)
434 retStruct = self.run(inputData['skyInfo'], inputs.tempExpRefList, inputs.imageScalerList,
435 inputs.weightList, supplementaryData=supplementaryData)
436
437 inputData.setdefault('brightObjectMask', None)
438 if self.config.doMaskBrightObjects and inputData["brightObjectMask"] is None:
439 log.warning("doMaskBrightObjects is set to True, but brightObjectMask not loaded")
440 self.processResults(retStruct.coaddExposure, inputData['brightObjectMask'], outputDataId)
441
442 if self.config.doWrite:
443 butlerQC.put(retStruct, outputRefs)
444 return retStruct
445
446 def processResults(self, coaddExposure, brightObjectMasks=None, dataId=None):
447 """Interpolate over missing data and mask bright stars.
448
449 Parameters
450 ----------
451 coaddExposure : `lsst.afw.image.Exposure`
452 The coadded exposure to process.
453 brightObjectMasks : `lsst.afw.table` or `None`, optional
454 Table of bright objects to mask.
455 dataId : `lsst.daf.butler.DataId` or `None`, optional
456 Data identification.
457 """
458 if self.config.doInterp:
459 self.interpImage.run(coaddExposure.getMaskedImage(), planeName="NO_DATA")
460 # The variance must be positive; work around for DM-3201.
461 varArray = coaddExposure.variance.array
462 with numpy.errstate(invalid="ignore"):
463 varArray[:] = numpy.where(varArray > 0, varArray, numpy.inf)
464
465 if self.config.doMaskBrightObjects:
466 self.setBrightObjectMasks(coaddExposure, brightObjectMasks, dataId)
467
468 def _makeSupplementaryData(self, butlerQC, inputRefs, outputRefs):
469 """Make additional inputs to run() specific to subclasses (Gen3).
470
471 Duplicates interface of `runQuantum` method.
472 Available to be implemented by subclasses only if they need the
473 coadd dataRef for performing preliminary processing before
474 assembling the coadd.
475
476 Parameters
477 ----------
478 butlerQC : `~lsst.pipe.base.ButlerQuantumContext`
479 Gen3 Butler object for fetching additional data products before
480 running the Task specialized for quantum being processed.
481 inputRefs : `~lsst.pipe.base.InputQuantizedConnection`
482 Attributes are the names of the connections describing input dataset types.
483 Values are DatasetRefs that task consumes for corresponding dataset type.
484 DataIds are guaranteed to match data objects in ``inputData``.
485 outputRefs : `~lsst.pipe.base.OutputQuantizedConnection`
486 Attributes are the names of the connections describing output dataset types.
487 Values are DatasetRefs that task is to produce
488 for corresponding dataset type.
489 """
490 return pipeBase.Struct()
491
492 @deprecated(
493 reason="makeSupplementaryDataGen3 is deprecated in favor of _makeSupplementaryData",
494 version="v25.0",
495 category=FutureWarning
496 )
497 def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs):
498 return self._makeSupplementaryData(butlerQC, inputRefs, outputRefs)
499
500 def prepareInputs(self, refList):
501 """Prepare the input warps for coaddition by measuring the weight for
502 each warp and the scaling for the photometric zero point.
503
504 Each Warp has its own photometric zeropoint and background variance.
505 Before coadding these Warps together, compute a scale factor to
506 normalize the photometric zeropoint and compute the weight for each Warp.
507
508 Parameters
509 ----------
510 refList : `list`
511 List of data references to tempExp.
512
513 Returns
514 -------
515 result : `~lsst.pipe.base.Struct`
516 Results as a struct with attributes:
517
518 ``tempExprefList``
519 `list` of data references to tempExp.
520 ``weightList``
521 `list` of weightings.
522 ``imageScalerList``
523 `list` of image scalers.
524 """
525 statsCtrl = afwMath.StatisticsControl()
526 statsCtrl.setNumSigmaClip(self.config.sigmaClip)
527 statsCtrl.setNumIter(self.config.clipIter)
528 statsCtrl.setAndMask(self.getBadPixelMask())
529 statsCtrl.setNanSafe(True)
530 # compute tempExpRefList: a list of tempExpRef that actually exist
531 # and weightList: a list of the weight of the associated coadd tempExp
532 # and imageScalerList: a list of scale factors for the associated coadd tempExp
533 tempExpRefList = []
534 weightList = []
535 imageScalerList = []
536 tempExpName = self.getTempExpDatasetName(self.warpType)
537 for tempExpRef in refList:
538 tempExp = tempExpRef.get()
539 # Ignore any input warp that is empty of data
540 if numpy.isnan(tempExp.image.array).all():
541 continue
542 maskedImage = tempExp.getMaskedImage()
543 imageScaler = self.scaleZeroPoint.computeImageScaler(
544 exposure=tempExp,
545 dataRef=tempExpRef, # FIXME
546 )
547 try:
548 imageScaler.scaleMaskedImage(maskedImage)
549 except Exception as e:
550 self.log.warning("Scaling failed for %s (skipping it): %s", tempExpRef.dataId, e)
551 continue
552 statObj = afwMath.makeStatistics(maskedImage.getVariance(), maskedImage.getMask(),
553 afwMath.MEANCLIP, statsCtrl)
554 meanVar, meanVarErr = statObj.getResult(afwMath.MEANCLIP)
555 weight = 1.0 / float(meanVar)
556 if not numpy.isfinite(weight):
557 self.log.warning("Non-finite weight for %s: skipping", tempExpRef.dataId)
558 continue
559 self.log.info("Weight of %s %s = %0.3f", tempExpName, tempExpRef.dataId, weight)
560
561 del maskedImage
562 del tempExp
563
564 tempExpRefList.append(tempExpRef)
565 weightList.append(weight)
566 imageScalerList.append(imageScaler)
567
568 return pipeBase.Struct(tempExpRefList=tempExpRefList, weightList=weightList,
569 imageScalerList=imageScalerList)
570
571 def prepareStats(self, mask=None):
572 """Prepare the statistics for coadding images.
573
574 Parameters
575 ----------
576 mask : `int`, optional
577 Bit mask value to exclude from coaddition.
578
579 Returns
580 -------
581 stats : `~lsst.pipe.base.Struct`
582 Statistics as a struct with attributes:
583
584 ``statsCtrl``
585 Statistics control object for coadd (`~lsst.afw.math.StatisticsControl`).
586 ``statsFlags``
587 Statistic for coadd (`~lsst.afw.math.Property`).
588 """
589 if mask is None:
590 mask = self.getBadPixelMask()
591 statsCtrl = afwMath.StatisticsControl()
592 statsCtrl.setNumSigmaClip(self.config.sigmaClip)
593 statsCtrl.setNumIter(self.config.clipIter)
594 statsCtrl.setAndMask(mask)
595 statsCtrl.setNanSafe(True)
596 statsCtrl.setWeighted(True)
597 statsCtrl.setCalcErrorFromInputVariance(self.config.calcErrorFromInputVariance)
598 for plane, threshold in self.config.maskPropagationThresholds.items():
599 bit = afwImage.Mask.getMaskPlane(plane)
600 statsCtrl.setMaskPropagationThreshold(bit, threshold)
601 statsFlags = afwMath.stringToStatisticsProperty(self.config.statistic)
602 return pipeBase.Struct(ctrl=statsCtrl, flags=statsFlags)
603
604 @timeMethod
605 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
606 altMaskList=None, mask=None, supplementaryData=None):
607 """Assemble a coadd from input warps.
608
609 Assemble the coadd using the provided list of coaddTempExps. Since
610 the full coadd covers a patch (a large area), the assembly is
611 performed over small areas on the image at a time in order to
612 conserve memory usage. Iterate over subregions within the outer
613 bbox of the patch using `assembleSubregion` to stack the corresponding
614 subregions from the coaddTempExps with the statistic specified.
615 Set the edge bits the coadd mask based on the weight map.
616
617 Parameters
618 ----------
619 skyInfo : `~lsst.pipe.base.Struct`
620 Struct with geometric information about the patch.
621 tempExpRefList : `list`
622 List of data references to Warps (previously called CoaddTempExps).
623 imageScalerList : `list`
624 List of image scalers.
625 weightList : `list`
626 List of weights.
627 altMaskList : `list`, optional
628 List of alternate masks to use rather than those stored with
629 tempExp.
630 mask : `int`, optional
631 Bit mask value to exclude from coaddition.
632 supplementaryData : `~lsst.pipe.base.Struct`, optional
633 Struct with additional data products needed to assemble coadd.
634 Only used by subclasses that implement ``_makeSupplementaryData``
635 and override `run`.
636
637 Returns
638 -------
639 result : `~lsst.pipe.base.Struct`
640 Results as a struct with attributes:
641
642 ``coaddExposure``
643 Coadded exposure (``lsst.afw.image.Exposure``).
644 ``nImage``
645 Exposure count image (``lsst.afw.image.Image``), if requested.
646 ``inputMap``
647 Bit-wise map of inputs, if requested.
648 ``warpRefList``
649 Input list of refs to the warps (``lsst.daf.butler.DeferredDatasetHandle``)
650 (unmodified).
651 ``imageScalerList``
652 Input list of image scalers (`list`) (unmodified).
653 ``weightList``
654 Input list of weights (`list`) (unmodified).
655 """
656 tempExpName = self.getTempExpDatasetName(self.warpType)
657 self.log.info("Assembling %s %s", len(tempExpRefList), tempExpName)
658 stats = self.prepareStats(mask=mask)
659
660 if altMaskList is None:
661 altMaskList = [None]*len(tempExpRefList)
662
663 coaddExposure = afwImage.ExposureF(skyInfo.bbox, skyInfo.wcs)
664 coaddExposure.setPhotoCalib(self.scaleZeroPoint.getPhotoCalib())
665 coaddExposure.getInfo().setCoaddInputs(self.inputRecorder.makeCoaddInputs())
666 self.assembleMetadata(coaddExposure, tempExpRefList, weightList)
667 coaddMaskedImage = coaddExposure.getMaskedImage()
668 subregionSizeArr = self.config.subregionSize
669 subregionSize = geom.Extent2I(subregionSizeArr[0], subregionSizeArr[1])
670 # if nImage is requested, create a zero one which can be passed to assembleSubregion
671 if self.config.doNImage:
672 nImage = afwImage.ImageU(skyInfo.bbox)
673 else:
674 nImage = None
675 # If inputMap is requested, create the initial version that can be masked in
676 # assembleSubregion.
677 if self.config.doInputMap:
678 self.inputMapper.build_ccd_input_map(skyInfo.bbox,
679 skyInfo.wcs,
680 coaddExposure.getInfo().getCoaddInputs().ccds)
681
682 if self.config.doOnlineForMean and self.config.statistic == "MEAN":
683 try:
684 self.assembleOnlineMeanCoadd(coaddExposure, tempExpRefList, imageScalerList,
685 weightList, altMaskList, stats.ctrl,
686 nImage=nImage)
687 except Exception as e:
688 self.log.exception("Cannot compute online coadd %s", e)
689 raise
690 else:
691 for subBBox in self._subBBoxIter(skyInfo.bbox, subregionSize):
692 try:
693 self.assembleSubregion(coaddExposure, subBBox, tempExpRefList, imageScalerList,
694 weightList, altMaskList, stats.flags, stats.ctrl,
695 nImage=nImage)
696 except Exception as e:
697 self.log.exception("Cannot compute coadd %s: %s", subBBox, e)
698 raise
699
700 # If inputMap is requested, we must finalize the map after the accumulation.
701 if self.config.doInputMap:
702 self.inputMapper.finalize_ccd_input_map_mask()
703 inputMap = self.inputMapper.ccd_input_map
704 else:
705 inputMap = None
706
707 self.setInexactPsf(coaddMaskedImage.getMask())
708 # Despite the name, the following doesn't really deal with "EDGE" pixels: it identifies
709 # pixels that didn't receive any unmasked inputs (as occurs around the edge of the field).
710 coaddUtils.setCoaddEdgeBits(coaddMaskedImage.getMask(), coaddMaskedImage.getVariance())
711 return pipeBase.Struct(coaddExposure=coaddExposure, nImage=nImage,
712 warpRefList=tempExpRefList, imageScalerList=imageScalerList,
713 weightList=weightList, inputMap=inputMap)
714
715 def assembleMetadata(self, coaddExposure, tempExpRefList, weightList):
716 """Set the metadata for the coadd.
717
718 This basic implementation sets the filter from the first input.
719
720 Parameters
721 ----------
722 coaddExposure : `lsst.afw.image.Exposure`
723 The target exposure for the coadd.
724 tempExpRefList : `list`
725 List of data references to tempExp.
726 weightList : `list`
727 List of weights.
728
729 Raises
730 ------
731 AssertionError
732 Raised if there is a length mismatch.
733 """
734 assert len(tempExpRefList) == len(weightList), "Length mismatch"
735
736 # We load a single pixel of each coaddTempExp, because we just want to get at the metadata
737 # (and we need more than just the PropertySet that contains the header), which is not possible
738 # with the current butler (see #2777).
739 bbox = geom.Box2I(coaddExposure.getBBox().getMin(), geom.Extent2I(1, 1))
740
741 tempExpList = [tempExpRef.get(parameters={'bbox': bbox}) for tempExpRef in tempExpRefList]
742
743 numCcds = sum(len(tempExp.getInfo().getCoaddInputs().ccds) for tempExp in tempExpList)
744
745 # Set the coadd FilterLabel to the band of the first input exposure:
746 # Coadds are calibrated, so the physical label is now meaningless.
747 coaddExposure.setFilter(afwImage.FilterLabel(tempExpList[0].getFilter().bandLabel))
748 coaddInputs = coaddExposure.getInfo().getCoaddInputs()
749 coaddInputs.ccds.reserve(numCcds)
750 coaddInputs.visits.reserve(len(tempExpList))
751
752 for tempExp, weight in zip(tempExpList, weightList):
753 self.inputRecorder.addVisitToCoadd(coaddInputs, tempExp, weight)
754
755 if self.config.doUsePsfMatchedPolygons:
756 self.shrinkValidPolygons(coaddInputs)
757
758 coaddInputs.visits.sort()
759 coaddInputs.ccds.sort()
760 if self.warpType == "psfMatched":
761 # The modelPsf BBox for a psfMatchedWarp/coaddTempExp was dynamically defined by
762 # ModelPsfMatchTask as the square box bounding its spatially-variable, pre-matched WarpedPsf.
763 # Likewise, set the PSF of a PSF-Matched Coadd to the modelPsf
764 # having the maximum width (sufficient because square)
765 modelPsfList = [tempExp.getPsf() for tempExp in tempExpList]
766 modelPsfWidthList = [modelPsf.computeBBox(modelPsf.getAveragePosition()).getWidth()
767 for modelPsf in modelPsfList]
768 psf = modelPsfList[modelPsfWidthList.index(max(modelPsfWidthList))]
769 else:
770 psf = measAlg.CoaddPsf(coaddInputs.ccds, coaddExposure.getWcs(),
771 self.config.coaddPsf.makeControl())
772 coaddExposure.setPsf(psf)
773 apCorrMap = measAlg.makeCoaddApCorrMap(coaddInputs.ccds, coaddExposure.getBBox(afwImage.PARENT),
774 coaddExposure.getWcs())
775 coaddExposure.getInfo().setApCorrMap(apCorrMap)
776 if self.config.doAttachTransmissionCurve:
777 transmissionCurve = measAlg.makeCoaddTransmissionCurve(coaddExposure.getWcs(), coaddInputs.ccds)
778 coaddExposure.getInfo().setTransmissionCurve(transmissionCurve)
779
780 def assembleSubregion(self, coaddExposure, bbox, tempExpRefList, imageScalerList, weightList,
781 altMaskList, statsFlags, statsCtrl, nImage=None):
782 """Assemble the coadd for a sub-region.
783
784 For each coaddTempExp, check for (and swap in) an alternative mask
785 if one is passed. Remove mask planes listed in
786 `config.removeMaskPlanes`. Finally, stack the actual exposures using
787 `lsst.afw.math.statisticsStack` with the statistic specified by
788 statsFlags. Typically, the statsFlag will be one of lsst.afw.math.MEAN for
789 a mean-stack or `lsst.afw.math.MEANCLIP` for outlier rejection using
790 an N-sigma clipped mean where N and iterations are specified by
791 statsCtrl. Assign the stacked subregion back to the coadd.
792
793 Parameters
794 ----------
795 coaddExposure : `lsst.afw.image.Exposure`
796 The target exposure for the coadd.
797 bbox : `lsst.geom.Box`
798 Sub-region to coadd.
799 tempExpRefList : `list`
800 List of data reference to tempExp.
801 imageScalerList : `list`
802 List of image scalers.
803 weightList : `list`
804 List of weights.
805 altMaskList : `list`
806 List of alternate masks to use rather than those stored with
807 tempExp, or None. Each element is dict with keys = mask plane
808 name to which to add the spans.
809 statsFlags : `lsst.afw.math.Property`
810 Property object for statistic for coadd.
812 Statistics control object for coadd.
813 nImage : `lsst.afw.image.ImageU`, optional
814 Keeps track of exposure count for each pixel.
815 """
816 self.log.debug("Computing coadd over %s", bbox)
817
818 coaddExposure.mask.addMaskPlane("REJECTED")
819 coaddExposure.mask.addMaskPlane("CLIPPED")
820 coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
821 maskMap = self.setRejectedMaskMapping(statsCtrl)
822 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
823 maskedImageList = []
824 if nImage is not None:
825 subNImage = afwImage.ImageU(bbox.getWidth(), bbox.getHeight())
826 for tempExpRef, imageScaler, altMask in zip(tempExpRefList, imageScalerList, altMaskList):
827
828 exposure = tempExpRef.get(parameters={'bbox': bbox})
829
830 maskedImage = exposure.getMaskedImage()
831 mask = maskedImage.getMask()
832 if altMask is not None:
833 self.applyAltMaskPlanes(mask, altMask)
834 imageScaler.scaleMaskedImage(maskedImage)
835
836 # Add 1 for each pixel which is not excluded by the exclude mask.
837 # In legacyCoadd, pixels may also be excluded by afwMath.statisticsStack.
838 if nImage is not None:
839 subNImage.getArray()[maskedImage.getMask().getArray() & statsCtrl.getAndMask() == 0] += 1
840 if self.config.removeMaskPlanes:
841 self.removeMaskPlanes(maskedImage)
842 maskedImageList.append(maskedImage)
843
844 if self.config.doInputMap:
845 visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
846 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
847
848 with self.timer("stack"):
849 coaddSubregion = afwMath.statisticsStack(maskedImageList, statsFlags, statsCtrl, weightList,
850 clipped, # also set output to CLIPPED if sigma-clipped
851 maskMap)
852 coaddExposure.maskedImage.assign(coaddSubregion, bbox)
853 if nImage is not None:
854 nImage.assign(subNImage, bbox)
855
856 def assembleOnlineMeanCoadd(self, coaddExposure, tempExpRefList, imageScalerList, weightList,
857 altMaskList, statsCtrl, nImage=None):
858 """Assemble the coadd using the "online" method.
859
860 This method takes a running sum of images and weights to save memory.
861 It only works for MEAN statistics.
862
863 Parameters
864 ----------
865 coaddExposure : `lsst.afw.image.Exposure`
866 The target exposure for the coadd.
867 tempExpRefList : `list`
868 List of data reference to tempExp.
869 imageScalerList : `list`
870 List of image scalers.
871 weightList : `list`
872 List of weights.
873 altMaskList : `list`
874 List of alternate masks to use rather than those stored with
875 tempExp, or None. Each element is dict with keys = mask plane
876 name to which to add the spans.
878 Statistics control object for coadd.
879 nImage : `lsst.afw.image.ImageU`, optional
880 Keeps track of exposure count for each pixel.
881 """
882 self.log.debug("Computing online coadd.")
883
884 coaddExposure.mask.addMaskPlane("REJECTED")
885 coaddExposure.mask.addMaskPlane("CLIPPED")
886 coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
887 maskMap = self.setRejectedMaskMapping(statsCtrl)
888 thresholdDict = AccumulatorMeanStack.stats_ctrl_to_threshold_dict(statsCtrl)
889
890 bbox = coaddExposure.maskedImage.getBBox()
891
892 stacker = AccumulatorMeanStack(
893 coaddExposure.image.array.shape,
894 statsCtrl.getAndMask(),
895 mask_threshold_dict=thresholdDict,
896 mask_map=maskMap,
897 no_good_pixels_mask=statsCtrl.getNoGoodPixelsMask(),
898 calc_error_from_input_variance=self.config.calcErrorFromInputVariance,
899 compute_n_image=(nImage is not None)
900 )
901
902 for tempExpRef, imageScaler, altMask, weight in zip(tempExpRefList,
903 imageScalerList,
904 altMaskList,
905 weightList):
906 exposure = tempExpRef.get()
907 maskedImage = exposure.getMaskedImage()
908 mask = maskedImage.getMask()
909 if altMask is not None:
910 self.applyAltMaskPlanes(mask, altMask)
911 imageScaler.scaleMaskedImage(maskedImage)
912 if self.config.removeMaskPlanes:
913 self.removeMaskPlanes(maskedImage)
914
915 stacker.add_masked_image(maskedImage, weight=weight)
916
917 if self.config.doInputMap:
918 visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
919 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
920
921 stacker.fill_stacked_masked_image(coaddExposure.maskedImage)
922
923 if nImage is not None:
924 nImage.array[:, :] = stacker.n_image
925
926 def removeMaskPlanes(self, maskedImage):
927 """Unset the mask of an image for mask planes specified in the config.
928
929 Parameters
930 ----------
931 maskedImage : `lsst.afw.image.MaskedImage`
932 The masked image to be modified.
933
934 Raises
935 ------
936 InvalidParameterError
937 Raised if no mask plane with that name was found.
938 """
939 mask = maskedImage.getMask()
940 for maskPlane in self.config.removeMaskPlanes:
941 try:
942 mask &= ~mask.getPlaneBitMask(maskPlane)
943 except pexExceptions.InvalidParameterError:
944 self.log.debug("Unable to remove mask plane %s: no mask plane with that name was found.",
945 maskPlane)
946
947 @staticmethod
949 """Map certain mask planes of the warps to new planes for the coadd.
950
951 If a pixel is rejected due to a mask value other than EDGE, NO_DATA,
952 or CLIPPED, set it to REJECTED on the coadd.
953 If a pixel is rejected due to EDGE, set the coadd pixel to SENSOR_EDGE.
954 If a pixel is rejected due to CLIPPED, set the coadd pixel to CLIPPED.
955
956 Parameters
957 ----------
959 Statistics control object for coadd.
960
961 Returns
962 -------
963 maskMap : `list` of `tuple` of `int`
964 A list of mappings of mask planes of the warped exposures to
965 mask planes of the coadd.
966 """
967 edge = afwImage.Mask.getPlaneBitMask("EDGE")
968 noData = afwImage.Mask.getPlaneBitMask("NO_DATA")
969 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
970 toReject = statsCtrl.getAndMask() & (~noData) & (~edge) & (~clipped)
971 maskMap = [(toReject, afwImage.Mask.getPlaneBitMask("REJECTED")),
972 (edge, afwImage.Mask.getPlaneBitMask("SENSOR_EDGE")),
973 (clipped, clipped)]
974 return maskMap
975
976 def applyAltMaskPlanes(self, mask, altMaskSpans):
977 """Apply in place alt mask formatted as SpanSets to a mask.
978
979 Parameters
980 ----------
981 mask : `lsst.afw.image.Mask`
982 Original mask.
983 altMaskSpans : `dict`
984 SpanSet lists to apply. Each element contains the new mask
985 plane name (e.g. "CLIPPED and/or "NO_DATA") as the key,
986 and list of SpanSets to apply to the mask.
987
988 Returns
989 -------
990 mask : `lsst.afw.image.Mask`
991 Updated mask.
992 """
993 if self.config.doUsePsfMatchedPolygons:
994 if ("NO_DATA" in altMaskSpans) and ("NO_DATA" in self.config.badMaskPlanes):
995 # Clear away any other masks outside the validPolygons. These pixels are no longer
996 # contributing to inexact PSFs, and will still be rejected because of NO_DATA
997 # self.config.doUsePsfMatchedPolygons should be True only in CompareWarpAssemble
998 # This mask-clearing step must only occur *before* applying the new masks below
999 for spanSet in altMaskSpans['NO_DATA']:
1000 spanSet.clippedTo(mask.getBBox()).clearMask(mask, self.getBadPixelMask())
1001
1002 for plane, spanSetList in altMaskSpans.items():
1003 maskClipValue = mask.addMaskPlane(plane)
1004 for spanSet in spanSetList:
1005 spanSet.clippedTo(mask.getBBox()).setMask(mask, 2**maskClipValue)
1006 return mask
1007
1008 def shrinkValidPolygons(self, coaddInputs):
1009 """Shrink coaddInputs' ccds' ValidPolygons in place.
1010
1011 Either modify each ccd's validPolygon in place, or if CoaddInputs
1012 does not have a validPolygon, create one from its bbox.
1013
1014 Parameters
1015 ----------
1016 coaddInputs : `lsst.afw.image.coaddInputs`
1017 Original mask.
1018 """
1019 for ccd in coaddInputs.ccds:
1020 polyOrig = ccd.getValidPolygon()
1021 validPolyBBox = polyOrig.getBBox() if polyOrig else ccd.getBBox()
1022 validPolyBBox.grow(-self.config.matchingKernelSize//2)
1023 if polyOrig:
1024 validPolygon = polyOrig.intersectionSingle(validPolyBBox)
1025 else:
1026 validPolygon = afwGeom.polygon.Polygon(geom.Box2D(validPolyBBox))
1027 ccd.setValidPolygon(validPolygon)
1028
1029 def setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None):
1030 """Set the bright object masks.
1031
1032 Parameters
1033 ----------
1034 exposure : `lsst.afw.image.Exposure`
1035 Exposure under consideration.
1036 brightObjectMasks : `lsst.afw.table`
1037 Table of bright objects to mask.
1038 dataId : `lsst.daf.butler.DataId`, optional
1039 Data identifier dict for patch.
1040 """
1041 if brightObjectMasks is None:
1042 self.log.warning("Unable to apply bright object mask: none supplied")
1043 return
1044 self.log.info("Applying %d bright object masks to %s", len(brightObjectMasks), dataId)
1045 mask = exposure.getMaskedImage().getMask()
1046 wcs = exposure.getWcs()
1047 plateScale = wcs.getPixelScale().asArcseconds()
1048
1049 for rec in brightObjectMasks:
1050 center = geom.PointI(wcs.skyToPixel(rec.getCoord()))
1051 if rec["type"] == "box":
1052 assert rec["angle"] == 0.0, ("Angle != 0 for mask object %s" % rec["id"])
1053 width = rec["width"].asArcseconds()/plateScale # convert to pixels
1054 height = rec["height"].asArcseconds()/plateScale # convert to pixels
1055
1056 halfSize = geom.ExtentI(0.5*width, 0.5*height)
1057 bbox = geom.Box2I(center - halfSize, center + halfSize)
1058
1059 bbox = geom.BoxI(geom.PointI(int(center[0] - 0.5*width), int(center[1] - 0.5*height)),
1060 geom.PointI(int(center[0] + 0.5*width), int(center[1] + 0.5*height)))
1061 spans = afwGeom.SpanSet(bbox)
1062 elif rec["type"] == "circle":
1063 radius = int(rec["radius"].asArcseconds()/plateScale) # convert to pixels
1064 spans = afwGeom.SpanSet.fromShape(radius, offset=center)
1065 else:
1066 self.log.warning("Unexpected region type %s at %s", rec["type"], center)
1067 continue
1068 spans.clippedTo(mask.getBBox()).setMask(mask, self.brightObjectBitmask)
1069
1070 def setInexactPsf(self, mask):
1071 """Set INEXACT_PSF mask plane.
1072
1073 If any of the input images isn't represented in the coadd (due to
1074 clipped pixels or chip gaps), the `CoaddPsf` will be inexact. Flag
1075 these pixels.
1076
1077 Parameters
1078 ----------
1079 mask : `lsst.afw.image.Mask`
1080 Coadded exposure's mask, modified in-place.
1081 """
1082 mask.addMaskPlane("INEXACT_PSF")
1083 inexactPsf = mask.getPlaneBitMask("INEXACT_PSF")
1084 sensorEdge = mask.getPlaneBitMask("SENSOR_EDGE") # chip edges (so PSF is discontinuous)
1085 clipped = mask.getPlaneBitMask("CLIPPED") # pixels clipped from coadd
1086 rejected = mask.getPlaneBitMask("REJECTED") # pixels rejected from coadd due to masks
1087 array = mask.getArray()
1088 selected = array & (sensorEdge | clipped | rejected) > 0
1089 array[selected] |= inexactPsf
1090
1091 @staticmethod
1092 def _subBBoxIter(bbox, subregionSize):
1093 """Iterate over subregions of a bbox.
1094
1095 Parameters
1096 ----------
1097 bbox : `lsst.geom.Box2I`
1098 Bounding box over which to iterate.
1099 subregionSize : `lsst.geom.Extent2I`
1100 Size of sub-bboxes.
1101
1102 Yields
1103 ------
1104 subBBox : `lsst.geom.Box2I`
1105 Next sub-bounding box of size ``subregionSize`` or smaller; each ``subBBox``
1106 is contained within ``bbox``, so it may be smaller than ``subregionSize`` at
1107 the edges of ``bbox``, but it will never be empty.
1108
1109 Raises
1110 ------
1111 RuntimeError
1112 Raised if any of the following occur:
1113 - The given bbox is empty.
1114 - The subregionSize is 0.
1115 """
1116 if bbox.isEmpty():
1117 raise RuntimeError("bbox %s is empty" % (bbox,))
1118 if subregionSize[0] < 1 or subregionSize[1] < 1:
1119 raise RuntimeError("subregionSize %s must be nonzero" % (subregionSize,))
1120
1121 for rowShift in range(0, bbox.getHeight(), subregionSize[1]):
1122 for colShift in range(0, bbox.getWidth(), subregionSize[0]):
1123 subBBox = geom.Box2I(bbox.getMin() + geom.Extent2I(colShift, rowShift), subregionSize)
1124 subBBox.clip(bbox)
1125 if subBBox.isEmpty():
1126 raise RuntimeError("Bug: empty bbox! bbox=%s, subregionSize=%s, "
1127 "colShift=%s, rowShift=%s" %
1128 (bbox, subregionSize, colShift, rowShift))
1129 yield subBBox
1130
1131 def filterWarps(self, inputs, goodVisits):
1132 """Return list of only inputRefs with visitId in goodVisits ordered by goodVisit.
1133
1134 Parameters
1135 ----------
1136 inputs : `list` of `~lsst.pipe.base.connections.DeferredDatasetRef`
1137 List of `lsst.pipe.base.connections.DeferredDatasetRef` with dataId containing visit.
1138 goodVisit : `dict`
1139 Dictionary with good visitIds as the keys. Value ignored.
1140
1141 Returns
1142 -------
1143 filteredInputs : `list` of `~lsst.pipe.base.connections.DeferredDatasetRef`
1144 Filtered and sorted list of inputRefs with visitId in goodVisits ordered by goodVisit.
1145 """
1146 inputWarpDict = {inputRef.ref.dataId['visit']: inputRef for inputRef in inputs}
1147 filteredInputs = []
1148 for visit in goodVisits.keys():
1149 if visit in inputWarpDict:
1150 filteredInputs.append(inputWarpDict[visit])
1151 return filteredInputs
1152
1153
1154def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask):
1155 """Function to count the number of pixels with a specific mask in a
1156 footprint.
1157
1158 Find the intersection of mask & footprint. Count all pixels in the mask
1159 that are in the intersection that have bitmask set but do not have
1160 ignoreMask set. Return the count.
1161
1162 Parameters
1163 ----------
1164 mask : `lsst.afw.image.Mask`
1165 Mask to define intersection region by.
1166 footprint : `lsst.afw.detection.Footprint`
1167 Footprint to define the intersection region by.
1168 bitmask : `Unknown`
1169 Specific mask that we wish to count the number of occurances of.
1170 ignoreMask : `Unknown`
1171 Pixels to not consider.
1172
1173 Returns
1174 -------
1175 result : `int`
1176 Number of pixels in footprint with specified mask.
1177 """
1178 bbox = footprint.getBBox()
1179 bbox.clip(mask.getBBox(afwImage.PARENT))
1180 fp = afwImage.Mask(bbox)
1181 subMask = mask.Factory(mask, bbox, afwImage.PARENT)
1182 footprint.spans.setMask(fp, bitmask)
1183 return numpy.logical_and((subMask.getArray() & fp.getArray()) > 0,
1184 (subMask.getArray() & ignoreMask) == 0).sum()
1185
1186
1188 psfMatchedWarps = pipeBase.connectionTypes.Input(
1189 doc=("PSF-Matched Warps are required by CompareWarp regardless of the coadd type requested. "
1190 "Only PSF-Matched Warps make sense for image subtraction. "
1191 "Therefore, they must be an additional declared input."),
1192 name="{inputCoaddName}Coadd_psfMatchedWarp",
1193 storageClass="ExposureF",
1194 dimensions=("tract", "patch", "skymap", "visit"),
1195 deferLoad=True,
1196 multiple=True
1197 )
1198 templateCoadd = pipeBase.connectionTypes.Output(
1199 doc=("Model of the static sky, used to find temporal artifacts. Typically a PSF-Matched, "
1200 "sigma-clipped coadd. Written if and only if assembleStaticSkyModel.doWrite=True"),
1201 name="{outputCoaddName}CoaddPsfMatched",
1202 storageClass="ExposureF",
1203 dimensions=("tract", "patch", "skymap", "band"),
1204 )
1205
1206 def __init__(self, *, config=None):
1207 super().__init__(config=config)
1208 if not config.assembleStaticSkyModel.doWrite:
1209 self.outputs.remove("templateCoadd")
1210 config.validate()
1211
1212
1213class CompareWarpAssembleCoaddConfig(AssembleCoaddConfig,
1214 pipelineConnections=CompareWarpAssembleCoaddConnections):
1215 assembleStaticSkyModel = pexConfig.ConfigurableField(
1216 target=AssembleCoaddTask,
1217 doc="Task to assemble an artifact-free, PSF-matched Coadd to serve as a"
1218 " naive/first-iteration model of the static sky.",
1219 )
1220 detect = pexConfig.ConfigurableField(
1221 target=SourceDetectionTask,
1222 doc="Detect outlier sources on difference between each psfMatched warp and static sky model"
1223 )
1224 detectTemplate = pexConfig.ConfigurableField(
1225 target=SourceDetectionTask,
1226 doc="Detect sources on static sky model. Only used if doPreserveContainedBySource is True"
1227 )
1228 maskStreaks = pexConfig.ConfigurableField(
1229 target=MaskStreaksTask,
1230 doc="Detect streaks on difference between each psfMatched warp and static sky model. Only used if "
1231 "doFilterMorphological is True. Adds a mask plane to an exposure, with the mask plane name set by"
1232 "streakMaskName"
1233 )
1234 streakMaskName = pexConfig.Field(
1235 dtype=str,
1236 default="STREAK",
1237 doc="Name of mask bit used for streaks"
1238 )
1239 maxNumEpochs = pexConfig.Field(
1240 doc="Charactistic maximum local number of epochs/visits in which an artifact candidate can appear "
1241 "and still be masked. The effective maxNumEpochs is a broken linear function of local "
1242 "number of epochs (N): min(maxFractionEpochsLow*N, maxNumEpochs + maxFractionEpochsHigh*N). "
1243 "For each footprint detected on the image difference between the psfMatched warp and static sky "
1244 "model, if a significant fraction of pixels (defined by spatialThreshold) are residuals in more "
1245 "than the computed effective maxNumEpochs, the artifact candidate is deemed persistant rather "
1246 "than transient and not masked.",
1247 dtype=int,
1248 default=2
1249 )
1250 maxFractionEpochsLow = pexConfig.RangeField(
1251 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for low N. "
1252 "Effective maxNumEpochs = "
1253 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1254 dtype=float,
1255 default=0.4,
1256 min=0., max=1.,
1257 )
1258 maxFractionEpochsHigh = pexConfig.RangeField(
1259 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for high N. "
1260 "Effective maxNumEpochs = "
1261 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1262 dtype=float,
1263 default=0.03,
1264 min=0., max=1.,
1265 )
1266 spatialThreshold = pexConfig.RangeField(
1267 doc="Unitless fraction of pixels defining how much of the outlier region has to meet the "
1268 "temporal criteria. If 0, clip all. If 1, clip none.",
1269 dtype=float,
1270 default=0.5,
1271 min=0., max=1.,
1272 inclusiveMin=True, inclusiveMax=True
1273 )
1274 doScaleWarpVariance = pexConfig.Field(
1275 doc="Rescale Warp variance plane using empirical noise?",
1276 dtype=bool,
1277 default=True,
1278 )
1279 scaleWarpVariance = pexConfig.ConfigurableField(
1280 target=ScaleVarianceTask,
1281 doc="Rescale variance on warps",
1282 )
1283 doPreserveContainedBySource = pexConfig.Field(
1284 doc="Rescue artifacts from clipping that completely lie within a footprint detected"
1285 "on the PsfMatched Template Coadd. Replicates a behavior of SafeClip.",
1286 dtype=bool,
1287 default=True,
1288 )
1289 doPrefilterArtifacts = pexConfig.Field(
1290 doc="Ignore artifact candidates that are mostly covered by the bad pixel mask, "
1291 "because they will be excluded anyway. This prevents them from contributing "
1292 "to the outlier epoch count image and potentially being labeled as persistant."
1293 "'Mostly' is defined by the config 'prefilterArtifactsRatio'.",
1294 dtype=bool,
1295 default=True
1296 )
1297 prefilterArtifactsMaskPlanes = pexConfig.ListField(
1298 doc="Prefilter artifact candidates that are mostly covered by these bad mask planes.",
1299 dtype=str,
1300 default=('NO_DATA', 'BAD', 'SAT', 'SUSPECT'),
1301 )
1302 prefilterArtifactsRatio = pexConfig.Field(
1303 doc="Prefilter artifact candidates with less than this fraction overlapping good pixels",
1304 dtype=float,
1305 default=0.05
1306 )
1307 doFilterMorphological = pexConfig.Field(
1308 doc="Filter artifact candidates based on morphological criteria, i.g. those that appear to "
1309 "be streaks.",
1310 dtype=bool,
1311 default=False
1312 )
1313 growStreakFp = pexConfig.Field(
1314 doc="Grow streak footprints by this number multiplied by the PSF width",
1315 dtype=float,
1316 default=5
1317 )
1318
1319 def setDefaults(self):
1320 AssembleCoaddConfig.setDefaults(self)
1321 self.statistic = 'MEAN'
1323
1324 # Real EDGE removed by psfMatched NO_DATA border half the width of the matching kernel
1325 # CompareWarp applies psfMatched EDGE pixels to directWarps before assembling
1326 if "EDGE" in self.badMaskPlanes:
1327 self.badMaskPlanes.remove('EDGE')
1328 self.removeMaskPlanes.append('EDGE')
1329 self.assembleStaticSkyModel.badMaskPlanes = ["NO_DATA", ]
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("No dataset type exists for a PSF-Matched Template N Image."
1357 "Please set assembleStaticSkyModel.doNImage=False")
1358
1359 if self.assembleStaticSkyModel.doWrite and (self.warpType == self.assembleStaticSkyModel.warpType):
1360 raise ValueError("warpType (%s) == assembleStaticSkyModel.warpType (%s) and will compete for "
1361 "the same dataset name. Please set assembleStaticSkyModel.doWrite to False "
1362 "or warpType to 'direct'. assembleStaticSkyModel.warpType should ways be "
1363 "'PsfMatched'" % (self.warpType, self.assembleStaticSkyModel.warpType))
1364
1365
1367 """Assemble a compareWarp coadded image from a set of warps
1368 by masking artifacts detected by comparing PSF-matched warps.
1369
1370 In ``AssembleCoaddTask``, we compute the coadd as an clipped mean (i.e.,
1371 we clip outliers). The problem with doing this is that when computing the
1372 coadd PSF at a given location, individual visit PSFs from visits with
1373 outlier pixels contribute to the coadd PSF and cannot be treated correctly.
1374 In this task, we correct for this behavior by creating a new badMaskPlane
1375 'CLIPPED' which marks pixels in the individual warps suspected to contain
1376 an artifact. We populate this plane on the input warps by comparing
1377 PSF-matched warps with a PSF-matched median coadd which serves as a
1378 model of the static sky. Any group of pixels that deviates from the
1379 PSF-matched template coadd by more than config.detect.threshold sigma,
1380 is an artifact candidate. The candidates are then filtered to remove
1381 variable sources and sources that are difficult to subtract such as
1382 bright stars. This filter is configured using the config parameters
1383 ``temporalThreshold`` and ``spatialThreshold``. The temporalThreshold is
1384 the maximum fraction of epochs that the deviation can appear in and still
1385 be considered an artifact. The spatialThreshold is the maximum fraction of
1386 pixels in the footprint of the deviation that appear in other epochs
1387 (where other epochs is defined by the temporalThreshold). If the deviant
1388 region meets this criteria of having a significant percentage of pixels
1389 that deviate in only a few epochs, these pixels have the 'CLIPPED' bit
1390 set in the mask. These regions will not contribute to the final coadd.
1391 Furthermore, any routine to determine the coadd PSF can now be cognizant
1392 of clipped regions. Note that the algorithm implemented by this task is
1393 preliminary and works correctly for HSC data. Parameter modifications and
1394 or considerable redesigning of the algorithm is likley required for other
1395 surveys.
1396
1397 ``CompareWarpAssembleCoaddTask`` sub-classes
1398 ``AssembleCoaddTask`` and instantiates ``AssembleCoaddTask``
1399 as a subtask to generate the TemplateCoadd (the model of the static sky).
1400
1401 Notes
1402 -----
1403 Debugging:
1404 This task supports the following debug variables:
1405 - ``saveCountIm``
1406 If True then save the Epoch Count Image as a fits file in the `figPath`
1407 - ``figPath``
1408 Path to save the debug fits images and figures
1409 """
1410
1411 ConfigClass = CompareWarpAssembleCoaddConfig
1412 _DefaultName = "compareWarpAssembleCoadd"
1413
1414 def __init__(self, *args, **kwargs):
1415 AssembleCoaddTask.__init__(self, *args, **kwargs)
1416 self.makeSubtask("assembleStaticSkyModel")
1417 detectionSchema = afwTable.SourceTable.makeMinimalSchema()
1418 self.makeSubtask("detect", schema=detectionSchema)
1419 if self.config.doPreserveContainedBySource:
1420 self.makeSubtask("detectTemplate", schema=afwTable.SourceTable.makeMinimalSchema())
1421 if self.config.doScaleWarpVariance:
1422 self.makeSubtask("scaleWarpVariance")
1423 if self.config.doFilterMorphological:
1424 self.makeSubtask("maskStreaks")
1425
1426 @utils.inheritDoc(AssembleCoaddTask)
1427 def _makeSupplementaryData(self, butlerQC, inputRefs, outputRefs):
1428 """Generate a templateCoadd to use as a naive model of static sky to
1429 subtract from PSF-Matched warps.
1430
1431 Returns
1432 -------
1433 result : `~lsst.pipe.base.Struct`
1434 Results as a struct with attributes:
1435
1436 ``templateCoadd``
1437 Coadded exposure (`lsst.afw.image.Exposure`).
1438 ``nImage``
1439 Keeps track of exposure count for each pixel (`lsst.afw.image.ImageU`).
1440
1441 Raises
1442 ------
1443 RuntimeError
1444 Raised if ``templateCoadd`` is `None`.
1445 """
1446 # Ensure that psfMatchedWarps are used as input warps for template generation
1447 staticSkyModelInputRefs = copy.deepcopy(inputRefs)
1448 staticSkyModelInputRefs.inputWarps = inputRefs.psfMatchedWarps
1449
1450 # Because subtasks don't have connections we have to make one.
1451 # The main task's `templateCoadd` is the subtask's `coaddExposure`
1452 staticSkyModelOutputRefs = copy.deepcopy(outputRefs)
1453 if self.config.assembleStaticSkyModel.doWrite:
1454 staticSkyModelOutputRefs.coaddExposure = staticSkyModelOutputRefs.templateCoadd
1455 # Remove template coadd from both subtask's and main tasks outputs,
1456 # because it is handled by the subtask as `coaddExposure`
1457 del outputRefs.templateCoadd
1458 del staticSkyModelOutputRefs.templateCoadd
1459
1460 # A PSF-Matched nImage does not exist as a dataset type
1461 if 'nImage' in staticSkyModelOutputRefs.keys():
1462 del staticSkyModelOutputRefs.nImage
1463
1464 templateCoadd = self.assembleStaticSkyModel.runQuantum(butlerQC, staticSkyModelInputRefs,
1465 staticSkyModelOutputRefs)
1466 if templateCoadd is None:
1467 raise RuntimeError(self._noTemplateMessage(self.assembleStaticSkyModel.warpType))
1468
1469 return pipeBase.Struct(templateCoadd=templateCoadd.coaddExposure,
1470 nImage=templateCoadd.nImage,
1471 warpRefList=templateCoadd.warpRefList,
1472 imageScalerList=templateCoadd.imageScalerList,
1473 weightList=templateCoadd.weightList)
1474
1475 def _noTemplateMessage(self, warpType):
1476 warpName = (warpType[0].upper() + warpType[1:])
1477 message = """No %(warpName)s warps were found to build the template coadd which is
1478 required to run CompareWarpAssembleCoaddTask. To continue assembling this type of coadd,
1479 first either rerun makeCoaddTempExp with config.make%(warpName)s=True or
1480 coaddDriver with config.makeCoadTempExp.make%(warpName)s=True, before assembleCoadd.
1481
1482 Alternatively, to use another algorithm with existing warps, retarget the CoaddDriverConfig to
1483 another algorithm like:
1484
1485 from lsst.pipe.tasks.assembleCoadd import SafeClipAssembleCoaddTask
1486 config.assemble.retarget(SafeClipAssembleCoaddTask)
1487 """ % {"warpName": warpName}
1488 return message
1489
1490 @utils.inheritDoc(AssembleCoaddTask)
1491 @timeMethod
1492 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
1493 supplementaryData):
1494 """Assemble the coadd.
1495
1496 Find artifacts and apply them to the warps' masks creating a list of
1497 alternative masks with a new "CLIPPED" plane and updated "NO_DATA"
1498 plane. Then pass these alternative masks to the base class's ``run``
1499 method.
1500 """
1501 # Check and match the order of the supplementaryData
1502 # (PSF-matched) inputs to the order of the direct inputs,
1503 # so that the artifact mask is applied to the right warp
1504 dataIds = [ref.dataId for ref in tempExpRefList]
1505 psfMatchedDataIds = [ref.dataId for ref in supplementaryData.warpRefList]
1506
1507 if dataIds != psfMatchedDataIds:
1508 self.log.info("Reordering and or/padding PSF-matched visit input list")
1509 supplementaryData.warpRefList = reorderAndPadList(supplementaryData.warpRefList,
1510 psfMatchedDataIds, dataIds)
1511 supplementaryData.imageScalerList = reorderAndPadList(supplementaryData.imageScalerList,
1512 psfMatchedDataIds, dataIds)
1513
1514 # Use PSF-Matched Warps (and corresponding scalers) and coadd to find artifacts
1515 spanSetMaskList = self.findArtifacts(supplementaryData.templateCoadd,
1516 supplementaryData.warpRefList,
1517 supplementaryData.imageScalerList)
1518
1519 badMaskPlanes = self.config.badMaskPlanes[:]
1520 badMaskPlanes.append("CLIPPED")
1521 badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes)
1522
1523 result = AssembleCoaddTask.run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
1524 spanSetMaskList, mask=badPixelMask)
1525
1526 # Propagate PSF-matched EDGE pixels to coadd SENSOR_EDGE and INEXACT_PSF
1527 # Psf-Matching moves the real edge inwards
1528 self.applyAltEdgeMask(result.coaddExposure.maskedImage.mask, spanSetMaskList)
1529 return result
1530
1531 def applyAltEdgeMask(self, mask, altMaskList):
1532 """Propagate alt EDGE mask to SENSOR_EDGE AND INEXACT_PSF planes.
1533
1534 Parameters
1535 ----------
1536 mask : `lsst.afw.image.Mask`
1537 Original mask.
1538 altMaskList : `list` of `dict`
1539 List of Dicts containing ``spanSet`` lists.
1540 Each element contains the new mask plane name (e.g. "CLIPPED
1541 and/or "NO_DATA") as the key, and list of ``SpanSets`` to apply to
1542 the mask.
1543 """
1544 maskValue = mask.getPlaneBitMask(["SENSOR_EDGE", "INEXACT_PSF"])
1545 for visitMask in altMaskList:
1546 if "EDGE" in visitMask:
1547 for spanSet in visitMask['EDGE']:
1548 spanSet.clippedTo(mask.getBBox()).setMask(mask, maskValue)
1549
1550 def findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList):
1551 """Find artifacts.
1552
1553 Loop through warps twice. The first loop builds a map with the count
1554 of how many epochs each pixel deviates from the templateCoadd by more
1555 than ``config.chiThreshold`` sigma. The second loop takes each
1556 difference image and filters the artifacts detected in each using
1557 count map to filter out variable sources and sources that are
1558 difficult to subtract cleanly.
1559
1560 Parameters
1561 ----------
1562 templateCoadd : `lsst.afw.image.Exposure`
1563 Exposure to serve as model of static sky.
1564 tempExpRefList : `list`
1565 List of data references to warps.
1566 imageScalerList : `list`
1567 List of image scalers.
1568
1569 Returns
1570 -------
1571 altMasks : `list` of `dict`
1572 List of dicts containing information about CLIPPED
1573 (i.e., artifacts), NO_DATA, and EDGE pixels.
1574 """
1575 self.log.debug("Generating Count Image, and mask lists.")
1576 coaddBBox = templateCoadd.getBBox()
1577 slateIm = afwImage.ImageU(coaddBBox)
1578 epochCountImage = afwImage.ImageU(coaddBBox)
1579 nImage = afwImage.ImageU(coaddBBox)
1580 spanSetArtifactList = []
1581 spanSetNoDataMaskList = []
1582 spanSetEdgeList = []
1583 spanSetBadMorphoList = []
1584 badPixelMask = self.getBadPixelMask()
1585
1586 # mask of the warp diffs should = that of only the warp
1587 templateCoadd.mask.clearAllMaskPlanes()
1588
1589 if self.config.doPreserveContainedBySource:
1590 templateFootprints = self.detectTemplate.detectFootprints(templateCoadd)
1591 else:
1592 templateFootprints = None
1593
1594 for warpRef, imageScaler in zip(tempExpRefList, imageScalerList):
1595 warpDiffExp = self._readAndComputeWarpDiff(warpRef, imageScaler, templateCoadd)
1596 if warpDiffExp is not None:
1597 # This nImage only approximates the final nImage because it uses the PSF-matched mask
1598 nImage.array += (numpy.isfinite(warpDiffExp.image.array)
1599 * ((warpDiffExp.mask.array & badPixelMask) == 0)).astype(numpy.uint16)
1600 fpSet = self.detect.detectFootprints(warpDiffExp, doSmooth=False, clearMask=True)
1601 fpSet.positive.merge(fpSet.negative)
1602 footprints = fpSet.positive
1603 slateIm.set(0)
1604 spanSetList = [footprint.spans for footprint in footprints.getFootprints()]
1605
1606 # Remove artifacts due to defects before they contribute to the epochCountImage
1607 if self.config.doPrefilterArtifacts:
1608 spanSetList = self.prefilterArtifacts(spanSetList, warpDiffExp)
1609
1610 # Clear mask before adding prefiltered spanSets
1611 self.detect.clearMask(warpDiffExp.mask)
1612 for spans in spanSetList:
1613 spans.setImage(slateIm, 1, doClip=True)
1614 spans.setMask(warpDiffExp.mask, warpDiffExp.mask.getPlaneBitMask("DETECTED"))
1615 epochCountImage += slateIm
1616
1617 if self.config.doFilterMorphological:
1618 maskName = self.config.streakMaskName
1619 _ = self.maskStreaks.run(warpDiffExp)
1620 streakMask = warpDiffExp.mask
1621 spanSetStreak = afwGeom.SpanSet.fromMask(streakMask,
1622 streakMask.getPlaneBitMask(maskName)).split()
1623 # Pad the streaks to account for low-surface brightness wings
1624 psf = warpDiffExp.getPsf()
1625 for s, sset in enumerate(spanSetStreak):
1626 psfShape = psf.computeShape(sset.computeCentroid())
1627 dilation = self.config.growStreakFp * psfShape.getDeterminantRadius()
1628 sset_dilated = sset.dilated(int(dilation))
1629 spanSetStreak[s] = sset_dilated
1630
1631 # PSF-Matched warps have less available area (~the matching kernel) because the calexps
1632 # undergo a second convolution. Pixels with data in the direct warp
1633 # but not in the PSF-matched warp will not have their artifacts detected.
1634 # NaNs from the PSF-matched warp therefore must be masked in the direct warp
1635 nans = numpy.where(numpy.isnan(warpDiffExp.maskedImage.image.array), 1, 0)
1636 nansMask = afwImage.makeMaskFromArray(nans.astype(afwImage.MaskPixel))
1637 nansMask.setXY0(warpDiffExp.getXY0())
1638 edgeMask = warpDiffExp.mask
1639 spanSetEdgeMask = afwGeom.SpanSet.fromMask(edgeMask,
1640 edgeMask.getPlaneBitMask("EDGE")).split()
1641 else:
1642 # If the directWarp has <1% coverage, the psfMatchedWarp can have 0% and not exist
1643 # In this case, mask the whole epoch
1644 nansMask = afwImage.MaskX(coaddBBox, 1)
1645 spanSetList = []
1646 spanSetEdgeMask = []
1647 spanSetStreak = []
1648
1649 spanSetNoDataMask = afwGeom.SpanSet.fromMask(nansMask).split()
1650
1651 spanSetNoDataMaskList.append(spanSetNoDataMask)
1652 spanSetArtifactList.append(spanSetList)
1653 spanSetEdgeList.append(spanSetEdgeMask)
1654 if self.config.doFilterMorphological:
1655 spanSetBadMorphoList.append(spanSetStreak)
1656
1657 if lsstDebug.Info(__name__).saveCountIm:
1658 path = self._dataRef2DebugPath("epochCountIm", tempExpRefList[0], coaddLevel=True)
1659 epochCountImage.writeFits(path)
1660
1661 for i, spanSetList in enumerate(spanSetArtifactList):
1662 if spanSetList:
1663 filteredSpanSetList = self.filterArtifacts(spanSetList, epochCountImage, nImage,
1664 templateFootprints)
1665 spanSetArtifactList[i] = filteredSpanSetList
1666 if self.config.doFilterMorphological:
1667 spanSetArtifactList[i] += spanSetBadMorphoList[i]
1668
1669 altMasks = []
1670 for artifacts, noData, edge in zip(spanSetArtifactList, spanSetNoDataMaskList, spanSetEdgeList):
1671 altMasks.append({'CLIPPED': artifacts,
1672 'NO_DATA': noData,
1673 'EDGE': edge})
1674 return altMasks
1675
1676 def prefilterArtifacts(self, spanSetList, exp):
1677 """Remove artifact candidates covered by bad mask plane.
1678
1679 Any future editing of the candidate list that does not depend on
1680 temporal information should go in this method.
1681
1682 Parameters
1683 ----------
1684 spanSetList : `list` of `lsst.afw.geom.SpanSet`
1685 List of SpanSets representing artifact candidates.
1687 Exposure containing mask planes used to prefilter.
1688
1689 Returns
1690 -------
1691 returnSpanSetList : `list` of `lsst.afw.geom.SpanSet`
1692 List of SpanSets with artifacts.
1693 """
1694 badPixelMask = exp.mask.getPlaneBitMask(self.config.prefilterArtifactsMaskPlanes)
1695 goodArr = (exp.mask.array & badPixelMask) == 0
1696 returnSpanSetList = []
1697 bbox = exp.getBBox()
1698 x0, y0 = exp.getXY0()
1699 for i, span in enumerate(spanSetList):
1700 y, x = span.clippedTo(bbox).indices()
1701 yIndexLocal = numpy.array(y) - y0
1702 xIndexLocal = numpy.array(x) - x0
1703 goodRatio = numpy.count_nonzero(goodArr[yIndexLocal, xIndexLocal])/span.getArea()
1704 if goodRatio > self.config.prefilterArtifactsRatio:
1705 returnSpanSetList.append(span)
1706 return returnSpanSetList
1707
1708 def filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None):
1709 """Filter artifact candidates.
1710
1711 Parameters
1712 ----------
1713 spanSetList : `list` of `lsst.afw.geom.SpanSet`
1714 List of SpanSets representing artifact candidates.
1715 epochCountImage : `lsst.afw.image.Image`
1716 Image of accumulated number of warpDiff detections.
1717 nImage : `lsst.afw.image.ImageU`
1718 Image of the accumulated number of total epochs contributing.
1719
1720 Returns
1721 -------
1722 maskSpanSetList : `list`
1723 List of SpanSets with artifacts.
1724 """
1725 maskSpanSetList = []
1726 x0, y0 = epochCountImage.getXY0()
1727 for i, span in enumerate(spanSetList):
1728 y, x = span.indices()
1729 yIdxLocal = [y1 - y0 for y1 in y]
1730 xIdxLocal = [x1 - x0 for x1 in x]
1731 outlierN = epochCountImage.array[yIdxLocal, xIdxLocal]
1732 totalN = nImage.array[yIdxLocal, xIdxLocal]
1733
1734 # effectiveMaxNumEpochs is broken line (fraction of N) with characteristic config.maxNumEpochs
1735 effMaxNumEpochsHighN = (self.config.maxNumEpochs
1736 + self.config.maxFractionEpochsHigh*numpy.mean(totalN))
1737 effMaxNumEpochsLowN = self.config.maxFractionEpochsLow * numpy.mean(totalN)
1738 effectiveMaxNumEpochs = int(min(effMaxNumEpochsLowN, effMaxNumEpochsHighN))
1739 nPixelsBelowThreshold = numpy.count_nonzero((outlierN > 0)
1740 & (outlierN <= effectiveMaxNumEpochs))
1741 percentBelowThreshold = nPixelsBelowThreshold / len(outlierN)
1742 if percentBelowThreshold > self.config.spatialThreshold:
1743 maskSpanSetList.append(span)
1744
1745 if self.config.doPreserveContainedBySource and footprintsToExclude is not None:
1746 # If a candidate is contained by a footprint on the template coadd, do not clip
1747 filteredMaskSpanSetList = []
1748 for span in maskSpanSetList:
1749 doKeep = True
1750 for footprint in footprintsToExclude.positive.getFootprints():
1751 if footprint.spans.contains(span):
1752 doKeep = False
1753 break
1754 if doKeep:
1755 filteredMaskSpanSetList.append(span)
1756 maskSpanSetList = filteredMaskSpanSetList
1757
1758 return maskSpanSetList
1759
1760 def _readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd):
1761 """Fetch a warp from the butler and return a warpDiff.
1762
1763 Parameters
1764 ----------
1765 warpRef : `lsst.daf.butler.DeferredDatasetHandle`
1766 Handle for the warp.
1768 An image scaler object.
1769 templateCoadd : `lsst.afw.image.Exposure`
1770 Exposure to be substracted from the scaled warp.
1771
1772 Returns
1773 -------
1775 Exposure of the image difference between the warp and template.
1776 """
1777 # If the PSF-Matched warp did not exist for this direct warp
1778 # None is holding its place to maintain order in Gen 3
1779 if warpRef is None:
1780 return None
1781
1782 warp = warpRef.get()
1783 # direct image scaler OK for PSF-matched Warp
1784 imageScaler.scaleMaskedImage(warp.getMaskedImage())
1785 mi = warp.getMaskedImage()
1786 if self.config.doScaleWarpVariance:
1787 try:
1788 self.scaleWarpVariance.run(mi)
1789 except Exception as exc:
1790 self.log.warning("Unable to rescale variance of warp (%s); leaving it as-is", exc)
1791 mi -= templateCoadd.getMaskedImage()
1792 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:170
def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask)