LSSTApplications  16.0-10-g0ee56ad+5,16.0-11-ga33d1f2+5,16.0-12-g3ef5c14+3,16.0-12-g71e5ef5+18,16.0-12-gbdf3636+3,16.0-13-g118c103+3,16.0-13-g8f68b0a+3,16.0-15-gbf5c1cb+4,16.0-16-gfd17674+3,16.0-17-g7c01f5c+3,16.0-18-g0a50484+1,16.0-20-ga20f992+8,16.0-21-g0e05fd4+6,16.0-21-g15e2d33+4,16.0-22-g62d8060+4,16.0-22-g847a80f+4,16.0-25-gf00d9b8+1,16.0-28-g3990c221+4,16.0-3-gf928089+3,16.0-32-g88a4f23+5,16.0-34-gd7987ad+3,16.0-37-gc7333cb+2,16.0-4-g10fc685+2,16.0-4-g18f3627+26,16.0-4-g5f3a788+26,16.0-5-gaf5c3d7+4,16.0-5-gcc1f4bb+1,16.0-6-g3b92700+4,16.0-6-g4412fcd+3,16.0-6-g7235603+4,16.0-69-g2562ce1b+2,16.0-8-g14ebd58+4,16.0-8-g2df868b+1,16.0-8-g4cec79c+6,16.0-8-gadf6c7a+1,16.0-8-gfc7ad86,16.0-82-g59ec2a54a+1,16.0-9-g5400cdc+2,16.0-9-ge6233d7+5,master-g2880f2d8cf+3,v17.0.rc1
LSSTDataManagementBasePackage
CModel.h
Go to the documentation of this file.
1 // -*- lsst-c++ -*-
2 /*
3  * LSST Data Management System
4  * Copyright 2008-2013 LSST Corporation.
5  *
6  * This product includes software developed by the
7  * LSST Project (http://www.lsst.org/).
8  *
9  * This program is free software: you can redistribute it and/or modify
10  * it under the terms of the GNU General Public License as published by
11  * the Free Software Foundation, either version 3 of the License, or
12  * (at your option) any later version.
13  *
14  * This program is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17  * GNU General Public License for more details.
18  *
19  * You should have received a copy of the LSST License Statement and
20  * the GNU General Public License along with this program. If not,
21  * see <http://www.lsstcorp.org/LegalNotices/>.
22  */
23 
24 #ifndef LSST_MEAS_MODELFIT_CModelFit_h_INCLUDED
25 #define LSST_MEAS_MODELFIT_CModelFit_h_INCLUDED
26 
27 #include "ndarray.h"
28 
29 #include "lsst/pex/config.h"
31 #include "lsst/afw/table/Source.h"
41 
42 namespace lsst { namespace meas { namespace modelfit {
43 
139 
141  profileName("lux"),
142  priorSource("EMPIRICAL"),
143  priorName(),
144  nComponents(8),
145  maxRadius(0),
146  usePixelWeights(false),
147  weightsMultiplier(1.0),
148  doRecordHistory(true),
149  doRecordTime(true)
150  {}
151 
154  }
155 
156  PTR(Model) getModel() const;
157 
158  PTR(Prior) getPrior() const;
159 
162  "Name of the shapelet.RadialProfile that defines the model to fit"
163  );
164 
167  "One of 'FILE', 'LINEAR', 'EMPIRICAL', or 'NONE', indicating whether the prior should be loaded "
168  "from disk, created from one of the nested prior config/control objects, or None"
169  );
170 
173  "Name of the Prior that defines the model to fit (a filename in $MEAS_MODELFIT_DIR/data, "
174  "with no extension), if priorSource='FILE'. Ignored for forced fitting."
175  );
176 
179  "Configuration for a linear prior, used if priorSource='LINEAR'."
180  );
181 
184  "Configuration for an empirical prior, used if priorSource='EMPIRICAL'."
185  );
186 
187  LSST_CONTROL_FIELD(nComponents, int, "Number of Gaussian used to approximate the profile");
188 
190  maxRadius,
191  int,
192  "Maximum radius used in approximating profile with Gaussians (0=default for this profile)"
193  );
194 
197  bool,
198  "Use per-pixel variances as weights in the nonlinear fit (the final linear fit for"
199  " flux never uses per-pixel variances)"
200  );
201 
204  double,
205  "Scale the likelihood by this factor to artificially reweight it w.r.t. the prior."
206  );
207 
210  "Configuration for how the objective surface is explored. Ignored for forced fitting"
211  );
212 
214  doRecordHistory, bool,
215  "Whether to record the steps the optimizer takes (or just the number, if running as a plugin)"
216  );
217 
219  doRecordTime, bool,
220  "Whether to record the time spent in this stage"
221  );
222 
223 };
224 
230 
232  psfName("modelfit_DoubleShapeletPsfApprox"),
233  minInitialRadius(0.1),
234  fallbackInitialMomentsPsfFactor(1.5)
235  {
236  initial.nComponents = 3; // use very rough model in initial fit
237  initial.optimizer.gradientThreshold = 1E-3; // with slightly coarser convergence criteria
238  initial.optimizer.minTrustRadiusThreshold = 1E-2;
239  initial.usePixelWeights = true;
240  dev.profileName = "luv";
241  exp.nComponents = 6;
242  exp.optimizer.maxOuterIterations = 250;
243  }
244 
246  psfName,
247  std::string,
248  "Field name prefix of the Shapelet PSF approximation used to convolve the galaxy model; "
249  "must contain a set of fields matching the schema defined by shapelet.MultiShapeletFunctionKey."
250  );
251 
254  "Configuration parameters related to the determination of the pixels to include in the fit."
255  );
256 
259  "An initial fit (usually with a fast, approximate model) used to warm-start the exp and dev fits, "
260  "convolved with only the zeroth-order terms in the multi-shapelet PSF approximation."
261  );
262 
265  "Independent fit of the exponential component"
266  );
267 
270  "Independent fit of the de Vaucouleur component"
271  );
272 
274  minInitialRadius, double,
275  "Minimum initial radius in pixels (used to regularize initial moments-based PSF deconvolution)"
276  );
277 
279  fallbackInitialMomentsPsfFactor, double,
280  "If the 2nd-moments shape used to initialize the fit failed, use the PSF moments multiplied by this."
281  " If <= 0.0, abort the fit early instead."
282  );
283 
284 };
285 
290 
292  enum FlagBit {
293  FAILED=0,
294  TR_SMALL,
295  MAX_ITERATIONS,
297  NUMERIC_ERROR,
299  BAD_REFERENCE,
301  N_FLAGS
302  };
303 
305 
316 
317  ndarray::Array<Scalar const,1,1> nonlinear;
318  ndarray::Array<Scalar const,1,1> amplitudes;
319  ndarray::Array<Scalar const,1,1> fixed;
320 
323 };
324 
329 struct CModelResult {
330 
332  enum FlagBit {
333  FAILED=0,
334  REGION_MAX_AREA,
336  REGION_MAX_BAD_PIXEL_FRACTION,
337  REGION_USED_FOOTPRINT_AREA,
339  REGION_USED_PSF_AREA,
342  REGION_USED_INITIAL_ELLIPSE_MIN,
344  REGION_USED_INITIAL_ELLIPSE_MAX,
346  NO_SHAPE,
348  SMALL_SHAPE,
350  NO_SHAPELET_PSF,
354  N_FLAGS
355  };
356 
357  CModelResult();
358 
363  Scalar objective;
365 
369 
372 
375 };
376 
386 public:
387 
390 
404  std::string const & name,
405  Control const & ctrl,
407  );
408 
425  std::string const & name,
426  Control const & ctrl,
427  afw::table::SchemaMapper & schemaMapper
428  );
429 
437  explicit CModelAlgorithm(Control const & ctrl);
438 
440  Control const & getControl() const { return _ctrl; }
441 
458  Result apply(
459  afw::image::Exposure<Pixel> const & exposure,
461  afw::geom::Point2D const & center,
463  Scalar approxFlux=-1,
464  Scalar kronRadius=-1,
465  int footprintArea=-1
466  ) const;
467 
480  Result applyForced(
481  afw::image::Exposure<Pixel> const & exposure,
483  afw::geom::Point2D const & center,
484  Result const & reference,
485  Scalar approxFlux=-1
486  ) const;
487 
500  void measure(
501  afw::table::SourceRecord & measRecord,
502  afw::image::Exposure<Pixel> const & exposure
503  ) const;
504 
520  void measure(
521  afw::table::SourceRecord & measRecord,
522  afw::image::Exposure<Pixel> const & exposure,
523  afw::table::SourceRecord const & refRecord
524  ) const;
525 
534  void fail(
535  afw::table::SourceRecord & measRecord,
537  ) const;
538 
540  void writeResultToRecord(Result const & result, afw::table::BaseRecord & record) const;
541 
542 private:
543 
544  friend class CModelAlgorithmControl;
545 
546  // Actual implementations go here; we use an output argument for the result so we can get partial
547  // results to the plugin version when we throw.
548  void _applyImpl(
549  Result & result,
550  afw::image::Exposure<Pixel> const & exposure,
552  afw::geom::Point2D const & center,
553  afw::geom::ellipses::Quadrupole const & moments,
554  Scalar approxFlux,
555  Scalar kronRadius=-1,
556  int footprintArea=-1
557  ) const;
558 
559  // Actual implementations go here; we use an output argument for the result so we can get partial
560  // results to the SourceRecord version when we throw.
561  void _applyForcedImpl(
562  Result & result,
563  afw::image::Exposure<Pixel> const & exposure,
565  afw::geom::Point2D const & center,
566  Result const & reference,
567  Scalar approxFlux
568  ) const;
569 
570  // gets/checks inputs from SourceRecord that are needed by both apply and applyForced
571  template <typename PixelT>
572  shapelet::MultiShapeletFunction _processInputs(
574  afw::image::Exposure<PixelT> const & exposure
575  ) const;
576 
577  class Impl;
578 
579  Control _ctrl;
580  PTR(Impl) _impl;
581 };
582 
583 }}} // namespace lsst::meas::modelfit
584 
585 #endif // !LSST_MEAS_MODELFIT_CModelFit_h_INCLUDED
An ellipse core with quadrupole moments as parameters.
Definition: Quadrupole.h:47
Defines the fields and offsets for a table.
Definition: Schema.h:50
shapelet::RadialProfile const & getProfile() const
Definition: CModel.h:152
LocalUnitTransform fitSysToMeasSys
Transforms to the coordinate system where parameters are defined.
Definition: CModel.h:373
boost::shared_ptr< OptimizerObjective > objfunc
Objective class used by the optimizer.
Definition: CModel.h:308
Control const & getControl() const
Return the control object the algorithm was constructed with.
Definition: CModel.h:440
Scalar objective
Value of the objective function at the best fit point: chisq/2 - ln(prior)
Definition: CModel.h:313
A mapping between the keys of two Schemas, used to copy data between them.
Definition: SchemaMapper.h:21
Scalar instFluxInner
Flux measured strictly within the fit region (no extrapolation).
Definition: CModel.h:361
std::bitset< N_FLAGS > flags
Array of flags.
Definition: CModel.h:374
ndarray::Array< Scalar const, 1, 1 > amplitudes
Opaque linear parameters in specialized units.
Definition: CModel.h:318
afw::geom::ellipses::Quadrupole ellipse
Best fit half-light ellipse in pixel coordinates.
Definition: CModel.h:315
py::object result
Definition: schema.cc:284
Scalar instFluxErr
Flux uncertainty from the final linear fit.
Definition: CModel.h:360
double Scalar
Typedefs to be used for probability and parameter values.
Definition: common.h:44
boost::shared_ptr< Model > model
Model object that defines the parametrization (defined fully by Control struct)
Definition: CModel.h:306
ndarray::Array< Scalar const, 1, 1 > nonlinear
Opaque nonlinear parameters in specialized units.
Definition: CModel.h:317
CModelResult Result
Typedef to the master Result struct.
Definition: CModel.h:389
Exception to be thrown when a measurement algorithm experiences a known failure mode.
Definition: exceptions.h:48
#define LSST_CONTROL_FIELD(NAME, TYPE, DOC)
A preprocessor macro used to define fields in C++ "control object" structs.
Definition: config.h:36
double weightsMultiplier
"Scale the likelihood by this factor to artificially reweight it w.r.t. the prior." ;
Definition: CModel.h:206
std::string profileName
"Name of the shapelet.RadialProfile that defines the model to fit" ;
Definition: CModel.h:163
bool doRecordTime
"Whether to record the time spent in this stage" ;
Definition: CModel.h:221
Nested control object for CModel that configures one of the three ("initial", "exp", "dev") nonlinear fitting stages.
Definition: CModel.h:138
CModelStageResult dev
Results from the de Vaucouleur (Sersic n=4) fit.
Definition: CModel.h:368
static RadialProfile & get(std::string const &name)
Return a predefined radial profile given its name.
STL class.
Registry and utility class for multi-Gaussian approximations to radial profiles.
Definition: RadialProfile.h:60
CModelStageResult initial
Results from the initial approximate nonlinear fit that feeds the others.
Definition: CModel.h:366
Abstract base class and concrete factories that define multi-shapelet galaxy models.
Definition: Model.h:56
boost::shared_ptr< Prior > prior
Bayesian priors on the parameters (defined fully by Control struct)
Definition: CModel.h:307
std::string priorSource
"One of &#39;FILE&#39;, &#39;LINEAR&#39;, &#39;EMPIRICAL&#39;, or &#39;NONE&#39;, indicating whether the prior should be loaded " "fr...
Definition: CModel.h:169
Main public interface class for CModel algorithm.
Definition: CModel.h:385
std::bitset< N_FLAGS > flags
Array of flags.
Definition: CModel.h:322
A base class for image defects.
A multi-scale shapelet function.
std::string priorName
"Name of the Prior that defines the model to fit (a filename in $MEAS_MODELFIT_DIR/data, " "with no extension), if priorSource=&#39;FILE&#39;. Ignored for forced fitting." ;
Definition: CModel.h:175
#define PTR(...)
Definition: base.h:41
The main control object for CModel, containing parameters for the final linear fit and aggregating th...
Definition: CModel.h:229
Scalar instFluxInner
Flux measured strictly within the fit region (no extrapolation).
Definition: CModel.h:312
afw::geom::ellipses::Quadrupole initialFitRegion
Pixels used in the initial fit.
Definition: CModel.h:370
afw::geom::ellipses::Quadrupole finalFitRegion
Pixels used in the exp, dev, and linear fits.
Definition: CModel.h:371
Scalar instFluxErr
Flux uncertainty from just this stage fit.
Definition: CModel.h:311
table::Schema schema
Definition: Camera.cc:161
Scalar time
Time spent in this fit in seconds.
Definition: CModel.h:314
#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...
Definition: config.h:62
const char * source()
Source function that allows astChannel to source from a Stream.
Definition: Stream.h:224
int nComponents
"Number of Gaussian used to approximate the profile" ;
Definition: CModel.h:187
bool doRecordHistory
"Whether to record the steps the optimizer takes (or just the number, if running as a plugin)" ; ...
Definition: CModel.h:216
SoftenedLinearPriorControl linearPriorConfig
"Configuration for a linear prior, used if priorSource=&#39;LINEAR&#39;." ;
Definition: CModel.h:180
VectorQ moments
Definition: simpleShape.cc:151
Configuration object for Optimizer.
Definition: optimizer.h:187
Scalar fracDev
Fraction of flux from the final linear fit in the de Vaucouleur component (always between 0 and 1)...
Definition: CModel.h:362
Base class for Bayesian priors.
Definition: Prior.h:36
FlagBit
Flags that apply to all four CModel fits or just the last one.
Definition: CModel.h:332
A local mapping between two UnitSystems.
Definition: UnitSystem.h:79
boost::shared_ptr< UnitTransformedLikelihood > likelihood
Object used to evaluate models and compare to data.
Definition: CModel.h:309
SemiEmpiricalPriorControl empiricalPriorConfig
"Configuration for an empirical prior, used if priorSource=&#39;EMPIRICAL&#39;." ;
Definition: CModel.h:185
Scalar instFlux
Flux from the final linear fit.
Definition: CModel.h:359
Base class for all records.
Definition: BaseRecord.h:31
CModelControl Control
Typedef to the master Control struct.
Definition: CModel.h:388
bool usePixelWeights
"Use per-pixel variances as weights in the nonlinear fit (the final linear fit for" " flux never uses...
Definition: CModel.h:200
int maxRadius
"Maximum radius used in approximating profile with Gaussians (0=default for this profile)" ; ...
Definition: CModel.h:193
def measure(mi, x, y, size, statistic, stats)
Definition: fringe.py:357
Reference fit failed, so forced fit will fail as well.
Definition: CModel.h:353
Record class that contains measurements made on a single exposure.
Definition: Source.h:82
FlagBit
Flags for a single CModel stage (note that there are additional flags for the full multi-stage fit) ...
Definition: CModel.h:292
Base class for objective functions for Optimizer.
Definition: optimizer.h:44
ndarray::Array< Scalar const, 1, 1 > fixed
Opaque fixed parameters in specialized units.
Definition: CModel.h:319
afw::table::BaseCatalog history
Trace of the optimizer&#39;s path, if enabled by diagnostic options.
Definition: CModel.h:321
Master result object for CModel, containing results for the final linear fit and three nested CModelS...
Definition: CModel.h:329
boost::shared_ptr< Prior > getPrior() const
Result object for a single nonlinear fitting stage of the CModel algorithm.
Definition: CModel.h:289
OptimizerControl optimizer
"Configuration for how the objective surface is explored. Ignored for forced fitting" ; ...
Definition: CModel.h:211
boost::shared_ptr< Model > getModel() const
A concrete Likelihood class that does not require its parameters and data to be in the same UnitSyste...
Input centroid did not land within the fit region.
Definition: CModel.h:352
Scalar instFlux
Flux measured from just this stage fit.
Definition: CModel.h:310
CModelStageResult exp
Results from the exponential (Sersic n=1) fit.
Definition: CModel.h:367