Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxana / geminate.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2011,2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #include "gmxpre.h"
36
37 #include "geminate.h"
38
39 #include <math.h>
40 #include <stdlib.h>
41
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/gmxomp.h"
46 #include "gromacs/utility/smalloc.h"
47
48 static void missing_code_message()
49 {
50     fprintf(stderr, "You have requested code to run that is deprecated.\n");
51     fprintf(stderr, "Revert to an older GROMACS version or help in porting the code.\n");
52 }
53
54 /* The first few sections of this file contain functions that were adopted,
55  * and to some extent modified, by Erik Marklund (erikm[aT]xray.bmc.uu.se,
56  * http://folding.bmc.uu.se) from code written by Omer Markovitch (email, url).
57  * This is also the case with the function eq10v2().
58  *
59  * The parts menetioned in the previous paragraph were contributed under the BSD license.
60  */
61
62
63 /* This first part is from complex.c which I recieved from Omer Markowitch.
64  * - Erik Marklund
65  *
66  * ------------- from complex.c ------------- */
67
68 /* Complexation of a paired number (r,i)                                     */
69 static gem_complex gem_cmplx(double r, double i)
70 {
71     gem_complex value;
72     value.r = r;
73     value.i = i;
74     return value;
75 }
76
77 /* Complexation of a real number, x */
78 static gem_complex gem_c(double x)
79 {
80     gem_complex value;
81     value.r = x;
82     value.i = 0;
83     return value;
84 }
85
86 /* Magnitude of a complex number z                                           */
87 static double gem_cx_abs(gem_complex z)
88 {
89     return (sqrt(z.r*z.r+z.i*z.i));
90 }
91
92 /* Addition of two complex numbers z1 and z2                                 */
93 static gem_complex gem_cxadd(gem_complex z1, gem_complex z2)
94 {
95     gem_complex value;
96     value.r = z1.r+z2.r;
97     value.i = z1.i+z2.i;
98     return value;
99 }
100
101 /* Addition of a complex number z1 and a real number r */
102 static gem_complex gem_cxradd(gem_complex z, double r)
103 {
104     gem_complex value;
105     value.r = z.r + r;
106     value.i = z.i;
107     return value;
108 }
109
110 /* Subtraction of two complex numbers z1 and z2                              */
111 static gem_complex gem_cxsub(gem_complex z1, gem_complex z2)
112 {
113     gem_complex value;
114     value.r = z1.r-z2.r;
115     value.i = z1.i-z2.i;
116     return value;
117 }
118
119 /* Multiplication of two complex numbers z1 and z2                           */
120 static gem_complex gem_cxmul(gem_complex z1, gem_complex z2)
121 {
122     gem_complex value;
123     value.r = z1.r*z2.r-z1.i*z2.i;
124     value.i = z1.r*z2.i+z1.i*z2.r;
125     return value;
126 }
127
128 /* Square of a complex number z                                              */
129 static gem_complex gem_cxsq(gem_complex z)
130 {
131     gem_complex value;
132     value.r = z.r*z.r-z.i*z.i;
133     value.i = z.r*z.i*2.;
134     return value;
135 }
136
137 /* multiplication of a complex number z and a real number r */
138 static gem_complex gem_cxrmul(gem_complex z, double r)
139 {
140     gem_complex value;
141     value.r = z.r*r;
142     value.i = z.i*r;
143     return value;
144 }
145
146 /* Division of two complex numbers z1 and z2                                 */
147 static gem_complex gem_cxdiv(gem_complex z1, gem_complex z2)
148 {
149     gem_complex value;
150     double      num;
151     num = z2.r*z2.r+z2.i*z2.i;
152     if (num == 0.)
153     {
154         fprintf(stderr, "ERROR in gem_cxdiv function\n");
155     }
156     value.r = (z1.r*z2.r+z1.i*z2.i)/num; value.i = (z1.i*z2.r-z1.r*z2.i)/num;
157     return value;
158 }
159
160 /* division of a complex z number by a real number x */
161 static gem_complex gem_cxrdiv(gem_complex z, double r)
162 {
163     gem_complex value;
164     value.r = z.r/r;
165     value.i = z.i/r;
166     return value;
167 }
168
169 /* division of a real number r by a complex number x */
170 static gem_complex gem_rxcdiv(double r, gem_complex z)
171 {
172     gem_complex value;
173     double      f;
174     f       = r/(z.r*z.r+z.i*z.i);
175     value.r = f*z.r;
176     value.i = -f*z.i;
177     return value;
178 }
179
180 /* Exponential of a complex number-- exp (z)=|exp(z.r)|*{cos(z.i)+I*sin(z.i)}*/
181 static gem_complex gem_cxdexp(gem_complex z)
182 {
183     gem_complex value;
184     double      exp_z_r;
185     exp_z_r = exp(z.r);
186     value.r = exp_z_r*cos(z.i);
187     value.i = exp_z_r*sin(z.i);
188     return value;
189 }
190
191 /* Logarithm of a complex number -- log(z)=log|z|+I*Arg(z),                  */
192 /*                                  where -PI < Arg(z) < PI                  */
193 static gem_complex gem_cxlog(gem_complex z)
194 {
195     gem_complex value;
196     double      mag2;
197     mag2 = z.r*z.r+z.i*z.i;
198     if (mag2 < 0.)
199     {
200         fprintf(stderr, "ERROR in gem_cxlog func\n");
201     }
202     value.r = log(sqrt(mag2));
203     if (z.r == 0.)
204     {
205         value.i = PI/2.;
206         if (z.i < 0.)
207         {
208             value.i = -value.i;
209         }
210     }
211     else
212     {
213         value.i = atan2(z.i, z.r);
214     }
215     return value;
216 }
217
218 /* Square root of a complex number z = |z| exp(I*the) -- z^(1/2)             */
219 /*                               z^(1/2)=|z|^(1/2)*[cos(the/2)+I*sin(the/2)] */
220 /*                               where 0 < the < 2*PI                        */
221 static gem_complex gem_cxdsqrt(gem_complex z)
222 {
223     gem_complex value;
224     double      sq;
225     sq      = gem_cx_abs(z);
226     value.r = sqrt(fabs((sq+z.r)*0.5)); /* z'.r={|z|*[1+cos(the)]/2}^(1/2) */
227     value.i = sqrt(fabs((sq-z.r)*0.5)); /* z'.i={|z|*[1-cos(the)]/2}^(1/2) */
228     if (z.i < 0.)
229     {
230         value.r = -value.r;
231     }
232     return value;
233 }
234
235 /* Complex power of a complex number z1^z2                                   */
236 static gem_complex gem_cxdpow(gem_complex z1, gem_complex z2)
237 {
238     gem_complex value;
239     value = gem_cxdexp(gem_cxmul(gem_cxlog(z1), z2));
240     return value;
241 }
242
243 /* ------------ end of complex.c ------------ */
244
245 /* This next part was derived from cubic.c, also received from Omer Markovitch.
246  * ------------- from cubic.c ------------- */
247
248 /* Solver for a cubic equation: x^3-a*x^2+b*x-c=0                            */
249 static void gem_solve(gem_complex *al, gem_complex *be, gem_complex *gam,
250                       double a, double b, double c)
251 {
252     double      t1, t2, two_3, temp;
253     gem_complex ctemp, ct3;
254
255     two_3 = pow(2., 1./3.); t1 = -a*a+3.*b; t2 = 2.*a*a*a-9.*a*b+27.*c;
256     temp  = 4.*t1*t1*t1+t2*t2;
257
258     ctemp = gem_cmplx(temp, 0.);   ctemp = gem_cxadd(gem_cmplx(t2, 0.), gem_cxdsqrt(ctemp));
259     ct3   = gem_cxdpow(ctemp, gem_cmplx(1./3., 0.));
260
261     ctemp = gem_rxcdiv(-two_3*t1/3., ct3);
262     ctemp = gem_cxadd(ctemp, gem_cxrdiv(ct3, 3.*two_3));
263
264     *gam = gem_cxadd(gem_cmplx(a/3., 0.), ctemp);
265
266     ctemp = gem_cxmul(gem_cxsq(*gam), gem_cxsq(gem_cxsub(*gam, gem_cmplx(a, 0.))));
267     ctemp = gem_cxadd(ctemp, gem_cxmul(gem_cmplx(-4.*c, 0.), *gam));
268     ctemp = gem_cxdiv(gem_cxdsqrt(ctemp), *gam);
269     *al   = gem_cxrmul(gem_cxsub(gem_cxsub(gem_cmplx(a, 0.), *gam), ctemp), 0.5);
270     *be   = gem_cxrmul(gem_cxadd(gem_cxsub(gem_cmplx(a, 0.), *gam), ctemp), 0.5);
271 }
272
273 /* ------------ end of cubic.c ------------ */
274
275 /* This next part was derived from cerror.c and rerror.c, also received from Omer Markovitch.
276  * ------------- from [cr]error.c ------------- */
277
278 /************************************************************/
279 /* Real valued error function and related functions         */
280 /* x, y     : real variables                                */
281 /* erf(x)   : error function                                */
282 /* erfc(x)  : complementary error function                  */
283 /* omega(x) : exp(x*x)*erfc(x)                              */
284 /* W(x,y)   : exp(-x*x)*omega(x+y)=exp(2*x*y+y^2)*erfc(x+y) */
285 /************************************************************/
286
287 /*---------------------------------------------------------------------------*/
288 /* Utilzed the series approximation of erf(x)                                */
289 /* Relative error=|err(x)|/erf(x)<EPS                                        */
290 /* Handbook of Mathematical functions, Abramowitz, p 297                     */
291 /* Note: When x>=6 series sum deteriorates -> Asymptotic series used instead */
292 /*---------------------------------------------------------------------------*/
293 /*  This gives erfc(z) correctly only upto >10-15 */
294
295 static double gem_erf(double x)
296 {
297     double n, sum, temp, exp2, xm, x2, x4, x6, x8, x10, x12;
298     temp = x;
299     sum  = temp;
300     xm   = 26.;
301     x2   = x*x;
302     x4   = x2*x2;
303     x6   = x4*x2;
304     x8   = x6*x2;
305     x10  = x8*x2;
306     x12  = x10*x2;
307     exp2 = exp(-x2);
308     if (x <= xm)
309     {
310         for (n = 1.; n <= 2000.; n += 1.)
311         {
312             temp *= 2.*x2/(2.*n+1.);
313             sum  += temp;
314             if (fabs(temp/sum) < 1.E-16)
315             {
316                 break;
317             }
318         }
319
320         if (n >= 2000.)
321         {
322             fprintf(stderr, "In Erf calc - iteration exceeds %lg\n", n);
323         }
324         sum *= 2./sPI*exp2;
325     }
326     else
327     {
328         /* from the asymptotic expansion of experfc(x) */
329         sum = (1. - 0.5/x2 + 0.75/x4
330                - 1.875/x6 + 6.5625/x8
331                - 29.53125/x10 + 162.421875/x12)
332             / sPI/x;
333         sum *= exp2; /* now sum is erfc(x) */
334         sum  = -sum+1.;
335     }
336     return x >= 0.0 ? sum : -sum;
337 }
338
339 /* Result --> Alex's code for erfc and experfc behaves better than this      */
340 /* complementray error function.                Returns 1.-erf(x)            */
341 static double gem_erfc(double x)
342 {
343     double t, z, ans;
344     z = fabs(x);
345     t = 1.0/(1.0+0.5*z);
346
347     ans = t * exp(-z*z - 1.26551223 +
348                   t*(1.00002368 +
349                      t*(0.37409196 +
350                         t*(0.09678418 +
351                            t*(-0.18628806 +
352                               t*(0.27886807 +
353                                  t*(-1.13520398 +
354                                     t*(1.48851587 +
355                                        t*(-0.82215223 +
356                                           t*0.17087277)))))))));
357
358     return x >= 0.0 ? ans : 2.0-ans;
359 }
360
361 /* omega(x)=exp(x^2)erfc(x)                                                  */
362 static double gem_omega(double x)
363 {
364     double xm, ans, xx, x4, x6, x8, x10, x12;
365     xm  = 26;
366     xx  = x*x;
367     x4  = xx*xx;
368     x6  = x4*xx;
369     x8  = x6*xx;
370     x10 = x8*xx;
371     x12 = x10*xx;
372
373     if (x <= xm)
374     {
375         ans = exp(xx)*gem_erfc(x);
376     }
377     else
378     {
379         /* Asymptotic expansion */
380         ans = (1. - 0.5/xx + 0.75/x4 - 1.875/x6 + 6.5625/x8 - 29.53125/x10 + 162.421875/x12) / sPI/x;
381     }
382     return ans;
383 }
384
385 /*---------------------------------------------------------------------------*/
386 /* Utilzed the series approximation of erf(z=x+iy)                           */
387 /* Relative error=|err(z)|/|erf(z)|<EPS                                      */
388 /* Handbook of Mathematical functions, Abramowitz, p 299                     */
389 /* comega(z=x+iy)=cexp(z^2)*cerfc(z)                                         */
390 /*---------------------------------------------------------------------------*/
391 static gem_complex gem_comega(gem_complex z)
392 {
393     gem_complex value;
394     double      x, y;
395     double      sumr, sumi, n, n2, f, temp, temp1;
396     double      x2, y2, cos_2xy, sin_2xy, cosh_2xy, sinh_2xy, cosh_ny, sinh_ny, exp_y2;
397
398     x        = z.r;
399     y        = z.i;
400     x2       = x*x;
401     y2       = y*y;
402     sumr     = 0.;
403     sumi     = 0.;
404     cos_2xy  = cos(2.*x*y);
405     sin_2xy  = sin(2.*x*y);
406     cosh_2xy = cosh(2.*x*y);
407     sinh_2xy = sinh(2.*x*y);
408     exp_y2   = exp(-y2);
409
410     for (n = 1.0, temp = 0.; n <= 2000.; n += 1.0)
411     {
412         n2      = n*n;
413         cosh_ny = cosh(n*y);
414         sinh_ny = sinh(n*y);
415         f       = exp(-n2/4.)/(n2+4.*x2);
416         /* if(f<1.E-200) break; */
417         sumr += (2.*x - 2.*x*cosh_ny*cos_2xy + n*sinh_ny*sin_2xy)*f;
418         sumi += (2.*x*cosh_ny*sin_2xy + n*sinh_ny*cos_2xy)*f;
419         temp1 = sqrt(sumr*sumr+sumi*sumi);
420         if (fabs((temp1-temp)/temp1) < 1.E-16)
421         {
422             break;
423         }
424         temp = temp1;
425     }
426     if (n == 2000.)
427     {
428         fprintf(stderr, "iteration exceeds %lg\n", n);
429     }
430     sumr *= 2./PI;
431     sumi *= 2./PI;
432
433     if (x != 0.)
434     {
435         f = 1./2./PI/x;
436     }
437     else
438     {
439         f = 0.;
440     }
441     value.r = gem_omega(x)-(f*(1.-cos_2xy)+sumr);
442     value.i = -(f*sin_2xy+sumi);
443     value   = gem_cxmul(value, gem_cmplx(exp_y2*cos_2xy, exp_y2*sin_2xy));
444     return (value);
445 }
446
447 /* ------------ end of [cr]error.c ------------ */
448
449 /*_ REVERSIBLE GEMINATE RECOMBINATION
450  *
451  * Here are the functions for reversible geminate recombination. */
452
453 /* Changes the unit from square cm per s to square Ã…ngström per ps,
454  * since Omers code uses the latter units while g_mds outputs the former.
455  * g_hbond expects a diffusion coefficent given in square cm per s. */
456 static double sqcm_per_s_to_sqA_per_ps (real D)
457 {
458     fprintf(stdout, "Diffusion coefficient is %f A^2/ps\n", D*1e4);
459     return (double)(D*1e4);
460 }
461
462
463 static double eq10v2(double theoryCt[], double time[], int manytimes,
464                      double ka, double kd, t_gemParams *params)
465 {
466     /* Finding the 3 roots */
467     int    i;
468     double kD, D, r, a, b, c, tsqrt, sumimaginary;
469     gem_complex
470            alpha, beta, gamma,
471         c1, c2, c3, c4,
472         oma, omb, omc,
473         part1, part2, part3, part4;
474
475     kD = params->kD;
476     D  = params->D;
477     r  = params->sigma;
478     a  = (1.0 + ka/kD) * sqrt(D)/r;
479     b  = kd;
480     c  = kd * sqrt(D)/r;
481
482     gem_solve(&alpha, &beta, &gamma, a, b, c);
483     /* Finding the 3 roots */
484
485     sumimaginary = 0;
486     part1        = gem_cxmul(alpha, gem_cxmul(gem_cxadd(beta,  gamma), gem_cxsub(beta,  gamma)));                 /* 1(2+3)(2-3) */
487     part2        = gem_cxmul(beta,  gem_cxmul(gem_cxadd(gamma, alpha), gem_cxsub(gamma, alpha)));                 /* 2(3+1)(3-1) */
488     part3        = gem_cxmul(gamma, gem_cxmul(gem_cxadd(alpha, beta), gem_cxsub(alpha, beta)));                   /* 3(1+2)(1-2) */
489     part4        = gem_cxmul(gem_cxsub(gamma, alpha), gem_cxmul(gem_cxsub(alpha, beta), gem_cxsub(beta, gamma))); /* (3-1)(1-2)(2-3) */
490
491 #pragma omp parallel for \
492     private(i, tsqrt, oma, omb, omc, c1, c2, c3, c4) \
493     reduction(+:sumimaginary) \
494     default(shared) \
495     schedule(guided)
496     for (i = 0; i < manytimes; i++)
497     {
498         tsqrt         = sqrt(time[i]);
499         oma           = gem_comega(gem_cxrmul(alpha, tsqrt));
500         c1            = gem_cxmul(oma, gem_cxdiv(part1, part4));
501         omb           = gem_comega(gem_cxrmul(beta, tsqrt));
502         c2            = gem_cxmul(omb, gem_cxdiv(part2, part4));
503         omc           = gem_comega(gem_cxrmul(gamma, tsqrt));
504         c3            = gem_cxmul(omc, gem_cxdiv(part3, part4));
505         c4.r          = c1.r+c2.r+c3.r;
506         c4.i          = c1.i+c2.i+c3.i;
507         theoryCt[i]   = c4.r;
508         sumimaginary += c4.i * c4.i;
509     }
510
511     return sumimaginary;
512
513 } /* eq10v2 */
514
515 /* This returns the real-valued index(!) to an ACF, equidistant on a log scale. */
516 static double getLogIndex(const int i, const t_gemParams *params)
517 {
518     return (exp(((double)(i)) * params->logQuota) -1);
519 }
520
521 extern t_gemParams *init_gemParams(const double sigma, const double D,
522                                    const real *t, const int len, const int nFitPoints,
523                                    const real begFit, const real endFit,
524                                    const real ballistic, const int nBalExp)
525 {
526     double       tDelta;
527     t_gemParams *p;
528
529     snew(p, 1);
530
531     /* A few hardcoded things here. For now anyway. */
532 /*   p->ka_min   = 0; */
533 /*   p->ka_max   = 100; */
534 /*   p->dka      = 10; */
535 /*   p->kd_min   = 0; */
536 /*   p->kd_max   = 2; */
537 /*   p->dkd      = 0.2; */
538     p->ka       = 0;
539     p->kd       = 0;
540 /*   p->lsq      = -1; */
541 /*   p->lifetime = 0; */
542     p->sigma    = sigma*10; /* Omer uses Ã…, not nm */
543 /*   p->lsq_old  = 99999; */
544     p->D        = sqcm_per_s_to_sqA_per_ps(D);
545     p->kD       = 4*acos(-1.0)*sigma*p->D;
546
547
548     /* Parameters used by calcsquare(). Better to calculate them
549      * here than in calcsquare every time it's called. */
550     p->len = len;
551 /*   p->logAfterTime = logAfterTime; */
552     tDelta       = (t[len-1]-t[0]) / len;
553     if (tDelta <= 0)
554     {
555         gmx_fatal(FARGS, "Time between frames is non-positive!");
556     }
557     else
558     {
559         p->tDelta = tDelta;
560     }
561
562     p->nExpFit      = nBalExp;
563 /*   p->nLin         = logAfterTime / tDelta; */
564     p->nFitPoints   = nFitPoints;
565     p->begFit       = begFit;
566     p->endFit       = endFit;
567     p->logQuota     = (double)(log(p->len))/(p->nFitPoints-1);
568 /*   if (p->nLin <= 0) { */
569 /*     fprintf(stderr, "Number of data points in the linear regime is non-positive!\n"); */
570 /*     sfree(p); */
571 /*     return NULL; */
572 /*   } */
573 /* We want the same number of data points in the log regime. Not crucial, but seems like a good idea. */
574 /* p->logDelta = log(((float)len)/p->nFitPoints) / p->nFitPoints;/\* log(((float)len)/p->nLin) / p->nLin; *\/ */
575 /*   p->logPF    = p->nFitPoints*p->nFitPoints/(float)len; /\* p->nLin*p->nLin/(float)len; *\/ */
576 /* logPF and logDelta are stitched together with the macro GETLOGINDEX defined in geminate.h */
577
578 /*   p->logMult      = pow((float)len, 1.0/nLin);/\* pow(t[len-1]-t[0], 1.0/p->nLin); *\/ */
579     p->ballistic    =  ballistic;
580     return p;
581 }
582
583 /* There was a misunderstanding regarding the fitting. From our
584  * recent correspondence it appears that Omer's code require
585  * the ACF data on a log-scale and does not operate on the raw data.
586  * This needs to be redone in gemFunc_residual() as well as in the
587  * t_gemParams structure. */
588
589 static real* d2r(const double *d, const int nn);
590
591 extern real fitGemRecomb(double gmx_unused      *ct,
592                          double gmx_unused      *time,
593                          double gmx_unused     **ctFit,
594                          const int gmx_unused    nData,
595                          t_gemParams gmx_unused *params)
596 {
597
598     int         nThreads, i, iter, status, maxiter;
599     real        size, d2, tol, *dumpdata;
600     size_t      p, n;
601     gemFitData *GD;
602     char       *dumpstr, dumpname[128];
603
604     missing_code_message();
605     return -1;
606
607 }
608
609
610 /* Removes the ballistic term from the beginning of the ACF,
611  * just like in Omer's paper.
612  */
613 extern void takeAwayBallistic(double gmx_unused *ct, double *t, int len, real tMax, int nexp, gmx_bool gmx_unused bDerivative)
614 {
615
616     /* Fit with 4 exponentials and one constant term,
617      * subtract the fatest exponential. */
618
619     int      nData, i, status, iter;
620     balData *BD;
621     double  *guess,           /* Initial guess. */
622     *A,                       /* The fitted parameters. (A1, B1, A2, B2,... C) */
623              a[2],
624              ddt[2];
625     gmx_bool sorted;
626     size_t   n;
627     size_t   p;
628
629     nData = 0;
630     do
631     {
632         nData++;
633     }
634     while (t[nData] < tMax+t[0] && nData < len);
635
636     p = nexp*2+1;            /* Number of parameters. */
637
638     missing_code_message();
639     return;
640 }
641
642 extern void dumpN(const real *e, const int nn, char *fn)
643 {
644     /* For debugging only */
645     int   i;
646     FILE *f;
647     char  standardName[] = "Nt.xvg";
648     if (fn == NULL)
649     {
650         fn = standardName;
651     }
652
653     f = fopen(fn, "w");
654     fprintf(f,
655             "@ type XY\n"
656             "@ xaxis label \"Frame\"\n"
657             "@ yaxis label \"N\"\n"
658             "@ s0 line type 3\n");
659
660     for (i = 0; i < nn; i++)
661     {
662         fprintf(f, "%-10i %-g\n", i, e[i]);
663     }
664
665     fclose(f);
666 }
667
668 static real* d2r(const double *d, const int nn)
669 {
670     real *r;
671     int   i;
672
673     snew(r, nn);
674     for (i = 0; i < nn; i++)
675     {
676         r[i] = (real)d[i];
677     }
678
679     return r;
680 }
681
682 static void _patchBad(double *ct, int n, double dy)
683 {
684     /* Just do lin. interpolation for now. */
685     int i;
686
687     for (i = 1; i < n; i++)
688     {
689         ct[i] = ct[0]+i*dy;
690     }
691 }
692
693 static void patchBadPart(double *ct, int n)
694 {
695     _patchBad(ct, n, (ct[n] - ct[0])/n);
696 }
697
698 static void patchBadTail(double *ct, int n)
699 {
700     _patchBad(ct+1, n-1, ct[1]-ct[0]);
701
702 }
703
704 extern void fixGemACF(double *ct, int len)
705 {
706     int      i, j, b, e;
707     gmx_bool bBad;
708
709     /* Let's separate two things:
710      * - identification of bad parts
711      * - patching of bad parts.
712      */
713
714     b    = 0; /* Start of a bad stretch */
715     e    = 0; /* End of a bad stretch */
716     bBad = FALSE;
717
718     /* An acf of binary data must be one at t=0. */
719     if (fabs(ct[0]-1.0) > 1e-6)
720     {
721         ct[0] = 1.0;
722         fprintf(stderr, "|ct[0]-1.0| = %1.6f. Setting ct[0] to 1.0.\n", fabs(ct[0]-1.0));
723     }
724
725     for (i = 0; i < len; i++)
726     {
727
728 #ifdef HAS_ISFINITE
729         if (isfinite(ct[i]))
730 #elif defined(HAS__ISFINITE)
731         if (_isfinite(ct[i]))
732 #else
733         if (1)
734 #endif
735         {
736             if (!bBad)
737             {
738                 /* Still on a good stretch. Proceed.*/
739                 continue;
740             }
741
742             /* Patch up preceding bad stretch. */
743             if (i == (len-1))
744             {
745                 /* It's the tail */
746                 if (b <= 1)
747                 {
748                     gmx_fatal(FARGS, "The ACF is mostly NaN or Inf. Aborting.");
749                 }
750                 patchBadTail(&(ct[b-2]), (len-b)+1);
751             }
752
753             e = i;
754             patchBadPart(&(ct[b-1]), (e-b)+1);
755             bBad = FALSE;
756         }
757         else
758         {
759             if (!bBad)
760             {
761                 b = i;
762
763                 bBad = TRUE;
764             }
765         }
766     }
767 }