Coverage for python/lsst/meas/extensions/scarlet/deblend.py : 15%

Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
1# This file is part of meas_extensions_scarlet.
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/>.
22from functools import partial
24import numpy as np
25import scarlet
26from scarlet.psf import PSF, gaussian
27from scarlet import PointSource, ExtendedSource, MultiComponentSource
29import lsst.log
30import lsst.pex.config as pexConfig
31import lsst.pipe.base as pipeBase
32from lsst.geom import Point2I, Box2I, Point2D
33import lsst.afw.math as afwMath
34import lsst.afw.geom as afwGeom
35import lsst.afw.geom.ellipses as afwEll
36import lsst.afw.image.utils
37import lsst.afw.image as afwImage
38import lsst.afw.detection as afwDet
39import lsst.afw.table as afwTable
41from .source import initSource, modelToHeavy
42from .blend import LsstBlend, checkBlendConvergence
43from .observation import LsstFrame, LsstObservation
45__all__ = ["deblend", "ScarletDeblendConfig", "ScarletDeblendTask"]
47logger = lsst.log.Log.getLogger("meas.deblender.deblend")
50def _getPsfFwhm(psf):
51 """Calculate the FWHM of the `psf`
52 """
53 return psf.computeShape().getDeterminantRadius() * 2.35
56def _estimateRMS(exposure, statsMask):
57 """Estimate the standard dev. of an image
59 Calculate the RMS of the `exposure`.
60 """
61 mi = exposure.getMaskedImage()
62 statsCtrl = afwMath.StatisticsControl()
63 statsCtrl.setAndMask(mi.getMask().getPlaneBitMask(statsMask))
64 stats = afwMath.makeStatistics(mi.variance, mi.mask, afwMath.STDEV | afwMath.MEAN, statsCtrl)
65 rms = np.sqrt(stats.getValue(afwMath.MEAN)**2 + stats.getValue(afwMath.STDEV)**2)
66 return rms
69def _computePsfImage(self, position=None):
70 """Get a multiband PSF image
71 The PSF Kernel Image is computed for each band
72 and combined into a (filter, y, x) array and stored
73 as `self._psfImage`.
74 The result is not cached, so if the same PSF is expected
75 to be used multiple times it is a good idea to store the
76 result in another variable.
77 Note: this is a temporary fix during the deblender sprint.
78 In the future this function will replace the current method
79 in `afw.MultibandExposure.computePsfImage` (DM-19789).
80 Parameters
81 ----------
82 position : `Point2D` or `tuple`
83 Coordinates to evaluate the PSF. If `position` is `None`
84 then `Psf.getAveragePosition()` is used.
85 Returns
86 -------
87 self._psfImage: array
88 The multiband PSF image.
89 """
90 psfs = []
91 # Make the coordinates into a Point2D (if necessary)
92 if not isinstance(position, Point2D) and position is not None:
93 position = Point2D(position[0], position[1])
95 for single in self.singles:
96 if position is None:
97 psf = single.getPsf().computeImage()
98 psfs.append(psf)
99 else:
100 psf = single.getPsf().computeImage(position)
101 psfs.append(psf)
102 left = np.min([psf.getBBox().getMinX() for psf in psfs])
103 bottom = np.min([psf.getBBox().getMinY() for psf in psfs])
104 right = np.max([psf.getBBox().getMaxX() for psf in psfs])
105 top = np.max([psf.getBBox().getMaxY() for psf in psfs])
106 bbox = Box2I(Point2I(left, bottom), Point2I(right, top))
107 psfs = [afwImage.utils.projectImage(psf, bbox) for psf in psfs]
108 psfImage = afwImage.MultibandImage.fromImages(self.filters, psfs)
109 return psfImage
112def getFootprintMask(footprint, mExposure):
113 """Mask pixels outside the footprint
115 Parameters
116 ----------
117 mExposure : `lsst.image.MultibandExposure`
118 - The multiband exposure containing the image,
119 mask, and variance data
120 footprint : `lsst.detection.Footprint`
121 - The footprint of the parent to deblend
123 Returns
124 -------
125 footprintMask : array
126 Boolean array with pixels not in the footprint set to one.
127 """
128 bbox = footprint.getBBox()
129 fpMask = afwImage.Mask(bbox)
130 footprint.spans.setMask(fpMask, 1)
131 fpMask = ~fpMask.getArray().astype(bool)
132 return fpMask
135def deblend(mExposure, footprint, config):
136 """Deblend a parent footprint
138 Parameters
139 ----------
140 mExposure : `lsst.image.MultibandExposure`
141 - The multiband exposure containing the image,
142 mask, and variance data
143 footprint : `lsst.detection.Footprint`
144 - The footprint of the parent to deblend
145 config : `ScarletDeblendConfig`
146 - Configuration of the deblending task
147 """
148 # Extract coordinates from each MultiColorPeak
149 bbox = footprint.getBBox()
151 # Create the data array from the masked images
152 images = mExposure.image[:, bbox].array
154 # Use the inverse variance as the weights
155 if config.useWeights:
156 weights = 1/mExposure.variance[:, bbox].array
157 else:
158 weights = np.ones_like(images)
160 # Mask out the pixels outside the footprint
161 mask = getFootprintMask(footprint, mExposure)
162 weights *= ~mask
164 psfs = _computePsfImage(mExposure, footprint.getCentroid()).array.astype(np.float32)
165 psfShape = (config.modelPsfSize, config.modelPsfSize)
166 model_psf = PSF(partial(gaussian, sigma=config.modelPsfSigma), shape=(None,)+psfShape)
168 frame = LsstFrame(images.shape, psfs=model_psf, channels=mExposure.filters)
169 observation = LsstObservation(images, psfs=psfs, weights=weights, channels=mExposure.filters)
170 observation.match(frame)
172 assert(config.sourceModel in ["single", "double", "point", "fit"])
174 # Only deblend sources that can be initialized
175 sources = []
176 skipped = []
177 for k, center in enumerate(footprint.peaks):
178 if config.sourceModel == "single":
179 components = 1
180 elif config.sourceModel == "double":
181 components = 2
182 elif config.sourceModel == "point":
183 components = 0
184 elif config.sourceModel == "fit":
185 # It is likely in the future that there will be some heuristic
186 # used to determine what type of model to use for each source,
187 # but that has not yet been implemented (see DM-22551)
188 raise NotImplementedError("sourceModel 'fit' has not been implemented yet")
189 else:
190 raise ValueError("Unrecognized sourceModel")
192 source = initSource(frame=frame, peak=center, observation=observation, bbox=bbox,
193 symmetric=config.symmetric, monotonic=config.monotonic,
194 thresh=config.morphThresh, components=components,
195 edgeDistance=config.edgeDistance, shifting=False)
196 if source is not None:
197 sources.append(source)
198 else:
199 skipped.append(k)
201 blend = LsstBlend(sources, observation)
202 blend.fit(max_iter=config.maxIter, e_rel=config.relativeError)
204 return blend, skipped
207class ScarletDeblendConfig(pexConfig.Config):
208 """MultibandDeblendConfig
210 Configuration for the multiband deblender.
211 The parameters are organized by the parameter types, which are
212 - Stopping Criteria: Used to determine if the fit has converged
213 - Position Fitting Criteria: Used to fit the positions of the peaks
214 - Constraints: Used to apply constraints to the peaks and their components
215 - Other: Parameters that don't fit into the above categories
216 """
217 # Stopping Criteria
218 maxIter = pexConfig.Field(dtype=int, default=300,
219 doc=("Maximum number of iterations to deblend a single parent"))
220 relativeError = pexConfig.Field(dtype=float, default=1e-4,
221 doc=("Change in the loss function between"
222 "iterations to exit fitter"))
224 # Blend Configuration options
225 edgeDistance = pexConfig.Field(dtype=int, default=1,
226 doc="All sources with flux within `edgeDistance` from the edge "
227 "will be considered edge sources.")
229 # Constraints
230 morphThresh = pexConfig.Field(dtype=float, default=1,
231 doc="Fraction of background RMS a pixel must have"
232 "to be included in the initial morphology")
233 monotonic = pexConfig.Field(dtype=bool, default=True, doc="Make models monotonic")
234 symmetric = pexConfig.Field(dtype=bool, default=False, doc="Make models symmetric")
236 # Other scarlet paremeters
237 useWeights = pexConfig.Field(
238 dtype=bool, default=True,
239 doc=("Whether or not use use inverse variance weighting."
240 "If `useWeights` is `False` then flat weights are used"))
241 modelPsfSize = pexConfig.Field(
242 dtype=int, default=11,
243 doc="Model PSF side length in pixels")
244 modelPsfSigma = pexConfig.Field(
245 dtype=float, default=0.8,
246 doc="Define sigma for the model frame PSF")
247 saveTemplates = pexConfig.Field(
248 dtype=bool, default=True,
249 doc="Whether or not to save the SEDs and templates")
250 processSingles = pexConfig.Field(
251 dtype=bool, default=False,
252 doc="Whether or not to process isolated sources in the deblender")
253 sourceModel = pexConfig.Field(
254 dtype=str, default="single",
255 doc=("How to determine which model to use for sources, from\n"
256 "- 'single': use a single component for all sources\n"
257 "- 'double': use a bulge disk model for all sources\n"
258 "- 'point: use a point-source model for all sources\n"
259 "- 'fit: use a PSF fitting model to determine the number of components (not yet implemented)")
260 )
262 # Mask-plane restrictions
263 badMask = pexConfig.ListField(
264 dtype=str, default=["BAD", "CR", "NO_DATA", "SAT", "SUSPECT"],
265 doc="Whether or not to process isolated sources in the deblender")
266 statsMask = pexConfig.ListField(dtype=str, default=["SAT", "INTRP", "NO_DATA"],
267 doc="Mask planes to ignore when performing statistics")
268 maskLimits = pexConfig.DictField(
269 keytype=str,
270 itemtype=float,
271 default={},
272 doc=("Mask planes with the corresponding limit on the fraction of masked pixels. "
273 "Sources violating this limit will not be deblended."),
274 )
276 # Size restrictions
277 maxNumberOfPeaks = pexConfig.Field(
278 dtype=int, default=0,
279 doc=("Only deblend the brightest maxNumberOfPeaks peaks in the parent"
280 " (<= 0: unlimited)"))
281 maxFootprintArea = pexConfig.Field(
282 dtype=int, default=1000000,
283 doc=("Maximum area for footprints before they are ignored as large; "
284 "non-positive means no threshold applied"))
285 maxFootprintSize = pexConfig.Field(
286 dtype=int, default=0,
287 doc=("Maximum linear dimension for footprints before they are ignored "
288 "as large; non-positive means no threshold applied"))
289 minFootprintAxisRatio = pexConfig.Field(
290 dtype=float, default=0.0,
291 doc=("Minimum axis ratio for footprints before they are ignored "
292 "as large; non-positive means no threshold applied"))
294 # Failure modes
295 notDeblendedMask = pexConfig.Field(
296 dtype=str, default="NOT_DEBLENDED", optional=True,
297 doc="Mask name for footprints not deblended, or None")
298 catchFailures = pexConfig.Field(
299 dtype=bool, default=False,
300 doc=("If True, catch exceptions thrown by the deblender, log them, "
301 "and set a flag on the parent, instead of letting them propagate up"))
302 propagateAllPeaks = pexConfig.Field(dtype=bool, default=False,
303 doc=('Guarantee that all peaks produce a child source.'))
306class ScarletDeblendTask(pipeBase.Task):
307 """ScarletDeblendTask
309 Split blended sources into individual sources.
311 This task has no return value; it only modifies the SourceCatalog in-place.
312 """
313 ConfigClass = ScarletDeblendConfig
314 _DefaultName = "scarletDeblend"
316 def __init__(self, schema, peakSchema=None, **kwargs):
317 """Create the task, adding necessary fields to the given schema.
319 Parameters
320 ----------
321 schema : `lsst.afw.table.schema.schema.Schema`
322 Schema object for measurement fields; will be modified in-place.
323 peakSchema : `lsst.afw.table.schema.schema.Schema`
324 Schema of Footprint Peaks that will be passed to the deblender.
325 Any fields beyond the PeakTable minimal schema will be transferred
326 to the main source Schema. If None, no fields will be transferred
327 from the Peaks.
328 filters : list of str
329 Names of the filters used for the eposures. This is needed to store
330 the SED as a field
331 **kwargs
332 Passed to Task.__init__.
333 """
334 pipeBase.Task.__init__(self, **kwargs)
336 peakMinimalSchema = afwDet.PeakTable.makeMinimalSchema()
337 if peakSchema is None:
338 # In this case, the peakSchemaMapper will transfer nothing, but
339 # we'll still have one
340 # to simplify downstream code
341 self.peakSchemaMapper = afwTable.SchemaMapper(peakMinimalSchema, schema)
342 else:
343 self.peakSchemaMapper = afwTable.SchemaMapper(peakSchema, schema)
344 for item in peakSchema:
345 if item.key not in peakMinimalSchema:
346 self.peakSchemaMapper.addMapping(item.key, item.field)
347 # Because SchemaMapper makes a copy of the output schema
348 # you give its ctor, it isn't updating this Schema in
349 # place. That's probably a design flaw, but in the
350 # meantime, we'll keep that schema in sync with the
351 # peakSchemaMapper.getOutputSchema() manually, by adding
352 # the same fields to both.
353 schema.addField(item.field)
354 assert schema == self.peakSchemaMapper.getOutputSchema(), "Logic bug mapping schemas"
355 self._addSchemaKeys(schema)
356 self.schema = schema
358 def _addSchemaKeys(self, schema):
359 """Add deblender specific keys to the schema
360 """
361 self.runtimeKey = schema.addField('runtime', type=np.float32, doc='runtime in ms')
363 self.iterKey = schema.addField('iterations', type=np.int32, doc='iterations to converge')
365 self.nChildKey = schema.addField('deblend_nChild', type=np.int32,
366 doc='Number of children this object has (defaults to 0)')
367 self.psfKey = schema.addField('deblend_deblendedAsPsf', type='Flag',
368 doc='Deblender thought this source looked like a PSF')
369 self.tooManyPeaksKey = schema.addField('deblend_tooManyPeaks', type='Flag',
370 doc='Source had too many peaks; '
371 'only the brightest were included')
372 self.tooBigKey = schema.addField('deblend_parentTooBig', type='Flag',
373 doc='Parent footprint covered too many pixels')
374 self.maskedKey = schema.addField('deblend_masked', type='Flag',
375 doc='Parent footprint was predominantly masked')
376 self.sedNotConvergedKey = schema.addField('deblend_sedConvergenceFailed', type='Flag',
377 doc='scarlet sed optimization did not converge before'
378 'config.maxIter')
379 self.morphNotConvergedKey = schema.addField('deblend_morphConvergenceFailed', type='Flag',
380 doc='scarlet morph optimization did not converge before'
381 'config.maxIter')
382 self.blendConvergenceFailedFlagKey = schema.addField('deblend_blendConvergenceFailedFlag',
383 type='Flag',
384 doc='at least one source in the blend'
385 'failed to converge')
386 self.edgePixelsKey = schema.addField('deblend_edgePixels', type='Flag',
387 doc='Source had flux on the edge of the parent footprint')
388 self.deblendFailedKey = schema.addField('deblend_failed', type='Flag',
389 doc="Deblending failed on source")
391 self.deblendSkippedKey = schema.addField('deblend_skipped', type='Flag',
392 doc="Deblender skipped this source")
393 self.modelCenter = afwTable.Point2DKey.addFields(schema, name="deblend_peak_center",
394 doc="Center used to apply constraints in scarlet",
395 unit="pixel")
396 self.modelCenterFlux = schema.addField('deblend_peak_instFlux', type=float, units='count',
397 doc="The instFlux at the peak position of deblended mode")
398 self.modelTypeKey = schema.addField("deblend_modelType", type="String", size=20,
399 doc="The type of model used, for example "
400 "MultiComponentSource, ExtendedSource, PointSource")
401 self.edgeFluxFlagKey = schema.addField("deblend_edgeFluxFlag", type="Flag",
402 doc="Source has flux on the edge of the image")
403 self.scarletFluxKey = schema.addField("deblend_scarletFlux", type=np.float32,
404 doc="Flux measurement from scarlet")
405 # self.log.trace('Added keys to schema: %s', ", ".join(str(x) for x in
406 # (self.nChildKey, self.tooManyPeaksKey, self.tooBigKey))
407 # )
409 @pipeBase.timeMethod
410 def run(self, mExposure, mergedSources):
411 """Get the psf from each exposure and then run deblend().
413 Parameters
414 ----------
415 mExposure : `MultibandExposure`
416 The exposures should be co-added images of the same
417 shape and region of the sky.
418 mergedSources : `SourceCatalog`
419 The merged `SourceCatalog` that contains parent footprints
420 to (potentially) deblend.
422 Returns
423 -------
424 fluxCatalogs: dict or None
425 Keys are the names of the filters and the values are
426 `lsst.afw.table.source.source.SourceCatalog`'s.
427 These are the flux-conserved catalogs with heavy footprints with
428 the image data weighted by the multiband templates.
429 If `self.config.conserveFlux` is `False`, then this item will be
430 None
431 templateCatalogs: dict or None
432 Keys are the names of the filters and the values are
433 `lsst.afw.table.source.source.SourceCatalog`'s.
434 These are catalogs with heavy footprints that are the templates
435 created by the multiband templates.
436 If `self.config.saveTemplates` is `False`, then this item will be
437 None
438 """
439 return self.deblend(mExposure, mergedSources)
441 @pipeBase.timeMethod
442 def deblend(self, mExposure, sources):
443 """Deblend a data cube of multiband images
445 Parameters
446 ----------
447 mExposure : `MultibandExposure`
448 The exposures should be co-added images of the same
449 shape and region of the sky.
450 sources : `SourceCatalog`
451 The merged `SourceCatalog` that contains parent footprints
452 to (potentially) deblend.
454 Returns
455 -------
456 fluxCatalogs : dict or None
457 Keys are the names of the filters and the values are
458 `lsst.afw.table.source.source.SourceCatalog`'s.
459 These are the flux-conserved catalogs with heavy footprints with
460 the image data weighted by the multiband templates.
461 If `self.config.conserveFlux` is `False`, then this item will be
462 None
463 templateCatalogs : dict or None
464 Keys are the names of the filters and the values are
465 `lsst.afw.table.source.source.SourceCatalog`'s.
466 These are catalogs with heavy footprints that are the templates
467 created by the multiband templates.
468 If `self.config.saveTemplates` is `False`, then this item will be
469 None
470 """
471 import time
473 filters = mExposure.filters
474 self.log.info("Deblending {0} sources in {1} exposure bands".format(len(sources), len(mExposure)))
476 # Create the output catalogs
477 templateCatalogs = {}
478 # This must be returned but is not calculated right now, setting it to
479 # None to be consistent with doc string
480 fluxCatalogs = None
481 for f in filters:
482 _catalog = afwTable.SourceCatalog(sources.table.clone())
483 _catalog.extend(sources)
484 templateCatalogs[f] = _catalog
486 n0 = len(sources)
487 nparents = 0
488 for pk, src in enumerate(sources):
489 foot = src.getFootprint()
490 bbox = foot.getBBox()
491 logger.info("id: {0}".format(src["id"]))
492 peaks = foot.getPeaks()
494 # Since we use the first peak for the parent object, we should
495 # propagate its flags to the parent source.
496 src.assign(peaks[0], self.peakSchemaMapper)
498 # Block of Skipping conditions
499 if len(peaks) < 2 and not self.config.processSingles:
500 for f in filters:
501 templateCatalogs[f][pk].set(self.runtimeKey, 0)
502 continue
503 if self._isLargeFootprint(foot):
504 src.set(self.tooBigKey, True)
505 self._skipParent(src, mExposure.mask)
506 self.log.trace('Parent %i: skipping large footprint', int(src.getId()))
507 continue
508 if self._isMasked(foot, mExposure):
509 src.set(self.maskedKey, True)
510 mask = np.bitwise_or.reduce(mExposure.mask[:, bbox].array, axis=0)
511 mask = afwImage.MaskX(mask, xy0=bbox.getMin())
512 self._skipParent(src, mask)
513 self.log.trace('Parent %i: skipping masked footprint', int(src.getId()))
514 continue
515 if len(peaks) > self.config.maxNumberOfPeaks:
516 src.set(self.tooManyPeaksKey, True)
517 msg = 'Parent {0}: Too many peaks, using the first {1} peaks'
518 self.log.trace(msg.format(int(src.getId()), self.config.maxNumberOfPeaks))
520 nparents += 1
521 self.log.trace('Parent %i: deblending %i peaks', int(src.getId()), len(peaks))
522 # Run the deblender
523 try:
524 t0 = time.time()
525 # Build the parameter lists with the same ordering
526 blend, skipped = deblend(mExposure, foot, self.config)
527 tf = time.time()
528 runtime = (tf-t0)*1000
529 src.set(self.deblendFailedKey, False)
530 src.set(self.runtimeKey, runtime)
531 converged = checkBlendConvergence(blend, self.config.relativeError)
532 src.set(self.blendConvergenceFailedFlagKey, converged)
533 sources = [src for src in blend.sources]
534 # Re-insert place holders for skipped sources
535 # to propagate them in the catalog so
536 # that the peaks stay consistent
537 for k in skipped:
538 sources.insert(k, None)
539 except Exception as e:
540 if self.config.catchFailures:
541 self.log.warn("Unable to deblend source %d: %s" % (src.getId(), e))
542 src.set(self.deblendFailedKey, True)
543 src.set(self.runtimeKey, 0)
544 import traceback
545 traceback.print_exc()
546 continue
547 else:
548 raise
550 # Add the merged source as a parent in the catalog for each band
551 templateParents = {}
552 parentId = src.getId()
553 for f in filters:
554 templateParents[f] = templateCatalogs[f][pk]
555 templateParents[f].set(self.runtimeKey, runtime)
556 templateParents[f].set(self.iterKey, len(blend.loss))
558 # Add each source to the catalogs in each band
559 templateSpans = {f: afwGeom.SpanSet() for f in filters}
560 nchild = 0
561 for k, source in enumerate(sources):
562 # Skip any sources with no flux or that scarlet skipped because
563 # it could not initialize
564 if k in skipped:
565 if not self.config.propagateAllPeaks:
566 # We don't care
567 continue
568 # We need to preserve the peak: make sure we have enough
569 # info to create a minimal child src
570 msg = "Peak at {0} failed deblending. Using minimal default info for child."
571 self.log.trace(msg.format(src.getFootprint().peaks[k]))
572 # copy the full footprint and strip out extra peaks
573 foot = afwDet.Footprint(src.getFootprint())
574 peakList = foot.getPeaks()
575 peakList.clear()
576 peakList.append(src.peaks[k])
577 zeroMimg = afwImage.MaskedImageF(foot.getBBox())
578 heavy = afwDet.makeHeavyFootprint(foot, zeroMimg)
579 models = afwDet.MultibandFootprint(mExposure.filters, [heavy]*len(mExposure.filters))
580 else:
581 src.set(self.deblendSkippedKey, False)
582 models = modelToHeavy(source, filters, xy0=bbox.getMin(),
583 observation=blend.observations[0])
584 # TODO: We should eventually write the morphology and SED to
585 # the catalog
586 # morph = source.morphToHeavy(xy0=bbox.getMin())
587 # sed = source.sed / source.sed.sum()
589 flux = scarlet.measure.flux(source)
590 for fidx, f in enumerate(filters):
591 if len(models[f].getPeaks()) != 1:
592 err = "Heavy footprint should have a single peak, got {0}"
593 raise ValueError(err.format(len(models[f].peaks)))
594 cat = templateCatalogs[f]
595 child = self._addChild(parentId, cat, models[f], source, converged,
596 xy0=bbox.getMin(), flux=flux[fidx])
597 if parentId == 0:
598 child.setId(src.getId())
599 child.set(self.runtimeKey, runtime)
600 else:
601 templateSpans[f] = templateSpans[f].union(models[f].getSpans())
602 nchild += 1
604 # Child footprints may extend beyond the full extent of their
605 # parent's which results in a failure of the replace-by-noise code
606 # to reinstate these pixels to their original values. The
607 # following updates the parent footprint in-place to ensure it
608 # contains the full union of itself and all of its
609 # children's footprints.
610 for f in filters:
611 templateParents[f].set(self.nChildKey, nchild)
612 templateParents[f].getFootprint().setSpans(templateSpans[f])
614 K = len(list(templateCatalogs.values())[0])
615 self.log.info('Deblended: of %i sources, %i were deblended, creating %i children, total %i sources'
616 % (n0, nparents, K-n0, K))
617 return fluxCatalogs, templateCatalogs
619 def _isLargeFootprint(self, footprint):
620 """Returns whether a Footprint is large
622 'Large' is defined by thresholds on the area, size and axis ratio.
623 These may be disabled independently by configuring them to be
624 non-positive.
626 This is principally intended to get rid of satellite streaks, which the
627 deblender or other downstream processing can have trouble dealing with
628 (e.g., multiple large HeavyFootprints can chew up memory).
629 """
630 if self.config.maxFootprintArea > 0 and footprint.getArea() > self.config.maxFootprintArea:
631 return True
632 if self.config.maxFootprintSize > 0:
633 bbox = footprint.getBBox()
634 if max(bbox.getWidth(), bbox.getHeight()) > self.config.maxFootprintSize:
635 return True
636 if self.config.minFootprintAxisRatio > 0:
637 axes = afwEll.Axes(footprint.getShape())
638 if axes.getB() < self.config.minFootprintAxisRatio*axes.getA():
639 return True
640 return False
642 def _isMasked(self, footprint, mExposure):
643 """Returns whether the footprint violates the mask limits"""
644 bbox = footprint.getBBox()
645 mask = np.bitwise_or.reduce(mExposure.mask[:, bbox].array, axis=0)
646 size = float(footprint.getArea())
647 for maskName, limit in self.config.maskLimits.items():
648 maskVal = mExposure.mask.getPlaneBitMask(maskName)
649 _mask = afwImage.MaskX(mask & maskVal, xy0=bbox.getMin())
650 unmaskedSpan = footprint.spans.intersectNot(_mask) # spanset of unmasked pixels
651 if (size - unmaskedSpan.getArea())/size > limit:
652 return True
653 return False
655 def _skipParent(self, source, masks):
656 """Indicate that the parent source is not being deblended
658 We set the appropriate flags and masks for each exposure.
660 Parameters
661 ----------
662 source : `lsst.afw.table.source.source.SourceRecord`
663 The source to flag as skipped
664 masks : list of `lsst.afw.image.MaskX`
665 The mask in each band to update with the non-detection
666 """
667 fp = source.getFootprint()
668 source.set(self.deblendSkippedKey, True)
669 source.set(self.nChildKey, len(fp.getPeaks())) # It would have this many if we deblended them all
670 if self.config.notDeblendedMask:
671 for mask in masks:
672 mask.addMaskPlane(self.config.notDeblendedMask)
673 fp.spans.setMask(mask, mask.getPlaneBitMask(self.config.notDeblendedMask))
675 def _addChild(self, parentId, sources, heavy, scarlet_source, blend_converged, xy0, flux):
676 """Add a child to a catalog
678 This creates a new child in the source catalog,
679 assigning it a parent id, adding a footprint,
680 and setting all appropriate flags based on the
681 deblender result.
682 """
683 assert len(heavy.getPeaks()) == 1
684 src = sources.addNew()
685 src.assign(heavy.getPeaks()[0], self.peakSchemaMapper)
686 src.setParent(parentId)
687 src.setFootprint(heavy)
688 src.set(self.psfKey, False)
689 src.set(self.runtimeKey, 0)
690 src.set(self.blendConvergenceFailedFlagKey, not blend_converged)
691 if isinstance(scarlet_source, ExtendedSource):
692 cy, cx = scarlet_source.pixel_center
693 morph = scarlet_source.morph
694 elif isinstance(scarlet_source, MultiComponentSource):
695 cy, cx = scarlet_source.components[0].pixel_center
696 morph = scarlet_source.components[0].morph
697 elif isinstance(scarlet_source, PointSource):
698 cy, cx = scarlet_source.parameters[1]
699 morph = scarlet_source.morph
700 else:
701 msg = "Did not recognize source type of `{0}`, could not write coordinates or center flux. "
702 msg += "Add `{0}` to meas_extensions_scarlet to properly persist this information."
703 logger.warning(msg.format(type(scarlet_source)))
704 return src
705 xmin, ymin = xy0
706 src.set(self.modelCenter, Point2D(cx+xmin, cy+ymin))
707 cy = np.max([np.min([int(np.round(cy)), morph.shape[0]-1]), 0])
708 cx = np.max([np.min([int(np.round(cx)), morph.shape[1]-1]), 0])
709 src.set(self.modelCenterFlux, morph[cy, cx])
710 src.set(self.modelTypeKey, scarlet_source.__class__.__name__)
711 src.set(self.edgeFluxFlagKey, scarlet_source.isEdge)
712 # Include the source flux in the model space in the catalog.
713 # This uses the narrower model PSF, which ensures that all sources
714 # not located on an edge have all of their flux included in the
715 # measurement.
716 src.set(self.scarletFluxKey, flux)
717 return src