24 #ifndef LSST_MEAS_MODELFIT_CModelFit_h_INCLUDED
25 #define LSST_MEAS_MODELFIT_CModelFit_h_INCLUDED
46 namespace lsst {
namespace meas {
namespace modelfit {
166 "Name of the shapelet.RadialProfile that defines the model to fit"
171 "One of 'FILE', 'LINEAR', 'EMPIRICAL', or 'NONE', indicating whether the prior should be loaded "
172 "from disk, created from one of the nested prior config/control objects, or None"
177 "Name of the Prior that defines the model to fit (a filename in $MEAS_MODELFIT_DIR/data, "
178 "with no extension), if priorSource='FILE'. Ignored for forced fitting."
183 "Configuration for a linear prior, used if priorSource='LINEAR'."
188 "Configuration for an empirical prior, used if priorSource='EMPIRICAL'."
196 "Maximum radius used in approximating profile with Gaussians (0=default for this profile)"
202 "Use per-pixel variances as weights in the nonlinear fit (the final linear fit for"
203 " flux never uses per-pixel variances)"
209 "Scale the likelihood by this factor to artificially reweight it w.r.t. the prior."
214 "Configuration for how the objective surface is explored. Ignored for forced fitting"
219 "Whether to record the steps the optimizer takes (or just the number, if running as a plugin)"
224 "Whether to record the time spent in this stage"
236 psfName(
"modelfit_DoubleShapeletPsfApprox"),
252 "Field name prefix of the Shapelet PSF approximation used to convolve the galaxy model; "
253 "must contain a set of fields matching the schema defined by shapelet.MultiShapeletFunctionKey."
258 "Configuration parameters related to the determination of the pixels to include in the fit."
263 "An initial fit (usually with a fast, approximate model) used to warm-start the exp and dev fits, "
264 "convolved with only the zeroth-order terms in the multi-shapelet PSF approximation."
269 "Independent fit of the exponential component"
274 "Independent fit of the de Vaucouleur component"
279 "Minimum initial radius in pixels (used to regularize initial moments-based PSF deconvolution)"
284 "If the 2nd-moments shape used to initialize the fit failed, use the PSF moments multiplied by this."
285 " If <= 0.0, abort the fit early instead."
324 ndarray::Array<Scalar const,1,1>
fixed;
567 void _applyForcedImpl(
577 template <
typename PixelT>
table::Key< std::string > name
An ellipse core with quadrupole moments as parameters.
Base class for all records.
Defines the fields and offsets for a table.
A mapping between the keys of two Schemas, used to copy data between them.
Record class that contains measurements made on a single exposure.
Exception to be thrown when a measurement algorithm experiences a known failure mode.
Main public interface class for CModel algorithm.
void writeResultToRecord(Result const &result, afw::table::BaseRecord &record) const
Copy values from a Result struct to a BaseRecord object.
friend class CModelAlgorithmControl
Result applyForced(afw::image::Exposure< Pixel > const &exposure, shapelet::MultiShapeletFunction const &psf, geom::Point2D const ¢er, Result const &reference, Scalar approxFlux=-1) const
Run the CModel algorithm in forced mode on an image, supplying inputs directly and returning outputs ...
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...
CModelResult Result
Typedef to the master Result struct.
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.
Control const & getControl() const
Return the control object the algorithm was constructed with.
CModelAlgorithm(std::string const &name, Control const &ctrl, afw::table::Schema &schema)
Construct an algorithm instance and add its fields to the Schema.
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.
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.
CModelControl Control
Typedef to the master Control struct.
CModelAlgorithm(Control const &ctrl)
Construct an algorithm instance that cannot use SourceRecords for input/output.
Result apply(afw::image::Exposure< Pixel > const &exposure, shapelet::MultiShapeletFunction const &psf, geom::Point2D const ¢er, 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.
Configuration object for Optimizer.
double gradientThreshold
"If the maximum of the gradient falls below this threshold, consider the algorithm converged" ;
double minTrustRadiusThreshold
"If the trust radius falls below this threshold, consider the algorithm converged" ;
int maxOuterIterations
"maximum number of steps" ;
A multi-scale shapelet function.
Registry and utility class for multi-Gaussian approximations to radial profiles.
static RadialProfile & get(std::string const &name)
Return a predefined radial profile given its name.
#define LSST_NESTED_CONTROL_FIELD(NAME, MODULE, TYPE, DOC)
A preprocessor macro used to define fields in C++ "control object" structs, for nested control object...
#define LSST_CONTROL_FIELD(NAME, TYPE, DOC)
A preprocessor macro used to define fields in C++ "control object" structs.
const char * source()
Source function that allows astChannel to source from a Stream.
double Scalar
Typedefs to be used for probability and parameter values.
A base class for image defects.
The main control object for CModel, containing parameters for the final linear fit and aggregating th...
CModelStageControl exp
"Independent fit of the exponential component" ;
PixelFitRegionControl region
"Configuration parameters related to the determination of the pixels to include in the fit....
double fallbackInitialMomentsPsfFactor
"If the 2nd-moments shape used to initialize the fit failed, use the PSF moments multiplied by this....
double minInitialRadius
"Minimum initial radius in pixels (used to regularize initial moments-based PSF deconvolution)" ;
std::string psfName
"Field name prefix of the Shapelet PSF approximation used to convolve the galaxy model; " "must conta...
CModelStageControl initial
"An initial fit (usually with a fast, approximate model) used to warm-start the exp and dev fits,...
CModelStageControl dev
"Independent fit of the de Vaucouleur component" ;
Master result object for CModel, containing results for the final linear fit and three nested CModelS...
Scalar instFlux
Flux from the final linear fit.
Scalar objective
Objective value at the best-fit point (chisq/2)
Scalar fracDev
Fraction of flux from the final linear fit in the de Vaucouleur component (always between 0 and 1).
FlagBit
Flags that apply to all four CModel fits or just the last one.
@ REGION_USED_INITIAL_ELLIPSE_MAX
Fit region implied by the best-fit ellipse of the initial was too large, so we used the configuration...
@ BAD_REFERENCE
Reference fit failed, so forced fit will fail as well.
@ REGION_USED_PSF_AREA
Kron radius was unavailable or outside bounds, so the second-moment ellipse scaled to the PSF area wa...
@ NO_SHAPELET_PSF
Set if the Psf shapelet approximation failed.
@ FAILED
General failure flag for the linear fit flux; set if any other CModel flag is set,...
@ SMALL_SHAPE
Initial moments were sufficiently small that we used minInitialRadius to set the initial parameters.
@ REGION_USED_INITIAL_ELLIPSE_MIN
Fit region implied by the best-fit ellipse of the initial was too small, so we used the configuration...
@ N_FLAGS
Non-flag counter to indicate the number of flags.
@ REGION_USED_FOOTPRINT_AREA
Kron radius was unavailable or outside bounds, so the second-moment ellipse scaled to the footprint a...
@ NO_SHAPE
Set if the input SourceRecord had no valid shape slot with which to start the fit.
@ REGION_MAX_BAD_PIXEL_FRACTION
Set if we aborted early because the fit region had too many bad pixels.
@ BAD_CENTROID
Input centroid did not land within the fit region.
@ NO_FLUX
No flux was measured.
@ REGION_MAX_AREA
Set if we aborted early because the fit region was too large.
Scalar instFluxErr
Flux uncertainty from the final linear fit.
afw::geom::ellipses::Quadrupole finalFitRegion
Pixels used in the exp, dev, and linear fits.
CModelStageResult dev
Results from the de Vaucouleur (Sersic n=4) fit.
CModelStageResult initial
Results from the initial approximate nonlinear fit that feeds the others.
CModelStageResult exp
Results from the exponential (Sersic n=1) fit.
Scalar instFluxInner
Flux measured strictly within the fit region (no extrapolation).
std::bitset< N_FLAGS > flags
Array of flags.
afw::geom::ellipses::Quadrupole initialFitRegion
Pixels used in the initial fit.
LocalUnitTransform fitSysToMeasSys
Transforms to the coordinate system where parameters are defined.
Nested control object for CModel that configures one of the three ("initial", "exp",...
std::string profileName
"Name of the shapelet.RadialProfile that defines the model to fit" ;
int nComponents
"Number of Gaussian used to approximate the profile" ;
shapelet::RadialProfile const & getProfile() const
std::string priorName
"Name of the Prior that defines the model to fit (a filename in $MEAS_MODELFIT_DIR/data,...
OptimizerControl optimizer
"Configuration for how the objective surface is explored. Ignored for forced fitting" ;
bool doRecordTime
"Whether to record the time spent in this stage" ;
bool usePixelWeights
"Use per-pixel variances as weights in the nonlinear fit (the final linear fit for" " flux never uses...
SemiEmpiricalPriorControl empiricalPriorConfig
"Configuration for an empirical prior, used if priorSource='EMPIRICAL'." ;
std::shared_ptr< Model > getModel() const
bool doRecordHistory
"Whether to record the steps the optimizer takes (or just the number, if running as a plugin)" ;
int maxRadius
"Maximum radius used in approximating profile with Gaussians (0=default for this profile)" ;
std::shared_ptr< Prior > getPrior() const
double weightsMultiplier
"Scale the likelihood by this factor to artificially reweight it w.r.t. the prior....
std::string priorSource
"One of 'FILE', 'LINEAR', 'EMPIRICAL', or 'NONE', indicating whether the prior should be loaded " "fr...
SoftenedLinearPriorControl linearPriorConfig
"Configuration for a linear prior, used if priorSource='LINEAR'." ;
Result object for a single nonlinear fitting stage of the CModel algorithm.
ndarray::Array< Scalar const, 1, 1 > fixed
Opaque fixed parameters in specialized units.
std::shared_ptr< OptimizerObjective > objfunc
Objective class used by the optimizer.
Scalar instFluxErr
Flux uncertainty from just this stage fit.
std::shared_ptr< Model > model
Model object that defines the parametrization (defined fully by Control struct)
Scalar time
Time spent in this fit in seconds.
std::shared_ptr< Prior > prior
Bayesian priors on the parameters (defined fully by Control struct)
afw::table::BaseCatalog history
Trace of the optimizer's path, if enabled by diagnostic options.
FlagBit
Flags for a single CModel stage (note that there are additional flags for the full multi-stage fit)
@ TR_SMALL
Whether convergence was due to the optimizer trust region getting too small (not a failure!...
@ N_FLAGS
Non-flag counter to indicate the number of flags.
@ MAX_ITERATIONS
Whether the optimizer exceeded the maximum number of iterations.
@ NUMERIC_ERROR
Optimizer encountered a numerical error (something likely went to infinity).
@ BAD_REFERENCE
Reference fit failed, so forced fit will fail as well.
@ FAILED
General flag, indicating whether the flux for this stage can be trusted.
@ NO_FLUX
No flux was measured.
std::shared_ptr< UnitTransformedLikelihood > likelihood
Object used to evaluate models and compare to data.
Scalar objective
Value of the objective function at the best fit point: chisq/2 - ln(prior)
ndarray::Array< Scalar const, 1, 1 > nonlinear
Opaque nonlinear parameters in specialized units.
std::bitset< N_FLAGS > flags
Array of flags.
Scalar instFluxInner
Flux measured strictly within the fit region (no extrapolation).
Scalar instFlux
Flux measured from just this stage fit.
afw::geom::ellipses::Quadrupole ellipse
Best fit half-light ellipse in pixel coordinates.
ndarray::Array< Scalar const, 1, 1 > amplitudes
Opaque linear parameters in specialized units.