2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "types/simple.h"
48 static void nrerror(const char error_text[], gmx_bool bExit)
50 fprintf(stderr,"Numerical Recipes run-time error...\n");
51 fprintf(stderr,"%s\n",error_text);
53 fprintf(stderr,"...now exiting to system...\n");
58 /* dont use the keyword vector - it will clash with the
59 * altivec extensions used for powerpc processors.
62 static real *rvector(int nl,int nh)
66 v=(real *)malloc((unsigned) (nh-nl+1)*sizeof(real));
67 if (!v) nrerror("allocation failure in rvector()", TRUE);
71 static int *ivector(int nl, int nh)
75 v=(int *)malloc((unsigned) (nh-nl+1)*sizeof(int));
76 if (!v) nrerror("allocation failure in ivector()", TRUE);
81 static real **matrix1(int nrl, int nrh, int ncl, int nch)
86 m=(real **) malloc((unsigned) (nrh-nrl+1)*sizeof(real*));
87 if (!m) nrerror("allocation failure 1 in matrix1()", TRUE);
90 for(i=nrl;i<=nrh;i++) {
91 m[i]=(real *) malloc((unsigned) (nch-ncl+1)*sizeof(real));
92 if (!m[i]) nrerror("allocation failure 2 in matrix1()", TRUE);
98 static double **dmatrix(int nrl, int nrh, int ncl, int nch)
103 m=(double **) malloc((unsigned) (nrh-nrl+1)*sizeof(double*));
104 if (!m) nrerror("allocation failure 1 in dmatrix()", TRUE);
107 for(i=nrl;i<=nrh;i++) {
108 m[i]=(double *) malloc((unsigned) (nch-ncl+1)*sizeof(double));
109 if (!m[i]) nrerror("allocation failure 2 in dmatrix()", TRUE);
115 static int **imatrix1(int nrl, int nrh, int ncl, int nch)
119 m=(int **)malloc((unsigned) (nrh-nrl+1)*sizeof(int*));
120 if (!m) nrerror("allocation failure 1 in imatrix1()", TRUE);
123 for(i=nrl;i<=nrh;i++) {
124 m[i]=(int *)malloc((unsigned) (nch-ncl+1)*sizeof(int));
125 if (!m[i]) nrerror("allocation failure 2 in imatrix1()", TRUE);
133 static real **submatrix(real **a, int oldrl, int oldrh, int oldcl,
134 int newrl, int newcl)
139 m=(real **) malloc((unsigned) (oldrh-oldrl+1)*sizeof(real*));
140 if (!m) nrerror("allocation failure in submatrix()", TRUE);
143 for(i=oldrl,j=newrl;i<=oldrh;i++,j++) m[j]=a[i]+oldcl-newcl;
150 static void free_vector(real *v, int nl)
152 free((char*) (v+nl));
155 static void free_ivector(int *v, int nl)
157 free((char*) (v+nl));
160 static void free_dvector(int *v, int nl)
162 free((char*) (v+nl));
167 static void free_matrix(real **m, int nrl, int nrh, int ncl)
171 for(i=nrh;i>=nrl;i--) free((char*) (m[i]+ncl));
172 free((char*) (m+nrl));
175 static real **convert_matrix(real *a, int nrl, int nrh, int ncl, int nch)
182 m = (real **) malloc((unsigned) (nrow)*sizeof(real*));
183 if (!m) nrerror("allocation failure in convert_matrix()", TRUE);
185 for(i=0,j=nrl;i<=nrow-1;i++,j++) m[j]=a+ncol*i-ncl;
191 static void free_convert_matrix(real **b, int nrl)
193 free((char*) (b+nrl));
196 #define SWAP(a,b) {real temp=(a);(a)=(b);(b)=temp;}
198 static void dump_mat(int n,real **a)
202 for(i=1; (i<=n); i++) {
203 for(j=1; (j<=n); j++)
204 fprintf(stderr," %10.3f",a[i][j]);
205 fprintf(stderr,"\n");
209 gmx_bool gaussj(real **a, int n, real **b, int m)
211 int *indxc,*indxr,*ipiv;
212 int i,icol=0,irow=0,j,k,l,ll;
218 for (j=1;j<=n;j++) ipiv[j]=0;
225 if (fabs(a[j][k]) >= big) {
230 } else if (ipiv[k] > 1) {
231 nrerror("GAUSSJ: Singular Matrix-1", FALSE);
237 for (l=1;l<=n;l++) SWAP(a[irow][l],a[icol][l])
238 for (l=1;l<=m;l++) SWAP(b[irow][l],b[icol][l])
242 if (a[icol][icol] == 0.0) {
243 fprintf(stderr,"irow = %d, icol = %d\n",irow,icol);
245 nrerror("GAUSSJ: Singular Matrix-2", FALSE);
248 pivinv=1.0/a[icol][icol];
250 for (l=1;l<=n;l++) a[icol][l] *= pivinv;
251 for (l=1;l<=m;l++) b[icol][l] *= pivinv;
252 for (ll=1;ll<=n;ll++)
256 for (l=1;l<=n;l++) a[ll][l] -= a[icol][l]*dum;
257 for (l=1;l<=m;l++) b[ll][l] -= b[icol][l]*dum;
261 if (indxr[l] != indxc[l])
263 SWAP(a[k][indxr[l]],a[k][indxc[l]]);
265 free_ivector(ipiv,1);
266 free_ivector(indxr,1);
267 free_ivector(indxc,1);
275 static void covsrt(real **covar, int ma, int lista[], int mfit)
281 for (i=j+1;i<=ma;i++) covar[i][j]=0.0;
283 for (j=i+1;j<=mfit;j++) {
284 if (lista[j] > lista[i])
285 covar[lista[j]][lista[i]]=covar[i][j];
287 covar[lista[i]][lista[j]]=covar[i][j];
290 for (j=1;j<=ma;j++) {
291 covar[1][j]=covar[j][j];
294 covar[lista[1]][lista[1]]=swap;
295 for (j=2;j<=mfit;j++) covar[lista[j]][lista[j]]=covar[1][j];
297 for (i=1;i<=j-1;i++) covar[i][j]=covar[j][i];
300 #define SWAP(a,b) {swap=(a);(a)=(b);(b)=swap;}
302 static void covsrt_new(real **covar,int ma, int ia[], int mfit)
303 /* Expand in storage the covariance matrix covar, so as to take
304 * into account parameters that are being held fixed. (For the
305 * latter, return zero covariances.)
310 for (i=mfit+1;i<=ma;i++)
311 for (j=1;j<=i;j++) covar[i][j]=covar[j][i]=0.0;
313 for (j=ma;j>=1;j--) {
315 for (i=1;i<=ma;i++) SWAP(covar[i][k],covar[i][j])
316 for (i=1;i<=ma;i++) SWAP(covar[k][i],covar[j][i])
323 static void mrqcof(real x[], real y[], real sig[], int ndata, real a[],
324 int ma, int lista[], int mfit,
325 real **alpha, real beta[], real *chisq,
326 void (*funcs)(real,real *,real *,real *))
329 real ymod,wt,sig2i,dy,*dyda;
332 for (j=1;j<=mfit;j++) {
333 for (k=1;k<=j;k++) alpha[j][k]=0.0;
337 for (i=1;i<=ndata;i++) {
338 (*funcs)(x[i],a,&ymod,dyda);
339 sig2i=1.0/(sig[i]*sig[i]);
341 for (j=1;j<=mfit;j++) {
342 wt=dyda[lista[j]]*sig2i;
344 alpha[j][k] += wt*dyda[lista[k]];
347 (*chisq) += dy*dy*sig2i;
349 for (j=2;j<=mfit;j++)
350 for (k=1;k<=j-1;k++) alpha[k][j]=alpha[j][k];
355 gmx_bool mrqmin(real x[], real y[], real sig[], int ndata, real a[],
356 int ma, int lista[], int mfit,
357 real **covar, real **alpha, real *chisq,
358 void (*funcs)(real,real *,real *,real *),
362 static real *da,*atry,**oneda,*beta,ochisq;
365 oneda=matrix1(1,mfit,1,1);
370 for (j=1;j<=ma;j++) {
372 for (k=1;k<=mfit;k++)
373 if (lista[k] == j) ihit++;
377 nrerror("Bad LISTA permutation in MRQMIN-1", FALSE);
382 nrerror("Bad LISTA permutation in MRQMIN-2", FALSE);
386 mrqcof(x,y,sig,ndata,a,ma,lista,mfit,alpha,beta,chisq,funcs);
389 for (j=1;j<=mfit;j++) {
390 for (k=1;k<=mfit;k++) covar[j][k]=alpha[j][k];
391 covar[j][j]=alpha[j][j]*(1.0+(*alamda));
394 if (!gaussj(covar,mfit,oneda,1))
396 for (j=1;j<=mfit;j++)
398 if (*alamda == 0.0) {
399 covsrt(covar,ma,lista,mfit);
403 free_matrix(oneda,1,mfit,1);
406 for (j=1;j<=ma;j++) atry[j]=a[j];
407 for (j=1;j<=mfit;j++)
408 atry[lista[j]] = a[lista[j]]+da[j];
409 mrqcof(x,y,sig,ndata,atry,ma,lista,mfit,covar,da,chisq,funcs);
410 if (*chisq < ochisq) {
413 for (j=1;j<=mfit;j++) {
414 for (k=1;k<=mfit;k++) alpha[j][k]=covar[j][k];
416 a[lista[j]]=atry[lista[j]];
426 gmx_bool mrqmin_new(real x[],real y[],real sig[],int ndata,real a[],
427 int ia[],int ma,real **covar,real **alpha,real *chisq,
428 void (*funcs)(real, real [], real *, real []),
430 /* Levenberg-Marquardt method, attempting to reduce the value Chi^2
431 * of a fit between a set of data points x[1..ndata], y[1..ndata]
432 * with individual standard deviations sig[1..ndata], and a nonlinear
433 * function dependent on ma coefficients a[1..ma]. The input array
434 * ia[1..ma] indicates by nonzero entries those components of a that
435 * should be fitted for, and by zero entries those components that
436 * should be held fixed at their input values. The program returns
437 * current best-fit values for the parameters a[1..ma], and
438 * Chi^2 = chisq. The arrays covar[1..ma][1..ma], alpha[1..ma][1..ma]
439 * are used as working space during most iterations. Supply a routine
440 * funcs(x,a,yfit,dyda,ma) that evaluates the fitting function yfit,
441 * and its derivatives dyda[1..ma] with respect to the fitting
442 * parameters a at x. On the first call provide an initial guess for
443 * the parameters a, and set alamda < 0 for initialization (which then
444 * sets alamda=.001). If a step succeeds chisq becomes smaller and
445 * alamda de-creases by a factor of 10. If a step fails alamda grows by
446 * a factor of 10. You must call this routine repeatedly until
447 * convergence is achieved. Then, make one final call with alamda=0,
448 * so that covar[1..ma][1..ma] returns the covariance matrix, and alpha
449 * the curvature matrix.
450 * (Parameters held fixed will return zero covariances.)
453 void covsrt(real **covar, int ma, int ia[], int mfit);
454 gmx_bool gaussj(real **a, int n, real **b,int m);
455 void mrqcof_new(real x[], real y[], real sig[], int ndata, real a[],
456 int ia[], int ma, real **alpha, real beta[], real *chisq,
457 void (*funcs)(real, real [], real *, real []));
460 static real ochisq,*atry,*beta,*da,**oneda;
462 if (*alamda < 0.0) { /* Initialization. */
466 for (mfit=0,j=1;j<=ma;j++)
468 oneda=matrix1(1,mfit,1,1);
470 mrqcof_new(x,y,sig,ndata,a,ia,ma,alpha,beta,chisq,funcs);
475 for (j=1;j<=mfit;j++) { /* Alter linearized fitting matrix, by augmenting. */
476 for (k=1;k<=mfit;k++)
477 covar[j][k]=alpha[j][k]; /* diagonal elements. */
478 covar[j][j]=alpha[j][j]*(1.0+(*alamda));
481 if (!gaussj(covar,mfit,oneda,1)) /* Matrix solution. */
483 for (j=1;j<=mfit;j++)
485 if (*alamda == 0.0) { /* Once converged, evaluate covariance matrix. */
486 covsrt_new(covar,ma,ia,mfit);
487 free_matrix(oneda,1,mfit,1);
493 for (j=0,l=1;l<=ma;l++) /* Did the trial succeed? */
494 if (ia[l]) atry[l]=a[l]+da[++j];
495 mrqcof_new(x,y,sig,ndata,atry,ia,ma,covar,da,chisq,funcs);
496 if (*chisq < ochisq) {
497 /* Success, accept the new solution. */
500 for (j=1;j<=mfit;j++) {
501 for (k=1;k<=mfit;k++) alpha[j][k]=covar[j][k];
504 for (l=1;l<=ma;l++) a[l]=atry[l];
505 } else { /* Failure, increase alamda and return. */
512 void mrqcof_new(real x[], real y[], real sig[], int ndata, real a[],
513 int ia[], int ma, real **alpha, real beta[], real *chisq,
514 void (*funcs)(real, real [], real *, real[]))
515 /* Used by mrqmin to evaluate the linearized fitting matrix alpha, and
516 * vector beta as in (15.5.8), and calculate Chi^2.
519 int i,j,k,l,m,mfit=0;
520 real ymod,wt,sig2i,dy,*dyda;
525 for (j=1;j<=mfit;j++) { /* Initialize (symmetric) alpha), beta. */
526 for (k=1;k<=j;k++) alpha[j][k]=0.0;
530 for (i=1;i<=ndata;i++) { /* Summation loop over all data. */
531 (*funcs)(x[i],a,&ymod,dyda);
532 sig2i=1.0/(sig[i]*sig[i]);
534 for (j=0,l=1;l<=ma;l++) {
537 for (j++,k=0,m=1;m<=l;m++)
538 if (ia[m]) alpha[j][++k] += wt*dyda[m];
542 *chisq += dy*dy*sig2i; /* And find Chi^2. */
544 for (j=2;j<=mfit;j++) /* Fill in the symmetric side. */
545 for (k=1;k<j;k++) alpha[k][j]=alpha[j][k];