Isis 3 Programmer Reference
NonLinearLSQ.cpp
1
7/* SPDX-License-Identifier: CC0-1.0 */
8
9#include <string>
10#include <vector>
11#include <numeric>
12#include <iostream>
13#include <sstream>
14
15#include <gsl/gsl_errno.h>
16#include <gsl/gsl_blas.h>
17#include <gsl/gsl_multifit_nlin.h>
18
19#include "NonLinearLSQ.h"
20#include "IException.h"
21#include "IString.h"
22
23using namespace std;
24
25namespace Isis {
26
27
28int NonLinearLSQ::curvefit() {
29
30 size_t n(nSize());
31 size_t p(nParms());
32
33 // Initialize the solver function information
34 _nlsqPointer d = { this };
35 gsl_multifit_function_fdf mf;
36 mf.f = &f;
37 mf.df = &df;
38 mf.fdf = &fdf;
39 mf.n = n;
40 mf.p = p;
41 mf.params = &d;
42
43 const gsl_multifit_fdfsolver_type *T = gsl_multifit_fdfsolver_lmsder;
44 gsl_multifit_fdfsolver *s = gsl_multifit_fdfsolver_alloc(T, n, p);
45
46 _fitParms = guess();
47 gsl_vector *x = NlsqTogsl(_fitParms);
48 gsl_matrix *covar = gsl_matrix_alloc(p, p);
49 gsl_multifit_fdfsolver_set(s, &mf, x);
50
51 _nIters = 0;
52 checkIteration(_nIters, gslToNlsq(s->x), NLVector(p,999.0),
53 gsl_blas_dnrm2(s->f), GSL_CONTINUE);
54
55
56 gsl_matrix *J = gsl_matrix_alloc(n, p);
57 do {
58 _nIters++;
59
60 _status = gsl_multifit_fdfsolver_iterate(s);
61 _fitParms = gslToNlsq(s->x);
62
63 gsl_multifit_fdfsolver_jac(s, J);
64 gsl_multifit_covar(J, 0.0, covar);
65 _uncert = getUncertainty(covar);
66
67 _status = checkIteration(_nIters, _fitParms, _uncert, gsl_blas_dnrm2(s->f),
68 _status);
69 if ( _status ) { break; }
70 if(!doContinue()) { break; }
71
72 _status = gsl_multifit_test_delta(s->dx, s->x, absErr(), relErr());
73 } while ((_status == GSL_CONTINUE) && (_nIters < _maxIters));
74
75 // Clean up
76 gsl_multifit_fdfsolver_free(s);
77 gsl_matrix_free(covar);
78 gsl_matrix_free(J);
79
80 return (_status);
81}
82
83
84void NonLinearLSQ::Terminate(const std::string &message) {
85 _userMessage = message;
86 _userTerminated = true;
87 _status = GSL_SUCCESS;
88 return;
89}
90
91void NonLinearLSQ::Abort(const std::string &reason) {
92 _userMessage = reason;
93 _userTerminated = true;
94 _status = GSL_FAILURE;
95 return;
96}
97
98
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]);
105 }
106 return (GSL_SUCCESS);
107}
108
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();
113
114 NLMatrix m = nlsq->df_x(nlsq->gslToNlsq(x));
115
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]);
119 }
120 }
121 return (GSL_SUCCESS);
122}
123
124int NonLinearLSQ::fdf(const gsl_vector *x, void *params, gsl_vector *fx,
125 gsl_matrix *J) {
126 f(x,params,fx);
127 df(x,params,J);
128 return (GSL_SUCCESS);
129}
130
131
132NonLinearLSQ::NLVector NonLinearLSQ::getUncertainty(const gsl_matrix *covar)
133 const {
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));
137 }
138 return (unc);
139}
140
141NonLinearLSQ::NLVector NonLinearLSQ::gslToNlsq(const gsl_vector *v) const {
142 size_t n = v->size;
143 NLVector Nv(n);
144 for (size_t i = 0 ; i < n ; i++) {
145 Nv[i] = gsl_vector_get(v, i);
146 }
147 return (Nv);
148}
149
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);
157 }
158 }
159 return (Nm);
160}
161
162gsl_vector *NonLinearLSQ::NlsqTogsl(const NonLinearLSQ::NLVector &v,
163 gsl_vector *gv) const {
164 if (gv == 0) {
165 gv = gsl_vector_alloc(v.dim());
166 }
167 else if (gv->size != (size_t) v.dim()) {
168 ostringstream mess;
169 mess << "Size of NL vector (" << v.dim() << ") not same as GSL vector ("
170 << gv->size << ")";
171 throw IException(IException::Programmer, mess.str(), _FILEINFO_);
172 }
173
174 for (int i = 0 ; i < v.dim() ; i++) {
175 gsl_vector_set(gv, i, v[i]);
176 }
177 return (gv);
178}
179
180gsl_matrix *NonLinearLSQ::NlsqTogsl(const NonLinearLSQ::NLMatrix &m,
181 gsl_matrix *gm) const {
182 if (gm == 0) {
183 gm = gsl_matrix_alloc(m.dim1(), m.dim2());
184 }
185 else if ((gm->size1 != (size_t) m.dim1()) &&
186 (gm->size2 != (size_t) m.dim2()) ) {
187 ostringstream mess;
188 mess << "Size of NL matrix (" << m.dim1() << "," << m.dim2()
189 << ") not same as GSL matrix (" << gm->size1 << "," << gm->size2
190 << ")";
191 throw IException(IException::Programmer, mess.str(), _FILEINFO_);
192 }
193
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]);
197 }
198 }
199 return (gm);
200}
201
202} // namespace ISIS
@ Programmer
This error is for when a programmer made an API call that was illegal.
Definition IException.h:146
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.
Definition Apollo.h:16
Namespace for the standard library.