1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
34 * Green Red Orange Magenta Azure Cyan Skyblue
46 #include "gmx_fatal.h"
54 #include "gmx_matrix.h"
55 #include "gmx_statistics.h"
60 /* must correspond to char *avbar_opt[] declared in main() */
61 enum { avbarSEL
, avbarNONE
, avbarSTDDEV
, avbarERROR
, avbar90
, avbarNR
};
63 static void power_fit(int n
,int nset
,real
**val
,real
*t
)
65 real
*x
,*y
,quality
,a
,b
,r
;
76 fprintf(stdout
,"First time is not larger than 0, using index number as time for power fit\n");
81 for(s
=0; s
<nset
; s
++) {
83 for(i
=0; i
<n
&& val
[s
][i
]>=0; i
++)
84 y
[i
] = log(val
[s
][i
]);
86 fprintf(stdout
,"Will power fit up to point %d, since it is not larger than 0\n",i
);
87 lsq_y_ax_b(i
,x
,y
,&a
,&b
,&r
,&quality
);
88 fprintf(stdout
,"Power fit set %3d: error %.3f a %g b %g\n",
89 s
+1,quality
,a
,exp(b
));
96 static real
cosine_content(int nhp
,int n
,real
*y
)
97 /* Assumes n equidistant points */
99 double fac
,cosyint
,yyint
;
105 fac
= M_PI
*nhp
/(n
-1);
110 cosyint
+= cos(fac
*i
)*y
[i
];
114 return 2*cosyint
*cosyint
/(n
*yyint
);
117 static void plot_coscont(const char *ccfile
,int n
,int nset
,real
**val
,
118 const output_env_t oenv
)
124 fp
= xvgropen(ccfile
,"Cosine content","set / half periods","cosine content",
127 for(s
=0; s
<nset
; s
++) {
128 cc
= cosine_content(s
+1,n
,val
[s
]);
129 fprintf(fp
," %d %g\n",s
+1,cc
);
130 fprintf(stdout
,"Cosine content of set %d with %.1f periods: %g\n",
133 fprintf(stdout
,"\n");
138 static void regression_analysis(int n
,bool bXYdy
,
139 real
*x
,int nset
,real
**val
)
141 real S
,chi2
,a
,b
,da
,db
,r
=0;
143 if (bXYdy
|| (nset
== 1))
145 printf("Fitting data to a function f(x) = ax + b\n");
146 printf("Minimizing residual chi2 = Sum_i w_i [f(x_i) - y_i]2\n");
147 printf("Error estimates will be given if w_i (sigma) values are given\n");
148 printf("(use option -xydy).\n\n");
150 lsq_y_ax_b_error(n
,x
,val
[0],val
[1],&a
,&b
,&da
,&db
,&r
,&S
);
152 lsq_y_ax_b(n
,x
,val
[0],&a
,&b
,&r
,&S
);
154 printf("Chi2 = %g\n",chi2
);
155 printf("S (Sqrt(Chi2/(n-2)) = %g\n",S
);
156 printf("Correlation coefficient = %.1f%%\n",100*r
);
159 printf("a = %g +/- %g\n",a
,da
);
160 printf("b = %g +/- %g\n",b
,db
);
163 printf("a = %g\n",a
);
164 printf("b = %g\n",b
);
169 double chi2
,*a
,**xx
,*y
;
174 for(j
=0; (j
<nset
-1); j
++)
179 for(j
=1; (j
<nset
); j
++)
180 xx
[j
-1][i
] = val
[j
][i
];
183 chi2
= multi_regression(NULL
,n
,y
,nset
-1,xx
,a
);
184 printf("Fitting %d data points in %d sets\n",n
,nset
-1);
185 printf("chi2 = %g\n",chi2
);
187 for(i
=0; (i
<nset
-1); i
++)
199 void histogram(const char *distfile
,real binwidth
,int n
, int nset
, real
**val
,
200 const output_env_t oenv
)
206 #if (defined SIZEOF_LONG_LONG_INT) && (SIZEOF_LONG_LONG_INT >= 8)
207 long long int *histo
;
214 for(s
=0; s
<nset
; s
++)
218 else if (val
[s
][i
] > max
)
221 min
= binwidth
*floor(min
/binwidth
);
222 max
= binwidth
*ceil(max
/binwidth
);
227 nbin
= (int)((max
- min
)/binwidth
+ 0.5) + 1;
228 fprintf(stderr
,"Making distributions with %d bins\n",nbin
);
230 fp
= xvgropen(distfile
,"Distribution","","",oenv
);
231 for(s
=0; s
<nset
; s
++) {
232 for(i
=0; i
<nbin
; i
++)
235 histo
[(int)((val
[s
][i
] - min
)/binwidth
+ 0.5)]++;
236 for(i
=0; i
<nbin
; i
++)
237 fprintf(fp
," %g %g\n",min
+i
*binwidth
,(double)histo
[i
]/(n
*binwidth
));
244 static int real_comp(const void *a
,const void *b
)
246 real dif
= *(real
*)a
- *(real
*)b
;
256 static void average(const char *avfile
,int avbar_opt
,
257 int n
, int nset
,real
**val
,real
*t
)
264 fp
= ffopen(avfile
,"w");
265 if ((avbar_opt
== avbarERROR
) && (nset
== 1))
266 avbar_opt
= avbarNONE
;
267 if (avbar_opt
!= avbarNONE
) {
268 if (avbar_opt
== avbar90
) {
270 fprintf(fp
,"@TYPE xydydy\n");
271 edge
= (int)(nset
*0.05+0.5);
272 fprintf(stdout
,"Errorbars: discarding %d points on both sides: %d%%"
273 " interval\n",edge
,(int)(100*(nset
-2*edge
)/nset
+0.5));
275 fprintf(fp
,"@TYPE xydy\n");
280 for(s
=0; s
<nset
; s
++)
283 fprintf(fp
," %g %g",t
[i
],av
);
285 if (avbar_opt
!= avbarNONE
) {
286 if (avbar_opt
== avbar90
) {
287 for(s
=0; s
<nset
; s
++)
289 qsort(tmp
,nset
,sizeof(tmp
[0]),real_comp
);
290 fprintf(fp
," %g %g",tmp
[nset
-1-edge
]-av
,av
-tmp
[edge
]);
292 for(s
=0; s
<nset
; s
++)
293 var
+= sqr(val
[s
][i
]-av
);
294 if (avbar_opt
== avbarSTDDEV
)
295 err
= sqrt(var
/nset
);
297 err
= sqrt(var
/(nset
*(nset
-1)));
298 fprintf(fp
," %g",err
);
305 if (avbar_opt
== avbar90
)
309 static real
anal_ee_inf(real
*parm
,real T
)
311 return sqrt(parm
[1]*2*parm
[0]/T
+parm
[3]*2*parm
[2]/T
);
314 static real
anal_ee(real
*parm
,real T
,real t
)
319 e1
= exp(-t
/parm
[0]);
323 e2
= exp(-t
/parm
[2]);
327 return sqrt(parm
[1]*2*parm
[0]/T
*((e1
- 1)*parm
[0]/t
+ 1) +
328 parm
[3]*2*parm
[2]/T
*((e2
- 1)*parm
[2]/t
+ 1));
331 static void estimate_error(const char *eefile
,int nb_min
,int resol
,int n
,
332 int nset
, double *av
,double *sig
,real
**val
,real dt
,
333 bool bFitAc
,bool bSingleExpFit
,bool bAllowNegLTCorr
,
334 const output_env_t oenv
)
337 int bs
,prev_bs
,nbs
,nb
;
342 real
*tbs
,*ybs
,rtmp
,dens
,*fitsig
,twooe
,tau1_est
,tau_sig
;
348 fprintf(stdout
,"The number of points is smaller than 4, can not make an error estimate\n");
353 fp
= xvgropen(eefile
,"Error estimates",
354 "Block size (time)","Error estimate", oenv
);
355 if (output_env_get_print_xvgr_codes(oenv
))
358 "@ subtitle \"using block averaging, total time %g (%d points)\"\n",
362 xvgr_legend(fp
,2*nset
,leg
,oenv
);
365 spacing
= pow(2,1.0/resol
);
369 for(s
=0; s
<nset
; s
++)
386 blav
+= val
[s
][bs
*i
+j
];
388 var
+= sqr(av
[s
] - blav
/bs
);
397 ybs
[nbs
] = var
/(nb
*(nb
-1.0))*(n
*dt
)/(sig
[s
]*sig
[s
]);
414 for(i
=0; i
<nbs
/2; i
++)
417 tbs
[i
] = tbs
[nbs
-1-i
];
420 ybs
[i
] = ybs
[nbs
-1-i
];
423 /* The initial slope of the normalized ybs^2 is 1.
424 * For a single exponential autocorrelation: ybs(tau1) = 2/e tau1
425 * From this we take our initial guess for tau1.
433 } while (i
< nbs
- 1 &&
434 (ybs
[i
] > ybs
[i
+1] || ybs
[i
] > twooe
*tau1_est
));
438 fprintf(stdout
,"Data set %d has strange time correlations:\n"
439 "the std. error using single points is larger than that of blocks of 2 points\n"
440 "The error estimate might be inaccurate, check the fit\n",
442 /* Use the total time as tau for the fitting weights */
443 tau_sig
= (n
- 1)*dt
;
452 fprintf(debug
,"set %d tau1 estimate %f\n",s
+1,tau1_est
);
455 /* Generate more or less appropriate sigma's,
456 * also taking the density of points into account.
462 dens
= tbs
[1]/tbs
[0] - 1;
466 dens
= tbs
[nbs
-1]/tbs
[nbs
-2] - 1;
470 dens
= 0.5*(tbs
[i
+1]/tbs
[i
-1] - 1);
472 fitsig
[i
] = sqrt((tau_sig
+ tbs
[i
])/dens
);
477 fitparm
[0] = tau1_est
;
479 /* We set the initial guess for tau2
480 * to halfway between tau1_est and the total time (on log scale).
482 fitparm
[2] = sqrt(tau1_est
*(n
-1)*dt
);
483 do_lmfit(nbs
,ybs
,fitsig
,0,tbs
,0,dt
*n
,oenv
,
484 bDebugMode(),effnERREST
,fitparm
,0);
485 fitparm
[3] = 1-fitparm
[1];
487 if (bSingleExpFit
|| fitparm
[0]<0 || fitparm
[2]<0 || fitparm
[1]<0
488 || (fitparm
[1]>1 && !bAllowNegLTCorr
) || fitparm
[2]>(n
-1)*dt
)
492 if (fitparm
[2] > (n
-1)*dt
)
495 "Warning: tau2 is longer than the length of the data (%g)\n"
496 " the statistics might be bad\n",
501 fprintf(stdout
,"a fitted parameter is negative\n");
503 fprintf(stdout
,"invalid fit: e.e. %g a %g tau1 %g tau2 %g\n",
504 sig
[s
]*anal_ee_inf(fitparm
,n
*dt
),
505 fitparm
[1],fitparm
[0],fitparm
[2]);
506 /* Do a fit with tau2 fixed at the total time.
507 * One could also choose any other large value for tau2.
509 fitparm
[0] = tau1_est
;
511 fitparm
[2] = (n
-1)*dt
;
512 fprintf(stderr
,"Will fix tau2 at the total time: %g\n",fitparm
[2]);
513 do_lmfit(nbs
,ybs
,fitsig
,0,tbs
,0,dt
*n
,oenv
,bDebugMode(),
514 effnERREST
,fitparm
,4);
515 fitparm
[3] = 1-fitparm
[1];
517 if (bSingleExpFit
|| fitparm
[0]<0 || fitparm
[1]<0
518 || (fitparm
[1]>1 && !bAllowNegLTCorr
))
520 if (!bSingleExpFit
) {
521 fprintf(stdout
,"a fitted parameter is negative\n");
522 fprintf(stdout
,"invalid fit: e.e. %g a %g tau1 %g tau2 %g\n",
523 sig
[s
]*anal_ee_inf(fitparm
,n
*dt
),
524 fitparm
[1],fitparm
[0],fitparm
[2]);
526 /* Do a single exponential fit */
527 fprintf(stderr
,"Will use a single exponential fit for set %d\n",s
+1);
528 fitparm
[0] = tau1_est
;
531 do_lmfit(nbs
,ybs
,fitsig
,0,tbs
,0,dt
*n
,oenv
,bDebugMode(),
532 effnERREST
,fitparm
,6);
533 fitparm
[3] = 1-fitparm
[1];
536 ee
= sig
[s
]*anal_ee_inf(fitparm
,n
*dt
);
541 fprintf(stdout
,"Set %3d: err.est. %g a %g tau1 %g tau2 %g\n",
543 fprintf(fp
,"@ legend string %d \"av %f\"\n",2*s
,av
[s
]);
544 fprintf(fp
,"@ legend string %d \"ee %6g\"\n",
545 2*s
+1,sig
[s
]*anal_ee_inf(fitparm
,n
*dt
));
548 fprintf(fp
,"%g %g %g\n",tbs
[i
],sig
[s
]*sqrt(ybs
[i
]/(n
*dt
)),
549 sig
[s
]*sqrt(fit_function(effnERREST
,fitparm
,tbs
[i
])/(n
*dt
)));
555 real
*ac
,acint
,ac_fit
[4];
559 ac
[i
] = val
[s
][i
] - av
[s
];
565 low_do_autocorr(NULL
,oenv
,NULL
,n
,1,-1,&ac
,
566 dt
,eacNormal
,1,FALSE
,TRUE
,
567 FALSE
,0,0,effnNONE
,0);
571 /* Integrate ACF only up to fitlen/2 to avoid integrating noise */
573 for(i
=1; i
<=fitlen
/2; i
++)
579 /* Generate more or less appropriate sigma's */
580 for(i
=0; i
<=fitlen
; i
++)
582 fitsig
[i
] = sqrt(acint
+ dt
*i
);
585 ac_fit
[0] = 0.5*acint
;
587 ac_fit
[2] = 10*acint
;
588 do_lmfit(n
/nb_min
,ac
,fitsig
,dt
,0,0,fitlen
*dt
,oenv
,
589 bDebugMode(),effnEXP3
,ac_fit
,0);
590 ac_fit
[3] = 1 - ac_fit
[1];
592 fprintf(stdout
,"Set %3d: ac erest %g a %g tau1 %g tau2 %g\n",
593 s
+1,sig
[s
]*anal_ee_inf(ac_fit
,n
*dt
),
594 ac_fit
[1],ac_fit
[0],ac_fit
[2]);
599 fprintf(fp
,"%g %g\n",tbs
[i
],
600 sig
[s
]*sqrt(fit_function(effnERREST
,ac_fit
,tbs
[i
]))/(n
*dt
));
616 static void luzar_correl(int nn
,real
*time
,int nset
,real
**val
,real temp
,
617 bool bError
,real fit_start
,real smooth_tail_start
,
618 const output_env_t oenv
)
620 const real tol
= 1e-8;
625 please_cite(stdout
,"Spoel2006b");
627 /* Compute negative derivative k(t) = -dc(t)/dt */
630 compute_derivative(nn
,time
,val
[0],kt
);
631 for(j
=0; (j
<nn
); j
++)
635 for(j
=0; (j
<nn
); j
++)
636 d2
+= sqr(kt
[j
] - val
[3][j
]);
637 fprintf(debug
,"RMS difference in derivatives is %g\n",sqrt(d2
/nn
));
639 analyse_corr(nn
,time
,val
[0],val
[2],kt
,NULL
,NULL
,NULL
,fit_start
,
640 temp
,smooth_tail_start
,oenv
);
643 else if (nset
== 6) {
644 analyse_corr(nn
,time
,val
[0],val
[2],val
[4],
645 val
[1],val
[3],val
[5],fit_start
,temp
,smooth_tail_start
,oenv
);
648 printf("Inconsistent input. I need c(t) sigma_c(t) n(t) sigma_n(t) K(t) sigma_K(t)\n");
649 printf("Not doing anything. Sorry.\n");
653 static void filter(real flen
,int n
,int nset
,real
**val
,real dt
,
654 const output_env_t oenv
)
657 double *filt
,sum
,vf
,fluc
,fluctot
;
659 f
= (int)(flen
/(2*dt
));
663 for(i
=1; i
<=f
; i
++) {
664 filt
[i
] = cos(M_PI
*dt
*i
/flen
);
669 fprintf(stdout
,"Will calculate the fluctuation over %d points\n",n
-2*f
);
670 fprintf(stdout
," using a filter of length %g of %d points\n",flen
,2*f
+1);
672 for(s
=0; s
<nset
; s
++) {
674 for(i
=f
; i
<n
-f
; i
++) {
675 vf
= filt
[0]*val
[s
][i
];
677 vf
+= filt
[j
]*(val
[s
][i
-f
]+val
[s
][i
+f
]);
678 fluc
+= sqr(val
[s
][i
] - vf
);
682 fprintf(stdout
,"Set %3d filtered fluctuation: %12.6e\n",s
+1,sqrt(fluc
));
684 fprintf(stdout
,"Overall filtered fluctuation: %12.6e\n",sqrt(fluctot
/nset
));
685 fprintf(stdout
,"\n");
690 static void do_fit(FILE *out
,int n
,bool bYdy
,int ny
,real
*x0
,real
**val
,
691 int npargs
,t_pargs
*ppa
,const output_env_t oenv
)
693 real
*c1
=NULL
,*sig
=NULL
,*fitparm
;
694 real dt
=0,tendfit
,tbeginfit
;
697 efitfn
= get_acffitfn();
698 nparm
= nfp_ffn
[efitfn
];
699 fprintf(out
,"Will fit to the following function:\n");
700 fprintf(out
,"%s\n",longs_ffn
[efitfn
]);
705 fprintf(out
,"Using two columns as y and sigma values\n");
709 if (opt2parg_bSet("-beginfit",npargs
,ppa
)) {
710 tbeginfit
= opt2parg_real("-beginfit",npargs
,ppa
);
712 tbeginfit
= x0
? x0
[0] : 0;
714 if (opt2parg_bSet("-endfit",npargs
,ppa
)) {
715 tendfit
= opt2parg_real("-endfit",npargs
,ppa
);
717 tendfit
= x0
? x0
[ny
-1] : (ny
-1)*dt
;
731 fitparm
[1] = 0.5*c1
[0];
735 fitparm
[0] = fitparm
[2] = 0.5*c1
[0];
741 fitparm
[0] = fitparm
[2] = fitparm
[4] = 0.33*c1
[0];
748 fitparm
[0] = fitparm
[2] = fitparm
[4] = fitparm
[6] = 0.25*c1
[0];
756 fprintf(out
,"Warning: don't know how to initialize the parameters\n");
757 for(i
=0; (i
<nparm
); i
++)
760 fprintf(out
,"Starting parameters:\n");
761 for(i
=0; (i
<nparm
); i
++)
762 fprintf(out
,"a%-2d = %12.5e\n",i
+1,fitparm
[i
]);
763 if (do_lmfit(ny
,c1
,sig
,dt
,x0
,tbeginfit
,tendfit
,
764 oenv
,bDebugMode(),efitfn
,fitparm
,0)) {
765 for(i
=0; (i
<nparm
); i
++)
766 fprintf(out
,"a%-2d = %12.5e\n",i
+1,fitparm
[i
]);
769 fprintf(out
,"No solution was found\n");
773 static void do_ballistic(const char *balFile
, int nData
,
774 real
*t
, real
**val
, int nSet
,
775 real balTime
, int nBalExp
,
777 const output_env_t oenv
)
779 double **ctd
=NULL
, *td
=NULL
;
780 t_gemParams
*GP
= init_gemParams(0, 0, t
, 0, nData
, balTime
, nBalExp
, bDerivative
);
781 static char *leg
[] = {"Ac'(t)"};
785 if (GP
->ballistic
/GP
->tDelta
>= GP
->nExpFit
*2+1)
790 fp
= xvgropen(balFile
, "Hydrogen Bond Autocorrelation","Time (ps)","C'(t)", oenv
);
791 xvgr_legend(fp
,asize(leg
),leg
,oenv
);
793 for (set
=0; set
<nSet
; set
++)
795 snew(ctd
[set
], nData
);
796 for (i
=0; i
<nData
; i
++) {
797 ctd
[set
][i
] = (double)val
[set
][i
];
799 td
[i
] = (double)t
[i
];
802 takeAwayBallistic(ctd
[set
], td
, nData
, GP
->ballistic
, GP
->nExpFit
, GP
->bDt
);
805 for (i
=0; i
<nData
; i
++)
807 fprintf(fp
, " %g",t
[i
]);
808 for (set
=0; set
<nSet
; set
++)
810 fprintf(fp
, " %g", ctd
[set
][i
]);
816 for (set
=0; set
<nSet
; set
++)
822 printf("Number of data points is less than the number of parameters to fit\n."
823 "The system is underdetermined, hence no ballistic term can be found.\n\n");
826 static void do_geminate(const char *gemFile
, int nData
,
827 real
*t
, real
**val
, int nSet
,
828 const real D
, const real logAfterTime
,
829 real rcut
, real balTime
,
830 const output_env_t oenv
)
832 double **ctd
=NULL
, **ctdGem
=NULL
, *td
=NULL
;
833 t_gemParams
*GP
= init_gemParams(rcut
, D
, t
, logAfterTime
,
834 nData
, balTime
, 1, FALSE
);
835 static char *leg
[] = {"Ac\\sgem\\N(t)"};
843 fp
= xvgropen(gemFile
, "Hydrogen Bond Autocorrelation","Time (ps)","C'(t)", oenv
);
844 xvgr_legend(fp
,asize(leg
),leg
,oenv
);
846 for (set
=0; set
<nSet
; set
++)
848 snew(ctd
[set
], nData
);
849 snew(ctdGem
[set
], nData
);
850 for (i
=0; i
<nData
; i
++) {
851 ctd
[set
][i
] = (double)val
[set
][i
];
853 td
[i
] = (double)t
[i
];
855 fitGemRecomb(ctd
[set
], td
, &(ctd
[set
]), nData
, GP
);
858 for (i
=0; i
<nData
; i
++)
860 fprintf(fp
, " %g",t
[i
]);
861 for (set
=0; set
<nSet
; set
++)
863 fprintf(fp
, " %g", ctdGem
[set
][i
]);
868 for (set
=0; set
<nSet
; set
++)
878 int gmx_analyze(int argc
,char *argv
[])
880 static const char *desc
[] = {
881 "g_analyze reads an ascii file and analyzes data sets.",
882 "A line in the input file may start with a time",
883 "(see option [TT]-time[tt]) and any number of y values may follow.",
884 "Multiple sets can also be",
885 "read when they are seperated by & (option [TT]-n[tt]),",
886 "in this case only one y value is read from each line.",
887 "All lines starting with # and @ are skipped.",
888 "All analyses can also be done for the derivative of a set",
889 "(option [TT]-d[tt]).[PAR]",
891 "All options, except for [TT]-av[tt] and [TT]-power[tt] assume that the",
892 "points are equidistant in time.[PAR]",
894 "g_analyze always shows the average and standard deviation of each",
895 "set. For each set it also shows the relative deviation of the third",
896 "and forth cumulant from those of a Gaussian distribution with the same",
897 "standard deviation.[PAR]",
899 "Option [TT]-ac[tt] produces the autocorrelation function(s).[PAR]",
901 "Option [TT]-cc[tt] plots the resemblance of set i with a cosine of",
902 "i/2 periods. The formula is:[BR]"
903 "2 (int0-T y(t) cos(i pi t) dt)^2 / int0-T y(t) y(t) dt[BR]",
904 "This is useful for principal components obtained from covariance",
905 "analysis, since the principal components of random diffusion are",
906 "pure cosines.[PAR]",
908 "Option [TT]-msd[tt] produces the mean square displacement(s).[PAR]",
910 "Option [TT]-dist[tt] produces distribution plot(s).[PAR]",
912 "Option [TT]-av[tt] produces the average over the sets.",
913 "Error bars can be added with the option [TT]-errbar[tt].",
914 "The errorbars can represent the standard deviation, the error",
915 "(assuming the points are independent) or the interval containing",
916 "90% of the points, by discarding 5% of the points at the top and",
919 "Option [TT]-ee[tt] produces error estimates using block averaging.",
920 "A set is divided in a number of blocks and averages are calculated for",
921 "each block. The error for the total average is calculated from",
922 "the variance between averages of the m blocks B_i as follows:",
923 "error^2 = Sum (B_i - <B>)^2 / (m*(m-1)).",
924 "These errors are plotted as a function of the block size.",
925 "Also an analytical block average curve is plotted, assuming",
926 "that the autocorrelation is a sum of two exponentials.",
927 "The analytical curve for the block average is:[BR]",
928 "f(t) = sigma sqrt(2/T ( a (tau1 ((exp(-t/tau1) - 1) tau1/t + 1)) +[BR]",
929 " (1-a) (tau2 ((exp(-t/tau2) - 1) tau2/t + 1)))),[BR]"
930 "where T is the total time.",
931 "a, tau1 and tau2 are obtained by fitting f^2(t) to error^2.",
932 "When the actual block average is very close to the analytical curve,",
933 "the error is sigma*sqrt(2/T (a tau1 + (1-a) tau2)).",
934 "The complete derivation is given in",
935 "B. Hess, J. Chem. Phys. 116:209-217, 2002.[PAR]",
937 "Option [TT]-bal[TT] finds and subtracts the ultrafast \"ballistic\"",
938 "component from a hydrogen bond autocorrelation function by the fitting",
939 "of a sum of exponentials, as described in e.g.",
940 "O. Markovitch, J. Chem. Phys. 129:084505, 2008. The fastest term",
941 "is the one with the most negative coefficient in the exponential,",
942 "or with [TT]-d[TT], the one with most negative time derivative at time 0.",
943 "[]TT-nbalexp[TT] sets the number of exponentials to fit.[PAR]",
945 "Option [TT]-gem[TT] fits bimolecular rate constants ka and kb",
946 "(and optionally kD) to the hydrogen bond autocorrelation function",
947 "according to the reversible geminate recombination model. Removal of",
948 "the ballistic component first is strongly adviced. The model is presented in",
949 "O. Markovitch, J. Chem. Phys. 129:084505, 2008.[PAR]",
951 "Option [TT]-filter[tt] prints the RMS high-frequency fluctuation",
952 "of each set and over all sets with respect to a filtered average.",
953 "The filter is proportional to cos(pi t/len) where t goes from -len/2",
954 "to len/2. len is supplied with the option [TT]-filter[tt].",
955 "This filter reduces oscillations with period len/2 and len by a factor",
956 "of 0.79 and 0.33 respectively.[PAR]",
958 "Option [TT]-g[tt] fits the data to the function given with option",
959 "[TT]-fitfn[tt].[PAR]",
961 "Option [TT]-power[tt] fits the data to b t^a, which is accomplished",
962 "by fitting to a t + b on log-log scale. All points after the first",
963 "zero or negative value are ignored.[PAR]"
965 "Option [TT]-luzar[tt] performs a Luzar & Chandler kinetics analysis",
966 "on output from [TT]g_hbond[tt]. The input file can be taken directly",
967 "from [TT]g_hbond -ac[tt], and then the same result should be produced."
969 static real tb
=-1,te
=-1,frac
=0.5,filtlen
=0,binwidth
=0.1,aver_start
=0;
970 static bool bHaveT
=TRUE
,bDer
=FALSE
,bSubAv
=TRUE
,bAverCorr
=FALSE
,bXYdy
=FALSE
;
971 static bool bEESEF
=FALSE
,bEENLC
=FALSE
,bEeFitAc
=FALSE
,bPower
=FALSE
;
972 static bool bIntegrate
=FALSE
,bRegression
=FALSE
,bLuzar
=FALSE
,bLuzarError
=FALSE
;
973 static int nsets_in
=1,d
=1,nb_min
=4,resol
=10, nBalExp
=4;
974 static real temp
=298.15,fit_start
=1,smooth_tail_start
=-1, logAfterTime
=10, balTime
=0.2, diffusion
=0,rcut
=0.35;
976 /* must correspond to enum avbar* declared at beginning of file */
977 static const char *avbar_opt
[avbarNR
+1] = {
978 NULL
, "none", "stddev", "error", "90", NULL
982 { "-time", FALSE
, etBOOL
, {&bHaveT
},
983 "Expect a time in the input" },
984 { "-b", FALSE
, etREAL
, {&tb
},
985 "First time to read from set" },
986 { "-e", FALSE
, etREAL
, {&te
},
987 "Last time to read from set" },
988 { "-n", FALSE
, etINT
, {&nsets_in
},
989 "Read # sets seperated by &" },
990 { "-d", FALSE
, etBOOL
, {&bDer
},
991 "Use the derivative" },
992 { "-dp", FALSE
, etINT
, {&d
},
993 "HIDDENThe derivative is the difference over # points" },
994 { "-bw", FALSE
, etREAL
, {&binwidth
},
995 "Binwidth for the distribution" },
996 { "-errbar", FALSE
, etENUM
, {avbar_opt
},
997 "Error bars for -av" },
998 { "-integrate",FALSE
,etBOOL
, {&bIntegrate
},
999 "Integrate data function(s) numerically using trapezium rule" },
1000 { "-aver_start",FALSE
, etREAL
, {&aver_start
},
1001 "Start averaging the integral from here" },
1002 { "-xydy", FALSE
, etBOOL
, {&bXYdy
},
1003 "Interpret second data set as error in the y values for integrating" },
1004 { "-regression",FALSE
,etBOOL
,{&bRegression
},
1005 "Perform a linear regression analysis on the data. If -xydy is set a second set will be interpreted as the error bar in the Y value. Otherwise, if multiple data sets are present a multilinear regression will be performed yielding the constant A that minimize chi^2 = (y - A0 x0 - A1 x1 - ... - AN xN)^2 where now Y is the first data set in the input file and xi the others. Do read the information at the option [TT]-time[tt]." },
1006 { "-luzar", FALSE
, etBOOL
, {&bLuzar
},
1007 "Do a Luzar and Chandler analysis on a correlation function and related as produced by g_hbond. When in addition the -xydy flag is given the second and fourth column will be interpreted as errors in c(t) and n(t)." },
1008 { "-temp", FALSE
, etREAL
, {&temp
},
1009 "Temperature for the Luzar hydrogen bonding kinetics analysis" },
1010 { "-fitstart", FALSE
, etREAL
, {&fit_start
},
1011 "Time (ps) from which to start fitting the correlation functions in order to obtain the forward and backward rate constants for HB breaking and formation" },
1012 { "-smooth",FALSE
, etREAL
, {&smooth_tail_start
},
1013 "If >= 0, the tail of the ACF will be smoothed by fitting it to an exponential function: y = A exp(-x/tau)" },
1014 { "-nbmin", FALSE
, etINT
, {&nb_min
},
1015 "HIDDENMinimum number of blocks for block averaging" },
1016 { "-resol", FALSE
, etINT
, {&resol
},
1017 "HIDDENResolution for the block averaging, block size increases with"
1018 " a factor 2^(1/#)" },
1019 { "-eeexpfit", FALSE
, etBOOL
, {&bEESEF
},
1020 "HIDDENAlways use a single exponential fit for the error estimate" },
1021 { "-eenlc", FALSE
, etBOOL
, {&bEENLC
},
1022 "HIDDENAllow a negative long-time correlation" },
1023 { "-eefitac", FALSE
, etBOOL
, {&bEeFitAc
},
1024 "HIDDENAlso plot analytical block average using a autocorrelation fit" },
1025 { "-filter", FALSE
, etREAL
, {&filtlen
},
1026 "Print the high-frequency fluctuation after filtering with a cosine filter of length #" },
1027 { "-power", FALSE
, etBOOL
, {&bPower
},
1028 "Fit data to: b t^a" },
1029 { "-subav", FALSE
, etBOOL
, {&bSubAv
},
1030 "Subtract the average before autocorrelating" },
1031 { "-oneacf", FALSE
, etBOOL
, {&bAverCorr
},
1032 "Calculate one ACF over all sets" },
1033 { "-nbalexp", FALSE
, etINT
, {&nBalExp
},
1034 "HIDDENNumber of exponentials to fit to the ultrafast component" },
1035 { "-baltime", FALSE
, etREAL
, {&balTime
},
1036 "HIDDENTime up to which the ballistic component will be fitted" },
1037 { "-rcut", FALSE
, etREAL
, {&rcut
},
1038 "Cut-off for hydrogen bonds in geminate algorithms" },
1039 { "-gemtype", FALSE
, etENUM
, {gemType
},
1040 "What type of gminate recombination to use"},
1041 { "-D", FALSE
, etREAL
, {&diffusion
},
1042 "The self diffusion coefficient which is used for the reversible geminate recombination model."
1043 "If non-positive, then the diffusion coefficient will be one of the parameters fitted"}
1045 #define NPA asize(pa)
1048 int n
,nlast
,s
,nset
,i
,j
=0;
1049 real
**val
,*t
,dt
,tot
,error
;
1050 double *av
,*sig
,cum1
,cum2
,cum3
,cum4
,db
;
1051 const char *acfile
,*msdfile
,*ccfile
,*distfile
,*avfile
,*eefile
,*balfile
,*gemfile
,*fitfile
;
1055 { efXVG
, "-f", "graph", ffREAD
},
1056 { efXVG
, "-ac", "autocorr", ffOPTWR
},
1057 { efXVG
, "-msd", "msd", ffOPTWR
},
1058 { efXVG
, "-cc", "coscont", ffOPTWR
},
1059 { efXVG
, "-dist", "distr", ffOPTWR
},
1060 { efXVG
, "-av", "average", ffOPTWR
},
1061 { efXVG
, "-ee", "errest", ffOPTWR
},
1062 { efXVG
, "-bal", "ballisitc",ffOPTWR
},
1063 { efXVG
, "-gem", "geminate", ffOPTWR
},
1064 { efLOG
, "-g", "fitlog", ffOPTWR
}
1066 #define NFILE asize(fnm)
1072 ppa
= add_acf_pargs(&npargs
,pa
);
1074 CopyRight(stderr
,argv
[0]);
1075 parse_common_args(&argc
,argv
,PCA_CAN_VIEW
,
1076 NFILE
,fnm
,npargs
,ppa
,asize(desc
),desc
,0,NULL
,&oenv
);
1078 acfile
= opt2fn_null("-ac",NFILE
,fnm
);
1079 msdfile
= opt2fn_null("-msd",NFILE
,fnm
);
1080 ccfile
= opt2fn_null("-cc",NFILE
,fnm
);
1081 distfile
= opt2fn_null("-dist",NFILE
,fnm
);
1082 avfile
= opt2fn_null("-av",NFILE
,fnm
);
1083 eefile
= opt2fn_null("-ee",NFILE
,fnm
);
1084 balfile
= opt2fn_null("-bal",NFILE
,fnm
);
1085 gemfile
= opt2fn_null("-gem",NFILE
,fnm
);
1086 if (opt2parg_bSet("-fitfn",npargs
,ppa
))
1087 fitfile
= opt2fn("-g",NFILE
,fnm
);
1089 fitfile
= opt2fn_null("-g",NFILE
,fnm
);
1091 val
=read_xvg_time(opt2fn("-f",NFILE
,fnm
),bHaveT
,
1092 opt2parg_bSet("-b",npargs
,ppa
),tb
,
1093 opt2parg_bSet("-e",npargs
,ppa
),te
,
1094 nsets_in
,&nset
,&n
,&dt
,&t
);
1095 printf("Read %d sets of %d points, dt = %g\n\n",nset
,n
,dt
);
1098 printf("Calculating the derivative as (f[i+%d]-f[i])/(%d*dt)\n\n",
1101 for(s
=0; s
<nset
; s
++)
1102 for(i
=0; (i
<n
); i
++)
1103 val
[s
][i
] = (val
[s
][i
+d
]-val
[s
][i
])/(d
*dt
);
1107 printf("Calculating the integral using the trapezium rule\n");
1110 sum
= evaluate_integral(n
,t
,val
[0],val
[1],aver_start
,&stddev
);
1111 printf("Integral %10.3f +/- %10.5f\n",sum
,stddev
);
1114 for(s
=0; s
<nset
; s
++) {
1115 sum
= evaluate_integral(n
,t
,val
[s
],NULL
,aver_start
,&stddev
);
1116 printf("Integral %d %10.5f +/- %10.5f\n",s
+1,sum
,stddev
);
1121 out_fit
= ffopen(fitfile
,"w");
1122 if (bXYdy
&& nset
>=2) {
1123 do_fit(out_fit
,0,TRUE
,n
,t
,val
,npargs
,ppa
,oenv
);
1125 for(s
=0; s
<nset
; s
++)
1126 do_fit(out_fit
,s
,FALSE
,n
,t
,val
,npargs
,ppa
,oenv
);
1131 printf(" std. dev. relative deviation of\n");
1132 printf(" standard --------- cumulants from those of\n");
1133 printf("set average deviation sqrt(n-1) a Gaussian distribition\n");
1134 printf(" cum. 3 cum. 4\n");
1137 for(s
=0; (s
<nset
); s
++) {
1142 for(i
=0; (i
<n
); i
++)
1145 for(i
=0; (i
<n
); i
++) {
1146 db
= val
[s
][i
]-cum1
;
1149 cum4
+= db
*db
*db
*db
;
1155 sig
[s
] = sqrt(cum2
);
1157 error
= sqrt(cum2
/(n
-1));
1160 printf("SS%d %13.6e %12.6e %12.6e %6.3f %6.3f\n",
1161 s
+1,av
[s
],sig
[s
],error
,
1162 sig
[s
] ? cum3
/(sig
[s
]*sig
[s
]*sig
[s
]*sqrt(8/M_PI
)) : 0,
1163 sig
[s
] ? cum4
/(sig
[s
]*sig
[s
]*sig
[s
]*sig
[s
]*3)-1 : 0);
1168 filter(filtlen
,n
,nset
,val
,dt
,oenv
);
1171 out
=xvgropen(msdfile
,"Mean square displacement",
1172 "time","MSD (nm\\S2\\N)",oenv
);
1173 nlast
= (int)(n
*frac
);
1174 for(s
=0; s
<nset
; s
++) {
1175 for(j
=0; j
<=nlast
; j
++) {
1177 fprintf(stderr
,"\r%d",j
);
1179 for(i
=0; i
<n
-j
; i
++)
1180 tot
+= sqr(val
[s
][i
]-val
[s
][i
+j
]);
1182 fprintf(out
," %g %8g\n",dt
*j
,tot
);
1188 fprintf(stderr
,"\r%d, time=%g\n",j
-1,(j
-1)*dt
);
1191 plot_coscont(ccfile
,n
,nset
,val
,oenv
);
1194 histogram(distfile
,binwidth
,n
,nset
,val
,oenv
);
1196 average(avfile
,nenum(avbar_opt
),n
,nset
,val
,t
);
1198 estimate_error(eefile
,nb_min
,resol
,n
,nset
,av
,sig
,val
,dt
,
1199 bEeFitAc
,bEESEF
,bEENLC
,oenv
);
1201 do_ballistic(balfile
,n
,t
,val
,nset
,balTime
,nBalExp
,bDer
,oenv
);
1203 do_geminate(gemfile
,n
,t
,val
,nset
,diffusion
,logAfterTime
,rcut
,balTime
,oenv
);
1205 power_fit(n
,nset
,val
,t
);
1208 for(s
=0; s
<nset
; s
++)
1211 do_autocorr(acfile
,oenv
,"Autocorrelation",n
,nset
,val
,dt
,
1212 eacNormal
,bAverCorr
);
1215 regression_analysis(n
,bXYdy
,t
,nset
,val
);
1218 luzar_correl(n
,t
,nset
,val
,temp
,bXYdy
,fit_start
,smooth_tail_start
,oenv
);
1220 view_all(oenv
,NFILE
, fnm
);