LSSTApplications  20.0.0
LSSTDataManagementBasePackage
optimizer.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_optimizer_h_INCLUDED
25 #define LSST_MEAS_MODELFIT_optimizer_h_INCLUDED
26 
27 #include "ndarray.h"
28 
29 #include "lsst/base.h"
30 #include "lsst/pex/config.h"
31 #include "lsst/afw/table/Schema.h"
34 
35 namespace lsst { namespace meas { namespace modelfit {
36 
37 class Likelihood;
38 class Prior;
39 class Optimizer;
40 
45 public:
46 
47  int const dataSize;
48  int const parameterSize;
49 
60  PTR(Likelihood) likelihood,
61  PTR(Prior) prior = PTR(Prior)()
62  );
63 
67  OptimizerObjective(int dataSize_, int parameterSize_) :
68  dataSize(dataSize_), parameterSize(parameterSize_)
69  {}
70 
91  ndarray::Array<Scalar const,2,1> const & parameters,
92  ndarray::Array<Scalar,1,1> const & output
93  ) const;
94 
103  virtual void computeResiduals(
104  ndarray::Array<Scalar const,1,1> const & parameters,
105  ndarray::Array<Scalar,1,1> const & residuals
106  ) const = 0;
107 
122  ndarray::Array<Scalar const,1,1> const & parameters,
123  ndarray::Array<Scalar,2,-2> const & derivatives
124  ) const {
125  return false;
126  }
127 
128 
134  virtual bool hasPrior() const { return false; }
135 
142  virtual Scalar computePrior(ndarray::Array<Scalar const,1,1> const & parameters) const { return 1.0; }
143 
158  virtual void differentiatePrior(
159  ndarray::Array<Scalar const,1,1> const & parameters,
160  ndarray::Array<Scalar,1,1> const & gradient,
161  ndarray::Array<Scalar,2,1> const & hessian
162  ) const {
163  gradient.deep() = 0.0;
164  hessian.deep() = 0.0;
165  }
166 
167  virtual ~OptimizerObjective() {}
168 };
169 
188 public:
190  noSR1Term, bool,
191  "If true, ignore the SR1 update term in the Hessian, resulting in a Levenberg-Marquardt-like method"
192  );
193 
195  skipSR1UpdateThreshold, double,
196  "Skip the SR1 update if |v||s| / (|v||s|) is less than this threshold"
197  );
198 
200  minTrustRadiusThreshold, double,
201  "If the trust radius falls below this threshold, consider the algorithm converged"
202  );
203 
205  gradientThreshold, double,
206  "If the maximum of the gradient falls below this threshold, consider the algorithm converged"
207  );
208 
210  numDiffRelStep, double,
211  "relative step size used for numerical derivatives (added to other steps)"
212  );
213 
215  numDiffAbsStep, double,
216  "absolute step size used for numerical derivatives (added to other steps)"
217  );
218 
220  numDiffTrustRadiusStep, double,
221  "step size (in units of trust radius) used for numerical derivatives (added to relative step)"
222  );
223 
225  stepAcceptThreshold, double,
226  "steps with reduction ratio greater than this are accepted"
227  );
228 
230  trustRegionInitialSize, double,
231  "the initial trust region will be set to this value"
232  );
233 
236  "steps with reduction radio greater than this may increase the trust radius"
237  );
238 
241  "steps with length this fraction of the trust radius may increase the trust radius"
242  );
243 
245  trustRegionGrowFactor, double,
246  "when increase the trust region size, multiply the radius by this factor"
247  );
248 
251  "steps with reduction radio less than this will decrease the trust radius"
252  );
253 
255  trustRegionShrinkFactor, double,
256  "when reducing the trust region size, multiply the radius by this factor"
257  );
258 
261  "value passed as the tolerance to solveTrustRegion"
262  );
263 
265  maxInnerIterations, int,
266  "maximum number of iterations (i.e. function evaluations and trust region subproblems) per step"
267  );
268 
270  maxOuterIterations, int,
271  "maximum number of steps"
272  );
273 
275  doSaveIterations, bool,
276  "whether to save all iterations for debugging purposes"
277  );
278 
280  noSR1Term(false), skipSR1UpdateThreshold(1E-8),
282  gradientThreshold(1E-5),
284  stepAcceptThreshold(0.0),
290  trustRegionShrinkFactor(1.0/3.0),
292  maxInnerIterations(20),
293  maxOuterIterations(500),
294  doSaveIterations(false)
295  {}
296 };
297 
299 public:
300 
303  PTR(Model) model,
304  bool doRecordDerivatives
305  );
306 
308 
309  void apply(
310  int outerIterCount,
311  int innerIterCount,
312  afw::table::BaseCatalog & history,
313  Optimizer const & optimizer
314  ) const;
315 
317  ndarray::Array<Scalar const,1,1> const & nested,
318  Vector & gradient,
319  Matrix & hessian
320  ) const;
321 
323  afw::table::BaseRecord const & record,
324  Vector & gradient,
325  Matrix & hessian
326  ) const;
327 
329  ndarray::Array<Scalar const,1,1> const & nested,
330  ndarray::Array<Scalar,1,1> const & gradient,
331  ndarray::Array<Scalar,2,2> const & hessian
332  ) const;
333 
335  afw::table::BaseRecord const & record,
336  ndarray::Array<Scalar,1,1> const & gradient,
337  ndarray::Array<Scalar,2,2> const & hessian
338  ) const;
339 
341  afw::table::BaseRecord const & record,
342  ndarray::Array<Scalar const,2,1> const & parameters,
343  ndarray::Array<Scalar,1,1> const & output
344  ) const;
345 
354 };
355 
399 class Optimizer {
400 public:
401 
405 
406  enum StateFlags {
414  FAILED_NAN = 0x0080,
424  };
425 
427  PTR(Objective const) objective,
428  ndarray::Array<Scalar const,1,1> const & parameters,
429  Control const & ctrl
430  );
431 
432  PTR(Objective const) getObjective() const { return _objective; }
433 
434  Control const & getControl() const { return _ctrl; }
435 
436  bool step() { return _stepImpl(0); }
437 
438  bool step(HistoryRecorder const & recorder, afw::table::BaseCatalog & history) {
439  return _stepImpl(0, &recorder, &history);
440  }
441 
442  int run() { return _runImpl(); }
443 
444  int run(HistoryRecorder const & recorder, afw::table::BaseCatalog & history) {
445  return _runImpl(&recorder, &history);
446  }
447 
448  int getState() const { return _state; }
449 
450  Scalar getObjectiveValue() const { return _current.objectiveValue; }
451 
452  ndarray::Array<Scalar const,1,1> getParameters() const { return _current.parameters; }
453 
454  ndarray::Array<Scalar const,1,1> getResiduals() const { return _current.residuals; }
455 
456  ndarray::Array<Scalar const,1,1> getGradient() const { return _gradient; }
457 
458  ndarray::Array<Scalar const,2,2> getHessian() const { return _hessian; }
459 
462 
463 private:
464 
465  struct IterationData {
466  Scalar objectiveValue;
467  Scalar priorValue;
468  ndarray::Array<Scalar,1,1> parameters;
469  ndarray::Array<Scalar,1,1> residuals;
470 
471  IterationData(int dataSize, int parameterSize);
472 
473  void swap(IterationData & other);
474  };
475 
477 
478  bool _stepImpl(
479  int outerIterCount,
480  HistoryRecorder const * recorder=NULL,
481  afw::table::BaseCatalog * history=NULL
482  );
483 
484  int _runImpl(HistoryRecorder const * recorder=NULL, afw::table::BaseCatalog * history=NULL);
485 
486  void _computeDerivatives();
487 
488  int _state;
489  PTR(Objective const) _objective;
490  Control _ctrl;
491  double _trustRadius;
492  IterationData _current;
493  IterationData _next;
494  ndarray::Array<Scalar,1,1> _step;
495  ndarray::Array<Scalar,1,1> _gradient;
496  ndarray::Array<Scalar,2,2> _hessian;
497  ndarray::Array<Scalar,2,-2> _residualDerivative;
498  Matrix _sr1b;
499  Vector _sr1v;
500  Vector _sr1jtr;
501 };
502 
519  ndarray::Array<Scalar,1,1> const & x,
520  ndarray::Array<Scalar const,2,1> const & F, ndarray::Array<Scalar const,1,1> const & g,
521  double r, double tolerance
522 );
523 
524 }}} // namespace lsst::meas::modelfit
525 
526 #endif // !LSST_MEAS_MODELFIT_optimizer_h_INCLUDED
schema
table::Schema schema
Definition: Amplifier.cc:115
lsst::meas::modelfit::OptimizerHistoryRecorder::fillObjectiveModelGrid
void fillObjectiveModelGrid(afw::table::BaseRecord const &record, ndarray::Array< Scalar const, 2, 1 > const &parameters, ndarray::Array< Scalar, 1, 1 > const &output) const
nested
table::Key< int > nested
Definition: TransmissionCurve.cc:555
lsst::meas::modelfit::OptimizerControl::skipSR1UpdateThreshold
double skipSR1UpdateThreshold
"Skip the SR1 update if |v||s| / (|v||s|) is less than this threshold" ;
Definition: optimizer.h:197
lsst::meas::modelfit::OptimizerHistoryRecorder::OptimizerHistoryRecorder
OptimizerHistoryRecorder(afw::table::Schema const &schema)
lsst::meas::modelfit::OptimizerControl::doSaveIterations
bool doSaveIterations
"whether to save all iterations for debugging purposes" ;
Definition: optimizer.h:277
lsst::meas::modelfit::Optimizer::getObjective
boost::shared_ptr< Objective const > getObjective() const
Definition: optimizer.h:432
lsst::meas::modelfit::OptimizerHistoryRecorder::OptimizerHistoryRecorder
OptimizerHistoryRecorder(afw::table::Schema &schema, boost::shared_ptr< Model > model, bool doRecordDerivatives)
lsst::meas::modelfit::Optimizer::CONVERGED_TR_SMALL
@ CONVERGED_TR_SMALL
Definition: optimizer.h:408
lsst::meas::modelfit::Optimizer::STATUS_TR_INCREASED
@ STATUS_TR_INCREASED
Definition: optimizer.h:421
lsst::meas::modelfit::OptimizerHistoryRecorder::unpackDerivatives
void unpackDerivatives(afw::table::BaseRecord const &record, Vector &gradient, Matrix &hessian) const
lsst::meas::modelfit::Optimizer::FAILED_NAN
@ FAILED_NAN
Definition: optimizer.h:414
lsst::meas::modelfit::OptimizerHistoryRecorder::unpackDerivatives
void unpackDerivatives(afw::table::BaseRecord const &record, ndarray::Array< Scalar, 1, 1 > const &gradient, ndarray::Array< Scalar, 2, 2 > const &hessian) const
lsst::meas::modelfit::Optimizer::Control
OptimizerControl Control
Definition: optimizer.h:403
lsst::meas::modelfit::OptimizerControl::trustRegionInitialSize
double trustRegionInitialSize
"the initial trust region will be set to this value" ;
Definition: optimizer.h:232
lsst::meas::modelfit::Optimizer::STATUS_STEP_REJECTED
@ STATUS_STEP_REJECTED
Definition: optimizer.h:416
lsst::meas::modelfit::Optimizer::STATUS_STEP_ACCEPTED
@ STATUS_STEP_ACCEPTED
Definition: optimizer.h:417
lsst::meas::modelfit::Prior
Base class for Bayesian priors.
Definition: Prior.h:36
lsst::meas::modelfit::Optimizer::HistoryRecorder
OptimizerHistoryRecorder HistoryRecorder
Definition: optimizer.h:404
lsst::meas::modelfit::Optimizer::removeSR1Term
void removeSR1Term()
Remove the symmetric-rank-1 secant term from the Hessian, making it just (J^T J)
lsst::meas::modelfit::OptimizerControl::numDiffRelStep
double numDiffRelStep
"relative step size used for numerical derivatives (added to other steps)" ;
Definition: optimizer.h:212
lsst::meas::modelfit::Optimizer::getResiduals
ndarray::Array< Scalar const, 1, 1 > getResiduals() const
Definition: optimizer.h:454
lsst::meas::modelfit::OptimizerObjective::fillObjectiveValueGrid
void fillObjectiveValueGrid(ndarray::Array< Scalar const, 2, 1 > const &parameters, ndarray::Array< Scalar, 1, 1 > const &output) const
Evaluate the Objective on a 1-d grid.
lsst::meas::modelfit::Optimizer::getState
int getState() const
Definition: optimizer.h:448
lsst::meas::modelfit::Optimizer::FAILED_EXCEPTION
@ FAILED_EXCEPTION
Definition: optimizer.h:413
lsst::meas::modelfit::Scalar
double Scalar
Typedefs to be used for probability and parameter values.
Definition: common.h:44
lsst::meas::modelfit::OptimizerObjective::OptimizerObjective
OptimizerObjective(int dataSize_, int parameterSize_)
Base class constructor; must be called by all subclasses.
Definition: optimizer.h:67
lsst::meas::modelfit::Model
Abstract base class and concrete factories that define multi-shapelet galaxy models.
Definition: Model.h:56
lsst::afw::table::Schema
Defines the fields and offsets for a table.
Definition: Schema.h:50
lsst::meas::modelfit::Optimizer::getObjectiveValue
Scalar getObjectiveValue() const
Definition: optimizer.h:450
lsst::meas::modelfit::Optimizer::STATUS_TR
@ STATUS_TR
Definition: optimizer.h:422
lsst::meas::modelfit::Optimizer::step
bool step(HistoryRecorder const &recorder, afw::table::BaseCatalog &history)
Definition: optimizer.h:438
lsst::meas::modelfit::Optimizer::FAILED_MAX_OUTER_ITERATIONS
@ FAILED_MAX_OUTER_ITERATIONS
Definition: optimizer.h:411
lsst::meas::modelfit::OptimizerHistoryRecorder::outer
afw::table::Key< int > outer
Definition: optimizer.h:346
lsst::meas::modelfit::Optimizer::STATUS
@ STATUS
Definition: optimizer.h:423
lsst::meas::modelfit::Likelihood
Base class for optimizer/sampler likelihood functions that compute likelihood at a point.
Definition: Likelihood.h:70
lsst::meas::modelfit::OptimizerHistoryRecorder
Definition: optimizer.h:298
lsst::meas::modelfit::OptimizerControl::trustRegionShrinkReductionRatio
double trustRegionShrinkReductionRatio
"steps with reduction radio less than this will decrease the trust radius" ;
Definition: optimizer.h:252
lsst::meas::modelfit::Optimizer::FAILED
@ FAILED
Definition: optimizer.h:415
lsst::meas::modelfit::OptimizerObjective
Base class for objective functions for Optimizer.
Definition: optimizer.h:44
lsst::meas::modelfit::Optimizer::getControl
Control const & getControl() const
Definition: optimizer.h:434
lsst::meas::modelfit::OptimizerObjective::~OptimizerObjective
virtual ~OptimizerObjective()
Definition: optimizer.h:167
lsst::meas::modelfit::OptimizerHistoryRecorder::parameters
ArrayKey parameters
Definition: optimizer.h:352
lsst::meas::modelfit::OptimizerControl
Configuration object for Optimizer.
Definition: optimizer.h:187
lsst::meas::modelfit::OptimizerHistoryRecorder::derivatives
ArrayKey derivatives
Definition: optimizer.h:353
lsst::meas::modelfit::OptimizerControl::maxInnerIterations
int maxInnerIterations
"maximum number of iterations (i.e. function evaluations and trust region subproblems) per step" ;
Definition: optimizer.h:267
lsst::meas::modelfit::Optimizer::Optimizer
Optimizer(boost::shared_ptr< Objective const > objective, ndarray::Array< Scalar const, 1, 1 > const &parameters, Control const &ctrl)
lsst::meas::modelfit::Optimizer::run
int run()
Definition: optimizer.h:442
lsst::meas::modelfit::OptimizerControl::trustRegionShrinkFactor
double trustRegionShrinkFactor
"when reducing the trust region size, multiply the radius by this factor" ;
Definition: optimizer.h:257
lsst::meas::modelfit::OptimizerObjective::computePrior
virtual Scalar computePrior(ndarray::Array< Scalar const, 1, 1 > const &parameters) const
Compute the value of the Bayesian prior for the given parameter vector.
Definition: optimizer.h:142
x
double x
Definition: ChebyshevBoundedField.cc:277
other
ItemVariant const * other
Definition: Schema.cc:56
lsst::afw::table::BaseRecord
Base class for all records.
Definition: BaseRecord.h:31
lsst::meas::modelfit::Optimizer::step
bool step()
Definition: optimizer.h:436
lsst::meas::modelfit::OptimizerControl::stepAcceptThreshold
double stepAcceptThreshold
"steps with reduction ratio greater than this are accepted" ;
Definition: optimizer.h:227
lsst::meas::modelfit::Optimizer::CONVERGED_GRADZERO
@ CONVERGED_GRADZERO
Definition: optimizer.h:407
lsst::meas::modelfit::solveTrustRegion
void solveTrustRegion(ndarray::Array< Scalar, 1, 1 > const &x, ndarray::Array< Scalar const, 2, 1 > const &F, ndarray::Array< Scalar const, 1, 1 > const &g, double r, double tolerance)
Solve a symmetric quadratic matrix equation with a ball constraint.
lsst::meas::modelfit::OptimizerObjective::differentiatePrior
virtual void differentiatePrior(ndarray::Array< Scalar const, 1, 1 > const &parameters, ndarray::Array< Scalar, 1, 1 > const &gradient, ndarray::Array< Scalar, 2, 1 > const &hessian) const
Compute the first and second derivatives of the Bayesian prior with respect to the parameters.
Definition: optimizer.h:158
lsst::afw::table::Key< int >
lsst::meas::modelfit::Vector
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > Vector
Definition: common.h:46
Model.h
lsst::meas::modelfit::OptimizerControl::noSR1Term
bool noSR1Term
"If true, ignore the SR1 update term in the Hessian, resulting in a Levenberg-Marquardt-like method" ...
Definition: optimizer.h:192
LSST_CONTROL_FIELD
#define LSST_CONTROL_FIELD(NAME, TYPE, DOC)
A preprocessor macro used to define fields in C++ "control object" structs.
Definition: config.h:43
lsst::meas::modelfit::OptimizerHistoryRecorder::inner
afw::table::Key< int > inner
Definition: optimizer.h:347
lsst::meas::modelfit::OptimizerObjective::dataSize
int const dataSize
Definition: optimizer.h:47
base.h
lsst::meas::modelfit::Optimizer
A numerical optimizer customized for least-squares problems with Bayesian priors.
Definition: optimizer.h:399
lsst
A base class for image defects.
Definition: imageAlgorithm.dox:1
lsst::meas::modelfit::Optimizer::FAILED_MAX_INNER_ITERATIONS
@ FAILED_MAX_INNER_ITERATIONS
Definition: optimizer.h:410
lsst::meas::modelfit::OptimizerHistoryRecorder::unpackDerivatives
void unpackDerivatives(ndarray::Array< Scalar const, 1, 1 > const &nested, Vector &gradient, Matrix &hessian) const
lsst::meas::modelfit::OptimizerObjective::differentiateResiduals
virtual bool differentiateResiduals(ndarray::Array< Scalar const, 1, 1 > const &parameters, ndarray::Array< Scalar, 2,-2 > const &derivatives) const
Evaluate analytic derivatives of the model or signal that they are not available.
Definition: optimizer.h:121
lsst::meas::modelfit::Matrix
Eigen::Matrix< Scalar, Eigen::Dynamic, Eigen::Dynamic > Matrix
Definition: common.h:45
lsst::meas::modelfit::OptimizerObjective::hasPrior
virtual bool hasPrior() const
Return true if the Objective has a Bayesian prior as well as a likelihood.
Definition: optimizer.h:134
lsst::meas::modelfit::OptimizerControl::OptimizerControl
OptimizerControl()
Definition: optimizer.h:279
lsst::meas::modelfit::Optimizer::getParameters
ndarray::Array< Scalar const, 1, 1 > getParameters() const
Definition: optimizer.h:452
lsst::meas::modelfit::Optimizer::STATUS_STEP
@ STATUS_STEP
Definition: optimizer.h:418
lsst::meas::modelfit::OptimizerObjective::parameterSize
int const parameterSize
Definition: optimizer.h:48
lsst::meas::modelfit::Optimizer::run
int run(HistoryRecorder const &recorder, afw::table::BaseCatalog &history)
Definition: optimizer.h:444
common.h
lsst::meas::modelfit::OptimizerControl::numDiffTrustRadiusStep
double numDiffTrustRadiusStep
"step size (in units of trust radius) used for numerical derivatives (added to relative step)" ;
Definition: optimizer.h:222
lsst::meas::modelfit::Optimizer::getHessian
ndarray::Array< Scalar const, 2, 2 > getHessian() const
Definition: optimizer.h:458
lsst::meas::modelfit::OptimizerControl::maxOuterIterations
int maxOuterIterations
"maximum number of steps" ;
Definition: optimizer.h:272
lsst::meas::modelfit::OptimizerObjective::computeResiduals
virtual void computeResiduals(ndarray::Array< Scalar const, 1, 1 > const &parameters, ndarray::Array< Scalar, 1, 1 > const &residuals) const =0
Evaluate the residuals of the model for a given parameter vector.
lsst::meas::modelfit::OptimizerHistoryRecorder::state
afw::table::Key< int > state
Definition: optimizer.h:348
lsst::meas::modelfit::OptimizerHistoryRecorder::apply
void apply(int outerIterCount, int innerIterCount, afw::table::BaseCatalog &history, Optimizer const &optimizer) const
PTR
#define PTR(...)
Definition: base.h:41
lsst::meas::modelfit::Optimizer::StateFlags
StateFlags
Definition: optimizer.h:406
lsst::meas::modelfit::Optimizer::STATUS_TR_DECREASED
@ STATUS_TR_DECREASED
Definition: optimizer.h:420
lsst::meas::modelfit::Optimizer::Objective
OptimizerObjective Objective
Definition: optimizer.h:402
lsst::meas::modelfit::OptimizerControl::trustRegionGrowStepFraction
double trustRegionGrowStepFraction
"steps with length this fraction of the trust radius may increase the trust radius" ;
Definition: optimizer.h:242
lsst::afw::cameraGeom::swap
void swap(CameraSys &a, CameraSys &b)
Definition: CameraSys.h:157
lsst::meas::modelfit::OptimizerControl::trustRegionGrowFactor
double trustRegionGrowFactor
"when increase the trust region size, multiply the radius by this factor" ;
Definition: optimizer.h:247
lsst::meas::modelfit::OptimizerHistoryRecorder::trust
ScalarKey trust
Definition: optimizer.h:351
lsst::meas::modelfit::OptimizerControl::trustRegionSolverTolerance
double trustRegionSolverTolerance
"value passed as the tolerance to solveTrustRegion" ;
Definition: optimizer.h:262
config.h
lsst::meas::modelfit::OptimizerControl::trustRegionGrowReductionRatio
double trustRegionGrowReductionRatio
"steps with reduction radio greater than this may increase the trust radius" ;
Definition: optimizer.h:237
Schema.h
lsst::meas::modelfit::OptimizerHistoryRecorder::objective
ScalarKey objective
Definition: optimizer.h:349
lsst::meas::modelfit::OptimizerHistoryRecorder::unpackDerivatives
void unpackDerivatives(ndarray::Array< Scalar const, 1, 1 > const &nested, ndarray::Array< Scalar, 1, 1 > const &gradient, ndarray::Array< Scalar, 2, 2 > const &hessian) const
lsst::meas::modelfit::Optimizer::getGradient
ndarray::Array< Scalar const, 1, 1 > getGradient() const
Definition: optimizer.h:456
lsst::meas::modelfit::OptimizerObjective::makeFromLikelihood
static boost::shared_ptr< OptimizerObjective > makeFromLikelihood(boost::shared_ptr< Likelihood > likelihood, boost::shared_ptr< Prior > prior=boost::shared_ptr< Prior >())
Return a concrete Objective object built from a Likelihood and Prior.
lsst::meas::modelfit::Optimizer::STATUS_TR_UNCHANGED
@ STATUS_TR_UNCHANGED
Definition: optimizer.h:419
lsst::meas::modelfit::OptimizerControl::minTrustRadiusThreshold
double minTrustRadiusThreshold
"If the trust radius falls below this threshold, consider the algorithm converged" ;
Definition: optimizer.h:202
lsst::meas::modelfit::OptimizerControl::gradientThreshold
double gradientThreshold
"If the maximum of the gradient falls below this threshold, consider the algorithm converged" ;
Definition: optimizer.h:207
lsst::afw::table::CatalogT< BaseRecord >
lsst::meas::modelfit::OptimizerControl::numDiffAbsStep
double numDiffAbsStep
"absolute step size used for numerical derivatives (added to other steps)" ;
Definition: optimizer.h:217
lsst::meas::modelfit::OptimizerHistoryRecorder::prior
ScalarKey prior
Definition: optimizer.h:350
lsst::meas::modelfit::Optimizer::CONVERGED
@ CONVERGED
Definition: optimizer.h:409
lsst::meas::modelfit::Optimizer::FAILED_MAX_ITERATIONS
@ FAILED_MAX_ITERATIONS
Definition: optimizer.h:412