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