LSST Applications g0b6bd0c080+a72a5dd7e6,g1182afd7b4+2a019aa3bb,g17e5ecfddb+2b8207f7de,g1d67935e3f+06cf436103,g38293774b4+ac198e9f13,g396055baef+6a2097e274,g3b44f30a73+6611e0205b,g480783c3b1+98f8679e14,g48ccf36440+89c08d0516,g4b93dc025c+98f8679e14,g5c4744a4d9+a302e8c7f0,g613e996a0d+e1c447f2e0,g6c8d09e9e7+25247a063c,g7271f0639c+98f8679e14,g7a9cd813b8+124095ede6,g9d27549199+a302e8c7f0,ga1cf026fa3+ac198e9f13,ga32aa97882+7403ac30ac,ga786bb30fb+7a139211af,gaa63f70f4e+9994eb9896,gabf319e997+ade567573c,gba47b54d5d+94dc90c3ea,gbec6a3398f+06cf436103,gc6308e37c7+07dd123edb,gc655b1545f+ade567573c,gcc9029db3c+ab229f5caf,gd01420fc67+06cf436103,gd877ba84e5+06cf436103,gdb4cecd868+6f279b5b48,ge2d134c3d5+cc4dbb2e3f,ge448b5faa6+86d1ceac1d,gecc7e12556+98f8679e14,gf3ee170dca+25247a063c,gf4ac96e456+ade567573c,gf9f5ea5b4d+ac198e9f13,gff490e6085+8c2580be5c,w.2022.27
LSST Data Management Base Package
assembleCoadd.py
Go to the documentation of this file.
1# This file is part of pipe_tasks.
2#
3# LSST Data Management System
4# This product includes software developed by the
5# LSST Project (http://www.lsst.org/).
6# See COPYRIGHT file at the top of the source tree.
7#
8# This program is free software: you can redistribute it and/or modify
9# it under the terms of the GNU General Public License as published by
10# the Free Software Foundation, either version 3 of the License, or
11# (at your option) any later version.
12#
13# This program is distributed in the hope that it will be useful,
14# but WITHOUT ANY WARRANTY; without even the implied warranty of
15# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16# GNU General Public License for more details.
17#
18# You should have received a copy of the LSST License Statement and
19# the GNU General Public License along with this program. If not,
20# see <https://www.lsstcorp.org/LegalNotices/>.
21#
22import os
23import copy
24import numpy
25import warnings
26import logging
27import lsst.pex.config as pexConfig
28import lsst.pex.exceptions as pexExceptions
29import lsst.geom as geom
30import lsst.afw.geom as afwGeom
31import lsst.afw.image as afwImage
32import lsst.afw.math as afwMath
33import lsst.afw.table as afwTable
34import lsst.afw.detection as afwDet
35import lsst.coadd.utils as coaddUtils
36import lsst.pipe.base as pipeBase
37import lsst.meas.algorithms as measAlg
38import lsstDebug
39import lsst.utils as utils
40from lsst.skymap import BaseSkyMap
41from .coaddBase import CoaddBaseTask, SelectDataIdContainer, makeSkyInfo, makeCoaddSuffix, reorderAndPadList
42from .interpImage import InterpImageTask
43from .scaleZeroPoint import ScaleZeroPointTask
44from .coaddHelpers import groupPatchExposures, getGroupDataRef
45from .maskStreaks import MaskStreaksTask
46from .healSparseMapping import HealSparseInputMapTask
47from lsst.meas.algorithms import SourceDetectionTask, AccumulatorMeanStack, ScaleVarianceTask
48from lsst.daf.butler import DeferredDatasetHandle
49from lsst.utils.timer import timeMethod
50
51__all__ = ["AssembleCoaddTask", "AssembleCoaddConnections", "AssembleCoaddConfig",
52 "SafeClipAssembleCoaddTask", "SafeClipAssembleCoaddConfig",
53 "CompareWarpAssembleCoaddTask", "CompareWarpAssembleCoaddConfig"]
54
55log = logging.getLogger(__name__)
56
57
58class AssembleCoaddConnections(pipeBase.PipelineTaskConnections,
59 dimensions=("tract", "patch", "band", "skymap"),
60 defaultTemplates={"inputCoaddName": "deep",
61 "outputCoaddName": "deep",
62 "warpType": "direct",
63 "warpTypeSuffix": ""}):
64
65 inputWarps = pipeBase.connectionTypes.Input(
66 doc=("Input list of warps to be assemebled i.e. stacked."
67 "WarpType (e.g. direct, psfMatched) is controlled by the warpType config parameter"),
68 name="{inputCoaddName}Coadd_{warpType}Warp",
69 storageClass="ExposureF",
70 dimensions=("tract", "patch", "skymap", "visit", "instrument"),
71 deferLoad=True,
72 multiple=True
73 )
74 skyMap = pipeBase.connectionTypes.Input(
75 doc="Input definition of geometry/bbox and projection/wcs for coadded exposures",
76 name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME,
77 storageClass="SkyMap",
78 dimensions=("skymap", ),
79 )
80 selectedVisits = pipeBase.connectionTypes.Input(
81 doc="Selected visits to be coadded.",
82 name="{outputCoaddName}Visits",
83 storageClass="StructuredDataDict",
84 dimensions=("instrument", "tract", "patch", "skymap", "band")
85 )
86 brightObjectMask = pipeBase.connectionTypes.PrerequisiteInput(
87 doc=("Input Bright Object Mask mask produced with external catalogs to be applied to the mask plane"
88 " BRIGHT_OBJECT."),
89 name="brightObjectMask",
90 storageClass="ObjectMaskCatalog",
91 dimensions=("tract", "patch", "skymap", "band"),
92 )
93 coaddExposure = pipeBase.connectionTypes.Output(
94 doc="Output coadded exposure, produced by stacking input warps",
95 name="{outputCoaddName}Coadd{warpTypeSuffix}",
96 storageClass="ExposureF",
97 dimensions=("tract", "patch", "skymap", "band"),
98 )
99 nImage = pipeBase.connectionTypes.Output(
100 doc="Output image of number of input images per pixel",
101 name="{outputCoaddName}Coadd_nImage",
102 storageClass="ImageU",
103 dimensions=("tract", "patch", "skymap", "band"),
104 )
105 inputMap = pipeBase.connectionTypes.Output(
106 doc="Output healsparse map of input images",
107 name="{outputCoaddName}Coadd_inputMap",
108 storageClass="HealSparseMap",
109 dimensions=("tract", "patch", "skymap", "band"),
110 )
111
112 def __init__(self, *, config=None):
113 super().__init__(config=config)
114
115 # Override the connection's name template with config to replicate Gen2 behavior
116 # This duplicates some of the logic in the base class, due to wanting Gen2 and
117 # Gen3 configs to stay in sync. This should be removed when gen2 is deprecated
118 templateValues = {name: getattr(config.connections, name) for name in self.defaultTemplates}
119 templateValues['warpType'] = config.warpType
120 templateValues['warpTypeSuffix'] = makeCoaddSuffix(config.warpType)
121 self._nameOverrides = {name: getattr(config.connections, name).format(**templateValues)
122 for name in self.allConnections}
123 self._typeNameToVarName = {v: k for k, v in self._nameOverrides.items()}
124 # End code to remove after deprecation
125
126 if not config.doMaskBrightObjects:
127 self.prerequisiteInputs.remove("brightObjectMask")
128
129 if not config.doSelectVisits:
130 self.inputs.remove("selectedVisits")
131
132 if not config.doNImage:
133 self.outputs.remove("nImage")
134
135 if not self.config.doInputMap:
136 self.outputs.remove("inputMap")
137
138
139class AssembleCoaddConfig(CoaddBaseTask.ConfigClass, pipeBase.PipelineTaskConfig,
140 pipelineConnections=AssembleCoaddConnections):
141 """Configuration parameters for the `AssembleCoaddTask`.
142
143 Notes
144 -----
145 The `doMaskBrightObjects` and `brightObjectMaskName` configuration options
146 only set the bitplane config.brightObjectMaskName. To make this useful you
147 *must* also configure the flags.pixel algorithm, for example by adding
148
149 .. code-block:: none
150
151 config.measurement.plugins["base_PixelFlags"].masksFpCenter.append("BRIGHT_OBJECT")
152 config.measurement.plugins["base_PixelFlags"].masksFpAnywhere.append("BRIGHT_OBJECT")
153
154 to your measureCoaddSources.py and forcedPhotCoadd.py config overrides.
155 """
156 warpType = pexConfig.Field(
157 doc="Warp name: one of 'direct' or 'psfMatched'",
158 dtype=str,
159 default="direct",
160 )
161 subregionSize = pexConfig.ListField(
162 dtype=int,
163 doc="Width, height of stack subregion size; "
164 "make small enough that a full stack of images will fit into memory at once.",
165 length=2,
166 default=(2000, 2000),
167 )
168 statistic = pexConfig.Field(
169 dtype=str,
170 doc="Main stacking statistic for aggregating over the epochs.",
171 default="MEANCLIP",
172 )
173 doOnlineForMean = pexConfig.Field(
174 dtype=bool,
175 doc="Perform online coaddition when statistic=\"MEAN\" to save memory?",
176 default=False,
177 )
178 doSigmaClip = pexConfig.Field(
179 dtype=bool,
180 doc="Perform sigma clipped outlier rejection with MEANCLIP statistic? (DEPRECATED)",
181 default=False,
182 )
183 sigmaClip = pexConfig.Field(
184 dtype=float,
185 doc="Sigma for outlier rejection; ignored if non-clipping statistic selected.",
186 default=3.0,
187 )
188 clipIter = pexConfig.Field(
189 dtype=int,
190 doc="Number of iterations of outlier rejection; ignored if non-clipping statistic selected.",
191 default=2,
192 )
193 calcErrorFromInputVariance = pexConfig.Field(
194 dtype=bool,
195 doc="Calculate coadd variance from input variance by stacking statistic."
196 "Passed to StatisticsControl.setCalcErrorFromInputVariance()",
197 default=True,
198 )
199 scaleZeroPoint = pexConfig.ConfigurableField(
200 target=ScaleZeroPointTask,
201 doc="Task to adjust the photometric zero point of the coadd temp exposures",
202 )
203 doInterp = pexConfig.Field(
204 doc="Interpolate over NaN pixels? Also extrapolate, if necessary, but the results are ugly.",
205 dtype=bool,
206 default=True,
207 )
208 interpImage = pexConfig.ConfigurableField(
209 target=InterpImageTask,
210 doc="Task to interpolate (and extrapolate) over NaN pixels",
211 )
212 doWrite = pexConfig.Field(
213 doc="Persist coadd?",
214 dtype=bool,
215 default=True,
216 )
217 doNImage = pexConfig.Field(
218 doc="Create image of number of contributing exposures for each pixel",
219 dtype=bool,
220 default=False,
221 )
222 doUsePsfMatchedPolygons = pexConfig.Field(
223 doc="Use ValidPolygons from shrunk Psf-Matched Calexps? Should be set to True by CompareWarp only.",
224 dtype=bool,
225 default=False,
226 )
227 maskPropagationThresholds = pexConfig.DictField(
228 keytype=str,
229 itemtype=float,
230 doc=("Threshold (in fractional weight) of rejection at which we propagate a mask plane to "
231 "the coadd; that is, we set the mask bit on the coadd if the fraction the rejected frames "
232 "would have contributed exceeds this value."),
233 default={"SAT": 0.1},
234 )
235 removeMaskPlanes = pexConfig.ListField(dtype=str, default=["NOT_DEBLENDED"],
236 doc="Mask planes to remove before coadding")
237 doMaskBrightObjects = pexConfig.Field(dtype=bool, default=False,
238 doc="Set mask and flag bits for bright objects?")
239 brightObjectMaskName = pexConfig.Field(dtype=str, default="BRIGHT_OBJECT",
240 doc="Name of mask bit used for bright objects")
241 coaddPsf = pexConfig.ConfigField(
242 doc="Configuration for CoaddPsf",
243 dtype=measAlg.CoaddPsfConfig,
244 )
245 doAttachTransmissionCurve = pexConfig.Field(
246 dtype=bool, default=False, optional=False,
247 doc=("Attach a piecewise TransmissionCurve for the coadd? "
248 "(requires all input Exposures to have TransmissionCurves).")
249 )
250 hasFakes = pexConfig.Field(
251 dtype=bool,
252 default=False,
253 doc="Should be set to True if fake sources have been inserted into the input data."
254 )
255 doSelectVisits = pexConfig.Field(
256 doc="Coadd only visits selected by a SelectVisitsTask",
257 dtype=bool,
258 default=False,
259 )
260 doInputMap = pexConfig.Field(
261 doc="Create a bitwise map of coadd inputs",
262 dtype=bool,
263 default=False,
264 )
265 inputMapper = pexConfig.ConfigurableField(
266 doc="Input map creation subtask.",
267 target=HealSparseInputMapTask,
268 )
269
270 def setDefaults(self):
271 super().setDefaults()
272 self.badMaskPlanes = ["NO_DATA", "BAD", "SAT", "EDGE"]
273
274 def validate(self):
275 super().validate()
276 if self.doPsfMatch:
277 # Backwards compatibility.
278 # Configs do not have loggers
279 log.warning("Config doPsfMatch deprecated. Setting warpType='psfMatched'")
280 self.warpType = 'psfMatched'
281 if self.doSigmaClip and self.statistic != "MEANCLIP":
282 log.warning('doSigmaClip deprecated. To replicate behavior, setting statistic to "MEANCLIP"')
283 self.statistic = "MEANCLIP"
284 if self.doInterp and self.statistic not in ['MEAN', 'MEDIAN', 'MEANCLIP', 'VARIANCE', 'VARIANCECLIP']:
285 raise ValueError("Must set doInterp=False for statistic=%s, which does not "
286 "compute and set a non-zero coadd variance estimate." % (self.statistic))
287
288 unstackableStats = ['NOTHING', 'ERROR', 'ORMASK']
289 if not hasattr(afwMath.Property, self.statistic) or self.statistic in unstackableStats:
290 stackableStats = [str(k) for k in afwMath.Property.__members__.keys()
291 if str(k) not in unstackableStats]
292 raise ValueError("statistic %s is not allowed. Please choose one of %s."
293 % (self.statistic, stackableStats))
294
295
296class AssembleCoaddTask(CoaddBaseTask, pipeBase.PipelineTask):
297 """Assemble a coadded image from a set of warps (coadded temporary exposures).
298
299 We want to assemble a coadded image from a set of Warps (also called
300 coadded temporary exposures or ``coaddTempExps``).
301 Each input Warp covers a patch on the sky and corresponds to a single
302 run/visit/exposure of the covered patch. We provide the task with a list
303 of Warps (``selectDataList``) from which it selects Warps that cover the
304 specified patch (pointed at by ``dataRef``).
305 Each Warp that goes into a coadd will typically have an independent
306 photometric zero-point. Therefore, we must scale each Warp to set it to
307 a common photometric zeropoint. WarpType may be one of 'direct' or
308 'psfMatched', and the boolean configs `config.makeDirect` and
309 `config.makePsfMatched` set which of the warp types will be coadded.
310 The coadd is computed as a mean with optional outlier rejection.
311 Criteria for outlier rejection are set in `AssembleCoaddConfig`.
312 Finally, Warps can have bad 'NaN' pixels which received no input from the
313 source calExps. We interpolate over these bad (NaN) pixels.
314
315 `AssembleCoaddTask` uses several sub-tasks. These are
316
317 - `ScaleZeroPointTask`
318 - create and use an ``imageScaler`` object to scale the photometric zeropoint for each Warp
319 - `InterpImageTask`
320 - interpolate across bad pixels (NaN) in the final coadd
321
322 You can retarget these subtasks if you wish.
323
324 Notes
325 -----
326 The `lsst.pipe.base.cmdLineTask.CmdLineTask` interface supports a
327 flag ``-d`` to import ``debug.py`` from your ``PYTHONPATH``; see
328 `baseDebug` for more about ``debug.py`` files. `AssembleCoaddTask` has
329 no debug variables of its own. Some of the subtasks may support debug
330 variables. See the documentation for the subtasks for further information.
331
332 Examples
333 --------
334 `AssembleCoaddTask` assembles a set of warped images into a coadded image.
335 The `AssembleCoaddTask` can be invoked by running ``assembleCoadd.py``
336 with the flag '--legacyCoadd'. Usage of assembleCoadd.py expects two
337 inputs: a data reference to the tract patch and filter to be coadded, and
338 a list of Warps to attempt to coadd. These are specified using ``--id`` and
339 ``--selectId``, respectively:
340
341 .. code-block:: none
342
343 --id = [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]
344 --selectId [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]
345
346 Only the Warps that cover the specified tract and patch will be coadded.
347 A list of the available optional arguments can be obtained by calling
348 ``assembleCoadd.py`` with the ``--help`` command line argument:
349
350 .. code-block:: none
351
352 assembleCoadd.py --help
353
354 To demonstrate usage of the `AssembleCoaddTask` in the larger context of
355 multi-band processing, we will generate the HSC-I & -R band coadds from
356 HSC engineering test data provided in the ``ci_hsc`` package. To begin,
357 assuming that the lsst stack has been already set up, we must set up the
358 obs_subaru and ``ci_hsc`` packages. This defines the environment variable
359 ``$CI_HSC_DIR`` and points at the location of the package. The raw HSC
360 data live in the ``$CI_HSC_DIR/raw directory``. To begin assembling the
361 coadds, we must first
362
363 - processCcd
364 - process the individual ccds in $CI_HSC_RAW to produce calibrated exposures
365 - makeSkyMap
366 - create a skymap that covers the area of the sky present in the raw exposures
367 - makeCoaddTempExp
368 - warp the individual calibrated exposures to the tangent plane of the coadd
369
370 We can perform all of these steps by running
371
372 .. code-block:: none
373
374 $CI_HSC_DIR scons warp-903986 warp-904014 warp-903990 warp-904010 warp-903988
375
376 This will produce warped exposures for each visit. To coadd the warped
377 data, we call assembleCoadd.py as follows:
378
379 .. code-block:: none
380
381 assembleCoadd.py --legacyCoadd $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-I \
382 --selectId visit=903986 ccd=16 --selectId visit=903986 ccd=22 --selectId visit=903986 ccd=23 \
383 --selectId visit=903986 ccd=100 --selectId visit=904014 ccd=1 --selectId visit=904014 ccd=6 \
384 --selectId visit=904014 ccd=12 --selectId visit=903990 ccd=18 --selectId visit=903990 ccd=25 \
385 --selectId visit=904010 ccd=4 --selectId visit=904010 ccd=10 --selectId visit=904010 ccd=100 \
386 --selectId visit=903988 ccd=16 --selectId visit=903988 ccd=17 --selectId visit=903988 ccd=23 \
387 --selectId visit=903988 ccd=24
388
389 that will process the HSC-I band data. The results are written in
390 ``$CI_HSC_DIR/DATA/deepCoadd-results/HSC-I``.
391
392 You may also choose to run:
393
394 .. code-block:: none
395
396 scons warp-903334 warp-903336 warp-903338 warp-903342 warp-903344 warp-903346
397 assembleCoadd.py --legacyCoadd $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-R \
398 --selectId visit=903334 ccd=16 --selectId visit=903334 ccd=22 --selectId visit=903334 ccd=23 \
399 --selectId visit=903334 ccd=100 --selectId visit=903336 ccd=17 --selectId visit=903336 ccd=24 \
400 --selectId visit=903338 ccd=18 --selectId visit=903338 ccd=25 --selectId visit=903342 ccd=4 \
401 --selectId visit=903342 ccd=10 --selectId visit=903342 ccd=100 --selectId visit=903344 ccd=0 \
402 --selectId visit=903344 ccd=5 --selectId visit=903344 ccd=11 --selectId visit=903346 ccd=1 \
403 --selectId visit=903346 ccd=6 --selectId visit=903346 ccd=12
404
405 to generate the coadd for the HSC-R band if you are interested in
406 following multiBand Coadd processing as discussed in `pipeTasks_multiBand`
407 (but note that normally, one would use the `SafeClipAssembleCoaddTask`
408 rather than `AssembleCoaddTask` to make the coadd.
409 """
410 ConfigClass = AssembleCoaddConfig
411 _DefaultName = "assembleCoadd"
412
413 def __init__(self, *args, **kwargs):
414 # TODO: DM-17415 better way to handle previously allowed passed args e.g.`AssembleCoaddTask(config)`
415 if args:
416 argNames = ["config", "name", "parentTask", "log"]
417 kwargs.update({k: v for k, v in zip(argNames, args)})
418 warnings.warn("AssembleCoadd received positional args, and casting them as kwargs: %s. "
419 "PipelineTask will not take positional args" % argNames, FutureWarning)
420
421 super().__init__(**kwargs)
422 self.makeSubtask("interpImage")
423 self.makeSubtask("scaleZeroPoint")
424
425 if self.config.doMaskBrightObjects:
426 mask = afwImage.Mask()
427 try:
428 self.brightObjectBitmask = 1 << mask.addMaskPlane(self.config.brightObjectMaskName)
429 except pexExceptions.LsstCppException:
430 raise RuntimeError("Unable to define mask plane for bright objects; planes used are %s" %
431 mask.getMaskPlaneDict().keys())
432 del mask
433
434 if self.config.doInputMap:
435 self.makeSubtask("inputMapper")
436
437 self.warpType = self.config.warpType
438
439 @utils.inheritDoc(pipeBase.PipelineTask)
440 def runQuantum(self, butlerQC, inputRefs, outputRefs):
441 # Docstring to be formatted with info from PipelineTask.runQuantum
442 """
443 Notes
444 -----
445 Assemble a coadd from a set of Warps.
446
447 PipelineTask (Gen3) entry point to Coadd a set of Warps.
448 Analogous to `runDataRef`, it prepares all the data products to be
449 passed to `run`, and processes the results before returning a struct
450 of results to be written out. AssembleCoadd cannot fit all Warps in memory.
451 Therefore, its inputs are accessed subregion by subregion
452 by the Gen3 `DeferredDatasetHandle` that is analagous to the Gen2
453 `lsst.daf.persistence.ButlerDataRef`. Any updates to this method should
454 correspond to an update in `runDataRef` while both entry points
455 are used.
456 """
457 inputData = butlerQC.get(inputRefs)
458
459 # Construct skyInfo expected by run
460 # Do not remove skyMap from inputData in case makeSupplementaryDataGen3 needs it
461 skyMap = inputData["skyMap"]
462 outputDataId = butlerQC.quantum.dataId
463
464 inputData['skyInfo'] = makeSkyInfo(skyMap,
465 tractId=outputDataId['tract'],
466 patchId=outputDataId['patch'])
467
468 if self.config.doSelectVisits:
469 warpRefList = self.filterWarps(inputData['inputWarps'], inputData['selectedVisits'])
470 else:
471 warpRefList = inputData['inputWarps']
472
473 # Perform same middle steps as `runDataRef` does
474 inputs = self.prepareInputs(warpRefList)
475 self.log.info("Found %d %s", len(inputs.tempExpRefList),
476 self.getTempExpDatasetName(self.warpType))
477 if len(inputs.tempExpRefList) == 0:
478 raise pipeBase.NoWorkFound("No coadd temporary exposures found")
479
480 supplementaryData = self.makeSupplementaryDataGen3(butlerQC, inputRefs, outputRefs)
481 retStruct = self.run(inputData['skyInfo'], inputs.tempExpRefList, inputs.imageScalerList,
482 inputs.weightList, supplementaryData=supplementaryData)
483
484 inputData.setdefault('brightObjectMask', None)
485 self.processResults(retStruct.coaddExposure, inputData['brightObjectMask'], outputDataId)
486
487 if self.config.doWrite:
488 butlerQC.put(retStruct, outputRefs)
489 return retStruct
490
491 @timeMethod
492 def runDataRef(self, dataRef, selectDataList=None, warpRefList=None):
493 """Assemble a coadd from a set of Warps.
494
495 Pipebase.CmdlineTask entry point to Coadd a set of Warps.
496 Compute weights to be applied to each Warp and
497 find scalings to match the photometric zeropoint to a reference Warp.
498 Assemble the Warps using `run`. Interpolate over NaNs and
499 optionally write the coadd to disk. Return the coadded exposure.
500
501 Parameters
502 ----------
504 Data reference defining the patch for coaddition and the
505 reference Warp (if ``config.autoReference=False``).
506 Used to access the following data products:
507 - ``self.config.coaddName + "Coadd_skyMap"``
508 - ``self.config.coaddName + "Coadd_ + <warpType> + "Warp"`` (optionally)
509 - ``self.config.coaddName + "Coadd"``
510 selectDataList : `list`
511 List of data references to Calexps. Data to be coadded will be
512 selected from this list based on overlap with the patch defined
513 by dataRef, grouped by visit, and converted to a list of data
514 references to warps.
515 warpRefList : `list`
516 List of data references to Warps to be coadded.
517 Note: `warpRefList` is just the new name for `tempExpRefList`.
518
519 Returns
520 -------
521 retStruct : `lsst.pipe.base.Struct`
522 Result struct with components:
523
524 - ``coaddExposure``: coadded exposure (``Exposure``).
525 - ``nImage``: exposure count image (``Image``).
526 """
527 if selectDataList and warpRefList:
528 raise RuntimeError("runDataRef received both a selectDataList and warpRefList, "
529 "and which to use is ambiguous. Please pass only one.")
530
531 skyInfo = self.getSkyInfo(dataRef)
532 if warpRefList is None:
533 calExpRefList = self.selectExposures(dataRef, skyInfo, selectDataList=selectDataList)
534 if len(calExpRefList) == 0:
535 self.log.warning("No exposures to coadd")
536 return
537 self.log.info("Coadding %d exposures", len(calExpRefList))
538
539 warpRefList = self.getTempExpRefList(dataRef, calExpRefList)
540
541 inputData = self.prepareInputs(warpRefList)
542 self.log.info("Found %d %s", len(inputData.tempExpRefList),
543 self.getTempExpDatasetName(self.warpType))
544 if len(inputData.tempExpRefList) == 0:
545 self.log.warning("No coadd temporary exposures found")
546 return
547
548 supplementaryData = self.makeSupplementaryData(dataRef, warpRefList=inputData.tempExpRefList)
549
550 retStruct = self.run(skyInfo, inputData.tempExpRefList, inputData.imageScalerList,
551 inputData.weightList, supplementaryData=supplementaryData)
552
553 brightObjects = self.readBrightObjectMasks(dataRef) if self.config.doMaskBrightObjects else None
554 self.processResults(retStruct.coaddExposure, brightObjectMasks=brightObjects, dataId=dataRef.dataId)
555
556 if self.config.doWrite:
557 if self.getCoaddDatasetName(self.warpType) == "deepCoadd" and self.config.hasFakes:
558 coaddDatasetName = "fakes_" + self.getCoaddDatasetName(self.warpType)
559 else:
560 coaddDatasetName = self.getCoaddDatasetName(self.warpType)
561 self.log.info("Persisting %s", coaddDatasetName)
562 dataRef.put(retStruct.coaddExposure, coaddDatasetName)
563 if self.config.doNImage and retStruct.nImage is not None:
564 dataRef.put(retStruct.nImage, self.getCoaddDatasetName(self.warpType) + '_nImage')
565
566 return retStruct
567
568 def processResults(self, coaddExposure, brightObjectMasks=None, dataId=None):
569 """Interpolate over missing data and mask bright stars.
570
571 Parameters
572 ----------
573 coaddExposure : `lsst.afw.image.Exposure`
574 The coadded exposure to process.
576 Butler data reference for supplementary data.
577 """
578 if self.config.doInterp:
579 self.interpImage.run(coaddExposure.getMaskedImage(), planeName="NO_DATA")
580 # The variance must be positive; work around for DM-3201.
581 varArray = coaddExposure.variance.array
582 with numpy.errstate(invalid="ignore"):
583 varArray[:] = numpy.where(varArray > 0, varArray, numpy.inf)
584
585 if self.config.doMaskBrightObjects:
586 self.setBrightObjectMasks(coaddExposure, brightObjectMasks, dataId)
587
588 def makeSupplementaryData(self, dataRef, selectDataList=None, warpRefList=None):
589 """Make additional inputs to run() specific to subclasses (Gen2)
590
591 Duplicates interface of `runDataRef` method
592 Available to be implemented by subclasses only if they need the
593 coadd dataRef for performing preliminary processing before
594 assembling the coadd.
595
596 Parameters
597 ----------
599 Butler data reference for supplementary data.
600 selectDataList : `list` (optional)
601 Optional List of data references to Calexps.
602 warpRefList : `list` (optional)
603 Optional List of data references to Warps.
604 """
605 return pipeBase.Struct()
606
607 def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs):
608 """Make additional inputs to run() specific to subclasses (Gen3)
609
610 Duplicates interface of `runQuantum` method.
611 Available to be implemented by subclasses only if they need the
612 coadd dataRef for performing preliminary processing before
613 assembling the coadd.
614
615 Parameters
616 ----------
617 butlerQC : `lsst.pipe.base.ButlerQuantumContext`
618 Gen3 Butler object for fetching additional data products before
619 running the Task specialized for quantum being processed
620 inputRefs : `lsst.pipe.base.InputQuantizedConnection`
621 Attributes are the names of the connections describing input dataset types.
622 Values are DatasetRefs that task consumes for corresponding dataset type.
623 DataIds are guaranteed to match data objects in ``inputData``.
624 outputRefs : `lsst.pipe.base.OutputQuantizedConnection`
625 Attributes are the names of the connections describing output dataset types.
626 Values are DatasetRefs that task is to produce
627 for corresponding dataset type.
628 """
629 return pipeBase.Struct()
630
631 def getTempExpRefList(self, patchRef, calExpRefList):
632 """Generate list data references corresponding to warped exposures
633 that lie within the patch to be coadded.
634
635 Parameters
636 ----------
637 patchRef : `dataRef`
638 Data reference for patch.
639 calExpRefList : `list`
640 List of data references for input calexps.
641
642 Returns
643 -------
644 tempExpRefList : `list`
645 List of Warp/CoaddTempExp data references.
646 """
647 butler = patchRef.getButler()
648 groupData = groupPatchExposures(patchRef, calExpRefList, self.getCoaddDatasetName(self.warpType),
649 self.getTempExpDatasetName(self.warpType))
650 tempExpRefList = [getGroupDataRef(butler, self.getTempExpDatasetName(self.warpType),
651 g, groupData.keys) for
652 g in groupData.groups.keys()]
653 return tempExpRefList
654
655 def prepareInputs(self, refList):
656 """Prepare the input warps for coaddition by measuring the weight for
657 each warp and the scaling for the photometric zero point.
658
659 Each Warp has its own photometric zeropoint and background variance.
660 Before coadding these Warps together, compute a scale factor to
661 normalize the photometric zeropoint and compute the weight for each Warp.
662
663 Parameters
664 ----------
665 refList : `list`
666 List of data references to tempExp
667
668 Returns
669 -------
670 result : `lsst.pipe.base.Struct`
671 Result struct with components:
672
673 - ``tempExprefList``: `list` of data references to tempExp.
674 - ``weightList``: `list` of weightings.
675 - ``imageScalerList``: `list` of image scalers.
676 """
677 statsCtrl = afwMath.StatisticsControl()
678 statsCtrl.setNumSigmaClip(self.config.sigmaClip)
679 statsCtrl.setNumIter(self.config.clipIter)
680 statsCtrl.setAndMask(self.getBadPixelMask())
681 statsCtrl.setNanSafe(True)
682 # compute tempExpRefList: a list of tempExpRef that actually exist
683 # and weightList: a list of the weight of the associated coadd tempExp
684 # and imageScalerList: a list of scale factors for the associated coadd tempExp
685 tempExpRefList = []
686 weightList = []
687 imageScalerList = []
688 tempExpName = self.getTempExpDatasetName(self.warpType)
689 for tempExpRef in refList:
690 # Gen3's DeferredDatasetHandles are guaranteed to exist and
691 # therefore have no datasetExists() method
692 if not isinstance(tempExpRef, DeferredDatasetHandle):
693 if not tempExpRef.datasetExists(tempExpName):
694 self.log.warning("Could not find %s %s; skipping it", tempExpName, tempExpRef.dataId)
695 continue
696
697 tempExp = tempExpRef.get(datasetType=tempExpName, immediate=True)
698 # Ignore any input warp that is empty of data
699 if numpy.isnan(tempExp.image.array).all():
700 continue
701 maskedImage = tempExp.getMaskedImage()
702 imageScaler = self.scaleZeroPoint.computeImageScaler(
703 exposure=tempExp,
704 dataRef=tempExpRef,
705 )
706 try:
707 imageScaler.scaleMaskedImage(maskedImage)
708 except Exception as e:
709 self.log.warning("Scaling failed for %s (skipping it): %s", tempExpRef.dataId, e)
710 continue
711 statObj = afwMath.makeStatistics(maskedImage.getVariance(), maskedImage.getMask(),
712 afwMath.MEANCLIP, statsCtrl)
713 meanVar, meanVarErr = statObj.getResult(afwMath.MEANCLIP)
714 weight = 1.0 / float(meanVar)
715 if not numpy.isfinite(weight):
716 self.log.warning("Non-finite weight for %s: skipping", tempExpRef.dataId)
717 continue
718 self.log.info("Weight of %s %s = %0.3f", tempExpName, tempExpRef.dataId, weight)
719
720 del maskedImage
721 del tempExp
722
723 tempExpRefList.append(tempExpRef)
724 weightList.append(weight)
725 imageScalerList.append(imageScaler)
726
727 return pipeBase.Struct(tempExpRefList=tempExpRefList, weightList=weightList,
728 imageScalerList=imageScalerList)
729
730 def prepareStats(self, mask=None):
731 """Prepare the statistics for coadding images.
732
733 Parameters
734 ----------
735 mask : `int`, optional
736 Bit mask value to exclude from coaddition.
737
738 Returns
739 -------
740 stats : `lsst.pipe.base.Struct`
741 Statistics structure with the following fields:
742
743 - ``statsCtrl``: Statistics control object for coadd
745 - ``statsFlags``: Statistic for coadd (`lsst.afw.math.Property`)
746 """
747 if mask is None:
748 mask = self.getBadPixelMask()
749 statsCtrl = afwMath.StatisticsControl()
750 statsCtrl.setNumSigmaClip(self.config.sigmaClip)
751 statsCtrl.setNumIter(self.config.clipIter)
752 statsCtrl.setAndMask(mask)
753 statsCtrl.setNanSafe(True)
754 statsCtrl.setWeighted(True)
755 statsCtrl.setCalcErrorFromInputVariance(self.config.calcErrorFromInputVariance)
756 for plane, threshold in self.config.maskPropagationThresholds.items():
757 bit = afwImage.Mask.getMaskPlane(plane)
758 statsCtrl.setMaskPropagationThreshold(bit, threshold)
759 statsFlags = afwMath.stringToStatisticsProperty(self.config.statistic)
760 return pipeBase.Struct(ctrl=statsCtrl, flags=statsFlags)
761
762 @timeMethod
763 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
764 altMaskList=None, mask=None, supplementaryData=None):
765 """Assemble a coadd from input warps
766
767 Assemble the coadd using the provided list of coaddTempExps. Since
768 the full coadd covers a patch (a large area), the assembly is
769 performed over small areas on the image at a time in order to
770 conserve memory usage. Iterate over subregions within the outer
771 bbox of the patch using `assembleSubregion` to stack the corresponding
772 subregions from the coaddTempExps with the statistic specified.
773 Set the edge bits the coadd mask based on the weight map.
774
775 Parameters
776 ----------
777 skyInfo : `lsst.pipe.base.Struct`
778 Struct with geometric information about the patch.
779 tempExpRefList : `list`
780 List of data references to Warps (previously called CoaddTempExps).
781 imageScalerList : `list`
782 List of image scalers.
783 weightList : `list`
784 List of weights
785 altMaskList : `list`, optional
786 List of alternate masks to use rather than those stored with
787 tempExp.
788 mask : `int`, optional
789 Bit mask value to exclude from coaddition.
790 supplementaryData : lsst.pipe.base.Struct, optional
791 Struct with additional data products needed to assemble coadd.
792 Only used by subclasses that implement `makeSupplementaryData`
793 and override `run`.
794
795 Returns
796 -------
797 result : `lsst.pipe.base.Struct`
798 Result struct with components:
799
800 - ``coaddExposure``: coadded exposure (``lsst.afw.image.Exposure``).
801 - ``nImage``: exposure count image (``lsst.afw.image.Image``), if requested.
802 - ``inputMap``: bit-wise map of inputs, if requested.
803 - ``warpRefList``: input list of refs to the warps (
804 ``lsst.daf.butler.DeferredDatasetHandle`` or
806 (unmodified)
807 - ``imageScalerList``: input list of image scalers (unmodified)
808 - ``weightList``: input list of weights (unmodified)
809 """
810 tempExpName = self.getTempExpDatasetName(self.warpType)
811 self.log.info("Assembling %s %s", len(tempExpRefList), tempExpName)
812 stats = self.prepareStats(mask=mask)
813
814 if altMaskList is None:
815 altMaskList = [None]*len(tempExpRefList)
816
817 coaddExposure = afwImage.ExposureF(skyInfo.bbox, skyInfo.wcs)
818 coaddExposure.setPhotoCalib(self.scaleZeroPoint.getPhotoCalib())
819 coaddExposure.getInfo().setCoaddInputs(self.inputRecorder.makeCoaddInputs())
820 self.assembleMetadata(coaddExposure, tempExpRefList, weightList)
821 coaddMaskedImage = coaddExposure.getMaskedImage()
822 subregionSizeArr = self.config.subregionSize
823 subregionSize = geom.Extent2I(subregionSizeArr[0], subregionSizeArr[1])
824 # if nImage is requested, create a zero one which can be passed to assembleSubregion
825 if self.config.doNImage:
826 nImage = afwImage.ImageU(skyInfo.bbox)
827 else:
828 nImage = None
829 # If inputMap is requested, create the initial version that can be masked in
830 # assembleSubregion.
831 if self.config.doInputMap:
832 self.inputMapper.build_ccd_input_map(skyInfo.bbox,
833 skyInfo.wcs,
834 coaddExposure.getInfo().getCoaddInputs().ccds)
835
836 if self.config.doOnlineForMean and self.config.statistic == "MEAN":
837 try:
838 self.assembleOnlineMeanCoadd(coaddExposure, tempExpRefList, imageScalerList,
839 weightList, altMaskList, stats.ctrl,
840 nImage=nImage)
841 except Exception as e:
842 self.log.exception("Cannot compute online coadd %s", e)
843 raise
844 else:
845 for subBBox in self._subBBoxIter(skyInfo.bbox, subregionSize):
846 try:
847 self.assembleSubregion(coaddExposure, subBBox, tempExpRefList, imageScalerList,
848 weightList, altMaskList, stats.flags, stats.ctrl,
849 nImage=nImage)
850 except Exception as e:
851 self.log.exception("Cannot compute coadd %s: %s", subBBox, e)
852 raise
853
854 # If inputMap is requested, we must finalize the map after the accumulation.
855 if self.config.doInputMap:
856 self.inputMapper.finalize_ccd_input_map_mask()
857 inputMap = self.inputMapper.ccd_input_map
858 else:
859 inputMap = None
860
861 self.setInexactPsf(coaddMaskedImage.getMask())
862 # Despite the name, the following doesn't really deal with "EDGE" pixels: it identifies
863 # pixels that didn't receive any unmasked inputs (as occurs around the edge of the field).
864 coaddUtils.setCoaddEdgeBits(coaddMaskedImage.getMask(), coaddMaskedImage.getVariance())
865 return pipeBase.Struct(coaddExposure=coaddExposure, nImage=nImage,
866 warpRefList=tempExpRefList, imageScalerList=imageScalerList,
867 weightList=weightList, inputMap=inputMap)
868
869 def assembleMetadata(self, coaddExposure, tempExpRefList, weightList):
870 """Set the metadata for the coadd.
871
872 This basic implementation sets the filter from the first input.
873
874 Parameters
875 ----------
876 coaddExposure : `lsst.afw.image.Exposure`
877 The target exposure for the coadd.
878 tempExpRefList : `list`
879 List of data references to tempExp.
880 weightList : `list`
881 List of weights.
882 """
883 assert len(tempExpRefList) == len(weightList), "Length mismatch"
884 tempExpName = self.getTempExpDatasetName(self.warpType)
885 # We load a single pixel of each coaddTempExp, because we just want to get at the metadata
886 # (and we need more than just the PropertySet that contains the header), which is not possible
887 # with the current butler (see #2777).
888 bbox = geom.Box2I(coaddExposure.getBBox().getMin(), geom.Extent2I(1, 1))
889
890 if isinstance(tempExpRefList[0], DeferredDatasetHandle):
891 # Gen 3 API
892 tempExpList = [tempExpRef.get(parameters={'bbox': bbox}) for tempExpRef in tempExpRefList]
893 else:
894 # Gen 2 API. Delete this when Gen 2 retired
895 tempExpList = [tempExpRef.get(tempExpName + "_sub", bbox=bbox, immediate=True)
896 for tempExpRef in tempExpRefList]
897 numCcds = sum(len(tempExp.getInfo().getCoaddInputs().ccds) for tempExp in tempExpList)
898
899 # Set the coadd FilterLabel to the band of the first input exposure:
900 # Coadds are calibrated, so the physical label is now meaningless.
901 coaddExposure.setFilter(afwImage.FilterLabel(tempExpList[0].getFilter().bandLabel))
902 coaddInputs = coaddExposure.getInfo().getCoaddInputs()
903 coaddInputs.ccds.reserve(numCcds)
904 coaddInputs.visits.reserve(len(tempExpList))
905
906 for tempExp, weight in zip(tempExpList, weightList):
907 self.inputRecorder.addVisitToCoadd(coaddInputs, tempExp, weight)
908
909 if self.config.doUsePsfMatchedPolygons:
910 self.shrinkValidPolygons(coaddInputs)
911
912 coaddInputs.visits.sort()
913 coaddInputs.ccds.sort()
914 if self.warpType == "psfMatched":
915 # The modelPsf BBox for a psfMatchedWarp/coaddTempExp was dynamically defined by
916 # ModelPsfMatchTask as the square box bounding its spatially-variable, pre-matched WarpedPsf.
917 # Likewise, set the PSF of a PSF-Matched Coadd to the modelPsf
918 # having the maximum width (sufficient because square)
919 modelPsfList = [tempExp.getPsf() for tempExp in tempExpList]
920 modelPsfWidthList = [modelPsf.computeBBox(modelPsf.getAveragePosition()).getWidth()
921 for modelPsf in modelPsfList]
922 psf = modelPsfList[modelPsfWidthList.index(max(modelPsfWidthList))]
923 else:
924 psf = measAlg.CoaddPsf(coaddInputs.ccds, coaddExposure.getWcs(),
925 self.config.coaddPsf.makeControl())
926 coaddExposure.setPsf(psf)
927 apCorrMap = measAlg.makeCoaddApCorrMap(coaddInputs.ccds, coaddExposure.getBBox(afwImage.PARENT),
928 coaddExposure.getWcs())
929 coaddExposure.getInfo().setApCorrMap(apCorrMap)
930 if self.config.doAttachTransmissionCurve:
931 transmissionCurve = measAlg.makeCoaddTransmissionCurve(coaddExposure.getWcs(), coaddInputs.ccds)
932 coaddExposure.getInfo().setTransmissionCurve(transmissionCurve)
933
934 def assembleSubregion(self, coaddExposure, bbox, tempExpRefList, imageScalerList, weightList,
935 altMaskList, statsFlags, statsCtrl, nImage=None):
936 """Assemble the coadd for a sub-region.
937
938 For each coaddTempExp, check for (and swap in) an alternative mask
939 if one is passed. Remove mask planes listed in
940 `config.removeMaskPlanes`. Finally, stack the actual exposures using
941 `lsst.afw.math.statisticsStack` with the statistic specified by
942 statsFlags. Typically, the statsFlag will be one of lsst.afw.math.MEAN for
943 a mean-stack or `lsst.afw.math.MEANCLIP` for outlier rejection using
944 an N-sigma clipped mean where N and iterations are specified by
945 statsCtrl. Assign the stacked subregion back to the coadd.
946
947 Parameters
948 ----------
949 coaddExposure : `lsst.afw.image.Exposure`
950 The target exposure for the coadd.
951 bbox : `lsst.geom.Box`
952 Sub-region to coadd.
953 tempExpRefList : `list`
954 List of data reference to tempExp.
955 imageScalerList : `list`
956 List of image scalers.
957 weightList : `list`
958 List of weights.
959 altMaskList : `list`
960 List of alternate masks to use rather than those stored with
961 tempExp, or None. Each element is dict with keys = mask plane
962 name to which to add the spans.
963 statsFlags : `lsst.afw.math.Property`
964 Property object for statistic for coadd.
966 Statistics control object for coadd.
967 nImage : `lsst.afw.image.ImageU`, optional
968 Keeps track of exposure count for each pixel.
969 """
970 self.log.debug("Computing coadd over %s", bbox)
971 tempExpName = self.getTempExpDatasetName(self.warpType)
972 coaddExposure.mask.addMaskPlane("REJECTED")
973 coaddExposure.mask.addMaskPlane("CLIPPED")
974 coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
975 maskMap = self.setRejectedMaskMapping(statsCtrl)
976 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
977 maskedImageList = []
978 if nImage is not None:
979 subNImage = afwImage.ImageU(bbox.getWidth(), bbox.getHeight())
980 for tempExpRef, imageScaler, altMask in zip(tempExpRefList, imageScalerList, altMaskList):
981
982 if isinstance(tempExpRef, DeferredDatasetHandle):
983 # Gen 3 API
984 exposure = tempExpRef.get(parameters={'bbox': bbox})
985 else:
986 # Gen 2 API. Delete this when Gen 2 retired
987 exposure = tempExpRef.get(tempExpName + "_sub", bbox=bbox)
988
989 maskedImage = exposure.getMaskedImage()
990 mask = maskedImage.getMask()
991 if altMask is not None:
992 self.applyAltMaskPlanes(mask, altMask)
993 imageScaler.scaleMaskedImage(maskedImage)
994
995 # Add 1 for each pixel which is not excluded by the exclude mask.
996 # In legacyCoadd, pixels may also be excluded by afwMath.statisticsStack.
997 if nImage is not None:
998 subNImage.getArray()[maskedImage.getMask().getArray() & statsCtrl.getAndMask() == 0] += 1
999 if self.config.removeMaskPlanes:
1000 self.removeMaskPlanes(maskedImage)
1001 maskedImageList.append(maskedImage)
1002
1003 if self.config.doInputMap:
1004 visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
1005 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
1006
1007 with self.timer("stack"):
1008 coaddSubregion = afwMath.statisticsStack(maskedImageList, statsFlags, statsCtrl, weightList,
1009 clipped, # also set output to CLIPPED if sigma-clipped
1010 maskMap)
1011 coaddExposure.maskedImage.assign(coaddSubregion, bbox)
1012 if nImage is not None:
1013 nImage.assign(subNImage, bbox)
1014
1015 def assembleOnlineMeanCoadd(self, coaddExposure, tempExpRefList, imageScalerList, weightList,
1016 altMaskList, statsCtrl, nImage=None):
1017 """Assemble the coadd using the "online" method.
1018
1019 This method takes a running sum of images and weights to save memory.
1020 It only works for MEAN statistics.
1021
1022 Parameters
1023 ----------
1024 coaddExposure : `lsst.afw.image.Exposure`
1025 The target exposure for the coadd.
1026 tempExpRefList : `list`
1027 List of data reference to tempExp.
1028 imageScalerList : `list`
1029 List of image scalers.
1030 weightList : `list`
1031 List of weights.
1032 altMaskList : `list`
1033 List of alternate masks to use rather than those stored with
1034 tempExp, or None. Each element is dict with keys = mask plane
1035 name to which to add the spans.
1037 Statistics control object for coadd
1038 nImage : `lsst.afw.image.ImageU`, optional
1039 Keeps track of exposure count for each pixel.
1040 """
1041 self.log.debug("Computing online coadd.")
1042 tempExpName = self.getTempExpDatasetName(self.warpType)
1043 coaddExposure.mask.addMaskPlane("REJECTED")
1044 coaddExposure.mask.addMaskPlane("CLIPPED")
1045 coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
1046 maskMap = self.setRejectedMaskMapping(statsCtrl)
1047 thresholdDict = AccumulatorMeanStack.stats_ctrl_to_threshold_dict(statsCtrl)
1048
1049 bbox = coaddExposure.maskedImage.getBBox()
1050
1051 stacker = AccumulatorMeanStack(
1052 coaddExposure.image.array.shape,
1053 statsCtrl.getAndMask(),
1054 mask_threshold_dict=thresholdDict,
1055 mask_map=maskMap,
1056 no_good_pixels_mask=statsCtrl.getNoGoodPixelsMask(),
1057 calc_error_from_input_variance=self.config.calcErrorFromInputVariance,
1058 compute_n_image=(nImage is not None)
1059 )
1060
1061 for tempExpRef, imageScaler, altMask, weight in zip(tempExpRefList,
1062 imageScalerList,
1063 altMaskList,
1064 weightList):
1065 if isinstance(tempExpRef, DeferredDatasetHandle):
1066 # Gen 3 API
1067 exposure = tempExpRef.get()
1068 else:
1069 # Gen 2 API. Delete this when Gen 2 retired
1070 exposure = tempExpRef.get(tempExpName)
1071
1072 maskedImage = exposure.getMaskedImage()
1073 mask = maskedImage.getMask()
1074 if altMask is not None:
1075 self.applyAltMaskPlanes(mask, altMask)
1076 imageScaler.scaleMaskedImage(maskedImage)
1077 if self.config.removeMaskPlanes:
1078 self.removeMaskPlanes(maskedImage)
1079
1080 stacker.add_masked_image(maskedImage, weight=weight)
1081
1082 if self.config.doInputMap:
1083 visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
1084 self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
1085
1086 stacker.fill_stacked_masked_image(coaddExposure.maskedImage)
1087
1088 if nImage is not None:
1089 nImage.array[:, :] = stacker.n_image
1090
1091 def removeMaskPlanes(self, maskedImage):
1092 """Unset the mask of an image for mask planes specified in the config.
1093
1094 Parameters
1095 ----------
1096 maskedImage : `lsst.afw.image.MaskedImage`
1097 The masked image to be modified.
1098 """
1099 mask = maskedImage.getMask()
1100 for maskPlane in self.config.removeMaskPlanes:
1101 try:
1102 mask &= ~mask.getPlaneBitMask(maskPlane)
1104 self.log.debug("Unable to remove mask plane %s: no mask plane with that name was found.",
1105 maskPlane)
1106
1107 @staticmethod
1108 def setRejectedMaskMapping(statsCtrl):
1109 """Map certain mask planes of the warps to new planes for the coadd.
1110
1111 If a pixel is rejected due to a mask value other than EDGE, NO_DATA,
1112 or CLIPPED, set it to REJECTED on the coadd.
1113 If a pixel is rejected due to EDGE, set the coadd pixel to SENSOR_EDGE.
1114 If a pixel is rejected due to CLIPPED, set the coadd pixel to CLIPPED.
1115
1116 Parameters
1117 ----------
1119 Statistics control object for coadd
1120
1121 Returns
1122 -------
1123 maskMap : `list` of `tuple` of `int`
1124 A list of mappings of mask planes of the warped exposures to
1125 mask planes of the coadd.
1126 """
1127 edge = afwImage.Mask.getPlaneBitMask("EDGE")
1128 noData = afwImage.Mask.getPlaneBitMask("NO_DATA")
1129 clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
1130 toReject = statsCtrl.getAndMask() & (~noData) & (~edge) & (~clipped)
1131 maskMap = [(toReject, afwImage.Mask.getPlaneBitMask("REJECTED")),
1132 (edge, afwImage.Mask.getPlaneBitMask("SENSOR_EDGE")),
1133 (clipped, clipped)]
1134 return maskMap
1135
1136 def applyAltMaskPlanes(self, mask, altMaskSpans):
1137 """Apply in place alt mask formatted as SpanSets to a mask.
1138
1139 Parameters
1140 ----------
1141 mask : `lsst.afw.image.Mask`
1142 Original mask.
1143 altMaskSpans : `dict`
1144 SpanSet lists to apply. Each element contains the new mask
1145 plane name (e.g. "CLIPPED and/or "NO_DATA") as the key,
1146 and list of SpanSets to apply to the mask.
1147
1148 Returns
1149 -------
1150 mask : `lsst.afw.image.Mask`
1151 Updated mask.
1152 """
1153 if self.config.doUsePsfMatchedPolygons:
1154 if ("NO_DATA" in altMaskSpans) and ("NO_DATA" in self.config.badMaskPlanes):
1155 # Clear away any other masks outside the validPolygons. These pixels are no longer
1156 # contributing to inexact PSFs, and will still be rejected because of NO_DATA
1157 # self.config.doUsePsfMatchedPolygons should be True only in CompareWarpAssemble
1158 # This mask-clearing step must only occur *before* applying the new masks below
1159 for spanSet in altMaskSpans['NO_DATA']:
1160 spanSet.clippedTo(mask.getBBox()).clearMask(mask, self.getBadPixelMask())
1161
1162 for plane, spanSetList in altMaskSpans.items():
1163 maskClipValue = mask.addMaskPlane(plane)
1164 for spanSet in spanSetList:
1165 spanSet.clippedTo(mask.getBBox()).setMask(mask, 2**maskClipValue)
1166 return mask
1167
1168 def shrinkValidPolygons(self, coaddInputs):
1169 """Shrink coaddInputs' ccds' ValidPolygons in place.
1170
1171 Either modify each ccd's validPolygon in place, or if CoaddInputs
1172 does not have a validPolygon, create one from its bbox.
1173
1174 Parameters
1175 ----------
1176 coaddInputs : `lsst.afw.image.coaddInputs`
1177 Original mask.
1178
1179 """
1180 for ccd in coaddInputs.ccds:
1181 polyOrig = ccd.getValidPolygon()
1182 validPolyBBox = polyOrig.getBBox() if polyOrig else ccd.getBBox()
1183 validPolyBBox.grow(-self.config.matchingKernelSize//2)
1184 if polyOrig:
1185 validPolygon = polyOrig.intersectionSingle(validPolyBBox)
1186 else:
1187 validPolygon = afwGeom.polygon.Polygon(geom.Box2D(validPolyBBox))
1188 ccd.setValidPolygon(validPolygon)
1189
1190 def readBrightObjectMasks(self, dataRef):
1191 """Retrieve the bright object masks.
1192
1193 Returns None on failure.
1194
1195 Parameters
1196 ----------
1198 A Butler dataRef.
1199
1200 Returns
1201 -------
1203 Bright object mask from the Butler object, or None if it cannot
1204 be retrieved.
1205 """
1206 try:
1207 return dataRef.get(datasetType="brightObjectMask", immediate=True)
1208 except Exception as e:
1209 self.log.warning("Unable to read brightObjectMask for %s: %s", dataRef.dataId, e)
1210 return None
1211
1212 def setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None):
1213 """Set the bright object masks.
1214
1215 Parameters
1216 ----------
1217 exposure : `lsst.afw.image.Exposure`
1218 Exposure under consideration.
1220 Data identifier dict for patch.
1221 brightObjectMasks : `lsst.afw.table`
1222 Table of bright objects to mask.
1223 """
1224
1225 if brightObjectMasks is None:
1226 self.log.warning("Unable to apply bright object mask: none supplied")
1227 return
1228 self.log.info("Applying %d bright object masks to %s", len(brightObjectMasks), dataId)
1229 mask = exposure.getMaskedImage().getMask()
1230 wcs = exposure.getWcs()
1231 plateScale = wcs.getPixelScale().asArcseconds()
1232
1233 for rec in brightObjectMasks:
1234 center = geom.PointI(wcs.skyToPixel(rec.getCoord()))
1235 if rec["type"] == "box":
1236 assert rec["angle"] == 0.0, ("Angle != 0 for mask object %s" % rec["id"])
1237 width = rec["width"].asArcseconds()/plateScale # convert to pixels
1238 height = rec["height"].asArcseconds()/plateScale # convert to pixels
1239
1240 halfSize = geom.ExtentI(0.5*width, 0.5*height)
1241 bbox = geom.Box2I(center - halfSize, center + halfSize)
1242
1243 bbox = geom.BoxI(geom.PointI(int(center[0] - 0.5*width), int(center[1] - 0.5*height)),
1244 geom.PointI(int(center[0] + 0.5*width), int(center[1] + 0.5*height)))
1245 spans = afwGeom.SpanSet(bbox)
1246 elif rec["type"] == "circle":
1247 radius = int(rec["radius"].asArcseconds()/plateScale) # convert to pixels
1248 spans = afwGeom.SpanSet.fromShape(radius, offset=center)
1249 else:
1250 self.log.warning("Unexpected region type %s at %s", rec["type"], center)
1251 continue
1252 spans.clippedTo(mask.getBBox()).setMask(mask, self.brightObjectBitmask)
1253
1254 def setInexactPsf(self, mask):
1255 """Set INEXACT_PSF mask plane.
1256
1257 If any of the input images isn't represented in the coadd (due to
1258 clipped pixels or chip gaps), the `CoaddPsf` will be inexact. Flag
1259 these pixels.
1260
1261 Parameters
1262 ----------
1263 mask : `lsst.afw.image.Mask`
1264 Coadded exposure's mask, modified in-place.
1265 """
1266 mask.addMaskPlane("INEXACT_PSF")
1267 inexactPsf = mask.getPlaneBitMask("INEXACT_PSF")
1268 sensorEdge = mask.getPlaneBitMask("SENSOR_EDGE") # chip edges (so PSF is discontinuous)
1269 clipped = mask.getPlaneBitMask("CLIPPED") # pixels clipped from coadd
1270 rejected = mask.getPlaneBitMask("REJECTED") # pixels rejected from coadd due to masks
1271 array = mask.getArray()
1272 selected = array & (sensorEdge | clipped | rejected) > 0
1273 array[selected] |= inexactPsf
1274
1275 @classmethod
1276 def _makeArgumentParser(cls):
1277 """Create an argument parser.
1278 """
1279 parser = pipeBase.ArgumentParser(name=cls._DefaultName)
1280 parser.add_id_argument("--id", cls.ConfigClass().coaddName + "Coadd_"
1281 + cls.ConfigClass().warpType + "Warp",
1282 help="data ID, e.g. --id tract=12345 patch=1,2",
1283 ContainerClass=AssembleCoaddDataIdContainer)
1284 parser.add_id_argument("--selectId", "calexp", help="data ID, e.g. --selectId visit=6789 ccd=0..9",
1285 ContainerClass=SelectDataIdContainer)
1286 return parser
1287
1288 @staticmethod
1289 def _subBBoxIter(bbox, subregionSize):
1290 """Iterate over subregions of a bbox.
1291
1292 Parameters
1293 ----------
1294 bbox : `lsst.geom.Box2I`
1295 Bounding box over which to iterate.
1296 subregionSize: `lsst.geom.Extent2I`
1297 Size of sub-bboxes.
1298
1299 Yields
1300 ------
1301 subBBox : `lsst.geom.Box2I`
1302 Next sub-bounding box of size ``subregionSize`` or smaller; each ``subBBox``
1303 is contained within ``bbox``, so it may be smaller than ``subregionSize`` at
1304 the edges of ``bbox``, but it will never be empty.
1305 """
1306 if bbox.isEmpty():
1307 raise RuntimeError("bbox %s is empty" % (bbox,))
1308 if subregionSize[0] < 1 or subregionSize[1] < 1:
1309 raise RuntimeError("subregionSize %s must be nonzero" % (subregionSize,))
1310
1311 for rowShift in range(0, bbox.getHeight(), subregionSize[1]):
1312 for colShift in range(0, bbox.getWidth(), subregionSize[0]):
1313 subBBox = geom.Box2I(bbox.getMin() + geom.Extent2I(colShift, rowShift), subregionSize)
1314 subBBox.clip(bbox)
1315 if subBBox.isEmpty():
1316 raise RuntimeError("Bug: empty bbox! bbox=%s, subregionSize=%s, "
1317 "colShift=%s, rowShift=%s" %
1318 (bbox, subregionSize, colShift, rowShift))
1319 yield subBBox
1320
1321 def filterWarps(self, inputs, goodVisits):
1322 """Return list of only inputRefs with visitId in goodVisits ordered by goodVisit
1323
1324 Parameters
1325 ----------
1326 inputs : list
1327 List of `lsst.pipe.base.connections.DeferredDatasetRef` with dataId containing visit
1328 goodVisit : `dict`
1329 Dictionary with good visitIds as the keys. Value ignored.
1330
1331 Returns:
1332 --------
1333 filteredInputs : `list`
1334 Filtered and sorted list of `lsst.pipe.base.connections.DeferredDatasetRef`
1335 """
1336 inputWarpDict = {inputRef.ref.dataId['visit']: inputRef for inputRef in inputs}
1337 filteredInputs = []
1338 for visit in goodVisits.keys():
1339 if visit in inputWarpDict:
1340 filteredInputs.append(inputWarpDict[visit])
1341 return filteredInputs
1342
1343
1344class AssembleCoaddDataIdContainer(pipeBase.DataIdContainer):
1345 """A version of `lsst.pipe.base.DataIdContainer` specialized for assembleCoadd.
1346 """
1347
1348 def makeDataRefList(self, namespace):
1349 """Make self.refList from self.idList.
1350
1351 Parameters
1352 ----------
1353 namespace
1354 Results of parsing command-line (with ``butler`` and ``log`` elements).
1355 """
1356 datasetType = namespace.config.coaddName + "Coadd"
1357 keysCoadd = namespace.butler.getKeys(datasetType=datasetType, level=self.level)
1358
1359 for dataId in self.idList:
1360 # tract and patch are required
1361 for key in keysCoadd:
1362 if key not in dataId:
1363 raise RuntimeError("--id must include " + key)
1364
1365 dataRef = namespace.butler.dataRef(
1366 datasetType=datasetType,
1367 dataId=dataId,
1368 )
1369 self.refList.append(dataRef)
1370
1371
1372def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask):
1373 """Function to count the number of pixels with a specific mask in a
1374 footprint.
1375
1376 Find the intersection of mask & footprint. Count all pixels in the mask
1377 that are in the intersection that have bitmask set but do not have
1378 ignoreMask set. Return the count.
1379
1380 Parameters
1381 ----------
1382 mask : `lsst.afw.image.Mask`
1383 Mask to define intersection region by.
1384 footprint : `lsst.afw.detection.Footprint`
1385 Footprint to define the intersection region by.
1386 bitmask
1387 Specific mask that we wish to count the number of occurances of.
1388 ignoreMask
1389 Pixels to not consider.
1390
1391 Returns
1392 -------
1393 result : `int`
1394 Count of number of pixels in footprint with specified mask.
1395 """
1396 bbox = footprint.getBBox()
1397 bbox.clip(mask.getBBox(afwImage.PARENT))
1398 fp = afwImage.Mask(bbox)
1399 subMask = mask.Factory(mask, bbox, afwImage.PARENT)
1400 footprint.spans.setMask(fp, bitmask)
1401 return numpy.logical_and((subMask.getArray() & fp.getArray()) > 0,
1402 (subMask.getArray() & ignoreMask) == 0).sum()
1403
1404
1405class SafeClipAssembleCoaddConfig(AssembleCoaddConfig, pipelineConnections=AssembleCoaddConnections):
1406 """Configuration parameters for the SafeClipAssembleCoaddTask.
1407 """
1408 clipDetection = pexConfig.ConfigurableField(
1409 target=SourceDetectionTask,
1410 doc="Detect sources on difference between unclipped and clipped coadd")
1411 minClipFootOverlap = pexConfig.Field(
1412 doc="Minimum fractional overlap of clipped footprint with visit DETECTED to be clipped",
1413 dtype=float,
1414 default=0.6
1415 )
1416 minClipFootOverlapSingle = pexConfig.Field(
1417 doc="Minimum fractional overlap of clipped footprint with visit DETECTED to be "
1418 "clipped when only one visit overlaps",
1419 dtype=float,
1420 default=0.5
1421 )
1422 minClipFootOverlapDouble = pexConfig.Field(
1423 doc="Minimum fractional overlap of clipped footprints with visit DETECTED to be "
1424 "clipped when two visits overlap",
1425 dtype=float,
1426 default=0.45
1427 )
1428 maxClipFootOverlapDouble = pexConfig.Field(
1429 doc="Maximum fractional overlap of clipped footprints with visit DETECTED when "
1430 "considering two visits",
1431 dtype=float,
1432 default=0.15
1433 )
1434 minBigOverlap = pexConfig.Field(
1435 doc="Minimum number of pixels in footprint to use DETECTED mask from the single visits "
1436 "when labeling clipped footprints",
1437 dtype=int,
1438 default=100
1439 )
1440
1441 def setDefaults(self):
1442 """Set default values for clipDetection.
1443
1444 Notes
1445 -----
1446 The numeric values for these configuration parameters were
1447 empirically determined, future work may further refine them.
1448 """
1449 AssembleCoaddConfig.setDefaults(self)
1450 self.clipDetectionclipDetection.doTempLocalBackground = False
1451 self.clipDetectionclipDetection.reEstimateBackground = False
1452 self.clipDetectionclipDetection.returnOriginalFootprints = False
1453 self.clipDetectionclipDetection.thresholdPolarity = "both"
1454 self.clipDetectionclipDetection.thresholdValue = 2
1455 self.clipDetectionclipDetection.nSigmaToGrow = 2
1456 self.clipDetectionclipDetection.minPixels = 4
1457 self.clipDetectionclipDetection.isotropicGrow = True
1458 self.clipDetectionclipDetection.thresholdType = "pixel_stdev"
1459 self.sigmaClipsigmaClip = 1.5
1460 self.clipIterclipIter = 3
1461 self.statisticstatistic = "MEAN"
1462
1463 def validate(self):
1464 if self.doSigmaClipdoSigmaClip:
1465 log.warning("Additional Sigma-clipping not allowed in Safe-clipped Coadds. "
1466 "Ignoring doSigmaClip.")
1467 self.doSigmaClipdoSigmaClip = False
1468 if self.statisticstatistic != "MEAN":
1469 raise ValueError("Only MEAN statistic allowed for final stacking in SafeClipAssembleCoadd "
1470 "(%s chosen). Please set statistic to MEAN."
1471 % (self.statisticstatistic))
1472 AssembleCoaddTask.ConfigClass.validate(self)
1473
1474
1475class SafeClipAssembleCoaddTask(AssembleCoaddTask):
1476 """Assemble a coadded image from a set of coadded temporary exposures,
1477 being careful to clip & flag areas with potential artifacts.
1478
1479 In ``AssembleCoaddTask``, we compute the coadd as an clipped mean (i.e.,
1480 we clip outliers). The problem with doing this is that when computing the
1481 coadd PSF at a given location, individual visit PSFs from visits with
1482 outlier pixels contribute to the coadd PSF and cannot be treated correctly.
1483 In this task, we correct for this behavior by creating a new
1484 ``badMaskPlane`` 'CLIPPED'. We populate this plane on the input
1485 coaddTempExps and the final coadd where
1486
1487 i. difference imaging suggests that there is an outlier and
1488 ii. this outlier appears on only one or two images.
1489
1490 Such regions will not contribute to the final coadd. Furthermore, any
1491 routine to determine the coadd PSF can now be cognizant of clipped regions.
1492 Note that the algorithm implemented by this task is preliminary and works
1493 correctly for HSC data. Parameter modifications and or considerable
1494 redesigning of the algorithm is likley required for other surveys.
1495
1496 ``SafeClipAssembleCoaddTask`` uses a ``SourceDetectionTask``
1497 "clipDetection" subtask and also sub-classes ``AssembleCoaddTask``.
1498 You can retarget the ``SourceDetectionTask`` "clipDetection" subtask
1499 if you wish.
1500
1501 Notes
1502 -----
1503 The `lsst.pipe.base.cmdLineTask.CmdLineTask` interface supports a
1504 flag ``-d`` to import ``debug.py`` from your ``PYTHONPATH``;
1505 see `baseDebug` for more about ``debug.py`` files.
1506 `SafeClipAssembleCoaddTask` has no debug variables of its own.
1507 The ``SourceDetectionTask`` "clipDetection" subtasks may support debug
1508 variables. See the documetation for `SourceDetectionTask` "clipDetection"
1509 for further information.
1510
1511 Examples
1512 --------
1513 `SafeClipAssembleCoaddTask` assembles a set of warped ``coaddTempExp``
1514 images into a coadded image. The `SafeClipAssembleCoaddTask` is invoked by
1515 running assembleCoadd.py *without* the flag '--legacyCoadd'.
1516
1517 Usage of ``assembleCoadd.py`` expects a data reference to the tract patch
1518 and filter to be coadded (specified using
1519 '--id = [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]')
1520 along with a list of coaddTempExps to attempt to coadd (specified using
1521 '--selectId [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]').
1522 Only the coaddTempExps that cover the specified tract and patch will be
1523 coadded. A list of the available optional arguments can be obtained by
1524 calling assembleCoadd.py with the --help command line argument:
1525
1526 .. code-block:: none
1527
1528 assembleCoadd.py --help
1529
1530 To demonstrate usage of the `SafeClipAssembleCoaddTask` in the larger
1531 context of multi-band processing, we will generate the HSC-I & -R band
1532 coadds from HSC engineering test data provided in the ci_hsc package.
1533 To begin, assuming that the lsst stack has been already set up, we must
1534 set up the obs_subaru and ci_hsc packages. This defines the environment
1535 variable $CI_HSC_DIR and points at the location of the package. The raw
1536 HSC data live in the ``$CI_HSC_DIR/raw`` directory. To begin assembling
1537 the coadds, we must first
1538
1539 - ``processCcd``
1540 process the individual ccds in $CI_HSC_RAW to produce calibrated exposures
1541 - ``makeSkyMap``
1542 create a skymap that covers the area of the sky present in the raw exposures
1543 - ``makeCoaddTempExp``
1544 warp the individual calibrated exposures to the tangent plane of the coadd</DD>
1545
1546 We can perform all of these steps by running
1547
1548 .. code-block:: none
1549
1550 $CI_HSC_DIR scons warp-903986 warp-904014 warp-903990 warp-904010 warp-903988
1551
1552 This will produce warped coaddTempExps for each visit. To coadd the
1553 warped data, we call ``assembleCoadd.py`` as follows:
1554
1555 .. code-block:: none
1556
1557 assembleCoadd.py $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-I \
1558 --selectId visit=903986 ccd=16 --selectId visit=903986 ccd=22 --selectId visit=903986 ccd=23 \
1559 --selectId visit=903986 ccd=100--selectId visit=904014 ccd=1 --selectId visit=904014 ccd=6 \
1560 --selectId visit=904014 ccd=12 --selectId visit=903990 ccd=18 --selectId visit=903990 ccd=25 \
1561 --selectId visit=904010 ccd=4 --selectId visit=904010 ccd=10 --selectId visit=904010 ccd=100 \
1562 --selectId visit=903988 ccd=16 --selectId visit=903988 ccd=17 --selectId visit=903988 ccd=23 \
1563 --selectId visit=903988 ccd=24
1564
1565 This will process the HSC-I band data. The results are written in
1566 ``$CI_HSC_DIR/DATA/deepCoadd-results/HSC-I``.
1567
1568 You may also choose to run:
1569
1570 .. code-block:: none
1571
1572 scons warp-903334 warp-903336 warp-903338 warp-903342 warp-903344 warp-903346 nnn
1573 assembleCoadd.py $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-R --selectId visit=903334 ccd=16 \
1574 --selectId visit=903334 ccd=22 --selectId visit=903334 ccd=23 --selectId visit=903334 ccd=100 \
1575 --selectId visit=903336 ccd=17 --selectId visit=903336 ccd=24 --selectId visit=903338 ccd=18 \
1576 --selectId visit=903338 ccd=25 --selectId visit=903342 ccd=4 --selectId visit=903342 ccd=10 \
1577 --selectId visit=903342 ccd=100 --selectId visit=903344 ccd=0 --selectId visit=903344 ccd=5 \
1578 --selectId visit=903344 ccd=11 --selectId visit=903346 ccd=1 --selectId visit=903346 ccd=6 \
1579 --selectId visit=903346 ccd=12
1580
1581 to generate the coadd for the HSC-R band if you are interested in following
1582 multiBand Coadd processing as discussed in ``pipeTasks_multiBand``.
1583 """
1584 ConfigClass = SafeClipAssembleCoaddConfig
1585 _DefaultName = "safeClipAssembleCoadd"
1586
1587 def __init__(self, *args, **kwargs):
1588 AssembleCoaddTask.__init__(self, *args, **kwargs)
1589 schema = afwTable.SourceTable.makeMinimalSchema()
1590 self.makeSubtask("clipDetection", schema=schema)
1591
1592 @utils.inheritDoc(AssembleCoaddTask)
1593 @timeMethod
1594 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, *args, **kwargs):
1595 """Assemble the coadd for a region.
1596
1597 Compute the difference of coadds created with and without outlier
1598 rejection to identify coadd pixels that have outlier values in some
1599 individual visits.
1600 Detect clipped regions on the difference image and mark these regions
1601 on the one or two individual coaddTempExps where they occur if there
1602 is significant overlap between the clipped region and a source. This
1603 leaves us with a set of footprints from the difference image that have
1604 been identified as having occured on just one or two individual visits.
1605 However, these footprints were generated from a difference image. It
1606 is conceivable for a large diffuse source to have become broken up
1607 into multiple footprints acrosss the coadd difference in this process.
1608 Determine the clipped region from all overlapping footprints from the
1609 detected sources in each visit - these are big footprints.
1610 Combine the small and big clipped footprints and mark them on a new
1611 bad mask plane.
1612 Generate the coadd using `AssembleCoaddTask.run` without outlier
1613 removal. Clipped footprints will no longer make it into the coadd
1614 because they are marked in the new bad mask plane.
1615
1616 Notes
1617 -----
1618 args and kwargs are passed but ignored in order to match the call
1619 signature expected by the parent task.
1620 """
1621 exp = self.buildDifferenceImagebuildDifferenceImage(skyInfo, tempExpRefList, imageScalerList, weightList)
1622 mask = exp.getMaskedImage().getMask()
1623 mask.addMaskPlane("CLIPPED")
1624
1625 result = self.detectClipdetectClip(exp, tempExpRefList)
1626
1627 self.log.info('Found %d clipped objects', len(result.clipFootprints))
1628
1629 maskClipValue = mask.getPlaneBitMask("CLIPPED")
1630 maskDetValue = mask.getPlaneBitMask("DETECTED") | mask.getPlaneBitMask("DETECTED_NEGATIVE")
1631 # Append big footprints from individual Warps to result.clipSpans
1632 bigFootprints = self.detectClipBigdetectClipBig(result.clipSpans, result.clipFootprints, result.clipIndices,
1633 result.detectionFootprints, maskClipValue, maskDetValue,
1634 exp.getBBox())
1635 # Create mask of the current clipped footprints
1636 maskClip = mask.Factory(mask.getBBox(afwImage.PARENT))
1637 afwDet.setMaskFromFootprintList(maskClip, result.clipFootprints, maskClipValue)
1638
1639 maskClipBig = maskClip.Factory(mask.getBBox(afwImage.PARENT))
1640 afwDet.setMaskFromFootprintList(maskClipBig, bigFootprints, maskClipValue)
1641 maskClip |= maskClipBig
1642
1643 # Assemble coadd from base class, but ignoring CLIPPED pixels
1644 badMaskPlanes = self.config.badMaskPlanes[:]
1645 badMaskPlanes.append("CLIPPED")
1646 badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes)
1647 return AssembleCoaddTask.run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
1648 result.clipSpans, mask=badPixelMask)
1649
1650 def buildDifferenceImage(self, skyInfo, tempExpRefList, imageScalerList, weightList):
1651 """Return an exposure that contains the difference between unclipped
1652 and clipped coadds.
1653
1654 Generate a difference image between clipped and unclipped coadds.
1655 Compute the difference image by subtracting an outlier-clipped coadd
1656 from an outlier-unclipped coadd. Return the difference image.
1657
1658 Parameters
1659 ----------
1660 skyInfo : `lsst.pipe.base.Struct`
1661 Patch geometry information, from getSkyInfo
1662 tempExpRefList : `list`
1663 List of data reference to tempExp
1664 imageScalerList : `list`
1665 List of image scalers
1666 weightList : `list`
1667 List of weights
1668
1669 Returns
1670 -------
1672 Difference image of unclipped and clipped coadd wrapped in an Exposure
1673 """
1674 config = AssembleCoaddConfig()
1675 # getattr necessary because subtasks do not survive Config.toDict()
1676 # exclude connections because the class of self.config.connections is not
1677 # the same as AssembleCoaddConfig.connections, and the connections are not
1678 # needed to run this task anyway.
1679 configIntersection = {k: getattr(self.config, k)
1680 for k, v in self.config.toDict().items()
1681 if (k in config.keys() and k != "connections")}
1682 configIntersection['doInputMap'] = False
1683 configIntersection['doNImage'] = False
1684 config.update(**configIntersection)
1685
1686 # statistic MEAN copied from self.config.statistic, but for clarity explicitly assign
1687 config.statistic = 'MEAN'
1688 task = AssembleCoaddTask(config=config)
1689 coaddMean = task.run(skyInfo, tempExpRefList, imageScalerList, weightList).coaddExposure
1690
1691 config.statistic = 'MEANCLIP'
1692 task = AssembleCoaddTask(config=config)
1693 coaddClip = task.run(skyInfo, tempExpRefList, imageScalerList, weightList).coaddExposure
1694
1695 coaddDiff = coaddMean.getMaskedImage().Factory(coaddMean.getMaskedImage())
1696 coaddDiff -= coaddClip.getMaskedImage()
1697 exp = afwImage.ExposureF(coaddDiff)
1698 exp.setPsf(coaddMean.getPsf())
1699 return exp
1700
1701 def detectClip(self, exp, tempExpRefList):
1702 """Detect clipped regions on an exposure and set the mask on the
1703 individual tempExp masks.
1704
1705 Detect footprints in the difference image after smoothing the
1706 difference image with a Gaussian kernal. Identify footprints that
1707 overlap with one or two input ``coaddTempExps`` by comparing the
1708 computed overlap fraction to thresholds set in the config. A different
1709 threshold is applied depending on the number of overlapping visits
1710 (restricted to one or two). If the overlap exceeds the thresholds,
1711 the footprint is considered "CLIPPED" and is marked as such on the
1712 coaddTempExp. Return a struct with the clipped footprints, the indices
1713 of the ``coaddTempExps`` that end up overlapping with the clipped
1714 footprints, and a list of new masks for the ``coaddTempExps``.
1715
1716 Parameters
1717 ----------
1719 Exposure to run detection on.
1720 tempExpRefList : `list`
1721 List of data reference to tempExp.
1722
1723 Returns
1724 -------
1725 result : `lsst.pipe.base.Struct`
1726 Result struct with components:
1727
1728 - ``clipFootprints``: list of clipped footprints.
1729 - ``clipIndices``: indices for each ``clippedFootprint`` in
1730 ``tempExpRefList``.
1731 - ``clipSpans``: List of dictionaries containing spanSet lists
1732 to clip. Each element contains the new maskplane name
1733 ("CLIPPED") as the key and list of ``SpanSets`` as the value.
1734 - ``detectionFootprints``: List of DETECTED/DETECTED_NEGATIVE plane
1735 compressed into footprints.
1736 """
1737 mask = exp.getMaskedImage().getMask()
1738 maskDetValue = mask.getPlaneBitMask("DETECTED") | mask.getPlaneBitMask("DETECTED_NEGATIVE")
1739 fpSet = self.clipDetection.detectFootprints(exp, doSmooth=True, clearMask=True)
1740 # Merge positive and negative together footprints together
1741 fpSet.positive.merge(fpSet.negative)
1742 footprints = fpSet.positive
1743 self.log.info('Found %d potential clipped objects', len(footprints.getFootprints()))
1744 ignoreMask = self.getBadPixelMask()
1745
1746 clipFootprints = []
1747 clipIndices = []
1748 artifactSpanSets = [{'CLIPPED': list()} for _ in tempExpRefList]
1749
1750 # for use by detectClipBig
1751 visitDetectionFootprints = []
1752
1753 dims = [len(tempExpRefList), len(footprints.getFootprints())]
1754 overlapDetArr = numpy.zeros(dims, dtype=numpy.uint16)
1755 ignoreArr = numpy.zeros(dims, dtype=numpy.uint16)
1756
1757 # Loop over masks once and extract/store only relevant overlap metrics and detection footprints
1758 for i, warpRef in enumerate(tempExpRefList):
1759 tmpExpMask = warpRef.get(datasetType=self.getTempExpDatasetName(self.warpType),
1760 immediate=True).getMaskedImage().getMask()
1761 maskVisitDet = tmpExpMask.Factory(tmpExpMask, tmpExpMask.getBBox(afwImage.PARENT),
1762 afwImage.PARENT, True)
1763 maskVisitDet &= maskDetValue
1764 visitFootprints = afwDet.FootprintSet(maskVisitDet, afwDet.Threshold(1))
1765 visitDetectionFootprints.append(visitFootprints)
1766
1767 for j, footprint in enumerate(footprints.getFootprints()):
1768 ignoreArr[i, j] = countMaskFromFootprint(tmpExpMask, footprint, ignoreMask, 0x0)
1769 overlapDetArr[i, j] = countMaskFromFootprint(tmpExpMask, footprint, maskDetValue, ignoreMask)
1770
1771 # build a list of clipped spans for each visit
1772 for j, footprint in enumerate(footprints.getFootprints()):
1773 nPixel = footprint.getArea()
1774 overlap = [] # hold the overlap with each visit
1775 indexList = [] # index of visit in global list
1776 for i in range(len(tempExpRefList)):
1777 ignore = ignoreArr[i, j]
1778 overlapDet = overlapDetArr[i, j]
1779 totPixel = nPixel - ignore
1780
1781 # If we have more bad pixels than detection skip
1782 if ignore > overlapDet or totPixel <= 0.5*nPixel or overlapDet == 0:
1783 continue
1784 overlap.append(overlapDet/float(totPixel))
1785 indexList.append(i)
1786
1787 overlap = numpy.array(overlap)
1788 if not len(overlap):
1789 continue
1790
1791 keep = False # Should this footprint be marked as clipped?
1792 keepIndex = [] # Which tempExps does the clipped footprint belong to
1793
1794 # If footprint only has one overlap use a lower threshold
1795 if len(overlap) == 1:
1796 if overlap[0] > self.config.minClipFootOverlapSingle:
1797 keep = True
1798 keepIndex = [0]
1799 else:
1800 # This is the general case where only visit should be clipped
1801 clipIndex = numpy.where(overlap > self.config.minClipFootOverlap)[0]
1802 if len(clipIndex) == 1:
1803 keep = True
1804 keepIndex = [clipIndex[0]]
1805
1806 # Test if there are clipped objects that overlap two different visits
1807 clipIndex = numpy.where(overlap > self.config.minClipFootOverlapDouble)[0]
1808 if len(clipIndex) == 2 and len(overlap) > 3:
1809 clipIndexComp = numpy.where(overlap <= self.config.minClipFootOverlapDouble)[0]
1810 if numpy.max(overlap[clipIndexComp]) <= self.config.maxClipFootOverlapDouble:
1811 keep = True
1812 keepIndex = clipIndex
1813
1814 if not keep:
1815 continue
1816
1817 for index in keepIndex:
1818 globalIndex = indexList[index]
1819 artifactSpanSets[globalIndex]['CLIPPED'].append(footprint.spans)
1820
1821 clipIndices.append(numpy.array(indexList)[keepIndex])
1822 clipFootprints.append(footprint)
1823
1824 return pipeBase.Struct(clipFootprints=clipFootprints, clipIndices=clipIndices,
1825 clipSpans=artifactSpanSets, detectionFootprints=visitDetectionFootprints)
1826
1827 def detectClipBig(self, clipList, clipFootprints, clipIndices, detectionFootprints,
1828 maskClipValue, maskDetValue, coaddBBox):
1829 """Return individual warp footprints for large artifacts and append
1830 them to ``clipList`` in place.
1831
1832 Identify big footprints composed of many sources in the coadd
1833 difference that may have originated in a large diffuse source in the
1834 coadd. We do this by indentifying all clipped footprints that overlap
1835 significantly with each source in all the coaddTempExps.
1836
1837 Parameters
1838 ----------
1839 clipList : `list`
1840 List of alt mask SpanSets with clipping information. Modified.
1841 clipFootprints : `list`
1842 List of clipped footprints.
1843 clipIndices : `list`
1844 List of which entries in tempExpClipList each footprint belongs to.
1845 maskClipValue
1846 Mask value of clipped pixels.
1847 maskDetValue
1848 Mask value of detected pixels.
1849 coaddBBox : `lsst.geom.Box`
1850 BBox of the coadd and warps.
1851
1852 Returns
1853 -------
1854 bigFootprintsCoadd : `list`
1855 List of big footprints
1856 """
1857 bigFootprintsCoadd = []
1858 ignoreMask = self.getBadPixelMask()
1859 for index, (clippedSpans, visitFootprints) in enumerate(zip(clipList, detectionFootprints)):
1860 maskVisitDet = afwImage.MaskX(coaddBBox, 0x0)
1861 for footprint in visitFootprints.getFootprints():
1862 footprint.spans.setMask(maskVisitDet, maskDetValue)
1863
1864 # build a mask of clipped footprints that are in this visit
1865 clippedFootprintsVisit = []
1866 for foot, clipIndex in zip(clipFootprints, clipIndices):
1867 if index not in clipIndex:
1868 continue
1869 clippedFootprintsVisit.append(foot)
1870 maskVisitClip = maskVisitDet.Factory(maskVisitDet.getBBox(afwImage.PARENT))
1871 afwDet.setMaskFromFootprintList(maskVisitClip, clippedFootprintsVisit, maskClipValue)
1872
1873 bigFootprintsVisit = []
1874 for foot in visitFootprints.getFootprints():
1875 if foot.getArea() < self.config.minBigOverlap:
1876 continue
1877 nCount = countMaskFromFootprint(maskVisitClip, foot, maskClipValue, ignoreMask)
1878 if nCount > self.config.minBigOverlap:
1879 bigFootprintsVisit.append(foot)
1880 bigFootprintsCoadd.append(foot)
1881
1882 for footprint in bigFootprintsVisit:
1883 clippedSpans["CLIPPED"].append(footprint.spans)
1884
1885 return bigFootprintsCoadd
1886
1887
1889 psfMatchedWarps = pipeBase.connectionTypes.Input(
1890 doc=("PSF-Matched Warps are required by CompareWarp regardless of the coadd type requested. "
1891 "Only PSF-Matched Warps make sense for image subtraction. "
1892 "Therefore, they must be an additional declared input."),
1893 name="{inputCoaddName}Coadd_psfMatchedWarp",
1894 storageClass="ExposureF",
1895 dimensions=("tract", "patch", "skymap", "visit"),
1896 deferLoad=True,
1897 multiple=True
1898 )
1899 templateCoadd = pipeBase.connectionTypes.Output(
1900 doc=("Model of the static sky, used to find temporal artifacts. Typically a PSF-Matched, "
1901 "sigma-clipped coadd. Written if and only if assembleStaticSkyModel.doWrite=True"),
1902 name="{outputCoaddName}CoaddPsfMatched",
1903 storageClass="ExposureF",
1904 dimensions=("tract", "patch", "skymap", "band"),
1905 )
1906
1907 def __init__(self, *, config=None):
1908 super().__init__(config=config)
1909 if not config.assembleStaticSkyModel.doWrite:
1910 self.outputs.remove("templateCoadd")
1911 config.validate()
1912
1913
1914class CompareWarpAssembleCoaddConfig(AssembleCoaddConfig,
1915 pipelineConnections=CompareWarpAssembleCoaddConnections):
1916 assembleStaticSkyModel = pexConfig.ConfigurableField(
1917 target=AssembleCoaddTask,
1918 doc="Task to assemble an artifact-free, PSF-matched Coadd to serve as a"
1919 " naive/first-iteration model of the static sky.",
1920 )
1921 detect = pexConfig.ConfigurableField(
1922 target=SourceDetectionTask,
1923 doc="Detect outlier sources on difference between each psfMatched warp and static sky model"
1924 )
1925 detectTemplate = pexConfig.ConfigurableField(
1926 target=SourceDetectionTask,
1927 doc="Detect sources on static sky model. Only used if doPreserveContainedBySource is True"
1928 )
1929 maskStreaks = pexConfig.ConfigurableField(
1930 target=MaskStreaksTask,
1931 doc="Detect streaks on difference between each psfMatched warp and static sky model. Only used if "
1932 "doFilterMorphological is True. Adds a mask plane to an exposure, with the mask plane name set by"
1933 "streakMaskName"
1934 )
1935 streakMaskName = pexConfig.Field(
1936 dtype=str,
1937 default="STREAK",
1938 doc="Name of mask bit used for streaks"
1939 )
1940 maxNumEpochs = pexConfig.Field(
1941 doc="Charactistic maximum local number of epochs/visits in which an artifact candidate can appear "
1942 "and still be masked. The effective maxNumEpochs is a broken linear function of local "
1943 "number of epochs (N): min(maxFractionEpochsLow*N, maxNumEpochs + maxFractionEpochsHigh*N). "
1944 "For each footprint detected on the image difference between the psfMatched warp and static sky "
1945 "model, if a significant fraction of pixels (defined by spatialThreshold) are residuals in more "
1946 "than the computed effective maxNumEpochs, the artifact candidate is deemed persistant rather "
1947 "than transient and not masked.",
1948 dtype=int,
1949 default=2
1950 )
1951 maxFractionEpochsLow = pexConfig.RangeField(
1952 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for low N. "
1953 "Effective maxNumEpochs = "
1954 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1955 dtype=float,
1956 default=0.4,
1957 min=0., max=1.,
1958 )
1959 maxFractionEpochsHigh = pexConfig.RangeField(
1960 doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for high N. "
1961 "Effective maxNumEpochs = "
1962 "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1963 dtype=float,
1964 default=0.03,
1965 min=0., max=1.,
1966 )
1967 spatialThreshold = pexConfig.RangeField(
1968 doc="Unitless fraction of pixels defining how much of the outlier region has to meet the "
1969 "temporal criteria. If 0, clip all. If 1, clip none.",
1970 dtype=float,
1971 default=0.5,
1972 min=0., max=1.,
1973 inclusiveMin=True, inclusiveMax=True
1974 )
1975 doScaleWarpVariance = pexConfig.Field(
1976 doc="Rescale Warp variance plane using empirical noise?",
1977 dtype=bool,
1978 default=True,
1979 )
1980 scaleWarpVariance = pexConfig.ConfigurableField(
1981 target=ScaleVarianceTask,
1982 doc="Rescale variance on warps",
1983 )
1984 doPreserveContainedBySource = pexConfig.Field(
1985 doc="Rescue artifacts from clipping that completely lie within a footprint detected"
1986 "on the PsfMatched Template Coadd. Replicates a behavior of SafeClip.",
1987 dtype=bool,
1988 default=True,
1989 )
1990 doPrefilterArtifacts = pexConfig.Field(
1991 doc="Ignore artifact candidates that are mostly covered by the bad pixel mask, "
1992 "because they will be excluded anyway. This prevents them from contributing "
1993 "to the outlier epoch count image and potentially being labeled as persistant."
1994 "'Mostly' is defined by the config 'prefilterArtifactsRatio'.",
1995 dtype=bool,
1996 default=True
1997 )
1998 prefilterArtifactsMaskPlanes = pexConfig.ListField(
1999 doc="Prefilter artifact candidates that are mostly covered by these bad mask planes.",
2000 dtype=str,
2001 default=('NO_DATA', 'BAD', 'SAT', 'SUSPECT'),
2002 )
2003 prefilterArtifactsRatio = pexConfig.Field(
2004 doc="Prefilter artifact candidates with less than this fraction overlapping good pixels",
2005 dtype=float,
2006 default=0.05
2007 )
2008 doFilterMorphological = pexConfig.Field(
2009 doc="Filter artifact candidates based on morphological criteria, i.g. those that appear to "
2010 "be streaks.",
2011 dtype=bool,
2012 default=False
2013 )
2014 growStreakFp = pexConfig.Field(
2015 doc="Grow streak footprints by this number multiplied by the PSF width",
2016 dtype=float,
2017 default=5
2018 )
2019
2020 def setDefaults(self):
2021 AssembleCoaddConfig.setDefaults(self)
2022 self.statisticstatistic = 'MEAN'
2023 self.doUsePsfMatchedPolygonsdoUsePsfMatchedPolygons = True
2024
2025 # Real EDGE removed by psfMatched NO_DATA border half the width of the matching kernel
2026 # CompareWarp applies psfMatched EDGE pixels to directWarps before assembling
2027 if "EDGE" in self.badMaskPlanes:
2028 self.badMaskPlanes.remove('EDGE')
2029 self.removeMaskPlanes.append('EDGE')
2030 self.assembleStaticSkyModelassembleStaticSkyModel.badMaskPlanes = ["NO_DATA", ]
2031 self.assembleStaticSkyModelassembleStaticSkyModel.warpType = 'psfMatched'
2032 self.assembleStaticSkyModelassembleStaticSkyModel.connections.warpType = 'psfMatched'
2033 self.assembleStaticSkyModelassembleStaticSkyModel.statistic = 'MEANCLIP'
2034 self.assembleStaticSkyModelassembleStaticSkyModel.sigmaClip = 2.5
2035 self.assembleStaticSkyModelassembleStaticSkyModel.clipIter = 3
2036 self.assembleStaticSkyModelassembleStaticSkyModel.calcErrorFromInputVariance = False
2037 self.assembleStaticSkyModelassembleStaticSkyModel.doWrite = False
2038 self.detectdetect.doTempLocalBackground = False
2039 self.detectdetect.reEstimateBackground = False
2040 self.detectdetect.returnOriginalFootprints = False
2041 self.detectdetect.thresholdPolarity = "both"
2042 self.detectdetect.thresholdValue = 5
2043 self.detectdetect.minPixels = 4
2044 self.detectdetect.isotropicGrow = True
2045 self.detectdetect.thresholdType = "pixel_stdev"
2046 self.detectdetect.nSigmaToGrow = 0.4
2047 # The default nSigmaToGrow for SourceDetectionTask is already 2.4,
2048 # Explicitly restating because ratio with detect.nSigmaToGrow matters
2049 self.detectTemplatedetectTemplate.nSigmaToGrow = 2.4
2050 self.detectTemplatedetectTemplate.doTempLocalBackground = False
2051 self.detectTemplatedetectTemplate.reEstimateBackground = False
2052 self.detectTemplatedetectTemplate.returnOriginalFootprints = False
2053
2054 def validate(self):
2055 super().validate()
2056 if self.assembleStaticSkyModelassembleStaticSkyModel.doNImage:
2057 raise ValueError("No dataset type exists for a PSF-Matched Template N Image."
2058 "Please set assembleStaticSkyModel.doNImage=False")
2059
2060 if self.assembleStaticSkyModelassembleStaticSkyModel.doWrite and (self.warpTypewarpType == self.assembleStaticSkyModelassembleStaticSkyModel.warpType):
2061 raise ValueError("warpType (%s) == assembleStaticSkyModel.warpType (%s) and will compete for "
2062 "the same dataset name. Please set assembleStaticSkyModel.doWrite to False "
2063 "or warpType to 'direct'. assembleStaticSkyModel.warpType should ways be "
2064 "'PsfMatched'" % (self.warpTypewarpType, self.assembleStaticSkyModelassembleStaticSkyModel.warpType))
2065
2066
2067class CompareWarpAssembleCoaddTask(AssembleCoaddTask):
2068 """Assemble a compareWarp coadded image from a set of warps
2069 by masking artifacts detected by comparing PSF-matched warps.
2070
2071 In ``AssembleCoaddTask``, we compute the coadd as an clipped mean (i.e.,
2072 we clip outliers). The problem with doing this is that when computing the
2073 coadd PSF at a given location, individual visit PSFs from visits with
2074 outlier pixels contribute to the coadd PSF and cannot be treated correctly.
2075 In this task, we correct for this behavior by creating a new badMaskPlane
2076 'CLIPPED' which marks pixels in the individual warps suspected to contain
2077 an artifact. We populate this plane on the input warps by comparing
2078 PSF-matched warps with a PSF-matched median coadd which serves as a
2079 model of the static sky. Any group of pixels that deviates from the
2080 PSF-matched template coadd by more than config.detect.threshold sigma,
2081 is an artifact candidate. The candidates are then filtered to remove
2082 variable sources and sources that are difficult to subtract such as
2083 bright stars. This filter is configured using the config parameters
2084 ``temporalThreshold`` and ``spatialThreshold``. The temporalThreshold is
2085 the maximum fraction of epochs that the deviation can appear in and still
2086 be considered an artifact. The spatialThreshold is the maximum fraction of
2087 pixels in the footprint of the deviation that appear in other epochs
2088 (where other epochs is defined by the temporalThreshold). If the deviant
2089 region meets this criteria of having a significant percentage of pixels
2090 that deviate in only a few epochs, these pixels have the 'CLIPPED' bit
2091 set in the mask. These regions will not contribute to the final coadd.
2092 Furthermore, any routine to determine the coadd PSF can now be cognizant
2093 of clipped regions. Note that the algorithm implemented by this task is
2094 preliminary and works correctly for HSC data. Parameter modifications and
2095 or considerable redesigning of the algorithm is likley required for other
2096 surveys.
2097
2098 ``CompareWarpAssembleCoaddTask`` sub-classes
2099 ``AssembleCoaddTask`` and instantiates ``AssembleCoaddTask``
2100 as a subtask to generate the TemplateCoadd (the model of the static sky).
2101
2102 Notes
2103 -----
2104 The `lsst.pipe.base.cmdLineTask.CmdLineTask` interface supports a
2105 flag ``-d`` to import ``debug.py`` from your ``PYTHONPATH``; see
2106 ``baseDebug`` for more about ``debug.py`` files.
2107
2108 This task supports the following debug variables:
2109
2110 - ``saveCountIm``
2111 If True then save the Epoch Count Image as a fits file in the `figPath`
2112 - ``figPath``
2113 Path to save the debug fits images and figures
2114
2115 For example, put something like:
2116
2117 .. code-block:: python
2118
2119 import lsstDebug
2120 def DebugInfo(name):
2121 di = lsstDebug.getInfo(name)
2122 if name == "lsst.pipe.tasks.assembleCoadd":
2123 di.saveCountIm = True
2124 di.figPath = "/desired/path/to/debugging/output/images"
2125 return di
2126 lsstDebug.Info = DebugInfo
2127
2128 into your ``debug.py`` file and run ``assemebleCoadd.py`` with the
2129 ``--debug`` flag. Some subtasks may have their own debug variables;
2130 see individual Task documentation.
2131
2132 Examples
2133 --------
2134 ``CompareWarpAssembleCoaddTask`` assembles a set of warped images into a
2135 coadded image. The ``CompareWarpAssembleCoaddTask`` is invoked by running
2136 ``assembleCoadd.py`` with the flag ``--compareWarpCoadd``.
2137 Usage of ``assembleCoadd.py`` expects a data reference to the tract patch
2138 and filter to be coadded (specified using
2139 '--id = [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]')
2140 along with a list of coaddTempExps to attempt to coadd (specified using
2141 '--selectId [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]').
2142 Only the warps that cover the specified tract and patch will be coadded.
2143 A list of the available optional arguments can be obtained by calling
2144 ``assembleCoadd.py`` with the ``--help`` command line argument:
2145
2146 .. code-block:: none
2147
2148 assembleCoadd.py --help
2149
2150 To demonstrate usage of the ``CompareWarpAssembleCoaddTask`` in the larger
2151 context of multi-band processing, we will generate the HSC-I & -R band
2152 oadds from HSC engineering test data provided in the ``ci_hsc`` package.
2153 To begin, assuming that the lsst stack has been already set up, we must
2154 set up the ``obs_subaru`` and ``ci_hsc`` packages.
2155 This defines the environment variable ``$CI_HSC_DIR`` and points at the
2156 location of the package. The raw HSC data live in the ``$CI_HSC_DIR/raw``
2157 directory. To begin assembling the coadds, we must first
2158
2159 - processCcd
2160 process the individual ccds in $CI_HSC_RAW to produce calibrated exposures
2161 - makeSkyMap
2162 create a skymap that covers the area of the sky present in the raw exposures
2163 - makeCoaddTempExp
2164 warp the individual calibrated exposures to the tangent plane of the coadd
2165
2166 We can perform all of these steps by running
2167
2168 .. code-block:: none
2169
2170 $CI_HSC_DIR scons warp-903986 warp-904014 warp-903990 warp-904010 warp-903988
2171
2172 This will produce warped ``coaddTempExps`` for each visit. To coadd the
2173 warped data, we call ``assembleCoadd.py`` as follows:
2174
2175 .. code-block:: none
2176
2177 assembleCoadd.py --compareWarpCoadd $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-I \
2178 --selectId visit=903986 ccd=16 --selectId visit=903986 ccd=22 --selectId visit=903986 ccd=23 \
2179 --selectId visit=903986 ccd=100 --selectId visit=904014 ccd=1 --selectId visit=904014 ccd=6 \
2180 --selectId visit=904014 ccd=12 --selectId visit=903990 ccd=18 --selectId visit=903990 ccd=25 \
2181 --selectId visit=904010 ccd=4 --selectId visit=904010 ccd=10 --selectId visit=904010 ccd=100 \
2182 --selectId visit=903988 ccd=16 --selectId visit=903988 ccd=17 --selectId visit=903988 ccd=23 \
2183 --selectId visit=903988 ccd=24
2184
2185 This will process the HSC-I band data. The results are written in
2186 ``$CI_HSC_DIR/DATA/deepCoadd-results/HSC-I``.
2187 """
2188 ConfigClass = CompareWarpAssembleCoaddConfig
2189 _DefaultName = "compareWarpAssembleCoadd"
2190
2191 def __init__(self, *args, **kwargs):
2192 AssembleCoaddTask.__init__(self, *args, **kwargs)
2193 self.makeSubtask("assembleStaticSkyModel")
2194 detectionSchema = afwTable.SourceTable.makeMinimalSchema()
2195 self.makeSubtask("detect", schema=detectionSchema)
2196 if self.config.doPreserveContainedBySource:
2197 self.makeSubtask("detectTemplate", schema=afwTable.SourceTable.makeMinimalSchema())
2198 if self.config.doScaleWarpVariance:
2199 self.makeSubtask("scaleWarpVariance")
2200 if self.config.doFilterMorphological:
2201 self.makeSubtask("maskStreaks")
2202
2203 @utils.inheritDoc(AssembleCoaddTask)
2204 def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs):
2205 """
2206 Generate a templateCoadd to use as a naive model of static sky to
2207 subtract from PSF-Matched warps.
2208
2209 Returns
2210 -------
2211 result : `lsst.pipe.base.Struct`
2212 Result struct with components:
2213
2214 - ``templateCoadd`` : coadded exposure (``lsst.afw.image.Exposure``)
2215 - ``nImage`` : N Image (``lsst.afw.image.Image``)
2216 """
2217 # Ensure that psfMatchedWarps are used as input warps for template generation
2218 staticSkyModelInputRefs = copy.deepcopy(inputRefs)
2219 staticSkyModelInputRefs.inputWarps = inputRefs.psfMatchedWarps
2220
2221 # Because subtasks don't have connections we have to make one.
2222 # The main task's `templateCoadd` is the subtask's `coaddExposure`
2223 staticSkyModelOutputRefs = copy.deepcopy(outputRefs)
2224 if self.config.assembleStaticSkyModel.doWrite:
2225 staticSkyModelOutputRefs.coaddExposure = staticSkyModelOutputRefs.templateCoadd
2226 # Remove template coadd from both subtask's and main tasks outputs,
2227 # because it is handled by the subtask as `coaddExposure`
2228 del outputRefs.templateCoadd
2229 del staticSkyModelOutputRefs.templateCoadd
2230
2231 # A PSF-Matched nImage does not exist as a dataset type
2232 if 'nImage' in staticSkyModelOutputRefs.keys():
2233 del staticSkyModelOutputRefs.nImage
2234
2235 templateCoadd = self.assembleStaticSkyModel.runQuantum(butlerQC, staticSkyModelInputRefs,
2236 staticSkyModelOutputRefs)
2237 if templateCoadd is None:
2238 raise RuntimeError(self._noTemplateMessage_noTemplateMessage(self.assembleStaticSkyModel.warpType))
2239
2240 return pipeBase.Struct(templateCoadd=templateCoadd.coaddExposure,
2241 nImage=templateCoadd.nImage,
2242 warpRefList=templateCoadd.warpRefList,
2243 imageScalerList=templateCoadd.imageScalerList,
2244 weightList=templateCoadd.weightList)
2245
2246 @utils.inheritDoc(AssembleCoaddTask)
2247 def makeSupplementaryData(self, dataRef, selectDataList=None, warpRefList=None):
2248 """
2249 Generate a templateCoadd to use as a naive model of static sky to
2250 subtract from PSF-Matched warps.
2251
2252 Returns
2253 -------
2254 result : `lsst.pipe.base.Struct`
2255 Result struct with components:
2256
2257 - ``templateCoadd``: coadded exposure (``lsst.afw.image.Exposure``)
2258 - ``nImage``: N Image (``lsst.afw.image.Image``)
2259 """
2260 templateCoadd = self.assembleStaticSkyModel.runDataRef(dataRef, selectDataList, warpRefList)
2261 if templateCoadd is None:
2262 raise RuntimeError(self._noTemplateMessage_noTemplateMessage(self.assembleStaticSkyModel.warpType))
2263
2264 return pipeBase.Struct(templateCoadd=templateCoadd.coaddExposure,
2265 nImage=templateCoadd.nImage,
2266 warpRefList=templateCoadd.warpRefList,
2267 imageScalerList=templateCoadd.imageScalerList,
2268 weightList=templateCoadd.weightList)
2269
2270 def _noTemplateMessage(self, warpType):
2271 warpName = (warpType[0].upper() + warpType[1:])
2272 message = """No %(warpName)s warps were found to build the template coadd which is
2273 required to run CompareWarpAssembleCoaddTask. To continue assembling this type of coadd,
2274 first either rerun makeCoaddTempExp with config.make%(warpName)s=True or
2275 coaddDriver with config.makeCoadTempExp.make%(warpName)s=True, before assembleCoadd.
2276
2277 Alternatively, to use another algorithm with existing warps, retarget the CoaddDriverConfig to
2278 another algorithm like:
2279
2280 from lsst.pipe.tasks.assembleCoadd import SafeClipAssembleCoaddTask
2281 config.assemble.retarget(SafeClipAssembleCoaddTask)
2282 """ % {"warpName": warpName}
2283 return message
2284
2285 @utils.inheritDoc(AssembleCoaddTask)
2286 @timeMethod
2287 def run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
2288 supplementaryData, *args, **kwargs):
2289 """Assemble the coadd.
2290
2291 Find artifacts and apply them to the warps' masks creating a list of
2292 alternative masks with a new "CLIPPED" plane and updated "NO_DATA"
2293 plane. Then pass these alternative masks to the base class's `run`
2294 method.
2295
2296 The input parameters ``supplementaryData`` is a `lsst.pipe.base.Struct`
2297 that must contain a ``templateCoadd`` that serves as the
2298 model of the static sky.
2299 """
2300
2301 # Check and match the order of the supplementaryData
2302 # (PSF-matched) inputs to the order of the direct inputs,
2303 # so that the artifact mask is applied to the right warp
2304 dataIds = [ref.dataId for ref in tempExpRefList]
2305 psfMatchedDataIds = [ref.dataId for ref in supplementaryData.warpRefList]
2306
2307 if dataIds != psfMatchedDataIds:
2308 self.log.info("Reordering and or/padding PSF-matched visit input list")
2309 supplementaryData.warpRefList = reorderAndPadList(supplementaryData.warpRefList,
2310 psfMatchedDataIds, dataIds)
2311 supplementaryData.imageScalerList = reorderAndPadList(supplementaryData.imageScalerList,
2312 psfMatchedDataIds, dataIds)
2313
2314 # Use PSF-Matched Warps (and corresponding scalers) and coadd to find artifacts
2315 spanSetMaskList = self.findArtifactsfindArtifacts(supplementaryData.templateCoadd,
2316 supplementaryData.warpRefList,
2317 supplementaryData.imageScalerList)
2318
2319 badMaskPlanes = self.config.badMaskPlanes[:]
2320 badMaskPlanes.append("CLIPPED")
2321 badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes)
2322
2323 result = AssembleCoaddTask.run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
2324 spanSetMaskList, mask=badPixelMask)
2325
2326 # Propagate PSF-matched EDGE pixels to coadd SENSOR_EDGE and INEXACT_PSF
2327 # Psf-Matching moves the real edge inwards
2328 self.applyAltEdgeMaskapplyAltEdgeMask(result.coaddExposure.maskedImage.mask, spanSetMaskList)
2329 return result
2330
2331 def applyAltEdgeMask(self, mask, altMaskList):
2332 """Propagate alt EDGE mask to SENSOR_EDGE AND INEXACT_PSF planes.
2333
2334 Parameters
2335 ----------
2336 mask : `lsst.afw.image.Mask`
2337 Original mask.
2338 altMaskList : `list`
2339 List of Dicts containing ``spanSet`` lists.
2340 Each element contains the new mask plane name (e.g. "CLIPPED
2341 and/or "NO_DATA") as the key, and list of ``SpanSets`` to apply to
2342 the mask.
2343 """
2344 maskValue = mask.getPlaneBitMask(["SENSOR_EDGE", "INEXACT_PSF"])
2345 for visitMask in altMaskList:
2346 if "EDGE" in visitMask:
2347 for spanSet in visitMask['EDGE']:
2348 spanSet.clippedTo(mask.getBBox()).setMask(mask, maskValue)
2349
2350 def findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList):
2351 """Find artifacts.
2352
2353 Loop through warps twice. The first loop builds a map with the count
2354 of how many epochs each pixel deviates from the templateCoadd by more
2355 than ``config.chiThreshold`` sigma. The second loop takes each
2356 difference image and filters the artifacts detected in each using
2357 count map to filter out variable sources and sources that are
2358 difficult to subtract cleanly.
2359
2360 Parameters
2361 ----------
2362 templateCoadd : `lsst.afw.image.Exposure`
2363 Exposure to serve as model of static sky.
2364 tempExpRefList : `list`
2365 List of data references to warps.
2366 imageScalerList : `list`
2367 List of image scalers.
2368
2369 Returns
2370 -------
2371 altMasks : `list`
2372 List of dicts containing information about CLIPPED
2373 (i.e., artifacts), NO_DATA, and EDGE pixels.
2374 """
2375
2376 self.log.debug("Generating Count Image, and mask lists.")
2377 coaddBBox = templateCoadd.getBBox()
2378 slateIm = afwImage.ImageU(coaddBBox)
2379 epochCountImage = afwImage.ImageU(coaddBBox)
2380 nImage = afwImage.ImageU(coaddBBox)
2381 spanSetArtifactList = []
2382 spanSetNoDataMaskList = []
2383 spanSetEdgeList = []
2384 spanSetBadMorphoList = []
2385 badPixelMask = self.getBadPixelMask()
2386
2387 # mask of the warp diffs should = that of only the warp
2388 templateCoadd.mask.clearAllMaskPlanes()
2389
2390 if self.config.doPreserveContainedBySource:
2391 templateFootprints = self.detectTemplate.detectFootprints(templateCoadd)
2392 else:
2393 templateFootprints = None
2394
2395 for warpRef, imageScaler in zip(tempExpRefList, imageScalerList):
2396 warpDiffExp = self._readAndComputeWarpDiff_readAndComputeWarpDiff(warpRef, imageScaler, templateCoadd)
2397 if warpDiffExp is not None:
2398 # This nImage only approximates the final nImage because it uses the PSF-matched mask
2399 nImage.array += (numpy.isfinite(warpDiffExp.image.array)
2400 * ((warpDiffExp.mask.array & badPixelMask) == 0)).astype(numpy.uint16)
2401 fpSet = self.detect.detectFootprints(warpDiffExp, doSmooth=False, clearMask=True)
2402 fpSet.positive.merge(fpSet.negative)
2403 footprints = fpSet.positive
2404 slateIm.set(0)
2405 spanSetList = [footprint.spans for footprint in footprints.getFootprints()]
2406
2407 # Remove artifacts due to defects before they contribute to the epochCountImage
2408 if self.config.doPrefilterArtifacts:
2409 spanSetList = self.prefilterArtifactsprefilterArtifacts(spanSetList, warpDiffExp)
2410
2411 # Clear mask before adding prefiltered spanSets
2412 self.detect.clearMask(warpDiffExp.mask)
2413 for spans in spanSetList:
2414 spans.setImage(slateIm, 1, doClip=True)
2415 spans.setMask(warpDiffExp.mask, warpDiffExp.mask.getPlaneBitMask("DETECTED"))
2416 epochCountImage += slateIm
2417
2418 if self.config.doFilterMorphological:
2419 maskName = self.config.streakMaskName
2420 _ = self.maskStreaks.run(warpDiffExp)
2421 streakMask = warpDiffExp.mask
2422 spanSetStreak = afwGeom.SpanSet.fromMask(streakMask,
2423 streakMask.getPlaneBitMask(maskName)).split()
2424 # Pad the streaks to account for low-surface brightness wings
2425 psf = warpDiffExp.getPsf()
2426 for s, sset in enumerate(spanSetStreak):
2427 psfShape = psf.computeShape(sset.computeCentroid())
2428 dilation = self.config.growStreakFp * psfShape.getDeterminantRadius()
2429 sset_dilated = sset.dilated(int(dilation))
2430 spanSetStreak[s] = sset_dilated
2431
2432 # PSF-Matched warps have less available area (~the matching kernel) because the calexps
2433 # undergo a second convolution. Pixels with data in the direct warp
2434 # but not in the PSF-matched warp will not have their artifacts detected.
2435 # NaNs from the PSF-matched warp therefore must be masked in the direct warp
2436 nans = numpy.where(numpy.isnan(warpDiffExp.maskedImage.image.array), 1, 0)
2437 nansMask = afwImage.makeMaskFromArray(nans.astype(afwImage.MaskPixel))
2438 nansMask.setXY0(warpDiffExp.getXY0())
2439 edgeMask = warpDiffExp.mask
2440 spanSetEdgeMask = afwGeom.SpanSet.fromMask(edgeMask,
2441 edgeMask.getPlaneBitMask("EDGE")).split()
2442 else:
2443 # If the directWarp has <1% coverage, the psfMatchedWarp can have 0% and not exist
2444 # In this case, mask the whole epoch
2445 nansMask = afwImage.MaskX(coaddBBox, 1)
2446 spanSetList = []
2447 spanSetEdgeMask = []
2448 spanSetStreak = []
2449
2450 spanSetNoDataMask = afwGeom.SpanSet.fromMask(nansMask).split()
2451
2452 spanSetNoDataMaskList.append(spanSetNoDataMask)
2453 spanSetArtifactList.append(spanSetList)
2454 spanSetEdgeList.append(spanSetEdgeMask)
2455 if self.config.doFilterMorphological:
2456 spanSetBadMorphoList.append(spanSetStreak)
2457
2458 if lsstDebug.Info(__name__).saveCountIm:
2459 path = self._dataRef2DebugPath_dataRef2DebugPath("epochCountIm", tempExpRefList[0], coaddLevel=True)
2460 epochCountImage.writeFits(path)
2461
2462 for i, spanSetList in enumerate(spanSetArtifactList):
2463 if spanSetList:
2464 filteredSpanSetList = self.filterArtifactsfilterArtifacts(spanSetList, epochCountImage, nImage,
2465 templateFootprints)
2466 spanSetArtifactList[i] = filteredSpanSetList
2467 if self.config.doFilterMorphological:
2468 spanSetArtifactList[i] += spanSetBadMorphoList[i]
2469
2470 altMasks = []
2471 for artifacts, noData, edge in zip(spanSetArtifactList, spanSetNoDataMaskList, spanSetEdgeList):
2472 altMasks.append({'CLIPPED': artifacts,
2473 'NO_DATA': noData,
2474 'EDGE': edge})
2475 return altMasks
2476
2477 def prefilterArtifacts(self, spanSetList, exp):
2478 """Remove artifact candidates covered by bad mask plane.
2479
2480 Any future editing of the candidate list that does not depend on
2481 temporal information should go in this method.
2482
2483 Parameters
2484 ----------
2485 spanSetList : `list`
2486 List of SpanSets representing artifact candidates.
2488 Exposure containing mask planes used to prefilter.
2489
2490 Returns
2491 -------
2492 returnSpanSetList : `list`
2493 List of SpanSets with artifacts.
2494 """
2495 badPixelMask = exp.mask.getPlaneBitMask(self.config.prefilterArtifactsMaskPlanes)
2496 goodArr = (exp.mask.array & badPixelMask) == 0
2497 returnSpanSetList = []
2498 bbox = exp.getBBox()
2499 x0, y0 = exp.getXY0()
2500 for i, span in enumerate(spanSetList):
2501 y, x = span.clippedTo(bbox).indices()
2502 yIndexLocal = numpy.array(y) - y0
2503 xIndexLocal = numpy.array(x) - x0
2504 goodRatio = numpy.count_nonzero(goodArr[yIndexLocal, xIndexLocal])/span.getArea()
2505 if goodRatio > self.config.prefilterArtifactsRatio:
2506 returnSpanSetList.append(span)
2507 return returnSpanSetList
2508
2509 def filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None):
2510 """Filter artifact candidates.
2511
2512 Parameters
2513 ----------
2514 spanSetList : `list`
2515 List of SpanSets representing artifact candidates.
2516 epochCountImage : `lsst.afw.image.Image`
2517 Image of accumulated number of warpDiff detections.
2518 nImage : `lsst.afw.image.Image`
2519 Image of the accumulated number of total epochs contributing.
2520
2521 Returns
2522 -------
2523 maskSpanSetList : `list`
2524 List of SpanSets with artifacts.
2525 """
2526
2527 maskSpanSetList = []
2528 x0, y0 = epochCountImage.getXY0()
2529 for i, span in enumerate(spanSetList):
2530 y, x = span.indices()
2531 yIdxLocal = [y1 - y0 for y1 in y]
2532 xIdxLocal = [x1 - x0 for x1 in x]
2533 outlierN = epochCountImage.array[yIdxLocal, xIdxLocal]
2534 totalN = nImage.array[yIdxLocal, xIdxLocal]
2535
2536 # effectiveMaxNumEpochs is broken line (fraction of N) with characteristic config.maxNumEpochs
2537 effMaxNumEpochsHighN = (self.config.maxNumEpochs
2538 + self.config.maxFractionEpochsHigh*numpy.mean(totalN))
2539 effMaxNumEpochsLowN = self.config.maxFractionEpochsLow * numpy.mean(totalN)
2540 effectiveMaxNumEpochs = int(min(effMaxNumEpochsLowN, effMaxNumEpochsHighN))
2541 nPixelsBelowThreshold = numpy.count_nonzero((outlierN > 0)
2542 & (outlierN <= effectiveMaxNumEpochs))
2543 percentBelowThreshold = nPixelsBelowThreshold / len(outlierN)
2544 if percentBelowThreshold > self.config.spatialThreshold:
2545 maskSpanSetList.append(span)
2546
2547 if self.config.doPreserveContainedBySource and footprintsToExclude is not None:
2548 # If a candidate is contained by a footprint on the template coadd, do not clip
2549 filteredMaskSpanSetList = []
2550 for span in maskSpanSetList:
2551 doKeep = True
2552 for footprint in footprintsToExclude.positive.getFootprints():
2553 if footprint.spans.contains(span):
2554 doKeep = False
2555 break
2556 if doKeep:
2557 filteredMaskSpanSetList.append(span)
2558 maskSpanSetList = filteredMaskSpanSetList
2559
2560 return maskSpanSetList
2561
2562 def _readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd):
2563 """Fetch a warp from the butler and return a warpDiff.
2564
2565 Parameters
2566 ----------
2568 Butler dataRef for the warp.
2570 An image scaler object.
2571 templateCoadd : `lsst.afw.image.Exposure`
2572 Exposure to be substracted from the scaled warp.
2573
2574 Returns
2575 -------
2577 Exposure of the image difference between the warp and template.
2578 """
2579
2580 # If the PSF-Matched warp did not exist for this direct warp
2581 # None is holding its place to maintain order in Gen 3
2582 if warpRef is None:
2583 return None
2584 # Warp comparison must use PSF-Matched Warps regardless of requested coadd warp type
2585 warpName = self.getTempExpDatasetName('psfMatched')
2586 if not isinstance(warpRef, DeferredDatasetHandle):
2587 if not warpRef.datasetExists(warpName):
2588 self.log.warning("Could not find %s %s; skipping it", warpName, warpRef.dataId)
2589 return None
2590 warp = warpRef.get(datasetType=warpName, immediate=True)
2591 # direct image scaler OK for PSF-matched Warp
2592 imageScaler.scaleMaskedImage(warp.getMaskedImage())
2593 mi = warp.getMaskedImage()
2594 if self.config.doScaleWarpVariance:
2595 try:
2596 self.scaleWarpVariance.run(mi)
2597 except Exception as exc:
2598 self.log.warning("Unable to rescale variance of warp (%s); leaving it as-is", exc)
2599 mi -= templateCoadd.getMaskedImage()
2600 return warp
2601
2602 def _dataRef2DebugPath(self, prefix, warpRef, coaddLevel=False):
2603 """Return a path to which to write debugging output.
2604
2605 Creates a hyphen-delimited string of dataId values for simple filenames.
2606
2607 Parameters
2608 ----------
2609 prefix : `str`
2610 Prefix for filename.
2612 Butler dataRef to make the path from.
2613 coaddLevel : `bool`, optional.
2614 If True, include only coadd-level keys (e.g., 'tract', 'patch',
2615 'filter', but no 'visit').
2616
2617 Returns
2618 -------
2619 result : `str`
2620 Path for debugging output.
2621 """
2622 if coaddLevel:
2623 keys = warpRef.getButler().getKeys(self.getCoaddDatasetName(self.warpType))
2624 else:
2625 keys = warpRef.dataId.keys()
2626 keyList = sorted(keys, reverse=True)
2627 directory = lsstDebug.Info(__name__).figPath if lsstDebug.Info(__name__).figPath else "."
2628 filename = "%s-%s.fits" % (prefix, '-'.join([str(warpRef.dataId[k]) for k in keyList]))
2629 return os.path.join(directory, filename)
std::vector< SchemaItem< Flag > > * items
int min
int max
Class to describe the properties of a detected object from an image.
Definition: Footprint.h:63
A set of Footprints, associated with a MaskedImage.
Definition: FootprintSet.h:55
A Threshold is used to pass a threshold value to detection algorithms.
Definition: Threshold.h:43
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:73
Pass parameters to a Statistics object.
Definition: Statistics.h:92
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
def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs)
def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, supplementaryData, *args, **kwargs)
def findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList)
def _readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd)
def _dataRef2DebugPath(self, prefix, warpRef, coaddLevel=False)
def makeSupplementaryData(self, dataRef, selectDataList=None, warpRefList=None)
def filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None)
def buildDifferenceImage(self, skyInfo, tempExpRefList, imageScalerList, weightList)
def detectClipBig(self, clipList, clipFootprints, clipIndices, detectionFootprints, maskClipValue, maskDetValue, coaddBBox)
def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, *args, **kwargs)
Base class for coaddition.
Definition: coaddBase.py:141
daf::base::PropertyList * list
Definition: fits.cc:913
std::shared_ptr< FrameSet > append(FrameSet const &first, FrameSet const &second)
Construct a FrameSet that performs two transformations in series.
Definition: functional.cc:33
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects.
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:359
Property stringToStatisticsProperty(std::string const property)
Conversion function to switch a string to a Property (see Statistics.h)
Definition: Statistics.cc:738
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
bool all(CoordinateExpr< N > const &expr) noexcept
Return true if all elements are true.
Fit spatial kernel using approximate fluxes for candidates, and solving a linear system of equations.
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174
def reorderAndPadList(inputList, inputKeys, outputKeys, padWith=None)
def prepareStats(self, mask=None)
def readBrightObjectMasks(self, dataRef)
def makeSupplementaryData(self, dataRef, selectDataList=None, warpRefList=None)
def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask)
def assembleOnlineMeanCoadd(self, coaddExposure, tempExpRefList, imageScalerList, weightList, altMaskList, statsCtrl, nImage=None)
def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs)
def applyAltMaskPlanes(self, mask, altMaskSpans)
def shrinkValidPolygons(self, coaddInputs)
def setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None)
def getTempExpRefList(self, patchRef, calExpRefList)
def assembleMetadata(self, coaddExposure, tempExpRefList, weightList)
def removeMaskPlanes(self, maskedImage)
def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, altMaskList=None, mask=None, supplementaryData=None)
def filterWarps(self, inputs, goodVisits)
def prepareInputs(self, refList)
def assembleSubregion(self, coaddExposure, bbox, tempExpRefList, imageScalerList, weightList, altMaskList, statsFlags, statsCtrl, nImage=None)
def processResults(self, coaddExposure, brightObjectMasks=None, dataId=None)
def makeCoaddSuffix(warpType="direct")
Definition: coaddBase.py:350
def makeSkyInfo(skyMap, tractId, patchId)
Definition: coaddBase.py:293
def getGroupDataRef(butler, datasetType, groupTuple, keys)
Definition: coaddHelpers.py:99
def groupPatchExposures(patchDataRef, calexpDataRefList, coaddDatasetType="deepCoadd", tempExpDatasetType="deepCoadd_directWarp")
Definition: coaddHelpers.py:60