ReFRACtor
nlls_solver_lm.h
Go to the documentation of this file.
1 #ifndef NLLS_SOLVER_LM_H
2 #define NLLS_SOLVER_LM_H
3 
4 #include <Eigen/Dense>
5 #include <nlls_solver.h>
6 
7 
8 namespace FullPhysics {
9 
10 /******************************************************************
11  \brief Homemade NLLS solver based on Levenberg-Marquardt
12  algorithm.
13 
14  This class is the implementation of J. J. More's
15  version of the Levenberg-Marquardt NLLS solver.
16 *******************************************************************/
17 class NLLSSolverLM :
18  public NLLSSolver {
19 
20 public:
21 
22  /******************************************************************
23  * \brief Provides NLLSSolverLM with very algorithm specific
24  * parameters.
25  *
26  * This class is used to collect and provide a solver of NLLSSolverLM
27  * type with some of the parameters that is not commonly used by
28  * other types of solvers. The constructor sets each parameter to
29  * a reasonable value; however, the trust region radius "tr_rad" is
30  * problem dependent and could optionally be reset to a value more
31  * appropriate for the problem to be solved.
32  *******************************************************************/
33  class Options
34  {
35  public:
37  : min_W(1.0e-9), tr_rad_tol(0.1), tr_rad(1000.0), cr_ratio_tol(0.0001)
38  {}
39 
40  double min_W;
41 
43  double tr_rad_tol;
44 
47  double tr_rad;
48 
49  double cr_ratio_tol;
50  };
55 
56 
57 //-----------------------------------------------------------------------
84 //-----------------------------------------------------------------------
85  NLLSSolverLM( const boost::shared_ptr<NLLSProblem>& p, int max_cost_function_calls,
87  double dx_tol_abs=0.000001, double dx_tol_rel=0.000001,
88  double g_tol_abs=0.000001, double g_tol_rel=6.0555e-06,
89  bool vrbs=false );
90 
91 
92 //-----------------------------------------------------------------------
94 //-----------------------------------------------------------------------
95 
96  virtual ~NLLSSolverLM()
97  {}
98 
99 
100 //-----------------------------------------------------------------------
104 //-----------------------------------------------------------------------
105 
106  virtual void solve();
107 
108 
109 //-----------------------------------------------------------------------
117 //-----------------------------------------------------------------------
118 
119  virtual void print_state(std::ostream &ostr=std::cout);
120 
121 
122 //-----------------------------------------------------------------------
124 //-----------------------------------------------------------------------
125 
126  virtual void print(std::ostream& Os) const
127  { Os << "NLLSSolverLM"; }
128 
129 
130 protected:
131 
132  virtual status_t iterate();
133 
134  static status_t test_dx_rel( const Eigen::VectorXd& dx, const Eigen::VectorXd& x, double dx_tol_rel );
135 
136  static status_t test_dx_abs( const Eigen::VectorXd& dx, double dx_tol_abs );
137 
138  static status_t test_dx( const Eigen::VectorXd& dx, const Eigen::VectorXd& x, double dx_tol_rel, double dx_tol_abs );
139 
140  static status_t test_grad_rel( const Eigen::VectorXd& g, const Eigen::VectorXd& x, double cost, double g_tol_rel );
141 
142  static status_t test_grad_abs( const Eigen::VectorXd& g, double g_tol_abs );
143 
144  double Dx_tol_abs;
145  double Dx_tol_rel;
146  double G_tol_abs;
147  double G_tol_rel;
148 
150 
151  double CR_ratio;
152  double Lambda;
153 
154  Eigen::VectorXd W;
155 
156  Eigen::VectorXd Dx;
157 
158 }; // class NLLSProblemLM
159 
160 } // namespace FullPhysics
161 
162 #endif
virtual void print_state(std::ostream &ostr=std::cout)
Prints solver state.
virtual void solve()
The method that solves the optimization problem.
static status_t test_dx(const Eigen::VectorXd &dx, const Eigen::VectorXd &x, double dx_tol_rel, double dx_tol_abs)
The base class for the solvers of the Nonlinear-Least-Squares Problem.
Definition: nlls_solver.h:24
static status_t test_dx_abs(const Eigen::VectorXd &dx, double dx_tol_abs)
double tr_rad
A positive number and the initial trust region radius.
status_t
Enum type for the status of the iterative solver.
virtual status_t iterate()
static status_t test_grad_rel(const Eigen::VectorXd &g, const Eigen::VectorXd &x, double cost, double g_tol_rel)
static status_t test_dx_rel(const Eigen::VectorXd &dx, const Eigen::VectorXd &x, double dx_tol_rel)
double tr_rad_tol
A positive number and a tolerance for comparing to the trust region radius.
double min_W
A positive number and the smallest value for the diagonal entries of the diagonal weight matrix...
virtual void print(std::ostream &Os) const
Print description of object.
const Unit g("g", 1e-3 *kg)
Contains classes to abstract away details in various Spurr Radiative Transfer software.
Definition: doxygen_python.h:1
static status_t test_grad_abs(const Eigen::VectorXd &g, double g_tol_abs)
virtual ~NLLSSolverLM()
Destructor.
NLLSSolverLM(const boost::shared_ptr< NLLSProblem > &p, int max_cost_function_calls, const NLLSSolverLM::Options &opt=NLLSSolverLM::Options(), double dx_tol_abs=0.000001, double dx_tol_rel=0.000001, double g_tol_abs=0.000001, double g_tol_rel=6.0555e-06, bool vrbs=false)
Constructor.
double cr_ratio_tol
A positive number and tolerance on the ratio of the actual to the predicted reduction of the value of...

Copyright © 2017, California Institute of Technology.
ALL RIGHTS RESERVED.
U.S. Government Sponsorship acknowledged.
Generated Fri Aug 24 2018 15:44:10