LSST Applications 26.0.0,g0265f82a02+6660c170cc,g07994bdeae+30b05a742e,g0a0026dc87+17526d298f,g0a60f58ba1+17526d298f,g0e4bf8285c+96dd2c2ea9,g0ecae5effc+c266a536c8,g1e7d6db67d+6f7cb1f4bb,g26482f50c6+6346c0633c,g2bbee38e9b+6660c170cc,g2cc88a2952+0a4e78cd49,g3273194fdb+f6908454ef,g337abbeb29+6660c170cc,g337c41fc51+9a8f8f0815,g37c6e7c3d5+7bbafe9d37,g44018dc512+6660c170cc,g4a941329ef+4f7594a38e,g4c90b7bd52+5145c320d2,g58be5f913a+bea990ba40,g635b316a6c+8d6b3a3e56,g67924a670a+bfead8c487,g6ae5381d9b+81bc2a20b4,g93c4d6e787+26b17396bd,g98cecbdb62+ed2cb6d659,g98ffbb4407+81bc2a20b4,g9ddcbc5298+7f7571301f,ga1e77700b3+99e9273977,gae46bcf261+6660c170cc,gb2715bf1a1+17526d298f,gc86a011abf+17526d298f,gcf0d15dbbd+96dd2c2ea9,gdaeeff99f8+0d8dbea60f,gdb4ec4c597+6660c170cc,ge23793e450+96dd2c2ea9,gf041782ebf+171108ac67
LSST Data Management Base Package
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: # TODO: Remove this in DM-39841
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 Raises
658 ------
659 lsst.pipe.base.NoWorkFound
660 Raised if no data references are provided.
661 """
662 tempExpName = self.getTempExpDatasetName(self.warpType)
663 self.log.info("Assembling %s %s", len(tempExpRefList), tempExpName)
664 if not tempExpRefList:
665 raise pipeBase.NoWorkFound("No exposures provided for co-addition.")
666
667 stats = self.prepareStats(mask=mask)
668
669 if altMaskList is None:
670 altMaskList = [None]*len(tempExpRefList)
671
672 coaddExposure = afwImage.ExposureF(skyInfo.bbox, skyInfo.wcs)
673 coaddExposure.setPhotoCalib(self.scaleZeroPoint.getPhotoCalib())
674 coaddExposure.getInfo().setCoaddInputs(self.inputRecorder.makeCoaddInputs())
675 self.assembleMetadata(coaddExposure, tempExpRefList, weightList)
676 coaddMaskedImage = coaddExposure.getMaskedImage()
677 subregionSizeArr = self.config.subregionSize
678 subregionSize = geom.Extent2I(subregionSizeArr[0], subregionSizeArr[1])
679 # if nImage is requested, create a zero one which can be passed to assembleSubregion
680 if self.config.doNImage:
681 nImage = afwImage.ImageU(skyInfo.bbox)
682 else:
683 nImage = None
684 # If inputMap is requested, create the initial version that can be masked in
685 # assembleSubregion.
686 if self.config.doInputMap:
687 self.inputMapper.build_ccd_input_map(skyInfo.bbox,
688 skyInfo.wcs,
689 coaddExposure.getInfo().getCoaddInputs().ccds)
690
691 if self.config.doOnlineForMean and self.config.statistic == "MEAN":
692 try:
693 self.assembleOnlineMeanCoadd(coaddExposure, tempExpRefList, imageScalerList,
694 weightList, altMaskList, stats.ctrl,
695 nImage=nImage)
696 except Exception as e:
697 self.log.exception("Cannot compute online coadd %s", e)
698 raise
699 else:
700 for subBBox in self._subBBoxIter(skyInfo.bbox, subregionSize):
701 try:
702 self.assembleSubregion(coaddExposure, subBBox, tempExpRefList, imageScalerList,
703 weightList, altMaskList, stats.flags, stats.ctrl,
704 nImage=nImage)
705 except Exception as e:
706 self.log.exception("Cannot compute coadd %s: %s", subBBox, e)
707 raise
708
709 # If inputMap is requested, we must finalize the map after the accumulation.
710 if self.config.doInputMap:
711 self.inputMapper.finalize_ccd_input_map_mask()
712 inputMap = self.inputMapper.ccd_input_map
713 else:
714 inputMap = None
715
716 self.setInexactPsf(coaddMaskedImage.getMask())
717 # Despite the name, the following doesn't really deal with "EDGE" pixels: it identifies
718 # pixels that didn't receive any unmasked inputs (as occurs around the edge of the field).
719 coaddUtils.setCoaddEdgeBits(coaddMaskedImage.getMask(), coaddMaskedImage.getVariance())
720 return pipeBase.Struct(coaddExposure=coaddExposure, nImage=nImage,
721 warpRefList=tempExpRefList, imageScalerList=imageScalerList,
722 weightList=weightList, inputMap=inputMap)
723
724 def assembleMetadata(self, coaddExposure, tempExpRefList, weightList):
725 """Set the metadata for the coadd.
726
727 This basic implementation sets the filter from the first input.
728
729 Parameters
730 ----------
731 coaddExposure : `lsst.afw.image.Exposure`
732 The target exposure for the coadd.
733 tempExpRefList : `list`
734 List of data references to tempExp.
735 weightList : `list`
736 List of weights.
737
738 Raises
739 ------
740 AssertionError
741 Raised if there is a length mismatch.
742 """
743 assert len(tempExpRefList) == len(weightList), "Length mismatch"
744
745 # We load a single pixel of each coaddTempExp, because we just want to get at the metadata
746 # (and we need more than just the PropertySet that contains the header), which is not possible
747 # with the current butler (see #2777).
748 bbox = geom.Box2I(coaddExposure.getBBox().getMin(), geom.Extent2I(1, 1))
749
750 tempExpList = [tempExpRef.get(parameters={'bbox': bbox}) for tempExpRef in tempExpRefList]
751
752 numCcds = sum(len(tempExp.getInfo().getCoaddInputs().ccds) for tempExp in tempExpList)
753
754 # Set the coadd FilterLabel to the band of the first input exposure:
755 # Coadds are calibrated, so the physical label is now meaningless.
756 coaddExposure.setFilter(afwImage.FilterLabel(tempExpList[0].getFilter().bandLabel))
757 coaddInputs = coaddExposure.getInfo().getCoaddInputs()
758 coaddInputs.ccds.reserve(numCcds)
759 coaddInputs.visits.reserve(len(tempExpList))
760
761 for tempExp, weight in zip(tempExpList, weightList):
762 self.inputRecorder.addVisitToCoadd(coaddInputs, tempExp, weight)
763
764 if self.config.doUsePsfMatchedPolygons:
765 self.shrinkValidPolygons(coaddInputs)
766
767 coaddInputs.visits.sort()
768 coaddInputs.ccds.sort()
769 if self.warpType == "psfMatched":
770 # The modelPsf BBox for a psfMatchedWarp/coaddTempExp was dynamically defined by
771 # ModelPsfMatchTask as the square box bounding its spatially-variable, pre-matched WarpedPsf.
772 # Likewise, set the PSF of a PSF-Matched Coadd to the modelPsf
773 # having the maximum width (sufficient because square)
774 modelPsfList = [tempExp.getPsf() for tempExp in tempExpList]
775 modelPsfWidthList = [modelPsf.computeBBox(modelPsf.getAveragePosition()).getWidth()
776 for modelPsf in modelPsfList]
777 psf = modelPsfList[modelPsfWidthList.index(max(modelPsfWidthList))]
778 else:
779 psf = measAlg.CoaddPsf(coaddInputs.ccds, coaddExposure.getWcs(),
780 self.config.coaddPsf.makeControl())
781 coaddExposure.setPsf(psf)
782 apCorrMap = measAlg.makeCoaddApCorrMap(coaddInputs.ccds, coaddExposure.getBBox(afwImage.PARENT),
783 coaddExposure.getWcs())
784 coaddExposure.getInfo().setApCorrMap(apCorrMap)
785 if self.config.doAttachTransmissionCurve:
786 transmissionCurve = measAlg.makeCoaddTransmissionCurve(coaddExposure.getWcs(), coaddInputs.ccds)
787 coaddExposure.getInfo().setTransmissionCurve(transmissionCurve)
788
789 def assembleSubregion(self, coaddExposure, bbox, tempExpRefList, imageScalerList, weightList,
790 altMaskList, statsFlags, statsCtrl, nImage=None):
791 """Assemble the coadd for a sub-region.
792
793 For each coaddTempExp, check for (and swap in) an alternative mask
794 if one is passed. Remove mask planes listed in
795 `config.removeMaskPlanes`. Finally, stack the actual exposures using
796 `lsst.afw.math.statisticsStack` with the statistic specified by
797 statsFlags. Typically, the statsFlag will be one of lsst.afw.math.MEAN for
798 a mean-stack or `lsst.afw.math.MEANCLIP` for outlier rejection using
799 an N-sigma clipped mean where N and iterations are specified by
800 statsCtrl. Assign the stacked subregion back to the coadd.
801
802 Parameters
803 ----------
804 coaddExposure : `lsst.afw.image.Exposure`
805 The target exposure for the coadd.
806 bbox : `lsst.geom.Box`
807 Sub-region to coadd.
808 tempExpRefList : `list`
809 List of data reference to tempExp.
810 imageScalerList : `list`
811 List of image scalers.
812 weightList : `list`
813 List of weights.
814 altMaskList : `list`
815 List of alternate masks to use rather than those stored with
816 tempExp, or None. Each element is dict with keys = mask plane
817 name to which to add the spans.
818 statsFlags : `lsst.afw.math.Property`
819 Property object for statistic for coadd.
821 Statistics control object for coadd.
822 nImage : `lsst.afw.image.ImageU`, optional
823 Keeps track of exposure count for each pixel.
824 """
825 self.log.debug("Computing coadd over %s", bbox)
826
827 coaddExposure.mask.addMaskPlane("REJECTED")
828 coaddExposure.mask.addMaskPlane("CLIPPED")
829 coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
830 maskMap = self.setRejectedMaskMapping(statsCtrl)
831 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
832 maskedImageList = []
833 if nImage is not None:
834 subNImage = afwImage.ImageU(bbox.getWidth(), bbox.getHeight())
835 for tempExpRef, imageScaler, altMask in zip(tempExpRefList, imageScalerList, altMaskList):
836
837 exposure = tempExpRef.get(parameters={'bbox': bbox})
838
839 maskedImage = exposure.getMaskedImage()
840 mask = maskedImage.getMask()
841 if altMask is not None:
842 self.applyAltMaskPlanes(mask, altMask)
843 imageScaler.scaleMaskedImage(maskedImage)
844
845 # Add 1 for each pixel which is not excluded by the exclude mask.
846 # In legacyCoadd, pixels may also be excluded by afwMath.statisticsStack.
847 if nImage is not None:
848 subNImage.getArray()[maskedImage.getMask().getArray() & statsCtrl.getAndMask() == 0] += 1
849 if self.config.removeMaskPlanes:
850 self.removeMaskPlanes(maskedImage)
851 maskedImageList.append(maskedImage)
852
853 if self.config.doInputMap:
854 visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
855 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
856
857 with self.timer("stack"):
858 coaddSubregion = afwMath.statisticsStack(maskedImageList, statsFlags, statsCtrl, weightList,
859 clipped, # also set output to CLIPPED if sigma-clipped
860 maskMap)
861 coaddExposure.maskedImage.assign(coaddSubregion, bbox)
862 if nImage is not None:
863 nImage.assign(subNImage, bbox)
864
865 def assembleOnlineMeanCoadd(self, coaddExposure, tempExpRefList, imageScalerList, weightList,
866 altMaskList, statsCtrl, nImage=None):
867 """Assemble the coadd using the "online" method.
868
869 This method takes a running sum of images and weights to save memory.
870 It only works for MEAN statistics.
871
872 Parameters
873 ----------
874 coaddExposure : `lsst.afw.image.Exposure`
875 The target exposure for the coadd.
876 tempExpRefList : `list`
877 List of data reference to tempExp.
878 imageScalerList : `list`
879 List of image scalers.
880 weightList : `list`
881 List of weights.
882 altMaskList : `list`
883 List of alternate masks to use rather than those stored with
884 tempExp, or None. Each element is dict with keys = mask plane
885 name to which to add the spans.
887 Statistics control object for coadd.
888 nImage : `lsst.afw.image.ImageU`, optional
889 Keeps track of exposure count for each pixel.
890 """
891 self.log.debug("Computing online coadd.")
892
893 coaddExposure.mask.addMaskPlane("REJECTED")
894 coaddExposure.mask.addMaskPlane("CLIPPED")
895 coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
896 maskMap = self.setRejectedMaskMapping(statsCtrl)
897 thresholdDict = AccumulatorMeanStack.stats_ctrl_to_threshold_dict(statsCtrl)
898
899 bbox = coaddExposure.maskedImage.getBBox()
900
901 stacker = AccumulatorMeanStack(
902 coaddExposure.image.array.shape,
903 statsCtrl.getAndMask(),
904 mask_threshold_dict=thresholdDict,
905 mask_map=maskMap,
906 no_good_pixels_mask=statsCtrl.getNoGoodPixelsMask(),
907 calc_error_from_input_variance=self.config.calcErrorFromInputVariance,
908 compute_n_image=(nImage is not None)
909 )
910
911 for tempExpRef, imageScaler, altMask, weight in zip(tempExpRefList,
912 imageScalerList,
913 altMaskList,
914 weightList):
915 exposure = tempExpRef.get()
916 maskedImage = exposure.getMaskedImage()
917 mask = maskedImage.getMask()
918 if altMask is not None:
919 self.applyAltMaskPlanes(mask, altMask)
920 imageScaler.scaleMaskedImage(maskedImage)
921 if self.config.removeMaskPlanes:
922 self.removeMaskPlanes(maskedImage)
923
924 stacker.add_masked_image(maskedImage, weight=weight)
925
926 if self.config.doInputMap:
927 visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
928 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
929
930 stacker.fill_stacked_masked_image(coaddExposure.maskedImage)
931
932 if nImage is not None:
933 nImage.array[:, :] = stacker.n_image
934
935 def removeMaskPlanes(self, maskedImage):
936 """Unset the mask of an image for mask planes specified in the config.
937
938 Parameters
939 ----------
940 maskedImage : `lsst.afw.image.MaskedImage`
941 The masked image to be modified.
942
943 Raises
944 ------
945 InvalidParameterError
946 Raised if no mask plane with that name was found.
947 """
948 mask = maskedImage.getMask()
949 for maskPlane in self.config.removeMaskPlanes:
950 try:
951 mask &= ~mask.getPlaneBitMask(maskPlane)
953 self.log.debug("Unable to remove mask plane %s: no mask plane with that name was found.",
954 maskPlane)
955
956 @staticmethod
958 """Map certain mask planes of the warps to new planes for the coadd.
959
960 If a pixel is rejected due to a mask value other than EDGE, NO_DATA,
961 or CLIPPED, set it to REJECTED on the coadd.
962 If a pixel is rejected due to EDGE, set the coadd pixel to SENSOR_EDGE.
963 If a pixel is rejected due to CLIPPED, set the coadd pixel to CLIPPED.
964
965 Parameters
966 ----------
968 Statistics control object for coadd.
969
970 Returns
971 -------
972 maskMap : `list` of `tuple` of `int`
973 A list of mappings of mask planes of the warped exposures to
974 mask planes of the coadd.
975 """
976 edge = afwImage.Mask.getPlaneBitMask("EDGE")
977 noData = afwImage.Mask.getPlaneBitMask("NO_DATA")
978 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
979 toReject = statsCtrl.getAndMask() & (~noData) & (~edge) & (~clipped)
980 maskMap = [(toReject, afwImage.Mask.getPlaneBitMask("REJECTED")),
981 (edge, afwImage.Mask.getPlaneBitMask("SENSOR_EDGE")),
982 (clipped, clipped)]
983 return maskMap
984
985 def applyAltMaskPlanes(self, mask, altMaskSpans):
986 """Apply in place alt mask formatted as SpanSets to a mask.
987
988 Parameters
989 ----------
990 mask : `lsst.afw.image.Mask`
991 Original mask.
992 altMaskSpans : `dict`
993 SpanSet lists to apply. Each element contains the new mask
994 plane name (e.g. "CLIPPED and/or "NO_DATA") as the key,
995 and list of SpanSets to apply to the mask.
996
997 Returns
998 -------
999 mask : `lsst.afw.image.Mask`
1000 Updated mask.
1001 """
1002 if self.config.doUsePsfMatchedPolygons:
1003 if ("NO_DATA" in altMaskSpans) and ("NO_DATA" in self.config.badMaskPlanes):
1004 # Clear away any other masks outside the validPolygons. These pixels are no longer
1005 # contributing to inexact PSFs, and will still be rejected because of NO_DATA
1006 # self.config.doUsePsfMatchedPolygons should be True only in CompareWarpAssemble
1007 # This mask-clearing step must only occur *before* applying the new masks below
1008 for spanSet in altMaskSpans['NO_DATA']:
1009 spanSet.clippedTo(mask.getBBox()).clearMask(mask, self.getBadPixelMask())
1010
1011 for plane, spanSetList in altMaskSpans.items():
1012 maskClipValue = mask.addMaskPlane(plane)
1013 for spanSet in spanSetList:
1014 spanSet.clippedTo(mask.getBBox()).setMask(mask, 2**maskClipValue)
1015 return mask
1016
1017 def shrinkValidPolygons(self, coaddInputs):
1018 """Shrink coaddInputs' ccds' ValidPolygons in place.
1019
1020 Either modify each ccd's validPolygon in place, or if CoaddInputs
1021 does not have a validPolygon, create one from its bbox.
1022
1023 Parameters
1024 ----------
1025 coaddInputs : `lsst.afw.image.coaddInputs`
1026 Original mask.
1027 """
1028 for ccd in coaddInputs.ccds:
1029 polyOrig = ccd.getValidPolygon()
1030 validPolyBBox = polyOrig.getBBox() if polyOrig else ccd.getBBox()
1031 validPolyBBox.grow(-self.config.matchingKernelSize//2)
1032 if polyOrig:
1033 validPolygon = polyOrig.intersectionSingle(validPolyBBox)
1034 else:
1035 validPolygon = afwGeom.polygon.Polygon(geom.Box2D(validPolyBBox))
1036 ccd.setValidPolygon(validPolygon)
1037
1038 def setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None):
1039 """Set the bright object masks.
1040
1041 Parameters
1042 ----------
1043 exposure : `lsst.afw.image.Exposure`
1044 Exposure under consideration.
1045 brightObjectMasks : `lsst.afw.table`
1046 Table of bright objects to mask.
1047 dataId : `lsst.daf.butler.DataId`, optional
1048 Data identifier dict for patch.
1049 """
1050 if brightObjectMasks is None:
1051 self.log.warning("Unable to apply bright object mask: none supplied")
1052 return
1053 self.log.info("Applying %d bright object masks to %s", len(brightObjectMasks), dataId)
1054 mask = exposure.getMaskedImage().getMask()
1055 wcs = exposure.getWcs()
1056 plateScale = wcs.getPixelScale().asArcseconds()
1057
1058 for rec in brightObjectMasks:
1059 center = geom.PointI(wcs.skyToPixel(rec.getCoord()))
1060 if rec["type"] == "box":
1061 assert rec["angle"] == 0.0, ("Angle != 0 for mask object %s" % rec["id"])
1062 width = rec["width"].asArcseconds()/plateScale # convert to pixels
1063 height = rec["height"].asArcseconds()/plateScale # convert to pixels
1064
1065 halfSize = geom.ExtentI(0.5*width, 0.5*height)
1066 bbox = geom.Box2I(center - halfSize, center + halfSize)
1067
1068 bbox = geom.BoxI(geom.PointI(int(center[0] - 0.5*width), int(center[1] - 0.5*height)),
1069 geom.PointI(int(center[0] + 0.5*width), int(center[1] + 0.5*height)))
1070 spans = afwGeom.SpanSet(bbox)
1071 elif rec["type"] == "circle":
1072 radius = int(rec["radius"].asArcseconds()/plateScale) # convert to pixels
1073 spans = afwGeom.SpanSet.fromShape(radius, offset=center)
1074 else:
1075 self.log.warning("Unexpected region type %s at %s", rec["type"], center)
1076 continue
1077 spans.clippedTo(mask.getBBox()).setMask(mask, self.brightObjectBitmask)
1078
1079 def setInexactPsf(self, mask):
1080 """Set INEXACT_PSF mask plane.
1081
1082 If any of the input images isn't represented in the coadd (due to
1083 clipped pixels or chip gaps), the `CoaddPsf` will be inexact. Flag
1084 these pixels.
1085
1086 Parameters
1087 ----------
1088 mask : `lsst.afw.image.Mask`
1089 Coadded exposure's mask, modified in-place.
1090 """
1091 mask.addMaskPlane("INEXACT_PSF")
1092 inexactPsf = mask.getPlaneBitMask("INEXACT_PSF")
1093 sensorEdge = mask.getPlaneBitMask("SENSOR_EDGE") # chip edges (so PSF is discontinuous)
1094 clipped = mask.getPlaneBitMask("CLIPPED") # pixels clipped from coadd
1095 rejected = mask.getPlaneBitMask("REJECTED") # pixels rejected from coadd due to masks
1096 array = mask.getArray()
1097 selected = array & (sensorEdge | clipped | rejected) > 0
1098 array[selected] |= inexactPsf
1099
1100 @staticmethod
1101 def _subBBoxIter(bbox, subregionSize):
1102 """Iterate over subregions of a bbox.
1103
1104 Parameters
1105 ----------
1106 bbox : `lsst.geom.Box2I`
1107 Bounding box over which to iterate.
1108 subregionSize : `lsst.geom.Extent2I`
1109 Size of sub-bboxes.
1110
1111 Yields
1112 ------
1113 subBBox : `lsst.geom.Box2I`
1114 Next sub-bounding box of size ``subregionSize`` or smaller; each ``subBBox``
1115 is contained within ``bbox``, so it may be smaller than ``subregionSize`` at
1116 the edges of ``bbox``, but it will never be empty.
1117
1118 Raises
1119 ------
1120 RuntimeError
1121 Raised if any of the following occur:
1122 - The given bbox is empty.
1123 - The subregionSize is 0.
1124 """
1125 if bbox.isEmpty():
1126 raise RuntimeError("bbox %s is empty" % (bbox,))
1127 if subregionSize[0] < 1 or subregionSize[1] < 1:
1128 raise RuntimeError("subregionSize %s must be nonzero" % (subregionSize,))
1129
1130 for rowShift in range(0, bbox.getHeight(), subregionSize[1]):
1131 for colShift in range(0, bbox.getWidth(), subregionSize[0]):
1132 subBBox = geom.Box2I(bbox.getMin() + geom.Extent2I(colShift, rowShift), subregionSize)
1133 subBBox.clip(bbox)
1134 if subBBox.isEmpty():
1135 raise RuntimeError("Bug: empty bbox! bbox=%s, subregionSize=%s, "
1136 "colShift=%s, rowShift=%s" %
1137 (bbox, subregionSize, colShift, rowShift))
1138 yield subBBox
1139
1140 def filterWarps(self, inputs, goodVisits):
1141 """Return list of only inputRefs with visitId in goodVisits ordered by goodVisit.
1142
1143 Parameters
1144 ----------
1145 inputs : `list` of `~lsst.pipe.base.connections.DeferredDatasetRef`
1146 List of `lsst.pipe.base.connections.DeferredDatasetRef` with dataId containing visit.
1147 goodVisit : `dict`
1148 Dictionary with good visitIds as the keys. Value ignored.
1149
1150 Returns
1151 -------
1152 filteredInputs : `list` of `~lsst.pipe.base.connections.DeferredDatasetRef`
1153 Filtered and sorted list of inputRefs with visitId in goodVisits ordered by goodVisit.
1154 """
1155 inputWarpDict = {inputRef.ref.dataId['visit']: inputRef for inputRef in inputs}
1156 filteredInputs = []
1157 for visit in goodVisits.keys():
1158 if visit in inputWarpDict:
1159 filteredInputs.append(inputWarpDict[visit])
1160 return filteredInputs
1161
1162
1163def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask):
1164 """Function to count the number of pixels with a specific mask in a
1165 footprint.
1166
1167 Find the intersection of mask & footprint. Count all pixels in the mask
1168 that are in the intersection that have bitmask set but do not have
1169 ignoreMask set. Return the count.
1170
1171 Parameters
1172 ----------
1173 mask : `lsst.afw.image.Mask`
1174 Mask to define intersection region by.
1175 footprint : `lsst.afw.detection.Footprint`
1176 Footprint to define the intersection region by.
1177 bitmask : `Unknown`
1178 Specific mask that we wish to count the number of occurances of.
1179 ignoreMask : `Unknown`
1180 Pixels to not consider.
1181
1182 Returns
1183 -------
1184 result : `int`
1185 Number of pixels in footprint with specified mask.
1186 """
1187 bbox = footprint.getBBox()
1188 bbox.clip(mask.getBBox(afwImage.PARENT))
1189 fp = afwImage.Mask(bbox)
1190 subMask = mask.Factory(mask, bbox, afwImage.PARENT)
1191 footprint.spans.setMask(fp, bitmask)
1192 return numpy.logical_and((subMask.getArray() & fp.getArray()) > 0,
1193 (subMask.getArray() & ignoreMask) == 0).sum()
1194
1195
1197 psfMatchedWarps = pipeBase.connectionTypes.Input(
1198 doc=("PSF-Matched Warps are required by CompareWarp regardless of the coadd type requested. "
1199 "Only PSF-Matched Warps make sense for image subtraction. "
1200 "Therefore, they must be an additional declared input."),
1201 name="{inputCoaddName}Coadd_psfMatchedWarp",
1202 storageClass="ExposureF",
1203 dimensions=("tract", "patch", "skymap", "visit"),
1204 deferLoad=True,
1205 multiple=True
1206 )
1207 templateCoadd = pipeBase.connectionTypes.Output(
1208 doc=("Model of the static sky, used to find temporal artifacts. Typically a PSF-Matched, "
1209 "sigma-clipped coadd. Written if and only if assembleStaticSkyModel.doWrite=True"),
1210 name="{outputCoaddName}CoaddPsfMatched",
1211 storageClass="ExposureF",
1212 dimensions=("tract", "patch", "skymap", "band"),
1213 )
1214
1215 def __init__(self, *, config=None):
1216 super().__init__(config=config)
1217 if not config.assembleStaticSkyModel.doWrite:
1218 self.outputs.remove("templateCoadd")
1219 config.validate()
1220
1221
1222class CompareWarpAssembleCoaddConfig(AssembleCoaddConfig,
1223 pipelineConnections=CompareWarpAssembleCoaddConnections):
1224 assembleStaticSkyModel = pexConfig.ConfigurableField(
1225 target=AssembleCoaddTask,
1226 doc="Task to assemble an artifact-free, PSF-matched Coadd to serve as a"
1227 " naive/first-iteration model of the static sky.",
1228 )
1229 detect = pexConfig.ConfigurableField(
1230 target=SourceDetectionTask,
1231 doc="Detect outlier sources on difference between each psfMatched warp and static sky model"
1232 )
1233 detectTemplate = pexConfig.ConfigurableField(
1234 target=SourceDetectionTask,
1235 doc="Detect sources on static sky model. Only used if doPreserveContainedBySource is True"
1236 )
1237 maskStreaks = pexConfig.ConfigurableField(
1238 target=MaskStreaksTask,
1239 doc="Detect streaks on difference between each psfMatched warp and static sky model. Only used if "
1240 "doFilterMorphological is True. Adds a mask plane to an exposure, with the mask plane name set by"
1241 "streakMaskName"
1242 )
1243 streakMaskName = pexConfig.Field(
1244 dtype=str,
1245 default="STREAK",
1246 doc="Name of mask bit used for streaks"
1247 )
1248 maxNumEpochs = pexConfig.Field(
1249 doc="Charactistic maximum local number of epochs/visits in which an artifact candidate can appear "
1250 "and still be masked. The effective maxNumEpochs is a broken linear function of local "
1251 "number of epochs (N): min(maxFractionEpochsLow*N, maxNumEpochs + maxFractionEpochsHigh*N). "
1252 "For each footprint detected on the image difference between the psfMatched warp and static sky "
1253 "model, if a significant fraction of pixels (defined by spatialThreshold) are residuals in more "
1254 "than the computed effective maxNumEpochs, the artifact candidate is deemed persistant rather "
1255 "than transient and not masked.",
1256 dtype=int,
1257 default=2
1258 )
1259 maxFractionEpochsLow = pexConfig.RangeField(
1260 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for low N. "
1261 "Effective maxNumEpochs = "
1262 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1263 dtype=float,
1264 default=0.4,
1265 min=0., max=1.,
1266 )
1267 maxFractionEpochsHigh = pexConfig.RangeField(
1268 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for high N. "
1269 "Effective maxNumEpochs = "
1270 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1271 dtype=float,
1272 default=0.03,
1273 min=0., max=1.,
1274 )
1275 spatialThreshold = pexConfig.RangeField(
1276 doc="Unitless fraction of pixels defining how much of the outlier region has to meet the "
1277 "temporal criteria. If 0, clip all. If 1, clip none.",
1278 dtype=float,
1279 default=0.5,
1280 min=0., max=1.,
1281 inclusiveMin=True, inclusiveMax=True
1282 )
1283 doScaleWarpVariance = pexConfig.Field(
1284 doc="Rescale Warp variance plane using empirical noise?",
1285 dtype=bool,
1286 default=True,
1287 )
1288 scaleWarpVariance = pexConfig.ConfigurableField(
1289 target=ScaleVarianceTask,
1290 doc="Rescale variance on warps",
1291 )
1292 doPreserveContainedBySource = pexConfig.Field(
1293 doc="Rescue artifacts from clipping that completely lie within a footprint detected"
1294 "on the PsfMatched Template Coadd. Replicates a behavior of SafeClip.",
1295 dtype=bool,
1296 default=True,
1297 )
1298 doPrefilterArtifacts = pexConfig.Field(
1299 doc="Ignore artifact candidates that are mostly covered by the bad pixel mask, "
1300 "because they will be excluded anyway. This prevents them from contributing "
1301 "to the outlier epoch count image and potentially being labeled as persistant."
1302 "'Mostly' is defined by the config 'prefilterArtifactsRatio'.",
1303 dtype=bool,
1304 default=True
1305 )
1306 prefilterArtifactsMaskPlanes = pexConfig.ListField(
1307 doc="Prefilter artifact candidates that are mostly covered by these bad mask planes.",
1308 dtype=str,
1309 default=('NO_DATA', 'BAD', 'SAT', 'SUSPECT'),
1310 )
1311 prefilterArtifactsRatio = pexConfig.Field(
1312 doc="Prefilter artifact candidates with less than this fraction overlapping good pixels",
1313 dtype=float,
1314 default=0.05
1315 )
1316 doFilterMorphological = pexConfig.Field(
1317 doc="Filter artifact candidates based on morphological criteria, i.g. those that appear to "
1318 "be streaks.",
1319 dtype=bool,
1320 default=False
1321 )
1322 growStreakFp = pexConfig.Field(
1323 doc="Grow streak footprints by this number multiplied by the PSF width",
1324 dtype=float,
1325 default=5
1326 )
1327
1328 def setDefaults(self):
1329 AssembleCoaddConfig.setDefaults(self)
1330 self.statistic = 'MEAN'
1332
1333 # Real EDGE removed by psfMatched NO_DATA border half the width of the matching kernel
1334 # CompareWarp applies psfMatched EDGE pixels to directWarps before assembling
1335 if "EDGE" in self.badMaskPlanes:
1336 self.badMaskPlanes.remove('EDGE')
1337 self.removeMaskPlanes.append('EDGE')
1338 self.assembleStaticSkyModel.badMaskPlanes = ["NO_DATA", ]
1339 self.assembleStaticSkyModel.warpType = 'psfMatched'
1340 self.assembleStaticSkyModel.connections.warpType = 'psfMatched'
1341 self.assembleStaticSkyModel.statistic = 'MEANCLIP'
1342 self.assembleStaticSkyModel.sigmaClip = 2.5
1343 self.assembleStaticSkyModel.clipIter = 3
1344 self.assembleStaticSkyModel.calcErrorFromInputVariance = False
1345 self.assembleStaticSkyModel.doWrite = False
1346 self.detect.doTempLocalBackground = False
1347 self.detect.reEstimateBackground = False
1348 self.detect.returnOriginalFootprints = False
1349 self.detect.thresholdPolarity = "both"
1350 self.detect.thresholdValue = 5
1351 self.detect.minPixels = 4
1352 self.detect.isotropicGrow = True
1353 self.detect.thresholdType = "pixel_stdev"
1354 self.detect.nSigmaToGrow = 0.4
1355 # The default nSigmaToGrow for SourceDetectionTask is already 2.4,
1356 # Explicitly restating because ratio with detect.nSigmaToGrow matters
1357 self.detectTemplate.nSigmaToGrow = 2.4
1358 self.detectTemplate.doTempLocalBackground = False
1359 self.detectTemplate.reEstimateBackground = False
1360 self.detectTemplate.returnOriginalFootprints = False
1361
1362 def validate(self):
1363 super().validate()
1364 if self.assembleStaticSkyModel.doNImage:
1365 raise ValueError("No dataset type exists for a PSF-Matched Template N Image."
1366 "Please set assembleStaticSkyModel.doNImage=False")
1367
1368 if self.assembleStaticSkyModel.doWrite and (self.warpType == self.assembleStaticSkyModel.warpType):
1369 raise ValueError("warpType (%s) == assembleStaticSkyModel.warpType (%s) and will compete for "
1370 "the same dataset name. Please set assembleStaticSkyModel.doWrite to False "
1371 "or warpType to 'direct'. assembleStaticSkyModel.warpType should ways be "
1372 "'PsfMatched'" % (self.warpType, self.assembleStaticSkyModel.warpType))
1373
1374
1376 """Assemble a compareWarp coadded image from a set of warps
1377 by masking artifacts detected by comparing PSF-matched warps.
1378
1379 In ``AssembleCoaddTask``, we compute the coadd as an clipped mean (i.e.,
1380 we clip outliers). The problem with doing this is that when computing the
1381 coadd PSF at a given location, individual visit PSFs from visits with
1382 outlier pixels contribute to the coadd PSF and cannot be treated correctly.
1383 In this task, we correct for this behavior by creating a new badMaskPlane
1384 'CLIPPED' which marks pixels in the individual warps suspected to contain
1385 an artifact. We populate this plane on the input warps by comparing
1386 PSF-matched warps with a PSF-matched median coadd which serves as a
1387 model of the static sky. Any group of pixels that deviates from the
1388 PSF-matched template coadd by more than config.detect.threshold sigma,
1389 is an artifact candidate. The candidates are then filtered to remove
1390 variable sources and sources that are difficult to subtract such as
1391 bright stars. This filter is configured using the config parameters
1392 ``temporalThreshold`` and ``spatialThreshold``. The temporalThreshold is
1393 the maximum fraction of epochs that the deviation can appear in and still
1394 be considered an artifact. The spatialThreshold is the maximum fraction of
1395 pixels in the footprint of the deviation that appear in other epochs
1396 (where other epochs is defined by the temporalThreshold). If the deviant
1397 region meets this criteria of having a significant percentage of pixels
1398 that deviate in only a few epochs, these pixels have the 'CLIPPED' bit
1399 set in the mask. These regions will not contribute to the final coadd.
1400 Furthermore, any routine to determine the coadd PSF can now be cognizant
1401 of clipped regions. Note that the algorithm implemented by this task is
1402 preliminary and works correctly for HSC data. Parameter modifications and
1403 or considerable redesigning of the algorithm is likley required for other
1404 surveys.
1405
1406 ``CompareWarpAssembleCoaddTask`` sub-classes
1407 ``AssembleCoaddTask`` and instantiates ``AssembleCoaddTask``
1408 as a subtask to generate the TemplateCoadd (the model of the static sky).
1409
1410 Notes
1411 -----
1412 Debugging:
1413 This task supports the following debug variables:
1414 - ``saveCountIm``
1415 If True then save the Epoch Count Image as a fits file in the `figPath`
1416 - ``figPath``
1417 Path to save the debug fits images and figures
1418 """
1419
1420 ConfigClass = CompareWarpAssembleCoaddConfig
1421 _DefaultName = "compareWarpAssembleCoadd"
1422
1423 def __init__(self, *args, **kwargs):
1424 AssembleCoaddTask.__init__(self, *args, **kwargs)
1425 self.makeSubtask("assembleStaticSkyModel")
1426 detectionSchema = afwTable.SourceTable.makeMinimalSchema()
1427 self.makeSubtask("detect", schema=detectionSchema)
1428 if self.config.doPreserveContainedBySource:
1429 self.makeSubtask("detectTemplate", schema=afwTable.SourceTable.makeMinimalSchema())
1430 if self.config.doScaleWarpVariance:
1431 self.makeSubtask("scaleWarpVariance")
1432 if self.config.doFilterMorphological:
1433 self.makeSubtask("maskStreaks")
1434
1435 @utils.inheritDoc(AssembleCoaddTask)
1436 def _makeSupplementaryData(self, butlerQC, inputRefs, outputRefs):
1437 """Generate a templateCoadd to use as a naive model of static sky to
1438 subtract from PSF-Matched warps.
1439
1440 Returns
1441 -------
1442 result : `~lsst.pipe.base.Struct`
1443 Results as a struct with attributes:
1444
1445 ``templateCoadd``
1446 Coadded exposure (`lsst.afw.image.Exposure`).
1447 ``nImage``
1448 Keeps track of exposure count for each pixel (`lsst.afw.image.ImageU`).
1449
1450 Raises
1451 ------
1452 RuntimeError
1453 Raised if ``templateCoadd`` is `None`.
1454 """
1455 # Ensure that psfMatchedWarps are used as input warps for template generation
1456 staticSkyModelInputRefs = copy.deepcopy(inputRefs)
1457 staticSkyModelInputRefs.inputWarps = inputRefs.psfMatchedWarps
1458
1459 # Because subtasks don't have connections we have to make one.
1460 # The main task's `templateCoadd` is the subtask's `coaddExposure`
1461 staticSkyModelOutputRefs = copy.deepcopy(outputRefs)
1462 if self.config.assembleStaticSkyModel.doWrite:
1463 staticSkyModelOutputRefs.coaddExposure = staticSkyModelOutputRefs.templateCoadd
1464 # Remove template coadd from both subtask's and main tasks outputs,
1465 # because it is handled by the subtask as `coaddExposure`
1466 del outputRefs.templateCoadd
1467 del staticSkyModelOutputRefs.templateCoadd
1468
1469 # A PSF-Matched nImage does not exist as a dataset type
1470 if 'nImage' in staticSkyModelOutputRefs.keys():
1471 del staticSkyModelOutputRefs.nImage
1472
1473 templateCoadd = self.assembleStaticSkyModel.runQuantum(butlerQC, staticSkyModelInputRefs,
1474 staticSkyModelOutputRefs)
1475 if templateCoadd is None:
1476 raise RuntimeError(self._noTemplateMessage(self.assembleStaticSkyModel.warpType))
1477
1478 return pipeBase.Struct(templateCoadd=templateCoadd.coaddExposure,
1479 nImage=templateCoadd.nImage,
1480 warpRefList=templateCoadd.warpRefList,
1481 imageScalerList=templateCoadd.imageScalerList,
1482 weightList=templateCoadd.weightList)
1483
1484 def _noTemplateMessage(self, warpType):
1485 warpName = (warpType[0].upper() + warpType[1:])
1486 message = """No %(warpName)s warps were found to build the template coadd which is
1487 required to run CompareWarpAssembleCoaddTask. To continue assembling this type of coadd,
1488 first either rerun makeCoaddTempExp with config.make%(warpName)s=True or
1489 coaddDriver with config.makeCoadTempExp.make%(warpName)s=True, before assembleCoadd.
1490
1491 Alternatively, to use another algorithm with existing warps, retarget the CoaddDriverConfig to
1492 another algorithm like:
1493
1494 from lsst.pipe.tasks.assembleCoadd import SafeClipAssembleCoaddTask
1495 config.assemble.retarget(SafeClipAssembleCoaddTask)
1496 """ % {"warpName": warpName}
1497 return message
1498
1499 @utils.inheritDoc(AssembleCoaddTask)
1500 @timeMethod
1501 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
1502 supplementaryData):
1503 """Assemble the coadd.
1504
1505 Find artifacts and apply them to the warps' masks creating a list of
1506 alternative masks with a new "CLIPPED" plane and updated "NO_DATA"
1507 plane. Then pass these alternative masks to the base class's ``run``
1508 method.
1509 """
1510 # Check and match the order of the supplementaryData
1511 # (PSF-matched) inputs to the order of the direct inputs,
1512 # so that the artifact mask is applied to the right warp
1513 dataIds = [ref.dataId for ref in tempExpRefList]
1514 psfMatchedDataIds = [ref.dataId for ref in supplementaryData.warpRefList]
1515
1516 if dataIds != psfMatchedDataIds:
1517 self.log.info("Reordering and or/padding PSF-matched visit input list")
1518 supplementaryData.warpRefList = reorderAndPadList(supplementaryData.warpRefList,
1519 psfMatchedDataIds, dataIds)
1520 supplementaryData.imageScalerList = reorderAndPadList(supplementaryData.imageScalerList,
1521 psfMatchedDataIds, dataIds)
1522
1523 # Use PSF-Matched Warps (and corresponding scalers) and coadd to find artifacts
1524 spanSetMaskList = self.findArtifacts(supplementaryData.templateCoadd,
1525 supplementaryData.warpRefList,
1526 supplementaryData.imageScalerList)
1527
1528 badMaskPlanes = self.config.badMaskPlanes[:]
1529 badMaskPlanes.append("CLIPPED")
1530 badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes)
1531
1532 result = AssembleCoaddTask.run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
1533 spanSetMaskList, mask=badPixelMask)
1534
1535 # Propagate PSF-matched EDGE pixels to coadd SENSOR_EDGE and INEXACT_PSF
1536 # Psf-Matching moves the real edge inwards
1537 self.applyAltEdgeMask(result.coaddExposure.maskedImage.mask, spanSetMaskList)
1538 return result
1539
1540 def applyAltEdgeMask(self, mask, altMaskList):
1541 """Propagate alt EDGE mask to SENSOR_EDGE AND INEXACT_PSF planes.
1542
1543 Parameters
1544 ----------
1545 mask : `lsst.afw.image.Mask`
1546 Original mask.
1547 altMaskList : `list` of `dict`
1548 List of Dicts containing ``spanSet`` lists.
1549 Each element contains the new mask plane name (e.g. "CLIPPED
1550 and/or "NO_DATA") as the key, and list of ``SpanSets`` to apply to
1551 the mask.
1552 """
1553 maskValue = mask.getPlaneBitMask(["SENSOR_EDGE", "INEXACT_PSF"])
1554 for visitMask in altMaskList:
1555 if "EDGE" in visitMask:
1556 for spanSet in visitMask['EDGE']:
1557 spanSet.clippedTo(mask.getBBox()).setMask(mask, maskValue)
1558
1559 def findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList):
1560 """Find artifacts.
1561
1562 Loop through warps twice. The first loop builds a map with the count
1563 of how many epochs each pixel deviates from the templateCoadd by more
1564 than ``config.chiThreshold`` sigma. The second loop takes each
1565 difference image and filters the artifacts detected in each using
1566 count map to filter out variable sources and sources that are
1567 difficult to subtract cleanly.
1568
1569 Parameters
1570 ----------
1571 templateCoadd : `lsst.afw.image.Exposure`
1572 Exposure to serve as model of static sky.
1573 tempExpRefList : `list`
1574 List of data references to warps.
1575 imageScalerList : `list`
1576 List of image scalers.
1577
1578 Returns
1579 -------
1580 altMasks : `list` of `dict`
1581 List of dicts containing information about CLIPPED
1582 (i.e., artifacts), NO_DATA, and EDGE pixels.
1583 """
1584 self.log.debug("Generating Count Image, and mask lists.")
1585 coaddBBox = templateCoadd.getBBox()
1586 slateIm = afwImage.ImageU(coaddBBox)
1587 epochCountImage = afwImage.ImageU(coaddBBox)
1588 nImage = afwImage.ImageU(coaddBBox)
1589 spanSetArtifactList = []
1590 spanSetNoDataMaskList = []
1591 spanSetEdgeList = []
1592 spanSetBadMorphoList = []
1593 badPixelMask = self.getBadPixelMask()
1594
1595 # mask of the warp diffs should = that of only the warp
1596 templateCoadd.mask.clearAllMaskPlanes()
1597
1598 if self.config.doPreserveContainedBySource:
1599 templateFootprints = self.detectTemplate.detectFootprints(templateCoadd)
1600 else:
1601 templateFootprints = None
1602
1603 for warpRef, imageScaler in zip(tempExpRefList, imageScalerList):
1604 warpDiffExp = self._readAndComputeWarpDiff(warpRef, imageScaler, templateCoadd)
1605 if warpDiffExp is not None:
1606 # This nImage only approximates the final nImage because it uses the PSF-matched mask
1607 nImage.array += (numpy.isfinite(warpDiffExp.image.array)
1608 * ((warpDiffExp.mask.array & badPixelMask) == 0)).astype(numpy.uint16)
1609 fpSet = self.detect.detectFootprints(warpDiffExp, doSmooth=False, clearMask=True)
1610 fpSet.positive.merge(fpSet.negative)
1611 footprints = fpSet.positive
1612 slateIm.set(0)
1613 spanSetList = [footprint.spans for footprint in footprints.getFootprints()]
1614
1615 # Remove artifacts due to defects before they contribute to the epochCountImage
1616 if self.config.doPrefilterArtifacts:
1617 spanSetList = self.prefilterArtifacts(spanSetList, warpDiffExp)
1618
1619 # Clear mask before adding prefiltered spanSets
1620 self.detect.clearMask(warpDiffExp.mask)
1621 for spans in spanSetList:
1622 spans.setImage(slateIm, 1, doClip=True)
1623 spans.setMask(warpDiffExp.mask, warpDiffExp.mask.getPlaneBitMask("DETECTED"))
1624 epochCountImage += slateIm
1625
1626 if self.config.doFilterMorphological:
1627 maskName = self.config.streakMaskName
1628 _ = self.maskStreaks.run(warpDiffExp)
1629 streakMask = warpDiffExp.mask
1630 spanSetStreak = afwGeom.SpanSet.fromMask(streakMask,
1631 streakMask.getPlaneBitMask(maskName)).split()
1632 # Pad the streaks to account for low-surface brightness wings
1633 psf = warpDiffExp.getPsf()
1634 for s, sset in enumerate(spanSetStreak):
1635 psfShape = psf.computeShape(sset.computeCentroid())
1636 dilation = self.config.growStreakFp * psfShape.getDeterminantRadius()
1637 sset_dilated = sset.dilated(int(dilation))
1638 spanSetStreak[s] = sset_dilated
1639
1640 # PSF-Matched warps have less available area (~the matching kernel) because the calexps
1641 # undergo a second convolution. Pixels with data in the direct warp
1642 # but not in the PSF-matched warp will not have their artifacts detected.
1643 # NaNs from the PSF-matched warp therefore must be masked in the direct warp
1644 nans = numpy.where(numpy.isnan(warpDiffExp.maskedImage.image.array), 1, 0)
1645 nansMask = afwImage.makeMaskFromArray(nans.astype(afwImage.MaskPixel))
1646 nansMask.setXY0(warpDiffExp.getXY0())
1647 edgeMask = warpDiffExp.mask
1648 spanSetEdgeMask = afwGeom.SpanSet.fromMask(edgeMask,
1649 edgeMask.getPlaneBitMask("EDGE")).split()
1650 else:
1651 # If the directWarp has <1% coverage, the psfMatchedWarp can have 0% and not exist
1652 # In this case, mask the whole epoch
1653 nansMask = afwImage.MaskX(coaddBBox, 1)
1654 spanSetList = []
1655 spanSetEdgeMask = []
1656 spanSetStreak = []
1657
1658 spanSetNoDataMask = afwGeom.SpanSet.fromMask(nansMask).split()
1659
1660 spanSetNoDataMaskList.append(spanSetNoDataMask)
1661 spanSetArtifactList.append(spanSetList)
1662 spanSetEdgeList.append(spanSetEdgeMask)
1663 if self.config.doFilterMorphological:
1664 spanSetBadMorphoList.append(spanSetStreak)
1665
1666 if lsstDebug.Info(__name__).saveCountIm:
1667 path = self._dataRef2DebugPath("epochCountIm", tempExpRefList[0], coaddLevel=True)
1668 epochCountImage.writeFits(path)
1669
1670 for i, spanSetList in enumerate(spanSetArtifactList):
1671 if spanSetList:
1672 filteredSpanSetList = self.filterArtifacts(spanSetList, epochCountImage, nImage,
1673 templateFootprints)
1674 spanSetArtifactList[i] = filteredSpanSetList
1675 if self.config.doFilterMorphological:
1676 spanSetArtifactList[i] += spanSetBadMorphoList[i]
1677
1678 altMasks = []
1679 for artifacts, noData, edge in zip(spanSetArtifactList, spanSetNoDataMaskList, spanSetEdgeList):
1680 altMasks.append({'CLIPPED': artifacts,
1681 'NO_DATA': noData,
1682 'EDGE': edge})
1683 return altMasks
1684
1685 def prefilterArtifacts(self, spanSetList, exp):
1686 """Remove artifact candidates covered by bad mask plane.
1687
1688 Any future editing of the candidate list that does not depend on
1689 temporal information should go in this method.
1690
1691 Parameters
1692 ----------
1693 spanSetList : `list` of `lsst.afw.geom.SpanSet`
1694 List of SpanSets representing artifact candidates.
1696 Exposure containing mask planes used to prefilter.
1697
1698 Returns
1699 -------
1700 returnSpanSetList : `list` of `lsst.afw.geom.SpanSet`
1701 List of SpanSets with artifacts.
1702 """
1703 badPixelMask = exp.mask.getPlaneBitMask(self.config.prefilterArtifactsMaskPlanes)
1704 goodArr = (exp.mask.array & badPixelMask) == 0
1705 returnSpanSetList = []
1706 bbox = exp.getBBox()
1707 x0, y0 = exp.getXY0()
1708 for i, span in enumerate(spanSetList):
1709 y, x = span.clippedTo(bbox).indices()
1710 yIndexLocal = numpy.array(y) - y0
1711 xIndexLocal = numpy.array(x) - x0
1712 goodRatio = numpy.count_nonzero(goodArr[yIndexLocal, xIndexLocal])/span.getArea()
1713 if goodRatio > self.config.prefilterArtifactsRatio:
1714 returnSpanSetList.append(span)
1715 return returnSpanSetList
1716
1717 def filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None):
1718 """Filter artifact candidates.
1719
1720 Parameters
1721 ----------
1722 spanSetList : `list` of `lsst.afw.geom.SpanSet`
1723 List of SpanSets representing artifact candidates.
1724 epochCountImage : `lsst.afw.image.Image`
1725 Image of accumulated number of warpDiff detections.
1726 nImage : `lsst.afw.image.ImageU`
1727 Image of the accumulated number of total epochs contributing.
1728
1729 Returns
1730 -------
1731 maskSpanSetList : `list`
1732 List of SpanSets with artifacts.
1733 """
1734 maskSpanSetList = []
1735 x0, y0 = epochCountImage.getXY0()
1736 for i, span in enumerate(spanSetList):
1737 y, x = span.indices()
1738 yIdxLocal = [y1 - y0 for y1 in y]
1739 xIdxLocal = [x1 - x0 for x1 in x]
1740 outlierN = epochCountImage.array[yIdxLocal, xIdxLocal]
1741 totalN = nImage.array[yIdxLocal, xIdxLocal]
1742
1743 # effectiveMaxNumEpochs is broken line (fraction of N) with characteristic config.maxNumEpochs
1744 effMaxNumEpochsHighN = (self.config.maxNumEpochs
1745 + self.config.maxFractionEpochsHigh*numpy.mean(totalN))
1746 effMaxNumEpochsLowN = self.config.maxFractionEpochsLow * numpy.mean(totalN)
1747 effectiveMaxNumEpochs = int(min(effMaxNumEpochsLowN, effMaxNumEpochsHighN))
1748 nPixelsBelowThreshold = numpy.count_nonzero((outlierN > 0)
1749 & (outlierN <= effectiveMaxNumEpochs))
1750 percentBelowThreshold = nPixelsBelowThreshold / len(outlierN)
1751 if percentBelowThreshold > self.config.spatialThreshold:
1752 maskSpanSetList.append(span)
1753
1754 if self.config.doPreserveContainedBySource and footprintsToExclude is not None:
1755 # If a candidate is contained by a footprint on the template coadd, do not clip
1756 filteredMaskSpanSetList = []
1757 for span in maskSpanSetList:
1758 doKeep = True
1759 for footprint in footprintsToExclude.positive.getFootprints():
1760 if footprint.spans.contains(span):
1761 doKeep = False
1762 break
1763 if doKeep:
1764 filteredMaskSpanSetList.append(span)
1765 maskSpanSetList = filteredMaskSpanSetList
1766
1767 return maskSpanSetList
1768
1769 def _readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd):
1770 """Fetch a warp from the butler and return a warpDiff.
1771
1772 Parameters
1773 ----------
1774 warpRef : `lsst.daf.butler.DeferredDatasetHandle`
1775 Handle for the warp.
1777 An image scaler object.
1778 templateCoadd : `lsst.afw.image.Exposure`
1779 Exposure to be substracted from the scaled warp.
1780
1781 Returns
1782 -------
1784 Exposure of the image difference between the warp and template.
1785 """
1786 # If the PSF-Matched warp did not exist for this direct warp
1787 # None is holding its place to maintain order in Gen 3
1788 if warpRef is None:
1789 return None
1790
1791 warp = warpRef.get()
1792 # direct image scaler OK for PSF-matched Warp
1793 imageScaler.scaleMaskedImage(warp.getMaskedImage())
1794 mi = warp.getMaskedImage()
1795 if self.config.doScaleWarpVariance:
1796 try:
1797 self.scaleWarpVariance.run(mi)
1798 except Exception as exc:
1799 self.log.warning("Unable to rescale variance of warp (%s); leaving it as-is", exc)
1800 mi -= templateCoadd.getMaskedImage()
1801 return warp
int min
int max
Class to describe the properties of a detected object from an image.
Definition Footprint.h:63
A compact representation of a collection of pixels.
Definition SpanSet.h:78
A class to contain the data, WCS, and other information needed to describe an image of the sky.
Definition Exposure.h:72
A group of labels for a filter in an exposure or coadd.
Definition FilterLabel.h:58
A class to represent a 2-dimensional array of pixels.
Definition Image.h:51
Represent a 2-dimensional array of bitmask pixels.
Definition Mask.h:77
A class to manipulate images, masks, and variance as a single object.
Definition MaskedImage.h:74
Pass parameters to a Statistics object.
Definition Statistics.h:83
A floating-point coordinate rectangle geometry.
Definition Box.h:413
An integer coordinate rectangle.
Definition Box.h:55
Reports invalid arguments.
Definition Runtime.h:66
runQuantum(self, butlerQC, inputRefs, outputRefs)
assembleMetadata(self, coaddExposure, tempExpRefList, weightList)
assembleOnlineMeanCoadd(self, coaddExposure, tempExpRefList, imageScalerList, weightList, altMaskList, statsCtrl, nImage=None)
processResults(self, coaddExposure, brightObjectMasks=None, dataId=None)
assembleSubregion(self, coaddExposure, bbox, tempExpRefList, imageScalerList, weightList, altMaskList, statsFlags, statsCtrl, nImage=None)
setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None)
_makeSupplementaryData(self, butlerQC, inputRefs, outputRefs)
makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs)
run(self, skyInfo, tempExpRefList, imageScalerList, weightList, altMaskList=None, mask=None, supplementaryData=None)
_readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd)
_makeSupplementaryData(self, butlerQC, inputRefs, outputRefs)
findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList)
filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None)
getTempExpDatasetName(self, warpType="direct")
Definition coaddBase.py:128
std::shared_ptr< lsst::afw::image::Image< PixelT > > statisticsStack(std::vector< std::shared_ptr< lsst::afw::image::Image< PixelT > > > &images, Property flags, StatisticsControl const &sctrl=StatisticsControl(), std::vector< lsst::afw::image::VariancePixel > const &wvector=std::vector< lsst::afw::image::VariancePixel >(0))
A function to compute some statistics of a stack of Images.
Statistics makeStatistics(lsst::afw::image::Image< Pixel > const &img, lsst::afw::image::Mask< image::MaskPixel > const &msk, int const flags, StatisticsControl const &sctrl=StatisticsControl())
Handle a watered-down front-end to the constructor (no variance)
Definition Statistics.h:361
Property stringToStatisticsProperty(std::string const property)
Conversion function to switch a string to a Property (see Statistics.h)
void setCoaddEdgeBits(lsst::afw::image::Mask< lsst::afw::image::MaskPixel > &coaddMask, lsst::afw::image::Image< WeightPixelT > const &weightMap)
set edge bits of coadd mask based on weight map
countMaskFromFootprint(mask, footprint, bitmask, ignoreMask)