lsst.meas.base  14.0-12-g233aa8e+6
tests.py
Go to the documentation of this file.
1 #!/usr/bin/env python
2 #
3 # LSST Data Management System
4 # Copyright 2008-2015 AURA/LSST.
5 #
6 # This product includes software developed by the
7 # LSST Project (http://www.lsst.org/).
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 LSST License Statement and
20 # the GNU General Public License along with this program. If not,
21 # see <http://www.lsstcorp.org/LegalNotices/>.
22 #
23 
24 from builtins import zip
25 from builtins import object
26 import numpy as np
27 
28 import lsst.afw.table
29 import lsst.afw.image
30 import lsst.afw.detection
31 import lsst.afw.geom
33 import lsst.afw.coord
35 
36 from .sfm import SingleFrameMeasurementTask
37 from .forcedMeasurement import ForcedMeasurementTask
38 from . import CentroidResultKey
39 
40 __all__ = ("BlendContext", "TestDataset", "AlgorithmTestCase", "TransformTestCase",
41  "SingleFramePluginTransformSetupHelper", "ForcedPluginTransformSetupHelper",
42  "FluxTransformTestCase", "CentroidTransformTestCase")
43 
44 
45 class BlendContext(object):
46  """!
47  A Python context manager used to add multiple overlapping sources along with a parent source
48  that represents all of them together.
49 
50  This is used as the return value for TestDataset.addBlend(), and this is the only way it should
51  be used. The only public method is addChild().
52  """
53 
54  def __init__(self, owner):
55  self.owner = owner
56  self.parentRecord = self.owner.catalog.addNew()
57  self.parentImage = lsst.afw.image.ImageF(self.owner.exposure.getBBox())
58  self.children = []
59 
60  def __enter__(self):
61  # BlendContext is its own context manager, so we just return self.
62  return self
63 
64  def addChild(self, flux, centroid, shape=None):
65  """!
66  Add a child source to the blend, and return the truth catalog record that corresponds to it.
67 
68  @param[in] flux Total flux of the source to be added.
69  @param[in] centroid Position of the source to be added (lsst.afw.geom.Point2D).
70  @param[in] shape 2nd moments of the source before PSF convolution
71  (lsst.afw.geom.ellipses.Quadrupole). Note that the truth catalog
72  records post-convolution moments)
73  """
74  record, image = self.owner.addSource(flux, centroid, shape)
75  record.set(self.owner.keys["parent"], self.parentRecord.getId())
76  self.parentImage += image
77  self.children.append((record, image))
78  return record
79 
80  def __exit__(self, type_, value, tb):
81  # We're not using the context manager for any kind of exception safety or guarantees;
82  # we just want the nice "with" statement syntax.
83  if type_ is not None: # exception was raised; just skip all this and let it propagate
84  return
85  # On exit, we need to compute and set the truth values for the parent object.
86  self.parentRecord.set(self.owner.keys["nChild"], len(self.children))
87  # Compute flux from sum of component fluxes
88  flux = 0.0
89  for record, image in self.children:
90  flux += record.get(self.owner.keys["flux"])
91  self.parentRecord.set(self.owner.keys["flux"], flux)
92  # Compute centroid from flux-weighted mean of component centroids
93  x = 0.0
94  y = 0.0
95  for record, image in self.children:
96  w = record.get(self.owner.keys["flux"])/flux
97  x += record.get(self.owner.keys["centroid"].getX())*w
98  y += record.get(self.owner.keys["centroid"].getY())*w
99  self.parentRecord.set(self.owner.keys["centroid"], lsst.afw.geom.Point2D(x, y))
100  # Compute shape from flux-weighted mean of offset component shapes
101  xx = 0.0
102  yy = 0.0
103  xy = 0.0
104  for record, image in self.children:
105  w = record.get(self.owner.keys["flux"])/flux
106  dx = record.get(self.owner.keys["centroid"].getX()) - x
107  dy = record.get(self.owner.keys["centroid"].getY()) - y
108  xx += (record.get(self.owner.keys["shape"].getIxx()) + dx**2)*w
109  yy += (record.get(self.owner.keys["shape"].getIyy()) + dy**2)*w
110  xy += (record.get(self.owner.keys["shape"].getIxy()) + dx*dy)*w
111  self.parentRecord.set(self.owner.keys["shape"], lsst.afw.geom.ellipses.Quadrupole(xx, yy, xy))
112  # Run detection on the parent image to get the parent Footprint.
113  self.owner._installFootprint(self.parentRecord, self.parentImage)
114  # Create perfect HeavyFootprints for all children; these will need to be modified later to account
115  # for the noise we'll add to the image.
116  deblend = lsst.afw.image.MaskedImageF(self.owner.exposure.getMaskedImage(), True)
117  for record, image in self.children:
118  deblend.getImage().getArray()[:, :] = image.getArray()
119  heavyFootprint = lsst.afw.detection.HeavyFootprintF(self.parentRecord.getFootprint(), deblend)
120  record.setFootprint(heavyFootprint)
121 
122 
123 class TestDataset(object):
124  """!
125  A simulated dataset consisting of a test image and an associated truth catalog.
126 
127  TestDataset creates an idealized image made of pure Gaussians (including a Gaussian PSF),
128  with simple noise and idealized Footprints/HeavyFootprints that simulated the outputs
129  of detection and deblending. Multiple noise realizations can be created from the same
130  underlying sources, allowing uncertainty estimates to be verified via Monte Carlo.
131 
132  Typical usage:
133  @code
134  bbox = lsst.afw.geom.Box2I(lsst.afw.geom.Point2I(0,0), lsst.afw.geom.Point2I(100, 100))
135  dataset = TestDataset(bbox)
136  dataset.addSource(flux=1E5, centroid=lsst.afw.geom.Point2D(25, 26))
137  dataset.addSource(flux=2E5, centroid=lsst.afw.geom.Point2D(75, 24),
138  shape=lsst.afw.geom.ellipses.Quadrupole(8, 7, 2))
139  with dataset.addBlend() as family:
140  family.addChild(flux=2E5, centroid=lsst.afw.geom.Point2D(50, 72))
141  family.addChild(flux=1.5E5, centroid=lsst.afw.geom.Point2D(51, 74))
142  exposure, catalog = dataset.realize(noise=100.0, schema=TestDataset.makeMinimalSchema())
143  @endcode
144  """
145 
146  @classmethod
148  """Return the minimal schema needed to hold truth catalog fields.
149 
150  When TestDataset.realize() is called, the schema must include at least these fields.
151  Usually it will include additional fields for measurement algorithm outputs, allowing
152  the same catalog to be used for both truth values (the fields from the minimal schema)
153  and the measurements.
154  """
155  if not hasattr(cls, "_schema"):
157  cls.keys = {}
158  cls.keys["parent"] = schema.find("parent").key
159  cls.keys["nChild"] = schema.addField("deblend_nChild", type=np.int32)
160  cls.keys["flux"] = schema.addField("truth_flux", type=np.float64, doc="true flux", units="count")
161  cls.keys["centroid"] = lsst.afw.table.Point2DKey.addFields(
162  schema, "truth", "true simulated centroid", "pixel"
163  )
164  cls.keys["centroid_flag"] = schema.addField("truth_flag", type="Flag",
165  doc="set if the object is a star")
167  schema, "truth", "true shape after PSF convolution", lsst.afw.table.CoordinateType.PIXEL
168  )
169  cls.keys["isStar"] = schema.addField("truth_isStar", type="Flag",
170  doc="set if the object is a star")
171  schema.getAliasMap().set("slot_Shape", "truth")
172  schema.getAliasMap().set("slot_Centroid", "truth")
173  schema.getAliasMap().set("slot_ModelFlux", "truth")
174  schema.getCitizen().markPersistent()
175  cls._schema = schema
176  schema = lsst.afw.table.Schema(cls._schema)
177  schema.disconnectAliases()
178  return schema
179 
180  @staticmethod
181  def makePerturbedWcs(oldWcs, minScaleFactor=1.2, maxScaleFactor=1.5,
182  minRotation=None, maxRotation=None,
183  minRefShift=None, maxRefShift=None,
184  minPixShift=2.0, maxPixShift=4.0, randomSeed=None):
185  """!
186  Create a new undistorted TanWcs that is similar but not identical to another, with random
187  scaling, rotation, and offset (in both pixel position and reference position).
188 
189  The maximum and minimum arguments are interpreted as absolute values for a split
190  range that covers both positive and negative values (as this method is used
191  in testing, it is typically most important to avoid perturbations near zero).
192  Scale factors are treated somewhat differently: the actual scale factor is chosen between
193  minScaleFactor and maxScaleFactor OR (1/maxScaleFactor) and (1/minScaleFactor).
194 
195  The default range for rotation is 30-60 degrees, and the default range for reference shift
196  is 0.5-1.0 arcseconds (these cannot be safely included directly as default values because Angle
197  objects are mutable).
198 
199  The random number generator is primed with the seed given. If ``None``, a seed is
200  automatically chosen.
201  """
202  random_state = np.random.RandomState(randomSeed)
203  if minRotation is None:
204  minRotation = 30.0*lsst.afw.geom.degrees
205  if maxRotation is None:
206  maxRotation = 60.0*lsst.afw.geom.degrees
207  if minRefShift is None:
208  minRefShift = 0.5*lsst.afw.geom.arcseconds
209  if maxRefShift is None:
210  maxRefShift = 1.0*lsst.afw.geom.arcseconds
211 
212  def splitRandom(min1, max1, min2=None, max2=None):
213  if min2 is None:
214  min2 = -max1
215  if max2 is None:
216  max2 = -min1
217  if random_state.uniform() > 0.5:
218  return float(random_state.uniform(min1, max1))
219  else:
220  return float(random_state.uniform(min2, max2))
221  # Generate random perturbations
222  scaleFactor = splitRandom(minScaleFactor, maxScaleFactor, 1.0/maxScaleFactor, 1.0/minScaleFactor)
223  rotation = splitRandom(minRotation.asRadians(), maxRotation.asRadians())*lsst.afw.geom.radians
224  refShiftRa = splitRandom(minRefShift.asRadians(), maxRefShift.asRadians())*lsst.afw.geom.radians
225  refShiftDec = splitRandom(minRefShift.asRadians(), maxRefShift.asRadians())*lsst.afw.geom.radians
226  pixShiftX = splitRandom(minPixShift, maxPixShift)
227  pixShiftY = splitRandom(minPixShift, maxPixShift)
228  # Compute new CD matrix
229  oldTransform = lsst.afw.geom.LinearTransform(oldWcs.getCDMatrix())
230  rTransform = lsst.afw.geom.LinearTransform.makeRotation(rotation)
231  sTransform = lsst.afw.geom.LinearTransform.makeScaling(scaleFactor)
232  newTransform = oldTransform*rTransform*sTransform
233  matrix = newTransform.getMatrix()
234  # Compute new coordinate reference pixel (CRVAL)
235  oldSkyOrigin = oldWcs.getSkyOrigin().toIcrs()
236  newSkyOrigin = lsst.afw.coord.IcrsCoord(oldSkyOrigin.getRa() + refShiftRa,
237  oldSkyOrigin.getDec() + refShiftDec)
238  # Compute new pixel reference pixel (CRPIX)
239  oldPixOrigin = oldWcs.getPixelOrigin()
240  newPixOrigin = lsst.afw.geom.Point2D(oldPixOrigin.getX() + pixShiftX,
241  oldPixOrigin.getY() + pixShiftY)
242  return lsst.afw.image.makeWcs(newSkyOrigin, newPixOrigin,
243  matrix[0, 0], matrix[0, 1], matrix[1, 0], matrix[1, 1])
244 
245  @staticmethod
246  def makeEmptyExposure(bbox, wcs=None, crval=None, cdelt=None, psfSigma=2.0, psfDim=17, fluxMag0=1E12):
247  """!
248  Create an Exposure, with a Calib, Wcs, and Psf, but no pixel values set.
249 
250  @param[in] bbox Bounding box of the image (image coordinates) as returned by makeCatalog.
251  @param[in] wcs New Wcs for the exposure (created from crval and cdelt if None).
252  @param[in] crval afw.coord.Coord: center of the TAN WCS attached to the image.
253  @param[in] cdelt afw.geom.Angle: pixel scale of the image
254  @param[in] psfSigma Radius (sigma) of the Gaussian PSF attached to the image
255  @param[in] psfDim Width and height of the image's Gaussian PSF attached to the image
256  @param[in] fluxMag0 Flux at magnitude zero (in e-) used to set the Calib of the exposure.
257  """
258  if wcs is None:
259  if crval is None:
260  crval = lsst.afw.coord.IcrsCoord(45.0*lsst.afw.geom.degrees, 45.0*lsst.afw.geom.degrees)
261  if cdelt is None:
262  cdelt = 0.2*lsst.afw.geom.arcseconds
263  crpix = lsst.afw.geom.Box2D(bbox).getCenter()
264  wcs = lsst.afw.image.makeWcs(crval, crpix, cdelt.asDegrees(), 0.0, 0.0, cdelt.asDegrees())
265  exposure = lsst.afw.image.ExposureF(bbox)
266  psf = lsst.afw.detection.GaussianPsf(psfDim, psfDim, psfSigma)
267  calib = lsst.afw.image.Calib()
268  calib.setFluxMag0(fluxMag0)
269  exposure.setWcs(wcs)
270  exposure.setPsf(psf)
271  exposure.setCalib(calib)
272  return exposure
273 
274  @staticmethod
275  def drawGaussian(bbox, flux, ellipse):
276  """!
277  Create an image of an elliptical Gaussian.
278 
279  @param[in,out] bbox Bounding box of image to create.
280  @param[in] flux Total flux of the Gaussian (normalized analytically, not using pixel
281  values)
282  @param[in] ellipse lsst.afw.geom.ellipses.Ellipse holding the centroid and shape.
283  """
284  x, y = np.meshgrid(np.arange(bbox.getBeginX(), bbox.getEndX()),
285  np.arange(bbox.getBeginY(), bbox.getEndY()))
286  t = ellipse.getGridTransform()
287  xt = t[t.XX] * x + t[t.XY] * y + t[t.X]
288  yt = t[t.YX] * x + t[t.YY] * y + t[t.Y]
289  image = lsst.afw.image.ImageF(bbox)
290  image.getArray()[:, :] = np.exp(-0.5*(xt**2 + yt**2))*flux/(2.0*ellipse.getCore().getArea())
291  return image
292 
293  def __init__(self, bbox, threshold=10.0, exposure=None, **kwds):
294  """!
295  Initialize the dataset.
296 
297  @param[in] bbox Bounding box of the test image.
298  @param[in] threshold Threshold absolute value used to determine footprints for
299  simulated sources. This thresholding will be applied before noise is
300  actually added to images (or before the noise level is even known), so
301  this will necessarily produce somewhat artificial footprints.
302  @param[in] exposure lsst.afw.image.ExposureF test sources should be added to. Ownership should
303  be considered transferred from the caller to the TestDataset.
304  Must have a Gaussian Psf for truth catalog shapes to be exact.
305  @param[in] **kwds Keyword arguments forwarded to makeEmptyExposure if exposure is None.
306  """
307  if exposure is None:
308  exposure = self.makeEmptyExposure(bbox, **kwds)
309  self.threshold = lsst.afw.detection.Threshold(threshold, lsst.afw.detection.Threshold.VALUE)
310  self.exposure = exposure
311  self.psfShape = self.exposure.getPsf().computeShape()
312  self.schema = self.makeMinimalSchema()
314 
315  def _installFootprint(self, record, image):
316  """Create a Footprint for a simulated source and add it to its truth catalog record.
317  """
318  # Run detection on the single-source image
319  fpSet = lsst.afw.detection.FootprintSet(image, self.threshold)
320  # the call below to the FootprintSet ctor is actually a grow operation
321  fpSet = lsst.afw.detection.FootprintSet(fpSet, int(self.psfShape.getDeterminantRadius() + 1.0), True)
322  # Update the full exposure's mask plane to indicate the detection
323  fpSet.setMask(self.exposure.getMaskedImage().getMask(), "DETECTED")
324  # Attach the new footprint to the exposure
325  if len(fpSet.getFootprints()) > 1:
326  raise RuntimeError("Threshold value results in multiple Footprints for a single object")
327  if len(fpSet.getFootprints()) == 0:
328  raise RuntimeError("Threshold value results in zero Footprints for object")
329  record.setFootprint(fpSet.getFootprints()[0])
330 
331  def addSource(self, flux, centroid, shape=None):
332  """!
333  Add a source to the simulation
334 
335  @param[in] flux Total flux of the source to be added.
336  @param[in] centroid Position of the source to be added (lsst.afw.geom.Point2D).
337  @param[in] shape 2nd moments of the source before PSF convolution
338  (lsst.afw.geom.ellipses.Quadrupole). Note that the truth catalog
339  records post-convolution moments). If None, a point source
340  will be added.
341 
342  @return a truth catalog record and single-source image corresponding to the new source.
343  """
344  # Create and set the truth catalog fields
345  record = self.catalog.addNew()
346  record.set(self.keys["flux"], flux)
347  record.set(self.keys["centroid"], centroid)
348  if shape is None:
349  record.set(self.keys["isStar"], True)
350  fullShape = self.psfShape
351  else:
352  record.set(self.keys["isStar"], False)
353  fullShape = shape.convolve(self.psfShape)
354  record.set(self.keys["shape"], fullShape)
355  # Create an image containing just this source
356  image = self.drawGaussian(self.exposure.getBBox(), flux,
357  lsst.afw.geom.ellipses.Ellipse(fullShape, centroid))
358  # Generate a footprint for this source
359  self._installFootprint(record, image)
360  # Actually add the source to the full exposure
361  self.exposure.getMaskedImage().getImage().getArray()[:, :] += image.getArray()
362  return record, image
363 
364  def addBlend(self):
365  """!
366  Return a context manager that allows a blend of multiple sources to be added.
367 
368  Example:
369  @code
370  d = TestDataset(...)
371  with d.addBlend() as b:
372  b.addChild(flux1, centroid1)
373  b.addChild(flux2, centroid2, shape2)
374  @endcode
375 
376  Note that nothing stops you from creating overlapping sources just using the addSource() method,
377  but addBlend() is necesssary to create a parent object and deblended HeavyFootprints of the type
378  produced by the detection and deblending pipelines.
379  """
380  return BlendContext(self)
381 
382  def transform(self, wcs, **kwds):
383  """!
384  Create a copy of the dataset transformed to a new WCS, with new Psf and Calib.
385 
386  @param[in] wcs Wcs for the new dataset.
387  @param[in] **kwds Additional keyword arguments passed on to makeEmptyExposure. If not
388  specified, these revert to the defaults for makeEmptyExposure, not the
389  values in the current dataset.
390  """
391  bboxD = lsst.afw.geom.Box2D()
392  xyt = lsst.afw.image.XYTransformFromWcsPair(wcs, self.exposure.getWcs())
393  for corner in lsst.afw.geom.Box2D(self.exposure.getBBox()).getCorners():
394  bboxD.include(xyt.forwardTransform(lsst.afw.geom.Point2D(corner)))
395  bboxI = lsst.afw.geom.Box2I(bboxD)
396  result = TestDataset(bbox=bboxI, wcs=wcs, **kwds)
397  oldCalib = self.exposure.getCalib()
398  newCalib = result.exposure.getCalib()
399  oldPsfShape = self.exposure.getPsf().computeShape()
400  for record in self.catalog:
401  if record.get(self.keys["nChild"]):
402  raise NotImplementedError("Transforming blended sources in TestDatasets is not supported")
403  magnitude = oldCalib.getMagnitude(record.get(self.keys["flux"]))
404  newFlux = newCalib.getFlux(magnitude)
405  oldCentroid = record.get(self.keys["centroid"])
406  newCentroid = xyt.forwardTransform(oldCentroid)
407  if record.get(self.keys["isStar"]):
408  newDeconvolvedShape = None
409  else:
410  affine = xyt.linearizeForwardTransform(oldCentroid)
411  oldFullShape = record.get(self.keys["shape"])
412  oldDeconvolvedShape = lsst.afw.geom.ellipses.Quadrupole(
413  oldFullShape.getIxx() - oldPsfShape.getIxx(),
414  oldFullShape.getIyy() - oldPsfShape.getIyy(),
415  oldFullShape.getIxy() - oldPsfShape.getIxy(),
416  False
417  )
418  newDeconvolvedShape = oldDeconvolvedShape.transform(affine.getLinear())
419  result.addSource(newFlux, newCentroid, newDeconvolvedShape)
420  return result
421 
422  def realize(self, noise, schema, randomSeed=None):
423  """!
424  Create a simulated with noise and a simulated post-detection catalog with (Heavy)Footprints.
425 
426  @param[in] noise Standard deviation of noise to be added to the exposure. The noise will be
427  Gaussian and constant, appropriate for the sky-limited regime.
428  @param[in] schema Schema of the new catalog to be created. Must start with self.schema (i.e.
429  schema.contains(self.schema) must be True), but typically contains fields for
430  already-configured measurement algorithms as well.
431  @param[in] randomSeed Seed for the random number generator. If None, a seed is chosen automatically.
432 
433  @return a tuple of (exposure, catalog)
434  """
435  random_state = np.random.RandomState(randomSeed)
436  assert schema.contains(self.schema)
437  mapper = lsst.afw.table.SchemaMapper(self.schema)
438  mapper.addMinimalSchema(self.schema, True)
439  exposure = self.exposure.clone()
440  exposure.getMaskedImage().getVariance().getArray()[:, :] = noise**2
441  exposure.getMaskedImage().getImage().getArray()[:, :] \
442  += random_state.randn(exposure.getHeight(), exposure.getWidth())*noise
443  catalog = lsst.afw.table.SourceCatalog(schema)
444  catalog.extend(self.catalog, mapper=mapper)
445  # Loop over sources and generate new HeavyFootprints that divide up the noisy pixels, not the
446  # ideal no-noise pixels.
447  for record in catalog:
448  # parent objects have non-Heavy Footprints, which don't need to be updated after adding noise.
449  if record.getParent() == 0:
450  continue
451  # get flattened arrays that correspond to the no-noise and noisy parent images
452  parent = catalog.find(record.getParent())
453  footprint = parent.getFootprint()
454  parentFluxArrayNoNoise = np.zeros(footprint.getArea(), dtype=np.float32)
455  footprint.spans.flatten(parentFluxArrayNoNoise,
456  self.exposure.getMaskedImage().getImage().getArray(),
457  self.exposure.getXY0())
458  parentFluxArrayNoisy = np.zeros(footprint.getArea(), dtype=np.float32)
459  footprint.spans.flatten(parentFluxArrayNoisy,
460  exposure.getMaskedImage().getImage().getArray(),
461  exposure.getXY0())
462  oldHeavy = record.getFootprint()
463  fraction = (oldHeavy.getImageArray() / parentFluxArrayNoNoise)
464  # n.b. this isn't a copy ctor - it's a copy from a vanilla Footprint, so it doesn't copy
465  # the arrays we don't want to change, and hence we have to do that ourselves below.
466  newHeavy = lsst.afw.detection.HeavyFootprintF(oldHeavy)
467  newHeavy.getImageArray()[:] = parentFluxArrayNoisy*fraction
468  newHeavy.getMaskArray()[:] = oldHeavy.getMaskArray()
469  newHeavy.getVarianceArray()[:] = oldHeavy.getVarianceArray()
470  record.setFootprint(newHeavy)
471  return exposure, catalog
472 
473 
474 class AlgorithmTestCase(object):
475 
476  def makeSingleFrameMeasurementConfig(self, plugin=None, dependencies=()):
477  """Convenience function to create a Config instance for SingleFrameMeasurementTask
478 
479  The plugin and its dependencies will be the only plugins run, while the Centroid, Shape,
480  and ModelFlux slots will be set to the truth fields generated by the TestDataset class.
481  """
482  config = SingleFrameMeasurementTask.ConfigClass()
483  config.slots.centroid = "truth"
484  config.slots.shape = "truth"
485  config.slots.modelFlux = None
486  config.slots.apFlux = None
487  config.slots.psfFlux = None
488  config.slots.instFlux = None
489  config.slots.calibFlux = None
490  config.plugins.names = (plugin,) + tuple(dependencies)
491  return config
492 
493  def makeSingleFrameMeasurementTask(self, plugin=None, dependencies=(), config=None, schema=None,
494  algMetadata=None):
495  """Convenience function to create a SingleFrameMeasurementTask with a simple configuration.
496  """
497  if config is None:
498  if plugin is None:
499  raise ValueError("Either plugin or config argument must not be None")
500  config = self.makeSingleFrameMeasurementConfig(plugin=plugin, dependencies=dependencies)
501  if schema is None:
502  schema = TestDataset.makeMinimalSchema()
503  # Clear all aliases so only those defined by config are set.
504  schema.setAliasMap(None)
505  if algMetadata is None:
506  algMetadata = lsst.daf.base.PropertyList()
507  return SingleFrameMeasurementTask(schema=schema, algMetadata=algMetadata, config=config)
508 
509  def makeForcedMeasurementConfig(self, plugin=None, dependencies=()):
510  """Convenience function to create a Config instance for ForcedMeasurementTask
511 
512  In addition to the plugins specified in the plugin and dependencies arguments,
513  the TransformedCentroid and TransformedShape plugins will be run and used as the
514  Centroid and Shape slots; these simply transform the reference catalog centroid
515  and shape to the measurement coordinate system.
516  """
517  config = ForcedMeasurementTask.ConfigClass()
518  config.slots.centroid = "base_TransformedCentroid"
519  config.slots.shape = "base_TransformedShape"
520  config.slots.modelFlux = None
521  config.slots.apFlux = None
522  config.slots.psfFlux = None
523  config.slots.instFlux = None
524  config.plugins.names = (plugin,) + tuple(dependencies) + ("base_TransformedCentroid",
525  "base_TransformedShape")
526  return config
527 
528  def makeForcedMeasurementTask(self, plugin=None, dependencies=(), config=None, refSchema=None,
529  algMetadata=None):
530  """Convenience function to create a ForcedMeasurementTask with a simple configuration.
531  """
532  if config is None:
533  if plugin is None:
534  raise ValueError("Either plugin or config argument must not be None")
535  config = self.makeForcedMeasurementConfig(plugin=plugin, dependencies=dependencies)
536  if refSchema is None:
537  refSchema = TestDataset.makeMinimalSchema()
538  if algMetadata is None:
539  algMetadata = lsst.daf.base.PropertyList()
540  return ForcedMeasurementTask(refSchema=refSchema, algMetadata=algMetadata, config=config)
541 
542 
543 class TransformTestCase(object):
544  """!
545  Base class for testing measurement transformations.
546 
547  We test both that the transform itself operates successfully (fluxes are
548  converted to magnitudes, flags are propagated properly) and that the
549  transform is registered as the default for the appropriate measurement
550  algorithms.
551 
552  In the simple case of one-measurement-per-transformation, the developer
553  need not directly write any tests themselves: simply customizing the class
554  variables is all that is required. More complex measurements (e.g.
555  multiple aperture fluxes) require extra effort.
556  """
557  # The name used for the measurement algorithm; determines the names of the
558  # fields in the resulting catalog. This default should generally be fine,
559  # but subclasses can override if required.
560  name = "MeasurementTransformTest"
561 
562  # These should be customized by subclassing.
563  controlClass = None
564  algorithmClass = None
565  transformClass = None
566 
567  # Flags which may be set by the algorithm being tested. Can be customized
568  # in subclasses.
569  flagNames = ("flag",)
570 
571  # The plugin being tested should be registered under these names for
572  # single frame and forced measurement. Should be customized by
573  # subclassing.
574  singleFramePlugins = ()
575  forcedPlugins = ()
576 
577  def setUp(self):
579  self.calexp = TestDataset.makeEmptyExposure(bbox)
580  self._setupTransform()
581 
582  def tearDown(self):
583  del self.calexp
584  del self.inputCat
585  del self.mapper
586  del self.transform
587  del self.outputCat
588 
589  def _populateCatalog(self, baseNames):
590  records = []
591  for flagValue in (True, False):
592  records.append(self.inputCat.addNew())
593  for baseName in baseNames:
594  for flagName in self.flagNames:
595  if records[-1].schema.join(baseName, flagName) in records[-1].schema:
596  records[-1].set(records[-1].schema.join(baseName, flagName), flagValue)
597  self._setFieldsInRecords(records, baseName)
598 
599  def _checkOutput(self, baseNames):
600  for inSrc, outSrc in zip(self.inputCat, self.outputCat):
601  for baseName in baseNames:
602  self._compareFieldsInRecords(inSrc, outSrc, baseName)
603  for flagName in self.flagNames:
604  keyName = outSrc.schema.join(baseName, flagName)
605  if keyName in inSrc.schema:
606  self.assertEqual(outSrc.get(keyName), inSrc.get(keyName))
607  else:
608  self.assertFalse(keyName in outSrc.schema)
609 
610  def _runTransform(self, doExtend=True):
611  if doExtend:
612  self.outputCat.extend(self.inputCat, mapper=self.mapper)
613  self.transform(self.inputCat, self.outputCat, self.calexp.getWcs(), self.calexp.getCalib())
614 
615  def testTransform(self, baseNames=None):
616  """
617  Test the operation of the transformation on a catalog containing random data.
618 
619  We check that:
620 
621  * An appropriate exception is raised on an attempt to transform between catalogs with different
622  numbers of rows;
623  * Otherwise, all appropriate conversions are properly appled and that flags have been propagated.
624 
625  The `baseNames` argument requires some explanation. This should be an iterable of the leading parts of
626  the field names for each measurement; that is, everything that appears before `_flux`, `_flag`, etc.
627  In the simple case of a single measurement per plugin, this is simply equal to `self.name` (thus
628  measurements are stored as `self.name + "_flux"`, etc). More generally, the developer may specify
629  whatever iterable they require. For example, to handle multiple apertures, we could have
630  `(self.name + "_0", self.name + "_1", ...)`.
631 
632  @param[in] baseNames Iterable of the initial parts of measurement field names.
633  """
634  baseNames = baseNames or [self.name]
635  self._populateCatalog(baseNames)
636  self.assertRaises(lsst.pex.exceptions.LengthError, self._runTransform, False)
637  self._runTransform()
638  self._checkOutput(baseNames)
639 
640  def _checkRegisteredTransform(self, registry, name):
641  # If this is a Python-based transform, we can compare directly; if
642  # it's wrapped C++, we need to compare the wrapped class.
643  self.assertEqual(registry[name].PluginClass.getTransformClass(), self.transformClass)
644 
645  def testRegistration(self):
646  """
647  Test that the transformation is appropriately registered with the relevant measurement algorithms.
648  """
649  for pluginName in self.singleFramePlugins:
650  self._checkRegisteredTransform(lsst.meas.base.SingleFramePlugin.registry, pluginName)
651  for pluginName in self.forcedPlugins:
652  self._checkRegisteredTransform(lsst.meas.base.ForcedPlugin.registry, pluginName)
653 
654 
656 
657  def _setupTransform(self):
658  self.control = self.controlClass()
660  # Trick algorithms that depend on the slot centroid or alias into thinking they've been defined;
661  # it doesn't matter for this test since we won't actually use the plugins for anything besides
662  # defining the schema.
663  inputSchema.getAliasMap().set("slot_Centroid", "dummy")
664  inputSchema.getAliasMap().set("slot_Shape", "dummy")
665  self.algorithmClass(self.control, self.name, inputSchema)
666  inputSchema.getAliasMap().erase("slot_Centroid")
667  inputSchema.getAliasMap().erase("slot_Shape")
670  self.transform = self.transformClass(self.control, self.name, self.mapper)
671  self.outputCat = lsst.afw.table.BaseCatalog(self.mapper.getOutputSchema())
672 
673 
675 
676  def _setupTransform(self):
677  self.control = self.controlClass()
680  # Trick algorithms that depend on the slot centroid or alias into thinking they've been defined;
681  # it doesn't matter for this test since we won't actually use the plugins for anything besides
682  # defining the schema.
683  inputMapper.editOutputSchema().getAliasMap().set("slot_Centroid", "dummy")
684  inputMapper.editOutputSchema().getAliasMap().set("slot_Shape", "dummy")
685  self.algorithmClass(self.control, self.name, inputMapper, lsst.daf.base.PropertyList())
686  inputMapper.editOutputSchema().getAliasMap().erase("slot_Centroid")
687  inputMapper.editOutputSchema().getAliasMap().erase("slot_Shape")
688  self.inputCat = lsst.afw.table.SourceCatalog(inputMapper.getOutputSchema())
689  self.mapper = lsst.afw.table.SchemaMapper(inputMapper.getOutputSchema())
690  self.transform = self.transformClass(self.control, self.name, self.mapper)
691  self.outputCat = lsst.afw.table.BaseCatalog(self.mapper.getOutputSchema())
692 
693 
695 
696  def _setFieldsInRecords(self, records, name):
697  for record in records:
698  record[record.schema.join(name, 'flux')] = np.random.random()
699  record[record.schema.join(name, 'fluxSigma')] = np.random.random()
700 
701  # Negative fluxes should be converted to NaNs.
702  assert len(records) > 1
703  records[0][record.schema.join(name, 'flux')] = -1
704 
705  def _compareFieldsInRecords(self, inSrc, outSrc, name):
706  fluxName, fluxSigmaName = inSrc.schema.join(name, 'flux'), inSrc.schema.join(name, 'fluxSigma')
707  if inSrc[fluxName] > 0:
708  mag, magErr = self.calexp.getCalib().getMagnitude(inSrc[fluxName], inSrc[fluxSigmaName])
709  self.assertEqual(outSrc[outSrc.schema.join(name, 'mag')], mag)
710  self.assertEqual(outSrc[outSrc.schema.join(name, 'magErr')], magErr)
711  else:
712  self.assertTrue(np.isnan(outSrc[outSrc.schema.join(name, 'mag')]))
713  self.assertTrue(np.isnan(outSrc[outSrc.schema.join(name, 'magErr')]))
714 
715 
717 
718  def _setFieldsInRecords(self, records, name):
719  for record in records:
720  record[record.schema.join(name, 'x')] = np.random.random()
721  record[record.schema.join(name, 'y')] = np.random.random()
722  # Some algorithms set no errors; some set only sigma on x & y; some provide
723  # a full covariance matrix. Set only those which exist in the schema.
724  for fieldSuffix in ('xSigma', 'ySigma', 'x_y_Cov'):
725  fieldName = record.schema.join(name, fieldSuffix)
726  if fieldName in record.schema:
727  record[fieldName] = np.random.random()
728 
729  def _compareFieldsInRecords(self, inSrc, outSrc, name):
730  centroidResultKey = CentroidResultKey(inSrc.schema[self.name])
731  centroidResult = centroidResultKey.get(inSrc)
732 
733  coord = lsst.afw.table.CoordKey(outSrc.schema[self.name]).get(outSrc)
734  coordTruth = self.calexp.getWcs().pixelToSky(centroidResult.getCentroid())
735  self.assertEqual(coordTruth, coord)
736 
737  # If the centroid has an associated uncertainty matrix, the coordinate
738  # must have one too, and vice versa.
739  try:
740  coordErr = lsst.afw.table.CovarianceMatrix2fKey(outSrc.schema[self.name],
741  ["ra", "dec"]).get(outSrc)
743  self.assertFalse(centroidResultKey.getCentroidErr().isValid())
744  else:
745  transform = self.calexp.getWcs().linearizePixelToSky(coordTruth, lsst.afw.geom.radians)
746  coordErrTruth = np.dot(np.dot(transform.getLinear().getMatrix(),
747  centroidResult.getCentroidErr()),
748  transform.getLinear().getMatrix().transpose())
749  np.testing.assert_array_almost_equal(np.array(coordErrTruth), coordErr)
def _runTransform(self, doExtend=True)
Definition: tests.py:610
def realize(self, noise, schema, randomSeed=None)
Create a simulated with noise and a simulated post-detection catalog with (Heavy)Footprints.
Definition: tests.py:422
def makeSingleFrameMeasurementConfig(self, plugin=None, dependencies=())
Definition: tests.py:476
def transform(self, wcs, kwds)
Create a copy of the dataset transformed to a new WCS, with new Psf and Calib.
Definition: tests.py:382
def makePerturbedWcs(oldWcs, minScaleFactor=1.2, maxScaleFactor=1.5, minRotation=None, maxRotation=None, minRefShift=None, maxRefShift=None, minPixShift=2.0, maxPixShift=4.0, randomSeed=None)
Create a new undistorted TanWcs that is similar but not identical to another, with random scaling...
Definition: tests.py:184
std::shared_ptr< Wcs > makeWcs(coord::Coord const &crval, geom::Point2D const &crpix, double CD11, double CD12, double CD21, double CD22)
def __init__(self, bbox, threshold=10.0, exposure=None, kwds)
Initialize the dataset.
Definition: tests.py:293
def _checkOutput(self, baseNames)
Definition: tests.py:599
A subtask for measuring the properties of sources on a single exposure, using an existing "reference"...
A subtask for measuring the properties of sources on a single exposure.
Definition: sfm.py:152
def addBlend(self)
Return a context manager that allows a blend of multiple sources to be added.
Definition: tests.py:364
def makeEmptyExposure(bbox, wcs=None, crval=None, cdelt=None, psfSigma=2.0, psfDim=17, fluxMag0=1E12)
Create an Exposure, with a Calib, Wcs, and Psf, but no pixel values set.
Definition: tests.py:246
static LinearTransform makeScaling(double s)
def makeSingleFrameMeasurementTask(self, plugin=None, dependencies=(), config=None, schema=None, algMetadata=None)
Definition: tests.py:494
def __init__(self, owner)
Definition: tests.py:54
def testTransform(self, baseNames=None)
Definition: tests.py:615
static QuadrupoleKey addFields(Schema &schema, std::string const &name, std::string const &doc, CoordinateType coordType=CoordinateType::PIXEL)
def drawGaussian(bbox, flux, ellipse)
Create an image of an elliptical Gaussian.
Definition: tests.py:275
static LinearTransform makeRotation(Angle t)
def _installFootprint(self, record, image)
Definition: tests.py:315
def _checkRegisteredTransform(self, registry, name)
Definition: tests.py:640
def _populateCatalog(self, baseNames)
Definition: tests.py:589
def addSource(self, flux, centroid, shape=None)
Add a source to the simulation.
Definition: tests.py:331
static Schema makeMinimalSchema()
A simulated dataset consisting of a test image and an associated truth catalog.
Definition: tests.py:123
def makeForcedMeasurementTask(self, plugin=None, dependencies=(), config=None, refSchema=None, algMetadata=None)
Definition: tests.py:529
def addChild(self, flux, centroid, shape=None)
Add a child source to the blend, and return the truth catalog record that corresponds to it...
Definition: tests.py:64
Base class for testing measurement transformations.
Definition: tests.py:543
def __exit__(self, type_, value, tb)
Definition: tests.py:80
def makeForcedMeasurementConfig(self, plugin=None, dependencies=())
Definition: tests.py:509
A Python context manager used to add multiple overlapping sources along with a parent source that rep...
Definition: tests.py:45