lsst.pipe.tasks  21.0.0-80-g62ad60b1+450980b591
imageDifference.py
Go to the documentation of this file.
1 # This file is part of pipe_tasks.
2 #
3 # Developed for the LSST Data Management System.
4 # This product includes software developed by the LSST Project
5 # (https://www.lsst.org).
6 # See the COPYRIGHT file at the top-level directory of this distribution
7 # for details of code ownership.
8 #
9 # This program is free software: you can redistribute it and/or modify
10 # it under the terms of the GNU General Public License as published by
11 # the Free Software Foundation, either version 3 of the License, or
12 # (at your option) any later version.
13 #
14 # This program is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 # GNU General Public License for more details.
18 #
19 # You should have received a copy of the GNU General Public License
20 # along with this program. If not, see <https://www.gnu.org/licenses/>.
21 
22 import math
23 import random
24 import numpy
25 
26 import lsst.utils
27 import lsst.pex.config as pexConfig
28 import lsst.pipe.base as pipeBase
29 import lsst.daf.base as dafBase
30 import lsst.geom as geom
31 import lsst.afw.math as afwMath
32 import lsst.afw.table as afwTable
33 from lsst.meas.astrom import AstrometryConfig, AstrometryTask
34 from lsst.meas.base import ForcedMeasurementTask, ApplyApCorrTask
35 from lsst.meas.algorithms import LoadIndexedReferenceObjectsTask
36 from lsst.pipe.tasks.registerImage import RegisterTask
37 from lsst.pipe.tasks.scaleVariance import ScaleVarianceTask
38 from lsst.meas.algorithms import SourceDetectionTask, SingleGaussianPsf, ObjectSizeStarSelectorTask
39 from lsst.ip.diffim import (DipoleAnalysis, SourceFlagChecker, KernelCandidateF, makeKernelBasisList,
40  KernelCandidateQa, DiaCatalogSourceSelectorTask, DiaCatalogSourceSelectorConfig,
41  GetCoaddAsTemplateTask, GetCalexpAsTemplateTask, DipoleFitTask,
42  DecorrelateALKernelSpatialTask, subtractAlgorithmRegistry)
43 import lsst.ip.diffim.diffimTools as diffimTools
44 import lsst.ip.diffim.utils as diUtils
45 import lsst.afw.display as afwDisplay
46 from lsst.skymap import BaseSkyMap
47 from lsst.obs.base import ExposureIdInfo
48 
49 __all__ = ["ImageDifferenceConfig", "ImageDifferenceTask"]
50 FwhmPerSigma = 2*math.sqrt(2*math.log(2))
51 IqrToSigma = 0.741
52 
53 
54 class ImageDifferenceTaskConnections(pipeBase.PipelineTaskConnections,
55  dimensions=("instrument", "visit", "detector", "skymap"),
56  defaultTemplates={"coaddName": "deep",
57  "skyMapName": "deep",
58  "warpTypeSuffix": "",
59  "fakesType": ""}):
60 
61  exposure = pipeBase.connectionTypes.Input(
62  doc="Input science exposure to subtract from.",
63  dimensions=("instrument", "visit", "detector"),
64  storageClass="ExposureF",
65  name="{fakesType}calexp"
66  )
67 
68  # TODO DM-22953
69  # kernelSources = pipeBase.connectionTypes.Input(
70  # doc="Source catalog produced in calibrate task for kernel candidate sources",
71  # name="src",
72  # storageClass="SourceCatalog",
73  # dimensions=("instrument", "visit", "detector"),
74  # )
75 
76  skyMap = pipeBase.connectionTypes.Input(
77  doc="Input definition of geometry/bbox and projection/wcs for template exposures",
78  name=BaseSkyMap.SKYMAP_DATASET_TYPE_NAME,
79  dimensions=("skymap", ),
80  storageClass="SkyMap",
81  )
82  coaddExposures = pipeBase.connectionTypes.Input(
83  doc="Input template to match and subtract from the exposure",
84  dimensions=("tract", "patch", "skymap", "band"),
85  storageClass="ExposureF",
86  name="{fakesType}{coaddName}Coadd{warpTypeSuffix}",
87  multiple=True,
88  deferLoad=True
89  )
90  dcrCoadds = pipeBase.connectionTypes.Input(
91  doc="Input DCR template to match and subtract from the exposure",
92  name="{fakesType}dcrCoadd{warpTypeSuffix}",
93  storageClass="ExposureF",
94  dimensions=("tract", "patch", "skymap", "band", "subfilter"),
95  multiple=True,
96  deferLoad=True
97  )
98  outputSchema = pipeBase.connectionTypes.InitOutput(
99  doc="Schema (as an example catalog) for output DIASource catalog.",
100  storageClass="SourceCatalog",
101  name="{fakesType}{coaddName}Diff_diaSrc_schema",
102  )
103  subtractedExposure = pipeBase.connectionTypes.Output(
104  doc="Output difference image",
105  dimensions=("instrument", "visit", "detector"),
106  storageClass="ExposureF",
107  name="{fakesType}{coaddName}Diff_differenceExp",
108  )
109  warpedExposure = pipeBase.connectionTypes.Output(
110  doc="Warped template used to create `subtractedExposure`.",
111  dimensions=("instrument", "visit", "detector"),
112  storageClass="ExposureF",
113  name="{fakesType}{coaddName}Diff_warpedExp",
114  )
115  diaSources = pipeBase.connectionTypes.Output(
116  doc="Output detected diaSources on the difference image",
117  dimensions=("instrument", "visit", "detector"),
118  storageClass="SourceCatalog",
119  name="{fakesType}{coaddName}Diff_diaSrc",
120  )
121 
122  def __init__(self, *, config=None):
123  super().__init__(config=config)
124  if config.coaddName == 'dcr':
125  self.inputs.remove("coaddExposures")
126  else:
127  self.inputs.remove("dcrCoadds")
128 
129  # TODO DM-22953: Add support for refObjLoader (kernelSourcesFromRef)
130  # Make kernelSources optional
131 
132 
133 class ImageDifferenceConfig(pipeBase.PipelineTaskConfig,
134  pipelineConnections=ImageDifferenceTaskConnections):
135  """Config for ImageDifferenceTask
136  """
137  doAddCalexpBackground = pexConfig.Field(dtype=bool, default=False,
138  doc="Add background to calexp before processing it. "
139  "Useful as ipDiffim does background matching.")
140  doUseRegister = pexConfig.Field(dtype=bool, default=True,
141  doc="Use image-to-image registration to align template with "
142  "science image")
143  doDebugRegister = pexConfig.Field(dtype=bool, default=False,
144  doc="Writing debugging data for doUseRegister")
145  doSelectSources = pexConfig.Field(dtype=bool, default=True,
146  doc="Select stars to use for kernel fitting")
147  doSelectDcrCatalog = pexConfig.Field(dtype=bool, default=False,
148  doc="Select stars of extreme color as part of the control sample")
149  doSelectVariableCatalog = pexConfig.Field(dtype=bool, default=False,
150  doc="Select stars that are variable to be part "
151  "of the control sample")
152  doSubtract = pexConfig.Field(dtype=bool, default=True, doc="Compute subtracted exposure?")
153  doPreConvolve = pexConfig.Field(dtype=bool, default=True,
154  doc="Convolve science image by its PSF before PSF-matching?")
155  doScaleTemplateVariance = pexConfig.Field(dtype=bool, default=False,
156  doc="Scale variance of the template before PSF matching")
157  doScaleDiffimVariance = pexConfig.Field(dtype=bool, default=True,
158  doc="Scale variance of the diffim before PSF matching. "
159  "You may do either this or template variance scaling, "
160  "or neither. (Doing both is a waste of CPU.)")
161  useGaussianForPreConvolution = pexConfig.Field(dtype=bool, default=True,
162  doc="Use a simple gaussian PSF model for pre-convolution "
163  "(else use fit PSF)? Ignored if doPreConvolve false.")
164  doDetection = pexConfig.Field(dtype=bool, default=True, doc="Detect sources?")
165  doDecorrelation = pexConfig.Field(dtype=bool, default=False,
166  doc="Perform diffim decorrelation to undo pixel correlation due to A&L "
167  "kernel convolution? If True, also update the diffim PSF.")
168  doMerge = pexConfig.Field(dtype=bool, default=True,
169  doc="Merge positive and negative diaSources with grow radius "
170  "set by growFootprint")
171  doMatchSources = pexConfig.Field(dtype=bool, default=True,
172  doc="Match diaSources with input calexp sources and ref catalog sources")
173  doMeasurement = pexConfig.Field(dtype=bool, default=True, doc="Measure diaSources?")
174  doDipoleFitting = pexConfig.Field(dtype=bool, default=True, doc="Measure dipoles using new algorithm?")
175  doForcedMeasurement = pexConfig.Field(
176  dtype=bool,
177  default=True,
178  doc="Force photometer diaSource locations on PVI?")
179  doWriteSubtractedExp = pexConfig.Field(dtype=bool, default=True, doc="Write difference exposure?")
180  doWriteWarpedExp = pexConfig.Field(dtype=bool, default=False,
181  doc="Write WCS, warped template coadd exposure?")
182  doWriteMatchedExp = pexConfig.Field(dtype=bool, default=False,
183  doc="Write warped and PSF-matched template coadd exposure?")
184  doWriteSources = pexConfig.Field(dtype=bool, default=True, doc="Write sources?")
185  doAddMetrics = pexConfig.Field(dtype=bool, default=True,
186  doc="Add columns to the source table to hold analysis metrics?")
187 
188  coaddName = pexConfig.Field(
189  doc="coadd name: typically one of deep, goodSeeing, or dcr",
190  dtype=str,
191  default="deep",
192  )
193  convolveTemplate = pexConfig.Field(
194  doc="Which image gets convolved (default = template)",
195  dtype=bool,
196  default=True
197  )
198  refObjLoader = pexConfig.ConfigurableField(
199  target=LoadIndexedReferenceObjectsTask,
200  doc="reference object loader",
201  )
202  astrometer = pexConfig.ConfigurableField(
203  target=AstrometryTask,
204  doc="astrometry task; used to match sources to reference objects, but not to fit a WCS",
205  )
206  sourceSelector = pexConfig.ConfigurableField(
207  target=ObjectSizeStarSelectorTask,
208  doc="Source selection algorithm",
209  )
210  subtract = subtractAlgorithmRegistry.makeField("Subtraction Algorithm", default="al")
211  decorrelate = pexConfig.ConfigurableField(
212  target=DecorrelateALKernelSpatialTask,
213  doc="Decorrelate effects of A&L kernel convolution on image difference, only if doSubtract is True. "
214  "If this option is enabled, then detection.thresholdValue should be set to 5.0 (rather than the "
215  "default of 5.5).",
216  )
217  # Old style ImageMapper grid. ZogyTask has its own grid option
218  doSpatiallyVarying = pexConfig.Field(
219  dtype=bool,
220  default=False,
221  doc="Perform A&L decorrelation on a grid across the "
222  "image in order to allow for spatial variations. Zogy does not use this option."
223  )
224  detection = pexConfig.ConfigurableField(
225  target=SourceDetectionTask,
226  doc="Low-threshold detection for final measurement",
227  )
228  measurement = pexConfig.ConfigurableField(
229  target=DipoleFitTask,
230  doc="Enable updated dipole fitting method",
231  )
232  doApCorr = lsst.pex.config.Field(
233  dtype=bool,
234  default=True,
235  doc="Run subtask to apply aperture corrections"
236  )
237  applyApCorr = lsst.pex.config.ConfigurableField(
238  target=ApplyApCorrTask,
239  doc="Subtask to apply aperture corrections"
240  )
241  forcedMeasurement = pexConfig.ConfigurableField(
242  target=ForcedMeasurementTask,
243  doc="Subtask to force photometer PVI at diaSource location.",
244  )
245  getTemplate = pexConfig.ConfigurableField(
246  target=GetCoaddAsTemplateTask,
247  doc="Subtask to retrieve template exposure and sources",
248  )
249  scaleVariance = pexConfig.ConfigurableField(
250  target=ScaleVarianceTask,
251  doc="Subtask to rescale the variance of the template "
252  "to the statistically expected level"
253  )
254  controlStepSize = pexConfig.Field(
255  doc="What step size (every Nth one) to select a control sample from the kernelSources",
256  dtype=int,
257  default=5
258  )
259  controlRandomSeed = pexConfig.Field(
260  doc="Random seed for shuffing the control sample",
261  dtype=int,
262  default=10
263  )
264  register = pexConfig.ConfigurableField(
265  target=RegisterTask,
266  doc="Task to enable image-to-image image registration (warping)",
267  )
268  kernelSourcesFromRef = pexConfig.Field(
269  doc="Select sources to measure kernel from reference catalog if True, template if false",
270  dtype=bool,
271  default=False
272  )
273  templateSipOrder = pexConfig.Field(
274  dtype=int, default=2,
275  doc="Sip Order for fitting the Template Wcs (default is too high, overfitting)"
276  )
277  growFootprint = pexConfig.Field(
278  dtype=int, default=2,
279  doc="Grow positive and negative footprints by this amount before merging"
280  )
281  diaSourceMatchRadius = pexConfig.Field(
282  dtype=float, default=0.5,
283  doc="Match radius (in arcseconds) for DiaSource to Source association"
284  )
285 
286  def setDefaults(self):
287  # defaults are OK for catalog and diacatalog
288 
289  self.subtract['al'].kernel.name = "AL"
290  self.subtract['al'].kernel.active.fitForBackground = True
291  self.subtract['al'].kernel.active.spatialKernelOrder = 1
292  self.subtract['al'].kernel.active.spatialBgOrder = 2
293  self.doPreConvolve = False
294  self.doMatchSources = False
295  self.doAddMetrics = False
296  self.doUseRegister = False
297 
298  # DiaSource Detection
299  self.detection.thresholdPolarity = "both"
300  self.detection.thresholdValue = 5.5
301  self.detection.reEstimateBackground = False
302  self.detection.thresholdType = "pixel_stdev"
303 
304  # Add filtered flux measurement, the correct measurement for pre-convolved images.
305  # Enable all measurements, regardless of doPreConvolve, as it makes data harvesting easier.
306  # To change that you must modify algorithms.names in the task's applyOverrides method,
307  # after the user has set doPreConvolve.
308  self.measurement.algorithms.names.add('base_PeakLikelihoodFlux')
309  self.measurement.plugins.names |= ['base_LocalPhotoCalib',
310  'base_LocalWcs']
311 
312  self.forcedMeasurement.plugins = ["base_TransformedCentroid", "base_PsfFlux"]
313  self.forcedMeasurement.copyColumns = {
314  "id": "objectId", "parent": "parentObjectId", "coord_ra": "coord_ra", "coord_dec": "coord_dec"}
315  self.forcedMeasurement.slots.centroid = "base_TransformedCentroid"
316  self.forcedMeasurement.slots.shape = None
317 
318  # For shuffling the control sample
319  random.seed(self.controlRandomSeed)
320 
321  def validate(self):
322  pexConfig.Config.validate(self)
323  if self.doAddMetrics and not self.doSubtract:
324  raise ValueError("Subtraction must be enabled for kernel metrics calculation.")
325  if not self.doSubtract and not self.doDetection:
326  raise ValueError("Either doSubtract or doDetection must be enabled.")
327  if self.subtract.name == 'zogy' and self.doAddMetrics:
328  raise ValueError("Kernel metrics does not exist in zogy subtraction.")
329  if self.doMeasurement and not self.doDetection:
330  raise ValueError("Cannot run source measurement without source detection.")
331  if self.doMerge and not self.doDetection:
332  raise ValueError("Cannot run source merging without source detection.")
333  if self.doUseRegister and not self.doSelectSources:
334  raise ValueError("doUseRegister=True and doSelectSources=False. "
335  "Cannot run RegisterTask without selecting sources.")
336  if self.doPreConvolve and self.doDecorrelation and not self.convolveTemplate:
337  raise ValueError("doPreConvolve=True and doDecorrelation=True and "
338  "convolveTemplate=False is not supported.")
339  if hasattr(self.getTemplate, "coaddName"):
340  if self.getTemplate.coaddName != self.coaddName:
341  raise ValueError("Mis-matched coaddName and getTemplate.coaddName in the config.")
342  if self.doScaleDiffimVariance and self.doScaleTemplateVariance:
343  raise ValueError("Scaling the diffim variance and scaling the template variance "
344  "are both set. Please choose one or the other.")
345 
346 
347 class ImageDifferenceTaskRunner(pipeBase.ButlerInitializedTaskRunner):
348 
349  @staticmethod
350  def getTargetList(parsedCmd, **kwargs):
351  return pipeBase.TaskRunner.getTargetList(parsedCmd, templateIdList=parsedCmd.templateId.idList,
352  **kwargs)
353 
354 
355 class ImageDifferenceTask(pipeBase.CmdLineTask, pipeBase.PipelineTask):
356  """Subtract an image from a template and measure the result
357  """
358  ConfigClass = ImageDifferenceConfig
359  RunnerClass = ImageDifferenceTaskRunner
360  _DefaultName = "imageDifference"
361 
362  def __init__(self, butler=None, **kwargs):
363  """!Construct an ImageDifference Task
364 
365  @param[in] butler Butler object to use in constructing reference object loaders
366  """
367  super().__init__(**kwargs)
368  self.makeSubtask("getTemplate")
369 
370  self.makeSubtask("subtract")
371 
372  if self.config.subtract.name == 'al' and self.config.doDecorrelation:
373  self.makeSubtask("decorrelate")
374 
375  if self.config.doScaleTemplateVariance or self.config.doScaleDiffimVariance:
376  self.makeSubtask("scaleVariance")
377 
378  if self.config.doUseRegister:
379  self.makeSubtask("register")
380  self.schema = afwTable.SourceTable.makeMinimalSchema()
381 
382  if self.config.doSelectSources:
383  self.makeSubtask("sourceSelector")
384  if self.config.kernelSourcesFromRef:
385  self.makeSubtask('refObjLoader', butler=butler)
386  self.makeSubtask("astrometer", refObjLoader=self.refObjLoader)
387 
388  self.algMetadata = dafBase.PropertyList()
389  if self.config.doDetection:
390  self.makeSubtask("detection", schema=self.schema)
391  if self.config.doMeasurement:
392  self.makeSubtask("measurement", schema=self.schema,
393  algMetadata=self.algMetadata)
394  if self.config.doApCorr:
395  self.makeSubtask("applyApCorr", schema=self.measurement.schema)
396  if self.config.doForcedMeasurement:
397  self.schema.addField(
398  "ip_diffim_forced_PsfFlux_instFlux", "D",
399  "Forced PSF flux measured on the direct image.",
400  units="count")
401  self.schema.addField(
402  "ip_diffim_forced_PsfFlux_instFluxErr", "D",
403  "Forced PSF flux error measured on the direct image.",
404  units="count")
405  self.schema.addField(
406  "ip_diffim_forced_PsfFlux_area", "F",
407  "Forced PSF flux effective area of PSF.",
408  units="pixel")
409  self.schema.addField(
410  "ip_diffim_forced_PsfFlux_flag", "Flag",
411  "Forced PSF flux general failure flag.")
412  self.schema.addField(
413  "ip_diffim_forced_PsfFlux_flag_noGoodPixels", "Flag",
414  "Forced PSF flux not enough non-rejected pixels in data to attempt the fit.")
415  self.schema.addField(
416  "ip_diffim_forced_PsfFlux_flag_edge", "Flag",
417  "Forced PSF flux object was too close to the edge of the image to use the full PSF model.")
418  self.makeSubtask("forcedMeasurement", refSchema=self.schema)
419  if self.config.doMatchSources:
420  self.schema.addField("refMatchId", "L", "unique id of reference catalog match")
421  self.schema.addField("srcMatchId", "L", "unique id of source match")
422 
423  # initialize InitOutputs
424  self.outputSchema = afwTable.SourceCatalog(self.schema)
425  self.outputSchema.getTable().setMetadata(self.algMetadata)
426 
427  @staticmethod
428  def makeIdFactory(expId, expBits):
429  """Create IdFactory instance for unique 64 bit diaSource id-s.
430 
431  Parameters
432  ----------
433  expId : `int`
434  Exposure id.
435 
436  expBits: `int`
437  Number of used bits in ``expId``.
438 
439  Note
440  ----
441  The diasource id-s consists of the ``expId`` stored fixed in the highest value
442  ``expBits`` of the 64-bit integer plus (bitwise or) a generated sequence number in the
443  low value end of the integer.
444 
445  Returns
446  -------
447  idFactory: `lsst.afw.table.IdFactory`
448  """
449  return ExposureIdInfo(expId, expBits).makeSourceIdFactory()
450 
451  @lsst.utils.inheritDoc(pipeBase.PipelineTask)
452  def runQuantum(self, butlerQC: pipeBase.ButlerQuantumContext,
453  inputRefs: pipeBase.InputQuantizedConnection,
454  outputRefs: pipeBase.OutputQuantizedConnection):
455  inputs = butlerQC.get(inputRefs)
456  self.log.info(f"Processing {butlerQC.quantum.dataId}")
457  expId, expBits = butlerQC.quantum.dataId.pack("visit_detector",
458  returnMaxBits=True)
459  idFactory = self.makeIdFactory(expId=expId, expBits=expBits)
460  if self.config.coaddName == 'dcr':
461  templateExposures = inputRefs.dcrCoadds
462  else:
463  templateExposures = inputRefs.coaddExposures
464  templateStruct = self.getTemplate.runQuantum(
465  inputs['exposure'], butlerQC, inputRefs.skyMap, templateExposures
466  )
467 
468  outputs = self.run(exposure=inputs['exposure'],
469  templateExposure=templateStruct.exposure,
470  idFactory=idFactory)
471  butlerQC.put(outputs, outputRefs)
472 
473  @pipeBase.timeMethod
474  def runDataRef(self, sensorRef, templateIdList=None):
475  """Subtract an image from a template coadd and measure the result.
476 
477  Data I/O wrapper around `run` using the butler in Gen2.
478 
479  Parameters
480  ----------
481  sensorRef : `lsst.daf.persistence.ButlerDataRef`
482  Sensor-level butler data reference, used for the following data products:
483 
484  Input only:
485  - calexp
486  - psf
487  - ccdExposureId
488  - ccdExposureId_bits
489  - self.config.coaddName + "Coadd_skyMap"
490  - self.config.coaddName + "Coadd"
491  Input or output, depending on config:
492  - self.config.coaddName + "Diff_subtractedExp"
493  Output, depending on config:
494  - self.config.coaddName + "Diff_matchedExp"
495  - self.config.coaddName + "Diff_src"
496 
497  Returns
498  -------
499  results : `lsst.pipe.base.Struct`
500  Returns the Struct by `run`.
501  """
502  subtractedExposureName = self.config.coaddName + "Diff_differenceExp"
503  subtractedExposure = None
504  selectSources = None
505  calexpBackgroundExposure = None
506  self.log.info("Processing %s" % (sensorRef.dataId))
507 
508  # We make one IdFactory that will be used by both icSrc and src datasets;
509  # I don't know if this is the way we ultimately want to do things, but at least
510  # this ensures the source IDs are fully unique.
511  idFactory = self.makeIdFactory(expId=int(sensorRef.get("ccdExposureId")),
512  expBits=sensorRef.get("ccdExposureId_bits"))
513  if self.config.doAddCalexpBackground:
514  calexpBackgroundExposure = sensorRef.get("calexpBackground")
515 
516  # Retrieve the science image we wish to analyze
517  exposure = sensorRef.get("calexp", immediate=True)
518 
519  # Retrieve the template image
520  template = self.getTemplate.runDataRef(exposure, sensorRef, templateIdList=templateIdList)
521 
522  if sensorRef.datasetExists("src"):
523  self.log.info("Source selection via src product")
524  # Sources already exist; for data release processing
525  selectSources = sensorRef.get("src")
526 
527  if not self.config.doSubtract and self.config.doDetection:
528  # If we don't do subtraction, we need the subtracted exposure from the repo
529  subtractedExposure = sensorRef.get(subtractedExposureName)
530  # Both doSubtract and doDetection cannot be False
531 
532  results = self.run(exposure=exposure,
533  selectSources=selectSources,
534  templateExposure=template.exposure,
535  templateSources=template.sources,
536  idFactory=idFactory,
537  calexpBackgroundExposure=calexpBackgroundExposure,
538  subtractedExposure=subtractedExposure)
539 
540  if self.config.doWriteSources and results.diaSources is not None:
541  sensorRef.put(results.diaSources, self.config.coaddName + "Diff_diaSrc")
542  if self.config.doWriteWarpedExp:
543  sensorRef.put(results.warpedExposure, self.config.coaddName + "Diff_warpedExp")
544  if self.config.doWriteMatchedExp:
545  sensorRef.put(results.matchedExposure, self.config.coaddName + "Diff_matchedExp")
546  if self.config.doAddMetrics and self.config.doSelectSources:
547  sensorRef.put(results.selectSources, self.config.coaddName + "Diff_kernelSrc")
548  if self.config.doWriteSubtractedExp:
549  sensorRef.put(results.subtractedExposure, subtractedExposureName)
550  return results
551 
552  @pipeBase.timeMethod
553  def run(self, exposure=None, selectSources=None, templateExposure=None, templateSources=None,
554  idFactory=None, calexpBackgroundExposure=None, subtractedExposure=None):
555  """PSF matches, subtract two images and perform detection on the difference image.
556 
557  Parameters
558  ----------
559  exposure : `lsst.afw.image.ExposureF`, optional
560  The science exposure, the minuend in the image subtraction.
561  Can be None only if ``config.doSubtract==False``.
562  selectSources : `lsst.afw.table.SourceCatalog`, optional
563  Identified sources on the science exposure. This catalog is used to
564  select sources in order to perform the AL PSF matching on stamp images
565  around them. The selection steps depend on config options and whether
566  ``templateSources`` and ``matchingSources`` specified.
567  templateExposure : `lsst.afw.image.ExposureF`, optional
568  The template to be subtracted from ``exposure`` in the image subtraction.
569  The template exposure should cover the same sky area as the science exposure.
570  It is either a stich of patches of a coadd skymap image or a calexp
571  of the same pointing as the science exposure. Can be None only
572  if ``config.doSubtract==False`` and ``subtractedExposure`` is not None.
573  templateSources : `lsst.afw.table.SourceCatalog`, optional
574  Identified sources on the template exposure.
575  idFactory : `lsst.afw.table.IdFactory`
576  Generator object to assign ids to detected sources in the difference image.
577  calexpBackgroundExposure : `lsst.afw.image.ExposureF`, optional
578  Background exposure to be added back to the science exposure
579  if ``config.doAddCalexpBackground==True``
580  subtractedExposure : `lsst.afw.image.ExposureF`, optional
581  If ``config.doSubtract==False`` and ``config.doDetection==True``,
582  performs the post subtraction source detection only on this exposure.
583  Otherwise should be None.
584 
585  Returns
586  -------
587  results : `lsst.pipe.base.Struct`
588  ``subtractedExposure`` : `lsst.afw.image.ExposureF`
589  Difference image.
590  ``matchedExposure`` : `lsst.afw.image.ExposureF`
591  The matched PSF exposure.
592  ``subtractRes`` : `lsst.pipe.base.Struct`
593  The returned result structure of the ImagePsfMatchTask subtask.
594  ``diaSources`` : `lsst.afw.table.SourceCatalog`
595  The catalog of detected sources.
596  ``selectSources`` : `lsst.afw.table.SourceCatalog`
597  The input source catalog with optionally added Qa information.
598 
599  Notes
600  -----
601  The following major steps are included:
602 
603  - warp template coadd to match WCS of image
604  - PSF match image to warped template
605  - subtract image from PSF-matched, warped template
606  - detect sources
607  - measure sources
608 
609  For details about the image subtraction configuration modes
610  see `lsst.ip.diffim`.
611  """
612  subtractRes = None
613  controlSources = None
614  diaSources = None
615  kernelSources = None
616 
617  if self.config.doAddCalexpBackground:
618  mi = exposure.getMaskedImage()
619  mi += calexpBackgroundExposure.getImage()
620 
621  if not exposure.hasPsf():
622  raise pipeBase.TaskError("Exposure has no psf")
623  sciencePsf = exposure.getPsf()
624 
625  if self.config.doSubtract:
626  if self.config.doScaleTemplateVariance:
627  self.log.info("Rescaling template variance")
628  templateVarFactor = self.scaleVariance.run(
629  templateExposure.getMaskedImage())
630  self.log.info("Template variance scaling factor: %.2f" % templateVarFactor)
631  self.metadata.add("scaleTemplateVarianceFactor", templateVarFactor)
632 
633  if self.config.subtract.name == 'zogy':
634  subtractRes = self.subtract.run(exposure, templateExposure, doWarping=True)
635  if self.config.doPreConvolve:
636  subtractedExposure = subtractRes.scoreExp
637  else:
638  subtractedExposure = subtractRes.diffExp
639  subtractRes.subtractedExposure = subtractedExposure
640  subtractRes.matchedExposure = None
641 
642  elif self.config.subtract.name == 'al':
643  # compute scienceSigmaOrig: sigma of PSF of science image before pre-convolution
644  scienceSigmaOrig = sciencePsf.computeShape().getDeterminantRadius()
645  templateSigma = templateExposure.getPsf().computeShape().getDeterminantRadius()
646 
647  # if requested, convolve the science exposure with its PSF
648  # (properly, this should be a cross-correlation, but our code does not yet support that)
649  # compute scienceSigmaPost: sigma of science exposure with pre-convolution, if done,
650  # else sigma of original science exposure
651  # TODO: DM-22762 This functional block should be moved into its own method
652  preConvPsf = None
653  if self.config.doPreConvolve:
654  convControl = afwMath.ConvolutionControl()
655  # cannot convolve in place, so make a new MI to receive convolved image
656  srcMI = exposure.getMaskedImage()
657  exposureOrig = exposure.clone()
658  destMI = srcMI.Factory(srcMI.getDimensions())
659  srcPsf = sciencePsf
660  if self.config.useGaussianForPreConvolution:
661  # convolve with a simplified PSF model: a double Gaussian
662  kWidth, kHeight = sciencePsf.getLocalKernel().getDimensions()
663  preConvPsf = SingleGaussianPsf(kWidth, kHeight, scienceSigmaOrig)
664  else:
665  # convolve with science exposure's PSF model
666  preConvPsf = srcPsf
667  afwMath.convolve(destMI, srcMI, preConvPsf.getLocalKernel(), convControl)
668  exposure.setMaskedImage(destMI)
669  scienceSigmaPost = scienceSigmaOrig*math.sqrt(2)
670  else:
671  scienceSigmaPost = scienceSigmaOrig
672  exposureOrig = exposure
673 
674  # If requested, find and select sources from the image
675  # else, AL subtraction will do its own source detection
676  # TODO: DM-22762 This functional block should be moved into its own method
677  if self.config.doSelectSources:
678  if selectSources is None:
679  self.log.warn("Src product does not exist; running detection, measurement, selection")
680  # Run own detection and measurement; necessary in nightly processing
681  selectSources = self.subtract.getSelectSources(
682  exposure,
683  sigma=scienceSigmaPost,
684  doSmooth=not self.config.doPreConvolve,
685  idFactory=idFactory,
686  )
687 
688  if self.config.doAddMetrics:
689  # Number of basis functions
690 
691  nparam = len(makeKernelBasisList(self.subtract.config.kernel.active,
692  referenceFwhmPix=scienceSigmaPost*FwhmPerSigma,
693  targetFwhmPix=templateSigma*FwhmPerSigma))
694  # Modify the schema of all Sources
695  # DEPRECATED: This is a data dependent (nparam) output product schema
696  # outside the task constructor.
697  # NOTE: The pre-determination of nparam at this point
698  # may be incorrect as the template psf is warped later in
699  # ImagePsfMatchTask.matchExposures()
700  kcQa = KernelCandidateQa(nparam)
701  selectSources = kcQa.addToSchema(selectSources)
702  if self.config.kernelSourcesFromRef:
703  # match exposure sources to reference catalog
704  astromRet = self.astrometer.loadAndMatch(exposure=exposure, sourceCat=selectSources)
705  matches = astromRet.matches
706  elif templateSources:
707  # match exposure sources to template sources
708  mc = afwTable.MatchControl()
709  mc.findOnlyClosest = False
710  matches = afwTable.matchRaDec(templateSources, selectSources, 1.0*geom.arcseconds,
711  mc)
712  else:
713  raise RuntimeError("doSelectSources=True and kernelSourcesFromRef=False,"
714  "but template sources not available. Cannot match science "
715  "sources with template sources. Run process* on data from "
716  "which templates are built.")
717 
718  kernelSources = self.sourceSelector.run(selectSources, exposure=exposure,
719  matches=matches).sourceCat
720  random.shuffle(kernelSources, random.random)
721  controlSources = kernelSources[::self.config.controlStepSize]
722  kernelSources = [k for i, k in enumerate(kernelSources)
723  if i % self.config.controlStepSize]
724 
725  if self.config.doSelectDcrCatalog:
726  redSelector = DiaCatalogSourceSelectorTask(
727  DiaCatalogSourceSelectorConfig(grMin=self.sourceSelector.config.grMax,
728  grMax=99.999))
729  redSources = redSelector.selectStars(exposure, selectSources, matches=matches).starCat
730  controlSources.extend(redSources)
731 
732  blueSelector = DiaCatalogSourceSelectorTask(
733  DiaCatalogSourceSelectorConfig(grMin=-99.999,
734  grMax=self.sourceSelector.config.grMin))
735  blueSources = blueSelector.selectStars(exposure, selectSources,
736  matches=matches).starCat
737  controlSources.extend(blueSources)
738 
739  if self.config.doSelectVariableCatalog:
740  varSelector = DiaCatalogSourceSelectorTask(
741  DiaCatalogSourceSelectorConfig(includeVariable=True))
742  varSources = varSelector.selectStars(exposure, selectSources, matches=matches).starCat
743  controlSources.extend(varSources)
744 
745  self.log.info("Selected %d / %d sources for Psf matching (%d for control sample)"
746  % (len(kernelSources), len(selectSources), len(controlSources)))
747 
748  allresids = {}
749  # TODO: DM-22762 This functional block should be moved into its own method
750  if self.config.doUseRegister:
751  self.log.info("Registering images")
752 
753  if templateSources is None:
754  # Run detection on the template, which is
755  # temporarily background-subtracted
756  # sigma of PSF of template image before warping
757  templateSigma = templateExposure.getPsf().computeShape().getDeterminantRadius()
758  templateSources = self.subtract.getSelectSources(
759  templateExposure,
760  sigma=templateSigma,
761  doSmooth=True,
762  idFactory=idFactory
763  )
764 
765  # Third step: we need to fit the relative astrometry.
766  #
767  wcsResults = self.fitAstrometry(templateSources, templateExposure, selectSources)
768  warpedExp = self.register.warpExposure(templateExposure, wcsResults.wcs,
769  exposure.getWcs(), exposure.getBBox())
770  templateExposure = warpedExp
771 
772  # Create debugging outputs on the astrometric
773  # residuals as a function of position. Persistence
774  # not yet implemented; expected on (I believe) #2636.
775  if self.config.doDebugRegister:
776  # Grab matches to reference catalog
777  srcToMatch = {x.second.getId(): x.first for x in matches}
778 
779  refCoordKey = wcsResults.matches[0].first.getTable().getCoordKey()
780  inCentroidKey = wcsResults.matches[0].second.getTable().getCentroidSlot().getMeasKey()
781  sids = [m.first.getId() for m in wcsResults.matches]
782  positions = [m.first.get(refCoordKey) for m in wcsResults.matches]
783  residuals = [m.first.get(refCoordKey).getOffsetFrom(wcsResults.wcs.pixelToSky(
784  m.second.get(inCentroidKey))) for m in wcsResults.matches]
785  allresids = dict(zip(sids, zip(positions, residuals)))
786 
787  cresiduals = [m.first.get(refCoordKey).getTangentPlaneOffset(
788  wcsResults.wcs.pixelToSky(
789  m.second.get(inCentroidKey))) for m in wcsResults.matches]
790  colors = numpy.array([-2.5*numpy.log10(srcToMatch[x].get("g"))
791  + 2.5*numpy.log10(srcToMatch[x].get("r"))
792  for x in sids if x in srcToMatch.keys()])
793  dlong = numpy.array([r[0].asArcseconds() for s, r in zip(sids, cresiduals)
794  if s in srcToMatch.keys()])
795  dlat = numpy.array([r[1].asArcseconds() for s, r in zip(sids, cresiduals)
796  if s in srcToMatch.keys()])
797  idx1 = numpy.where(colors < self.sourceSelector.config.grMin)
798  idx2 = numpy.where((colors >= self.sourceSelector.config.grMin)
799  & (colors <= self.sourceSelector.config.grMax))
800  idx3 = numpy.where(colors > self.sourceSelector.config.grMax)
801  rms1Long = IqrToSigma*(
802  (numpy.percentile(dlong[idx1], 75) - numpy.percentile(dlong[idx1], 25)))
803  rms1Lat = IqrToSigma*(numpy.percentile(dlat[idx1], 75)
804  - numpy.percentile(dlat[idx1], 25))
805  rms2Long = IqrToSigma*(
806  (numpy.percentile(dlong[idx2], 75) - numpy.percentile(dlong[idx2], 25)))
807  rms2Lat = IqrToSigma*(numpy.percentile(dlat[idx2], 75)
808  - numpy.percentile(dlat[idx2], 25))
809  rms3Long = IqrToSigma*(
810  (numpy.percentile(dlong[idx3], 75) - numpy.percentile(dlong[idx3], 25)))
811  rms3Lat = IqrToSigma*(numpy.percentile(dlat[idx3], 75)
812  - numpy.percentile(dlat[idx3], 25))
813  self.log.info("Blue star offsets'': %.3f %.3f, %.3f %.3f" %
814  (numpy.median(dlong[idx1]), rms1Long,
815  numpy.median(dlat[idx1]), rms1Lat))
816  self.log.info("Green star offsets'': %.3f %.3f, %.3f %.3f" %
817  (numpy.median(dlong[idx2]), rms2Long,
818  numpy.median(dlat[idx2]), rms2Lat))
819  self.log.info("Red star offsets'': %.3f %.3f, %.3f %.3f" %
820  (numpy.median(dlong[idx3]), rms3Long,
821  numpy.median(dlat[idx3]), rms3Lat))
822 
823  self.metadata.add("RegisterBlueLongOffsetMedian", numpy.median(dlong[idx1]))
824  self.metadata.add("RegisterGreenLongOffsetMedian", numpy.median(dlong[idx2]))
825  self.metadata.add("RegisterRedLongOffsetMedian", numpy.median(dlong[idx3]))
826  self.metadata.add("RegisterBlueLongOffsetStd", rms1Long)
827  self.metadata.add("RegisterGreenLongOffsetStd", rms2Long)
828  self.metadata.add("RegisterRedLongOffsetStd", rms3Long)
829 
830  self.metadata.add("RegisterBlueLatOffsetMedian", numpy.median(dlat[idx1]))
831  self.metadata.add("RegisterGreenLatOffsetMedian", numpy.median(dlat[idx2]))
832  self.metadata.add("RegisterRedLatOffsetMedian", numpy.median(dlat[idx3]))
833  self.metadata.add("RegisterBlueLatOffsetStd", rms1Lat)
834  self.metadata.add("RegisterGreenLatOffsetStd", rms2Lat)
835  self.metadata.add("RegisterRedLatOffsetStd", rms3Lat)
836 
837  # warp template exposure to match exposure,
838  # PSF match template exposure to exposure,
839  # then return the difference
840 
841  # Return warped template... Construct sourceKernelCand list after subtract
842  self.log.info("Subtracting images")
843  subtractRes = self.subtract.subtractExposures(
844  templateExposure=templateExposure,
845  scienceExposure=exposure,
846  candidateList=kernelSources,
847  convolveTemplate=self.config.convolveTemplate,
848  doWarping=not self.config.doUseRegister
849  )
850  subtractedExposure = subtractRes.subtractedExposure
851 
852  if self.config.doDetection:
853  self.log.info("Computing diffim PSF")
854 
855  # Get Psf from the appropriate input image if it doesn't exist
856  if not subtractedExposure.hasPsf():
857  if self.config.convolveTemplate:
858  subtractedExposure.setPsf(exposure.getPsf())
859  else:
860  subtractedExposure.setPsf(templateExposure.getPsf())
861 
862  # If doSubtract is False, then subtractedExposure was fetched from disk (above),
863  # thus it may have already been decorrelated. Thus, we do not decorrelate if
864  # doSubtract is False.
865 
866  # NOTE: At this point doSubtract == True
867  if self.config.doDecorrelation and self.config.doSubtract:
868  preConvKernel = None
869  if preConvPsf is not None:
870  preConvKernel = preConvPsf.getLocalKernel()
871  decorrResult = self.decorrelate.run(exposureOrig, subtractRes.warpedExposure,
872  subtractedExposure,
873  subtractRes.psfMatchingKernel,
874  spatiallyVarying=self.config.doSpatiallyVarying,
875  preConvKernel=preConvKernel,
876  templateMatched=self.config.convolveTemplate)
877  subtractedExposure = decorrResult.correctedExposure
878 
879  # END (if subtractAlgorithm == 'AL')
880  # END (if self.config.doSubtract)
881  if self.config.doDetection:
882  self.log.info("Running diaSource detection")
883 
884  # Rescale difference image variance plane
885  if self.config.doScaleDiffimVariance:
886  self.log.info("Rescaling diffim variance")
887  diffimVarFactor = self.scaleVariance.run(subtractedExposure.getMaskedImage())
888  self.log.info("Diffim variance scaling factor: %.2f" % diffimVarFactor)
889  self.metadata.add("scaleDiffimVarianceFactor", diffimVarFactor)
890 
891  # Erase existing detection mask planes
892  mask = subtractedExposure.getMaskedImage().getMask()
893  mask &= ~(mask.getPlaneBitMask("DETECTED") | mask.getPlaneBitMask("DETECTED_NEGATIVE"))
894 
895  table = afwTable.SourceTable.make(self.schema, idFactory)
896  table.setMetadata(self.algMetadata)
897  results = self.detection.run(
898  table=table,
899  exposure=subtractedExposure,
900  doSmooth=not self.config.doPreConvolve
901  )
902 
903  if self.config.doMerge:
904  fpSet = results.fpSets.positive
905  fpSet.merge(results.fpSets.negative, self.config.growFootprint,
906  self.config.growFootprint, False)
907  diaSources = afwTable.SourceCatalog(table)
908  fpSet.makeSources(diaSources)
909  self.log.info("Merging detections into %d sources" % (len(diaSources)))
910  else:
911  diaSources = results.sources
912 
913  if self.config.doMeasurement:
914  newDipoleFitting = self.config.doDipoleFitting
915  self.log.info("Running diaSource measurement: newDipoleFitting=%r", newDipoleFitting)
916  if not newDipoleFitting:
917  # Just fit dipole in diffim
918  self.measurement.run(diaSources, subtractedExposure)
919  else:
920  # Use (matched) template and science image (if avail.) to constrain dipole fitting
921  if self.config.doSubtract and 'matchedExposure' in subtractRes.getDict():
922  self.measurement.run(diaSources, subtractedExposure, exposure,
923  subtractRes.matchedExposure)
924  else:
925  self.measurement.run(diaSources, subtractedExposure, exposure)
926  if self.config.doApCorr:
927  self.applyApCorr.run(
928  catalog=diaSources,
929  apCorrMap=subtractedExposure.getInfo().getApCorrMap()
930  )
931 
932  if self.config.doForcedMeasurement:
933  # Run forced psf photometry on the PVI at the diaSource locations.
934  # Copy the measured flux and error into the diaSource.
935  forcedSources = self.forcedMeasurement.generateMeasCat(
936  exposure, diaSources, subtractedExposure.getWcs())
937  self.forcedMeasurement.run(forcedSources, exposure, diaSources, subtractedExposure.getWcs())
938  mapper = afwTable.SchemaMapper(forcedSources.schema, diaSources.schema)
939  mapper.addMapping(forcedSources.schema.find("base_PsfFlux_instFlux")[0],
940  "ip_diffim_forced_PsfFlux_instFlux", True)
941  mapper.addMapping(forcedSources.schema.find("base_PsfFlux_instFluxErr")[0],
942  "ip_diffim_forced_PsfFlux_instFluxErr", True)
943  mapper.addMapping(forcedSources.schema.find("base_PsfFlux_area")[0],
944  "ip_diffim_forced_PsfFlux_area", True)
945  mapper.addMapping(forcedSources.schema.find("base_PsfFlux_flag")[0],
946  "ip_diffim_forced_PsfFlux_flag", True)
947  mapper.addMapping(forcedSources.schema.find("base_PsfFlux_flag_noGoodPixels")[0],
948  "ip_diffim_forced_PsfFlux_flag_noGoodPixels", True)
949  mapper.addMapping(forcedSources.schema.find("base_PsfFlux_flag_edge")[0],
950  "ip_diffim_forced_PsfFlux_flag_edge", True)
951  for diaSource, forcedSource in zip(diaSources, forcedSources):
952  diaSource.assign(forcedSource, mapper)
953 
954  # Match with the calexp sources if possible
955  if self.config.doMatchSources:
956  if selectSources is not None:
957  # Create key,val pair where key=diaSourceId and val=sourceId
958  matchRadAsec = self.config.diaSourceMatchRadius
959  matchRadPixel = matchRadAsec/exposure.getWcs().getPixelScale().asArcseconds()
960 
961  srcMatches = afwTable.matchXy(selectSources, diaSources, matchRadPixel)
962  srcMatchDict = dict([(srcMatch.second.getId(), srcMatch.first.getId()) for
963  srcMatch in srcMatches])
964  self.log.info("Matched %d / %d diaSources to sources" % (len(srcMatchDict),
965  len(diaSources)))
966  else:
967  self.log.warn("Src product does not exist; cannot match with diaSources")
968  srcMatchDict = {}
969 
970  # Create key,val pair where key=diaSourceId and val=refId
971  refAstromConfig = AstrometryConfig()
972  refAstromConfig.matcher.maxMatchDistArcSec = matchRadAsec
973  refAstrometer = AstrometryTask(refAstromConfig)
974  astromRet = refAstrometer.run(exposure=exposure, sourceCat=diaSources)
975  refMatches = astromRet.matches
976  if refMatches is None:
977  self.log.warn("No diaSource matches with reference catalog")
978  refMatchDict = {}
979  else:
980  self.log.info("Matched %d / %d diaSources to reference catalog" % (len(refMatches),
981  len(diaSources)))
982  refMatchDict = dict([(refMatch.second.getId(), refMatch.first.getId()) for
983  refMatch in refMatches])
984 
985  # Assign source Ids
986  for diaSource in diaSources:
987  sid = diaSource.getId()
988  if sid in srcMatchDict:
989  diaSource.set("srcMatchId", srcMatchDict[sid])
990  if sid in refMatchDict:
991  diaSource.set("refMatchId", refMatchDict[sid])
992 
993  if self.config.doAddMetrics and self.config.doSelectSources:
994  self.log.info("Evaluating metrics and control sample")
995 
996  kernelCandList = []
997  for cell in subtractRes.kernelCellSet.getCellList():
998  for cand in cell.begin(False): # include bad candidates
999  kernelCandList.append(cand)
1000 
1001  # Get basis list to build control sample kernels
1002  basisList = kernelCandList[0].getKernel(KernelCandidateF.ORIG).getKernelList()
1003  nparam = len(kernelCandList[0].getKernel(KernelCandidateF.ORIG).getKernelParameters())
1004 
1005  controlCandList = (
1006  diffimTools.sourceTableToCandidateList(controlSources,
1007  subtractRes.warpedExposure, exposure,
1008  self.config.subtract.kernel.active,
1009  self.config.subtract.kernel.active.detectionConfig,
1010  self.log, doBuild=True, basisList=basisList))
1011 
1012  KernelCandidateQa.apply(kernelCandList, subtractRes.psfMatchingKernel,
1013  subtractRes.backgroundModel, dof=nparam)
1014  KernelCandidateQa.apply(controlCandList, subtractRes.psfMatchingKernel,
1015  subtractRes.backgroundModel)
1016 
1017  if self.config.doDetection:
1018  KernelCandidateQa.aggregate(selectSources, self.metadata, allresids, diaSources)
1019  else:
1020  KernelCandidateQa.aggregate(selectSources, self.metadata, allresids)
1021 
1022  self.runDebug(exposure, subtractRes, selectSources, kernelSources, diaSources)
1023  return pipeBase.Struct(
1024  subtractedExposure=subtractedExposure,
1025  warpedExposure=subtractRes.warpedExposure,
1026  matchedExposure=subtractRes.matchedExposure,
1027  subtractRes=subtractRes,
1028  diaSources=diaSources,
1029  selectSources=selectSources
1030  )
1031 
1032  def fitAstrometry(self, templateSources, templateExposure, selectSources):
1033  """Fit the relative astrometry between templateSources and selectSources
1034 
1035  Todo
1036  ----
1037 
1038  Remove this method. It originally fit a new WCS to the template before calling register.run
1039  because our TAN-SIP fitter behaved badly for points far from CRPIX, but that's been fixed.
1040  It remains because a subtask overrides it.
1041  """
1042  results = self.register.run(templateSources, templateExposure.getWcs(),
1043  templateExposure.getBBox(), selectSources)
1044  return results
1045 
1046  def runDebug(self, exposure, subtractRes, selectSources, kernelSources, diaSources):
1047  """Make debug plots and displays.
1048 
1049  Todo
1050  ----
1051  Test and update for current debug display and slot names
1052  """
1053  import lsstDebug
1054  display = lsstDebug.Info(__name__).display
1055  showSubtracted = lsstDebug.Info(__name__).showSubtracted
1056  showPixelResiduals = lsstDebug.Info(__name__).showPixelResiduals
1057  showDiaSources = lsstDebug.Info(__name__).showDiaSources
1058  showDipoles = lsstDebug.Info(__name__).showDipoles
1059  maskTransparency = lsstDebug.Info(__name__).maskTransparency
1060  if display:
1061  disp = afwDisplay.getDisplay(frame=lsstDebug.frame)
1062  if not maskTransparency:
1063  maskTransparency = 0
1064  disp.setMaskTransparency(maskTransparency)
1065 
1066  if display and showSubtracted:
1067  disp.mtv(subtractRes.subtractedExposure, title="Subtracted image")
1068  mi = subtractRes.subtractedExposure.getMaskedImage()
1069  x0, y0 = mi.getX0(), mi.getY0()
1070  with disp.Buffering():
1071  for s in diaSources:
1072  x, y = s.getX() - x0, s.getY() - y0
1073  ctype = "red" if s.get("flags_negative") else "yellow"
1074  if (s.get("base_PixelFlags_flag_interpolatedCenter")
1075  or s.get("base_PixelFlags_flag_saturatedCenter")
1076  or s.get("base_PixelFlags_flag_crCenter")):
1077  ptype = "x"
1078  elif (s.get("base_PixelFlags_flag_interpolated")
1079  or s.get("base_PixelFlags_flag_saturated")
1080  or s.get("base_PixelFlags_flag_cr")):
1081  ptype = "+"
1082  else:
1083  ptype = "o"
1084  disp.dot(ptype, x, y, size=4, ctype=ctype)
1085  lsstDebug.frame += 1
1086 
1087  if display and showPixelResiduals and selectSources:
1088  nonKernelSources = []
1089  for source in selectSources:
1090  if source not in kernelSources:
1091  nonKernelSources.append(source)
1092 
1093  diUtils.plotPixelResiduals(exposure,
1094  subtractRes.warpedExposure,
1095  subtractRes.subtractedExposure,
1096  subtractRes.kernelCellSet,
1097  subtractRes.psfMatchingKernel,
1098  subtractRes.backgroundModel,
1099  nonKernelSources,
1100  self.subtract.config.kernel.active.detectionConfig,
1101  origVariance=False)
1102  diUtils.plotPixelResiduals(exposure,
1103  subtractRes.warpedExposure,
1104  subtractRes.subtractedExposure,
1105  subtractRes.kernelCellSet,
1106  subtractRes.psfMatchingKernel,
1107  subtractRes.backgroundModel,
1108  nonKernelSources,
1109  self.subtract.config.kernel.active.detectionConfig,
1110  origVariance=True)
1111  if display and showDiaSources:
1112  flagChecker = SourceFlagChecker(diaSources)
1113  isFlagged = [flagChecker(x) for x in diaSources]
1114  isDipole = [x.get("ip_diffim_ClassificationDipole_value") for x in diaSources]
1115  diUtils.showDiaSources(diaSources, subtractRes.subtractedExposure, isFlagged, isDipole,
1116  frame=lsstDebug.frame)
1117  lsstDebug.frame += 1
1118 
1119  if display and showDipoles:
1120  DipoleAnalysis().displayDipoles(subtractRes.subtractedExposure, diaSources,
1121  frame=lsstDebug.frame)
1122  lsstDebug.frame += 1
1123 
1124  def _getConfigName(self):
1125  """Return the name of the config dataset
1126  """
1127  return "%sDiff_config" % (self.config.coaddName,)
1128 
1129  def _getMetadataName(self):
1130  """Return the name of the metadata dataset
1131  """
1132  return "%sDiff_metadata" % (self.config.coaddName,)
1133 
1134  def getSchemaCatalogs(self):
1135  """Return a dict of empty catalogs for each catalog dataset produced by this task."""
1136  return {self.config.coaddName + "Diff_diaSrc": self.outputSchema}
1137 
1138  @classmethod
1139  def _makeArgumentParser(cls):
1140  """Create an argument parser
1141  """
1142  parser = pipeBase.ArgumentParser(name=cls._DefaultName)
1143  parser.add_id_argument("--id", "calexp", help="data ID, e.g. --id visit=12345 ccd=1,2")
1144  parser.add_id_argument("--templateId", "calexp", doMakeDataRefList=True,
1145  help="Template data ID in case of calexp template,"
1146  " e.g. --templateId visit=6789")
1147  return parser
1148 
1149 
1150 class Winter2013ImageDifferenceConfig(ImageDifferenceConfig):
1151  winter2013WcsShift = pexConfig.Field(dtype=float, default=0.0,
1152  doc="Shift stars going into RegisterTask by this amount")
1153  winter2013WcsRms = pexConfig.Field(dtype=float, default=0.0,
1154  doc="Perturb stars going into RegisterTask by this amount")
1155 
1156  def setDefaults(self):
1157  ImageDifferenceConfig.setDefaults(self)
1158  self.getTemplate.retarget(GetCalexpAsTemplateTask)
1159 
1160 
1161 class Winter2013ImageDifferenceTask(ImageDifferenceTask):
1162  """!Image difference Task used in the Winter 2013 data challege.
1163  Enables testing the effects of registration shifts and scatter.
1164 
1165  For use with winter 2013 simulated images:
1166  Use --templateId visit=88868666 for sparse data
1167  --templateId visit=22222200 for dense data (g)
1168  --templateId visit=11111100 for dense data (i)
1169  """
1170  ConfigClass = Winter2013ImageDifferenceConfig
1171  _DefaultName = "winter2013ImageDifference"
1172 
1173  def __init__(self, **kwargs):
1174  ImageDifferenceTask.__init__(self, **kwargs)
1175 
1176  def fitAstrometry(self, templateSources, templateExposure, selectSources):
1177  """Fit the relative astrometry between templateSources and selectSources"""
1178  if self.config.winter2013WcsShift > 0.0:
1179  offset = geom.Extent2D(self.config.winter2013WcsShift,
1180  self.config.winter2013WcsShift)
1181  cKey = templateSources[0].getTable().getCentroidSlot().getMeasKey()
1182  for source in templateSources:
1183  centroid = source.get(cKey)
1184  source.set(cKey, centroid + offset)
1185  elif self.config.winter2013WcsRms > 0.0:
1186  cKey = templateSources[0].getTable().getCentroidSlot().getMeasKey()
1187  for source in templateSources:
1188  offset = geom.Extent2D(self.config.winter2013WcsRms*numpy.random.normal(),
1189  self.config.winter2013WcsRms*numpy.random.normal())
1190  centroid = source.get(cKey)
1191  source.set(cKey, centroid + offset)
1192 
1193  results = self.register.run(templateSources, templateExposure.getWcs(),
1194  templateExposure.getBBox(), selectSources)
1195  return results
def run(self, skyInfo, tempExpRefList, imageScalerList, weightList, altMaskList=None, mask=None, supplementaryData=None)