18#include "tnt_array1d.h" 
   19#include "tnt_array1d_utils.h" 
   20#include "tnt_array2d.h" 
   21#include "tnt_array2d_utils.h" 
   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> 
   48      NonLinearLSQ() : _fitParms(), _uncert(), _nIters(0), _maxIters(50),
 
   49                       _status(0), _userTerminated(false), _userMessage() {}
 
   76      virtual double absErr()
 const { 
return (1.0E-4); }
 
   77      virtual double relErr()
 const { 
return (1.0E-4); }
 
   81      inline int status()
 const { 
return (_status); }
 
   83      inline bool success()
 const { 
return (_status == GSL_SUCCESS); }
 
   88        return (std::string(gsl_strerror(_status)));
 
   92        return (std::string(gsl_strerror(
status)));
 
   97                                 const NLVector &uncerts, 
double cplxconj,
 
  110      void Terminate(
const std::string &message = 
"");
 
  111      void Abort(
const std::string &reason = 
"");
 
  122      bool        _userTerminated;
 
  123      std::string _userMessage;
 
  125      struct _nlsqPointer {
 
  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,
 
  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;
 
NonLinearLSQ Computes a fit using a Levenberg-Marquardt algorithm.
Definition: NonLinearLSQ.h:42
 
NLVector uncert() const
Return uncertainties from last fit processing.
Definition: NonLinearLSQ.h:105
 
virtual int nParms() const =0
 
int curvefit()
Definition: NonLinearLSQ.cpp:28
 
std::string statusstr(int status) const
Return error message given status condition.
Definition: NonLinearLSQ.h:91
 
virtual double absErr() const
Definition: NonLinearLSQ.h:76
 
int status() const
Return status of last fit processing.
Definition: NonLinearLSQ.h:81
 
TNT::Array1D< double > NLVector
Definition: NonLinearLSQ.h:44
 
virtual NLMatrix df_x(const NLVector &x)=0
 
virtual NLVector f_x(const NLVector &x)=0
 
virtual NLVector guess()=0
 
void setMaxIters(int m)
Sets the maximum number of iterations.
Definition: NonLinearLSQ.h:62
 
TNT::Array2D< double > NLMatrix
Definition: NonLinearLSQ.h:45
 
void Terminate(const std::string &message="")
Definition: NonLinearLSQ.cpp:84
 
int nIterations() const
Return number of iterations from last fit processing.
Definition: NonLinearLSQ.h:107
 
bool success() const
Determine success from last fit processing.
Definition: NonLinearLSQ.h:83
 
int maxIters() const
Maximum number iterations for valid solution.
Definition: NonLinearLSQ.h:70
 
std::string statusstr() const
Return error message pertaining to last fit procesing.
Definition: NonLinearLSQ.h:87
 
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
 
virtual double relErr() const
Definition: NonLinearLSQ.h:77
 
bool success(int status) const
Check status for success of the given condition
Definition: NonLinearLSQ.h:85
 
bool doContinue() const
Definition: NonLinearLSQ.h:113
 
virtual ~NonLinearLSQ()
Destructor.
Definition: NonLinearLSQ.h:52
 
NonLinearLSQ()
Definition: NonLinearLSQ.h:48
 
void Abort(const std::string &reason="")
Definition: NonLinearLSQ.cpp:91
 
virtual int nSize() const =0
 
NLVector coefs() const
Return coefficients from last fit processing.
Definition: NonLinearLSQ.h:103
 
This is free and unencumbered software released into the public domain.
Definition: Apollo.h:16