2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
40 #include "gromacs/legacyheaders/typedefs.h"
41 #include "gromacs/utility/smalloc.h"
42 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/fatalerror.h"
46 #include "gromacs/utility/gmxomp.h"
48 static void missing_code_message()
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");
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().
59 * The parts menetioned in the previous paragraph were contributed under the BSD license.
63 /* This first part is from complex.c which I recieved from Omer Markowitch.
66 * ------------- from complex.c ------------- */
68 /* Complexation of a paired number (r,i) */
69 static gem_complex gem_cmplx(double r, double i)
77 /* Complexation of a real number, x */
78 static gem_complex gem_c(double x)
86 /* Magnitude of a complex number z */
87 static double gem_cx_abs(gem_complex z)
89 return (sqrt(z.r*z.r+z.i*z.i));
92 /* Addition of two complex numbers z1 and z2 */
93 static gem_complex gem_cxadd(gem_complex z1, gem_complex z2)
101 /* Addition of a complex number z1 and a real number r */
102 static gem_complex gem_cxradd(gem_complex z, double r)
110 /* Subtraction of two complex numbers z1 and z2 */
111 static gem_complex gem_cxsub(gem_complex z1, gem_complex z2)
119 /* Multiplication of two complex numbers z1 and z2 */
120 static gem_complex gem_cxmul(gem_complex z1, gem_complex z2)
123 value.r = z1.r*z2.r-z1.i*z2.i;
124 value.i = z1.r*z2.i+z1.i*z2.r;
128 /* Square of a complex number z */
129 static gem_complex gem_cxsq(gem_complex z)
132 value.r = z.r*z.r-z.i*z.i;
133 value.i = z.r*z.i*2.;
137 /* multiplication of a complex number z and a real number r */
138 static gem_complex gem_cxrmul(gem_complex z, double r)
146 /* Division of two complex numbers z1 and z2 */
147 static gem_complex gem_cxdiv(gem_complex z1, gem_complex z2)
151 num = z2.r*z2.r+z2.i*z2.i;
154 fprintf(stderr, "ERROR in gem_cxdiv function\n");
156 value.r = (z1.r*z2.r+z1.i*z2.i)/num; value.i = (z1.i*z2.r-z1.r*z2.i)/num;
160 /* division of a complex z number by a real number x */
161 static gem_complex gem_cxrdiv(gem_complex z, double r)
169 /* division of a real number r by a complex number x */
170 static gem_complex gem_rxcdiv(double r, gem_complex z)
174 f = r/(z.r*z.r+z.i*z.i);
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)
186 value.r = exp_z_r*cos(z.i);
187 value.i = exp_z_r*sin(z.i);
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)
197 mag2 = z.r*z.r+z.i*z.i;
200 fprintf(stderr, "ERROR in gem_cxlog func\n");
202 value.r = log(sqrt(mag2));
213 value.i = atan2(z.i, z.r);
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)
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) */
235 /* Complex power of a complex number z1^z2 */
236 static gem_complex gem_cxdpow(gem_complex z1, gem_complex z2)
239 value = gem_cxdexp(gem_cxmul(gem_cxlog(z1), z2));
243 /* ------------ end of complex.c ------------ */
245 /* This next part was derived from cubic.c, also received from Omer Markovitch.
246 * ------------- from cubic.c ------------- */
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)
252 double t1, t2, two_3, temp;
253 gem_complex ctemp, ct3;
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;
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.));
261 ctemp = gem_rxcdiv(-two_3*t1/3., ct3);
262 ctemp = gem_cxadd(ctemp, gem_cxrdiv(ct3, 3.*two_3));
264 *gam = gem_cxadd(gem_cmplx(a/3., 0.), ctemp);
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);
273 /* ------------ end of cubic.c ------------ */
275 /* This next part was derived from cerror.c and rerror.c, also received from Omer Markovitch.
276 * ------------- from [cr]error.c ------------- */
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 /************************************************************/
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 */
295 static double gem_erf(double x)
297 double n, sum, temp, exp2, xm, x2, x4, x6, x8, x10, x12;
310 for (n = 1.; n <= 2000.; n += 1.)
312 temp *= 2.*x2/(2.*n+1.);
314 if (fabs(temp/sum) < 1.E-16)
322 fprintf(stderr, "In Erf calc - iteration exceeds %lg\n", n);
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)
333 sum *= exp2; /* now sum is erfc(x) */
336 return x >= 0.0 ? sum : -sum;
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)
347 ans = t * exp(-z*z - 1.26551223 +
356 t*0.17087277)))))))));
358 return x >= 0.0 ? ans : 2.0-ans;
361 /* omega(x)=exp(x^2)erfc(x) */
362 static double gem_omega(double x)
364 double xm, ans, xx, x4, x6, x8, x10, x12;
375 ans = exp(xx)*gem_erfc(x);
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;
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)
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;
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);
410 for (n = 1.0, temp = 0.; n <= 2000.; n += 1.0)
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)
428 fprintf(stderr, "iteration exceeds %lg\n", n);
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));
447 /* ------------ end of [cr]error.c ------------ */
449 /*_ REVERSIBLE GEMINATE RECOMBINATION
451 * Here are the functions for reversible geminate recombination. */
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)
458 fprintf(stdout, "Diffusion coefficient is %f A^2/ps\n", D*1e4);
459 return (double)(D*1e4);
463 static double eq10v2(double theoryCt[], double time[], int manytimes,
464 double ka, double kd, t_gemParams *params)
466 /* Finding the 3 roots */
468 double kD, D, r, a, b, c, tsqrt, sumimaginary;
473 part1, part2, part3, part4;
478 a = (1.0 + ka/kD) * sqrt(D)/r;
482 gem_solve(&alpha, &beta, &gamma, a, b, c);
483 /* Finding the 3 roots */
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) */
491 #pragma omp parallel for \
492 private(i, tsqrt, oma, omb, omc, c1, c2, c3, c4) \
493 reduction(+:sumimaginary) \
496 for (i = 0; i < manytimes; i++)
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;
508 sumimaginary += c4.i * c4.i;
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)
518 return (exp(((double)(i)) * params->logQuota) -1);
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)
531 /* A few hardcoded things here. For now anyway. */
533 /* p->ka_max = 100; */
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;
548 /* Parameters used by calcsquare(). Better to calculate them
549 * here than in calcsquare every time it's called. */
551 /* p->logAfterTime = logAfterTime; */
552 tDelta = (t[len-1]-t[0]) / len;
555 gmx_fatal(FARGS, "Time between frames is non-positive!");
562 p->nExpFit = nBalExp;
563 /* p->nLin = logAfterTime / tDelta; */
564 p->nFitPoints = nFitPoints;
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"); */
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 */
578 /* p->logMult = pow((float)len, 1.0/nLin);/\* pow(t[len-1]-t[0], 1.0/p->nLin); *\/ */
579 p->ballistic = ballistic;
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. */
589 static real* d2r(const double *d, const int nn);
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)
598 int nThreads, i, iter, status, maxiter;
599 real size, d2, tol, *dumpdata;
602 char *dumpstr, dumpname[128];
604 missing_code_message();
610 /* Removes the ballistic term from the beginning of the ACF,
611 * just like in Omer's paper.
613 extern void takeAwayBallistic(double gmx_unused *ct, double *t, int len, real tMax, int nexp, gmx_bool gmx_unused bDerivative)
616 /* Fit with 4 exponentials and one constant term,
617 * subtract the fatest exponential. */
619 int nData, i, status, iter;
621 double *guess, /* Initial guess. */
622 *A, /* The fitted parameters. (A1, B1, A2, B2,... C) */
634 while (t[nData] < tMax+t[0] && nData < len);
636 p = nexp*2+1; /* Number of parameters. */
638 missing_code_message();
642 extern void dumpN(const real *e, const int nn, char *fn)
644 /* For debugging only */
647 char standardName[] = "Nt.xvg";
656 "@ xaxis label \"Frame\"\n"
657 "@ yaxis label \"N\"\n"
658 "@ s0 line type 3\n");
660 for (i = 0; i < nn; i++)
662 fprintf(f, "%-10i %-g\n", i, e[i]);
668 static real* d2r(const double *d, const int nn)
674 for (i = 0; i < nn; i++)
682 static void _patchBad(double *ct, int n, double dy)
684 /* Just do lin. interpolation for now. */
687 for (i = 1; i < n; i++)
693 static void patchBadPart(double *ct, int n)
695 _patchBad(ct, n, (ct[n] - ct[0])/n);
698 static void patchBadTail(double *ct, int n)
700 _patchBad(ct+1, n-1, ct[1]-ct[0]);
704 extern void fixGemACF(double *ct, int len)
709 /* Let's separate two things:
710 * - identification of bad parts
711 * - patching of bad parts.
714 b = 0; /* Start of a bad stretch */
715 e = 0; /* End of a bad stretch */
718 /* An acf of binary data must be one at t=0. */
719 if (fabs(ct[0]-1.0) > 1e-6)
722 fprintf(stderr, "|ct[0]-1.0| = %1.6f. Setting ct[0] to 1.0.\n", fabs(ct[0]-1.0));
725 for (i = 0; i < len; i++)
730 #elif defined(HAS__ISFINITE)
731 if (_isfinite(ct[i]))
738 /* Still on a good stretch. Proceed.*/
742 /* Patch up preceding bad stretch. */
748 gmx_fatal(FARGS, "The ACF is mostly NaN or Inf. Aborting.");
750 patchBadTail(&(ct[b-2]), (len-b)+1);
754 patchBadPart(&(ct[b-1]), (e-b)+1);