ReFRACtor
nlls_problem.h
Go to the documentation of this file.
1 #ifndef NLLS_PROBLEM_H
2 #define NLLS_PROBLEM_H
3 #include <cost_func_diff.h>
4 
5 namespace FullPhysics {
6 
7 //-----------------------------------------------------------------------
52 //-----------------------------------------------------------------------
53 
54 class NLLSProblem :
55  public CostFuncDiff {
56 
57 public:
58 
59 
60 //-----------------------------------------------------------------------
62 //-----------------------------------------------------------------------
63 
65  : CostFuncDiff()
66  {}
67 
68 
69  virtual ~NLLSProblem() {}
70 
71 
72 //-----------------------------------------------------------------------
74 //-----------------------------------------------------------------------
75 
76  virtual double cost();
77 
78 
79 //-----------------------------------------------------------------------
81 //-----------------------------------------------------------------------
82 
83  virtual blitz::Array<double, 1> gradient();
84 
85 
86 //-----------------------------------------------------------------------
111 //-----------------------------------------------------------------------
112 
113  virtual blitz::Array<double, 1> residual() = 0;
114 
115 
116 //-----------------------------------------------------------------------
130 //-----------------------------------------------------------------------
131 
132  virtual blitz::Array<double, 1> residual_x(const blitz::Array<double, 1>& x)
133  { parameters(x); return residual(); }
134 
135 
136 //-----------------------------------------------------------------------
165 //-----------------------------------------------------------------------
166 
167  virtual blitz::Array<double, 2> jacobian() = 0;
168 
169 
170 //-----------------------------------------------------------------------
184 //-----------------------------------------------------------------------
185 
186  virtual blitz::Array<double, 2> jacobian_x(const blitz::Array<double, 1>& x)
187  { parameters(x); return jacobian(); }
188 
189 
190 //-----------------------------------------------------------------------
213 //-----------------------------------------------------------------------
214 
215  virtual void residual_jacobian(
216  blitz::Array<double, 1>& r, blitz::Array<double, 2>& j);
217 
218 
219 //-----------------------------------------------------------------------
233 //-----------------------------------------------------------------------
234 
235  virtual void residual_jacobian_x(const blitz::Array<double, 1>& x,
236  blitz::Array<double, 1>& r, blitz::Array<double, 2>& j)
237  { parameters(x); residual_jacobian(r,j); }
238 
239 
240 //-----------------------------------------------------------------------
244 //-----------------------------------------------------------------------
245 
246  virtual int num_residual_evaluations() const
247  { return num_cost_evaluations(); }
248 
249 
250 //-----------------------------------------------------------------------
254 //-----------------------------------------------------------------------
255 
256  virtual int num_jacobian_evaluations() const
257  { return num_der1_evaluations(); }
258 
259 
260 //-----------------------------------------------------------------------
267 //-----------------------------------------------------------------------
268 
269  virtual int residual_size() const = 0;
270 
271 
272 //-----------------------------------------------------------------------
274 //-----------------------------------------------------------------------
275 
276  virtual void print(std::ostream& Os) const
277  { Os << "NLLSProblem"; }
278 
279 
280 protected:
281 
282 };
283 }
284 #endif
virtual blitz::Array< double, 1 > parameters() const
Returns the current parameters.
NLLSProblem()
Default Constructor.
Definition: nlls_problem.h:64
virtual int num_der1_evaluations() const
Returns the number of the times gradient has been evaluated.
virtual void residual_jacobian_x(const blitz::Array< double, 1 > &x, blitz::Array< double, 1 > &r, blitz::Array< double, 2 > &j)
The residual and its Jacobian with parameters.
Definition: nlls_problem.h:235
virtual double cost()
Read comments on CostFunc::cost()
Definition: nlls_problem.cc:24
virtual blitz::Array< double, 1 > residual()=0
The residual vector function.
virtual blitz::Array< double, 1 > residual_x(const blitz::Array< double, 1 > &x)
The residual function with parameters.
Definition: nlls_problem.h:132
virtual int num_jacobian_evaluations() const
Returns the number of the times Jacobian has been evaluated.
Definition: nlls_problem.h:256
virtual blitz::Array< double, 2 > jacobian()=0
The Jacobian matrix function.
virtual int num_residual_evaluations() const
Returns the number of the times residual has been evaluated.
Definition: nlls_problem.h:246
The base class for the Non-Linear Least Squares problem.
Definition: nlls_problem.h:54
The base class for all problem classes that implement a cost function and its gradient.
virtual int num_cost_evaluations() const
Returns the number of the times cost has been evaluated.
Definition: cost_func.h:99
virtual blitz::Array< double, 1 > gradient()
Read comments on CostFuncDiff::gradient()
Definition: nlls_problem.cc:29
virtual blitz::Array< double, 2 > jacobian_x(const blitz::Array< double, 1 > &x)
The Jacobian function with parameters.
Definition: nlls_problem.h:186
Contains classes to abstract away details in various Spurr Radiative Transfer software.
Definition: doxygen_python.h:1
virtual void residual_jacobian(blitz::Array< double, 1 > &r, blitz::Array< double, 2 > &j)
The residual function and its Jacobian together.
Definition: nlls_problem.cc:18
virtual void print(std::ostream &Os) const
Prints description of object.
Definition: nlls_problem.h:276
virtual int residual_size() const =0
The size of the residual returned by residual()

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