242 def run(self, *, masked_image: MaskedImageF, catalog: SourceCatalog, seed: int) -> Struct:
243 """Compute raw moments, select stars, and fit a shapelet expansion to
249 Masked image to measure on. Must be background-subtracted.
251 Catalog of detections to extract footprints from and fill output
252 columns of. Its schema must be a superset of ``self.schema``.
254 A random-number generator seed, used for the robust covariance
259 `lsst.pipe.base.Struct`
260 A struct of results containing:
262 - ``shapelet`` (`lsst.shapelet.ShapeletFunction`): A
263 Gauss-Laguerre (polar shaplet) expansion of the PSF.
264 - ``radial`` (`list` [`float`]): the purely radial coefficients
265 of the shapelet expansion.
266 - all attributes returned by the `select_stars` method.
270 star_moments = [moments[star_id]
for star_id
in result.star_ids]
271 result.shapelet = SpanSetMoments.fit_shapelets(
274 self.config.shapelet_order,
275 self.config.shapelet_scale_factor,
277 result.shapelet.getEllipse().setCore(result.mean_shape)
278 result.shapelet.changeBasisType(LAGUERRE)
279 result.radial = result.shapelet.getCoefficients()[
280 [computeOffset(i)
for i
in range(0, self.config.shapelet_order + 1, 2)]
282 self.log.info(
"Rough PSF shapelet radial terms: %s.", result.radial)
286 self, *, masked_image: MaskedImageF, catalog: SourceCatalog
287 ) -> dict[int, SpanSetMoments]:
288 """Compute the unweighted moments of the footprints in a catalog.
293 Masked image to measure on. Must be background-subtracted.
295 Catalog of detections to extract footprints from and fill output
296 columns of. Its schema must be a superset of ``self.schema``.
300 `dict` [`int`, `SpanSetMoments`]
301 Objects used to construct and hold the unweighted moments and the
302 pixel region used to computed them, keyed by source ID.
304 bitmask = masked_image.mask.getPlaneBitMask(self.config.bad_mask_planes)
305 all_moments: dict[int, SpanSetMoments] = {}
306 for record
in catalog:
307 if record.getParent() != 0:
309 self.log.debug(
"Skipping child source %s", record.getId())
311 footprint_spans = record.getFootprint().getSpans()
312 if footprint_spans.getArea() > self.config.max_footprint_area:
315 "Skipping source %s with footprint area %d > %d.",
317 footprint_spans.getArea(),
318 self.config.max_footprint_area,
321 moments = SpanSetMoments.compute(
322 record.getFootprint().getSpans(),
323 masked_image=masked_image,
325 bad_pixel_max_fraction=self.config.bad_pixel_max_fraction,
326 bad_pixel_exclusion_radius=self.config.bad_pixel_exclusion_radius,
332 record.set(self.
_flag_key, moments.any_flags_set)
333 if not moments.any_flags_set:
334 all_moments[record.getId()] = moments
337 "Successfully measured raw moments for %d of %d sources.",
346 """Select probable stars from the distribution of second moments.
351 Catalog of detections to extract footprints from and fill output
352 columns of. Its schema must be a superset of ``self.schema``.
354 A random-number generator seed, used for the robust covariance
359 `lsst.pipe.base.Struct`
360 A struct of results containing:
362 - ``star_ids`` (`numpy.ndarray`): the source IDs that are expected
364 - ``mean_shape`` (`lsst.afw.geom.ellipses.BaseCore`): the mean of
365 the shape distribution.
366 - ``shape_covariance`` (`numpy.ndarray`): the covariance of the
367 distribution of shapes; a 3x3 matrix. This uses the same
368 parameterization of the shapes as ``mean_shape``.
369 - ``radius_cut`` (`float`): the indended radius cut (i.e. the
370 mode of the radius distribution multipled by the
371 ``radius_factor`` configuration option).
372 - ``radius_kde`` (`sklearn.neighbors.KernelDensity`): kernel
373 density estimator on the radius distribution, used to determine
377 indices = np.arange(len(catalog), dtype=int)[np.logical_not(catalog[self.
_flag_key])]
381 threshold=self.config.min_snr,
382 min_count=self.config.min_n_stars,
383 max_count=len(catalog),
389 radii = np.zeros(indices.size, dtype=np.float64)
390 for n, index
in enumerate(indices):
391 record = catalog[index]
393 radii[n] = shape.getTraceRadius()
395 radius_cut = self.config.max_radius_factor * radius_mode
399 threshold=radius_cut,
400 min_count=self.config.min_n_stars,
401 max_count=len(catalog),
406 shape_data = np.zeros((len(indices), 3), dtype=np.float64)
407 for n, index
in enumerate(indices):
408 record = catalog[index]
411 if self.config.logarithmic_shapes:
413 shape_data[n, :] = shape.getParameterVector()
414 shape_dist = MinCovDet(random_state=seed).fit(shape_data)
415 m_distances = shape_dist.mahalanobis(shape_data)
419 threshold=self.config.max_shape_distance,
420 min_count=self.config.min_n_stars,
421 max_count=self.config.max_n_stars,
422 name=
"Mahalanobis distance",
426 for index
in indices:
428 star_ids = catalog[
"id"][indices]
429 if self.config.logarithmic_shapes:
435 mean_shape=mean_shape,
436 shape_covariance=shape_dist.covariance_,
437 radius_cut=radius_cut,
438 radius_kde=radius_kde,
442 self, figure: matplotlib.figure.Figure, *, catalog: SourceCatalog, results: Struct
444 """Create plots of the shape distribution space used to select stars.
449 Matplotlib figure to plot to.
451 Catalog of sources with columns populated by the `run` method (at
452 least through the `select_stars` step).
454 Result struct returned by `run` or `select_stars`.
456 from matplotlib.lines
import Line2D
458 shape_data = np.zeros((len(catalog), 3), dtype=np.float64)
459 radii = np.zeros(len(catalog), dtype=np.float64)
460 for n, record
in enumerate(catalog):
464 if self.config.logarithmic_shapes:
466 shape_data[n, :] = shape.getParameterVector()
467 radii[n] = shape.getTraceRadius()
469 candidate_mask = np.logical_and(catalog[self.
_candidate_key], np.logical_not(used_mask))
470 measured_mask = np.logical_and(
474 axes = figure.subplot_mosaic(
476 [
"radius",
"radius",
"radius"],
478 [
"scatter01",
"hist1",
"."],
479 [
"scatter02",
"scatter12",
"hist2"],
481 gridspec_kw=dict(bottom=0.2, top=1.0),
483 axes[
"scatter01"].sharex(axes[
"hist0"])
484 axes[
"scatter02"].sharex(axes[
"hist0"])
485 axes[
"scatter12"].sharex(axes[
"hist1"])
486 axes[
"scatter12"].sharey(axes[
"scatter02"])
487 for tk
in itertools.chain(
488 axes[
"hist0"].get_xticklabels(),
489 axes[
"scatter01"].get_xticklabels(),
490 axes[
"hist1"].get_xticklabels(),
491 axes[
"scatter12"].get_yticklabels(),
493 tk.set_visible(
False)
495 for ax
in [axes[
"hist0"], axes[
"hist1"], axes[
"hist2"]]:
496 ax.yaxis.set_label_position(
"right")
497 ax.yaxis.tick_right()
499 names = [
"Ixx",
"Iyy",
"Ixy"]
if not self.config.logarithmic_shapes
else [
"𝜂1",
"𝜂2",
"ln(r)"]
500 axes[
"scatter02"].set_xlabel(names[0])
501 axes[
"scatter12"].set_xlabel(names[1])
502 axes[
"hist2"].set_xlabel(names[2])
503 axes[
"scatter01"].set_ylabel(names[1])
504 axes[
"scatter02"].set_ylabel(names[2])
506 mu = results.mean_shape.getParameterVector()
507 sigma = np.diagonal(results.shape_covariance) ** 0.5
508 lower_bounds = [max(mu[i] - 3 * sigma[i], min(shape_data[:, i]))
for i
in range(3)]
509 upper_bounds = [min(mu[i] + 3 * sigma[i], max(shape_data[:, i]))
for i
in range(3)]
510 grids = [np.linspace(lower_bounds[i], upper_bounds[i], 50)
for i
in range(3)]
511 axes[
"radius"].axvline(results.radius_cut, color=
"k")
512 for color, mask, alpha
in [
513 (
"grey", measured_mask, 0.5),
514 (
"blue", candidate_mask, 0.75),
515 (
"green", used_mask, 1.0),
523 range=(radii.min(), 15.0),
527 axes[f
"hist{i}"].hist(
530 range=(lower_bounds[i], upper_bounds[i]),
537 if (ax := axes.get(f
"scatter{i}{j}"))
is not None:
547 axes[f
"hist{i}"].
plot(
548 grids[i], scipy.stats.norm.pdf(grids[i], loc=mu[i], scale=sigma[i]),
"k", alpha=0.5
550 axes[f
"hist{i}"].set_xlim(lower_bounds[i], upper_bounds[i])
552 if (ax := axes.get(f
"scatter{i}{j}"))
is not None:
554 results.shape_covariance[i, i],
555 results.shape_covariance[j, j],
556 results.shape_covariance[i, j],
558 for factor
in [1, 2, 3]:
569 ax.set_xlim(lower_bounds[i], upper_bounds[i])
570 ax.set_ylim(lower_bounds[j], upper_bounds[j])
573 Line2D([], [], color=
"green", alpha=1.0),
574 Line2D([], [], color=
"blue", alpha=0.75),
575 Line2D([], [], color=
"gray", alpha=0.5),
578 f
"RoughPsfShapelets_used ({used_mask.sum()})",
579 f
"RoughPsfShapelets_candidate & ~RoughPsfShapelets_used ({candidate_mask.sum()})",
580 f
"~RoughPsfShapelets_flag & ~RoughPsfShapelets_candidate ({measured_mask.sum()})",
588 figure: matplotlib.figure.Figure,
591 catalog: SourceCatalog,
594 stamp_size: float = 2.0,
596 """Create data/model/residual plots of stars and the shapelet model.
601 Matplotlib figure to plot to.
603 The image the stars were measured on.
605 Catalog of sources with columns populated by the `run` method .
607 Result struct returned by `run`.
609 Number of stars to include.
611 Stamp size in inches.
613 from matplotlib.colors
import Normalize
615 width = stamp_size * 3 + 1.5
616 figure.set_size_inches(w=width, h=stamp_size * n_stars)
617 axes = figure.subplot_mosaic(
619 [
"image_cbar", f
"d{star_id}", f
"m{star_id}", f
"r{star_id}",
"res_cbar"]
620 for star_id
in results.star_ids[:n_stars]
623 wspace=0.01, hspace=0.01, left=0.5 / width, right=1.0 - 0.5 / width, bottom=0.01, top=0.99
625 width_ratios=[0.25, stamp_size, stamp_size, stamp_size, 0.25],
627 for name, ax
in axes.items():
628 if not name.endswith(
"cbar"):
630 norm: Normalize |
None =
None
631 res_norm: Normalize |
None =
None
632 for star_id
in results.star_ids[:n_stars]:
633 record = catalog.find(star_id)
634 star_bbox = record.getFootprint().getBBox()
635 star_model = ImageD(star_bbox)
639 star_shapelet.setEllipse(star_ellipse)
640 star_shapelet.evaluate().addToImage(star_model)
642 norm = Normalize(vmin=star_model.array.min(), vmax=star_model.array.max())
643 star_image = image[star_bbox].clone()
645 self.
_draw_image(axes[f
"d{star_id}"], star_image, norm=norm, cmap=
"YlGnBu")
646 self.
_draw_ellipse(axes[f
"d{star_id}"], star_ellipse, fill=
False, edgecolor=
"blue", alpha=0.5)
647 image_plot = self.
_draw_image(axes[f
"m{star_id}"], star_model, norm=norm, cmap=
"YlGnBu")
648 self.
_draw_ellipse(axes[f
"m{star_id}"], star_ellipse, fill=
False, edgecolor=
"blue", alpha=0.5)
649 star_image -= star_model.convertF()
650 amax = np.abs(star_image.array).max()
652 res_norm = Normalize(vmin=-amax, vmax=amax)
653 res_plot = self.
_draw_image(axes[f
"r{star_id}"], star_image, norm=res_norm, cmap=
"RdBu")
654 self.
_draw_ellipse(axes[f
"r{star_id}"], star_ellipse, fill=
False, edgecolor=
"blue", alpha=0.5)
655 figure.colorbar(image_plot, cax=axes[
"image_cbar"], location=
"left")
656 figure.colorbar(res_plot, cax=axes[
"res_cbar"], location=
"right")
666 kind: Literal[
"<",
">"],
668 """Return the indices of an array that satisfy an inequality
669 and/or lower and upper bounds on the number of indices returned.
674 Array of values to threshold on.
676 Threshold value that selected elements must be above or below.
678 The minimum number of indices returned. When thresholding would
679 yield fewer than this number, the threshold is ignored. Note that
680 the number of indices may still be less than this if the size of
681 ``values`` is less than this.
683 The maximum number of indices returned.
685 Name of the quantity being thresholded, for log messages.
687 Whether the threshold is a upper bound (``<``) or lower bound
688 (``>``). This also sets how values are ranked when they are added
689 or dropped to satisfy the count constraints.
694 Indices into ``values``.
696 if min_count > len(values):
698 f
"Not enough sources ({len(values)}) for {name} cut that must yield at least {min_count}."
700 sorter = values.argsort()
701 n = np.searchsorted(values[sorter], threshold)
703 sorter = sorter[::-1]
707 "Applying a %s %s %f cut yields only %d sources; keeping the top %d (%s %s %f) instead.",
715 values[sorter[min_count - 1]],
720 "%d sources have %s %s %f; keeping only the top %d (%s %s %f) instead.",
728 values[sorter[max_count - 1]],
732 self.log.verbose(
"Keeping %d sources with %s %s %f.", n, name, kind, threshold)
736 """Find the first peak in a 1-d distribution of radii."""
737 kde = KernelDensity(bandwidth=self.config.radius_kde_bandwidth).fit(radii.reshape(-1, 1))
738 sorted_radii = radii.copy()
740 sorted_radii = sorted_radii[sorted_radii.searchsorted(self.config.radius_mode_min):]
741 scores = kde.score_samples(sorted_radii.reshape(-1, 1))
742 peaks, _ = scipy.signal.find_peaks(scores)
745 return sorted_radii[peaks.min()], kde
749 axes: matplotlib.axes.Axes,
752 x: float |
None =
None,
753 y: float |
None =
None,
756 ) -> matplotlib.patches.Ellipse:
757 from matplotlib.patches
import Ellipse
as EllipsePatch
761 x = ellipse.getCenter().getX()
763 y = ellipse.getCenter().getY()
764 ellipse = ellipse.getCore()
771 patch = EllipsePatch(
773 ellipse.getA() * 2 * scale,
774 ellipse.getB() * 2 * scale,
775 angle=ellipse.getTheta() * 180.0 / np.pi,
778 axes.add_patch(patch)