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