File Mixture.h

namespace lsst

Class for a simple mapping implementing a generic AstrometryTransform.

Remove all non-astronomical counts from the Chunk Exposure’s pixels.

Forward declarations for lsst::utils::Cache

For details on the Cache class, see the Cache.h file.

It uses a template rather than a pointer so that the derived classes can use the specifics of the transform. The class simplePolyMapping overloads a few routines.

A base class for image defects

Numeric constants used by the Integrate.h integrator routines.

Compute Image Statistics

Note

Gauss-Kronrod-Patterson quadrature coefficients for use in quadpack routine qng. These coefficients were calculated with 101 decimal digit arithmetic by L. W. Fullerton, Bell Labs, Nov 1981.

Note

The Statistics class itself can only handle lsst::afw::image::MaskedImage() types. The philosophy has been to handle other types by making them look like lsst::afw::image::MaskedImage() and reusing that code. Users should have no need to instantiate a Statistics object directly, but should use the overloaded makeStatistics() factory functions.

namespace meas
namespace modelfit
class Mixture : public lsst::afw::table::io::PersistableFacade<Mixture>, public lsst::afw::table::io::Persistable

Unnamed Group

iterator begin()

Iterator and indexed access to components.

While mutable iterators and accessors are provided, any modifications to the component weights should be followed by a call to normalize(), as other member functions will not work properly if the mixture is not normalized.

iterator end()
const_iterator begin() const
const_iterator end() const
Component &operator[](std::size_t i)
Component const &operator[](std::size_t i) const

Public Types

typedef MixtureComponent Component
typedef MixtureUpdateRestriction UpdateRestriction
typedef std::vector<Component> ComponentList
typedef ComponentList::iterator iterator
typedef ComponentList::const_iterator const_iterator

Public Functions

std::size_t size() const

Return the number of components.

virtual int getComponentCount() const

Return the number of components.

PTR(Mixture) const

Project the distribution onto the given dimensions (marginalize over all others)

PTR(Mixture)

Project the distribution onto the given dimensions (marginalize over all others)

int getDimension() const

Return the number of dimensions.

void normalize()

Iterate over all components, rescaling their weights so they sum to one.

void shift(int dim, Scalar offset)

Shift the mixture in the given dimension, adding the given offset to all mu vectors.

std::size_t clip(Scalar threshold = 0.0)

Iterate over all components, removing those with weight less than or equal to threshold.

The weights will be normalized if any are removed.

Return

the number of components removed.

Scalar getDegreesOfFreedom() const

Get the number of degrees of freedom in the component Student’s T distributions (inf=Gaussian)

void setDegreesOfFreedom(Scalar df = std::numeric_limits<Scalar>::infinity())

Set the number of degrees of freedom in the component Student’s T distributions (inf=Gaussian)

template<typename Derived>
Scalar evaluate(Component const &component, Eigen::MatrixBase<Derived> const &x) const

Evaluate the probability density at the given point for the given component distribution.

This evaluates the probability of a single component, including the current weight of that component.

template<typename Derived>
Scalar evaluate(Eigen::MatrixBase<Derived> const &x) const

Evaluate the mixture distribution probability density function (PDF) at the given points.

Parameters
  • [in] x: point to evaluate, as an Eigen expression, shape=(dim,)

void evaluate(ndarray::Array<Scalar const, 2, 1> const &x, ndarray::Array<Scalar, 1, 0> const &p) const

Evaluate the distribution probability density function (PDF) at the given points.

Parameters
  • [in] x: array of points, shape=(numSamples, dim)

  • [out] p: array of probability values, shape=(numSamples,)

void evaluateComponents(ndarray::Array<Scalar const, 2, 1> const &x, ndarray::Array<Scalar, 2, 1> const &p) const

Evaluate the contributions of each component to the full probability at the given points.

Parameters
  • [in] x: points to evaluate at, with number of columns equal to the number of dimensions

  • [in] p: array to fill, with number of columns equal to the number of components

void evaluateDerivatives(ndarray::Array<Scalar const, 1, 1> const &x, ndarray::Array<Scalar, 1, 1> const &gradient, ndarray::Array<Scalar, 2, 1> const &hessian) const

Evaluate the derivative of the distribution at the given point.

Parameters
  • [in] x: point to evaluate the derivative, with size equal to the number of dimensions

  • [in] gradient: 1st derivative array to fill

  • [in] hessian: 2nd derivative array to fill

void evaluateDerivatives(Eigen::Matrix<Scalar, Eigen::Dynamic, 1> &x, Eigen::Matrix<Scalar, Eigen::Dynamic, 1> &gradient, Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> &hessian) const

Evaluate the derivative of the distribution at the given point.

Parameters
  • [in] x: point to evaluate the derivative, with size equal to the number of dimensions

  • [in] gradient: 1st derivative array to fill

  • [in] hessian: 2nd derivative array to fill

void evaluateDerivatives(Eigen::Matrix<Scalar, Eigen::Dynamic, 1> &x, Eigen::Matrix<Scalar, Eigen::Dynamic, 1> &gradient) const

Evaluate the derivative of the distribution at the given point.

Parameters
  • [in] x: point to evaluate the derivative, with size equal to the number of dimensions

  • [in] gradient: 1st derivative array to fill

void draw(afw::math::Random &rng, ndarray::Array<Scalar, 2, 1> const &x) const

Draw random variates from the distribution.

Parameters
  • [inout] rng: random number generator

  • [out] x: array of points, shape=(numSamples, dim)

void updateEM(ndarray::Array<Scalar const, 2, 1> const &x, ndarray::Array<Scalar const, 1, 0> const &w, Scalar tau1 = 0.0, Scalar tau2 = 0.5)

Perform an Expectation-Maximization step, updating the component parameters to match the given weighted samples.

The updates to the

\(\sigma\) matrices are damped according to:
\[ \sigma_d = \alpha\sigma_1 + (1-\alpha)\sigma_0 \]
Where \(\sigma_0\) is the previous matrix, \(\sigma_1\) is the undamped update, and \(\sigma_d\) is the damped update. The parameter \(\alpha\) is set by the ratio of the determinants:
\[ r \equiv \frac{|\sigma_1|}{|\sigma_0|} \]
When \(r \ge \tau_1\), \(\alpha=1\); when \(r \lt \tau_1\), it is rolled off quadratically to \(\tau_2\).
Parameters
  • [in] x: array of variables, shape=(numSamples, dim)

  • [in] w: array of weights, shape=(numSamples,)

  • [in] tau1: damping parameter (see below)

  • [in] tau2: damping parameter (see below)

void updateEM(ndarray::Array<Scalar const, 2, 1> const &x, ndarray::Array<Scalar const, 1, 0> const &w, UpdateRestriction const &restriction, Scalar tau1 = 0.0, Scalar tau2 = 0.5)

Perform an Expectation-Maximization step, updating the component parameters to match the given weighted samples.

Parameters
  • [in] x: array of variables, shape=(numSamples, dim)

  • [in] w: array of weights, shape=(numSamples,)

  • [in] restriction: Functor used to restrict the form of the updated mu and sigma

  • [in] tau1: damping parameter (see Mixture::updateEM)

  • [in] tau2: damping parameter (see Mixture::updateEM)

void updateEM(ndarray::Array<Scalar const, 2, 1> const &x, UpdateRestriction const &restriction, Scalar tau1 = 0.0, Scalar tau2 = 0.5)

Perform an Expectation-Maximization step, updating the component parameters to match the given unweighted samples.

Parameters
  • [in] x: array of variables, shape=(numSamples, dim)

  • [in] restriction: Functor used to restrict the form of the updated mu and sigma

  • [in] tau1: damping parameter (see Mixture::updateEM)

  • [in] tau2: damping parameter (see Mixture::updateEM)

virtual PTR(Mixture) const

Polymorphic deep copy.

Mixture(int dim, ComponentList &components, Scalar df = std::numeric_limits<Scalar>::infinity())

Construct a mixture model.

The components will be automatically normalized after construction.

Parameters
  • [in] dim: Dimensionality of the distribution

  • [in] df: Number of degrees of freedom for component Student’s T distributions (inf=Gaussian)

  • [in] components: List of components; will be emptied on return.

virtual bool isPersistable() const

Return true if this particular object can be persisted using afw::table::io.

Public Members

int dim2 lsst::meas::modelfit::Mixture::const

Protected Functions

std::string getPythonModule() const

Return the fully-qualified Python module that should be imported to guarantee that its factory is registered.

Must be less than ArchiveIndexSchema::MAX_MODULE_LENGTH characters.

Will be ignored if empty.

std::string getPersistenceName() const

Return the unique name used to persist this object and look up its factory.

Must be less than ArchiveIndexSchema::MAX_NAME_LENGTH characters.

void write(OutputArchiveHandle &handle) const

Write the object to one or more catalogs.

The handle object passed to this function provides an interface for adding new catalogs and adding nested objects to the same archive (while checking for duplicates). See OutputArchiveHandle for more information.

Private Functions

template<typename A, typename B, typename C>
void _evaluateDerivativesImpl(A const &x, B &gradient, C *hessian, bool computeHessian = true) const
template<typename Derived>
Scalar _computeZ(Component const &component, Eigen::MatrixBase<Derived> const &x) const
void updateDampedSigma(int k, Matrix const &sigma, double tau1, double tau2)
Scalar _evaluate(Scalar z) const
void _stream(std::ostream &os) const

Private Members

bool _isGaussian
int _dim
Scalar _df
Scalar _norm
Vector _workspace
ComponentList _components

Friends

std::ostream &operator<<(std::ostream &os, Mixture const &self)
class MixtureComponent
#include <Mixture.h>

A weighted Student’s T or Gaussian distribution used as a component in a Mixture.

Unnamed Group

Vector getMu() const

Get/set the location parameter (mean/median/mode) of this component.

void setMu(Vector const &mu)

Unnamed Group

Matrix getSigma() const

Get/set the shape/size parameter.

For the Gaussian distribution, this is simply the covariance matrix. For the Student’s T distribution with df > 2, covariance = sigma * df / (df - 2); for df <= 2, the Student’s T distribution has infinite variance, but is still a valid distribution.

void setSigma(Matrix const &sigma)

Public Functions

int getDimension() const

Return the number of dimensions.

MixtureComponent project(int dim) const

Project the distribution onto the given dimension (marginalize over all others)

MixtureComponent project(int dim1, int dim2) const

Project the distribution onto the given dimensions (marginalize over all others)

MixtureComponent(int dim)

Default-construct a mixture component with weight=1, mu=0, sigma=identity.

MixtureComponent(Scalar weight_, Vector const &mu, Matrix const &sigma)

Default-construct a mixture component with the given parameters.

MixtureComponent &operator=(MixtureComponent const &other)

Public Members

Scalar weight

Weight of this distribution in the mixture.

Private Functions

void _stream(std::ostream &os, int offset = 0) const

Private Members

Scalar _sqrtDet
Vector _mu
Eigen::LLT<Matrix> _sigmaLLT

Friends

friend lsst::meas::modelfit::Mixture
std::ostream &operator<<(std::ostream &os, MixtureComponent const &self)
class MixtureUpdateRestriction
#include <Mixture.h>

Helper class used to define restrictions to the form of the component parameters in Mixture::updateEM.

The base class implementation does not apply any restrictions.

Public Functions

int getDimension() const
virtual void restrictMu(Vector &mu) const
virtual void restrictSigma(Matrix &sigma) const
virtual ~MixtureUpdateRestriction()
MixtureUpdateRestriction(int dim)

Private Members

int _dim