added os-specific defines from cmake required by memtestG80
[gromacs/qmmm-gamess-us.git] / src / mdlib / tables.c
blob2f32c1a57dc429c3fa1aa739e1f57c67b99ab02a
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
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
32 * And Hey:
33 * GROwing Monsters And Cloning Shrimps
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
39 #include <math.h>
40 #include "maths.h"
41 #include "typedefs.h"
42 #include "names.h"
43 #include "smalloc.h"
44 #include "gmx_fatal.h"
45 #include "futil.h"
46 #include "xvgr.h"
47 #include "vec.h"
48 #include "main.h"
49 #include "network.h"
50 #include "physics.h"
51 #include "force.h"
52 #include "gmxfio.h"
54 /* All the possible (implemented) table functions */
55 enum {
56 etabLJ6,
57 etabLJ12,
58 etabLJ6Shift,
59 etabLJ12Shift,
60 etabShift,
61 etabRF,
62 etabRF_ZERO,
63 etabCOUL,
64 etabEwald,
65 etabEwaldSwitch,
66 etabEwaldUser,
67 etabEwaldUserSwitch,
68 etabLJ6Switch,
69 etabLJ12Switch,
70 etabCOULSwitch,
71 etabLJ6Encad,
72 etabLJ12Encad,
73 etabCOULEncad,
74 etabEXPMIN,
75 etabUSER,
76 etabNR
79 typedef struct {
80 const char *name;
81 bool bCoulomb;
82 } t_tab_props;
84 /* This structure holds name and a flag that tells whether
85 this is a Coulomb type funtion */
86 static const t_tab_props tprops[etabNR] = {
87 { "LJ6", FALSE },
88 { "LJ12", FALSE },
89 { "LJ6Shift", FALSE },
90 { "LJ12Shift", FALSE },
91 { "Shift", TRUE },
92 { "RF", TRUE },
93 { "RF-zero", TRUE },
94 { "COUL", TRUE },
95 { "Ewald", TRUE },
96 { "Ewald-Switch", TRUE },
97 { "Ewald-User", TRUE },
98 { "Ewald-User-Switch", TRUE },
99 { "LJ6Switch", FALSE },
100 { "LJ12Switch", FALSE },
101 { "COULSwitch", TRUE },
102 { "LJ6-Encad shift", FALSE },
103 { "LJ12-Encad shift", FALSE },
104 { "COUL-Encad shift", TRUE },
105 { "EXPMIN", FALSE },
106 { "USER", FALSE }
109 /* Index in the table that says which function to use */
110 enum { etiCOUL, etiLJ6, etiLJ12, etiNR };
112 typedef struct {
113 int nx,nx0;
114 double tabscale;
115 double *x,*v,*f;
116 } t_tabledata;
118 #define pow2(x) ((x)*(x))
119 #define pow3(x) ((x)*(x)*(x))
120 #define pow4(x) ((x)*(x)*(x)*(x))
121 #define pow5(x) ((x)*(x)*(x)*(x)*(x))
123 /* Calculate the potential and force for an r value
124 * in exactly the same way it is done in the inner loop.
125 * VFtab is a pointer to the table data, offset is
126 * the point where we should begin and stride is
127 * 4 if we have a buckingham table, 3 otherwise.
128 * If you want to evaluate table no N, set offset to 4*N.
130 * We use normal precision here, since that is what we
131 * will use in the inner loops.
133 static void evaluate_table(real VFtab[], int offset, int stride,
134 real tabscale, real r, real *y, real *yp)
136 int n;
137 real rt,eps,eps2;
138 real Y,F,Geps,Heps2,Fp;
140 rt = r*tabscale;
141 n = (int)rt;
142 eps = rt - n;
143 eps2 = eps*eps;
144 n = offset+stride*n;
145 Y = VFtab[n];
146 F = VFtab[n+1];
147 Geps = eps*VFtab[n+2];
148 Heps2 = eps2*VFtab[n+3];
149 Fp = F+Geps+Heps2;
150 *y = Y+eps*Fp;
151 *yp = (Fp+Geps+2.0*Heps2)*tabscale;
155 static void splint(real xa[],real ya[],real y2a[],
156 int n,real x,real *y,real *yp)
158 int klo,khi,k;
159 real h,b,a,eps;
160 real F,G,H;
162 klo=1;
163 khi=n;
165 while ((khi-klo) > 1) {
166 k=(khi+klo) >> 1;
167 if (xa[k] > x)
168 khi=k;
169 else
170 klo=k;
172 h = xa[khi]-xa[klo];
173 if (h == 0.0)
174 gmx_fatal(FARGS,"Bad XA input to function splint");
175 a = (xa[khi]-x)/h;
176 b = (x-xa[klo])/h;
177 *y = a*ya[klo]+b*ya[khi]+((a*a*a-a)*y2a[klo]+(b*b*b-b)*y2a[khi])*(h*h)/6.0;
178 *yp = (ya[khi]-ya[klo])/h+((3*a*a-1)*y2a[klo]-(3*b*b-1)*y2a[khi])*h/6.0;
180 eps = b;
181 F = (ya[khi]-ya[klo]-(h*h/6.0)*(2*y2a[klo]+y2a[khi]));
182 G = (h*h/2.0)*y2a[klo];
183 H = (h*h/6.0)*(y2a[khi]-y2a[klo]);
184 *y = ya[klo] + eps*F + eps*eps*G + eps*eps*eps*H;
185 *yp = (F + 2*eps*G + 3*eps*eps*H)/h;
189 static void copy2table(int n,int offset,int stride,
190 double x[],double Vtab[],double Ftab[],
191 real dest[])
193 /* Use double prec. for the intermediary variables
194 * and temporary x/vtab/vtab2 data to avoid unnecessary
195 * loss of precision.
197 int i,nn0;
198 double F,G,H,h;
200 h = 0;
201 for(i=0; (i<n); i++) {
202 if (i < n-1) {
203 h = x[i+1] - x[i];
204 F = -Ftab[i]*h;
205 G = 3*(Vtab[i+1] - Vtab[i]) + (Ftab[i+1] + 2*Ftab[i])*h;
206 H = -2*(Vtab[i+1] - Vtab[i]) - (Ftab[i+1] + Ftab[i])*h;
207 } else {
208 /* Fill the last entry with a linear potential,
209 * this is mainly for rounding issues with angle and dihedral potentials.
211 F = -Ftab[i]*h;
212 G = 0;
213 H = 0;
215 nn0 = offset + i*stride;
216 dest[nn0] = Vtab[i];
217 dest[nn0+1] = F;
218 dest[nn0+2] = G;
219 dest[nn0+3] = H;
223 static void init_table(FILE *fp,int n,int nx0,
224 double tabscale,t_tabledata *td,bool bAlloc)
226 int i;
228 td->nx = n;
229 td->nx0 = nx0;
230 td->tabscale = tabscale;
231 if (bAlloc) {
232 snew(td->x,td->nx);
233 snew(td->v,td->nx);
234 snew(td->f,td->nx);
236 for(i=td->nx0; (i<td->nx); i++)
237 td->x[i] = i/tabscale;
240 static void spline_forces(int nx,double h,double v[],bool bS3,bool bE3,
241 double f[])
243 int start,end,i;
244 double v3,b_s,b_e,b;
245 double beta,*gamma;
247 /* Formulas can be found in:
248 * H.J.C. Berendsen, Simulating the Physical World, Cambridge 2007
251 if (nx < 4 && (bS3 || bE3))
252 gmx_fatal(FARGS,"Can not generate splines with third derivative boundary conditions with less than 4 (%d) points",nx);
254 /* To make life easy we initially set the spacing to 1
255 * and correct for this at the end.
257 beta = 2;
258 if (bS3) {
259 /* Fit V''' at the start */
260 v3 = v[3] - 3*v[2] + 3*v[1] - v[0];
261 if (debug)
262 fprintf(debug,"The left third derivative is %g\n",v3/(h*h*h));
263 b_s = 2*(v[1] - v[0]) + v3/6;
264 start = 0;
266 if (FALSE) {
267 /* Fit V'' at the start */
268 real v2;
270 v2 = -v[3] + 4*v[2] - 5*v[1] + 2*v[0];
271 /* v2 = v[2] - 2*v[1] + v[0]; */
272 if (debug)
273 fprintf(debug,"The left second derivative is %g\n",v2/(h*h));
274 b_s = 3*(v[1] - v[0]) - v2/2;
275 start = 0;
277 } else {
278 b_s = 3*(v[2] - v[0]) + f[0]*h;
279 start = 1;
281 if (bE3) {
282 /* Fit V''' at the end */
283 v3 = v[nx-1] - 3*v[nx-2] + 3*v[nx-3] - v[nx-4];
284 if (debug)
285 fprintf(debug,"The right third derivative is %g\n",v3/(h*h*h));
286 b_e = 2*(v[nx-1] - v[nx-2]) + v3/6;
287 end = nx;
288 } else {
289 /* V'=0 at the end */
290 b_e = 3*(v[nx-1] - v[nx-3]) + f[nx-1]*h;
291 end = nx - 1;
294 snew(gamma,nx);
295 beta = (bS3 ? 1 : 4);
297 /* For V'' fitting */
298 /* beta = (bS3 ? 2 : 4); */
300 f[start] = b_s/beta;
301 for(i=start+1; i<end; i++) {
302 gamma[i] = 1/beta;
303 beta = 4 - gamma[i];
304 b = 3*(v[i+1] - v[i-1]);
305 f[i] = (b - f[i-1])/beta;
307 gamma[end-1] = 1/beta;
308 beta = (bE3 ? 1 : 4) - gamma[end-1];
309 f[end-1] = (b_e - f[end-2])/beta;
311 for(i=end-2; i>=start; i--)
312 f[i] -= gamma[i+1]*f[i+1];
313 sfree(gamma);
315 /* Correct for the minus sign and the spacing */
316 for(i=start; i<end; i++)
317 f[i] = -f[i]/h;
320 static void set_forces(FILE *fp,int angle,
321 int nx,double h,double v[],double f[],
322 int table)
324 int start,end;
326 if (angle == 2)
327 gmx_fatal(FARGS,
328 "Force generation for dihedral tables is not (yet) implemented");
330 start = 0;
331 while (v[start] == 0)
332 start++;
334 end = nx;
335 while(v[end-1] == 0)
336 end--;
337 if (end > nx - 2)
338 end = nx;
339 else
340 end++;
342 if (fp)
343 fprintf(fp,"Generating forces for table %d, boundary conditions: V''' at %g, %s at %g\n",
344 table+1,start*h,end==nx ? "V'''" : "V'=0",(end-1)*h);
345 spline_forces(end-start,h,v+start,TRUE,end==nx,f+start);
348 static void read_tables(FILE *fp,const char *fn,
349 int ntab,int angle,t_tabledata td[])
351 char *libfn;
352 char buf[STRLEN];
353 double **yy=NULL,start,end,dx0,dx1,ssd,vm,vp,f,numf;
354 int k,i,nx,nx0=0,ny,nny,ns;
355 bool bAllZero,bZeroV,bZeroF;
356 double tabscale;
358 nny = 2*ntab+1;
359 libfn = gmxlibfn(fn);
360 nx = read_xvg(libfn,&yy,&ny);
361 if (ny != nny)
362 gmx_fatal(FARGS,"Trying to read file %s, but nr columns = %d, should be %d",
363 libfn,ny,nny);
364 if (angle == 0) {
365 if (yy[0][0] != 0.0)
366 gmx_fatal(FARGS,
367 "The first distance in file %s is %f nm instead of %f nm",
368 libfn,yy[0][0],0.0);
369 } else {
370 if (angle == 1)
371 start = 0.0;
372 else
373 start = -180.0;
374 end = 180.0;
375 if (yy[0][0] != start || yy[0][nx-1] != end)
376 gmx_fatal(FARGS,"The angles in file %s should go from %f to %f instead of %f to %f\n",
377 libfn,start,end,yy[0][0],yy[0][nx-1]);
380 tabscale = (nx-1)/(yy[0][nx-1] - yy[0][0]);
382 if (fp) {
383 fprintf(fp,"Read user tables from %s with %d data points.\n",libfn,nx);
384 if (angle == 0)
385 fprintf(fp,"Tabscale = %g points/nm\n",tabscale);
388 bAllZero = TRUE;
389 for(k=0; k<ntab; k++) {
390 bZeroV = TRUE;
391 bZeroF = TRUE;
392 for(i=0; (i < nx); i++) {
393 if (i >= 2) {
394 dx0 = yy[0][i-1] - yy[0][i-2];
395 dx1 = yy[0][i] - yy[0][i-1];
396 /* Check for 1% deviation in spacing */
397 if (fabs(dx1 - dx0) >= 0.005*(fabs(dx0) + fabs(dx1))) {
398 gmx_fatal(FARGS,"In table file '%s' the x values are not equally spaced: %f %f %f",fn,yy[0][i-2],yy[0][i-1],yy[0][i]);
401 if (yy[1+k*2][i] != 0) {
402 bZeroV = FALSE;
403 if (bAllZero) {
404 bAllZero = FALSE;
405 nx0 = i;
407 if (yy[1+k*2][i] > 0.01*GMX_REAL_MAX ||
408 yy[1+k*2][i] < -0.01*GMX_REAL_MAX) {
409 gmx_fatal(FARGS,"Out of range potential value %g in file '%s'",
410 yy[1+k*2][i],fn);
413 if (yy[1+k*2+1][i] != 0) {
414 bZeroF = FALSE;
415 if (bAllZero) {
416 bAllZero = FALSE;
417 nx0 = i;
419 if (yy[1+k*2+1][i] > 0.01*GMX_REAL_MAX ||
420 yy[1+k*2+1][i] < -0.01*GMX_REAL_MAX) {
421 gmx_fatal(FARGS,"Out of range force value %g in file '%s'",
422 yy[1+k*2+1][i],fn);
427 if (!bZeroV && bZeroF) {
428 set_forces(fp,angle,nx,1/tabscale,yy[1+k*2],yy[1+k*2+1],k);
429 } else {
430 /* Check if the second column is close to minus the numerical
431 * derivative of the first column.
433 ssd = 0;
434 ns = 0;
435 for(i=1; (i < nx-1); i++) {
436 vm = yy[1+2*k][i-1];
437 vp = yy[1+2*k][i+1];
438 f = yy[1+2*k+1][i];
439 if (vm != 0 && vp != 0 && f != 0) {
440 /* Take the centered difference */
441 numf = -(vp - vm)*0.5*tabscale;
442 ssd += fabs(2*(f - numf)/(f + numf));
443 ns++;
446 if (ns > 0) {
447 ssd /= ns;
448 sprintf(buf,"For the %d non-zero entries for table %d in %s the forces deviate on average %d%% from minus the numerical derivative of the potential\n",ns,k,libfn,(int)(100*ssd+0.5));
449 if (debug)
450 fprintf(debug,"%s",buf);
451 if (ssd > 0.2) {
452 if (fp)
453 fprintf(fp,"\nWARNING: %s\n",buf);
454 fprintf(stderr,"\nWARNING: %s\n",buf);
459 if (bAllZero && fp) {
460 fprintf(fp,"\nNOTE: All elements in table %s are zero\n\n",libfn);
463 for(k=0; (k<ntab); k++) {
464 init_table(fp,nx,nx0,tabscale,&(td[k]),TRUE);
465 for(i=0; (i<nx); i++) {
466 td[k].x[i] = yy[0][i];
467 td[k].v[i] = yy[2*k+1][i];
468 td[k].f[i] = yy[2*k+2][i];
471 for(i=0; (i<ny); i++)
472 sfree(yy[i]);
473 sfree(yy);
474 sfree(libfn);
477 static void done_tabledata(t_tabledata *td)
479 int i;
481 if (!td)
482 return;
484 sfree(td->x);
485 sfree(td->v);
486 sfree(td->f);
489 static void fill_table(t_tabledata *td,int tp,const t_forcerec *fr)
491 /* Fill the table according to the formulas in the manual.
492 * In principle, we only need the potential and the second
493 * derivative, but then we would have to do lots of calculations
494 * in the inner loop. By precalculating some terms (see manual)
495 * we get better eventual performance, despite a larger table.
497 * Since some of these higher-order terms are very small,
498 * we always use double precision to calculate them here, in order
499 * to avoid unnecessary loss of precision.
501 #ifdef DEBUG_SWITCH
502 FILE *fp;
503 #endif
504 int i;
505 double reppow,p;
506 double r1,rc,r12,r13;
507 double r,r2,r6,rc6;
508 double expr,Vtab,Ftab;
509 /* Parameters for David's function */
510 double A=0,B=0,C=0,A_3=0,B_4=0;
511 /* Parameters for the switching function */
512 double ksw,swi,swi1;
513 /* Temporary parameters */
514 bool bSwitch,bShift;
515 double ewc=fr->ewaldcoeff;
516 double isp= 0.564189583547756;
518 bSwitch = ((tp == etabLJ6Switch) || (tp == etabLJ12Switch) ||
519 (tp == etabCOULSwitch) ||
520 (tp == etabEwaldSwitch) || (tp == etabEwaldUserSwitch));
521 bShift = ((tp == etabLJ6Shift) || (tp == etabLJ12Shift) ||
522 (tp == etabShift));
524 reppow = fr->reppow;
526 if (tprops[tp].bCoulomb) {
527 r1 = fr->rcoulomb_switch;
528 rc = fr->rcoulomb;
530 else {
531 r1 = fr->rvdw_switch;
532 rc = fr->rvdw;
534 if (bSwitch)
535 ksw = 1.0/(pow5(rc-r1));
536 else
537 ksw = 0.0;
538 if (bShift) {
539 if (tp == etabShift)
540 p = 1;
541 else if (tp == etabLJ6Shift)
542 p = 6;
543 else
544 p = reppow;
546 A = p * ((p+1)*r1-(p+4)*rc)/(pow(rc,p+2)*pow2(rc-r1));
547 B = -p * ((p+1)*r1-(p+3)*rc)/(pow(rc,p+2)*pow3(rc-r1));
548 C = 1.0/pow(rc,p)-A/3.0*pow3(rc-r1)-B/4.0*pow4(rc-r1);
549 if (tp == etabLJ6Shift) {
550 A=-A;
551 B=-B;
552 C=-C;
554 A_3=A/3.0;
555 B_4=B/4.0;
557 if (debug) { fprintf(debug,"Setting up tables\n"); fflush(debug); }
559 #ifdef DEBUG_SWITCH
560 fp=xvgropen("switch.xvg","switch","r","s");
561 #endif
563 for(i=td->nx0; (i<td->nx); i++) {
564 r = td->x[i];
565 r2 = r*r;
566 r6 = 1.0/(r2*r2*r2);
567 r12 = r6*r6;
568 Vtab = 0.0;
569 Ftab = 0.0;
570 if (bSwitch) {
571 /* swi is function, swi1 1st derivative and swi2 2nd derivative */
572 /* The switch function is 1 for r<r1, 0 for r>rc, and smooth for
573 * r1<=r<=rc. The 1st and 2nd derivatives are both zero at
574 * r1 and rc.
575 * ksw is just the constant 1/(rc-r1)^5, to save some calculations...
577 if(r<=r1) {
578 swi = 1.0;
579 swi1 = 0.0;
580 } else if (r>=rc) {
581 swi = 0.0;
582 swi1 = 0.0;
583 } else {
584 swi = 1 - 10*pow3(r-r1)*ksw*pow2(rc-r1)
585 + 15*pow4(r-r1)*ksw*(rc-r1) - 6*pow5(r-r1)*ksw;
586 swi1 = -30*pow2(r-r1)*ksw*pow2(rc-r1)
587 + 60*pow3(r-r1)*ksw*(rc-r1) - 30*pow4(r-r1)*ksw;
590 else { /* not really needed, but avoids compiler warnings... */
591 swi = 1.0;
592 swi1 = 0.0;
594 #ifdef DEBUG_SWITCH
595 fprintf(fp,"%10g %10g %10g %10g\n",r,swi,swi1,swi2);
596 #endif
598 rc6 = rc*rc*rc;
599 rc6 = 1.0/(rc6*rc6);
601 switch (tp) {
602 case etabLJ6:
603 /* Dispersion */
604 Vtab = -r6;
605 Ftab = 6.0*Vtab/r;
606 break;
607 case etabLJ6Switch:
608 case etabLJ6Shift:
609 /* Dispersion */
610 if (r < rc) {
611 Vtab = -r6;
612 Ftab = 6.0*Vtab/r;
614 break;
615 case etabLJ12:
616 /* Repulsion */
617 Vtab = r12;
618 Ftab = reppow*Vtab/r;
619 break;
620 case etabLJ12Switch:
621 case etabLJ12Shift:
622 /* Repulsion */
623 if (r < rc) {
624 Vtab = r12;
625 Ftab = reppow*Vtab/r;
627 break;
628 case etabLJ6Encad:
629 if(r < rc) {
630 Vtab = -(r6-6.0*(rc-r)*rc6/rc-rc6);
631 Ftab = -(6.0*r6/r-6.0*rc6/rc);
632 } else { /* r>rc */
633 Vtab = 0;
634 Ftab = 0;
636 break;
637 case etabLJ12Encad:
638 if(r < rc) {
639 Vtab = r12-12.0*(rc-r)*rc6*rc6/rc-1.0*rc6*rc6;
640 Ftab = 12.0*r12/r-12.0*rc6*rc6/rc;
641 } else { /* r>rc */
642 Vtab = 0;
643 Ftab = 0;
645 break;
646 case etabCOUL:
647 Vtab = 1.0/r;
648 Ftab = 1.0/r2;
649 break;
650 case etabCOULSwitch:
651 case etabShift:
652 if (r < rc) {
653 Vtab = 1.0/r;
654 Ftab = 1.0/r2;
656 break;
657 case etabEwald:
658 case etabEwaldSwitch:
659 Vtab = gmx_erfc(ewc*r)/r;
660 Ftab = gmx_erfc(ewc*r)/r2+2*exp(-(ewc*ewc*r2))*ewc*isp/r;
661 break;
662 case etabEwaldUser:
663 case etabEwaldUserSwitch:
664 /* Only calculate minus the reciprocal space contribution */
665 Vtab = -gmx_erf(ewc*r)/r;
666 Ftab = -gmx_erf(ewc*r)/r2+2*exp(-(ewc*ewc*r2))*ewc*isp/r;
667 break;
668 case etabRF:
669 case etabRF_ZERO:
670 Vtab = 1.0/r + fr->k_rf*r2 - fr->c_rf;
671 Ftab = 1.0/r2 - 2*fr->k_rf*r;
672 if (tp == etabRF_ZERO && r >= rc) {
673 Vtab = 0;
674 Ftab = 0;
676 break;
677 case etabEXPMIN:
678 expr = exp(-r);
679 Vtab = expr;
680 Ftab = expr;
681 break;
682 case etabCOULEncad:
683 if(r < rc) {
684 Vtab = 1.0/r-(rc-r)/(rc*rc)-1.0/rc;
685 Ftab = 1.0/r2-1.0/(rc*rc);
686 } else { /* r>rc */
687 Vtab = 0;
688 Ftab = 0;
690 break;
691 default:
692 gmx_fatal(FARGS,"Table type %d not implemented yet. (%s,%d)",
693 tp,__FILE__,__LINE__);
695 if (bShift) {
696 /* Normal coulomb with cut-off correction for potential */
697 if (r < rc) {
698 Vtab -= C;
699 /* If in Shifting range add something to it */
700 if (r > r1) {
701 r12 = (r-r1)*(r-r1);
702 r13 = (r-r1)*r12;
703 Vtab += - A_3*r13 - B_4*r12*r12;
704 Ftab += A*r12 + B*r13;
709 if (tp == etabEwaldUser) {
710 Vtab += td->v[i];
711 Ftab += td->f[i];
714 if ((r > r1) && bSwitch) {
715 Ftab = Ftab*swi - Vtab*swi1;
716 Vtab = Vtab*swi;
719 /* Convert to single precision when we store to mem */
720 td->v[i] = Vtab;
721 td->f[i] = Ftab;
724 #ifdef DEBUG_SWITCH
725 gmx_fio_fclose(fp);
726 #endif
729 static void set_table_type(int tabsel[],const t_forcerec *fr,bool b14only)
731 int eltype,vdwtype;
733 /* Set the different table indices.
734 * Coulomb first.
738 if (b14only) {
739 switch (fr->eeltype) {
740 case eelRF_NEC:
741 eltype = eelRF;
742 break;
743 case eelUSER:
744 case eelPMEUSER:
745 eltype = eelUSER;
746 break;
747 default:
748 eltype = eelCUT;
750 } else {
751 eltype = fr->eeltype;
754 switch (eltype) {
755 case eelCUT:
756 tabsel[etiCOUL] = etabCOUL;
757 break;
758 case eelPPPM:
759 case eelPOISSON:
760 tabsel[etiCOUL] = etabShift;
761 break;
762 case eelSHIFT:
763 if (fr->rcoulomb > fr->rcoulomb_switch)
764 tabsel[etiCOUL] = etabShift;
765 else
766 tabsel[etiCOUL] = etabCOUL;
767 break;
768 case eelEWALD:
769 case eelPME:
770 tabsel[etiCOUL] = etabEwald;
771 break;
772 case eelPMESWITCH:
773 tabsel[etiCOUL] = etabEwaldSwitch;
774 break;
775 case eelPMEUSER:
776 tabsel[etiCOUL] = etabEwaldUser;
777 break;
778 case eelRF:
779 case eelGRF:
780 case eelRF_NEC:
781 tabsel[etiCOUL] = etabRF;
782 break;
783 case eelRF_ZERO:
784 tabsel[etiCOUL] = etabRF_ZERO;
785 break;
786 case eelSWITCH:
787 tabsel[etiCOUL] = etabCOULSwitch;
788 break;
789 case eelUSER:
790 tabsel[etiCOUL] = etabUSER;
791 break;
792 case eelENCADSHIFT:
793 tabsel[etiCOUL] = etabCOULEncad;
794 break;
795 default:
796 gmx_fatal(FARGS,"Invalid eeltype %d",eltype);
799 /* Van der Waals time */
800 if (fr->bBHAM && !b14only) {
801 tabsel[etiLJ6] = etabLJ6;
802 tabsel[etiLJ12] = etabEXPMIN;
803 } else {
804 if (b14only && fr->vdwtype != evdwUSER)
805 vdwtype = evdwCUT;
806 else
807 vdwtype = fr->vdwtype;
809 switch (vdwtype) {
810 case evdwSWITCH:
811 tabsel[etiLJ6] = etabLJ6Switch;
812 tabsel[etiLJ12] = etabLJ12Switch;
813 break;
814 case evdwSHIFT:
815 tabsel[etiLJ6] = etabLJ6Shift;
816 tabsel[etiLJ12] = etabLJ12Shift;
817 break;
818 case evdwUSER:
819 tabsel[etiLJ6] = etabUSER;
820 tabsel[etiLJ12] = etabUSER;
821 break;
822 case evdwCUT:
823 tabsel[etiLJ6] = etabLJ6;
824 tabsel[etiLJ12] = etabLJ12;
825 break;
826 case evdwENCADSHIFT:
827 tabsel[etiLJ6] = etabLJ6Encad;
828 tabsel[etiLJ12] = etabLJ12Encad;
829 break;
830 default:
831 gmx_fatal(FARGS,"Invalid vdwtype %d in %s line %d",vdwtype,
832 __FILE__,__LINE__);
837 t_forcetable make_tables(FILE *out,const output_env_t oenv,
838 const t_forcerec *fr,
839 bool bVerbose,const char *fn,
840 real rtab,int flags)
842 const char *fns[3] = { "ctab.xvg", "dtab.xvg", "rtab.xvg" };
843 const char *fns14[3] = { "ctab14.xvg", "dtab14.xvg", "rtab14.xvg" };
844 FILE *fp;
845 t_tabledata *td;
846 bool b14only,bReadTab,bGenTab;
847 real x0,y0,yp;
848 int i,j,k,nx,nx0,tabsel[etiNR];
849 void * p_tmp;
851 t_forcetable table;
853 b14only = (flags & GMX_MAKETABLES_14ONLY);
855 if (flags & GMX_MAKETABLES_FORCEUSER) {
856 tabsel[etiCOUL] = etabUSER;
857 tabsel[etiLJ6] = etabUSER;
858 tabsel[etiLJ12] = etabUSER;
859 } else {
860 set_table_type(tabsel,fr,b14only);
862 snew(td,etiNR);
863 table.r = rtab;
864 table.scale = 0;
865 table.n = 0;
866 table.scale_exp = 0;
867 nx0 = 10;
868 nx = 0;
870 /* Check whether we have to read or generate */
871 bReadTab = FALSE;
872 bGenTab = FALSE;
873 for(i=0; (i<etiNR); i++) {
874 if (tabsel[i] == etabUSER || tabsel[i] == etabEwaldUser)
875 bReadTab = TRUE;
876 if (tabsel[i] != etabUSER)
877 bGenTab = TRUE;
879 if (bReadTab) {
880 read_tables(out,fn,etiNR,0,td);
881 if (rtab == 0 || (flags & GMX_MAKETABLES_14ONLY)) {
882 rtab = td[0].x[td[0].nx-1];
883 table.n = td[0].nx;
884 nx = table.n;
885 } else {
886 if (td[0].x[td[0].nx-1] < rtab)
887 gmx_fatal(FARGS,"Tables in file %s not long enough for cut-off:\n"
888 "\tshould be at least %f nm\n",fn,rtab);
889 nx = table.n = (int)(rtab*td[0].tabscale + 0.5);
891 table.scale = td[0].tabscale;
892 nx0 = td[0].nx0;
894 if (bGenTab) {
895 if (!bReadTab) {
896 #ifdef GMX_DOUBLE
897 table.scale = 2000.0;
898 #else
899 table.scale = 500.0;
900 #endif
901 nx = table.n = rtab*table.scale;
904 if (fr->bBHAM) {
905 if(fr->bham_b_max!=0)
906 table.scale_exp = table.scale/fr->bham_b_max;
907 else
908 table.scale_exp = table.scale;
911 /* Each table type (e.g. coul,lj6,lj12) requires four
912 * numbers per datapoint. For performance reasons we want
913 * the table data to be aligned to 16-byte. This is accomplished
914 * by allocating 16 bytes extra to a temporary pointer, and then
915 * calculating an aligned pointer. This new pointer must not be
916 * used in a free() call, but thankfully we're sloppy enough not
917 * to do this :-)
920 /* 12 fp entries per table point, nx+1 points, and 16 bytes extra to align it. */
921 p_tmp = malloc(12*(nx+1)*sizeof(real)+16);
923 /* align it - size_t has the same same as a pointer */
924 table.tab = (real *) (((size_t) p_tmp + 16) & (~((size_t) 15)));
927 for(k=0; (k<etiNR); k++) {
928 if (tabsel[k] != etabUSER) {
929 init_table(out,nx,nx0,
930 (tabsel[k] == etabEXPMIN) ? table.scale_exp : table.scale,
931 &(td[k]),!bReadTab);
932 fill_table(&(td[k]),tabsel[k],fr);
933 if (out)
934 fprintf(out,"%s table with %d data points for %s%s.\n"
935 "Tabscale = %g points/nm\n",
936 tabsel[k]==etabEwaldUser ? "Modified" : "Generated",
937 td[k].nx,b14only?"1-4 ":"",tprops[tabsel[k]].name,
938 td[k].tabscale);
940 copy2table(table.n,k*4,12,td[k].x,td[k].v,td[k].f,table.tab);
942 if (bDebugMode() && bVerbose) {
943 if (b14only)
944 fp=xvgropen(fns14[k],fns14[k],"r","V",oenv);
945 else
946 fp=xvgropen(fns[k],fns[k],"r","V",oenv);
947 /* plot the output 5 times denser than the table data */
948 for(i=5*((nx0+1)/2); i<5*table.n; i++) {
949 x0 = i*table.r/(5*(table.n-1));
950 evaluate_table(table.tab,4*k,12,table.scale,x0,&y0,&yp);
951 fprintf(fp,"%15.10e %15.10e %15.10e\n",x0,y0,yp);
953 gmx_fio_fclose(fp);
955 done_tabledata(&(td[k]));
957 sfree(td);
959 return table;
962 t_forcetable make_gb_table(FILE *out,const output_env_t oenv,
963 const t_forcerec *fr,
964 const char *fn,
965 real rtab)
967 const char *fns[3] = { "gbctab.xvg", "gbdtab.xvg", "gbrtab.xvg" };
968 const char *fns14[3] = { "gbctab14.xvg", "gbdtab14.xvg", "gbrtab14.xvg" };
969 FILE *fp;
970 t_tabledata *td;
971 bool bReadTab,bGenTab;
972 real x0,y0,yp;
973 int i,j,k,nx,nx0,tabsel[etiNR];
974 void * p_tmp;
975 double r,r2,Vtab,Ftab,expterm;
977 t_forcetable table;
979 double abs_error_r, abs_error_r2;
980 double rel_error_r, rel_error_r2;
981 double rel_error_r_old=0, rel_error_r2_old=0;
982 double x0_r_error, x0_r2_error;
985 /* Only set a Coulomb table for GB */
987 tabsel[0]=etabGB;
988 tabsel[1]=-1;
989 tabsel[2]=-1;
992 /* Set the table dimensions for GB, not really necessary to
993 * use etiNR (since we only have one table, but ...)
995 snew(td,1);
996 table.r = fr->gbtabr;
997 table.scale = fr->gbtabscale;
998 table.scale_exp = 0;
999 table.n = table.scale*table.r;
1000 nx0 = 0;
1001 nx = table.scale*table.r;
1003 /* Check whether we have to read or generate
1004 * We will always generate a table, so remove the read code
1005 * (Compare with original make_table function
1007 bReadTab = FALSE;
1008 bGenTab = TRUE;
1010 /* Each table type (e.g. coul,lj6,lj12) requires four
1011 * numbers per datapoint. For performance reasons we want
1012 * the table data to be aligned to 16-byte. This is accomplished
1013 * by allocating 16 bytes extra to a temporary pointer, and then
1014 * calculating an aligned pointer. This new pointer must not be
1015 * used in a free() call, but thankfully we're sloppy enough not
1016 * to do this :-)
1019 /* 4 fp entries per table point, nx+1 points, and 16 bytes extra
1020 to align it. */
1021 p_tmp = malloc(4*(nx+1)*sizeof(real)+16);
1023 /* align it - size_t has the same same as a pointer */
1024 table.tab = (real *) (((size_t) p_tmp + 16) & (~((size_t) 15)));
1026 init_table(out,nx,nx0,table.scale,&(td[0]),!bReadTab);
1028 /* Local implementation so we don't have to use the etabGB
1029 * enum above, which will cause problems later when
1030 * making the other tables (right now even though we are using
1031 * GB, the normal Coulomb tables will be created, but this
1032 * will cause a problem since fr->eeltype==etabGB which will not
1033 * be defined in fill_table and set_table_type
1036 for(i=nx0;i<nx;i++)
1038 Vtab = 0.0;
1039 Ftab = 0.0;
1040 r = td->x[i];
1041 r2 = r*r;
1042 expterm = exp(-0.25*r2);
1044 Vtab = 1/sqrt(r2+expterm);
1045 Ftab = (r-0.25*r*expterm)/((r2+expterm)*sqrt(r2+expterm));
1047 /* Convert to single precision when we store to mem */
1048 td->v[i] = Vtab;
1049 td->f[i] = Ftab;
1053 copy2table(table.n,0,4,td[0].x,td[0].v,td[0].f,table.tab);
1055 if(bDebugMode())
1057 fp=xvgropen(fns[0],fns[0],"r","V",oenv);
1058 /* plot the output 5 times denser than the table data */
1059 /* for(i=5*nx0;i<5*table.n;i++) */
1060 for(i=nx0;i<table.n;i++)
1062 /* x0=i*table.r/(5*table.n); */
1063 x0=i*table.r/table.n;
1064 evaluate_table(table.tab,0,4,table.scale,x0,&y0,&yp);
1065 fprintf(fp,"%15.10e %15.10e %15.10e\n",x0,y0,yp);
1068 ffclose(fp);
1072 for(i=100*nx0;i<99.81*table.n;i++)
1074 r = i*table.r/(100*table.n);
1075 r2 = r*r;
1076 expterm = exp(-0.25*r2);
1078 Vtab = 1/sqrt(r2+expterm);
1079 Ftab = (r-0.25*r*expterm)/((r2+expterm)*sqrt(r2+expterm));
1082 evaluate_table(table.tab,0,4,table.scale,r,&y0,&yp);
1083 printf("gb: i=%d, x0=%g, y0=%15.15f, Vtab=%15.15f, yp=%15.15f, Ftab=%15.15f\n",i,r, y0, Vtab, yp, Ftab);
1085 abs_error_r=fabs(y0-Vtab);
1086 abs_error_r2=fabs(yp-(-1)*Ftab);
1088 rel_error_r=abs_error_r/y0;
1089 rel_error_r2=fabs(abs_error_r2/yp);
1092 if(rel_error_r>rel_error_r_old)
1094 rel_error_r_old=rel_error_r;
1095 x0_r_error=x0;
1098 if(rel_error_r2>rel_error_r2_old)
1100 rel_error_r2_old=rel_error_r2;
1101 x0_r2_error=x0;
1105 printf("gb: MAX REL ERROR IN R=%15.15f, MAX REL ERROR IN R2=%15.15f\n",rel_error_r_old, rel_error_r2_old);
1106 printf("gb: XO_R=%g, X0_R2=%g\n",x0_r_error, x0_r2_error);
1108 exit(1); */
1109 done_tabledata(&(td[0]));
1110 sfree(td);
1112 return table;
1117 bondedtable_t make_bonded_table(FILE *fplog,char *fn,int angle)
1119 t_tabledata td;
1120 double start;
1121 int i;
1122 bondedtable_t tab;
1124 if (angle < 2)
1125 start = 0;
1126 else
1127 start = -180.0;
1128 read_tables(fplog,fn,1,angle,&td);
1129 if (angle > 0) {
1130 /* Convert the table from degrees to radians */
1131 for(i=0; i<td.nx; i++) {
1132 td.x[i] *= DEG2RAD;
1133 td.f[i] *= RAD2DEG;
1135 td.tabscale *= RAD2DEG;
1137 tab.n = td.nx;
1138 tab.scale = td.tabscale;
1139 snew(tab.tab,tab.n*4);
1140 copy2table(tab.n,0,4,td.x,td.v,td.f,tab.tab);
1141 done_tabledata(&td);
1143 return tab;