LSST Applications  21.0.0-147-g0e635eb1+1acddb5be5,22.0.0+052faf71bd,22.0.0+1ea9a8b2b2,22.0.0+6312710a6c,22.0.0+729191ecac,22.0.0+7589c3a021,22.0.0+9f079a9461,22.0.1-1-g7d6de66+b8044ec9de,22.0.1-1-g87000a6+536b1ee016,22.0.1-1-g8e32f31+6312710a6c,22.0.1-10-gd060f87+016f7cdc03,22.0.1-12-g9c3108e+df145f6f68,22.0.1-16-g314fa6d+c825727ab8,22.0.1-19-g93a5c75+d23f2fb6d8,22.0.1-19-gb93eaa13+aab3ef7709,22.0.1-2-g8ef0a89+b8044ec9de,22.0.1-2-g92698f7+9f079a9461,22.0.1-2-ga9b0f51+052faf71bd,22.0.1-2-gac51dbf+052faf71bd,22.0.1-2-gb66926d+6312710a6c,22.0.1-2-gcb770ba+09e3807989,22.0.1-20-g32debb5+b8044ec9de,22.0.1-23-gc2439a9a+fb0756638e,22.0.1-3-g496fd5d+09117f784f,22.0.1-3-g59f966b+1e6ba2c031,22.0.1-3-g849a1b8+f8b568069f,22.0.1-3-gaaec9c0+c5c846a8b1,22.0.1-32-g5ddfab5d3+60ce4897b0,22.0.1-4-g037fbe1+64e601228d,22.0.1-4-g8623105+b8044ec9de,22.0.1-5-g096abc9+d18c45d440,22.0.1-5-g15c806e+57f5c03693,22.0.1-7-gba73697+57f5c03693,master-g6e05de7fdc+c1283a92b8,master-g72cdda8301+729191ecac,w.2021.39
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 #
22 import os
23 import copy
24 import numpy
25 import warnings
26 import logging
27 import lsst.pex.config as pexConfig
28 import lsst.pex.exceptions as pexExceptions
29 import lsst.geom as geom
30 import lsst.afw.geom as afwGeom
31 import lsst.afw.image as afwImage
32 import lsst.afw.math as afwMath
33 import lsst.afw.table as afwTable
34 import lsst.afw.detection as afwDet
35 import lsst.coadd.utils as coaddUtils
36 import lsst.pipe.base as pipeBase
37 import lsst.meas.algorithms as measAlg
38 import lsstDebug
39 import lsst.utils as utils
40 from lsst.skymap import BaseSkyMap
41 from .coaddBase import CoaddBaseTask, SelectDataIdContainer, makeSkyInfo, makeCoaddSuffix, reorderAndPadList
42 from .interpImage import InterpImageTask
43 from .scaleZeroPoint import ScaleZeroPointTask
44 from .coaddHelpers import groupPatchExposures, getGroupDataRef
45 from .scaleVariance import ScaleVarianceTask
46 from .maskStreaks import MaskStreaksTask
47 from .healSparseMapping import HealSparseInputMapTask
48 from lsst.meas.algorithms import SourceDetectionTask, AccumulatorMeanStack
49 from lsst.daf.butler import DeferredDatasetHandle
50 
51 __all__ = ["AssembleCoaddTask", "AssembleCoaddConnections", "AssembleCoaddConfig",
52  "SafeClipAssembleCoaddTask", "SafeClipAssembleCoaddConfig",
53  "CompareWarpAssembleCoaddTask", "CompareWarpAssembleCoaddConfig"]
54 
55 log = logging.getLogger(__name__.partition(".")[2])
56 
57 
58 class 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 
139 class 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 
296 class 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  @pipeBase.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  ----------
503  dataRef : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
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.
575  dataRef : `lsst.daf.persistence.ButlerDataRef`
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  ----------
598  dataRef : `lsst.daf.persistence.ButlerDataRef`
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
744  (`lsst.afw.math.StatisticsControl`)
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  @pipeBase.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
805  ``lsst.daf.persistence.ButlerDataRef``)
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.fatal("Cannot compute online coadd %s", e)
843  else:
844  for subBBox in self._subBBoxIter(skyInfo.bbox, subregionSize):
845  try:
846  self.assembleSubregion(coaddExposure, subBBox, tempExpRefList, imageScalerList,
847  weightList, altMaskList, stats.flags, stats.ctrl,
848  nImage=nImage)
849  except Exception as e:
850  self.log.fatal("Cannot compute coadd %s: %s", subBBox, e)
851 
852  # If inputMap is requested, we must finalize the map after the accumulation.
853  if self.config.doInputMap:
854  self.inputMapper.finalize_ccd_input_map_mask()
855  inputMap = self.inputMapper.ccd_input_map
856  else:
857  inputMap = None
858 
859  self.setInexactPsf(coaddMaskedImage.getMask())
860  # Despite the name, the following doesn't really deal with "EDGE" pixels: it identifies
861  # pixels that didn't receive any unmasked inputs (as occurs around the edge of the field).
862  coaddUtils.setCoaddEdgeBits(coaddMaskedImage.getMask(), coaddMaskedImage.getVariance())
863  return pipeBase.Struct(coaddExposure=coaddExposure, nImage=nImage,
864  warpRefList=tempExpRefList, imageScalerList=imageScalerList,
865  weightList=weightList, inputMap=inputMap)
866 
867  def assembleMetadata(self, coaddExposure, tempExpRefList, weightList):
868  """Set the metadata for the coadd.
869 
870  This basic implementation sets the filter from the first input.
871 
872  Parameters
873  ----------
874  coaddExposure : `lsst.afw.image.Exposure`
875  The target exposure for the coadd.
876  tempExpRefList : `list`
877  List of data references to tempExp.
878  weightList : `list`
879  List of weights.
880  """
881  assert len(tempExpRefList) == len(weightList), "Length mismatch"
882  tempExpName = self.getTempExpDatasetName(self.warpType)
883  # We load a single pixel of each coaddTempExp, because we just want to get at the metadata
884  # (and we need more than just the PropertySet that contains the header), which is not possible
885  # with the current butler (see #2777).
886  bbox = geom.Box2I(coaddExposure.getBBox().getMin(), geom.Extent2I(1, 1))
887 
888  if isinstance(tempExpRefList[0], DeferredDatasetHandle):
889  # Gen 3 API
890  tempExpList = [tempExpRef.get(parameters={'bbox': bbox}) for tempExpRef in tempExpRefList]
891  else:
892  # Gen 2 API. Delete this when Gen 2 retired
893  tempExpList = [tempExpRef.get(tempExpName + "_sub", bbox=bbox, immediate=True)
894  for tempExpRef in tempExpRefList]
895  numCcds = sum(len(tempExp.getInfo().getCoaddInputs().ccds) for tempExp in tempExpList)
896 
897  # Set the coadd FilterLabel to the band of the first input exposure:
898  # Coadds are calibrated, so the physical label is now meaningless.
899  coaddExposure.setFilterLabel(afwImage.FilterLabel(tempExpList[0].getFilterLabel().bandLabel))
900  coaddInputs = coaddExposure.getInfo().getCoaddInputs()
901  coaddInputs.ccds.reserve(numCcds)
902  coaddInputs.visits.reserve(len(tempExpList))
903 
904  for tempExp, weight in zip(tempExpList, weightList):
905  self.inputRecorder.addVisitToCoadd(coaddInputs, tempExp, weight)
906 
907  if self.config.doUsePsfMatchedPolygons:
908  self.shrinkValidPolygons(coaddInputs)
909 
910  coaddInputs.visits.sort()
911  if self.warpType == "psfMatched":
912  # The modelPsf BBox for a psfMatchedWarp/coaddTempExp was dynamically defined by
913  # ModelPsfMatchTask as the square box bounding its spatially-variable, pre-matched WarpedPsf.
914  # Likewise, set the PSF of a PSF-Matched Coadd to the modelPsf
915  # having the maximum width (sufficient because square)
916  modelPsfList = [tempExp.getPsf() for tempExp in tempExpList]
917  modelPsfWidthList = [modelPsf.computeBBox().getWidth() for modelPsf in modelPsfList]
918  psf = modelPsfList[modelPsfWidthList.index(max(modelPsfWidthList))]
919  else:
920  psf = measAlg.CoaddPsf(coaddInputs.ccds, coaddExposure.getWcs(),
921  self.config.coaddPsf.makeControl())
922  coaddExposure.setPsf(psf)
923  apCorrMap = measAlg.makeCoaddApCorrMap(coaddInputs.ccds, coaddExposure.getBBox(afwImage.PARENT),
924  coaddExposure.getWcs())
925  coaddExposure.getInfo().setApCorrMap(apCorrMap)
926  if self.config.doAttachTransmissionCurve:
927  transmissionCurve = measAlg.makeCoaddTransmissionCurve(coaddExposure.getWcs(), coaddInputs.ccds)
928  coaddExposure.getInfo().setTransmissionCurve(transmissionCurve)
929 
930  def assembleSubregion(self, coaddExposure, bbox, tempExpRefList, imageScalerList, weightList,
931  altMaskList, statsFlags, statsCtrl, nImage=None):
932  """Assemble the coadd for a sub-region.
933 
934  For each coaddTempExp, check for (and swap in) an alternative mask
935  if one is passed. Remove mask planes listed in
936  `config.removeMaskPlanes`. Finally, stack the actual exposures using
937  `lsst.afw.math.statisticsStack` with the statistic specified by
938  statsFlags. Typically, the statsFlag will be one of lsst.afw.math.MEAN for
939  a mean-stack or `lsst.afw.math.MEANCLIP` for outlier rejection using
940  an N-sigma clipped mean where N and iterations are specified by
941  statsCtrl. Assign the stacked subregion back to the coadd.
942 
943  Parameters
944  ----------
945  coaddExposure : `lsst.afw.image.Exposure`
946  The target exposure for the coadd.
947  bbox : `lsst.geom.Box`
948  Sub-region to coadd.
949  tempExpRefList : `list`
950  List of data reference to tempExp.
951  imageScalerList : `list`
952  List of image scalers.
953  weightList : `list`
954  List of weights.
955  altMaskList : `list`
956  List of alternate masks to use rather than those stored with
957  tempExp, or None. Each element is dict with keys = mask plane
958  name to which to add the spans.
959  statsFlags : `lsst.afw.math.Property`
960  Property object for statistic for coadd.
961  statsCtrl : `lsst.afw.math.StatisticsControl`
962  Statistics control object for coadd.
963  nImage : `lsst.afw.image.ImageU`, optional
964  Keeps track of exposure count for each pixel.
965  """
966  self.log.debug("Computing coadd over %s", bbox)
967  tempExpName = self.getTempExpDatasetName(self.warpType)
968  coaddExposure.mask.addMaskPlane("REJECTED")
969  coaddExposure.mask.addMaskPlane("CLIPPED")
970  coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
971  maskMap = self.setRejectedMaskMapping(statsCtrl)
972  clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
973  maskedImageList = []
974  if nImage is not None:
975  subNImage = afwImage.ImageU(bbox.getWidth(), bbox.getHeight())
976  for tempExpRef, imageScaler, altMask in zip(tempExpRefList, imageScalerList, altMaskList):
977 
978  if isinstance(tempExpRef, DeferredDatasetHandle):
979  # Gen 3 API
980  exposure = tempExpRef.get(parameters={'bbox': bbox})
981  else:
982  # Gen 2 API. Delete this when Gen 2 retired
983  exposure = tempExpRef.get(tempExpName + "_sub", bbox=bbox)
984 
985  maskedImage = exposure.getMaskedImage()
986  mask = maskedImage.getMask()
987  if altMask is not None:
988  self.applyAltMaskPlanes(mask, altMask)
989  imageScaler.scaleMaskedImage(maskedImage)
990 
991  # Add 1 for each pixel which is not excluded by the exclude mask.
992  # In legacyCoadd, pixels may also be excluded by afwMath.statisticsStack.
993  if nImage is not None:
994  subNImage.getArray()[maskedImage.getMask().getArray() & statsCtrl.getAndMask() == 0] += 1
995  if self.config.removeMaskPlanes:
996  self.removeMaskPlanes(maskedImage)
997  maskedImageList.append(maskedImage)
998 
999  if self.config.doInputMap:
1000  visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
1001  self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
1002 
1003  with self.timer("stack"):
1004  coaddSubregion = afwMath.statisticsStack(maskedImageList, statsFlags, statsCtrl, weightList,
1005  clipped, # also set output to CLIPPED if sigma-clipped
1006  maskMap)
1007  coaddExposure.maskedImage.assign(coaddSubregion, bbox)
1008  if nImage is not None:
1009  nImage.assign(subNImage, bbox)
1010 
1011  def assembleOnlineMeanCoadd(self, coaddExposure, tempExpRefList, imageScalerList, weightList,
1012  altMaskList, statsCtrl, nImage=None):
1013  """Assemble the coadd using the "online" method.
1014 
1015  This method takes a running sum of images and weights to save memory.
1016  It only works for MEAN statistics.
1017 
1018  Parameters
1019  ----------
1020  coaddExposure : `lsst.afw.image.Exposure`
1021  The target exposure for the coadd.
1022  tempExpRefList : `list`
1023  List of data reference to tempExp.
1024  imageScalerList : `list`
1025  List of image scalers.
1026  weightList : `list`
1027  List of weights.
1028  altMaskList : `list`
1029  List of alternate masks to use rather than those stored with
1030  tempExp, or None. Each element is dict with keys = mask plane
1031  name to which to add the spans.
1032  statsCtrl : `lsst.afw.math.StatisticsControl`
1033  Statistics control object for coadd
1034  nImage : `lsst.afw.image.ImageU`, optional
1035  Keeps track of exposure count for each pixel.
1036  """
1037  self.log.debug("Computing online coadd.")
1038  tempExpName = self.getTempExpDatasetName(self.warpType)
1039  coaddExposure.mask.addMaskPlane("REJECTED")
1040  coaddExposure.mask.addMaskPlane("CLIPPED")
1041  coaddExposure.mask.addMaskPlane("SENSOR_EDGE")
1042  maskMap = self.setRejectedMaskMapping(statsCtrl)
1043  thresholdDict = AccumulatorMeanStack.stats_ctrl_to_threshold_dict(statsCtrl)
1044 
1045  bbox = coaddExposure.maskedImage.getBBox()
1046 
1047  stacker = AccumulatorMeanStack(
1048  coaddExposure.image.array.shape,
1049  statsCtrl.getAndMask(),
1050  mask_threshold_dict=thresholdDict,
1051  mask_map=maskMap,
1052  no_good_pixels_mask=statsCtrl.getNoGoodPixelsMask(),
1053  calc_error_from_input_variance=self.config.calcErrorFromInputVariance,
1054  compute_n_image=(nImage is not None)
1055  )
1056 
1057  for tempExpRef, imageScaler, altMask, weight in zip(tempExpRefList,
1058  imageScalerList,
1059  altMaskList,
1060  weightList):
1061  if isinstance(tempExpRef, DeferredDatasetHandle):
1062  # Gen 3 API
1063  exposure = tempExpRef.get()
1064  else:
1065  # Gen 2 API. Delete this when Gen 2 retired
1066  exposure = tempExpRef.get(tempExpName)
1067 
1068  maskedImage = exposure.getMaskedImage()
1069  mask = maskedImage.getMask()
1070  if altMask is not None:
1071  self.applyAltMaskPlanes(mask, altMask)
1072  imageScaler.scaleMaskedImage(maskedImage)
1073  if self.config.removeMaskPlanes:
1074  self.removeMaskPlanes(maskedImage)
1075 
1076  stacker.add_masked_image(maskedImage, weight=weight)
1077 
1078  if self.config.doInputMap:
1079  visit = exposure.getInfo().getCoaddInputs().visits[0].getId()
1080  self.inputMapper.mask_warp_bbox(bbox, visit, mask, statsCtrl.getAndMask())
1081 
1082  stacker.fill_stacked_masked_image(coaddExposure.maskedImage)
1083 
1084  if nImage is not None:
1085  nImage.array[:, :] = stacker.n_image
1086 
1087  def removeMaskPlanes(self, maskedImage):
1088  """Unset the mask of an image for mask planes specified in the config.
1089 
1090  Parameters
1091  ----------
1092  maskedImage : `lsst.afw.image.MaskedImage`
1093  The masked image to be modified.
1094  """
1095  mask = maskedImage.getMask()
1096  for maskPlane in self.config.removeMaskPlanes:
1097  try:
1098  mask &= ~mask.getPlaneBitMask(maskPlane)
1100  self.log.debug("Unable to remove mask plane %s: no mask plane with that name was found.",
1101  maskPlane)
1102 
1103  @staticmethod
1104  def setRejectedMaskMapping(statsCtrl):
1105  """Map certain mask planes of the warps to new planes for the coadd.
1106 
1107  If a pixel is rejected due to a mask value other than EDGE, NO_DATA,
1108  or CLIPPED, set it to REJECTED on the coadd.
1109  If a pixel is rejected due to EDGE, set the coadd pixel to SENSOR_EDGE.
1110  If a pixel is rejected due to CLIPPED, set the coadd pixel to CLIPPED.
1111 
1112  Parameters
1113  ----------
1114  statsCtrl : `lsst.afw.math.StatisticsControl`
1115  Statistics control object for coadd
1116 
1117  Returns
1118  -------
1119  maskMap : `list` of `tuple` of `int`
1120  A list of mappings of mask planes of the warped exposures to
1121  mask planes of the coadd.
1122  """
1123  edge = afwImage.Mask.getPlaneBitMask("EDGE")
1124  noData = afwImage.Mask.getPlaneBitMask("NO_DATA")
1125  clipped = afwImage.Mask.getPlaneBitMask("CLIPPED")
1126  toReject = statsCtrl.getAndMask() & (~noData) & (~edge) & (~clipped)
1127  maskMap = [(toReject, afwImage.Mask.getPlaneBitMask("REJECTED")),
1128  (edge, afwImage.Mask.getPlaneBitMask("SENSOR_EDGE")),
1129  (clipped, clipped)]
1130  return maskMap
1131 
1132  def applyAltMaskPlanes(self, mask, altMaskSpans):
1133  """Apply in place alt mask formatted as SpanSets to a mask.
1134 
1135  Parameters
1136  ----------
1137  mask : `lsst.afw.image.Mask`
1138  Original mask.
1139  altMaskSpans : `dict`
1140  SpanSet lists to apply. Each element contains the new mask
1141  plane name (e.g. "CLIPPED and/or "NO_DATA") as the key,
1142  and list of SpanSets to apply to the mask.
1143 
1144  Returns
1145  -------
1146  mask : `lsst.afw.image.Mask`
1147  Updated mask.
1148  """
1149  if self.config.doUsePsfMatchedPolygons:
1150  if ("NO_DATA" in altMaskSpans) and ("NO_DATA" in self.config.badMaskPlanes):
1151  # Clear away any other masks outside the validPolygons. These pixels are no longer
1152  # contributing to inexact PSFs, and will still be rejected because of NO_DATA
1153  # self.config.doUsePsfMatchedPolygons should be True only in CompareWarpAssemble
1154  # This mask-clearing step must only occur *before* applying the new masks below
1155  for spanSet in altMaskSpans['NO_DATA']:
1156  spanSet.clippedTo(mask.getBBox()).clearMask(mask, self.getBadPixelMask())
1157 
1158  for plane, spanSetList in altMaskSpans.items():
1159  maskClipValue = mask.addMaskPlane(plane)
1160  for spanSet in spanSetList:
1161  spanSet.clippedTo(mask.getBBox()).setMask(mask, 2**maskClipValue)
1162  return mask
1163 
1164  def shrinkValidPolygons(self, coaddInputs):
1165  """Shrink coaddInputs' ccds' ValidPolygons in place.
1166 
1167  Either modify each ccd's validPolygon in place, or if CoaddInputs
1168  does not have a validPolygon, create one from its bbox.
1169 
1170  Parameters
1171  ----------
1172  coaddInputs : `lsst.afw.image.coaddInputs`
1173  Original mask.
1174 
1175  """
1176  for ccd in coaddInputs.ccds:
1177  polyOrig = ccd.getValidPolygon()
1178  validPolyBBox = polyOrig.getBBox() if polyOrig else ccd.getBBox()
1179  validPolyBBox.grow(-self.config.matchingKernelSize//2)
1180  if polyOrig:
1181  validPolygon = polyOrig.intersectionSingle(validPolyBBox)
1182  else:
1183  validPolygon = afwGeom.polygon.Polygon(geom.Box2D(validPolyBBox))
1184  ccd.setValidPolygon(validPolygon)
1185 
1186  def readBrightObjectMasks(self, dataRef):
1187  """Retrieve the bright object masks.
1188 
1189  Returns None on failure.
1190 
1191  Parameters
1192  ----------
1193  dataRef : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
1194  A Butler dataRef.
1195 
1196  Returns
1197  -------
1198  result : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
1199  Bright object mask from the Butler object, or None if it cannot
1200  be retrieved.
1201  """
1202  try:
1203  return dataRef.get(datasetType="brightObjectMask", immediate=True)
1204  except Exception as e:
1205  self.log.warning("Unable to read brightObjectMask for %s: %s", dataRef.dataId, e)
1206  return None
1207 
1208  def setBrightObjectMasks(self, exposure, brightObjectMasks, dataId=None):
1209  """Set the bright object masks.
1210 
1211  Parameters
1212  ----------
1213  exposure : `lsst.afw.image.Exposure`
1214  Exposure under consideration.
1215  dataId : `lsst.daf.persistence.dataId`
1216  Data identifier dict for patch.
1217  brightObjectMasks : `lsst.afw.table`
1218  Table of bright objects to mask.
1219  """
1220 
1221  if brightObjectMasks is None:
1222  self.log.warning("Unable to apply bright object mask: none supplied")
1223  return
1224  self.log.info("Applying %d bright object masks to %s", len(brightObjectMasks), dataId)
1225  mask = exposure.getMaskedImage().getMask()
1226  wcs = exposure.getWcs()
1227  plateScale = wcs.getPixelScale().asArcseconds()
1228 
1229  for rec in brightObjectMasks:
1230  center = geom.PointI(wcs.skyToPixel(rec.getCoord()))
1231  if rec["type"] == "box":
1232  assert rec["angle"] == 0.0, ("Angle != 0 for mask object %s" % rec["id"])
1233  width = rec["width"].asArcseconds()/plateScale # convert to pixels
1234  height = rec["height"].asArcseconds()/plateScale # convert to pixels
1235 
1236  halfSize = geom.ExtentI(0.5*width, 0.5*height)
1237  bbox = geom.Box2I(center - halfSize, center + halfSize)
1238 
1239  bbox = geom.BoxI(geom.PointI(int(center[0] - 0.5*width), int(center[1] - 0.5*height)),
1240  geom.PointI(int(center[0] + 0.5*width), int(center[1] + 0.5*height)))
1241  spans = afwGeom.SpanSet(bbox)
1242  elif rec["type"] == "circle":
1243  radius = int(rec["radius"].asArcseconds()/plateScale) # convert to pixels
1244  spans = afwGeom.SpanSet.fromShape(radius, offset=center)
1245  else:
1246  self.log.warning("Unexpected region type %s at %s", rec["type"], center)
1247  continue
1248  spans.clippedTo(mask.getBBox()).setMask(mask, self.brightObjectBitmask)
1249 
1250  def setInexactPsf(self, mask):
1251  """Set INEXACT_PSF mask plane.
1252 
1253  If any of the input images isn't represented in the coadd (due to
1254  clipped pixels or chip gaps), the `CoaddPsf` will be inexact. Flag
1255  these pixels.
1256 
1257  Parameters
1258  ----------
1259  mask : `lsst.afw.image.Mask`
1260  Coadded exposure's mask, modified in-place.
1261  """
1262  mask.addMaskPlane("INEXACT_PSF")
1263  inexactPsf = mask.getPlaneBitMask("INEXACT_PSF")
1264  sensorEdge = mask.getPlaneBitMask("SENSOR_EDGE") # chip edges (so PSF is discontinuous)
1265  clipped = mask.getPlaneBitMask("CLIPPED") # pixels clipped from coadd
1266  rejected = mask.getPlaneBitMask("REJECTED") # pixels rejected from coadd due to masks
1267  array = mask.getArray()
1268  selected = array & (sensorEdge | clipped | rejected) > 0
1269  array[selected] |= inexactPsf
1270 
1271  @classmethod
1272  def _makeArgumentParser(cls):
1273  """Create an argument parser.
1274  """
1275  parser = pipeBase.ArgumentParser(name=cls._DefaultName)
1276  parser.add_id_argument("--id", cls.ConfigClass().coaddName + "Coadd_"
1277  + cls.ConfigClass().warpType + "Warp",
1278  help="data ID, e.g. --id tract=12345 patch=1,2",
1279  ContainerClass=AssembleCoaddDataIdContainer)
1280  parser.add_id_argument("--selectId", "calexp", help="data ID, e.g. --selectId visit=6789 ccd=0..9",
1281  ContainerClass=SelectDataIdContainer)
1282  return parser
1283 
1284  @staticmethod
1285  def _subBBoxIter(bbox, subregionSize):
1286  """Iterate over subregions of a bbox.
1287 
1288  Parameters
1289  ----------
1290  bbox : `lsst.geom.Box2I`
1291  Bounding box over which to iterate.
1292  subregionSize: `lsst.geom.Extent2I`
1293  Size of sub-bboxes.
1294 
1295  Yields
1296  ------
1297  subBBox : `lsst.geom.Box2I`
1298  Next sub-bounding box of size ``subregionSize`` or smaller; each ``subBBox``
1299  is contained within ``bbox``, so it may be smaller than ``subregionSize`` at
1300  the edges of ``bbox``, but it will never be empty.
1301  """
1302  if bbox.isEmpty():
1303  raise RuntimeError("bbox %s is empty" % (bbox,))
1304  if subregionSize[0] < 1 or subregionSize[1] < 1:
1305  raise RuntimeError("subregionSize %s must be nonzero" % (subregionSize,))
1306 
1307  for rowShift in range(0, bbox.getHeight(), subregionSize[1]):
1308  for colShift in range(0, bbox.getWidth(), subregionSize[0]):
1309  subBBox = geom.Box2I(bbox.getMin() + geom.Extent2I(colShift, rowShift), subregionSize)
1310  subBBox.clip(bbox)
1311  if subBBox.isEmpty():
1312  raise RuntimeError("Bug: empty bbox! bbox=%s, subregionSize=%s, "
1313  "colShift=%s, rowShift=%s" %
1314  (bbox, subregionSize, colShift, rowShift))
1315  yield subBBox
1316 
1317  def filterWarps(self, inputs, goodVisits):
1318  """Return list of only inputRefs with visitId in goodVisits ordered by goodVisit
1319 
1320  Parameters
1321  ----------
1322  inputs : list
1323  List of `lsst.pipe.base.connections.DeferredDatasetRef` with dataId containing visit
1324  goodVisit : `dict`
1325  Dictionary with good visitIds as the keys. Value ignored.
1326 
1327  Returns:
1328  --------
1329  filteredInputs : `list`
1330  Filtered and sorted list of `lsst.pipe.base.connections.DeferredDatasetRef`
1331  """
1332  inputWarpDict = {inputRef.ref.dataId['visit']: inputRef for inputRef in inputs}
1333  filteredInputs = []
1334  for visit in goodVisits.keys():
1335  if visit in inputWarpDict:
1336  filteredInputs.append(inputWarpDict[visit])
1337  return filteredInputs
1338 
1339 
1340 class AssembleCoaddDataIdContainer(pipeBase.DataIdContainer):
1341  """A version of `lsst.pipe.base.DataIdContainer` specialized for assembleCoadd.
1342  """
1343 
1344  def makeDataRefList(self, namespace):
1345  """Make self.refList from self.idList.
1346 
1347  Parameters
1348  ----------
1349  namespace
1350  Results of parsing command-line (with ``butler`` and ``log`` elements).
1351  """
1352  datasetType = namespace.config.coaddName + "Coadd"
1353  keysCoadd = namespace.butler.getKeys(datasetType=datasetType, level=self.level)
1354 
1355  for dataId in self.idList:
1356  # tract and patch are required
1357  for key in keysCoadd:
1358  if key not in dataId:
1359  raise RuntimeError("--id must include " + key)
1360 
1361  dataRef = namespace.butler.dataRef(
1362  datasetType=datasetType,
1363  dataId=dataId,
1364  )
1365  self.refList.append(dataRef)
1366 
1367 
1368 def countMaskFromFootprint(mask, footprint, bitmask, ignoreMask):
1369  """Function to count the number of pixels with a specific mask in a
1370  footprint.
1371 
1372  Find the intersection of mask & footprint. Count all pixels in the mask
1373  that are in the intersection that have bitmask set but do not have
1374  ignoreMask set. Return the count.
1375 
1376  Parameters
1377  ----------
1378  mask : `lsst.afw.image.Mask`
1379  Mask to define intersection region by.
1380  footprint : `lsst.afw.detection.Footprint`
1381  Footprint to define the intersection region by.
1382  bitmask
1383  Specific mask that we wish to count the number of occurances of.
1384  ignoreMask
1385  Pixels to not consider.
1386 
1387  Returns
1388  -------
1389  result : `int`
1390  Count of number of pixels in footprint with specified mask.
1391  """
1392  bbox = footprint.getBBox()
1393  bbox.clip(mask.getBBox(afwImage.PARENT))
1394  fp = afwImage.Mask(bbox)
1395  subMask = mask.Factory(mask, bbox, afwImage.PARENT)
1396  footprint.spans.setMask(fp, bitmask)
1397  return numpy.logical_and((subMask.getArray() & fp.getArray()) > 0,
1398  (subMask.getArray() & ignoreMask) == 0).sum()
1399 
1400 
1401 class SafeClipAssembleCoaddConfig(AssembleCoaddConfig, pipelineConnections=AssembleCoaddConnections):
1402  """Configuration parameters for the SafeClipAssembleCoaddTask.
1403  """
1404  clipDetection = pexConfig.ConfigurableField(
1405  target=SourceDetectionTask,
1406  doc="Detect sources on difference between unclipped and clipped coadd")
1407  minClipFootOverlap = pexConfig.Field(
1408  doc="Minimum fractional overlap of clipped footprint with visit DETECTED to be clipped",
1409  dtype=float,
1410  default=0.6
1411  )
1412  minClipFootOverlapSingle = pexConfig.Field(
1413  doc="Minimum fractional overlap of clipped footprint with visit DETECTED to be "
1414  "clipped when only one visit overlaps",
1415  dtype=float,
1416  default=0.5
1417  )
1418  minClipFootOverlapDouble = pexConfig.Field(
1419  doc="Minimum fractional overlap of clipped footprints with visit DETECTED to be "
1420  "clipped when two visits overlap",
1421  dtype=float,
1422  default=0.45
1423  )
1424  maxClipFootOverlapDouble = pexConfig.Field(
1425  doc="Maximum fractional overlap of clipped footprints with visit DETECTED when "
1426  "considering two visits",
1427  dtype=float,
1428  default=0.15
1429  )
1430  minBigOverlap = pexConfig.Field(
1431  doc="Minimum number of pixels in footprint to use DETECTED mask from the single visits "
1432  "when labeling clipped footprints",
1433  dtype=int,
1434  default=100
1435  )
1436 
1437  def setDefaults(self):
1438  """Set default values for clipDetection.
1439 
1440  Notes
1441  -----
1442  The numeric values for these configuration parameters were
1443  empirically determined, future work may further refine them.
1444  """
1445  AssembleCoaddConfig.setDefaults(self)
1446  self.clipDetectionclipDetection.doTempLocalBackground = False
1447  self.clipDetectionclipDetection.reEstimateBackground = False
1448  self.clipDetectionclipDetection.returnOriginalFootprints = False
1449  self.clipDetectionclipDetection.thresholdPolarity = "both"
1450  self.clipDetectionclipDetection.thresholdValue = 2
1451  self.clipDetectionclipDetection.nSigmaToGrow = 2
1452  self.clipDetectionclipDetection.minPixels = 4
1453  self.clipDetectionclipDetection.isotropicGrow = True
1454  self.clipDetectionclipDetection.thresholdType = "pixel_stdev"
1455  self.sigmaClipsigmaClip = 1.5
1456  self.clipIterclipIter = 3
1457  self.statisticstatistic = "MEAN"
1458 
1459  def validate(self):
1460  if self.doSigmaClipdoSigmaClip:
1461  log.warning("Additional Sigma-clipping not allowed in Safe-clipped Coadds. "
1462  "Ignoring doSigmaClip.")
1463  self.doSigmaClipdoSigmaClip = False
1464  if self.statisticstatistic != "MEAN":
1465  raise ValueError("Only MEAN statistic allowed for final stacking in SafeClipAssembleCoadd "
1466  "(%s chosen). Please set statistic to MEAN."
1467  % (self.statisticstatistic))
1468  AssembleCoaddTask.ConfigClass.validate(self)
1469 
1470 
1471 class SafeClipAssembleCoaddTask(AssembleCoaddTask):
1472  """Assemble a coadded image from a set of coadded temporary exposures,
1473  being careful to clip & flag areas with potential artifacts.
1474 
1475  In ``AssembleCoaddTask``, we compute the coadd as an clipped mean (i.e.,
1476  we clip outliers). The problem with doing this is that when computing the
1477  coadd PSF at a given location, individual visit PSFs from visits with
1478  outlier pixels contribute to the coadd PSF and cannot be treated correctly.
1479  In this task, we correct for this behavior by creating a new
1480  ``badMaskPlane`` 'CLIPPED'. We populate this plane on the input
1481  coaddTempExps and the final coadd where
1482 
1483  i. difference imaging suggests that there is an outlier and
1484  ii. this outlier appears on only one or two images.
1485 
1486  Such regions will not contribute to the final coadd. Furthermore, any
1487  routine to determine the coadd PSF can now be cognizant of clipped regions.
1488  Note that the algorithm implemented by this task is preliminary and works
1489  correctly for HSC data. Parameter modifications and or considerable
1490  redesigning of the algorithm is likley required for other surveys.
1491 
1492  ``SafeClipAssembleCoaddTask`` uses a ``SourceDetectionTask``
1493  "clipDetection" subtask and also sub-classes ``AssembleCoaddTask``.
1494  You can retarget the ``SourceDetectionTask`` "clipDetection" subtask
1495  if you wish.
1496 
1497  Notes
1498  -----
1499  The `lsst.pipe.base.cmdLineTask.CmdLineTask` interface supports a
1500  flag ``-d`` to import ``debug.py`` from your ``PYTHONPATH``;
1501  see `baseDebug` for more about ``debug.py`` files.
1502  `SafeClipAssembleCoaddTask` has no debug variables of its own.
1503  The ``SourceDetectionTask`` "clipDetection" subtasks may support debug
1504  variables. See the documetation for `SourceDetectionTask` "clipDetection"
1505  for further information.
1506 
1507  Examples
1508  --------
1509  `SafeClipAssembleCoaddTask` assembles a set of warped ``coaddTempExp``
1510  images into a coadded image. The `SafeClipAssembleCoaddTask` is invoked by
1511  running assembleCoadd.py *without* the flag '--legacyCoadd'.
1512 
1513  Usage of ``assembleCoadd.py`` expects a data reference to the tract patch
1514  and filter to be coadded (specified using
1515  '--id = [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]')
1516  along with a list of coaddTempExps to attempt to coadd (specified using
1517  '--selectId [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]').
1518  Only the coaddTempExps that cover the specified tract and patch will be
1519  coadded. A list of the available optional arguments can be obtained by
1520  calling assembleCoadd.py with the --help command line argument:
1521 
1522  .. code-block:: none
1523 
1524  assembleCoadd.py --help
1525 
1526  To demonstrate usage of the `SafeClipAssembleCoaddTask` in the larger
1527  context of multi-band processing, we will generate the HSC-I & -R band
1528  coadds from HSC engineering test data provided in the ci_hsc package.
1529  To begin, assuming that the lsst stack has been already set up, we must
1530  set up the obs_subaru and ci_hsc packages. This defines the environment
1531  variable $CI_HSC_DIR and points at the location of the package. The raw
1532  HSC data live in the ``$CI_HSC_DIR/raw`` directory. To begin assembling
1533  the coadds, we must first
1534 
1535  - ``processCcd``
1536  process the individual ccds in $CI_HSC_RAW to produce calibrated exposures
1537  - ``makeSkyMap``
1538  create a skymap that covers the area of the sky present in the raw exposures
1539  - ``makeCoaddTempExp``
1540  warp the individual calibrated exposures to the tangent plane of the coadd</DD>
1541 
1542  We can perform all of these steps by running
1543 
1544  .. code-block:: none
1545 
1546  $CI_HSC_DIR scons warp-903986 warp-904014 warp-903990 warp-904010 warp-903988
1547 
1548  This will produce warped coaddTempExps for each visit. To coadd the
1549  warped data, we call ``assembleCoadd.py`` as follows:
1550 
1551  .. code-block:: none
1552 
1553  assembleCoadd.py $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-I \
1554  --selectId visit=903986 ccd=16 --selectId visit=903986 ccd=22 --selectId visit=903986 ccd=23 \
1555  --selectId visit=903986 ccd=100--selectId visit=904014 ccd=1 --selectId visit=904014 ccd=6 \
1556  --selectId visit=904014 ccd=12 --selectId visit=903990 ccd=18 --selectId visit=903990 ccd=25 \
1557  --selectId visit=904010 ccd=4 --selectId visit=904010 ccd=10 --selectId visit=904010 ccd=100 \
1558  --selectId visit=903988 ccd=16 --selectId visit=903988 ccd=17 --selectId visit=903988 ccd=23 \
1559  --selectId visit=903988 ccd=24
1560 
1561  This will process the HSC-I band data. The results are written in
1562  ``$CI_HSC_DIR/DATA/deepCoadd-results/HSC-I``.
1563 
1564  You may also choose to run:
1565 
1566  .. code-block:: none
1567 
1568  scons warp-903334 warp-903336 warp-903338 warp-903342 warp-903344 warp-903346 nnn
1569  assembleCoadd.py $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-R --selectId visit=903334 ccd=16 \
1570  --selectId visit=903334 ccd=22 --selectId visit=903334 ccd=23 --selectId visit=903334 ccd=100 \
1571  --selectId visit=903336 ccd=17 --selectId visit=903336 ccd=24 --selectId visit=903338 ccd=18 \
1572  --selectId visit=903338 ccd=25 --selectId visit=903342 ccd=4 --selectId visit=903342 ccd=10 \
1573  --selectId visit=903342 ccd=100 --selectId visit=903344 ccd=0 --selectId visit=903344 ccd=5 \
1574  --selectId visit=903344 ccd=11 --selectId visit=903346 ccd=1 --selectId visit=903346 ccd=6 \
1575  --selectId visit=903346 ccd=12
1576 
1577  to generate the coadd for the HSC-R band if you are interested in following
1578  multiBand Coadd processing as discussed in ``pipeTasks_multiBand``.
1579  """
1580  ConfigClass = SafeClipAssembleCoaddConfig
1581  _DefaultName = "safeClipAssembleCoadd"
1582 
1583  def __init__(self, *args, **kwargs):
1584  AssembleCoaddTask.__init__(self, *args, **kwargs)
1585  schema = afwTable.SourceTable.makeMinimalSchema()
1586  self.makeSubtask("clipDetection", schema=schema)
1587 
1588  @utils.inheritDoc(AssembleCoaddTask)
1589  @pipeBase.timeMethod
1590  def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, *args, **kwargs):
1591  """Assemble the coadd for a region.
1592 
1593  Compute the difference of coadds created with and without outlier
1594  rejection to identify coadd pixels that have outlier values in some
1595  individual visits.
1596  Detect clipped regions on the difference image and mark these regions
1597  on the one or two individual coaddTempExps where they occur if there
1598  is significant overlap between the clipped region and a source. This
1599  leaves us with a set of footprints from the difference image that have
1600  been identified as having occured on just one or two individual visits.
1601  However, these footprints were generated from a difference image. It
1602  is conceivable for a large diffuse source to have become broken up
1603  into multiple footprints acrosss the coadd difference in this process.
1604  Determine the clipped region from all overlapping footprints from the
1605  detected sources in each visit - these are big footprints.
1606  Combine the small and big clipped footprints and mark them on a new
1607  bad mask plane.
1608  Generate the coadd using `AssembleCoaddTask.run` without outlier
1609  removal. Clipped footprints will no longer make it into the coadd
1610  because they are marked in the new bad mask plane.
1611 
1612  Notes
1613  -----
1614  args and kwargs are passed but ignored in order to match the call
1615  signature expected by the parent task.
1616  """
1617  exp = self.buildDifferenceImagebuildDifferenceImage(skyInfo, tempExpRefList, imageScalerList, weightList)
1618  mask = exp.getMaskedImage().getMask()
1619  mask.addMaskPlane("CLIPPED")
1620 
1621  result = self.detectClipdetectClip(exp, tempExpRefList)
1622 
1623  self.log.info('Found %d clipped objects', len(result.clipFootprints))
1624 
1625  maskClipValue = mask.getPlaneBitMask("CLIPPED")
1626  maskDetValue = mask.getPlaneBitMask("DETECTED") | mask.getPlaneBitMask("DETECTED_NEGATIVE")
1627  # Append big footprints from individual Warps to result.clipSpans
1628  bigFootprints = self.detectClipBigdetectClipBig(result.clipSpans, result.clipFootprints, result.clipIndices,
1629  result.detectionFootprints, maskClipValue, maskDetValue,
1630  exp.getBBox())
1631  # Create mask of the current clipped footprints
1632  maskClip = mask.Factory(mask.getBBox(afwImage.PARENT))
1633  afwDet.setMaskFromFootprintList(maskClip, result.clipFootprints, maskClipValue)
1634 
1635  maskClipBig = maskClip.Factory(mask.getBBox(afwImage.PARENT))
1636  afwDet.setMaskFromFootprintList(maskClipBig, bigFootprints, maskClipValue)
1637  maskClip |= maskClipBig
1638 
1639  # Assemble coadd from base class, but ignoring CLIPPED pixels
1640  badMaskPlanes = self.config.badMaskPlanes[:]
1641  badMaskPlanes.append("CLIPPED")
1642  badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes)
1643  return AssembleCoaddTask.run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
1644  result.clipSpans, mask=badPixelMask)
1645 
1646  def buildDifferenceImage(self, skyInfo, tempExpRefList, imageScalerList, weightList):
1647  """Return an exposure that contains the difference between unclipped
1648  and clipped coadds.
1649 
1650  Generate a difference image between clipped and unclipped coadds.
1651  Compute the difference image by subtracting an outlier-clipped coadd
1652  from an outlier-unclipped coadd. Return the difference image.
1653 
1654  Parameters
1655  ----------
1656  skyInfo : `lsst.pipe.base.Struct`
1657  Patch geometry information, from getSkyInfo
1658  tempExpRefList : `list`
1659  List of data reference to tempExp
1660  imageScalerList : `list`
1661  List of image scalers
1662  weightList : `list`
1663  List of weights
1664 
1665  Returns
1666  -------
1667  exp : `lsst.afw.image.Exposure`
1668  Difference image of unclipped and clipped coadd wrapped in an Exposure
1669  """
1670  config = AssembleCoaddConfig()
1671  # getattr necessary because subtasks do not survive Config.toDict()
1672  # exclude connections because the class of self.config.connections is not
1673  # the same as AssembleCoaddConfig.connections, and the connections are not
1674  # needed to run this task anyway.
1675  configIntersection = {k: getattr(self.config, k)
1676  for k, v in self.config.toDict().items()
1677  if (k in config.keys() and k != "connections")}
1678  configIntersection['doInputMap'] = False
1679  configIntersection['doNImage'] = False
1680  config.update(**configIntersection)
1681 
1682  # statistic MEAN copied from self.config.statistic, but for clarity explicitly assign
1683  config.statistic = 'MEAN'
1684  task = AssembleCoaddTask(config=config)
1685  coaddMean = task.run(skyInfo, tempExpRefList, imageScalerList, weightList).coaddExposure
1686 
1687  config.statistic = 'MEANCLIP'
1688  task = AssembleCoaddTask(config=config)
1689  coaddClip = task.run(skyInfo, tempExpRefList, imageScalerList, weightList).coaddExposure
1690 
1691  coaddDiff = coaddMean.getMaskedImage().Factory(coaddMean.getMaskedImage())
1692  coaddDiff -= coaddClip.getMaskedImage()
1693  exp = afwImage.ExposureF(coaddDiff)
1694  exp.setPsf(coaddMean.getPsf())
1695  return exp
1696 
1697  def detectClip(self, exp, tempExpRefList):
1698  """Detect clipped regions on an exposure and set the mask on the
1699  individual tempExp masks.
1700 
1701  Detect footprints in the difference image after smoothing the
1702  difference image with a Gaussian kernal. Identify footprints that
1703  overlap with one or two input ``coaddTempExps`` by comparing the
1704  computed overlap fraction to thresholds set in the config. A different
1705  threshold is applied depending on the number of overlapping visits
1706  (restricted to one or two). If the overlap exceeds the thresholds,
1707  the footprint is considered "CLIPPED" and is marked as such on the
1708  coaddTempExp. Return a struct with the clipped footprints, the indices
1709  of the ``coaddTempExps`` that end up overlapping with the clipped
1710  footprints, and a list of new masks for the ``coaddTempExps``.
1711 
1712  Parameters
1713  ----------
1714  exp : `lsst.afw.image.Exposure`
1715  Exposure to run detection on.
1716  tempExpRefList : `list`
1717  List of data reference to tempExp.
1718 
1719  Returns
1720  -------
1721  result : `lsst.pipe.base.Struct`
1722  Result struct with components:
1723 
1724  - ``clipFootprints``: list of clipped footprints.
1725  - ``clipIndices``: indices for each ``clippedFootprint`` in
1726  ``tempExpRefList``.
1727  - ``clipSpans``: List of dictionaries containing spanSet lists
1728  to clip. Each element contains the new maskplane name
1729  ("CLIPPED") as the key and list of ``SpanSets`` as the value.
1730  - ``detectionFootprints``: List of DETECTED/DETECTED_NEGATIVE plane
1731  compressed into footprints.
1732  """
1733  mask = exp.getMaskedImage().getMask()
1734  maskDetValue = mask.getPlaneBitMask("DETECTED") | mask.getPlaneBitMask("DETECTED_NEGATIVE")
1735  fpSet = self.clipDetection.detectFootprints(exp, doSmooth=True, clearMask=True)
1736  # Merge positive and negative together footprints together
1737  fpSet.positive.merge(fpSet.negative)
1738  footprints = fpSet.positive
1739  self.log.info('Found %d potential clipped objects', len(footprints.getFootprints()))
1740  ignoreMask = self.getBadPixelMask()
1741 
1742  clipFootprints = []
1743  clipIndices = []
1744  artifactSpanSets = [{'CLIPPED': list()} for _ in tempExpRefList]
1745 
1746  # for use by detectClipBig
1747  visitDetectionFootprints = []
1748 
1749  dims = [len(tempExpRefList), len(footprints.getFootprints())]
1750  overlapDetArr = numpy.zeros(dims, dtype=numpy.uint16)
1751  ignoreArr = numpy.zeros(dims, dtype=numpy.uint16)
1752 
1753  # Loop over masks once and extract/store only relevant overlap metrics and detection footprints
1754  for i, warpRef in enumerate(tempExpRefList):
1755  tmpExpMask = warpRef.get(datasetType=self.getTempExpDatasetName(self.warpType),
1756  immediate=True).getMaskedImage().getMask()
1757  maskVisitDet = tmpExpMask.Factory(tmpExpMask, tmpExpMask.getBBox(afwImage.PARENT),
1758  afwImage.PARENT, True)
1759  maskVisitDet &= maskDetValue
1760  visitFootprints = afwDet.FootprintSet(maskVisitDet, afwDet.Threshold(1))
1761  visitDetectionFootprints.append(visitFootprints)
1762 
1763  for j, footprint in enumerate(footprints.getFootprints()):
1764  ignoreArr[i, j] = countMaskFromFootprint(tmpExpMask, footprint, ignoreMask, 0x0)
1765  overlapDetArr[i, j] = countMaskFromFootprint(tmpExpMask, footprint, maskDetValue, ignoreMask)
1766 
1767  # build a list of clipped spans for each visit
1768  for j, footprint in enumerate(footprints.getFootprints()):
1769  nPixel = footprint.getArea()
1770  overlap = [] # hold the overlap with each visit
1771  indexList = [] # index of visit in global list
1772  for i in range(len(tempExpRefList)):
1773  ignore = ignoreArr[i, j]
1774  overlapDet = overlapDetArr[i, j]
1775  totPixel = nPixel - ignore
1776 
1777  # If we have more bad pixels than detection skip
1778  if ignore > overlapDet or totPixel <= 0.5*nPixel or overlapDet == 0:
1779  continue
1780  overlap.append(overlapDet/float(totPixel))
1781  indexList.append(i)
1782 
1783  overlap = numpy.array(overlap)
1784  if not len(overlap):
1785  continue
1786 
1787  keep = False # Should this footprint be marked as clipped?
1788  keepIndex = [] # Which tempExps does the clipped footprint belong to
1789 
1790  # If footprint only has one overlap use a lower threshold
1791  if len(overlap) == 1:
1792  if overlap[0] > self.config.minClipFootOverlapSingle:
1793  keep = True
1794  keepIndex = [0]
1795  else:
1796  # This is the general case where only visit should be clipped
1797  clipIndex = numpy.where(overlap > self.config.minClipFootOverlap)[0]
1798  if len(clipIndex) == 1:
1799  keep = True
1800  keepIndex = [clipIndex[0]]
1801 
1802  # Test if there are clipped objects that overlap two different visits
1803  clipIndex = numpy.where(overlap > self.config.minClipFootOverlapDouble)[0]
1804  if len(clipIndex) == 2 and len(overlap) > 3:
1805  clipIndexComp = numpy.where(overlap <= self.config.minClipFootOverlapDouble)[0]
1806  if numpy.max(overlap[clipIndexComp]) <= self.config.maxClipFootOverlapDouble:
1807  keep = True
1808  keepIndex = clipIndex
1809 
1810  if not keep:
1811  continue
1812 
1813  for index in keepIndex:
1814  globalIndex = indexList[index]
1815  artifactSpanSets[globalIndex]['CLIPPED'].append(footprint.spans)
1816 
1817  clipIndices.append(numpy.array(indexList)[keepIndex])
1818  clipFootprints.append(footprint)
1819 
1820  return pipeBase.Struct(clipFootprints=clipFootprints, clipIndices=clipIndices,
1821  clipSpans=artifactSpanSets, detectionFootprints=visitDetectionFootprints)
1822 
1823  def detectClipBig(self, clipList, clipFootprints, clipIndices, detectionFootprints,
1824  maskClipValue, maskDetValue, coaddBBox):
1825  """Return individual warp footprints for large artifacts and append
1826  them to ``clipList`` in place.
1827 
1828  Identify big footprints composed of many sources in the coadd
1829  difference that may have originated in a large diffuse source in the
1830  coadd. We do this by indentifying all clipped footprints that overlap
1831  significantly with each source in all the coaddTempExps.
1832 
1833  Parameters
1834  ----------
1835  clipList : `list`
1836  List of alt mask SpanSets with clipping information. Modified.
1837  clipFootprints : `list`
1838  List of clipped footprints.
1839  clipIndices : `list`
1840  List of which entries in tempExpClipList each footprint belongs to.
1841  maskClipValue
1842  Mask value of clipped pixels.
1843  maskDetValue
1844  Mask value of detected pixels.
1845  coaddBBox : `lsst.geom.Box`
1846  BBox of the coadd and warps.
1847 
1848  Returns
1849  -------
1850  bigFootprintsCoadd : `list`
1851  List of big footprints
1852  """
1853  bigFootprintsCoadd = []
1854  ignoreMask = self.getBadPixelMask()
1855  for index, (clippedSpans, visitFootprints) in enumerate(zip(clipList, detectionFootprints)):
1856  maskVisitDet = afwImage.MaskX(coaddBBox, 0x0)
1857  for footprint in visitFootprints.getFootprints():
1858  footprint.spans.setMask(maskVisitDet, maskDetValue)
1859 
1860  # build a mask of clipped footprints that are in this visit
1861  clippedFootprintsVisit = []
1862  for foot, clipIndex in zip(clipFootprints, clipIndices):
1863  if index not in clipIndex:
1864  continue
1865  clippedFootprintsVisit.append(foot)
1866  maskVisitClip = maskVisitDet.Factory(maskVisitDet.getBBox(afwImage.PARENT))
1867  afwDet.setMaskFromFootprintList(maskVisitClip, clippedFootprintsVisit, maskClipValue)
1868 
1869  bigFootprintsVisit = []
1870  for foot in visitFootprints.getFootprints():
1871  if foot.getArea() < self.config.minBigOverlap:
1872  continue
1873  nCount = countMaskFromFootprint(maskVisitClip, foot, maskClipValue, ignoreMask)
1874  if nCount > self.config.minBigOverlap:
1875  bigFootprintsVisit.append(foot)
1876  bigFootprintsCoadd.append(foot)
1877 
1878  for footprint in bigFootprintsVisit:
1879  clippedSpans["CLIPPED"].append(footprint.spans)
1880 
1881  return bigFootprintsCoadd
1882 
1883 
1885  psfMatchedWarps = pipeBase.connectionTypes.Input(
1886  doc=("PSF-Matched Warps are required by CompareWarp regardless of the coadd type requested. "
1887  "Only PSF-Matched Warps make sense for image subtraction. "
1888  "Therefore, they must be an additional declared input."),
1889  name="{inputCoaddName}Coadd_psfMatchedWarp",
1890  storageClass="ExposureF",
1891  dimensions=("tract", "patch", "skymap", "visit"),
1892  deferLoad=True,
1893  multiple=True
1894  )
1895  templateCoadd = pipeBase.connectionTypes.Output(
1896  doc=("Model of the static sky, used to find temporal artifacts. Typically a PSF-Matched, "
1897  "sigma-clipped coadd. Written if and only if assembleStaticSkyModel.doWrite=True"),
1898  name="{outputCoaddName}CoaddPsfMatched",
1899  storageClass="ExposureF",
1900  dimensions=("tract", "patch", "skymap", "band"),
1901  )
1902 
1903  def __init__(self, *, config=None):
1904  super().__init__(config=config)
1905  if not config.assembleStaticSkyModel.doWrite:
1906  self.outputs.remove("templateCoadd")
1907  config.validate()
1908 
1909 
1910 class CompareWarpAssembleCoaddConfig(AssembleCoaddConfig,
1911  pipelineConnections=CompareWarpAssembleCoaddConnections):
1912  assembleStaticSkyModel = pexConfig.ConfigurableField(
1913  target=AssembleCoaddTask,
1914  doc="Task to assemble an artifact-free, PSF-matched Coadd to serve as a"
1915  " naive/first-iteration model of the static sky.",
1916  )
1917  detect = pexConfig.ConfigurableField(
1918  target=SourceDetectionTask,
1919  doc="Detect outlier sources on difference between each psfMatched warp and static sky model"
1920  )
1921  detectTemplate = pexConfig.ConfigurableField(
1922  target=SourceDetectionTask,
1923  doc="Detect sources on static sky model. Only used if doPreserveContainedBySource is True"
1924  )
1925  maskStreaks = pexConfig.ConfigurableField(
1926  target=MaskStreaksTask,
1927  doc="Detect streaks on difference between each psfMatched warp and static sky model. Only used if "
1928  "doFilterMorphological is True. Adds a mask plane to an exposure, with the mask plane name set by"
1929  "streakMaskName"
1930  )
1931  streakMaskName = pexConfig.Field(
1932  dtype=str,
1933  default="STREAK",
1934  doc="Name of mask bit used for streaks"
1935  )
1936  maxNumEpochs = pexConfig.Field(
1937  doc="Charactistic maximum local number of epochs/visits in which an artifact candidate can appear "
1938  "and still be masked. The effective maxNumEpochs is a broken linear function of local "
1939  "number of epochs (N): min(maxFractionEpochsLow*N, maxNumEpochs + maxFractionEpochsHigh*N). "
1940  "For each footprint detected on the image difference between the psfMatched warp and static sky "
1941  "model, if a significant fraction of pixels (defined by spatialThreshold) are residuals in more "
1942  "than the computed effective maxNumEpochs, the artifact candidate is deemed persistant rather "
1943  "than transient and not masked.",
1944  dtype=int,
1945  default=2
1946  )
1947  maxFractionEpochsLow = pexConfig.RangeField(
1948  doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for low N. "
1949  "Effective maxNumEpochs = "
1950  "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1951  dtype=float,
1952  default=0.4,
1953  min=0., max=1.,
1954  )
1955  maxFractionEpochsHigh = pexConfig.RangeField(
1956  doc="Fraction of local number of epochs (N) to use as effective maxNumEpochs for high N. "
1957  "Effective maxNumEpochs = "
1958  "min(maxFractionEpochsLow * N, maxNumEpochs + maxFractionEpochsHigh * N)",
1959  dtype=float,
1960  default=0.03,
1961  min=0., max=1.,
1962  )
1963  spatialThreshold = pexConfig.RangeField(
1964  doc="Unitless fraction of pixels defining how much of the outlier region has to meet the "
1965  "temporal criteria. If 0, clip all. If 1, clip none.",
1966  dtype=float,
1967  default=0.5,
1968  min=0., max=1.,
1969  inclusiveMin=True, inclusiveMax=True
1970  )
1971  doScaleWarpVariance = pexConfig.Field(
1972  doc="Rescale Warp variance plane using empirical noise?",
1973  dtype=bool,
1974  default=True,
1975  )
1976  scaleWarpVariance = pexConfig.ConfigurableField(
1977  target=ScaleVarianceTask,
1978  doc="Rescale variance on warps",
1979  )
1980  doPreserveContainedBySource = pexConfig.Field(
1981  doc="Rescue artifacts from clipping that completely lie within a footprint detected"
1982  "on the PsfMatched Template Coadd. Replicates a behavior of SafeClip.",
1983  dtype=bool,
1984  default=True,
1985  )
1986  doPrefilterArtifacts = pexConfig.Field(
1987  doc="Ignore artifact candidates that are mostly covered by the bad pixel mask, "
1988  "because they will be excluded anyway. This prevents them from contributing "
1989  "to the outlier epoch count image and potentially being labeled as persistant."
1990  "'Mostly' is defined by the config 'prefilterArtifactsRatio'.",
1991  dtype=bool,
1992  default=True
1993  )
1994  prefilterArtifactsMaskPlanes = pexConfig.ListField(
1995  doc="Prefilter artifact candidates that are mostly covered by these bad mask planes.",
1996  dtype=str,
1997  default=('NO_DATA', 'BAD', 'SAT', 'SUSPECT'),
1998  )
1999  prefilterArtifactsRatio = pexConfig.Field(
2000  doc="Prefilter artifact candidates with less than this fraction overlapping good pixels",
2001  dtype=float,
2002  default=0.05
2003  )
2004  doFilterMorphological = pexConfig.Field(
2005  doc="Filter artifact candidates based on morphological criteria, i.g. those that appear to "
2006  "be streaks.",
2007  dtype=bool,
2008  default=False
2009  )
2010 
2011  def setDefaults(self):
2012  AssembleCoaddConfig.setDefaults(self)
2013  self.statisticstatistic = 'MEAN'
2014  self.doUsePsfMatchedPolygonsdoUsePsfMatchedPolygons = True
2015 
2016  # Real EDGE removed by psfMatched NO_DATA border half the width of the matching kernel
2017  # CompareWarp applies psfMatched EDGE pixels to directWarps before assembling
2018  if "EDGE" in self.badMaskPlanes:
2019  self.badMaskPlanes.remove('EDGE')
2020  self.removeMaskPlanes.append('EDGE')
2021  self.assembleStaticSkyModelassembleStaticSkyModel.badMaskPlanes = ["NO_DATA", ]
2022  self.assembleStaticSkyModelassembleStaticSkyModel.warpType = 'psfMatched'
2023  self.assembleStaticSkyModelassembleStaticSkyModel.connections.warpType = 'psfMatched'
2024  self.assembleStaticSkyModelassembleStaticSkyModel.statistic = 'MEANCLIP'
2025  self.assembleStaticSkyModelassembleStaticSkyModel.sigmaClip = 2.5
2026  self.assembleStaticSkyModelassembleStaticSkyModel.clipIter = 3
2027  self.assembleStaticSkyModelassembleStaticSkyModel.calcErrorFromInputVariance = False
2028  self.assembleStaticSkyModelassembleStaticSkyModel.doWrite = False
2029  self.detectdetect.doTempLocalBackground = False
2030  self.detectdetect.reEstimateBackground = False
2031  self.detectdetect.returnOriginalFootprints = False
2032  self.detectdetect.thresholdPolarity = "both"
2033  self.detectdetect.thresholdValue = 5
2034  self.detectdetect.minPixels = 4
2035  self.detectdetect.isotropicGrow = True
2036  self.detectdetect.thresholdType = "pixel_stdev"
2037  self.detectdetect.nSigmaToGrow = 0.4
2038  # The default nSigmaToGrow for SourceDetectionTask is already 2.4,
2039  # Explicitly restating because ratio with detect.nSigmaToGrow matters
2040  self.detectTemplatedetectTemplate.nSigmaToGrow = 2.4
2041  self.detectTemplatedetectTemplate.doTempLocalBackground = False
2042  self.detectTemplatedetectTemplate.reEstimateBackground = False
2043  self.detectTemplatedetectTemplate.returnOriginalFootprints = False
2044 
2045  def validate(self):
2046  super().validate()
2047  if self.assembleStaticSkyModelassembleStaticSkyModel.doNImage:
2048  raise ValueError("No dataset type exists for a PSF-Matched Template N Image."
2049  "Please set assembleStaticSkyModel.doNImage=False")
2050 
2051  if self.assembleStaticSkyModelassembleStaticSkyModel.doWrite and (self.warpTypewarpType == self.assembleStaticSkyModelassembleStaticSkyModel.warpType):
2052  raise ValueError("warpType (%s) == assembleStaticSkyModel.warpType (%s) and will compete for "
2053  "the same dataset name. Please set assembleStaticSkyModel.doWrite to False "
2054  "or warpType to 'direct'. assembleStaticSkyModel.warpType should ways be "
2055  "'PsfMatched'" % (self.warpTypewarpType, self.assembleStaticSkyModelassembleStaticSkyModel.warpType))
2056 
2057 
2058 class CompareWarpAssembleCoaddTask(AssembleCoaddTask):
2059  """Assemble a compareWarp coadded image from a set of warps
2060  by masking artifacts detected by comparing PSF-matched warps.
2061 
2062  In ``AssembleCoaddTask``, we compute the coadd as an clipped mean (i.e.,
2063  we clip outliers). The problem with doing this is that when computing the
2064  coadd PSF at a given location, individual visit PSFs from visits with
2065  outlier pixels contribute to the coadd PSF and cannot be treated correctly.
2066  In this task, we correct for this behavior by creating a new badMaskPlane
2067  'CLIPPED' which marks pixels in the individual warps suspected to contain
2068  an artifact. We populate this plane on the input warps by comparing
2069  PSF-matched warps with a PSF-matched median coadd which serves as a
2070  model of the static sky. Any group of pixels that deviates from the
2071  PSF-matched template coadd by more than config.detect.threshold sigma,
2072  is an artifact candidate. The candidates are then filtered to remove
2073  variable sources and sources that are difficult to subtract such as
2074  bright stars. This filter is configured using the config parameters
2075  ``temporalThreshold`` and ``spatialThreshold``. The temporalThreshold is
2076  the maximum fraction of epochs that the deviation can appear in and still
2077  be considered an artifact. The spatialThreshold is the maximum fraction of
2078  pixels in the footprint of the deviation that appear in other epochs
2079  (where other epochs is defined by the temporalThreshold). If the deviant
2080  region meets this criteria of having a significant percentage of pixels
2081  that deviate in only a few epochs, these pixels have the 'CLIPPED' bit
2082  set in the mask. These regions will not contribute to the final coadd.
2083  Furthermore, any routine to determine the coadd PSF can now be cognizant
2084  of clipped regions. Note that the algorithm implemented by this task is
2085  preliminary and works correctly for HSC data. Parameter modifications and
2086  or considerable redesigning of the algorithm is likley required for other
2087  surveys.
2088 
2089  ``CompareWarpAssembleCoaddTask`` sub-classes
2090  ``AssembleCoaddTask`` and instantiates ``AssembleCoaddTask``
2091  as a subtask to generate the TemplateCoadd (the model of the static sky).
2092 
2093  Notes
2094  -----
2095  The `lsst.pipe.base.cmdLineTask.CmdLineTask` interface supports a
2096  flag ``-d`` to import ``debug.py`` from your ``PYTHONPATH``; see
2097  ``baseDebug`` for more about ``debug.py`` files.
2098 
2099  This task supports the following debug variables:
2100 
2101  - ``saveCountIm``
2102  If True then save the Epoch Count Image as a fits file in the `figPath`
2103  - ``figPath``
2104  Path to save the debug fits images and figures
2105 
2106  For example, put something like:
2107 
2108  .. code-block:: python
2109 
2110  import lsstDebug
2111  def DebugInfo(name):
2112  di = lsstDebug.getInfo(name)
2113  if name == "lsst.pipe.tasks.assembleCoadd":
2114  di.saveCountIm = True
2115  di.figPath = "/desired/path/to/debugging/output/images"
2116  return di
2117  lsstDebug.Info = DebugInfo
2118 
2119  into your ``debug.py`` file and run ``assemebleCoadd.py`` with the
2120  ``--debug`` flag. Some subtasks may have their own debug variables;
2121  see individual Task documentation.
2122 
2123  Examples
2124  --------
2125  ``CompareWarpAssembleCoaddTask`` assembles a set of warped images into a
2126  coadded image. The ``CompareWarpAssembleCoaddTask`` is invoked by running
2127  ``assembleCoadd.py`` with the flag ``--compareWarpCoadd``.
2128  Usage of ``assembleCoadd.py`` expects a data reference to the tract patch
2129  and filter to be coadded (specified using
2130  '--id = [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]')
2131  along with a list of coaddTempExps to attempt to coadd (specified using
2132  '--selectId [KEY=VALUE1[^VALUE2[^VALUE3...] [KEY=VALUE1[^VALUE2[^VALUE3...] ...]]').
2133  Only the warps that cover the specified tract and patch will be coadded.
2134  A list of the available optional arguments can be obtained by calling
2135  ``assembleCoadd.py`` with the ``--help`` command line argument:
2136 
2137  .. code-block:: none
2138 
2139  assembleCoadd.py --help
2140 
2141  To demonstrate usage of the ``CompareWarpAssembleCoaddTask`` in the larger
2142  context of multi-band processing, we will generate the HSC-I & -R band
2143  oadds from HSC engineering test data provided in the ``ci_hsc`` package.
2144  To begin, assuming that the lsst stack has been already set up, we must
2145  set up the ``obs_subaru`` and ``ci_hsc`` packages.
2146  This defines the environment variable ``$CI_HSC_DIR`` and points at the
2147  location of the package. The raw HSC data live in the ``$CI_HSC_DIR/raw``
2148  directory. To begin assembling the coadds, we must first
2149 
2150  - processCcd
2151  process the individual ccds in $CI_HSC_RAW to produce calibrated exposures
2152  - makeSkyMap
2153  create a skymap that covers the area of the sky present in the raw exposures
2154  - makeCoaddTempExp
2155  warp the individual calibrated exposures to the tangent plane of the coadd
2156 
2157  We can perform all of these steps by running
2158 
2159  .. code-block:: none
2160 
2161  $CI_HSC_DIR scons warp-903986 warp-904014 warp-903990 warp-904010 warp-903988
2162 
2163  This will produce warped ``coaddTempExps`` for each visit. To coadd the
2164  warped data, we call ``assembleCoadd.py`` as follows:
2165 
2166  .. code-block:: none
2167 
2168  assembleCoadd.py --compareWarpCoadd $CI_HSC_DIR/DATA --id patch=5,4 tract=0 filter=HSC-I \
2169  --selectId visit=903986 ccd=16 --selectId visit=903986 ccd=22 --selectId visit=903986 ccd=23 \
2170  --selectId visit=903986 ccd=100 --selectId visit=904014 ccd=1 --selectId visit=904014 ccd=6 \
2171  --selectId visit=904014 ccd=12 --selectId visit=903990 ccd=18 --selectId visit=903990 ccd=25 \
2172  --selectId visit=904010 ccd=4 --selectId visit=904010 ccd=10 --selectId visit=904010 ccd=100 \
2173  --selectId visit=903988 ccd=16 --selectId visit=903988 ccd=17 --selectId visit=903988 ccd=23 \
2174  --selectId visit=903988 ccd=24
2175 
2176  This will process the HSC-I band data. The results are written in
2177  ``$CI_HSC_DIR/DATA/deepCoadd-results/HSC-I``.
2178  """
2179  ConfigClass = CompareWarpAssembleCoaddConfig
2180  _DefaultName = "compareWarpAssembleCoadd"
2181 
2182  def __init__(self, *args, **kwargs):
2183  AssembleCoaddTask.__init__(self, *args, **kwargs)
2184  self.makeSubtask("assembleStaticSkyModel")
2185  detectionSchema = afwTable.SourceTable.makeMinimalSchema()
2186  self.makeSubtask("detect", schema=detectionSchema)
2187  if self.config.doPreserveContainedBySource:
2188  self.makeSubtask("detectTemplate", schema=afwTable.SourceTable.makeMinimalSchema())
2189  if self.config.doScaleWarpVariance:
2190  self.makeSubtask("scaleWarpVariance")
2191  if self.config.doFilterMorphological:
2192  self.makeSubtask("maskStreaks")
2193 
2194  @utils.inheritDoc(AssembleCoaddTask)
2195  def makeSupplementaryDataGen3(self, butlerQC, inputRefs, outputRefs):
2196  """
2197  Generate a templateCoadd to use as a naive model of static sky to
2198  subtract from PSF-Matched warps.
2199 
2200  Returns
2201  -------
2202  result : `lsst.pipe.base.Struct`
2203  Result struct with components:
2204 
2205  - ``templateCoadd`` : coadded exposure (``lsst.afw.image.Exposure``)
2206  - ``nImage`` : N Image (``lsst.afw.image.Image``)
2207  """
2208  # Ensure that psfMatchedWarps are used as input warps for template generation
2209  staticSkyModelInputRefs = copy.deepcopy(inputRefs)
2210  staticSkyModelInputRefs.inputWarps = inputRefs.psfMatchedWarps
2211 
2212  # Because subtasks don't have connections we have to make one.
2213  # The main task's `templateCoadd` is the subtask's `coaddExposure`
2214  staticSkyModelOutputRefs = copy.deepcopy(outputRefs)
2215  if self.config.assembleStaticSkyModel.doWrite:
2216  staticSkyModelOutputRefs.coaddExposure = staticSkyModelOutputRefs.templateCoadd
2217  # Remove template coadd from both subtask's and main tasks outputs,
2218  # because it is handled by the subtask as `coaddExposure`
2219  del outputRefs.templateCoadd
2220  del staticSkyModelOutputRefs.templateCoadd
2221 
2222  # A PSF-Matched nImage does not exist as a dataset type
2223  if 'nImage' in staticSkyModelOutputRefs.keys():
2224  del staticSkyModelOutputRefs.nImage
2225 
2226  templateCoadd = self.assembleStaticSkyModel.runQuantum(butlerQC, staticSkyModelInputRefs,
2227  staticSkyModelOutputRefs)
2228  if templateCoadd is None:
2229  raise RuntimeError(self._noTemplateMessage_noTemplateMessage(self.assembleStaticSkyModel.warpType))
2230 
2231  return pipeBase.Struct(templateCoadd=templateCoadd.coaddExposure,
2232  nImage=templateCoadd.nImage,
2233  warpRefList=templateCoadd.warpRefList,
2234  imageScalerList=templateCoadd.imageScalerList,
2235  weightList=templateCoadd.weightList)
2236 
2237  @utils.inheritDoc(AssembleCoaddTask)
2238  def makeSupplementaryData(self, dataRef, selectDataList=None, warpRefList=None):
2239  """
2240  Generate a templateCoadd to use as a naive model of static sky to
2241  subtract from PSF-Matched warps.
2242 
2243  Returns
2244  -------
2245  result : `lsst.pipe.base.Struct`
2246  Result struct with components:
2247 
2248  - ``templateCoadd``: coadded exposure (``lsst.afw.image.Exposure``)
2249  - ``nImage``: N Image (``lsst.afw.image.Image``)
2250  """
2251  templateCoadd = self.assembleStaticSkyModel.runDataRef(dataRef, selectDataList, warpRefList)
2252  if templateCoadd is None:
2253  raise RuntimeError(self._noTemplateMessage_noTemplateMessage(self.assembleStaticSkyModel.warpType))
2254 
2255  return pipeBase.Struct(templateCoadd=templateCoadd.coaddExposure,
2256  nImage=templateCoadd.nImage,
2257  warpRefList=templateCoadd.warpRefList,
2258  imageScalerList=templateCoadd.imageScalerList,
2259  weightList=templateCoadd.weightList)
2260 
2261  def _noTemplateMessage(self, warpType):
2262  warpName = (warpType[0].upper() + warpType[1:])
2263  message = """No %(warpName)s warps were found to build the template coadd which is
2264  required to run CompareWarpAssembleCoaddTask. To continue assembling this type of coadd,
2265  first either rerun makeCoaddTempExp with config.make%(warpName)s=True or
2266  coaddDriver with config.makeCoadTempExp.make%(warpName)s=True, before assembleCoadd.
2267 
2268  Alternatively, to use another algorithm with existing warps, retarget the CoaddDriverConfig to
2269  another algorithm like:
2270 
2271  from lsst.pipe.tasks.assembleCoadd import SafeClipAssembleCoaddTask
2272  config.assemble.retarget(SafeClipAssembleCoaddTask)
2273  """ % {"warpName": warpName}
2274  return message
2275 
2276  @utils.inheritDoc(AssembleCoaddTask)
2277  @pipeBase.timeMethod
2278  def run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
2279  supplementaryData, *args, **kwargs):
2280  """Assemble the coadd.
2281 
2282  Find artifacts and apply them to the warps' masks creating a list of
2283  alternative masks with a new "CLIPPED" plane and updated "NO_DATA"
2284  plane. Then pass these alternative masks to the base class's `run`
2285  method.
2286 
2287  The input parameters ``supplementaryData`` is a `lsst.pipe.base.Struct`
2288  that must contain a ``templateCoadd`` that serves as the
2289  model of the static sky.
2290  """
2291 
2292  # Check and match the order of the supplementaryData
2293  # (PSF-matched) inputs to the order of the direct inputs,
2294  # so that the artifact mask is applied to the right warp
2295  dataIds = [ref.dataId for ref in tempExpRefList]
2296  psfMatchedDataIds = [ref.dataId for ref in supplementaryData.warpRefList]
2297 
2298  if dataIds != psfMatchedDataIds:
2299  self.log.info("Reordering and or/padding PSF-matched visit input list")
2300  supplementaryData.warpRefList = reorderAndPadList(supplementaryData.warpRefList,
2301  psfMatchedDataIds, dataIds)
2302  supplementaryData.imageScalerList = reorderAndPadList(supplementaryData.imageScalerList,
2303  psfMatchedDataIds, dataIds)
2304 
2305  # Use PSF-Matched Warps (and corresponding scalers) and coadd to find artifacts
2306  spanSetMaskList = self.findArtifactsfindArtifacts(supplementaryData.templateCoadd,
2307  supplementaryData.warpRefList,
2308  supplementaryData.imageScalerList)
2309 
2310  badMaskPlanes = self.config.badMaskPlanes[:]
2311  badMaskPlanes.append("CLIPPED")
2312  badPixelMask = afwImage.Mask.getPlaneBitMask(badMaskPlanes)
2313 
2314  result = AssembleCoaddTask.run(self, skyInfo, tempExpRefList, imageScalerList, weightList,
2315  spanSetMaskList, mask=badPixelMask)
2316 
2317  # Propagate PSF-matched EDGE pixels to coadd SENSOR_EDGE and INEXACT_PSF
2318  # Psf-Matching moves the real edge inwards
2319  self.applyAltEdgeMaskapplyAltEdgeMask(result.coaddExposure.maskedImage.mask, spanSetMaskList)
2320  return result
2321 
2322  def applyAltEdgeMask(self, mask, altMaskList):
2323  """Propagate alt EDGE mask to SENSOR_EDGE AND INEXACT_PSF planes.
2324 
2325  Parameters
2326  ----------
2327  mask : `lsst.afw.image.Mask`
2328  Original mask.
2329  altMaskList : `list`
2330  List of Dicts containing ``spanSet`` lists.
2331  Each element contains the new mask plane name (e.g. "CLIPPED
2332  and/or "NO_DATA") as the key, and list of ``SpanSets`` to apply to
2333  the mask.
2334  """
2335  maskValue = mask.getPlaneBitMask(["SENSOR_EDGE", "INEXACT_PSF"])
2336  for visitMask in altMaskList:
2337  if "EDGE" in visitMask:
2338  for spanSet in visitMask['EDGE']:
2339  spanSet.clippedTo(mask.getBBox()).setMask(mask, maskValue)
2340 
2341  def findArtifacts(self, templateCoadd, tempExpRefList, imageScalerList):
2342  """Find artifacts.
2343 
2344  Loop through warps twice. The first loop builds a map with the count
2345  of how many epochs each pixel deviates from the templateCoadd by more
2346  than ``config.chiThreshold`` sigma. The second loop takes each
2347  difference image and filters the artifacts detected in each using
2348  count map to filter out variable sources and sources that are
2349  difficult to subtract cleanly.
2350 
2351  Parameters
2352  ----------
2353  templateCoadd : `lsst.afw.image.Exposure`
2354  Exposure to serve as model of static sky.
2355  tempExpRefList : `list`
2356  List of data references to warps.
2357  imageScalerList : `list`
2358  List of image scalers.
2359 
2360  Returns
2361  -------
2362  altMasks : `list`
2363  List of dicts containing information about CLIPPED
2364  (i.e., artifacts), NO_DATA, and EDGE pixels.
2365  """
2366 
2367  self.log.debug("Generating Count Image, and mask lists.")
2368  coaddBBox = templateCoadd.getBBox()
2369  slateIm = afwImage.ImageU(coaddBBox)
2370  epochCountImage = afwImage.ImageU(coaddBBox)
2371  nImage = afwImage.ImageU(coaddBBox)
2372  spanSetArtifactList = []
2373  spanSetNoDataMaskList = []
2374  spanSetEdgeList = []
2375  spanSetBadMorphoList = []
2376  badPixelMask = self.getBadPixelMask()
2377 
2378  # mask of the warp diffs should = that of only the warp
2379  templateCoadd.mask.clearAllMaskPlanes()
2380 
2381  if self.config.doPreserveContainedBySource:
2382  templateFootprints = self.detectTemplate.detectFootprints(templateCoadd)
2383  else:
2384  templateFootprints = None
2385 
2386  for warpRef, imageScaler in zip(tempExpRefList, imageScalerList):
2387  warpDiffExp = self._readAndComputeWarpDiff_readAndComputeWarpDiff(warpRef, imageScaler, templateCoadd)
2388  if warpDiffExp is not None:
2389  # This nImage only approximates the final nImage because it uses the PSF-matched mask
2390  nImage.array += (numpy.isfinite(warpDiffExp.image.array)
2391  * ((warpDiffExp.mask.array & badPixelMask) == 0)).astype(numpy.uint16)
2392  fpSet = self.detect.detectFootprints(warpDiffExp, doSmooth=False, clearMask=True)
2393  fpSet.positive.merge(fpSet.negative)
2394  footprints = fpSet.positive
2395  slateIm.set(0)
2396  spanSetList = [footprint.spans for footprint in footprints.getFootprints()]
2397 
2398  # Remove artifacts due to defects before they contribute to the epochCountImage
2399  if self.config.doPrefilterArtifacts:
2400  spanSetList = self.prefilterArtifactsprefilterArtifacts(spanSetList, warpDiffExp)
2401 
2402  # Clear mask before adding prefiltered spanSets
2403  self.detect.clearMask(warpDiffExp.mask)
2404  for spans in spanSetList:
2405  spans.setImage(slateIm, 1, doClip=True)
2406  spans.setMask(warpDiffExp.mask, warpDiffExp.mask.getPlaneBitMask("DETECTED"))
2407  epochCountImage += slateIm
2408 
2409  if self.config.doFilterMorphological:
2410  maskName = self.config.streakMaskName
2411  _ = self.maskStreaks.run(warpDiffExp)
2412  streakMask = warpDiffExp.mask
2413  spanSetStreak = afwGeom.SpanSet.fromMask(streakMask,
2414  streakMask.getPlaneBitMask(maskName)).split()
2415 
2416  # PSF-Matched warps have less available area (~the matching kernel) because the calexps
2417  # undergo a second convolution. Pixels with data in the direct warp
2418  # but not in the PSF-matched warp will not have their artifacts detected.
2419  # NaNs from the PSF-matched warp therefore must be masked in the direct warp
2420  nans = numpy.where(numpy.isnan(warpDiffExp.maskedImage.image.array), 1, 0)
2421  nansMask = afwImage.makeMaskFromArray(nans.astype(afwImage.MaskPixel))
2422  nansMask.setXY0(warpDiffExp.getXY0())
2423  edgeMask = warpDiffExp.mask
2424  spanSetEdgeMask = afwGeom.SpanSet.fromMask(edgeMask,
2425  edgeMask.getPlaneBitMask("EDGE")).split()
2426  else:
2427  # If the directWarp has <1% coverage, the psfMatchedWarp can have 0% and not exist
2428  # In this case, mask the whole epoch
2429  nansMask = afwImage.MaskX(coaddBBox, 1)
2430  spanSetList = []
2431  spanSetEdgeMask = []
2432  spanSetStreak = []
2433 
2434  spanSetNoDataMask = afwGeom.SpanSet.fromMask(nansMask).split()
2435 
2436  spanSetNoDataMaskList.append(spanSetNoDataMask)
2437  spanSetArtifactList.append(spanSetList)
2438  spanSetEdgeList.append(spanSetEdgeMask)
2439  if self.config.doFilterMorphological:
2440  spanSetBadMorphoList.append(spanSetStreak)
2441 
2442  if lsstDebug.Info(__name__).saveCountIm:
2443  path = self._dataRef2DebugPath_dataRef2DebugPath("epochCountIm", tempExpRefList[0], coaddLevel=True)
2444  epochCountImage.writeFits(path)
2445 
2446  for i, spanSetList in enumerate(spanSetArtifactList):
2447  if spanSetList:
2448  filteredSpanSetList = self.filterArtifactsfilterArtifacts(spanSetList, epochCountImage, nImage,
2449  templateFootprints)
2450  spanSetArtifactList[i] = filteredSpanSetList
2451  if self.config.doFilterMorphological:
2452  spanSetArtifactList[i] += spanSetBadMorphoList[i]
2453 
2454  altMasks = []
2455  for artifacts, noData, edge in zip(spanSetArtifactList, spanSetNoDataMaskList, spanSetEdgeList):
2456  altMasks.append({'CLIPPED': artifacts,
2457  'NO_DATA': noData,
2458  'EDGE': edge})
2459  return altMasks
2460 
2461  def prefilterArtifacts(self, spanSetList, exp):
2462  """Remove artifact candidates covered by bad mask plane.
2463 
2464  Any future editing of the candidate list that does not depend on
2465  temporal information should go in this method.
2466 
2467  Parameters
2468  ----------
2469  spanSetList : `list`
2470  List of SpanSets representing artifact candidates.
2471  exp : `lsst.afw.image.Exposure`
2472  Exposure containing mask planes used to prefilter.
2473 
2474  Returns
2475  -------
2476  returnSpanSetList : `list`
2477  List of SpanSets with artifacts.
2478  """
2479  badPixelMask = exp.mask.getPlaneBitMask(self.config.prefilterArtifactsMaskPlanes)
2480  goodArr = (exp.mask.array & badPixelMask) == 0
2481  returnSpanSetList = []
2482  bbox = exp.getBBox()
2483  x0, y0 = exp.getXY0()
2484  for i, span in enumerate(spanSetList):
2485  y, x = span.clippedTo(bbox).indices()
2486  yIndexLocal = numpy.array(y) - y0
2487  xIndexLocal = numpy.array(x) - x0
2488  goodRatio = numpy.count_nonzero(goodArr[yIndexLocal, xIndexLocal])/span.getArea()
2489  if goodRatio > self.config.prefilterArtifactsRatio:
2490  returnSpanSetList.append(span)
2491  return returnSpanSetList
2492 
2493  def filterArtifacts(self, spanSetList, epochCountImage, nImage, footprintsToExclude=None):
2494  """Filter artifact candidates.
2495 
2496  Parameters
2497  ----------
2498  spanSetList : `list`
2499  List of SpanSets representing artifact candidates.
2500  epochCountImage : `lsst.afw.image.Image`
2501  Image of accumulated number of warpDiff detections.
2502  nImage : `lsst.afw.image.Image`
2503  Image of the accumulated number of total epochs contributing.
2504 
2505  Returns
2506  -------
2507  maskSpanSetList : `list`
2508  List of SpanSets with artifacts.
2509  """
2510 
2511  maskSpanSetList = []
2512  x0, y0 = epochCountImage.getXY0()
2513  for i, span in enumerate(spanSetList):
2514  y, x = span.indices()
2515  yIdxLocal = [y1 - y0 for y1 in y]
2516  xIdxLocal = [x1 - x0 for x1 in x]
2517  outlierN = epochCountImage.array[yIdxLocal, xIdxLocal]
2518  totalN = nImage.array[yIdxLocal, xIdxLocal]
2519 
2520  # effectiveMaxNumEpochs is broken line (fraction of N) with characteristic config.maxNumEpochs
2521  effMaxNumEpochsHighN = (self.config.maxNumEpochs
2522  + self.config.maxFractionEpochsHigh*numpy.mean(totalN))
2523  effMaxNumEpochsLowN = self.config.maxFractionEpochsLow * numpy.mean(totalN)
2524  effectiveMaxNumEpochs = int(min(effMaxNumEpochsLowN, effMaxNumEpochsHighN))
2525  nPixelsBelowThreshold = numpy.count_nonzero((outlierN > 0)
2526  & (outlierN <= effectiveMaxNumEpochs))
2527  percentBelowThreshold = nPixelsBelowThreshold / len(outlierN)
2528  if percentBelowThreshold > self.config.spatialThreshold:
2529  maskSpanSetList.append(span)
2530 
2531  if self.config.doPreserveContainedBySource and footprintsToExclude is not None:
2532  # If a candidate is contained by a footprint on the template coadd, do not clip
2533  filteredMaskSpanSetList = []
2534  for span in maskSpanSetList:
2535  doKeep = True
2536  for footprint in footprintsToExclude.positive.getFootprints():
2537  if footprint.spans.contains(span):
2538  doKeep = False
2539  break
2540  if doKeep:
2541  filteredMaskSpanSetList.append(span)
2542  maskSpanSetList = filteredMaskSpanSetList
2543 
2544  return maskSpanSetList
2545 
2546  def _readAndComputeWarpDiff(self, warpRef, imageScaler, templateCoadd):
2547  """Fetch a warp from the butler and return a warpDiff.
2548 
2549  Parameters
2550  ----------
2551  warpRef : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
2552  Butler dataRef for the warp.
2553  imageScaler : `lsst.pipe.tasks.scaleZeroPoint.ImageScaler`
2554  An image scaler object.
2555  templateCoadd : `lsst.afw.image.Exposure`
2556  Exposure to be substracted from the scaled warp.
2557 
2558  Returns
2559  -------
2560  warp : `lsst.afw.image.Exposure`
2561  Exposure of the image difference between the warp and template.
2562  """
2563 
2564  # If the PSF-Matched warp did not exist for this direct warp
2565  # None is holding its place to maintain order in Gen 3
2566  if warpRef is None:
2567  return None
2568  # Warp comparison must use PSF-Matched Warps regardless of requested coadd warp type
2569  warpName = self.getTempExpDatasetName('psfMatched')
2570  if not isinstance(warpRef, DeferredDatasetHandle):
2571  if not warpRef.datasetExists(warpName):
2572  self.log.warning("Could not find %s %s; skipping it", warpName, warpRef.dataId)
2573  return None
2574  warp = warpRef.get(datasetType=warpName, immediate=True)
2575  # direct image scaler OK for PSF-matched Warp
2576  imageScaler.scaleMaskedImage(warp.getMaskedImage())
2577  mi = warp.getMaskedImage()
2578  if self.config.doScaleWarpVariance:
2579  try:
2580  self.scaleWarpVariance.run(mi)
2581  except Exception as exc:
2582  self.log.warning("Unable to rescale variance of warp (%s); leaving it as-is", exc)
2583  mi -= templateCoadd.getMaskedImage()
2584  return warp
2585 
2586  def _dataRef2DebugPath(self, prefix, warpRef, coaddLevel=False):
2587  """Return a path to which to write debugging output.
2588 
2589  Creates a hyphen-delimited string of dataId values for simple filenames.
2590 
2591  Parameters
2592  ----------
2593  prefix : `str`
2594  Prefix for filename.
2595  warpRef : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
2596  Butler dataRef to make the path from.
2597  coaddLevel : `bool`, optional.
2598  If True, include only coadd-level keys (e.g., 'tract', 'patch',
2599  'filter', but no 'visit').
2600 
2601  Returns
2602  -------
2603  result : `str`
2604  Path for debugging output.
2605  """
2606  if coaddLevel:
2607  keys = warpRef.getButler().getKeys(self.getCoaddDatasetName(self.warpType))
2608  else:
2609  keys = warpRef.dataId.keys()
2610  keyList = sorted(keys, reverse=True)
2611  directory = lsstDebug.Info(__name__).figPath if lsstDebug.Info(__name__).figPath else "."
2612  filename = "%s-%s.fits" % (prefix, '-'.join([str(warpRef.dataId[k]) for k in keyList]))
2613  return os.path.join(directory, filename)
std::vector< SchemaItem< Flag > > * items
int min
int max
A set of Footprints, associated with a MaskedImage.
Definition: FootprintSet.h:53
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 group of labels for a filter in an exposure or coadd.
Definition: FilterLabel.h:58
Represent a 2-dimensional array of bitmask pixels.
Definition: Mask.h:77
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.
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
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.
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:346
def makeSkyInfo(skyMap, tractId, patchId)
Definition: coaddBase.py:289
def getGroupDataRef(butler, datasetType, groupTuple, keys)
Definition: coaddHelpers.py:99
def groupPatchExposures(patchDataRef, calexpDataRefList, coaddDatasetType="deepCoadd", tempExpDatasetType="deepCoadd_directWarp")
Definition: coaddHelpers.py:60