Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxana / autocorr.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include <math.h>
40 #include <stdio.h>
41 #include <string.h>
42
43 #include "gromacs/fileio/xvgr.h"
44 #include "gromacs/gmxana/correl.h"
45 #include "gromacs/gmxana/gstat.h"
46 #include "gromacs/legacyheaders/macros.h"
47 #include "gromacs/legacyheaders/names.h"
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/futil.h"
52 #include "gromacs/utility/smalloc.h"
53
54 #define MODE(x) ((mode & (x)) == (x))
55
56 typedef struct {
57     unsigned long mode;
58     int           nrestart, nout, P, fitfn;
59     gmx_bool      bFour, bNormalize;
60     real          tbeginfit, tendfit;
61 } t_acf;
62
63 static gmx_bool  bACFinit = FALSE;
64 static t_acf     acf;
65
66 enum {
67     enNorm, enCos, enSin
68 };
69
70 int sffn2effn(const char **sffn)
71 {
72     int eFitFn, i;
73
74     eFitFn = 0;
75     for (i = 0; i < effnNR; i++)
76     {
77         if (sffn[i+1] && strcmp(sffn[0], sffn[i+1]) == 0)
78         {
79             eFitFn = i;
80         }
81     }
82
83     return eFitFn;
84 }
85
86 static void low_do_four_core(int nfour, int nframes, real c1[], real cfour[],
87                              int nCos, gmx_bool bPadding)
88 {
89     int  i = 0;
90     real aver, *ans;
91
92     aver = 0.0;
93     switch (nCos)
94     {
95         case enNorm:
96             for (i = 0; (i < nframes); i++)
97             {
98                 aver    += c1[i];
99                 cfour[i] = c1[i];
100             }
101             break;
102         case enCos:
103             for (i = 0; (i < nframes); i++)
104             {
105                 cfour[i] = cos(c1[i]);
106             }
107             break;
108         case enSin:
109             for (i = 0; (i < nframes); i++)
110             {
111                 cfour[i] = sin(c1[i]);
112             }
113             break;
114         default:
115             gmx_fatal(FARGS, "nCos = %d, %s %d", nCos, __FILE__, __LINE__);
116     }
117     for (; (i < nfour); i++)
118     {
119         cfour[i] = 0.0;
120     }
121
122     if (bPadding)
123     {
124         aver /= nframes;
125         /* printf("aver = %g\n",aver); */
126         for (i = 0; (i < nframes); i++)
127         {
128             cfour[i] -= aver;
129         }
130     }
131
132     snew(ans, 2*nfour);
133     correl(cfour-1, cfour-1, nfour, ans-1);
134
135     if (bPadding)
136     {
137         for (i = 0; (i < nfour); i++)
138         {
139             cfour[i] = ans[i]+sqr(aver);
140         }
141     }
142     else
143     {
144         for (i = 0; (i < nfour); i++)
145         {
146             cfour[i] = ans[i];
147         }
148     }
149
150     sfree(ans);
151 }
152
153 static void do_ac_core(int nframes, int nout,
154                        real corr[], real c1[], int nrestart,
155                        unsigned long mode)
156 {
157     int     j, k, j3, jk3, m, n;
158     real    ccc, cth;
159     rvec    xj, xk, rr;
160
161     if (nrestart < 1)
162     {
163         printf("WARNING: setting number of restarts to 1\n");
164         nrestart = 1;
165     }
166     if (debug)
167     {
168         fprintf(debug,
169                 "Starting do_ac_core: nframes=%d, nout=%d, nrestart=%d,mode=%lu\n",
170                 nframes, nout, nrestart, mode);
171     }
172
173     for (j = 0; (j < nout); j++)
174     {
175         corr[j] = 0;
176     }
177
178     /* Loop over starting points. */
179     for (j = 0; (j < nframes); j += nrestart)
180     {
181         j3  = DIM*j;
182
183         /* Loop over the correlation length for this starting point */
184         for (k = 0; (k < nout) && (j+k < nframes); k++)
185         {
186             jk3 = DIM*(j+k);
187
188             /* Switch over possible ACF types.
189              * It might be more efficient to put the loops inside the switch,
190              * but this is more clear, and save development time!
191              */
192             if (MODE(eacNormal))
193             {
194                 corr[k] += c1[j]*c1[j+k];
195             }
196             else if (MODE(eacCos))
197             {
198                 /* Compute the cos (phi(t)-phi(t+dt)) */
199                 corr[k] += cos(c1[j]-c1[j+k]);
200             }
201             else if (MODE(eacIden))
202             {
203                 /* Check equality (phi(t)==phi(t+dt)) */
204                 corr[k] += (c1[j] == c1[j+k]) ? 1 : 0;
205             }
206             else if (MODE(eacP1) || MODE(eacP2) || MODE(eacP3))
207             {
208                 for (m = 0; (m < DIM); m++)
209                 {
210                     xj[m] = c1[j3+m];
211                     xk[m] = c1[jk3+m];
212                 }
213                 cth = cos_angle(xj, xk);
214
215                 if (cth-1.0 > 1.0e-15)
216                 {
217                     printf("j: %d, k: %d, xj:(%g,%g,%g), xk:(%g,%g,%g)\n",
218                            j, k, xj[XX], xj[YY], xj[ZZ], xk[XX], xk[YY], xk[ZZ]);
219                 }
220
221                 corr[k] += LegendreP(cth, mode); /* 1.5*cth*cth-0.5; */
222             }
223             else if (MODE(eacRcross))
224             {
225                 for (m = 0; (m < DIM); m++)
226                 {
227                     xj[m] = c1[j3+m];
228                     xk[m] = c1[jk3+m];
229                 }
230                 cprod(xj, xk, rr);
231
232                 corr[k] += iprod(rr, rr);
233             }
234             else if (MODE(eacVector))
235             {
236                 for (m = 0; (m < DIM); m++)
237                 {
238                     xj[m] = c1[j3+m];
239                     xk[m] = c1[jk3+m];
240                 }
241                 ccc = iprod(xj, xk);
242
243                 corr[k] += ccc;
244             }
245             else
246             {
247                 gmx_fatal(FARGS, "\nInvalid mode (%d) in do_ac_core", mode);
248             }
249         }
250     }
251     /* Correct for the number of points and copy results to the data array */
252     for (j = 0; (j < nout); j++)
253     {
254         n     = (nframes-j+(nrestart-1))/nrestart;
255         c1[j] = corr[j]/n;
256     }
257 }
258
259 void normalize_acf(int nout, real corr[])
260 {
261     int  j;
262     real c0;
263
264     if (debug)
265     {
266         fprintf(debug, "Before normalization\n");
267         for (j = 0; (j < nout); j++)
268         {
269             fprintf(debug, "%5d  %10f\n", j, corr[j]);
270         }
271     }
272
273     /* Normalisation makes that c[0] = 1.0 and that other points are scaled
274      * accordingly.
275      */
276     if (corr[0] == 0.0)
277     {
278         c0 = 1.0;
279     }
280     else
281     {
282         c0 = 1.0/corr[0];
283     }
284     for (j = 0; (j < nout); j++)
285     {
286         corr[j] *= c0;
287     }
288
289     if (debug)
290     {
291         fprintf(debug, "After normalization\n");
292         for (j = 0; (j < nout); j++)
293         {
294             fprintf(debug, "%5d  %10f\n", j, corr[j]);
295         }
296     }
297 }
298
299 void average_acf(gmx_bool bVerbose, int n, int nitem, real **c1)
300 {
301     real c0;
302     int  i, j;
303
304     if (bVerbose)
305     {
306         printf("Averaging correlation functions\n");
307     }
308
309     for (j = 0; (j < n); j++)
310     {
311         c0 = 0;
312         for (i = 0; (i < nitem); i++)
313         {
314             c0 += c1[i][j];
315         }
316         c1[0][j] = c0/nitem;
317     }
318 }
319
320 void norm_and_scale_vectors(int nframes, real c1[], real scale)
321 {
322     int   j, m;
323     real *rij;
324
325     for (j = 0; (j < nframes); j++)
326     {
327         rij = &(c1[j*DIM]);
328         unitv(rij, rij);
329         for (m = 0; (m < DIM); m++)
330         {
331             rij[m] *= scale;
332         }
333     }
334 }
335
336 void dump_tmp(char *s, int n, real c[])
337 {
338     FILE *fp;
339     int   i;
340
341     fp = gmx_ffopen(s, "w");
342     for (i = 0; (i < n); i++)
343     {
344         fprintf(fp, "%10d  %10g\n", i, c[i]);
345     }
346     gmx_ffclose(fp);
347 }
348
349 real print_and_integrate(FILE *fp, int n, real dt, real c[], real *fit, int nskip)
350 {
351     real c0, sum;
352     int  j;
353
354     /* Use trapezoidal rule for calculating integral */
355     sum = 0.0;
356     for (j = 0; (j < n); j++)
357     {
358         c0 = c[j];
359         if (fp && (nskip == 0 || j % nskip == 0))
360         {
361             fprintf(fp, "%10.3f  %10.5f\n", j*dt, c0);
362         }
363         if (j > 0)
364         {
365             sum += dt*(c0+c[j-1]);
366         }
367     }
368     if (fp)
369     {
370         fprintf(fp, "&\n");
371         if (fit)
372         {
373             for (j = 0; (j < n); j++)
374             {
375                 if (nskip == 0 || j % nskip == 0)
376                 {
377                     fprintf(fp, "%10.3f  %10.5f\n", j*dt, fit[j]);
378                 }
379             }
380             fprintf(fp, "&\n");
381         }
382     }
383     return sum*0.5;
384 }
385
386 real evaluate_integral(int n, real x[], real y[], real dy[], real aver_start,
387                        real *stddev)
388 {
389     double sum, sum_var, w;
390     double sum_tail = 0, sum2_tail = 0;
391     int    j, nsum_tail = 0;
392
393     /* Use trapezoidal rule for calculating integral */
394     if (n <= 0)
395     {
396         gmx_fatal(FARGS, "Evaluating integral: n = %d (file %s, line %d)",
397                   n, __FILE__, __LINE__);
398     }
399
400     sum     = 0;
401     sum_var = 0;
402     for (j = 0; (j < n); j++)
403     {
404         w = 0;
405         if (j > 0)
406         {
407             w += 0.5*(x[j] - x[j-1]);
408         }
409         if (j < n-1)
410         {
411             w += 0.5*(x[j+1] - x[j]);
412         }
413         sum += w*y[j];
414         if (dy)
415         {
416             /* Assume all errors are uncorrelated */
417             sum_var += sqr(w*dy[j]);
418         }
419
420         if ((aver_start > 0) && (x[j] >= aver_start))
421         {
422             sum_tail  += sum;
423             sum2_tail += sqrt(sum_var);
424             nsum_tail += 1;
425         }
426     }
427
428     if (nsum_tail > 0)
429     {
430         sum = sum_tail/nsum_tail;
431         /* This is a worst case estimate, assuming all stddev's are correlated. */
432         *stddev = sum2_tail/nsum_tail;
433     }
434     else
435     {
436         *stddev = sqrt(sum_var);
437     }
438
439     return sum;
440 }
441
442 void do_four_core(unsigned long mode, int nfour, int nf2, int nframes,
443                   real c1[], real csum[], real ctmp[])
444 {
445     real   *cfour;
446     char    buf[32];
447     real    fac;
448     int     j, m, m1;
449
450     snew(cfour, nfour);
451
452     if (MODE(eacNormal))
453     {
454         /********************************************
455          *  N O R M A L
456          ********************************************/
457         low_do_four_core(nfour, nf2, c1, csum, enNorm, FALSE);
458     }
459     else if (MODE(eacCos))
460     {
461         /***************************************************
462          * C O S I N E
463          ***************************************************/
464         /* Copy the data to temp array. Since we need it twice
465          * we can't overwrite original.
466          */
467         for (j = 0; (j < nf2); j++)
468         {
469             ctmp[j] = c1[j];
470         }
471
472         /* Cosine term of AC function */
473         low_do_four_core(nfour, nf2, ctmp, cfour, enCos, FALSE);
474         for (j = 0; (j < nf2); j++)
475         {
476             c1[j]  = cfour[j];
477         }
478
479         /* Sine term of AC function */
480         low_do_four_core(nfour, nf2, ctmp, cfour, enSin, FALSE);
481         for (j = 0; (j < nf2); j++)
482         {
483             c1[j]  += cfour[j];
484             csum[j] = c1[j];
485         }
486     }
487     else if (MODE(eacP2))
488     {
489         /***************************************************
490          * Legendre polynomials
491          ***************************************************/
492         /* First normalize the vectors */
493         norm_and_scale_vectors(nframes, c1, 1.0);
494
495         /* For P2 thingies we have to do six FFT based correls
496          * First for XX^2, then for YY^2, then for ZZ^2
497          * Then we have to do XY, YZ and XZ (counting these twice)
498          * After that we sum them and normalise
499          * P2(x) = (3 * cos^2 (x) - 1)/2
500          * for unit vectors u and v we compute the cosine as the inner product
501          * cos(u,v) = uX vX + uY vY + uZ vZ
502          *
503          *        oo
504          *        /
505          * C(t) = |  (3 cos^2(u(t'),u(t'+t)) - 1)/2 dt'
506          *        /
507          *        0
508          *
509          * For ACF we need:
510          * P2(u(0),u(t)) = [3 * (uX(0) uX(t) +
511          *                       uY(0) uY(t) +
512          *                       uZ(0) uZ(t))^2 - 1]/2
513          *               = [3 * ((uX(0) uX(t))^2 +
514          *                       (uY(0) uY(t))^2 +
515          *                       (uZ(0) uZ(t))^2 +
516          *                 2(uX(0) uY(0) uX(t) uY(t)) +
517          *                 2(uX(0) uZ(0) uX(t) uZ(t)) +
518          *                 2(uY(0) uZ(0) uY(t) uZ(t))) - 1]/2
519          *
520          *               = [(3/2) * (<uX^2> + <uY^2> + <uZ^2> +
521          *                         2<uXuY> + 2<uXuZ> + 2<uYuZ>) - 0.5]
522          *
523          */
524
525         /* Because of normalization the number of -0.5 to subtract
526          * depends on the number of data points!
527          */
528         for (j = 0; (j < nf2); j++)
529         {
530             csum[j]  = -0.5*(nf2-j);
531         }
532
533         /***** DIAGONAL ELEMENTS ************/
534         for (m = 0; (m < DIM); m++)
535         {
536             /* Copy the vector data in a linear array */
537             for (j = 0; (j < nf2); j++)
538             {
539                 ctmp[j]  = sqr(c1[DIM*j+m]);
540             }
541             if (debug)
542             {
543                 sprintf(buf, "c1diag%d.xvg", m);
544                 dump_tmp(buf, nf2, ctmp);
545             }
546
547             low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
548
549             if (debug)
550             {
551                 sprintf(buf, "c1dfout%d.xvg", m);
552                 dump_tmp(buf, nf2, cfour);
553             }
554             fac = 1.5;
555             for (j = 0; (j < nf2); j++)
556             {
557                 csum[j] += fac*(cfour[j]);
558             }
559         }
560         /******* OFF-DIAGONAL ELEMENTS **********/
561         for (m = 0; (m < DIM); m++)
562         {
563             /* Copy the vector data in a linear array */
564             m1 = (m+1) % DIM;
565             for (j = 0; (j < nf2); j++)
566             {
567                 ctmp[j] = c1[DIM*j+m]*c1[DIM*j+m1];
568             }
569
570             if (debug)
571             {
572                 sprintf(buf, "c1off%d.xvg", m);
573                 dump_tmp(buf, nf2, ctmp);
574             }
575             low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
576             if (debug)
577             {
578                 sprintf(buf, "c1ofout%d.xvg", m);
579                 dump_tmp(buf, nf2, cfour);
580             }
581             fac = 3.0;
582             for (j = 0; (j < nf2); j++)
583             {
584                 csum[j] += fac*cfour[j];
585             }
586         }
587     }
588     else if (MODE(eacP1) || MODE(eacVector))
589     {
590         /***************************************************
591          * V E C T O R & P1
592          ***************************************************/
593         if (MODE(eacP1))
594         {
595             /* First normalize the vectors */
596             norm_and_scale_vectors(nframes, c1, 1.0);
597         }
598
599         /* For vector thingies we have to do three FFT based correls
600          * First for XX, then for YY, then for ZZ
601          * After that we sum them and normalise
602          */
603         for (j = 0; (j < nf2); j++)
604         {
605             csum[j] = 0.0;
606         }
607         for (m = 0; (m < DIM); m++)
608         {
609             /* Copy the vector data in a linear array */
610             for (j = 0; (j < nf2); j++)
611             {
612                 ctmp[j] = c1[DIM*j+m];
613             }
614             low_do_four_core(nfour, nf2, ctmp, cfour, enNorm, FALSE);
615             for (j = 0; (j < nf2); j++)
616             {
617                 csum[j] += cfour[j];
618             }
619         }
620     }
621     else
622     {
623         gmx_fatal(FARGS, "\nUnknown mode in do_autocorr (%d)", mode);
624     }
625
626     sfree(cfour);
627     for (j = 0; (j < nf2); j++)
628     {
629         c1[j] = csum[j]/(real)(nframes-j);
630     }
631 }
632
633 real fit_acf(int ncorr, int fitfn, const output_env_t oenv, gmx_bool bVerbose,
634              real tbeginfit, real tendfit, real dt, real c1[], real *fit)
635 {
636     real        fitparm[3];
637     real        tStart, tail_corr, sum, sumtot = 0, c_start, ct_estimate, *sig;
638     int         i, j, jmax, nf_int;
639     gmx_bool    bPrint;
640
641     bPrint = bVerbose || bDebugMode();
642
643     if (bPrint)
644     {
645         printf("COR:\n");
646     }
647
648     if (tendfit <= 0)
649     {
650         tendfit = ncorr*dt;
651     }
652     nf_int = min(ncorr, (int)(tendfit/dt));
653     sum    = print_and_integrate(debug, nf_int, dt, c1, NULL, 1);
654
655     /* Estimate the correlation time for better fitting */
656     ct_estimate = 0.5*c1[0];
657     for (i = 1; (i < ncorr) && (c1[i] > 0); i++)
658     {
659         ct_estimate += c1[i];
660     }
661     ct_estimate *= dt/c1[0];
662
663     if (bPrint)
664     {
665         printf("COR: Correlation time (plain integral from %6.3f to %6.3f ps) = %8.5f ps\n",
666                0.0, dt*nf_int, sum);
667         printf("COR: Relaxation times are computed as fit to an exponential:\n");
668         printf("COR:   %s\n", longs_ffn[fitfn]);
669         printf("COR: Fit to correlation function from %6.3f ps to %6.3f ps, results in a\n", tbeginfit, min(ncorr*dt, tendfit));
670     }
671
672     tStart = 0;
673     if (bPrint)
674     {
675         printf("COR:%11s%11s%11s%11s%11s%11s%11s\n",
676                "Fit from", "Integral", "Tail Value", "Sum (ps)", " a1 (ps)",
677                (nfp_ffn[fitfn] >= 2) ? " a2 ()" : "",
678                (nfp_ffn[fitfn] >= 3) ? " a3 (ps)" : "");
679     }
680
681     snew(sig, ncorr);
682
683     if (tbeginfit > 0)
684     {
685         jmax = 3;
686     }
687     else
688     {
689         jmax = 1;
690     }
691     for (j = 0; ((j < jmax) && (tStart < tendfit) && (tStart < ncorr*dt)); j++)
692     {
693         /* Estimate the correlation time for better fitting */
694         c_start     = -1;
695         ct_estimate = 0;
696         for (i = 0; (i < ncorr) && (dt*i < tStart || c1[i] > 0); i++)
697         {
698             if (c_start < 0)
699             {
700                 if (dt*i >= tStart)
701                 {
702                     c_start     = c1[i];
703                     ct_estimate = 0.5*c1[i];
704                 }
705             }
706             else
707             {
708                 ct_estimate += c1[i];
709             }
710         }
711         if (c_start > 0)
712         {
713             ct_estimate *= dt/c_start;
714         }
715         else
716         {
717             /* The data is strange, so we need to choose somehting */
718             ct_estimate = tendfit;
719         }
720         if (debug)
721         {
722             fprintf(debug, "tStart %g ct_estimate: %g\n", tStart, ct_estimate);
723         }
724
725         if (fitfn == effnEXP3)
726         {
727             fitparm[0] = 0.002*ncorr*dt;
728             fitparm[1] = 0.95;
729             fitparm[2] = 0.2*ncorr*dt;
730         }
731         else
732         {
733             /* Good initial guess, this increases the probability of convergence */
734             fitparm[0] = ct_estimate;
735             fitparm[1] = 1.0;
736             fitparm[2] = 1.0;
737         }
738
739         /* Generate more or less appropriate sigma's */
740         for (i = 0; i < ncorr; i++)
741         {
742             sig[i] = sqrt(ct_estimate+dt*i);
743         }
744
745         nf_int    = min(ncorr, (int)((tStart+1e-4)/dt));
746         sum       = print_and_integrate(debug, nf_int, dt, c1, NULL, 1);
747         tail_corr = do_lmfit(ncorr, c1, sig, dt, NULL, tStart, tendfit, oenv,
748                              bDebugMode(), fitfn, fitparm, 0);
749         sumtot = sum+tail_corr;
750         if (fit && ((jmax == 1) || (j == 1)))
751         {
752             for (i = 0; (i < ncorr); i++)
753             {
754                 fit[i] = fit_function(fitfn, fitparm, i*dt);
755             }
756         }
757         if (bPrint)
758         {
759             printf("COR:%11.4e%11.4e%11.4e%11.4e", tStart, sum, tail_corr, sumtot);
760             for (i = 0; (i < nfp_ffn[fitfn]); i++)
761             {
762                 printf(" %11.4e", fitparm[i]);
763             }
764             printf("\n");
765         }
766         tStart += tbeginfit;
767     }
768     sfree(sig);
769
770     return sumtot;
771 }
772
773 void low_do_autocorr(const char *fn, const output_env_t oenv, const char *title,
774                      int nframes, int nitem, int nout, real **c1,
775                      real dt, unsigned long mode, int nrestart,
776                      gmx_bool bAver, gmx_bool bNormalize,
777                      gmx_bool bVerbose, real tbeginfit, real tendfit,
778                      int eFitFn)
779 {
780     FILE       *fp, *gp = NULL;
781     int         i, k, nfour;
782     real       *csum;
783     real       *ctmp, *fit;
784     real        c0, sum, Ct2av, Ctav;
785     gmx_bool    bFour = acf.bFour;
786
787     /* Check flags and parameters */
788     nout = get_acfnout();
789     if (nout == -1)
790     {
791         nout = acf.nout = (nframes+1)/2;
792     }
793     else if (nout > nframes)
794     {
795         nout = nframes;
796     }
797
798     if (MODE(eacCos) && MODE(eacVector))
799     {
800         gmx_fatal(FARGS, "Incompatible options bCos && bVector (%s, %d)",
801                   __FILE__, __LINE__);
802     }
803     if ((MODE(eacP3) || MODE(eacRcross)) && bFour)
804     {
805         fprintf(stderr, "Can't combine mode %lu with FFT, turning off FFT\n", mode);
806         bFour = FALSE;
807     }
808     if (MODE(eacNormal) && MODE(eacVector))
809     {
810         gmx_fatal(FARGS, "Incompatible mode bits: normal and vector (or Legendre)");
811     }
812
813     /* Print flags and parameters */
814     if (bVerbose)
815     {
816         printf("Will calculate %s of %d thingies for %d frames\n",
817                title ? title : "autocorrelation", nitem, nframes);
818         printf("bAver = %s, bFour = %s bNormalize= %s\n",
819                bool_names[bAver], bool_names[bFour], bool_names[bNormalize]);
820         printf("mode = %lu, dt = %g, nrestart = %d\n", mode, dt, nrestart);
821     }
822     if (bFour)
823     {
824         c0 = log((double)nframes)/log(2.0);
825         k  = c0;
826         if (k < c0)
827         {
828             k++;
829         }
830         k++;
831         nfour = 1<<k;
832         if (debug)
833         {
834             fprintf(debug, "Using FFT to calculate %s, #points for FFT = %d\n",
835                     title, nfour);
836         }
837
838         /* Allocate temp arrays */
839         snew(csum, nfour);
840         snew(ctmp, nfour);
841     }
842     else
843     {
844         nfour = 0; /* To keep the compiler happy */
845         snew(csum, nframes);
846         snew(ctmp, nframes);
847     }
848
849     /* Loop over items (e.g. molecules or dihedrals)
850      * In this loop the actual correlation functions are computed, but without
851      * normalizing them.
852      */
853     k = max(1, pow(10, (int)(log(nitem)/log(100))));
854     for (i = 0; i < nitem; i++)
855     {
856         if (bVerbose && ((i%k == 0 || i == nitem-1)))
857         {
858             fprintf(stderr, "\rThingie %d", i+1);
859         }
860
861         if (bFour)
862         {
863             do_four_core(mode, nfour, nframes, nframes, c1[i], csum, ctmp);
864         }
865         else
866         {
867             do_ac_core(nframes, nout, ctmp, c1[i], nrestart, mode);
868         }
869     }
870     if (bVerbose)
871     {
872         fprintf(stderr, "\n");
873     }
874     sfree(ctmp);
875     sfree(csum);
876
877     if (fn)
878     {
879         snew(fit, nout);
880         fp = xvgropen(fn, title, "Time (ps)", "C(t)", oenv);
881     }
882     else
883     {
884         fit = NULL;
885         fp  = NULL;
886     }
887     if (bAver)
888     {
889         if (nitem > 1)
890         {
891             average_acf(bVerbose, nframes, nitem, c1);
892         }
893
894         if (bNormalize)
895         {
896             normalize_acf(nout, c1[0]);
897         }
898
899         if (eFitFn != effnNONE)
900         {
901             fit_acf(nout, eFitFn, oenv, fn != NULL, tbeginfit, tendfit, dt, c1[0], fit);
902             sum = print_and_integrate(fp, nout, dt, c1[0], fit, 1);
903         }
904         else
905         {
906             sum = print_and_integrate(fp, nout, dt, c1[0], NULL, 1);
907             if (bVerbose)
908             {
909                 printf("Correlation time (integral over corrfn): %g (ps)\n", sum);
910             }
911         }
912     }
913     else
914     {
915         /* Not averaging. Normalize individual ACFs */
916         Ctav = Ct2av = 0;
917         if (debug)
918         {
919             gp = xvgropen("ct-distr.xvg", "Correlation times", "item", "time (ps)", oenv);
920         }
921         for (i = 0; i < nitem; i++)
922         {
923             if (bNormalize)
924             {
925                 normalize_acf(nout, c1[i]);
926             }
927             if (eFitFn != effnNONE)
928             {
929                 fit_acf(nout, eFitFn, oenv, fn != NULL, tbeginfit, tendfit, dt, c1[i], fit);
930                 sum = print_and_integrate(fp, nout, dt, c1[i], fit, 1);
931             }
932             else
933             {
934                 sum = print_and_integrate(fp, nout, dt, c1[i], NULL, 1);
935                 if (debug)
936                 {
937                     fprintf(debug,
938                             "CORRelation time (integral over corrfn %d): %g (ps)\n",
939                             i, sum);
940                 }
941             }
942             Ctav  += sum;
943             Ct2av += sum*sum;
944             if (debug)
945             {
946                 fprintf(gp, "%5d  %.3f\n", i, sum);
947             }
948         }
949         if (debug)
950         {
951             gmx_ffclose(gp);
952         }
953         if (nitem > 1)
954         {
955             Ctav  /= nitem;
956             Ct2av /= nitem;
957             printf("Average correlation time %.3f Std. Dev. %.3f Error %.3f (ps)\n",
958                    Ctav, sqrt((Ct2av - sqr(Ctav))),
959                    sqrt((Ct2av - sqr(Ctav))/(nitem-1)));
960         }
961     }
962     if (fp)
963     {
964         gmx_ffclose(fp);
965     }
966     sfree(fit);
967 }
968
969 static const char *Leg[]   = { NULL, "0", "1", "2", "3", NULL };
970
971 t_pargs *add_acf_pargs(int *npargs, t_pargs *pa)
972 {
973     t_pargs  acfpa[] = {
974         { "-acflen",     FALSE, etINT,  {&acf.nout},
975           "Length of the ACF, default is half the number of frames" },
976         { "-normalize", FALSE, etBOOL, {&acf.bNormalize},
977           "Normalize ACF" },
978         { "-fftcorr",  FALSE, etBOOL, {&acf.bFour},
979           "HIDDENUse fast fourier transform for correlation function" },
980         { "-nrestart", FALSE, etINT,  {&acf.nrestart},
981           "HIDDENNumber of frames between time origins for ACF when no FFT is used" },
982         { "-P",        FALSE, etENUM, {Leg},
983           "Order of Legendre polynomial for ACF (0 indicates none)" },
984         { "-fitfn",    FALSE, etENUM, {s_ffn},
985           "Fit function" },
986         { "-beginfit", FALSE, etREAL, {&acf.tbeginfit},
987           "Time where to begin the exponential fit of the correlation function" },
988         { "-endfit",   FALSE, etREAL, {&acf.tendfit},
989           "Time where to end the exponential fit of the correlation function, -1 is until the end" },
990     };
991 #define NPA asize(acfpa)
992     t_pargs *ppa;
993     int      i;
994
995     snew(ppa, *npargs+NPA);
996     for (i = 0; (i < *npargs); i++)
997     {
998         ppa[i] = pa[i];
999     }
1000     for (i = 0; (i < NPA); i++)
1001     {
1002         ppa[*npargs+i] = acfpa[i];
1003     }
1004     (*npargs) += NPA;
1005
1006     acf.mode       = 0;
1007     acf.nrestart   = 1;
1008     acf.nout       = -1;
1009     acf.P          = 0;
1010     acf.fitfn      = effnEXP1;
1011     acf.bFour      = TRUE;
1012     acf.bNormalize = TRUE;
1013     acf.tbeginfit  = 0.0;
1014     acf.tendfit    = -1;
1015
1016     bACFinit = TRUE;
1017
1018     return ppa;
1019 }
1020
1021 void do_autocorr(const char *fn, const output_env_t oenv, const char *title,
1022                  int nframes, int nitem, real **c1,
1023                  real dt, unsigned long mode, gmx_bool bAver)
1024 {
1025     if (!bACFinit)
1026     {
1027         printf("ACF data structures have not been initialised. Call add_acf_pargs\n");
1028     }
1029
1030     /* Handle enumerated types */
1031     sscanf(Leg[0], "%d", &acf.P);
1032     acf.fitfn = sffn2effn(s_ffn);
1033
1034     switch (acf.P)
1035     {
1036         case 1:
1037             mode = mode | eacP1;
1038             break;
1039         case 2:
1040             mode = mode | eacP2;
1041             break;
1042         case 3:
1043             mode = mode | eacP3;
1044             break;
1045         default:
1046             break;
1047     }
1048
1049     low_do_autocorr(fn, oenv, title, nframes, nitem, acf.nout, c1, dt, mode,
1050                     acf.nrestart, bAver, acf.bNormalize,
1051                     bDebugMode(), acf.tbeginfit, acf.tendfit,
1052                     acf.fitfn);
1053 }
1054
1055 int get_acfnout(void)
1056 {
1057     if (!bACFinit)
1058     {
1059         gmx_fatal(FARGS, "ACF data not initialized yet");
1060     }
1061
1062     return acf.nout;
1063 }
1064
1065 int get_acffitfn(void)
1066 {
1067     if (!bACFinit)
1068     {
1069         gmx_fatal(FARGS, "ACF data not initialized yet");
1070     }
1071
1072     return sffn2effn(s_ffn);
1073 }
1074
1075 void cross_corr(int n, real f[], real g[], real corr[])
1076 {
1077     int i, j;
1078
1079     if (acf.bFour)
1080     {
1081         correl(f-1, g-1, n, corr-1);
1082     }
1083     else
1084     {
1085         for (i = 0; (i < n); i++)
1086         {
1087             for (j = i; (j < n); j++)
1088             {
1089                 corr[j-i] += f[i]*g[j];
1090             }
1091             corr[i] /= (real)(n-i);
1092         }
1093     }
1094 }