File CModel.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 CModelAlgorithm
#include <CModel.h>

Main public interface class for CModel algorithm.

See CModel Magnitudes for a full description of the algorithm.

This class provides the methods that actually execute the algorithm, and (depending on how it is constructed) holds the Key objects necessary to use SourceRecords for input and output.

Public Types

typedef CModelControl Control

Typedef to the master Control struct.

typedef CModelResult Result

Typedef to the master Result struct.

Public Functions

CModelAlgorithm(std::string const &name, Control const &ctrl, afw::table::Schema &schema)

Construct an algorithm instance and add its fields to the Schema.

All fields needed to write the outputs of a regular, non-forced fit will be added to the given Schema. In addition, keys needed to retrieve the PSF shapelet approximation (assuming the ShapeletPsfApprox plugin has been run) will be extracted from the Schema.

Parameters
  • [in] name: Name of the algorithm used as a prefix for all fields added to the Schema.

  • [in] ctrl: Control object that configures the algorithm.

  • [inout] schema: Schema to which fields will be added, and from which keys for the PSF shapelet approximation will be extacted.

CModelAlgorithm(std::string const &name, Control const &ctrl, afw::table::SchemaMapper &schemaMapper)

Construct an algorithm instance suitable for forced photometry and add its fields to the Schema.

All fields needed to write the outputs of a forced fit will be added to the given SchemaMapper’s output schema. Keys needed to retrieve the reference ellipses for the exp and dev fits will be extracted from the SchemaMapper’s input schema. In addition, keys needed to retrieve the PSF shapelet approximation (assuming the ShapeletPsfApprox plugin has been run) will be extracted from the SchemaMapper’s output schema (note that the ShapeletPsfApprox plugin must be run in forced mode as well, to approximate the measurement image’s PSF rather than the reference image’s PSF, so its outputs are found in the output schema, not the input schema).

Parameters
  • [in] name: Name of the algorithm used as a prefix for all fields added to the Schema.

  • [in] ctrl: Control object that configures the algorithm.

  • [inout] schemaMapper: SchemaMapper containing input (reference) and output schemas.

CModelAlgorithm(Control const &ctrl)

Construct an algorithm instance that cannot use SourceRecords for input/output.

This constructor initializes the algorithm without initializing any of the keys necessary to operate on SourceRecords. As a result, only methods that take inputs directly and return Result objects may be called.

Control const &getControl() const

Return the control object the algorithm was constructed with.

Result apply(afw::image::Exposure<Pixel> const &exposure, shapelet::MultiShapeletFunction const &psf, geom::Point2D const &center, afw::geom::ellipses::Quadrupole const &moments, Scalar approxFlux = -1, Scalar kronRadius = -1, int footprintArea = -1) const

Run the CModel algorithm on an image, supplying inputs directly and returning outputs in a Result.

Parameters
  • [in] exposure: Image to measure. Must have a valid Psf, Wcs and PhotoCalib.

  • [in] psf: multi-shapelet approximation to the PSF at the position of the source

  • [in] center: Centroid of the source to be fit.

  • [in] moments: Non-PSF-corrected moments of the source, used to initialize the model parameters

  • [in] approxFlux: Rough estimate of the flux of the source, used to set the fit coordinate system and ensure internal parameters are of order unity. If less than or equal to zero, the sum of the flux within the footprint will be used.

  • [in] kronRadius: Estimate of the Kron radius (optional); used as the first choice when estimating the region of pixel to include in the fit.

  • [in] footprintArea: Area of the detection Fooptrint; used as the fallback when estimating the region of pixel to include in the fit.

Result applyForced(afw::image::Exposure<Pixel> const &exposure, shapelet::MultiShapeletFunction const &psf, geom::Point2D const &center, Result const &reference, Scalar approxFlux = -1) const

Run the CModel algorithm in forced mode on an image, supplying inputs directly and returning outputs in a Result.

Parameters
  • [in] exposure: Image to measure. Must have a valid Psf, Wcs and PhotoCalib.

  • [in] psf: multi-shapelet approximation to the PSF at the position of the source

  • [in] center: Centroid of the source to be fit.

  • [in] reference: Result object from a previous, non-forced run of CModelAlgorithm.

  • [in] approxFlux: Rough estimate of the flux of the source, used to set the fit coordinate system and ensure internal parameters are of order unity. If less than or equal to zero, the sum of the flux within the footprint will be used.

void measure(afw::table::SourceRecord &measRecord, afw::image::Exposure<Pixel> const &exposure) const

Run the CModel algorithm on an image, using a SourceRecord for inputs and outputs.

To run this method, the

CModelAlgorithm instance must have been created using the constructor that takes a Schema argument, and that Schema must match the Schema of the SourceRecord passed here.
Parameters
  • [inout] measRecord: A SourceRecord instance used to provide a Footprint, the centroid and shape of the source, a MultiShapeletFunction PSF, and an approximate estimate of the (via the PsfFlux slot), and to which all outputs will be written.

  • [in] exposure: Image to be measured. Must have a valid Psf, Wcs, and PhotoCalib.

void measure(afw::table::SourceRecord &measRecord, afw::image::Exposure<Pixel> const &exposure, afw::table::SourceRecord const &refRecord) const

Run the CModel algorithm in forced mode on an image, using a SourceRecord for inputs and outputs.

To run this method, the

CModelAlgorithm instance must have been created using the constructor that takes a Schema argument, and that Schema must match the Schema of the SourceRecord passed here.
Parameters
  • [inout] measRecord: A SourceRecord instance used to provide a Footprint, the centroid of the source, a MultiShapeletFunction PSF, and an approximate estimate of the (via the PsfFlux slot), and to which all outputs will be written.

  • [in] exposure: Image to be measured. Must have a valid Psf, Wcs, and PhotoCalib.

  • [in] refRecord: A SourceRecord that contains the outputs of a previous non-forced run of CModelAlgorithm (which may have taken place on an image with a different Wcs).

void fail(afw::table::SourceRecord &measRecord, meas::base::MeasurementError *error) const

Handle an exception thrown by one of the measure() methods, setting the appropriate flag in the given record.

Parameters
  • [out] measRecord: Record on which the flag should be set.

  • [in] error: Error containing the bit to be set. If null, only the general failure bit will be set.

void writeResultToRecord(Result const &result, afw::table::BaseRecord &record) const

Copy values from a Result struct to a BaseRecord object.

Private Functions

void _applyImpl(Result &result, afw::image::Exposure<Pixel> const &exposure, shapelet::MultiShapeletFunction const &psf, geom::Point2D const &center, afw::geom::ellipses::Quadrupole const &moments, Scalar approxFlux, Scalar kronRadius = -1, int footprintArea = -1) const
void _applyForcedImpl(Result &result, afw::image::Exposure<Pixel> const &exposure, shapelet::MultiShapeletFunction const &psf, geom::Point2D const &center, Result const &reference, Scalar approxFlux) const
template<typename PixelT>
shapelet::MultiShapeletFunction _processInputs(afw::table::SourceRecord &source, afw::image::Exposure<PixelT> const &exposure) const
PTR(Impl)

Private Members

Control _ctrl

Friends

friend lsst::meas::modelfit::CModelAlgorithmControl
struct CModelControl
#include <CModel.h>

The main control object for CModel, containing parameters for the final linear fit and aggregating the other control objects.

Public Functions

CModelControl()
lsst::meas::modelfit::CModelControl::LSST_CONTROL_FIELD(psfName, std::string, "Field name prefix of the Shapelet PSF approximation used to convolve the galaxy model; " "must contain a set of fields matching the schema defined by shapelet.MultiShapeletFunctionKey.")
lsst::meas::modelfit::CModelControl::LSST_NESTED_CONTROL_FIELD(region, lsst.meas.modelfit. pixelFitRegion, PixelFitRegionControl, "Configuration parameters related to the determination of the pixels to include in the fit.")
lsst::meas::modelfit::CModelControl::LSST_NESTED_CONTROL_FIELD(initial, lsst.meas.modelfit. cmodel, CModelStageControl, "An initial fit (usually with a fast, approximate model) used to warm-start the exp and dev fits, " "convolved with only the zeroth-order terms in the multi-shapelet PSF approximation.")
lsst::meas::modelfit::CModelControl::LSST_NESTED_CONTROL_FIELD(exp, lsst.meas.modelfit. cmodel, CModelStageControl, "Independent fit of the exponential component")
lsst::meas::modelfit::CModelControl::LSST_NESTED_CONTROL_FIELD(dev, lsst.meas.modelfit. cmodel, CModelStageControl, "Independent fit of the de Vaucouleur component")
lsst::meas::modelfit::CModelControl::LSST_CONTROL_FIELD(minInitialRadius, double, "Minimum initial radius in pixels (used to regularize initial moments-based PSF deconvolution)")
lsst::meas::modelfit::CModelControl::LSST_CONTROL_FIELD(fallbackInitialMomentsPsfFactor, double, "If the 2nd-moments shape used to initialize the fit failed, use the PSF moments multiplied by this." " If<=0. 0, abort the fit early instead.")
struct CModelResult
#include <CModel.h>

Master result object for CModel, containing results for the final linear fit and three nested CModelStageResult objects for the results of the previous stages.

Public Types

enum FlagBit

Flags that apply to all four CModel fits or just the last one.

Values:

FAILED = 0

General failure flag for the linear fit flux; set if any other CModel flag is set, or if any of the three previous stages failed.

REGION_MAX_AREA

Set if we aborted early because the fit region was too large.

REGION_MAX_BAD_PIXEL_FRACTION

Set if we aborted early because the fit region had too many bad pixels.

REGION_USED_FOOTPRINT_AREA

Kron radius was unavailable or outside bounds, so the second-moment ellipse scaled to the footprint area was used instead.

REGION_USED_PSF_AREA

Kron radius was unavailable or outside bounds, so the second-moment ellipse scaled to the PSF area was used instead.

REGION_USED_INITIAL_ELLIPSE_MIN

Fit region implied by the best-fit ellipse of the initial was too small, so we used the configuration minimum instead.

REGION_USED_INITIAL_ELLIPSE_MAX

Fit region implied by the best-fit ellipse of the initial was too large, so we used the configuration maximum instead.

NO_SHAPE

Set if the input SourceRecord had no valid shape slot with which to start the fit.

SMALL_SHAPE

Initial moments were sufficiently small that we used minInitialRadius to set the initial parameters.

NO_SHAPELET_PSF

Set if the Psf shapelet approximation failed.

BAD_CENTROID

Input centroid did not land within the fit region.

BAD_REFERENCE

Reference fit failed, so forced fit will fail as well.

N_FLAGS

Non-flag counter to indicate the number of flags.

Public Functions

CModelResult()

Public Members

Scalar instFlux

Flux from the final linear fit.

Scalar instFluxErr

Flux uncertainty from the final linear fit.

Scalar instFluxInner

Flux measured strictly within the fit region (no extrapolation).

Scalar fracDev

Fraction of flux from the final linear fit in the de Vaucouleur component (always between 0 and 1).

Scalar objective

Objective value at the best-fit point (chisq/2)

CModelStageResult initial

Results from the initial approximate nonlinear fit that feeds the others.

CModelStageResult exp

Results from the exponential (Sersic n=1) fit.

CModelStageResult dev

Results from the de Vaucouleur (Sersic n=4) fit.

afw::geom::ellipses::Quadrupole initialFitRegion

Pixels used in the initial fit.

afw::geom::ellipses::Quadrupole finalFitRegion

Pixels used in the exp, dev, and linear fits.

LocalUnitTransform fitSysToMeasSys

Transforms to the coordinate system where parameters are defined.

std::bitset<N_FLAGS> flags

Array of flags.

struct CModelStageControl
#include <CModel.h>

Nested control object for CModel that configures one of the three (“initial”, “exp”, “dev”) nonlinear fitting stages.

Public Functions

CModelStageControl()
shapelet::RadialProfile const &getProfile() const
PTR(Model) const
PTR(Prior) const
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(profileName, std::string, "Name of the shapelet.RadialProfile that defines the model to fit")
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(priorSource, std::string, "One of 'FILE', 'LINEAR', 'EMPIRICAL', or 'NONE', indicating whether the prior should be loaded " "from disk, created from one of the nested prior config/control objects, or None")
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(priorName, std::string, "Name of the Prior that defines the model to fit a filename in $MEAS_MODELFIT_DIR/data, " "with no extension)
lsst::meas::modelfit::CModelStageControl::LSST_NESTED_CONTROL_FIELD(linearPriorConfig, lsst.meas.modelfit. priors, SoftenedLinearPriorControl, "Configuration for a linear prior)
lsst::meas::modelfit::CModelStageControl::LSST_NESTED_CONTROL_FIELD(empiricalPriorConfig, lsst.meas.modelfit. priors, SemiEmpiricalPriorControl, "Configuration for an empirical prior)
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(nComponents, int, "Number of Gaussian used to approximate the profile")
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(maxRadius, int, "Maximum radius used in approximating profile with Gaussians (0=default for this profile)")
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(usePixelWeights, bool, "Use per-pixel variances as weights in the nonlinear fit (the final linear fit for" " flux never uses per-pixel variances)")
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(weightsMultiplier, double, "Scale the likelihood by this factor to artificially reweight it w.r.t. the prior.")
lsst::meas::modelfit::CModelStageControl::LSST_NESTED_CONTROL_FIELD(optimizer, lsst.meas.modelfit. optimizer, OptimizerControl, "Configuration for how the objective surface is explored. Ignored for forced fitting")
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(doRecordHistory, bool, "Whether to record the steps the optimizer takes (or just the number, if running as a plugin)")
lsst::meas::modelfit::CModelStageControl::LSST_CONTROL_FIELD(doRecordTime, bool, "Whether to record the time spent in this stage")
struct CModelStageResult
#include <CModel.h>

Result object for a single nonlinear fitting stage of the CModel algorithm

Public Types

enum FlagBit

Flags for a single CModel stage (note that there are additional flags for the full multi-stage fit)

Values:

FAILED = 0

General flag, indicating whether the flux for this stage can be trusted.

TR_SMALL

Whether convergence was due to the optimizer trust region getting too small (not a failure!)

MAX_ITERATIONS

Whether the optimizer exceeded the maximum number of iterations. Indicates a suspect fit, but not necessarily a bad one (implies FAILED).

NUMERIC_ERROR

Optimizer encountered a numerical error (something likely went to infinity). Result will be unusable; implies FAILED.

BAD_REFERENCE

Reference fit failed, so forced fit will fail as well.

N_FLAGS

Non-flag counter to indicate the number of flags.

Public Functions

CModelStageResult()
PTR(Model)

Model object that defines the parametrization (defined fully by Control struct)

PTR(Prior)

Bayesian priors on the parameters (defined fully by Control struct)

PTR(OptimizerObjective)

Objective class used by the optimizer.

PTR(UnitTransformedLikelihood)

Object used to evaluate models and compare to data.

Public Members

Scalar instFlux

Flux measured from just this stage fit.

Scalar instFluxErr

Flux uncertainty from just this stage fit.

Scalar instFluxInner

Flux measured strictly within the fit region (no extrapolation).

Scalar objective

Value of the objective function at the best fit point: chisq/2 - ln(prior)

Scalar time

Time spent in this fit in seconds.

afw::geom::ellipses::Quadrupole ellipse

Best fit half-light ellipse in pixel coordinates.

ndarray::Array<Scalar const, 1, 1> nonlinear

Opaque nonlinear parameters in specialized units.

ndarray::Array<Scalar const, 1, 1> amplitudes

Opaque linear parameters in specialized units.

ndarray::Array<Scalar const, 1, 1> fixed

Opaque fixed parameters in specialized units.

afw::table::BaseCatalog history

Trace of the optimizer’s path, if enabled by diagnostic options.

std::bitset<N_FLAGS> flags

Array of flags.