lsst.meas.modelfit  15.0-4-g535e784+9
Public Types | Public Member Functions | List of all members
lsst::meas::modelfit::MomentsModel Class Reference

A model for calculating observed shape moments and their derivatives. More...

#include <RegularizedMoments.h>

Public Types

using Element = double
 
using Moments = Eigen::Matrix< Element, 6, 1 >
 
using Jacobian = Eigen::Matrix< Element, 6, 6 >
 
using FirstMoment = Eigen::Matrix< Element, 2, 1 >
 
using SecondMoment = Eigen::Matrix< Element, 2, 2 >
 

Public Member Functions

 MomentsModel (Moments const &W)
 Construct a MomentsModel object with a given vector of weight function moments. More...
 
 MomentsModel (const MomentsModel &)=default
 
 MomentsModel (MomentsModel &&)=default
 
void setParameters (Moments const &Q)
 Sets the moments corresponding to the astronomical object's intrinsic moments. More...
 
Moments computeValues () const
 Returns the value of the biased weighted moments given the input weight and object's intrinsic moments. More...
 
Jacobian computeJacobian () const
 Computes and returns the gradient of the model biased weighted moments along each of the intrinsic moment dimensions. More...
 

Detailed Description

A model for calculating observed shape moments and their derivatives.

This class models the biased weighted moments of an astronomical source given a vector of moments intrinsic to an object, and a vector of moments representing a weighting function.

The constructor takes the moments of weight function, and the moments of the object are set by calling the setParameters method. This allows a single instance to be reused to evaluate multiple positions in object moment space. This structure is computationally more efficient as well since the weight function would remain constant for many evaluations of image moments

An instance of this class is designed to be used as an objective function in an optimizer which will find the maximally likely image moments given moments measured from an image. As such the class provides two methods in addition to the setParameters method already mentioned. One method computes the biased weighted moments given the specified weight and "true" moments, and the other computes the derivative of the objective function along each of the image moment variables.

Definition at line 52 of file RegularizedMoments.h.

Member Typedef Documentation

◆ Element

Definition at line 54 of file RegularizedMoments.h.

◆ FirstMoment

Definition at line 57 of file RegularizedMoments.h.

◆ Jacobian

Definition at line 56 of file RegularizedMoments.h.

◆ Moments

Definition at line 55 of file RegularizedMoments.h.

◆ SecondMoment

Definition at line 58 of file RegularizedMoments.h.

Constructor & Destructor Documentation

◆ MomentsModel() [1/3]

lsst::meas::modelfit::MomentsModel::MomentsModel ( Moments const &  W)
inlineexplicit

Construct a MomentsModel object with a given vector of weight function moments.

Parameters
[in]WAn eigen vector of the moments of the elliptical Gaussian weighting function. Most likely the first element (the zeroth moment) should be one, representing a normalized weight function.

Definition at line 67 of file RegularizedMoments.h.

◆ MomentsModel() [2/3]

lsst::meas::modelfit::MomentsModel::MomentsModel ( const MomentsModel )
default

◆ MomentsModel() [3/3]

lsst::meas::modelfit::MomentsModel::MomentsModel ( MomentsModel &&  )
default

Member Function Documentation

◆ computeJacobian()

Jacobian lsst::meas::modelfit::MomentsModel::computeJacobian ( ) const

Computes and returns the gradient of the model biased weighted moments along each of the intrinsic moment dimensions.

Parameters
[out]Asix by six Eigen matrix, with the biased weighted moments down the column, and the derivative of the column moment with respect to the intrinsic moments along the rows.

◆ computeValues()

Moments lsst::meas::modelfit::MomentsModel::computeValues ( ) const

Returns the value of the biased weighted moments given the input weight and object's intrinsic moments.

Parameters
[out]Eigenvector of the model biased weighted moments given the input weight and intrinsic moments

◆ setParameters()

void lsst::meas::modelfit::MomentsModel::setParameters ( Moments const &  Q)

Sets the moments corresponding to the astronomical object's intrinsic moments.

When this function is called, it sets the astronomical object's intrinsic moments for all subsequent calculations. Additionally after setting the moment, this function calculates and caches the value of the biased weighted moments given the weights and intrinsic moments.

Parameters
[in]QAn eigen vector of the intrinsic moments of the astronomical object

The documentation for this class was generated from the following file: