LSSTApplications  10.0+286,10.0+36,10.0+46,10.0-2-g4f67435,10.1+152,10.1+37,11.0,11.0+1,11.0-1-g47edd16,11.0-1-g60db491,11.0-1-g7418c06,11.0-2-g04d2804,11.0-2-g68503cd,11.0-2-g818369d,11.0-2-gb8b8ce7
LSSTDataManagementBasePackage
Public Member Functions | Static Public Member Functions | Public Attributes | Private Member Functions | List of all members
lsst.meas.base.tests.TestDataset Class Reference

A simulated dataset consisting of a test image and an associated truth catalog. More...

Inheritance diagram for lsst.meas.base.tests.TestDataset:

Public Member Functions

def makeMinimalSchema
 
def __init__
 Initialize the dataset. More...
 
def addSource
 Add a source to the simulation. More...
 
def addBlend
 Return a context manager that allows a blend of multiple sources to be added. More...
 
def transform
 Create a copy of the dataset transformed to a new WCS, with new Psf and Calib. More...
 
def realize
 Create a simulated with noise and a simulated post-detection catalog with (Heavy)Footprints. More...
 

Static Public Member Functions

def makePerturbedWcs
 Create a new undistorted TanWcs that is similar but not identical to another, with random scaling, rotation, and offset (in both pixel position and reference position). More...
 
def makeEmptyExposure
 Create an Exposure, with a Calib, Wcs, and Psf, but no pixel values set. More...
 
def drawGaussian
 Create an image of an elliptical Gaussian. More...
 

Public Attributes

 threshold
 
 exposure
 
 psfShape
 
 schema
 
 catalog
 

Private Member Functions

def _installFootprint
 

Detailed Description

A simulated dataset consisting of a test image and an associated truth catalog.

TestDataset creates an idealized image made of pure Gaussians (including a Gaussian PSF), with simple noise and idealized Footprints/HeavyFootprints that simulated the outputs of detection and deblending. Multiple noise realizations can be created from the same underlying sources, allowing uncertainty estimates to be verified via Monte Carlo.

Typical usage:

2 dataset = TestDataset(bbox)
3 dataset.addSource(flux=1E5, centroid=lsst.afw.geom.Point2D(25, 26))
4 dataset.addSource(flux=2E5, centroid=lsst.afw.geom.Point2D(75, 24),
6 with dataset.addBlend() as family:
7  family.addChild(flux=2E5, centroid=lsst.afw.geom.Point2D(50, 72))
8  family.addChild(flux=1.5E5, centroid=lsst.afw.geom.Point2D(51, 74))
9 exposure, catalog = dataset.realize(noise=100.0, schema=TestDataset.makeMinimalSchema())

Definition at line 121 of file tests.py.

Constructor & Destructor Documentation

def lsst.meas.base.tests.TestDataset.__init__ (   self,
  bbox,
  threshold = 10.0,
  exposure = None,
  kwds 
)

Initialize the dataset.

Parameters
[in]bboxBounding box of the test image.
[in]thresholdThreshold absolute value used to determine footprints for simulated sources. This thresholding will be applied before noise is actually added to images (or before the noise level is even known), so this will necessarily produce somewhat artificial footprints.
[in]exposurelsst.afw.image.ExposureF test sources should be added to. Ownership should be considered transferred from the caller to the TestDataset. Must have a Gaussian Psf for truth catalog shapes to be exact.
[in]**kwdsKeyword arguments forwarded to makeEmptyExposure if exposure is None.

Definition at line 278 of file tests.py.

279  def __init__(self, bbox, threshold=10.0, exposure=None, **kwds):
280  """!
281  Initialize the dataset.
282 
283  @param[in] bbox Bounding box of the test image.
284  @param[in] threshold Threshold absolute value used to determine footprints for
285  simulated sources. This thresholding will be applied before noise is
286  actually added to images (or before the noise level is even known), so
287  this will necessarily produce somewhat artificial footprints.
288  @param[in] exposure lsst.afw.image.ExposureF test sources should be added to. Ownership should
289  be considered transferred from the caller to the TestDataset.
290  Must have a Gaussian Psf for truth catalog shapes to be exact.
291  @param[in] **kwds Keyword arguments forwarded to makeEmptyExposure if exposure is None.
292  """
293  if exposure is None:
294  exposure = self.makeEmptyExposure(bbox, **kwds)
295  self.threshold = lsst.afw.detection.Threshold(threshold, lsst.afw.detection.Threshold.VALUE)
296  self.exposure = exposure
297  self.psfShape = self.exposure.getPsf().computeShape()
298  self.schema = self.makeMinimalSchema()
def makeEmptyExposure
Create an Exposure, with a Calib, Wcs, and Psf, but no pixel values set.
Definition: tests.py:231
A Threshold is used to pass a threshold value to detection algorithms.
Definition: Threshold.h:44
Custom catalog class for record/table subclasses that are guaranteed to have an ID, and should generally be sorted by that ID.
Definition: fwd.h:55
def __init__
Initialize the dataset.
Definition: tests.py:278

Member Function Documentation

def lsst.meas.base.tests.TestDataset._installFootprint (   self,
  record,
  image 
)
private
Create a Footprint for a simulated source and add it to its truth catalog record.

Definition at line 300 of file tests.py.

301  def _installFootprint(self, record, image):
302  """Create a Footprint for a simulated source and add it to its truth catalog record.
303  """
304  # Run detection on the single-source image
305  fpSet = lsst.afw.detection.FootprintSet(image, self.threshold)
306  # the call below to the FootprintSet ctor is actually a grow operation
307  fpSet = lsst.afw.detection.FootprintSet(fpSet, int(self.psfShape.getDeterminantRadius() + 1.0), True)
308  # Update the full exposure's mask plane to indicate the detection
309  fpSet.setMask(self.exposure.getMaskedImage().getMask(), "DETECTED")
310  # Attach the new footprint to the exposure
311  if len(fpSet.getFootprints()) > 1:
312  raise RuntimeError("Threshold value results in multiple Footprints for a single object")
313  if len(fpSet.getFootprints()) == 0:
314  raise RuntimeError("Threshold value results in zero Footprints for object")
315  record.setFootprint(fpSet.getFootprints()[0])
A set of Footprints, associated with a MaskedImage.
Definition: FootprintSet.h:53
def lsst.meas.base.tests.TestDataset.addBlend (   self)

Return a context manager that allows a blend of multiple sources to be added.

Example:

1 d = TestDataset(...)
2 with d.addBlend() as b:
3  b.addChild(flux1, centroid1)
4  b.addChild(flux2, centroid2, shape2)

Note that nothing stops you from creating overlapping sources just using the addSource() method, but addBlend() is necesssary to create a parent object and deblended HeavyFootprints of the type produced by the detection and deblending pipelines.

Definition at line 349 of file tests.py.

350  def addBlend(self):
351  """!
352  Return a context manager that allows a blend of multiple sources to be added.
353 
354  Example:
355  @code
356  d = TestDataset(...)
357  with d.addBlend() as b:
358  b.addChild(flux1, centroid1)
359  b.addChild(flux2, centroid2, shape2)
360  @endcode
361 
362  Note that nothing stops you from creating overlapping sources just using the addSource() method,
363  but addBlend() is necesssary to create a parent object and deblended HeavyFootprints of the type
364  produced by the detection and deblending pipelines.
365  """
366  return BlendContext(self)
def addBlend
Return a context manager that allows a blend of multiple sources to be added.
Definition: tests.py:349
A Python context manager used to add multiple overlapping sources along with a parent source that rep...
Definition: tests.py:43
def lsst.meas.base.tests.TestDataset.addSource (   self,
  flux,
  centroid,
  shape = None 
)

Add a source to the simulation.

Parameters
[in]fluxTotal flux of the source to be added.
[in]centroidPosition of the source to be added (lsst.afw.geom.Point2D).
[in]shape2nd moments of the source before PSF convolution (lsst.afw.geom.ellipses.Quadrupole). Note that the truth catalog records post-convolution moments). If None, a point source will be added.
Returns
a truth catalog record and single-source image corresponding to the new source.

Definition at line 316 of file tests.py.

317  def addSource(self, flux, centroid, shape=None):
318  """!
319  Add a source to the simulation
320 
321  @param[in] flux Total flux of the source to be added.
322  @param[in] centroid Position of the source to be added (lsst.afw.geom.Point2D).
323  @param[in] shape 2nd moments of the source before PSF convolution
324  (lsst.afw.geom.ellipses.Quadrupole). Note that the truth catalog
325  records post-convolution moments). If None, a point source
326  will be added.
327 
328  @return a truth catalog record and single-source image corresponding to the new source.
329  """
330  # Create and set the truth catalog fields
331  record = self.catalog.addNew()
332  record.set(self.keys["flux"], flux)
333  record.set(self.keys["centroid"], centroid)
334  if shape is None:
335  record.set(self.keys["isStar"], True)
336  fullShape = self.psfShape
337  else:
338  record.set(self.keys["isStar"], False)
339  fullShape = shape.convolve(self.psfShape)
340  record.set(self.keys["shape"], fullShape)
341  # Create an image containing just this source
342  image = self.drawGaussian(self.exposure.getBBox(), flux,
343  lsst.afw.geom.ellipses.Ellipse(fullShape, centroid))
344  # Generate a footprint for this source
345  self._installFootprint(record, image)
346  # Actually add the source to the full exposure
347  self.exposure.getMaskedImage().getImage().getArray()[:,:] += image.getArray()
348  return record, image
def addSource
Add a source to the simulation.
Definition: tests.py:316
def drawGaussian
Create an image of an elliptical Gaussian.
Definition: tests.py:260
An ellipse defined by an arbitrary BaseCore and a center point.
Definition: Ellipse.h:50
def lsst.meas.base.tests.TestDataset.drawGaussian (   bbox,
  flux,
  ellipse 
)
static

Create an image of an elliptical Gaussian.

Parameters
[in,out]bboxBounding box of image to create.
[in]fluxTotal flux of the Gaussian (normalized analytically, not using pixel values)
[in]ellipselsst.afw.geom.ellipses.Ellipse holding the centroid and shape.

Definition at line 260 of file tests.py.

261  def drawGaussian(bbox, flux, ellipse):
262  """!
263  Create an image of an elliptical Gaussian.
264 
265  @param[in,out] bbox Bounding box of image to create.
266  @param[in] flux Total flux of the Gaussian (normalized analytically, not using pixel
267  values)
268  @param[in] ellipse lsst.afw.geom.ellipses.Ellipse holding the centroid and shape.
269  """
270  x, y = numpy.meshgrid(numpy.arange(bbox.getBeginX(), bbox.getEndX()),
271  numpy.arange(bbox.getBeginY(), bbox.getEndY()))
272  t = ellipse.getGridTransform()
273  xt = t[t.XX] * x + t[t.XY] * y + t[t.X]
274  yt = t[t.YX] * x + t[t.YY] * y + t[t.Y]
275  image = lsst.afw.image.ImageF(bbox)
276  image.getArray()[:,:] = numpy.exp(-0.5*(xt**2 + yt**2))*flux/(2.0*ellipse.getCore().getArea())
277  return image
def drawGaussian
Create an image of an elliptical Gaussian.
Definition: tests.py:260
def lsst.meas.base.tests.TestDataset.makeEmptyExposure (   bbox,
  wcs = None,
  crval = None,
  cdelt = None,
  psfSigma = 2.0,
  psfDim = 17,
  fluxMag0 = 1E12 
)
static

Create an Exposure, with a Calib, Wcs, and Psf, but no pixel values set.

Parameters
[in]bboxBounding box of the image (image coordinates) as returned by makeCatalog.
[in]wcsNew Wcs for the exposure (created from crval and cdelt if None).
[in]crvalafw.coord.Coord: center of the TAN WCS attached to the image.
[in]cdeltafw.geom.Angle: pixel scale of the image
[in]psfSigmaRadius (sigma) of the Gaussian PSF attached to the image
[in]psfDimWidth and height of the image's Gaussian PSF attached to the image
[in]fluxMag0Flux at magnitude zero (in e-) used to set the Calib of the exposure.

Definition at line 231 of file tests.py.

232  def makeEmptyExposure(bbox, wcs=None, crval=None, cdelt=None, psfSigma=2.0, psfDim=17, fluxMag0=1E12):
233  """!
234  Create an Exposure, with a Calib, Wcs, and Psf, but no pixel values set.
235 
236  @param[in] bbox Bounding box of the image (image coordinates) as returned by makeCatalog.
237  @param[in] wcs New Wcs for the exposure (created from crval and cdelt if None).
238  @param[in] crval afw.coord.Coord: center of the TAN WCS attached to the image.
239  @param[in] cdelt afw.geom.Angle: pixel scale of the image
240  @param[in] psfSigma Radius (sigma) of the Gaussian PSF attached to the image
241  @param[in] psfDim Width and height of the image's Gaussian PSF attached to the image
242  @param[in] fluxMag0 Flux at magnitude zero (in e-) used to set the Calib of the exposure.
243  """
244  if wcs is None:
245  if crval is None:
246  crval = lsst.afw.coord.IcrsCoord(45.0*lsst.afw.geom.degrees, 45.0*lsst.afw.geom.degrees)
247  if cdelt is None:
248  cdelt = 0.2*lsst.afw.geom.arcseconds
249  crpix = lsst.afw.geom.Box2D(bbox).getCenter()
250  wcs = lsst.afw.image.makeWcs(crval, crpix, cdelt.asDegrees(), 0.0, 0.0, cdelt.asDegrees())
251  exposure = lsst.afw.image.ExposureF(bbox)
252  psf = lsst.afw.detection.GaussianPsf(psfDim, psfDim, psfSigma)
253  calib = lsst.afw.image.Calib()
254  calib.setFluxMag0(fluxMag0)
255  exposure.setWcs(wcs)
256  exposure.setPsf(psf)
257  exposure.setCalib(calib)
258  return exposure
def makeEmptyExposure
Create an Exposure, with a Calib, Wcs, and Psf, but no pixel values set.
Definition: tests.py:231
Wcs::Ptr makeWcs(coord::Coord const &crval, geom::Point2D const &crpix, double CD11, double CD12, double CD21, double CD22)
Create a Wcs object from crval, crpix, CD, using CD elements (useful from python) ...
Definition: makeWcs.cc:141
A floating-point coordinate rectangle geometry.
Definition: Box.h:271
A class to handle Icrs coordinates (inherits from Coord)
Definition: Coord.h:157
A circularly symmetric Gaussian Psf class with no spatial variation, intended mostly for testing purp...
Definition: GaussianPsf.h:41
def lsst.meas.base.tests.TestDataset.makeMinimalSchema (   cls)
Return the minimal schema needed to hold truth catalog fields.

When TestDataset.realize() is called, the schema must include at least these fields.
Usually it will include additional fields for measurement algorithm outputs, allowing
the same catalog to be used for both truth values (the fields from the minimal schema)
and the measurements.

Definition at line 145 of file tests.py.

146  def makeMinimalSchema(cls):
147  """Return the minimal schema needed to hold truth catalog fields.
148 
149  When TestDataset.realize() is called, the schema must include at least these fields.
150  Usually it will include additional fields for measurement algorithm outputs, allowing
151  the same catalog to be used for both truth values (the fields from the minimal schema)
152  and the measurements.
153  """
154  if not hasattr(cls, "_schema"):
156  cls.keys = {}
157  cls.keys["parent"] = schema.find("parent").key
158  cls.keys["nChild"] = schema.addField("deblend_nchild", type=int)
159  cls.keys["flux"] = schema.addField("truth_flux", type=float, doc="true flux", units="dn")
160  cls.keys["centroid"] = lsst.afw.table.Point2DKey.addFields(
161  schema, "truth", "true simulated centroid", "pixels"
162  )
163  cls.keys["shape"] = lsst.afw.table.QuadrupoleKey.addFields(
164  schema, "truth", "true shape after PSF convolution", lsst.afw.table.PIXEL
165  )
166  cls.keys["isStar"] = schema.addField("truth_isStar", type="Flag",
167  doc="set if the object is a star")
168  schema.getAliasMap().set("slot_Shape", "truth")
169  schema.getAliasMap().set("slot_Centroid", "truth")
170  schema.getAliasMap().set("slot_ModelFlux", "truth")
171  schema.getCitizen().markPersistent()
172  cls._schema = schema
173  schema = lsst.afw.table.Schema(cls._schema)
174  schema.disconnectAliases()
175  return schema
Defines the fields and offsets for a table.
Definition: Schema.h:46
static QuadrupoleKey addFields(Schema &schema, std::string const &name, std::string const &doc, CoordinateType coordType=CoordinateType::PIXEL)
static Schema makeMinimalSchema()
Return a minimal schema for Source tables and records.
Definition: Source.h:242
def lsst.meas.base.tests.TestDataset.makePerturbedWcs (   oldWcs,
  minScaleFactor = 1.2,
  maxScaleFactor = 1.5,
  minRotation = None,
  maxRotation = None,
  minRefShift = None,
  maxRefShift = None,
  minPixShift = 2.0,
  maxPixShift = 4.0 
)
static

Create a new undistorted TanWcs that is similar but not identical to another, with random scaling, rotation, and offset (in both pixel position and reference position).

The maximum and minimum arguments are interpreted as absolute values for a split range that covers both positive and negative values (as this method is used in testing, it is typically most important to avoid perturbations near zero). Scale factors are treated somewhat differently: the actual scale factor is chosen between minScaleFactor and maxScaleFactor OR (1/maxScaleFactor) and (1/minScaleFactor).

The default range for rotation is 30-60 degrees, and the default range for reference shift is 0.5-1.0 arcseconds (these cannot be safely included directly as default values because Angle objects are mutable).

Definition at line 180 of file tests.py.

181  minPixShift=2.0, maxPixShift=4.0):
182  """!
183  Create a new undistorted TanWcs that is similar but not identical to another, with random
184  scaling, rotation, and offset (in both pixel position and reference position).
185 
186  The maximum and minimum arguments are interpreted as absolute values for a split
187  range that covers both positive and negative values (as this method is used
188  in testing, it is typically most important to avoid perturbations near zero).
189  Scale factors are treated somewhat differently: the actual scale factor is chosen between
190  minScaleFactor and maxScaleFactor OR (1/maxScaleFactor) and (1/minScaleFactor).
191 
192  The default range for rotation is 30-60 degrees, and the default range for reference shift
193  is 0.5-1.0 arcseconds (these cannot be safely included directly as default values because Angle
194  objects are mutable).
195  """
196  if minRotation is None: minRotation = 30.0*lsst.afw.geom.degrees
197  if maxRotation is None: maxRotation = 60.0*lsst.afw.geom.degrees
198  if minRefShift is None: minRefShift = 0.5*lsst.afw.geom.arcseconds
199  if maxRefShift is None: maxRefShift = 1.0*lsst.afw.geom.arcseconds
200  def splitRandom(min1, max1, min2=None, max2=None):
201  if min2 is None: min2 = -max1
202  if max2 is None: max2 = -min1
203  if numpy.random.uniform() > 0.5:
204  return float(numpy.random.uniform(min1, max1))
205  else:
206  return float(numpy.random.uniform(min2, max2))
207  # Generate random perturbations
208  scaleFactor = splitRandom(minScaleFactor, maxScaleFactor, 1.0/maxScaleFactor, 1.0/minScaleFactor)
209  rotation = splitRandom(minRotation.asRadians(), maxRotation.asRadians())*lsst.afw.geom.radians
210  refShiftRa = splitRandom(minRefShift.asRadians(), maxRefShift.asRadians())*lsst.afw.geom.radians
211  refShiftDec = splitRandom(minRefShift.asRadians(), maxRefShift.asRadians())*lsst.afw.geom.radians
212  pixShiftX = splitRandom(minPixShift, maxPixShift)
213  pixShiftY = splitRandom(minPixShift, maxPixShift)
214  # Compute new CD matrix
215  oldTransform = lsst.afw.geom.LinearTransform(oldWcs.getCDMatrix())
216  rTransform = lsst.afw.geom.LinearTransform.makeRotation(rotation)
217  sTransform = lsst.afw.geom.LinearTransform.makeScaling(scaleFactor)
218  newTransform = oldTransform*rTransform*sTransform
219  matrix = newTransform.getMatrix()
220  # Compute new coordinate reference pixel (CRVAL)
221  oldSkyOrigin = oldWcs.getSkyOrigin().toIcrs()
222  newSkyOrigin = lsst.afw.coord.IcrsCoord(oldSkyOrigin.getRa() + refShiftRa,
223  oldSkyOrigin.getDec() + refShiftDec)
224  # Compute new pixel reference pixel (CRPIX)
225  oldPixOrigin = oldWcs.getPixelOrigin()
226  newPixOrigin = lsst.afw.geom.Point2D(oldPixOrigin.getX() + pixShiftX,
227  oldPixOrigin.getY() + pixShiftY)
228  return lsst.afw.image.makeWcs(newSkyOrigin, newPixOrigin,
229  matrix[0,0], matrix[0,1], matrix[1,0], matrix[1,1])
static LinearTransform makeScaling(double s)
A 2D linear coordinate transformation.
static LinearTransform makeRotation(Angle t)
Wcs::Ptr makeWcs(coord::Coord const &crval, geom::Point2D const &crpix, double CD11, double CD12, double CD21, double CD22)
Create a Wcs object from crval, crpix, CD, using CD elements (useful from python) ...
Definition: makeWcs.cc:141
A class to handle Icrs coordinates (inherits from Coord)
Definition: Coord.h:157
def lsst.meas.base.tests.TestDataset.realize (   self,
  noise,
  schema 
)

Create a simulated with noise and a simulated post-detection catalog with (Heavy)Footprints.

Parameters
[in]noiseStandard deviation of noise to be added to the exposure. The noise will be Gaussian and constant, appropriate for the sky-limited regime.
[in]schemaSchema of the new catalog to be created. Must start with self.schema (i.e. schema.contains(self.schema) must be True), but typically contains fields for already-configured measurement algorithms as well.
Returns
a tuple of (exposure, catalog)

Definition at line 407 of file tests.py.

408  def realize(self, noise, schema):
409  """!
410  Create a simulated with noise and a simulated post-detection catalog with (Heavy)Footprints.
411 
412  @param[in] noise Standard deviation of noise to be added to the exposure. The noise will be
413  Gaussian and constant, appropriate for the sky-limited regime.
414  @param[in] schema Schema of the new catalog to be created. Must start with self.schema (i.e.
415  schema.contains(self.schema) must be True), but typically contains fields for
416  already-configured measurement algorithms as well.
417 
418  @return a tuple of (exposure, catalog)
419  """
420  assert schema.contains(self.schema)
421  mapper = lsst.afw.table.SchemaMapper(self.schema)
422  mapper.addMinimalSchema(self.schema, True)
423  exposure = self.exposure.clone()
424  exposure.getMaskedImage().getVariance().getArray()[:,:] = noise**2
425  exposure.getMaskedImage().getImage().getArray()[:,:] \
426  += numpy.random.randn(exposure.getHeight(), exposure.getWidth())*noise
427  catalog = lsst.afw.table.SourceCatalog(schema)
428  catalog.extend(self.catalog, mapper=mapper)
429  # Loop over sources and generate new HeavyFootprints that divide up the noisy pixels, not the
430  # ideal no-noise pixels.
431  for record in catalog:
432  # parent objects have non-Heavy Footprints, which don't need to be updated after adding noise.
433  if record.getParent() == 0: continue
434  # get flattened arrays that correspond to the no-noise and noisy parent images
435  parent = catalog.find(record.getParent())
436  footprint = parent.getFootprint()
437  parentFluxArrayNoNoise = numpy.zeros(footprint.getArea(), dtype=numpy.float32)
439  footprint,
440  self.exposure.getMaskedImage().getImage().getArray(),
441  parentFluxArrayNoNoise,
442  self.exposure.getXY0()
443  )
444  parentFluxArrayNoisy = numpy.zeros(footprint.getArea(), dtype=numpy.float32)
446  footprint,
447  exposure.getMaskedImage().getImage().getArray(),
448  parentFluxArrayNoisy,
449  exposure.getXY0()
450  )
451  oldHeavy = lsst.afw.detection.HeavyFootprintF.cast(record.getFootprint())
452  fraction = (oldHeavy.getImageArray() / parentFluxArrayNoNoise)
453  # n.b. this isn't a copy ctor - it's a copy from a vanilla Footprint, so it doesn't copy
454  # the arrays we don't want to change, and hence we have to do that ourselves below.
455  newHeavy = lsst.afw.detection.HeavyFootprintF(oldHeavy)
456  newHeavy.getImageArray()[:] = parentFluxArrayNoisy*fraction
457  newHeavy.getMaskArray()[:] = oldHeavy.getMaskArray()
458  newHeavy.getVarianceArray()[:] = oldHeavy.getVarianceArray()
459  record.setFootprint(newHeavy)
460  return exposure, catalog
461 
A mapping between the keys of two Schemas, used to copy data between them.
Definition: SchemaMapper.h:19
Custom catalog class for record/table subclasses that are guaranteed to have an ID, and should generally be sorted by that ID.
Definition: fwd.h:55
ndarray::Array< typename boost::remove_const< T >::type, N-1, N-1 > flattenArray(Footprint const &fp, ndarray::Array< T, N, C > const &src, lsst::afw::geom::Point2I const &xy0=lsst::afw::geom::Point2I())
Flatten the first two dimensions of an array Use this footprint to map 2-D points in the source to 1-...
def realize
Create a simulated with noise and a simulated post-detection catalog with (Heavy)Footprints.
Definition: tests.py:407
def lsst.meas.base.tests.TestDataset.transform (   self,
  wcs,
  kwds 
)

Create a copy of the dataset transformed to a new WCS, with new Psf and Calib.

Parameters
[in]wcsWcs for the new dataset.
[in]**kwdsAdditional keyword arguments passed on to makeEmptyExposure. If not specified, these revert to the defaults for makeEmptyExposure, not the values in the current dataset.

Definition at line 367 of file tests.py.

368  def transform(self, wcs, **kwds):
369  """!
370  Create a copy of the dataset transformed to a new WCS, with new Psf and Calib.
371 
372  @param[in] wcs Wcs for the new dataset.
373  @param[in] **kwds Additional keyword arguments passed on to makeEmptyExposure. If not
374  specified, these revert to the defaults for makeEmptyExposure, not the
375  values in the current dataset.
376  """
377  bboxD = lsst.afw.geom.Box2D()
378  xyt = lsst.afw.image.XYTransformFromWcsPair(wcs, self.exposure.getWcs())
379  for corner in lsst.afw.geom.Box2D(self.exposure.getBBox()).getCorners():
380  bboxD.include(xyt.forwardTransform(lsst.afw.geom.Point2D(corner)))
381  bboxI = lsst.afw.geom.Box2I(bboxD)
382  result = TestDataset(bbox=bboxI, wcs=wcs, **kwds)
383  oldCalib = self.exposure.getCalib()
384  newCalib = result.exposure.getCalib()
385  oldPsfShape = self.exposure.getPsf().computeShape()
386  for record in self.catalog:
387  if record.get(self.keys["nChild"]):
388  raise NotImplementedError("Transforming blended sources in TestDatasets is not supported")
389  magnitude = oldCalib.getMagnitude(record.get(self.keys["flux"]))
390  newFlux = newCalib.getFlux(magnitude)
391  oldCentroid = record.get(self.keys["centroid"])
392  newCentroid = xyt.forwardTransform(oldCentroid)
393  if record.get(self.keys["isStar"]):
394  newDeconvolvedShape = None
395  else:
396  affine = xyt.linearizeForwardTransform(oldCentroid)
397  oldFullShape = record.get(self.keys["shape"])
398  oldDeconvolvedShape = lsst.afw.geom.ellipses.Quadrupole(
399  oldFullShape.getIxx() - oldPsfShape.getIxx(),
400  oldFullShape.getIyy() - oldPsfShape.getIyy(),
401  oldFullShape.getIxy() - oldPsfShape.getIxy(),
402  False
403  )
404  newDeconvolvedShape = oldDeconvolvedShape.transform(affine.getLinear())
405  result.addSource(newFlux, newCentroid, newDeconvolvedShape)
406  return result
An ellipse core with quadrupole moments as parameters.
Definition: Quadrupole.h:45
XYTransformFromWcsPair: An XYTransform obtained by putting two Wcs objects &quot;back to back&quot;...
Definition: Wcs.h:432
An integer coordinate rectangle.
Definition: Box.h:53
def transform
Create a copy of the dataset transformed to a new WCS, with new Psf and Calib.
Definition: tests.py:367
A simulated dataset consisting of a test image and an associated truth catalog.
Definition: tests.py:121
A floating-point coordinate rectangle geometry.
Definition: Box.h:271

Member Data Documentation

lsst.meas.base.tests.TestDataset.catalog

Definition at line 298 of file tests.py.

lsst.meas.base.tests.TestDataset.exposure

Definition at line 295 of file tests.py.

lsst.meas.base.tests.TestDataset.psfShape

Definition at line 296 of file tests.py.

lsst.meas.base.tests.TestDataset.schema

Definition at line 297 of file tests.py.

lsst.meas.base.tests.TestDataset.threshold

Definition at line 294 of file tests.py.


The documentation for this class was generated from the following file: