Namespace lsst::ip::diffim

namespace diffim

Functions

lsst::afw::math::KernelList makeDeltaFunctionBasisList(int width, int height)

Build a set of Delta Function basis kernels.

Note

Total number of basis functions is width*height

Parameters
  • width: Width of basis set (cols)

  • height: Height of basis set (rows)

Eigen::MatrixXd makeRegularizationMatrix(lsst::daf::base::PropertySet const &ps)

Build a regularization matrix for Delta function kernels.

Note

Calls either makeForwardDifferenceMatrix or makeCentralDifferenceMatrix based on the PropertySet config.

Parameters
  • ps: PropertySet dictating which type of matrix to make

Eigen::MatrixXd makeForwardDifferenceMatrix(int width, int height, std::vector<int> const &orders, float borderPenalty, bool fitForBackground)

Build a forward difference regularization matrix for Delta function kernels.

Parameters
  • width: Width of basis set you want to regularize

  • height: Height of basis set you want to regularize

  • orders: Which derivatives to penalize (1,2,3)

  • borderPenalty: Amount of penalty (if any) to apply to border pixels; > 0

  • fitForBackground: Fit for differential background?

Eigen::MatrixXd makeCentralDifferenceMatrix(int width, int height, int stencil, float borderPenalty, bool fitForBackground)

Build a central difference Laplacian regularization matrix for Delta function kernels.

Parameters
  • width: Width of basis set you want to regularize

  • height: Height of basis set you want to regularize

  • stencil: Which type of Laplacian approximation to use

  • borderPenalty: Amount of penalty (if any) to apply to border pixels; > 0

  • fitForBackground: Fit for differential background?

lsst::afw::math::KernelList renormalizeKernelList(lsst::afw::math::KernelList const &kernelListIn)

Renormalize a list of basis kernels.

Note

Renormalization means make Ksum_0 = 1.0, Ksum_i = 0.0, K_i.dot.K_i = 1.0

Note

Output list of shared pointers to FixedKernels

Note

Images are checked for their current kernel sum. If it is larger than std::numeric_limits<double>::epsilon(), the kernel is first divided by the kernel sum, giving it a kSum of 1.0, and then the first (normalized) component is subtracted from it, giving it a kSum of 0.0.

Parameters
  • kernelListIn: input list of basis kernels

lsst::afw::math::KernelList makeAlardLuptonBasisList(int halfWidth, int nGauss, std::vector<double> const &sigGauss, std::vector<int> const &degGauss)

Build a set of Alard/Lupton basis kernels.

Note

Should consider implementing as SeparableKernels for additional speed, but this will make the normalization a bit more complicated

Parameters
  • halfWidth: size is 2*N + 1

  • nGauss: number of gaussians

  • sigGauss: Widths of the Gaussian Kernels

  • degGauss: Local spatial variation of bases

template<typename PixelT, typename BackgroundT>
lsst::afw::image::MaskedImage<PixelT> convolveAndSubtract(lsst::afw::image::MaskedImage<PixelT> const &templateImage, lsst::afw::image::MaskedImage<PixelT> const &scienceMaskedImage, lsst::afw::math::Kernel const &convolutionKernel, BackgroundT background, bool invert = true)

Execute fundamental task of convolving template and subtracting it from science image.

Note

This version accepts a MaskedImage for the template

Parameters
  • templateImage: MaskedImage to apply convolutionKernel to

  • scienceMaskedImage: MaskedImage from which convolved templateImage is subtracted

  • convolutionKernel: Kernel to apply to templateImage

  • background: Background scalar or function to subtract after convolution

  • invert: Invert the output difference image

template<typename PixelT, typename BackgroundT>
lsst::afw::image::MaskedImage<PixelT> convolveAndSubtract(lsst::afw::image::Image<PixelT> const &templateImage, lsst::afw::image::MaskedImage<PixelT> const &scienceMaskedImage, lsst::afw::math::Kernel const &convolutionKernel, BackgroundT background, bool invert = true)

Execute fundamental task of convolving template and subtracting it from science image.

Note

This version accepts an Image for the template, and is thus faster during convolution

Parameters
  • templateImage: Image to apply convolutionKernel to

  • scienceMaskedImage: MaskedImage from which convolved templateImage is subtracted

  • convolutionKernel: Kernel to apply to templateImage

  • background: Background scalar or function to subtract after convolution

  • invert: Invert the output difference image

template<typename PixelT>
Eigen::MatrixXd imageToEigenMatrix(lsst::afw::image::Image<PixelT> const &img)

Turns a 2-d Image into a 2-d Eigen Matrix.

Parameters
  • img: Image whose pixel values are read into an Eigen::MatrixXd

Eigen::MatrixXi maskToEigenMatrix(lsst::afw::image::Mask<lsst::afw::image::MaskPixel> const &mask)
template<typename PixelT>
std::shared_ptr<KernelCandidate<PixelT>> makeKernelCandidate(float const xCenter, float const yCenter, std::shared_ptr<afw::image::MaskedImage<PixelT>> const &templateMaskedImage, std::shared_ptr<afw::image::MaskedImage<PixelT>> const &scienceMaskedImage, daf::base::PropertySet const &ps)

Return a KernelCandidate pointer of the right sort.

Parameters
  • xCenter: X-center of candidate

  • yCenter: Y-center of candidate

  • templateMaskedImage: Template subimage

  • scienceMaskedImage: Science image subimage

  • ps: PropertySet for creation of rating

template<typename PixelT>
std::shared_ptr<KernelCandidate<PixelT>> makeKernelCandidate(std::shared_ptr<afw::table::SourceRecord> const &source, std::shared_ptr<afw::image::MaskedImage<PixelT>> const &templateMaskedImage, std::shared_ptr<afw::image::MaskedImage<PixelT>> const &scienceMaskedImage, daf::base::PropertySet const &ps)

Return a KernelCandidate pointer of the right sort.

Parameters

class DipoleCentroidAlgorithm : public lsst::meas::base::SimpleAlgorithm
#include <DipoleAlgorithms.h>

Intermediate base class for algorithms that compute a centroid.

Subclassed by lsst::ip::diffim::NaiveDipoleCentroid

class DipoleFluxAlgorithm : public lsst::meas::base::SimpleAlgorithm
#include <DipoleAlgorithms.h>

Intermediate base class for algorithms that compute a flux.

Subclassed by lsst::ip::diffim::NaiveDipoleFlux, lsst::ip::diffim::PsfDipoleFlux

class DipoleFluxControl

Subclassed by lsst::ip::diffim::PsfDipoleFluxControl

template<typename MaskT>
class FindSetBits
#include <FindSetBits.h>

Class to accumulate Mask bits.

Note

Search through a Mask for any set bits.

template<typename PixelT>
class ImageStatistics
#include <ImageStatistics.h>

Class to calculate difference image statistics.

Note

Find mean and unbiased variance of pixel residuals in units of sqrt(variance)

template<typename _PixelT>
class KernelCandidate : public lsst::afw::math::SpatialCellImageCandidate
#include <KernelCandidate.h>

Class stored in SpatialCells for spatial Kernel fitting.

Note

KernelCandidate is a single Kernel derived around a source. We’ll assign them to sets of SpatialCells; these sets will then be used to fit a spatial model to the Kernel.

template<typename PixelT>
class KernelCandidateDetection
#include <KernelCandidateDetection.h>

Search through images for Footprints with no masked pixels.

Note

Runs detection on the template; searches through both images for masked pixels

Parameters
  • templateMaskedImage: MaskedImage that will be convolved with kernel

  • scienceMaskedImage: MaskedImage to subtract convolved template from

  • ps: PropertySet for operations; in particular object detection

class KernelSolution

Subclassed by lsst::ip::diffim::SpatialKernelSolution, lsst::ip::diffim::StaticKernelSolution< InputT >

class NaiveDipoleCentroid : public lsst::ip::diffim::DipoleCentroidAlgorithm
#include <DipoleAlgorithms.h>

Intermediate base class for algorithms that compute a centroid.

class PsfDipoleFlux : public lsst::ip::diffim::DipoleFluxAlgorithm
#include <DipoleAlgorithms.h>

Implementation of Psf dipole flux

class PsfDipoleFluxControl : public lsst::ip::diffim::DipoleFluxControl
#include <DipoleAlgorithms.h>

C++ control object for PSF dipole fluxes.

template<typename InputT>
class StaticKernelSolution : public lsst::ip::diffim::KernelSolution

Subclassed by lsst::ip::diffim::MaskedKernelSolution< InputT >, lsst::ip::diffim::RegularizedKernelSolution< InputT >

namespace detail

Functions

template<typename PixelT>
std::shared_ptr<AssessSpatialKernelVisitor<PixelT>> makeAssessSpatialKernelVisitor(std::shared_ptr<lsst::afw::math::LinearCombinationKernel> spatialKernel, lsst::afw::math::Kernel::SpatialFunctionPtr spatialBackground, lsst::daf::base::PropertySet const &ps)
template<typename PixelT>
std::shared_ptr<BuildSingleKernelVisitor<PixelT>> makeBuildSingleKernelVisitor(lsst::afw::math::KernelList const &basisList, lsst::daf::base::PropertySet const &ps)
template<typename PixelT>
std::shared_ptr<BuildSingleKernelVisitor<PixelT>> makeBuildSingleKernelVisitor(lsst::afw::math::KernelList const &basisList, lsst::daf::base::PropertySet const &ps, Eigen::MatrixXd const &hMat)
template<typename PixelT>
std::shared_ptr<BuildSpatialKernelVisitor<PixelT>> makeBuildSpatialKernelVisitor(lsst::afw::math::KernelList const &basisList, lsst::geom::Box2I const &regionBBox, lsst::daf::base::PropertySet const &ps)
template<typename PixelT>
std::shared_ptr<KernelPcaVisitor<PixelT>> makeKernelPcaVisitor(std::shared_ptr<KernelPca<typename KernelPcaVisitor<PixelT>::ImageT>> imagePca)
template<typename PixelT>
std::shared_ptr<KernelSumVisitor<PixelT>> makeKernelSumVisitor(lsst::daf::base::PropertySet const &ps)