3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
40 #include <gsl/gsl_rng.h>
41 #include <gsl/gsl_randist.h>
42 #include <gsl/gsl_vector.h>
43 #include <gsl/gsl_blas.h>
44 #include <gsl/gsl_multimin.h>
45 #include <gsl/gsl_multifit_nlin.h>
46 #include <gsl/gsl_sf.h>
47 #include <gsl/gsl_version.h>
56 #include "gromacs/utility/gmxomp.h"
58 /* The first few sections of this file contain functions that were adopted,
59 * and to some extent modified, by Erik Marklund (erikm[aT]xray.bmc.uu.se,
60 * http://folding.bmc.uu.se) from code written by Omer Markovitch (email, url).
61 * This is also the case with the function eq10v2().
63 * The parts menetioned in the previous paragraph were contributed under the BSD license.
67 /* This first part is from complex.c which I recieved from Omer Markowitch.
70 * ------------- from complex.c ------------- */
72 /* Complexation of a paired number (r,i) */
73 static gem_complex gem_cmplx(double r, double i)
81 /* Complexation of a real number, x */
82 static gem_complex gem_c(double x)
90 /* Magnitude of a complex number z */
91 static double gem_cx_abs(gem_complex z)
93 return (sqrt(z.r*z.r+z.i*z.i));
96 /* Addition of two complex numbers z1 and z2 */
97 static gem_complex gem_cxadd(gem_complex z1, gem_complex z2)
105 /* Addition of a complex number z1 and a real number r */
106 static gem_complex gem_cxradd(gem_complex z, double r)
114 /* Subtraction of two complex numbers z1 and z2 */
115 static gem_complex gem_cxsub(gem_complex z1, gem_complex z2)
123 /* Multiplication of two complex numbers z1 and z2 */
124 static gem_complex gem_cxmul(gem_complex z1, gem_complex z2)
127 value.r = z1.r*z2.r-z1.i*z2.i;
128 value.i = z1.r*z2.i+z1.i*z2.r;
132 /* Square of a complex number z */
133 static gem_complex gem_cxsq(gem_complex z)
136 value.r = z.r*z.r-z.i*z.i;
137 value.i = z.r*z.i*2.;
141 /* multiplication of a complex number z and a real number r */
142 static gem_complex gem_cxrmul(gem_complex z, double r)
150 /* Division of two complex numbers z1 and z2 */
151 static gem_complex gem_cxdiv(gem_complex z1, gem_complex z2)
155 num = z2.r*z2.r+z2.i*z2.i;
158 fprintf(stderr, "ERROR in gem_cxdiv function\n");
160 value.r = (z1.r*z2.r+z1.i*z2.i)/num; value.i = (z1.i*z2.r-z1.r*z2.i)/num;
164 /* division of a complex z number by a real number x */
165 static gem_complex gem_cxrdiv(gem_complex z, double r)
173 /* division of a real number r by a complex number x */
174 static gem_complex gem_rxcdiv(double r, gem_complex z)
178 f = r/(z.r*z.r+z.i*z.i);
184 /* Exponential of a complex number-- exp (z)=|exp(z.r)|*{cos(z.i)+I*sin(z.i)}*/
185 static gem_complex gem_cxdexp(gem_complex z)
190 value.r = exp_z_r*cos(z.i);
191 value.i = exp_z_r*sin(z.i);
195 /* Logarithm of a complex number -- log(z)=log|z|+I*Arg(z), */
196 /* where -PI < Arg(z) < PI */
197 static gem_complex gem_cxlog(gem_complex z)
201 mag2 = z.r*z.r+z.i*z.i;
204 fprintf(stderr, "ERROR in gem_cxlog func\n");
206 value.r = log(sqrt(mag2));
217 value.i = atan2(z.i, z.r);
222 /* Square root of a complex number z = |z| exp(I*the) -- z^(1/2) */
223 /* z^(1/2)=|z|^(1/2)*[cos(the/2)+I*sin(the/2)] */
224 /* where 0 < the < 2*PI */
225 static gem_complex gem_cxdsqrt(gem_complex z)
230 value.r = sqrt(fabs((sq+z.r)*0.5)); /* z'.r={|z|*[1+cos(the)]/2}^(1/2) */
231 value.i = sqrt(fabs((sq-z.r)*0.5)); /* z'.i={|z|*[1-cos(the)]/2}^(1/2) */
239 /* Complex power of a complex number z1^z2 */
240 static gem_complex gem_cxdpow(gem_complex z1, gem_complex z2)
243 value = gem_cxdexp(gem_cxmul(gem_cxlog(z1), z2));
247 /* ------------ end of complex.c ------------ */
249 /* This next part was derived from cubic.c, also received from Omer Markovitch.
250 * ------------- from cubic.c ------------- */
252 /* Solver for a cubic equation: x^3-a*x^2+b*x-c=0 */
253 static void gem_solve(gem_complex *al, gem_complex *be, gem_complex *gam,
254 double a, double b, double c)
256 double t1, t2, two_3, temp;
257 gem_complex ctemp, ct3;
259 two_3 = pow(2., 1./3.); t1 = -a*a+3.*b; t2 = 2.*a*a*a-9.*a*b+27.*c;
260 temp = 4.*t1*t1*t1+t2*t2;
262 ctemp = gem_cmplx(temp, 0.); ctemp = gem_cxadd(gem_cmplx(t2, 0.), gem_cxdsqrt(ctemp));
263 ct3 = gem_cxdpow(ctemp, gem_cmplx(1./3., 0.));
265 ctemp = gem_rxcdiv(-two_3*t1/3., ct3);
266 ctemp = gem_cxadd(ctemp, gem_cxrdiv(ct3, 3.*two_3));
268 *gam = gem_cxadd(gem_cmplx(a/3., 0.), ctemp);
270 ctemp = gem_cxmul(gem_cxsq(*gam), gem_cxsq(gem_cxsub(*gam, gem_cmplx(a, 0.))));
271 ctemp = gem_cxadd(ctemp, gem_cxmul(gem_cmplx(-4.*c, 0.), *gam));
272 ctemp = gem_cxdiv(gem_cxdsqrt(ctemp), *gam);
273 *al = gem_cxrmul(gem_cxsub(gem_cxsub(gem_cmplx(a, 0.), *gam), ctemp), 0.5);
274 *be = gem_cxrmul(gem_cxadd(gem_cxsub(gem_cmplx(a, 0.), *gam), ctemp), 0.5);
277 /* ------------ end of cubic.c ------------ */
279 /* This next part was derived from cerror.c and rerror.c, also received from Omer Markovitch.
280 * ------------- from [cr]error.c ------------- */
282 /************************************************************/
283 /* Real valued error function and related functions */
284 /* x, y : real variables */
285 /* erf(x) : error function */
286 /* erfc(x) : complementary error function */
287 /* omega(x) : exp(x*x)*erfc(x) */
288 /* W(x,y) : exp(-x*x)*omega(x+y)=exp(2*x*y+y^2)*erfc(x+y) */
289 /************************************************************/
291 /*---------------------------------------------------------------------------*/
292 /* Utilzed the series approximation of erf(x) */
293 /* Relative error=|err(x)|/erf(x)<EPS */
294 /* Handbook of Mathematical functions, Abramowitz, p 297 */
295 /* Note: When x>=6 series sum deteriorates -> Asymptotic series used instead */
296 /*---------------------------------------------------------------------------*/
297 /* This gives erfc(z) correctly only upto >10-15 */
299 static double gem_erf(double x)
301 double n, sum, temp, exp2, xm, x2, x4, x6, x8, x10, x12;
314 for (n = 1.; n <= 2000.; n += 1.)
316 temp *= 2.*x2/(2.*n+1.);
318 if (fabs(temp/sum) < 1.E-16)
326 fprintf(stderr, "In Erf calc - iteration exceeds %lg\n", n);
332 /* from the asymptotic expansion of experfc(x) */
333 sum = (1. - 0.5/x2 + 0.75/x4
334 - 1.875/x6 + 6.5625/x8
335 - 29.53125/x10 + 162.421875/x12)
337 sum *= exp2; /* now sum is erfc(x) */
340 return x >= 0.0 ? sum : -sum;
343 /* Result --> Alex's code for erfc and experfc behaves better than this */
344 /* complementray error function. Returns 1.-erf(x) */
345 static double gem_erfc(double x)
351 ans = t * exp(-z*z - 1.26551223 +
360 t*0.17087277)))))))));
362 return x >= 0.0 ? ans : 2.0-ans;
365 /* omega(x)=exp(x^2)erfc(x) */
366 static double gem_omega(double x)
368 double xm, ans, xx, x4, x6, x8, x10, x12;
379 ans = exp(xx)*gem_erfc(x);
383 /* Asymptotic expansion */
384 ans = (1. - 0.5/xx + 0.75/x4 - 1.875/x6 + 6.5625/x8 - 29.53125/x10 + 162.421875/x12) / sPI/x;
389 /*---------------------------------------------------------------------------*/
390 /* Utilzed the series approximation of erf(z=x+iy) */
391 /* Relative error=|err(z)|/|erf(z)|<EPS */
392 /* Handbook of Mathematical functions, Abramowitz, p 299 */
393 /* comega(z=x+iy)=cexp(z^2)*cerfc(z) */
394 /*---------------------------------------------------------------------------*/
395 static gem_complex gem_comega(gem_complex z)
399 double sumr, sumi, n, n2, f, temp, temp1;
400 double x2, y2, cos_2xy, sin_2xy, cosh_2xy, sinh_2xy, cosh_ny, sinh_ny, exp_y2;
408 cos_2xy = cos(2.*x*y);
409 sin_2xy = sin(2.*x*y);
410 cosh_2xy = cosh(2.*x*y);
411 sinh_2xy = sinh(2.*x*y);
414 for (n = 1.0, temp = 0.; n <= 2000.; n += 1.0)
419 f = exp(-n2/4.)/(n2+4.*x2);
420 /* if(f<1.E-200) break; */
421 sumr += (2.*x - 2.*x*cosh_ny*cos_2xy + n*sinh_ny*sin_2xy)*f;
422 sumi += (2.*x*cosh_ny*sin_2xy + n*sinh_ny*cos_2xy)*f;
423 temp1 = sqrt(sumr*sumr+sumi*sumi);
424 if (fabs((temp1-temp)/temp1) < 1.E-16)
432 fprintf(stderr, "iteration exceeds %lg\n", n);
445 value.r = gem_omega(x)-(f*(1.-cos_2xy)+sumr);
446 value.i = -(f*sin_2xy+sumi);
447 value = gem_cxmul(value, gem_cmplx(exp_y2*cos_2xy, exp_y2*sin_2xy));
451 /* ------------ end of [cr]error.c ------------ */
453 /*_ REVERSIBLE GEMINATE RECOMBINATION
455 * Here are the functions for reversible geminate recombination. */
457 /* Changes the unit from square cm per s to square Ångström per ps,
458 * since Omers code uses the latter units while g_mds outputs the former.
459 * g_hbond expects a diffusion coefficent given in square cm per s. */
460 static double sqcm_per_s_to_sqA_per_ps (real D)
462 fprintf(stdout, "Diffusion coefficient is %f A^2/ps\n", D*1e4);
463 return (double)(D*1e4);
467 static double eq10v2(double theoryCt[], double time[], int manytimes,
468 double ka, double kd, t_gemParams *params)
470 /* Finding the 3 roots */
472 double kD, D, r, a, b, c, tsqrt, sumimaginary;
477 part1, part2, part3, part4;
482 a = (1.0 + ka/kD) * sqrt(D)/r;
486 gem_solve(&alpha, &beta, &gamma, a, b, c);
487 /* Finding the 3 roots */
490 part1 = gem_cxmul(alpha, gem_cxmul(gem_cxadd(beta, gamma), gem_cxsub(beta, gamma))); /* 1(2+3)(2-3) */
491 part2 = gem_cxmul(beta, gem_cxmul(gem_cxadd(gamma, alpha), gem_cxsub(gamma, alpha))); /* 2(3+1)(3-1) */
492 part3 = gem_cxmul(gamma, gem_cxmul(gem_cxadd(alpha, beta), gem_cxsub(alpha, beta))); /* 3(1+2)(1-2) */
493 part4 = gem_cxmul(gem_cxsub(gamma, alpha), gem_cxmul(gem_cxsub(alpha, beta), gem_cxsub(beta, gamma))); /* (3-1)(1-2)(2-3) */
495 #pragma omp parallel for \
496 private(i, tsqrt, oma, omb, omc, c1, c2, c3, c4) \
497 reduction(+:sumimaginary) \
500 for (i = 0; i < manytimes; i++)
502 tsqrt = sqrt(time[i]);
503 oma = gem_comega(gem_cxrmul(alpha, tsqrt));
504 c1 = gem_cxmul(oma, gem_cxdiv(part1, part4));
505 omb = gem_comega(gem_cxrmul(beta, tsqrt));
506 c2 = gem_cxmul(omb, gem_cxdiv(part2, part4));
507 omc = gem_comega(gem_cxrmul(gamma, tsqrt));
508 c3 = gem_cxmul(omc, gem_cxdiv(part3, part4));
509 c4.r = c1.r+c2.r+c3.r;
510 c4.i = c1.i+c2.i+c3.i;
512 sumimaginary += c4.i * c4.i;
519 /* This returns the real-valued index(!) to an ACF, equidistant on a log scale. */
520 static double getLogIndex(const int i, const t_gemParams *params)
522 return (exp(((double)(i)) * params->logQuota) -1);
525 extern t_gemParams *init_gemParams(const double sigma, const double D,
526 const real *t, const int len, const int nFitPoints,
527 const real begFit, const real endFit,
528 const real ballistic, const int nBalExp)
535 /* A few hardcoded things here. For now anyway. */
537 /* p->ka_max = 100; */
545 /* p->lifetime = 0; */
546 p->sigma = sigma*10; /* Omer uses Ã…, not nm */
547 /* p->lsq_old = 99999; */
548 p->D = sqcm_per_s_to_sqA_per_ps(D);
549 p->kD = 4*acos(-1.0)*sigma*p->D;
552 /* Parameters used by calcsquare(). Better to calculate them
553 * here than in calcsquare every time it's called. */
555 /* p->logAfterTime = logAfterTime; */
556 tDelta = (t[len-1]-t[0]) / len;
559 gmx_fatal(FARGS, "Time between frames is non-positive!");
566 p->nExpFit = nBalExp;
567 /* p->nLin = logAfterTime / tDelta; */
568 p->nFitPoints = nFitPoints;
571 p->logQuota = (double)(log(p->len))/(p->nFitPoints-1);
572 /* if (p->nLin <= 0) { */
573 /* fprintf(stderr, "Number of data points in the linear regime is non-positive!\n"); */
577 /* We want the same number of data points in the log regime. Not crucial, but seems like a good idea. */
578 /* p->logDelta = log(((float)len)/p->nFitPoints) / p->nFitPoints;/\* log(((float)len)/p->nLin) / p->nLin; *\/ */
579 /* p->logPF = p->nFitPoints*p->nFitPoints/(float)len; /\* p->nLin*p->nLin/(float)len; *\/ */
580 /* logPF and logDelta are stitched together with the macro GETLOGINDEX defined in geminate.h */
582 /* p->logMult = pow((float)len, 1.0/nLin);/\* pow(t[len-1]-t[0], 1.0/p->nLin); *\/ */
583 p->ballistic = ballistic;
587 /* There was a misunderstanding regarding the fitting. From our
588 * recent correspondence it appears that Omer's code require
589 * the ACF data on a log-scale and does not operate on the raw data.
590 * This needs to be redone in gemFunc_residual() as well as in the
591 * t_gemParams structure. */
593 static double gemFunc_residual2(const gsl_vector *p, void *data)
596 int i, iLog, nLin, nFitPoints, nData;
597 double r, residual2, *ctTheory, *y;
599 GD = (gemFitData *)data;
600 nLin = GD->params->nLin;
601 nFitPoints = GD->params->nFitPoints;
604 ctTheory = GD->ctTheory;
607 /* Now, we'd save a lot of time by not calculating eq10v2 for all
608 * time points, but only those that we sample to calculate the mse.
611 eq10v2(GD->ctTheory, GD->doubleLogTime /* GD->time */, nFitPoints /* GD->nData */,
612 gsl_vector_get(p, 0), gsl_vector_get(p, 1),
615 fixGemACF(GD->ctTheory, nFitPoints);
617 /* Removing a bunch of points from the log-part. */
618 #pragma omp parallel for schedule(dynamic) \
619 firstprivate(nData, ctTheory, y, nFitPoints) \
620 private (i, iLog, r) \
621 reduction(+:residual2) \
623 for (i = 0; i < nFitPoints; i++)
625 iLog = GD->logtime[i];
626 r = log(ctTheory[i /* iLog */]);
627 residual2 += sqr(r-log(y[iLog]));
629 residual2 /= nFitPoints; /* Not really necessary for the fitting, but gives more meaning to the returned data. */
630 /* printf("ka=%3.5f\tkd=%3.5f\trmse=%3.5f\n", gsl_vector_get(p,0), gsl_vector_get(p,1), residual2); */
633 /* for (i=0; i<nLin; i++) { */
634 /* /\* Linear part ----------*\/ */
635 /* r = ctTheory[i]; */
636 /* residual2 += sqr(r-y[i]); */
637 /* /\* Log part -------------*\/ */
638 /* iLog = GETLOGINDEX(i, GD->params); */
639 /* if (iLog >= nData) */
640 /* gmx_fatal(FARGS, "log index out of bounds: %i", iLog); */
641 /* r = ctTheory[iLog]; */
642 /* residual2 += sqr(r-y[iLog]); */
645 /* residual2 /= GD->n; */
646 /* return residual2; */
650 static real* d2r(const double *d, const int nn);
652 extern real fitGemRecomb(double gmx_unused *ct,
653 double gmx_unused *time,
654 double gmx_unused **ctFit,
655 const int gmx_unused nData,
656 t_gemParams gmx_unused *params)
659 int nThreads, i, iter, status, maxiter;
660 real size, d2, tol, *dumpdata;
663 char *dumpstr, dumpname[128];
665 /* nmsimplex2 had convergence problems prior to gsl v1.14,
666 * but it's O(N) instead of O(N) operations, so let's use it if v >= 1.14 */
668 gsl_multimin_fminimizer *s;
669 gsl_vector *x, *dx; /* parameters and initial step size */
670 gsl_multimin_function fitFunc;
671 #ifdef GSL_MAJOR_VERSION
672 #ifdef GSL_MINOR_VERSION
673 #if ((GSL_MAJOR_VERSION == 1 && GSL_MINOR_VERSION >= 14) || \
674 (GSL_MAJOR_VERSION > 1))
675 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex2;
677 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
679 #endif /* GSL_MINOR_VERSION */
681 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
682 #endif /* GSL_MAJOR_VERSION */
683 fprintf(stdout, "Will fit ka and kd to the ACF according to the reversible geminate recombination model.\n");
684 #else /* HAVE_LIBGSL */
685 fprintf(stderr, "Sorry, can't do reversible geminate recombination without gsl. "
686 "Recompile using --with-gsl.\n");
688 #endif /* HAVE_LIBGSL */
692 nThreads = gmx_omp_get_max_threads();
693 gmx_omp_set_num_threads(nThreads);
694 fprintf(stdout, "We will be using %i threads.\n", nThreads);
702 p = 2; /* Number of parameters to fit. ka and kd. */
703 n = params->nFitPoints; /* params->nLin*2 */; /* Number of points in the reduced dataset */
707 fprintf(stderr, "Fitting of D is not implemented yet. It must be provided on the command line.\n");
712 /* fprintf(stderr, "Reduced data set larger than the complete data set!\n"); */
715 snew(dumpdata, nData);
721 snew(GD->ctTheory, nData);
727 GD->tDelta = time[1]-time[0];
730 snew(GD->logtime, params->nFitPoints);
731 snew(GD->doubleLogTime, params->nFitPoints);
733 for (i = 0; i < params->nFitPoints; i++)
735 GD->doubleLogTime[i] = (double)(getLogIndex(i, params));
736 GD->logtime[i] = (int)(GD->doubleLogTime[i]);
737 GD->doubleLogTime[i] *= GD->tDelta;
739 if (GD->logtime[i] >= nData)
741 fprintf(stderr, "Ayay. It seems we're indexing out of bounds.\n");
742 params->nFitPoints = i;
746 fitFunc.f = &gemFunc_residual2;
748 fitFunc.params = (void*)GD;
750 x = gsl_vector_alloc (fitFunc.n);
751 dx = gsl_vector_alloc (fitFunc.n);
752 gsl_vector_set (x, 0, 25);
753 gsl_vector_set (x, 1, 0.5);
754 gsl_vector_set (dx, 0, 0.1);
755 gsl_vector_set (dx, 1, 0.01);
758 s = gsl_multimin_fminimizer_alloc (T, fitFunc.n);
759 gsl_multimin_fminimizer_set (s, &fitFunc, x, dx);
761 gsl_vector_free (dx);
766 status = gsl_multimin_fminimizer_iterate (s);
770 gmx_fatal(FARGS, "Something went wrong in the iteration in minimizer %s:\n \"%s\"\n",
771 gsl_multimin_fminimizer_name(s), gsl_strerror(status));
774 d2 = gsl_multimin_fminimizer_minimum(s);
775 size = gsl_multimin_fminimizer_size(s);
776 params->ka = gsl_vector_get (s->x, 0);
777 params->kd = gsl_vector_get (s->x, 1);
781 fprintf(stderr, "%s\n", gsl_strerror(status));
785 status = gsl_multimin_test_size(size, tol);
787 if (status == GSL_SUCCESS)
789 fprintf(stdout, "Converged to minimum at\n");
792 printf ("iter %5d: ka = %2.5f kd = %2.5f f() = %7.3f size = %.3f chi2 = %2.5f\n",
800 eq10v2(GD->ctTheory, time, nData, params->ka, params->kd, params);
801 /* fixGemACF(GD->ctTheory, nFitPoints); */
802 sprintf(dumpname, "Iter_%i.xvg", iter);
803 for (i = 0; i < GD->nData; i++)
805 dumpdata[i] = (real)(GD->ctTheory[i]);
806 if (!gmx_isfinite(dumpdata[i]))
808 gmx_fatal(FARGS, "Non-finite value in acf.");
811 dumpN(dumpdata, GD->nData, dumpname);
814 while ((status == GSL_CONTINUE) && (iter < maxiter));
816 /* /\* Calculate the theoretical ACF from the parameters one last time. *\/ */
817 eq10v2(GD->ctTheory, time, nData, params->ka, params->kd, params);
818 *ctFit = GD->ctTheory;
821 gsl_multimin_fminimizer_free (s);
826 #endif /* HAVE_LIBGSL */
830 static int balFunc_f(const gsl_vector *x, void *data, gsl_vector *f)
832 /* C + sum{ A_i * exp(-B_i * t) }*/
836 double *y, *A, *B, C, /* There are the parameters to be optimized. */
839 BD = (balData *)data;
846 for (i = 0; i < nexp; i++)
848 A[i] = gsl_vector_get(x, i*2);
849 B[i] = gsl_vector_get(x, i*2+1);
851 C = gsl_vector_get(x, nexp*2);
853 for (i = 0; i < n; i++)
857 for (j = 0; j < nexp; j++)
859 ct += A[j] * exp(B[j] * t);
862 gsl_vector_set (f, i, ct - y[i]);
867 /* The derivative stored in jacobian form (J)*/
868 static int balFunc_df(const gsl_vector *params, void *data, gsl_matrix *J)
872 double *y, *A, *B, C, /* There are the parameters. */
884 for (i = 0; i < nexp; i++)
886 A[i] = gsl_vector_get(params, i*2);
887 B[i] = gsl_vector_get(params, i*2+1);
889 C = gsl_vector_get(params, nexp*2);
890 for (i = 0; i < n; i++)
893 for (j = 0; j < nexp; j++)
895 gsl_matrix_set (J, i, j*2, exp(B[j]*t)); /* df(t)/dA_j */
896 gsl_matrix_set (J, i, j*2+1, A[j]*t*exp(B[j]*t)); /* df(t)/dB_j */
898 gsl_matrix_set (J, i, nexp*2, 1); /* df(t)/dC */
903 /* Calculation of the function and its derivative */
904 static int balFunc_fdf(const gsl_vector *params, void *data,
905 gsl_vector *f, gsl_matrix *J)
907 balFunc_f(params, data, f);
908 balFunc_df(params, data, J);
911 #endif /* HAVE_LIBGSL */
913 /* Removes the ballistic term from the beginning of the ACF,
914 * just like in Omer's paper.
916 extern void takeAwayBallistic(double gmx_unused *ct, double *t, int len, real tMax, int nexp, gmx_bool gmx_unused bDerivative)
919 /* Use nonlinear regression with GSL instead.
920 * Fit with 4 exponentials and one constant term,
921 * subtract the fatest exponential. */
923 int nData, i, status, iter;
925 double *guess, /* Initial guess. */
926 *A, /* The fitted parameters. (A1, B1, A2, B2,... C) */
938 while (t[nData] < tMax+t[0] && nData < len);
940 p = nexp*2+1; /* Number of parameters. */
943 const gsl_multifit_fdfsolver_type *T
944 = gsl_multifit_fdfsolver_lmsder;
946 gsl_multifit_fdfsolver *s; /* The solver itself. */
947 gsl_multifit_function_fdf fitFunction; /* The function to be fitted. */
948 gsl_matrix *covar; /* Covariance matrix for the parameters.
949 * We'll not use the result, though. */
950 gsl_vector_view theParams;
957 while (t[nData] < tMax+t[0] && nData < len);
964 covar = gsl_matrix_alloc (p, p);
966 /* Set up an initial gess for the parameters.
967 * The solver is somewhat sensitive to the initial guess,
968 * but this worked fine for a TIP3P box with -geminate dd
969 * EDIT: In fact, this seems like a good starting pont for other watermodels too. */
970 for (i = 0; i < nexp; i++)
973 guess[i*2+1] = -0.5 + (((double)i)/nexp - 0.5)*0.3;
975 guess[nexp * 2] = 0.01;
977 theParams = gsl_vector_view_array(guess, p);
982 BD->tDelta = t[1]-t[0];
985 fitFunction.f = &balFunc_f;
986 fitFunction.df = &balFunc_df;
987 fitFunction.fdf = &balFunc_fdf;
988 fitFunction.n = nData;
990 fitFunction.params = BD;
992 s = gsl_multifit_fdfsolver_alloc (T, nData, p);
995 gmx_fatal(FARGS, "Could not set up the nonlinear solver.");
998 gsl_multifit_fdfsolver_set(s, &fitFunction, &theParams.vector);
1000 /* \=============================================/ */
1006 status = gsl_multifit_fdfsolver_iterate (s);
1012 status = gsl_multifit_test_delta (s->dx, s->x, 1e-4, 1e-4);
1014 while (iter < 5000 && status == GSL_CONTINUE);
1018 fprintf(stderr, "The non-linear fitting did not converge in 5000 steps.\n"
1019 "Check the quality of the fit!\n");
1023 fprintf(stderr, "Non-linear fitting of ballistic term converged in %d steps.\n\n", (int)iter);
1025 for (i = 0; i < nexp; i++)
1027 fprintf(stdout, "%c * exp(%c * t) + ", 'A'+(char)i*2, 'B'+(char)i*2);
1030 fprintf(stdout, "%c\n", 'A'+(char)nexp*2);
1031 fprintf(stdout, "Here are the actual numbers for A-%c:\n", 'A'+nexp*2);
1033 for (i = 0; i < nexp; i++)
1035 A[i*2] = gsl_vector_get(s->x, i*2);
1036 A[i*2+1] = gsl_vector_get(s->x, i*2+1);
1037 fprintf(stdout, " %g*exp(%g * x) +", A[i*2], A[i*2+1]);
1039 A[i*2] = gsl_vector_get(s->x, i*2); /* The last and constant term */
1040 fprintf(stdout, " %g\n", A[i*2]);
1044 /* Implement some check for parameter quality */
1045 for (i = 0; i < nexp; i++)
1047 if (A[i*2] < 0 || A[i*2] > 1)
1049 fprintf(stderr, "WARNING: ----------------------------------\n"
1050 " | A coefficient does not lie within [0,1].\n"
1051 " | This may or may not be a problem.\n"
1052 " | Double check the quality of the fit!\n");
1056 fprintf(stderr, "WARNING: ----------------------------------\n"
1057 " | One factor in the exponent is positive.\n"
1058 " | This could be a problem if the coefficient\n"
1059 " | is large. Double check the quality of the fit!\n");
1062 if (A[i*2] < 0 || A[i*2] > 1)
1064 fprintf(stderr, "WARNING: ----------------------------------\n"
1065 " | The constant term does not lie within [0,1].\n"
1066 " | This may or may not be a problem.\n"
1067 " | Double check the quality of the fit!\n");
1070 /* Sort the terms */
1071 sorted = (nexp > 1) ? FALSE : TRUE;
1075 for (i = 0; i < nexp-1; i++)
1077 ddt[0] = A[i*2] * A[i*2+1];
1078 ddt[1] = A[i*2+2] * A[i*2+3];
1080 if ((bDerivative && (ddt[0] < 0 && ddt[1] < 0 && ddt[0] > ddt[1])) || /* Compare derivative at t=0... */
1081 (!bDerivative && (A[i*2+1] > A[i*2+3]))) /* Or just the coefficient in the exponent */
1084 a[0] = A[i*2]; /* coefficient */
1085 a[1] = A[i*2+1]; /* parameter in the exponent */
1088 A[i*2+1] = A[i*2+3];
1096 /* Subtract the fastest component */
1097 fprintf(stdout, "Fastest component is %g * exp(%g * t)\n"
1098 "Subtracting fastest component from ACF.\n", A[0], A[1]);
1100 for (i = 0; i < len; i++)
1102 ct[i] = (ct[i] - A[0] * exp(A[1] * i*BD->tDelta)) / (1-A[0]);
1108 gsl_multifit_fdfsolver_free(s);
1109 gsl_matrix_free(covar);
1113 /* We have no gsl. */
1114 fprintf(stderr, "Sorry, can't take away ballistic component without gsl. "
1115 "Recompile using --with-gsl.\n");
1117 #endif /* HAVE_LIBGSL */
1121 extern void dumpN(const real *e, const int nn, char *fn)
1123 /* For debugging only */
1126 char standardName[] = "Nt.xvg";
1135 "@ xaxis label \"Frame\"\n"
1136 "@ yaxis label \"N\"\n"
1137 "@ s0 line type 3\n");
1139 for (i = 0; i < nn; i++)
1141 fprintf(f, "%-10i %-g\n", i, e[i]);
1147 static real* d2r(const double *d, const int nn)
1153 for (i = 0; i < nn; i++)
1161 static void _patchBad(double *ct, int n, double dy)
1163 /* Just do lin. interpolation for now. */
1166 for (i = 1; i < n; i++)
1172 static void patchBadPart(double *ct, int n)
1174 _patchBad(ct, n, (ct[n] - ct[0])/n);
1177 static void patchBadTail(double *ct, int n)
1179 _patchBad(ct+1, n-1, ct[1]-ct[0]);
1183 extern void fixGemACF(double *ct, int len)
1188 /* Let's separate two things:
1189 * - identification of bad parts
1190 * - patching of bad parts.
1193 b = 0; /* Start of a bad stretch */
1194 e = 0; /* End of a bad stretch */
1197 /* An acf of binary data must be one at t=0. */
1198 if (abs(ct[0]-1.0) > 1e-6)
1201 fprintf(stderr, "|ct[0]-1.0| = %1.6d. Setting ct[0] to 1.0.\n", abs(ct[0]-1.0));
1204 for (i = 0; i < len; i++)
1208 if (isfinite(ct[i]))
1209 #elif defined(HAS__ISFINITE)
1210 if (_isfinite(ct[i]))
1217 /* Still on a good stretch. Proceed.*/
1221 /* Patch up preceding bad stretch. */
1227 gmx_fatal(FARGS, "The ACF is mostly NaN or Inf. Aborting.");
1229 patchBadTail(&(ct[b-2]), (len-b)+1);
1233 patchBadPart(&(ct[b-1]), (e-b)+1);