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>
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) { return (sqrt(z.r*z.r+z.i*z.i)); }
93 /* Addition of two complex numbers z1 and z2 */
94 static gem_complex gem_cxadd(gem_complex z1, gem_complex z2)
102 /* Addition of a complex number z1 and a real number r */
103 static gem_complex gem_cxradd(gem_complex z, double r)
111 /* Subtraction of two complex numbers z1 and z2 */
112 static gem_complex gem_cxsub(gem_complex z1, gem_complex z2)
120 /* Multiplication of two complex numbers z1 and z2 */
121 static gem_complex gem_cxmul(gem_complex z1, gem_complex z2)
124 value.r = z1.r*z2.r-z1.i*z2.i;
125 value.i = z1.r*z2.i+z1.i*z2.r;
129 /* Square of a complex number z */
130 static gem_complex gem_cxsq(gem_complex z)
133 value.r = z.r*z.r-z.i*z.i;
134 value.i = z.r*z.i*2.;
138 /* multiplication of a complex number z and a real number r */
139 static gem_complex gem_cxrmul(gem_complex z, double r)
147 /* Division of two complex numbers z1 and z2 */
148 static gem_complex gem_cxdiv(gem_complex z1, gem_complex z2)
152 num = z2.r*z2.r+z2.i*z2.i;
155 fprintf(stderr, "ERROR in gem_cxdiv function\n");
157 value.r = (z1.r*z2.r+z1.i*z2.i)/num; value.i = (z1.i*z2.r-z1.r*z2.i)/num;
161 /* division of a complex z number by a real number x */
162 static gem_complex gem_cxrdiv(gem_complex z, double r)
170 /* division of a real number r by a complex number x */
171 static gem_complex gem_rxcdiv(double r, gem_complex z)
175 f = r/(z.r*z.r+z.i*z.i);
181 /* Exponential of a complex number-- exp (z)=|exp(z.r)|*{cos(z.i)+I*sin(z.i)}*/
182 static gem_complex gem_cxdexp(gem_complex z)
187 value.r = exp_z_r*cos(z.i);
188 value.i = exp_z_r*sin(z.i);
192 /* Logarithm of a complex number -- log(z)=log|z|+I*Arg(z), */
193 /* where -PI < Arg(z) < PI */
194 static gem_complex gem_cxlog(gem_complex z)
198 mag2 = z.r*z.r+z.i*z.i;
200 fprintf(stderr, "ERROR in gem_cxlog func\n");
202 value.r = log(sqrt(mag2));
209 value.i = atan2(z.i, z.r);
214 /* Square root of a complex number z = |z| exp(I*the) -- z^(1/2) */
215 /* z^(1/2)=|z|^(1/2)*[cos(the/2)+I*sin(the/2)] */
216 /* where 0 < the < 2*PI */
217 static gem_complex gem_cxdsqrt(gem_complex z)
222 value.r = sqrt(fabs((sq+z.r)*0.5)); /* z'.r={|z|*[1+cos(the)]/2}^(1/2) */
223 value.i = sqrt(fabs((sq-z.r)*0.5)); /* z'.i={|z|*[1-cos(the)]/2}^(1/2) */
230 /* Complex power of a complex number z1^z2 */
231 static gem_complex gem_cxdpow(gem_complex z1, gem_complex z2)
234 value=gem_cxdexp(gem_cxmul(gem_cxlog(z1),z2));
238 /* ------------ end of complex.c ------------ */
240 /* This next part was derived from cubic.c, also received from Omer Markovitch.
241 * ------------- from cubic.c ------------- */
243 /* Solver for a cubic equation: x^3-a*x^2+b*x-c=0 */
244 static void gem_solve(gem_complex *al, gem_complex *be, gem_complex *gam,
245 double a, double b, double c)
247 double t1, t2, two_3, temp;
248 gem_complex ctemp, ct3;
250 two_3=pow(2., 1./3.); t1 = -a*a+3.*b; t2 = 2.*a*a*a-9.*a*b+27.*c;
251 temp = 4.*t1*t1*t1+t2*t2;
253 ctemp = gem_cmplx(temp,0.); ctemp = gem_cxadd(gem_cmplx(t2,0.),gem_cxdsqrt(ctemp));
254 ct3 = gem_cxdpow(ctemp, gem_cmplx(1./3.,0.));
256 ctemp = gem_rxcdiv(-two_3*t1/3., ct3);
257 ctemp = gem_cxadd(ctemp, gem_cxrdiv(ct3, 3.*two_3));
259 *gam = gem_cxadd(gem_cmplx(a/3.,0.), ctemp);
261 ctemp=gem_cxmul(gem_cxsq(*gam), gem_cxsq(gem_cxsub(*gam, gem_cmplx(a,0.))));
262 ctemp=gem_cxadd(ctemp, gem_cxmul(gem_cmplx(-4.*c,0.), *gam));
263 ctemp = gem_cxdiv(gem_cxdsqrt(ctemp), *gam);
264 *al = gem_cxrmul(gem_cxsub(gem_cxsub(gem_cmplx(a,0.), *gam),ctemp),0.5);
265 *be = gem_cxrmul(gem_cxadd(gem_cxsub(gem_cmplx(a,0.), *gam),ctemp),0.5);
268 /* ------------ end of cubic.c ------------ */
270 /* This next part was derived from cerror.c and rerror.c, also received from Omer Markovitch.
271 * ------------- from [cr]error.c ------------- */
273 /************************************************************/
274 /* Real valued error function and related functions */
275 /* x, y : real variables */
276 /* erf(x) : error function */
277 /* erfc(x) : complementary error function */
278 /* omega(x) : exp(x*x)*erfc(x) */
279 /* W(x,y) : exp(-x*x)*omega(x+y)=exp(2*x*y+y^2)*erfc(x+y) */
280 /************************************************************/
282 /*---------------------------------------------------------------------------*/
283 /* Utilzed the series approximation of erf(x) */
284 /* Relative error=|err(x)|/erf(x)<EPS */
285 /* Handbook of Mathematical functions, Abramowitz, p 297 */
286 /* Note: When x>=6 series sum deteriorates -> Asymptotic series used instead */
287 /*---------------------------------------------------------------------------*/
288 /* This gives erfc(z) correctly only upto >10-15 */
290 static double gem_erf(double x)
292 double n,sum,temp,exp2,xm,x2,x4,x6,x8,x10,x12;
305 for(n=1.; n<=2000.; n+=1.){
306 temp *= 2.*x2/(2.*n+1.);
308 if(fabs(temp/sum)<1.E-16)
313 fprintf(stderr, "In Erf calc - iteration exceeds %lg\n",n);
318 /* from the asymptotic expansion of experfc(x) */
319 sum = (1. - 0.5/x2 + 0.75/x4
320 - 1.875/x6 + 6.5625/x8
321 - 29.53125/x10 + 162.421875/x12)
323 sum*=exp2; /* now sum is erfc(x) */
326 return x>=0.0 ? sum : -sum;
329 /* Result --> Alex's code for erfc and experfc behaves better than this */
330 /* complementray error function. Returns 1.-erf(x) */
331 static double gem_erfc(double x)
337 ans = t * exp(-z*z - 1.26551223 +
346 t*0.17087277)))))))));
348 return x>=0.0 ? ans : 2.0-ans;
351 /* omega(x)=exp(x^2)erfc(x) */
352 static double gem_omega(double x)
354 double xm, ans, xx, x4, x6, x8, x10, x12;
364 ans = exp(xx)*gem_erfc(x);
366 /* Asymptotic expansion */
367 ans = (1. - 0.5/xx + 0.75/x4 - 1.875/x6 + 6.5625/x8 - 29.53125/x10 + 162.421875/x12) / sPI/x;
372 /*---------------------------------------------------------------------------*/
373 /* Utilzed the series approximation of erf(z=x+iy) */
374 /* Relative error=|err(z)|/|erf(z)|<EPS */
375 /* Handbook of Mathematical functions, Abramowitz, p 299 */
376 /* comega(z=x+iy)=cexp(z^2)*cerfc(z) */
377 /*---------------------------------------------------------------------------*/
378 static gem_complex gem_comega(gem_complex z)
382 double sumr,sumi,n,n2,f,temp,temp1;
383 double x2, y2,cos_2xy,sin_2xy,cosh_2xy,sinh_2xy,cosh_ny,sinh_ny,exp_y2;
391 cos_2xy = cos(2.*x*y);
392 sin_2xy = sin(2.*x*y);
393 cosh_2xy = cosh(2.*x*y);
394 sinh_2xy = sinh(2.*x*y);
397 for(n=1.0, temp=0.; n<=2000.; n+=1.0){
401 f = exp(-n2/4.)/(n2+4.*x2);
402 /* if(f<1.E-200) break; */
403 sumr += (2.*x - 2.*x*cosh_ny*cos_2xy + n*sinh_ny*sin_2xy)*f;
404 sumi += (2.*x*cosh_ny*sin_2xy + n*sinh_ny*cos_2xy)*f;
405 temp1 = sqrt(sumr*sumr+sumi*sumi);
406 if(fabs((temp1-temp)/temp1)<1.E-16) {
412 fprintf(stderr, "iteration exceeds %lg\n",n);
422 value.r = gem_omega(x)-(f*(1.-cos_2xy)+sumr);
423 value.i = -(f*sin_2xy+sumi);
424 value = gem_cxmul(value,gem_cmplx(exp_y2*cos_2xy,exp_y2*sin_2xy));
428 /* ------------ end of [cr]error.c ------------ */
430 /*_ REVERSIBLE GEMINATE RECOMBINATION
432 * Here are the functions for reversible geminate recombination. */
434 /* Changes the unit from square cm per s to square Ångström per ps,
435 * since Omers code uses the latter units while g_mds outputs the former.
436 * g_hbond expects a diffusion coefficent given in square cm per s. */
437 static double sqcm_per_s_to_sqA_per_ps (real D) {
438 fprintf(stdout, "Diffusion coefficient is %f A^2/ps\n", D*1e4);
439 return (double)(D*1e4);
443 static double eq10v2(double theoryCt[], double time[], int manytimes,
444 double ka, double kd, t_gemParams *params)
446 /* Finding the 3 roots */
448 double kD, D, r, a, b, c, tsqrt, sumimaginary;
453 part1, part2, part3, part4;
458 a = (1.0 + ka/kD) * sqrt(D)/r;
462 gem_solve(&alpha, &beta, &gamma, a, b, c);
463 /* Finding the 3 roots */
466 part1 = gem_cxmul(alpha, gem_cxmul(gem_cxadd(beta, gamma), gem_cxsub(beta, gamma))); /* 1(2+3)(2-3) */
467 part2 = gem_cxmul(beta, gem_cxmul(gem_cxadd(gamma, alpha), gem_cxsub(gamma, alpha))); /* 2(3+1)(3-1) */
468 part3 = gem_cxmul(gamma, gem_cxmul(gem_cxadd(alpha, beta) , gem_cxsub(alpha, beta))); /* 3(1+2)(1-2) */
469 part4 = gem_cxmul(gem_cxsub(gamma, alpha), gem_cxmul(gem_cxsub(alpha, beta), gem_cxsub(beta, gamma))); /* (3-1)(1-2)(2-3) */
471 #pragma omp parallel for \
472 private(i, tsqrt, oma, omb, omc, c1, c2, c3, c4), \
473 reduction(+:sumimaginary), \
476 for (i=0; i<manytimes; i++){
477 tsqrt = sqrt(time[i]);
478 oma = gem_comega(gem_cxrmul(alpha, tsqrt));
479 c1 = gem_cxmul(oma, gem_cxdiv(part1, part4));
480 omb = gem_comega(gem_cxrmul(beta, tsqrt));
481 c2 = gem_cxmul(omb, gem_cxdiv(part2, part4));
482 omc = gem_comega(gem_cxrmul(gamma, tsqrt));
483 c3 = gem_cxmul(omc, gem_cxdiv(part3, part4));
484 c4.r = c1.r+c2.r+c3.r;
485 c4.i = c1.i+c2.i+c3.i;
487 sumimaginary += c4.i * c4.i;
494 /* This returns the real-valued index(!) to an ACF, equidistant on a log scale. */
495 static double getLogIndex(const int i, const t_gemParams *params)
497 return (exp(((double)(i)) * params->logQuota) -1);
500 extern t_gemParams *init_gemParams(const double sigma, const double D,
501 const real *t, const int len, const int nFitPoints,
502 const real begFit, const real endFit,
503 const real ballistic, const int nBalExp, const gmx_bool bDt)
510 /* A few hardcoded things here. For now anyway. */
512 /* p->ka_max = 100; */
520 /* p->lifetime = 0; */
521 p->sigma = sigma*10; /* Omer uses Ã…, not nm */
522 /* p->lsq_old = 99999; */
523 p->D = sqcm_per_s_to_sqA_per_ps(D);
524 p->kD = 4*acos(-1.0)*sigma*p->D;
527 /* Parameters used by calcsquare(). Better to calculate them
528 * here than in calcsquare every time it's called. */
530 /* p->logAfterTime = logAfterTime; */
531 tDelta = (t[len-1]-t[0]) / len;
533 gmx_fatal(FARGS, "Time between frames is non-positive!");
538 p->nExpFit = nBalExp;
539 /* p->nLin = logAfterTime / tDelta; */
540 p->nFitPoints = nFitPoints;
543 p->logQuota = (double)(log(p->len))/(p->nFitPoints-1);
544 /* if (p->nLin <= 0) { */
545 /* fprintf(stderr, "Number of data points in the linear regime is non-positive!\n"); */
549 /* We want the same number of data points in the log regime. Not crucial, but seems like a good idea. */
550 /* p->logDelta = log(((float)len)/p->nFitPoints) / p->nFitPoints;/\* log(((float)len)/p->nLin) / p->nLin; *\/ */
551 /* p->logPF = p->nFitPoints*p->nFitPoints/(float)len; /\* p->nLin*p->nLin/(float)len; *\/ */
552 /* logPF and logDelta are stitched together with the macro GETLOGINDEX defined in geminate.h */
554 /* p->logMult = pow((float)len, 1.0/nLin);/\* pow(t[len-1]-t[0], 1.0/p->nLin); *\/ */
555 p->ballistic = ballistic;
559 /* There was a misunderstanding regarding the fitting. From our
560 * recent correspondence it appears that Omer's code require
561 * the ACF data on a log-scale and does not operate on the raw data.
562 * This needs to be redone in gemFunc_residual() as well as in the
563 * t_gemParams structure. */
565 static double gemFunc_residual2(const gsl_vector *p, void *data)
568 int i, iLog, nLin, nFitPoints, nData;
569 double r, residual2, *ctTheory, *y;
571 GD = (gemFitData *)data;
572 nLin = GD->params->nLin;
573 nFitPoints = GD->params->nFitPoints;
576 ctTheory = GD->ctTheory;
579 /* Now, we'd save a lot of time by not calculating eq10v2 for all
580 * time points, but only those that we sample to calculate the mse.
583 eq10v2(GD->ctTheory, GD->doubleLogTime/* GD->time */, nFitPoints/* GD->nData */,
584 gsl_vector_get(p, 0), gsl_vector_get(p, 1),
587 fixGemACF(GD->ctTheory, nFitPoints);
589 /* Removing a bunch of points from the log-part. */
590 #pragma omp parallel for schedule(dynamic) \
591 firstprivate(nData, ctTheory, y, nFitPoints) \
592 private (i, iLog, r) \
593 reduction(+:residual2) \
595 for(i=0; i<nFitPoints; i++)
597 iLog = GD->logtime[i];
598 r = log(ctTheory[i /* iLog */]);
599 residual2 += sqr(r-log(y[iLog]));
601 residual2 /= nFitPoints; /* Not really necessary for the fitting, but gives more meaning to the returned data. */
602 /* printf("ka=%3.5f\tkd=%3.5f\trmse=%3.5f\n", gsl_vector_get(p,0), gsl_vector_get(p,1), residual2); */
605 /* for (i=0; i<nLin; i++) { */
606 /* /\* Linear part ----------*\/ */
607 /* r = ctTheory[i]; */
608 /* residual2 += sqr(r-y[i]); */
609 /* /\* Log part -------------*\/ */
610 /* iLog = GETLOGINDEX(i, GD->params); */
611 /* if (iLog >= nData) */
612 /* gmx_fatal(FARGS, "log index out of bounds: %i", iLog); */
613 /* r = ctTheory[iLog]; */
614 /* residual2 += sqr(r-y[iLog]); */
617 /* residual2 /= GD->n; */
618 /* return residual2; */
622 static real* d2r(const double *d, const int nn);
624 extern real fitGemRecomb(double *ct, double *time, double **ctFit,
625 const int nData, t_gemParams *params)
628 int nThreads, i, iter, status, maxiter;
629 real size, d2, tol, *dumpdata;
632 char *dumpstr, dumpname[128];
634 /* nmsimplex2 had convergence problems prior to gsl v1.14,
635 * but it's O(N) instead of O(N) operations, so let's use it if v >= 1.14 */
637 gsl_multimin_fminimizer *s;
638 gsl_vector *x,*dx; /* parameters and initial step size */
639 gsl_multimin_function fitFunc;
640 #ifdef GSL_MAJOR_VERSION
641 #ifdef GSL_MINOR_VERSION
642 #if ((GSL_MAJOR_VERSION == 1 && GSL_MINOR_VERSION >= 14) || \
643 (GSL_MAJOR_VERSION > 1))
644 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex2;
646 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
648 #endif /* GSL_MINOR_VERSION */
650 const gsl_multimin_fminimizer_type *T = gsl_multimin_fminimizer_nmsimplex;
651 #endif /* GSL_MAJOR_VERSION */
652 fprintf(stdout, "Will fit ka and kd to the ACF according to the reversible geminate recombination model.\n");
653 #else /* HAVE_LIBGSL */
654 fprintf(stderr, "Sorry, can't do reversible geminate recombination without gsl. "
655 "Recompile using --with-gsl.\n");
657 #endif /* HAVE_LIBGSL */
661 nThreads = gmx_omp_get_max_threads();
662 gmx_omp_set_num_threads(nThreads);
663 fprintf(stdout, "We will be using %i threads.\n", nThreads);
671 p = 2; /* Number of parameters to fit. ka and kd. */
672 n = params->nFitPoints; /* params->nLin*2 */; /* Number of points in the reduced dataset */
676 fprintf(stderr, "Fitting of D is not implemented yet. It must be provided on the command line.\n");
681 /* fprintf(stderr, "Reduced data set larger than the complete data set!\n"); */
684 snew(dumpdata, nData);
690 snew(GD->ctTheory, nData);
696 GD->tDelta = time[1]-time[0];
699 snew(GD->logtime,params->nFitPoints);
700 snew(GD->doubleLogTime,params->nFitPoints);
702 for (i=0; i<params->nFitPoints; i++)
704 GD->doubleLogTime[i] = (double)(getLogIndex(i, params));
705 GD->logtime[i] = (int)(GD->doubleLogTime[i]);
706 GD->doubleLogTime[i]*=GD->tDelta;
708 if (GD->logtime[i] >= nData)
710 fprintf(stderr, "Ayay. It seems we're indexing out of bounds.\n");
711 params->nFitPoints = i;
715 fitFunc.f = &gemFunc_residual2;
717 fitFunc.params = (void*)GD;
719 x = gsl_vector_alloc (fitFunc.n);
720 dx = gsl_vector_alloc (fitFunc.n);
721 gsl_vector_set (x, 0, 25);
722 gsl_vector_set (x, 1, 0.5);
723 gsl_vector_set (dx, 0, 0.1);
724 gsl_vector_set (dx, 1, 0.01);
727 s = gsl_multimin_fminimizer_alloc (T, fitFunc.n);
728 gsl_multimin_fminimizer_set (s, &fitFunc, x, dx);
730 gsl_vector_free (dx);
734 status = gsl_multimin_fminimizer_iterate (s);
737 gmx_fatal(FARGS,"Something went wrong in the iteration in minimizer %s:\n \"%s\"\n",
738 gsl_multimin_fminimizer_name(s), gsl_strerror(status));
740 d2 = gsl_multimin_fminimizer_minimum(s);
741 size = gsl_multimin_fminimizer_size(s);
742 params->ka = gsl_vector_get (s->x, 0);
743 params->kd = gsl_vector_get (s->x, 1);
747 fprintf(stderr, "%s\n", gsl_strerror(status));
751 status = gsl_multimin_test_size(size,tol);
753 if (status == GSL_SUCCESS) {
754 fprintf(stdout, "Converged to minimum at\n");
757 printf ("iter %5d: ka = %2.5f kd = %2.5f f() = %7.3f size = %.3f chi2 = %2.5f\n",
765 eq10v2(GD->ctTheory, time, nData, params->ka, params->kd, params);
766 /* fixGemACF(GD->ctTheory, nFitPoints); */
767 sprintf(dumpname, "Iter_%i.xvg", iter);
768 for(i=0; i<GD->nData; i++)
770 dumpdata[i] = (real)(GD->ctTheory[i]);
771 if (!gmx_isfinite(dumpdata[i]))
773 gmx_fatal(FARGS, "Non-finite value in acf.");
776 dumpN(dumpdata, GD->nData, dumpname);
779 while ((status == GSL_CONTINUE) && (iter < maxiter));
781 /* /\* Calculate the theoretical ACF from the parameters one last time. *\/ */
782 eq10v2(GD->ctTheory, time, nData, params->ka, params->kd, params);
783 *ctFit = GD->ctTheory;
786 gsl_multimin_fminimizer_free (s);
791 #endif /* HAVE_LIBGSL */
795 static int balFunc_f(const gsl_vector *x, void *data, gsl_vector *f)
797 /* C + sum{ A_i * exp(-B_i * t) }*/
801 double *y, *A, *B, C, /* There are the parameters to be optimized. */
804 BD = (balData *)data;
811 for (i = 0; i<nexp; i++)
813 A[i] = gsl_vector_get(x, i*2);
814 B[i] = gsl_vector_get(x, i*2+1);
816 C = gsl_vector_get(x, nexp*2);
822 for (j=0; j<nexp; j++) {
823 ct += A[j] * exp(B[j] * t);
826 gsl_vector_set (f, i, ct - y[i]);
831 /* The derivative stored in jacobian form (J)*/
832 static int balFunc_df(const gsl_vector *params, void *data, gsl_matrix *J)
836 double *y, *A, *B, C, /* There are the parameters. */
848 for (i=0; i<nexp; i++)
850 A[i] = gsl_vector_get(params, i*2);
851 B[i] = gsl_vector_get(params, i*2+1);
853 C = gsl_vector_get(params, nexp*2);
857 for (j=0; j<nexp; j++)
859 gsl_matrix_set (J, i, j*2, exp(B[j]*t)); /* df(t)/dA_j */
860 gsl_matrix_set (J, i, j*2+1, A[j]*t*exp(B[j]*t)); /* df(t)/dB_j */
862 gsl_matrix_set (J, i, nexp*2, 1); /* df(t)/dC */
867 /* Calculation of the function and its derivative */
868 static int balFunc_fdf(const gsl_vector *params, void *data,
869 gsl_vector *f, gsl_matrix *J)
871 balFunc_f(params, data, f);
872 balFunc_df(params, data, J);
875 #endif /* HAVE_LIBGSL */
877 /* Removes the ballistic term from the beginning of the ACF,
878 * just like in Omer's paper.
880 extern void takeAwayBallistic(double *ct, double *t, int len, real tMax, int nexp, gmx_bool bDerivative)
883 /* Use nonlinear regression with GSL instead.
884 * Fit with 4 exponentials and one constant term,
885 * subtract the fatest exponential. */
887 int nData,i,status, iter;
889 double *guess, /* Initial guess. */
890 *A, /* The fitted parameters. (A1, B1, A2, B2,... C) */
900 } while (t[nData]<tMax+t[0] && nData<len);
902 p = nexp*2+1; /* Number of parameters. */
905 const gsl_multifit_fdfsolver_type *T
906 = gsl_multifit_fdfsolver_lmsder;
908 gsl_multifit_fdfsolver *s; /* The solver itself. */
909 gsl_multifit_function_fdf fitFunction; /* The function to be fitted. */
910 gsl_matrix *covar; /* Covariance matrix for the parameters.
911 * We'll not use the result, though. */
912 gsl_vector_view theParams;
917 } while (t[nData]<tMax+t[0] && nData<len);
924 covar = gsl_matrix_alloc (p, p);
926 /* Set up an initial gess for the parameters.
927 * The solver is somewhat sensitive to the initial guess,
928 * but this worked fine for a TIP3P box with -geminate dd
929 * EDIT: In fact, this seems like a good starting pont for other watermodels too. */
930 for (i=0; i<nexp; i++)
933 guess[i*2+1] = -0.5 + (((double)i)/nexp - 0.5)*0.3;
935 guess[nexp * 2] = 0.01;
937 theParams = gsl_vector_view_array(guess, p);
942 BD->tDelta = t[1]-t[0];
945 fitFunction.f = &balFunc_f;
946 fitFunction.df = &balFunc_df;
947 fitFunction.fdf = &balFunc_fdf;
948 fitFunction.n = nData;
950 fitFunction.params = BD;
952 s = gsl_multifit_fdfsolver_alloc (T, nData, p);
954 gmx_fatal(FARGS, "Could not set up the nonlinear solver.");
956 gsl_multifit_fdfsolver_set(s, &fitFunction, &theParams.vector);
958 /* \=============================================/ */
964 status = gsl_multifit_fdfsolver_iterate (s);
968 status = gsl_multifit_test_delta (s->dx, s->x, 1e-4, 1e-4);
970 while (iter < 5000 && status == GSL_CONTINUE);
974 fprintf(stderr, "The non-linear fitting did not converge in 5000 steps.\n"
975 "Check the quality of the fit!\n");
979 fprintf(stderr, "Non-linear fitting of ballistic term converged in %d steps.\n\n", (int)iter);
981 for (i=0; i<nexp; i++) {
982 fprintf(stdout, "%c * exp(%c * t) + ", 'A'+(char)i*2, 'B'+(char)i*2);
985 fprintf(stdout, "%c\n", 'A'+(char)nexp*2);
986 fprintf(stdout, "Here are the actual numbers for A-%c:\n", 'A'+nexp*2);
988 for (i=0; i<nexp; i++)
990 A[i*2] = gsl_vector_get(s->x, i*2);
991 A[i*2+1] = gsl_vector_get(s->x, i*2+1);
992 fprintf(stdout, " %g*exp(%g * x) +", A[i*2], A[i*2+1]);
994 A[i*2] = gsl_vector_get(s->x, i*2); /* The last and constant term */
995 fprintf(stdout, " %g\n", A[i*2]);
999 /* Implement some check for parameter quality */
1000 for (i=0; i<nexp; i++)
1002 if (A[i*2]<0 || A[i*2]>1) {
1003 fprintf(stderr, "WARNING: ----------------------------------\n"
1004 " | A coefficient does not lie within [0,1].\n"
1005 " | This may or may not be a problem.\n"
1006 " | Double check the quality of the fit!\n");
1009 fprintf(stderr, "WARNING: ----------------------------------\n"
1010 " | One factor in the exponent is positive.\n"
1011 " | This could be a problem if the coefficient\n"
1012 " | is large. Double check the quality of the fit!\n");
1015 if (A[i*2]<0 || A[i*2]>1) {
1016 fprintf(stderr, "WARNING: ----------------------------------\n"
1017 " | The constant term does not lie within [0,1].\n"
1018 " | This may or may not be a problem.\n"
1019 " | Double check the quality of the fit!\n");
1022 /* Sort the terms */
1023 sorted = (nexp > 1) ? FALSE : TRUE;
1027 for (i=0;i<nexp-1;i++)
1029 ddt[0] = A[i*2] * A[i*2+1];
1030 ddt[1] =A[i*2+2] * A[i*2+3];
1032 if ((bDerivative && (ddt[0]<0 && ddt[1]<0 && ddt[0]>ddt[1])) || /* Compare derivative at t=0... */
1033 (!bDerivative && (A[i*2+1] > A[i*2+3]))) /* Or just the coefficient in the exponent */
1036 a[0] = A[i*2]; /* coefficient */
1037 a[1] = A[i*2+1]; /* parameter in the exponent */
1040 A[i*2+1] = A[i*2+3];
1048 /* Subtract the fastest component */
1049 fprintf(stdout, "Fastest component is %g * exp(%g * t)\n"
1050 "Subtracting fastest component from ACF.\n", A[0], A[1]);
1052 for (i=0; i<len; i++) {
1053 ct[i] = (ct[i] - A[0] * exp(A[1] * i*BD->tDelta)) / (1-A[0]);
1059 gsl_multifit_fdfsolver_free(s);
1060 gsl_matrix_free(covar);
1064 /* We have no gsl. */
1065 fprintf(stderr, "Sorry, can't take away ballistic component without gsl. "
1066 "Recompile using --with-gsl.\n");
1068 #endif /* HAVE_LIBGSL */
1072 extern void dumpN(const real *e, const int nn, char *fn)
1074 /* For debugging only */
1077 char standardName[] = "Nt.xvg";
1085 "@ xaxis label \"Frame\"\n"
1086 "@ yaxis label \"N\"\n"
1087 "@ s0 line type 3\n");
1089 for (i=0; i<nn; i++) {
1090 fprintf(f, "%-10i %-g\n", i, e[i]);
1096 static real* d2r(const double *d, const int nn)
1102 for (i=0; i<nn; i++)
1108 static void _patchBad(double *ct, int n, double dy)
1110 /* Just do lin. interpolation for now. */
1119 static void patchBadPart(double *ct, int n)
1121 _patchBad(ct,n,(ct[n] - ct[0])/n);
1124 static void patchBadTail(double *ct, int n)
1126 _patchBad(ct+1,n-1,ct[1]-ct[0]);
1130 extern void fixGemACF(double *ct, int len)
1135 /* Let's separate two things:
1136 * - identification of bad parts
1137 * - patching of bad parts.
1140 b = 0; /* Start of a bad stretch */
1141 e = 0; /* End of a bad stretch */
1144 /* An acf of binary data must be one at t=0. */
1145 if (abs(ct[0]-1.0) > 1e-6)
1148 fprintf(stderr, "|ct[0]-1.0| = %1.6d. Setting ct[0] to 1.0.\n", abs(ct[0]-1.0));
1151 for (i=0; i<len; i++)
1155 if (isfinite(ct[i]))
1156 #elif defined(HAS__ISFINITE)
1157 if (_isfinite(ct[i]))
1164 /* Still on a good stretch. Proceed.*/
1168 /* Patch up preceding bad stretch. */
1174 gmx_fatal(FARGS, "The ACF is mostly NaN or Inf. Aborting.");
1176 patchBadTail(&(ct[b-2]), (len-b)+1);
1180 patchBadPart(&(ct[b-1]), (e-b)+1);