15#include <gsl/gsl_errno.h> 
   16#include <gsl/gsl_blas.h> 
   17#include <gsl/gsl_multifit_nlin.h> 
   19#include "NonLinearLSQ.h" 
   20#include "IException.h" 
   28int NonLinearLSQ::curvefit() {
 
   34  _nlsqPointer d = { 
this };
 
   35  gsl_multifit_function_fdf mf;
 
   43  const gsl_multifit_fdfsolver_type *T = gsl_multifit_fdfsolver_lmsder;
 
   44  gsl_multifit_fdfsolver *s = gsl_multifit_fdfsolver_alloc(T, n, p);
 
   47  gsl_vector *x = NlsqTogsl(_fitParms);
 
   48  gsl_matrix *covar = gsl_matrix_alloc(p, p);
 
   49  gsl_multifit_fdfsolver_set(s, &mf, x);
 
   53                  gsl_blas_dnrm2(s->f), GSL_CONTINUE);
 
   56  gsl_matrix *J = gsl_matrix_alloc(n, p);
 
   60    _status = gsl_multifit_fdfsolver_iterate(s);
 
   61    _fitParms = gslToNlsq(s->x);
 
   63    gsl_multifit_fdfsolver_jac(s, J);
 
   64    gsl_multifit_covar(J, 0.0, covar);
 
   65    _uncert = getUncertainty(covar);
 
   67    _status = 
checkIteration(_nIters, _fitParms, _uncert, gsl_blas_dnrm2(s->f),
 
   69    if ( _status  ) { 
break; }
 
   70    if(!doContinue()) { 
break; }
 
   72    _status = gsl_multifit_test_delta(s->dx, s->x, absErr(), relErr());
 
   73  } 
while ((_status == GSL_CONTINUE) && (_nIters < _maxIters));
 
   76  gsl_multifit_fdfsolver_free(s);
 
   77  gsl_matrix_free(covar);
 
   84void NonLinearLSQ::Terminate(
const std::string &message) {
 
   85  _userMessage = message;
 
   86  _userTerminated = 
true;
 
   87  _status = GSL_SUCCESS;
 
   91void NonLinearLSQ::Abort(
const std::string &reason) {
 
   92  _userMessage = reason;
 
   93  _userTerminated = 
true;
 
   94  _status = GSL_FAILURE;
 
   99int NonLinearLSQ::f(
const gsl_vector *x, 
void *params, gsl_vector *fx) {
 
  100  NonLinearLSQ *nlsq = (
static_cast<_nlsqPointer *
> (params))->nlsq;
 
  101  int n = nlsq->nSize();
 
  102  NLVector fxv = nlsq->f_x(nlsq->gslToNlsq(x));
 
  103  for (
int i = 0 ; i < n ; i++ ) {
 
  104    gsl_vector_set(fx, i, fxv[i]);
 
  106  return  (GSL_SUCCESS);
 
  109int NonLinearLSQ::df(
const gsl_vector *x, 
void *params, gsl_matrix *J) {
 
  110  NonLinearLSQ *nlsq = (
static_cast<_nlsqPointer *
> (params))->nlsq;
 
  111  int n = nlsq->nSize();
 
  112  int p = nlsq->nParms();
 
  114  NLMatrix m = nlsq->df_x(nlsq->gslToNlsq(x));
 
  116  for (
int i = 0 ; i < n ; i++ ) {
 
  117    for (
int j = 0 ; j < p ; j++ ) {
 
  118      gsl_matrix_set(J, i, j, m[i][j]);
 
  121  return  (GSL_SUCCESS);
 
  124int NonLinearLSQ::fdf(
const gsl_vector *x, 
void *params, gsl_vector *fx,
 
  128  return  (GSL_SUCCESS);
 
  132NonLinearLSQ::NLVector NonLinearLSQ::getUncertainty(
const gsl_matrix *covar)
 
  134   NLVector unc(covar->size1);
 
  135   for (
size_t i = 0 ; i < covar->size1 ; i++ ) {
 
  136     unc[i] = sqrt(gsl_matrix_get(covar, i, i));
 
  141NonLinearLSQ::NLVector NonLinearLSQ::gslToNlsq(
const gsl_vector *v)
 const {
 
  144  for (
size_t i = 0 ; i < n ; i++) {
 
  145      Nv[i] = gsl_vector_get(v, i);
 
  150NonLinearLSQ::NLMatrix NonLinearLSQ::gslToNlsq(
const gsl_matrix *m)
 const {
 
  151  size_t nrows = m->size1;
 
  152  size_t ncols = m->size2;
 
  153  NLMatrix Nm(nrows, ncols);
 
  154  for (
size_t i = 0 ; i < nrows ; i++) {
 
  155    for (
size_t j = 0 ; j < ncols ; j++) {
 
  156      Nm[i][j] = gsl_matrix_get(m,i,j);
 
  162gsl_vector *NonLinearLSQ::NlsqTogsl(
const NonLinearLSQ::NLVector &v,
 
  163                                    gsl_vector *gv)
 const {
 
  165    gv = gsl_vector_alloc(v.dim());
 
  167  else if (gv->size != (
size_t) v.dim()) {
 
  169    mess << 
"Size of NL vector (" << v.dim() << 
") not same as GSL vector (" 
  174  for (
int i = 0 ; i < v.dim() ; i++) {
 
  175      gsl_vector_set(gv, i, v[i]);
 
  180gsl_matrix *NonLinearLSQ::NlsqTogsl(
const NonLinearLSQ::NLMatrix &m,
 
  181                                    gsl_matrix *gm)
 const {
 
  183    gm = gsl_matrix_alloc(m.dim1(), m.dim2());
 
  185  else if ((gm->size1 != (
size_t) m.dim1()) &&
 
  186           (gm->size2 != (
size_t) m.dim2()) ) {
 
  188    mess << 
"Size of NL matrix (" << m.dim1() << 
"," << m.dim2()
 
  189         << 
") not same as GSL matrix (" << gm->size1 << 
"," << gm->size2
 
  194  for (
int i = 0 ; i < m.dim1() ; i++) {
 
  195    for (
int j = 0 ; j < m.dim2() ; j++) {
 
  196      gsl_matrix_set(gm, i, j, m[i][j]);
 
@ Programmer
This error is for when a programmer made an API call that was illegal.
virtual int checkIteration(const int Iter, const NLVector &fitcoefs, const NLVector &uncerts, double cplxconj, int Istatus)
Default interation test simply returns input status.
This is free and unencumbered software released into the public domain.
Namespace for the standard library.