331980151ca7b9f145116ff724585fe82249180f
[alexxy/gromacs.git] / src / tools / expfit.c
1 /*
2  * $Id: expfit.c,v 1.33 2005/08/29 19:39:11 lindahl Exp $
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * Green Red Orange Magenta Azure Cyan Skyblue
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40
41 #include <sysstuff.h>
42 #include <string.h>
43 #include <math.h>
44 #include "typedefs.h"
45 #include "smalloc.h"
46 #include "xvgr.h"
47 #include "futil.h"
48 #include "gstat.h"
49 #include "vec.h"
50 #include "statutil.h"
51 #include "index.h"
52
53 const int  nfp_ffn[effnNR] = { 0, 1, 2, 3, 2, 5, 7, 9, 4, 3};
54
55 const char *s_ffn[effnNR+2] = { NULL, "none", "exp", "aexp", "exp_exp", "vac", 
56                           "exp5", "exp7", "exp9","erffit", NULL, NULL };
57 /* We don't allow errest as a choice on the command line */
58
59 const char *longs_ffn[effnNR] = {
60   "no fit",
61   "y = exp(-x/a1)",
62   "y = a2 exp(-x/a1)",
63   "y = a2 exp(-x/a1) + (1-a2) exp(-x/a3)",
64   "y = exp(-v) (cosh(wv) + 1/w sinh(wv)), v = x/(2 a1), w = sqrt(1 - a2)",
65   "y = a1 exp(-x/a2) +  a3 exp(-x/a4) + a5",
66   "y = a1 exp(-x/a2) +  a3 exp(-x/a4) + a5 exp(-x/a6) + a7",
67   "y = a1 exp(-x/a2) +  a3 exp(-x/a4) + a5 exp(-x/a6) + a7 exp(-x/a8) + a9",
68   "y = 1/2*(a1+a2) - 1/2*(a1-a2)*erf( (x-a3) /a4^2)",
69   "y = a2*ee(a1,x) + (1-a2)*ee(a2,x)"
70 };
71
72 extern gmx_bool mrqmin_new(real x[],real y[],real sig[],int ndata,real a[], 
73                        int ia[],int ma,real **covar,real **alpha,real *chisq, 
74                        void (*funcs)(real, real [], real *, real []), 
75                        real *alamda);
76                        
77 static real myexp(real x,real A,real tau)
78 {
79   if ((A == 0) || (tau == 0))
80     return 0;
81   return A*exp(-x/tau);
82 }
83
84 void erffit (real x, real a[], real *y, real dyda[])
85 {
86 /* Fuction 
87  *      y=(a1+a2)/2-(a1-a2)/2*erf((x-a3)/a4^2)
88  */
89
90 real erfarg;
91 real erfval;
92 real erfarg2;
93 real derf;
94
95  erfarg=(x-a[3])/(a[4]*a[4]);
96         erfarg2=erfarg*erfarg;
97         erfval=gmx_erf(erfarg)/2;
98         derf=M_2_SQRTPI*(a[1]-a[2])/2*exp(-erfarg2)/(a[4]*a[4]);
99         *y=(a[1]+a[2])/2-(a[1]-a[2])*erfval;
100         dyda[1]=1/2-erfval;
101         dyda[2]=1/2+erfval;
102         dyda[3]=derf;
103         dyda[4]=2*derf*erfarg;
104 }       
105         
106
107                    
108 static void exp_one_parm(real x,real a[],real *y,real dyda[])
109 {
110   /* Fit to function 
111    *
112    * y = exp(-x/a1)
113    *
114    */
115    
116   real e1;
117   
118   e1      = exp(-x/a[1]);
119   *y      = e1;
120   dyda[1] = x*e1/sqr(a[1]);
121 }
122
123 static void exp_two_parm(real x,real a[],real *y,real dyda[])
124 {
125   /* Fit to function 
126    *
127    * y = a2 exp(-x/a1)
128    *
129    */
130    
131   real e1;
132   
133   e1      = exp(-x/a[1]);
134   *y      = a[2]*e1;
135   dyda[1] = x*a[2]*e1/sqr(a[1]);
136   dyda[2] = e1;
137 }
138
139 static void exp_3_parm(real x,real a[],real *y,real dyda[])
140 {
141   /* Fit to function 
142    *
143    * y = a2 exp(-x/a1) + (1-a2) exp(-x/a3)
144    *
145    */
146    
147   real e1,e2;
148   
149   e1      = exp(-x/a[1]);
150   e2      = exp(-x/a[3]);
151   *y      = a[2]*e1 + (1-a[2])*e2;
152   dyda[1] = x*a[2]*e1/sqr(a[1]);
153   dyda[2] = e1-e2;
154   dyda[3] = x*(1-a[2])*e2/sqr(a[3]);
155 }
156
157 static void exp_5_parm(real x,real a[],real *y,real dyda[])
158 {
159   /* Fit to function 
160    *
161    * y = a1 exp(-x/a2) + a3 exp(-x/a4) + a5
162    *
163    */
164    
165   real e1,e2;
166
167   e1      = exp(-x/a[2]);
168   e2      = exp(-x/a[4]);
169   *y      = a[1]*e1 + a[3]*e2 + a[5];
170
171   if (debug)
172     fprintf(debug,"exp_5_parm called: x = %10.3f  y = %10.3f\n"
173             "a = ( %8.3f  %8.3f  %8.3f  %8.3f  %8.3f)\n",
174             x,*y,a[1],a[2],a[3],a[4],a[5]);
175   dyda[1] = e1;
176   dyda[2] = x*e1/sqr(a[2]);
177   dyda[3] = e2;
178   dyda[4] = x*e2/sqr(a[4]);
179   dyda[5] = 0;
180 }
181
182 static void exp_7_parm(real x,real a[],real *y,real dyda[])
183 {
184   /* Fit to function 
185    *
186    * y = a1 exp(-x/a2) + a3 exp(-x/a4) + a5 exp(-x/a6) + a7
187    *
188    */
189    
190   real e1,e2,e3;
191   
192   e1      = exp(-x/a[2]);
193   e2      = exp(-x/a[4]);
194   e3      = exp(-x/a[6]);
195   *y      = a[1]*e1 + a[3]*e2 + a[5]*e3 + a[7];
196
197   dyda[1] = e1;
198   dyda[2] = x*e1/sqr(a[2]);
199   dyda[3] = e2;
200   dyda[4] = x*e2/sqr(a[4]);
201   dyda[5] = e3;
202   dyda[6] = x*e3/sqr(a[6]);
203   dyda[7] = 0;
204 }
205
206 static void exp_9_parm(real x,real a[],real *y,real dyda[])
207 {
208   /* Fit to function 
209    *
210    * y = a1 exp(-x/a2) + a3 exp(-x/a4) + a5 exp(-x/a6) + a7
211    *
212    */
213    
214   real e1,e2,e3,e4;
215   
216   e1      = exp(-x/a[2]);
217   e2      = exp(-x/a[4]);
218   e3      = exp(-x/a[6]);
219   e4      = exp(-x/a[8]);
220   *y      = a[1]*e1 + a[3]*e2 + a[5]*e3 + a[7]*e4 + a[9];
221
222   dyda[1] = e1;
223   dyda[2] = x*e1/sqr(a[2]);
224   dyda[3] = e2;
225   dyda[4] = x*e2/sqr(a[4]);
226   dyda[5] = e3;
227   dyda[6] = x*e3/sqr(a[6]);
228   dyda[7] = e4;
229   dyda[8] = x*e4/sqr(a[8]);
230   dyda[9] = 0;
231 }
232
233 static void vac_2_parm(real x,real a[],real *y,real dyda[])
234 {
235   /* Fit to function 
236    *
237    * y = 1/2 (1 - 1/w) exp(-(1+w)v) + 1/2 (1 + 1/w) exp(-(1-w)v)
238    *
239    *   = exp(-v) (cosh(wv) + 1/w sinh(wv))
240    *
241    *    v = x/(2 a1)
242    *    w = sqrt(1 - a2)
243    *
244    *    For tranverse current autocorrelation functions:
245    *       a1 = tau
246    *       a2 = 4 tau (eta/rho) k^2
247    *
248    */
249    
250   double v,det,omega,omega2,em,ec,es;
251   
252   v   = x/(2*a[1]);
253   det = 1 - a[2];
254   em  = exp(-v);
255   if (det != 0) {
256     omega2  = fabs(det);
257     omega   = sqrt(omega2);
258     if (det > 0) {
259       ec = em*0.5*(exp(omega*v)+exp(-omega*v));
260       es = em*0.5*(exp(omega*v)-exp(-omega*v))/omega;
261     } else {
262       ec = em*cos(omega*v);
263       es = em*sin(omega*v)/omega;
264     }
265     *y      = ec + es;
266     dyda[2] = (v/det*ec+(v-1/det)*es)/(-2.0);
267     dyda[1] = (1-det)*v/a[1]*es;
268   } else {
269     *y      = (1+v)*em;
270     dyda[2] = -v*v*em*(0.5+v/6);
271     dyda[1] = v*v/a[1]*em;
272   }
273 }
274
275 static void errest_3_parm(real x,real a[],real *y,real dyda[])
276 {
277   real e1,e2,v1,v2;  
278
279   if (a[1])
280     e1 = exp(-x/a[1]) - 1;
281   else
282     e1 = 0;
283   if (a[3])
284     e2 = exp(-x/a[3]) - 1;
285   else
286     e2 = 0;
287
288   if (x > 0) {
289     v1 = 2*a[1]*(e1*a[1]/x + 1);
290     v2 = 2*a[3]*(e2*a[3]/x + 1);
291     *y      = a[2]*v1 + (1-a[2])*v2;
292     dyda[1] = 2*     a[2] *(v1/a[1] + e1);
293     dyda[3] = 2*(1 - a[2])*(v2/a[3] + e2);
294     dyda[2] = (v1 - v2);
295   } else {
296     *y      = 0;
297     dyda[1] = 0;
298     dyda[3] = 0;
299     dyda[2] = 0;
300   }
301 }
302
303 typedef void (*myfitfn)(real x,real a[],real *y,real dyda[]);
304 myfitfn mfitfn[effnNR] = { 
305   exp_one_parm, exp_one_parm, exp_two_parm, exp_3_parm, vac_2_parm,
306   exp_5_parm,   exp_7_parm,   exp_9_parm, erffit,  errest_3_parm
307 };
308
309 real fit_function(int eFitFn,real *parm,real x)
310 {
311   static real y,dum[8];
312
313   mfitfn[eFitFn](x,parm-1,&y,dum);
314
315   return y;
316 }
317
318 /* lmfit_exp supports up to 3 parameter fitting of exponential functions */
319 static gmx_bool lmfit_exp(int nfit,real x[],real y[],real dy[],real ftol,
320                       real parm[],real dparm[],gmx_bool bVerbose,
321                       int eFitFn,int fix)
322 {
323   real chisq,ochisq,alamda;
324   real *a,**covar,**alpha,*dum;
325   gmx_bool bCont;
326   int  i,j,ma,mfit,*lista,*ia;
327
328   if ((eFitFn < 0) || (eFitFn >= effnNR))
329     gmx_fatal(FARGS,"fitfn = %d, should be in 0..%d (%s,%d)",
330                 effnNR-1,eFitFn,__FILE__,__LINE__);
331
332   ma=mfit=nfp_ffn[eFitFn];         /* number of parameters to fit */
333   snew(a,ma+1);
334   snew(covar,ma+1);
335   snew(alpha,ma+1);
336   snew(lista,ma+1);
337   snew(ia,ma+1);
338   snew(dum,ma+1);
339   for(i=1; (i<ma+1); i++) {
340     lista[i] = i;
341     ia[i] = 1; /* fixed bug B.S.S 19/11  */
342     snew(covar[i],ma+1);
343     snew(alpha[i],ma+1);
344   }
345   if (fix) {
346     if (bVerbose)
347       fprintf(stderr,"Will keep parameters fixed during fit procedure: %d\n",
348               fix);
349     for(i=0; i<ma; i++)
350       if (fix & 1<<i)
351         ia[i+1] = 0;
352   }
353   if (debug)
354     fprintf(debug,"%d parameter fit\n",mfit);
355
356   /* Initial params */
357   alamda = -1;    /* Starting value   */
358   chisq  = 1e12;
359   for(i=0; (i<mfit); i++)
360     a[i+1] = parm[i];
361
362   j = 0;      
363   if (bVerbose)
364     fprintf(stderr,"%4s  %10s  %10s  %10s  %10s  %10s %10s\n",
365             "Step","chi^2","Lambda","A1","A2","A3","A4");
366   do {
367     ochisq = chisq;
368     /* mrqmin(x-1,y-1,dy-1,nfit,a,ma,lista,mfit,covar,alpha,
369      *   &chisq,expfn[mfit-1],&alamda)
370      */
371     if (!mrqmin_new(x-1,y-1,dy-1,nfit,a,ia,ma,covar,alpha,&chisq,
372                     mfitfn[eFitFn],&alamda))
373       return FALSE;
374      
375     if (bVerbose) {
376       fprintf(stderr,"%4d  %10.5e  %10.5e  %10.5e",
377               j,chisq,alamda,a[1]);
378       if (mfit > 1)
379         fprintf(stderr,"  %10.5e",a[2]);
380       if (mfit > 2)
381         fprintf(stderr,"  %10.5e",a[3]);
382       if (mfit > 3)
383          fprintf(stderr," %10.5e",a[4]);
384       fprintf(stderr,"\n");
385     }
386     j++;
387     bCont = ((fabs(ochisq - chisq) > fabs(ftol*chisq)) ||
388              ((ochisq == chisq)));
389   } while (bCont && (alamda != 0.0) && (j < 50));
390   if (bVerbose)
391     fprintf(stderr,"\n");
392     
393   /* Now get the covariance matrix out */
394   alamda = 0;
395
396   /*  mrqmin(x-1,y-1,dy-1,nfit,a,ma,lista,mfit,covar,alpha,
397    * &chisq,expfn[mfit-1],&alamda)
398    */
399   if ( !mrqmin_new(x-1,y-1,dy-1,nfit,a,ia,ma,covar,alpha,&chisq,
400                    mfitfn[eFitFn],&alamda))
401     return FALSE;
402
403   for(j=0; (j<mfit); j++) {
404     parm[j]  = a[j+1];
405     dparm[j] = covar[j+1][j+1];
406   }
407
408   for(i=0; (i<ma+1); i++) {
409     sfree(covar[i]);
410     sfree(alpha[i]);
411   }
412   sfree(a);
413   sfree(covar);
414   sfree(alpha);
415   sfree(lista);
416   sfree(dum);
417   
418   return TRUE;
419 }
420
421 real do_lmfit(int ndata,real c1[],real sig[],real dt,real x0[],
422               real begintimefit,real endtimefit,const output_env_t oenv,
423               gmx_bool bVerbose, int eFitFn,real fitparms[],int fix)
424 {
425   FILE *fp;
426   char buf[32];
427
428   int  i,j,nparm,nfitpnts;
429   real integral,ttt;
430   real *parm,*dparm;
431   real *x,*y,*dy;
432   real ftol = 1e-4;
433
434   nparm = nfp_ffn[eFitFn];
435   if (debug) {
436     fprintf(debug,"There are %d points to fit %d vars!\n",ndata,nparm);
437     fprintf(debug,"Fit to function %d from %g through %g, dt=%g\n",
438             eFitFn,begintimefit,endtimefit,dt);
439   }
440
441   snew(x,ndata);
442   snew(y,ndata);
443   snew(dy,ndata);
444
445   j=0;
446   for(i=0; i<ndata; i++) {
447     ttt = x0 ? x0[i] : dt*i;
448     if (ttt>=begintimefit && ttt<=endtimefit) {
449       x[j] = ttt;
450       y[j] = c1[i];
451
452       /* mrqmin does not like sig to be zero */
453       if (sig[i]<1.0e-7)
454         dy[j]=1.0e-7;
455       else
456         dy[j]=sig[i];
457       if (debug)
458         fprintf(debug,"j= %d, i= %d, x= %g, y= %g, dy= %g\n",
459                 j,i,x[j],y[j],dy[j]);
460       j++;
461     }
462   }
463   nfitpnts = j;
464   integral = 0;
465   if (nfitpnts < nparm) 
466     fprintf(stderr,"Not enough data points for fitting!\n");
467   else {
468     snew(parm,nparm);
469     snew(dparm,nparm);
470     if (fitparms)
471       for(i=0; (i < nparm); i++)
472         parm[i]=fitparms[i];
473     
474     if (!lmfit_exp(nfitpnts,x,y,dy,ftol,parm,dparm,bVerbose,eFitFn,fix))
475       fprintf(stderr,"Fit failed!\n");
476     else if (nparm <= 3) {
477       /* Compute the integral from begintimefit */
478       if (nparm == 3) 
479         integral=(parm[0]*myexp(begintimefit,parm[1],  parm[0]) +
480                   parm[2]*myexp(begintimefit,1-parm[1],parm[2]));
481       else if (nparm == 2)
482         integral=parm[0]*myexp(begintimefit,parm[1],  parm[0]);
483       else if (nparm == 1)
484         integral=parm[0]*myexp(begintimefit,1,  parm[0]);
485       else
486         gmx_fatal(FARGS,"nparm = %d in file %s, line %d",
487                     nparm,__FILE__,__LINE__);
488       
489       /* Generate THE output */
490       if (bVerbose) {
491         fprintf(stderr,"FIT: # points used in fit is: %d\n",nfitpnts);
492         fprintf(stderr,"FIT: %21s%21s%21s\n",
493                 "parm0     ","parm1 (ps)   ","parm2 (ps)    ");
494         fprintf(stderr,"FIT: ------------------------------------------------------------\n");
495         fprintf(stderr,"FIT: %8.3g +/- %8.3g%9.4g +/- %8.3g%8.3g +/- %8.3g\n",
496                 parm[0],dparm[0],parm[1],dparm[1],parm[2],dparm[2]);
497         fprintf(stderr,"FIT: Integral (calc with fitted function) from %g ps to inf. is: %g\n",
498                 begintimefit,integral);
499         
500         sprintf(buf,"test%d.xvg",nfitpnts);
501         fp = xvgropen(buf,"C(t) + Fit to C(t)","Time (ps)","C(t)",oenv);
502         fprintf(fp,"# parm0 = %g, parm1 = %g, parm2 = %g\n",
503                 parm[0],parm[1],parm[2]);
504         for(j=0; (j<nfitpnts); j++) {
505           ttt = x0 ? x0[j] : dt*j;
506           fprintf(fp,"%10.5e  %10.5e  %10.5e\n",
507                   ttt,c1[j],fit_function(eFitFn,parm,ttt));
508         }
509         xvgrclose(fp);
510       }
511     }
512     if (fitparms)
513     {
514         for(i=0;(i<nparm);i++)
515         {
516             fitparms[i] = parm[i];
517         }
518     }
519     sfree(parm);
520     sfree(dparm);
521   }
522   
523   sfree(x);
524   sfree(y);
525   sfree(dy);
526   
527   return integral;
528 }
529
530 void do_expfit(int ndata,real c1[],real dt,real begintimefit,real endtimefit)
531 {
532   int i,n;
533   real *x,*y,*Dy;
534   real aa,bb,saa,sbb,A,tau,dA,dtau;
535
536   fprintf(stderr,"Will fit data from %g (ps) to %g (ps).\n",
537           begintimefit,endtimefit);
538
539   snew(x,ndata);   /* allocate the maximum necessary space */
540   snew(y,ndata);
541   snew(Dy,ndata);
542   n=0;
543
544   for(i=0; (i<ndata); i++) {
545     if ( (dt*i >= begintimefit) && (dt*i <= endtimefit) ) {
546       x[n]=dt*i;
547       y[n]=c1[i];
548       Dy[n]=0.5;
549       fprintf(stderr,"n= %d, i= %d, x= %g, y= %g\n",n,i,x[n],y[n]);
550       n++;
551     }
552   }
553   fprintf(stderr,"# of data points used in the fit is : %d\n\n",n);
554   expfit(n,x,y,Dy,&aa,&saa,&bb,&sbb);
555
556   A=exp(aa);
557   dA=exp(aa)*saa;
558   tau=-1.0/bb;
559   dtau=sbb/sqr(bb);
560   fprintf(stderr,"Fitted to y=exp(a+bx):\n");
561   fprintf(stderr,"a = %10.5f\t b = %10.5f",aa,bb);
562   fprintf(stderr,"\n");
563   fprintf(stderr,"Fitted to y=Aexp(-x/tau):\n");
564   fprintf(stderr,"A  = %10.5f\t tau  = %10.5f\n",A,tau);
565   fprintf(stderr,"dA = %10.5f\t dtau = %10.5f\n",dA,dtau);
566 }
567
568
569 void expfit(int n, real *x, real *y, real *Dy, real *a, real *sa, 
570             real *b, real *sb)
571 {
572   real *w,*ly,A,SA,B,SB;
573   int  i;
574   real sum,xbar,ybar,Sxx,Sxy,wr2,chi2,gamma,Db;
575   
576 #define ZERO 0.0
577 #define ONE 1.0
578 #define ONEP5 1.5
579 #define TWO 2.0
580   
581 #define sqr(x) ((x)*(x))
582
583   /*allocate memory */
584   snew(w,n);
585   snew(ly,n);
586
587   /* Calculate weights and values of ln(y). */
588   for(i=0;(i<n); i++){
589     w[i]=sqr(y[i]/Dy[i]);
590     ly[i]=log(y[i]);
591   }
592   
593   /* The weighted averages of x and y: xbar and ybar. */
594   sum=ZERO;
595   xbar=ZERO;
596   ybar=ZERO;
597   for(i=0;(i<n);i++){
598     sum+=w[i];
599     xbar+=w[i]*x[i];
600     ybar+=w[i]*ly[i];
601   }
602   xbar/=sum;
603   ybar/=sum;
604   
605   /* The centered product sums Sxx and Sxy, and hence A and B. */
606   Sxx=ZERO;
607   Sxy=ZERO;
608   for(i=0;(i<n);i++){
609     Sxx+=w[i]*sqr(x[i]-xbar);
610     Sxy+=w[i]*(x[i]-xbar)*(ly[i]-ybar);
611   }
612   B=Sxy/Sxx;
613   A=ybar-B*xbar;
614   
615   /* Chi-squared (chi2) and gamma. */
616   chi2=ZERO;
617   gamma=ZERO;
618   for(i=0;(i<n);i++){
619     wr2=w[i]*sqr(ly[i]-A-B*x[i]);
620     chi2+=wr2;
621     gamma+=wr2*(x[i]-xbar);
622   }
623   
624   /* Refined values of A and B. Also SA and SB. */
625   Db=-ONEP5*gamma/Sxx;
626   B+=Db;
627   A-=ONEP5*chi2/sum-xbar*Db;
628   SB=sqrt((chi2/(n-2))/Sxx);
629   SA=SB*sqrt(Sxx/sum+sqr(xbar));
630   *a=A;
631   *b=B;
632   *sa=SA;
633   *sb=SB;
634 }
635
636