Isis Developer Reference
NonLinearLSQ.h
Go to the documentation of this file.
1 #ifndef NonLinearLSQ_h
2 #define NonLinearLSQ_h
3 
10 /* SPDX-License-Identifier: CC0-1.0 */
11 
12 #include <iostream>
13 #include <sstream>
14 #include <string>
15 #include <vector>
16 #include <algorithm>
17 
18 #include "tnt_array1d.h"
19 #include "tnt_array1d_utils.h"
20 #include "tnt_array2d.h"
21 #include "tnt_array2d_utils.h"
22 
23 #include <gsl/gsl_vector.h>
24 #include <gsl/gsl_matrix.h>
25 #include <gsl/gsl_blas.h>
26 #include <gsl/gsl_multifit_nlin.h>
27 
28 #include "IException.h"
29 
30 namespace Isis {
31 
42  class NonLinearLSQ {
43  public:
44  typedef TNT::Array1D<double> NLVector;
45  typedef TNT::Array2D<double> NLMatrix;
46 
47  // Constructors and Destructor
48  NonLinearLSQ() : _fitParms(), _uncert(), _nIters(0), _maxIters(50),
49  _status(0), _userTerminated(false), _userMessage() {}
50 
52  virtual ~NonLinearLSQ() { }
53 
54  virtual int nSize() const = 0;
55  virtual int nParms() const = 0;
56 
62  void setMaxIters(int m) { _maxIters = m; }
63 
70  int maxIters() const { return (_maxIters); }
71 
72  virtual NLVector guess() = 0;
73  virtual NLVector f_x(const NLVector &x) = 0;
74  virtual NLMatrix df_x(const NLVector &x) = 0;
75 
76  virtual double absErr() const { return (1.0E-4); }
77  virtual double relErr() const { return (1.0E-4); }
78 
79  int curvefit();
81  inline int status() const { return (_status); }
83  inline bool success() const { return (_status == GSL_SUCCESS); }
85  inline bool success(int status) const { return (status == GSL_SUCCESS); }
87  inline std::string statusstr() const {
88  return (std::string(gsl_strerror(_status)));
89  }
91  inline std::string statusstr(int status) const {
92  return (std::string(gsl_strerror(status)));
93  }
94 
96  virtual int checkIteration(const int Iter, const NLVector &fitcoefs,
97  const NLVector &uncerts, double cplxconj,
98  int Istatus) {
99  return (Istatus);
100  }
101 
103  inline NLVector coefs() const { return (_fitParms); }
105  inline NLVector uncert() const { return (_uncert); }
107  inline int nIterations() const { return (_nIters); }
108 
109  protected:
110  void Terminate(const std::string &message = "");
111  void Abort(const std::string &reason = "");
112 
113  bool doContinue() const { return (!_userTerminated); }
114 
115 
116  private:
117  NLVector _fitParms;
118  NLVector _uncert;
119  int _nIters;
120  int _maxIters;
121  int _status;
122  bool _userTerminated;
123  std::string _userMessage;
124 
125  struct _nlsqPointer {
126  NonLinearLSQ *nlsq;
127  };
128  static int f(const gsl_vector *x, void *params, gsl_vector *fx);
129  static int df(const gsl_vector *x, void *params, gsl_matrix *J);
130  static int fdf(const gsl_vector *x, void *params, gsl_vector *fx,
131  gsl_matrix *J);
132 
133  NLVector gslToNlsq(const gsl_vector *v) const;
134  NLMatrix gslToNlsq(const gsl_matrix *m) const;
135  gsl_vector *NlsqTogsl(const NLVector &v, gsl_vector *gv = 0) const;
136  gsl_matrix *NlsqTogsl(const NLMatrix &m, gsl_matrix *gm = 0) const;
137  NLVector getUncertainty(const gsl_matrix *m) const;
138  };
139 
140 } // namespace Isis
141 #endif
Isis::NonLinearLSQ::checkIteration
virtual int checkIteration(const int Iter, const NLVector &fitcoefs, const NLVector &uncerts, double cplxconj, int Istatus)
Default interation test simply returns input status.
Definition: NonLinearLSQ.h:96
Isis::NonLinearLSQ::NonLinearLSQ
NonLinearLSQ()
Definition: NonLinearLSQ.h:48
Abort
void Abort(int)
Isis::NonLinearLSQ::curvefit
int curvefit()
Definition: NonLinearLSQ.cpp:28
Isis::NonLinearLSQ::Abort
void Abort(const std::string &reason="")
Definition: NonLinearLSQ.cpp:91
Isis::NonLinearLSQ::status
int status() const
Return status of last fit processing.
Definition: NonLinearLSQ.h:81
Isis::NonLinearLSQ::uncert
NLVector uncert() const
Return uncertainties from last fit processing.
Definition: NonLinearLSQ.h:105
Isis::NonLinearLSQ::setMaxIters
void setMaxIters(int m)
Sets the maximum number of iterations.
Definition: NonLinearLSQ.h:62
IString.h
Isis::NonLinearLSQ::coefs
NLVector coefs() const
Return coefficients from last fit processing.
Definition: NonLinearLSQ.h:103
Isis::NonLinearLSQ::absErr
virtual double absErr() const
Definition: NonLinearLSQ.h:76
Isis::NonLinearLSQ::NLVector
TNT::Array1D< double > NLVector
Definition: NonLinearLSQ.h:44
Isis::NonLinearLSQ::nSize
virtual int nSize() const =0
Isis::NonLinearLSQ::success
bool success() const
Determine success from last fit processing.
Definition: NonLinearLSQ.h:83
Isis::NonLinearLSQ::df_x
virtual NLMatrix df_x(const NLVector &x)=0
_FILEINFO_
#define _FILEINFO_
Macro for the filename and line number.
Definition: IException.h:24
Isis::NonLinearLSQ
NonLinearLSQ Computes a fit using a Levenberg-Marquardt algorithm.
Definition: NonLinearLSQ.h:42
Isis::NonLinearLSQ::maxIters
int maxIters() const
Maximum number iterations for valid solution.
Definition: NonLinearLSQ.h:70
Isis::NonLinearLSQ::guess
virtual NLVector guess()=0
Isis::NonLinearLSQ::Terminate
void Terminate(const std::string &message="")
Definition: NonLinearLSQ.cpp:84
Isis::NonLinearLSQ::success
bool success(int status) const
Check status for success of the given condition
Definition: NonLinearLSQ.h:85
NonLinearLSQ.h
Isis::NonLinearLSQ::statusstr
std::string statusstr() const
Return error message pertaining to last fit procesing.
Definition: NonLinearLSQ.h:87
IException.h
std
Namespace for the standard library.
Isis::NonLinearLSQ::NLMatrix
TNT::Array2D< double > NLMatrix
Definition: NonLinearLSQ.h:45
Isis::NonLinearLSQ::statusstr
std::string statusstr(int status) const
Return error message given status condition.
Definition: NonLinearLSQ.h:91
Isis::NonLinearLSQ::doContinue
bool doContinue() const
Definition: NonLinearLSQ.h:113
Isis::E
const double E
Sets some basic constants for use in ISIS programming.
Definition: Constants.h:39
Isis::NonLinearLSQ::~NonLinearLSQ
virtual ~NonLinearLSQ()
Destructor.
Definition: NonLinearLSQ.h:52
Isis::NonLinearLSQ::nIterations
int nIterations() const
Return number of iterations from last fit processing.
Definition: NonLinearLSQ.h:107
Isis::NonLinearLSQ::f_x
virtual NLVector f_x(const NLVector &x)=0
Isis
This is free and unencumbered software released into the public domain.
Definition: Apollo.h:16
Isis::NonLinearLSQ::relErr
virtual double relErr() const
Definition: NonLinearLSQ.h:77
Isis::NonLinearLSQ::nParms
virtual int nParms() const =0