1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROwing Monsters And Cloning Shrimps
45 #include "gmx_fatal.h"
57 /* All the possible (implemented) table functions */
82 /** Evaluates to true if the table type contains user data. */
83 #define ETAB_USER(e) ((e) == etabUSER || \
84 (e) == etabEwaldUser || (e) == etabEwaldUserSwitch)
91 /* This structure holds name and a flag that tells whether
92 this is a Coulomb type funtion */
93 static const t_tab_props tprops[etabNR] = {
96 { "LJ6Shift", FALSE },
97 { "LJ12Shift", FALSE },
103 { "Ewald-Switch", TRUE },
104 { "Ewald-User", TRUE },
105 { "Ewald-User-Switch", TRUE },
106 { "LJ6Switch", FALSE },
107 { "LJ12Switch", FALSE },
108 { "COULSwitch", TRUE },
109 { "LJ6-Encad shift", FALSE },
110 { "LJ12-Encad shift", FALSE },
111 { "COUL-Encad shift", TRUE },
116 /* Index in the table that says which function to use */
118 etiCOUL, etiLJ6, etiLJ12, etiNR
127 #define pow2(x) ((x)*(x))
128 #define pow3(x) ((x)*(x)*(x))
129 #define pow4(x) ((x)*(x)*(x)*(x))
130 #define pow5(x) ((x)*(x)*(x)*(x)*(x))
133 static double v_ewald_lr(double beta, double r)
137 return beta*2/sqrt(M_PI);
141 return gmx_erfd(beta*r)/r;
145 void table_spline3_fill_ewald_lr(real *table_f,
155 gmx_bool bOutOfRange;
156 double v_r0, v_r1, v_inrange, vi, a0, a1, a2dx;
161 gmx_fatal(FARGS, "Can not make a spline table with less than 2 points");
164 /* We need some margin to be able to divide table values by r
165 * in the kernel and also to do the integration arithmetics
166 * without going out of range. Furthemore, we divide by dx below.
168 tab_max = GMX_REAL_MAX*0.0001;
170 /* This function produces a table with:
171 * maximum energy error: V'''/(6*12*sqrt(3))*dx^3
172 * maximum force error: V'''/(6*4)*dx^2
173 * The rms force error is the max error times 1/sqrt(5)=0.45.
180 for (i = ntab-1; i >= 0; i--)
184 v_r0 = v_ewald_lr(beta, x_r0);
195 /* Linear continuation for the last point in range */
196 vi = v_inrange - dc*(i - i_inrange)*dx;
209 /* Get the potential at table point i-1 */
210 v_r1 = v_ewald_lr(beta, (i-1)*dx);
212 if (v_r1 != v_r1 || v_r1 < -tab_max || v_r1 > tab_max)
219 /* Calculate the average second derivative times dx over interval i-1 to i.
220 * Using the function values at the end points and in the middle.
222 a2dx = (v_r0 + v_r1 - 2*v_ewald_lr(beta, x_r0-0.5*dx))/(0.25*dx);
223 /* Set the derivative of the spline to match the difference in potential
224 * over the interval plus the average effect of the quadratic term.
225 * This is the essential step for minimizing the error in the force.
227 dc = (v_r0 - v_r1)/dx + 0.5*a2dx;
232 /* Fill the table with the force, minus the derivative of the spline */
237 /* tab[i] will contain the average of the splines over the two intervals */
238 table_f[i] += -0.5*dc;
243 /* Make spline s(x) = a0 + a1*(x - xr) + 0.5*a2*(x - xr)^2
244 * matching the potential at the two end points
245 * and the derivative dc at the end point xr.
249 a2dx = (a1*dx + v_r1 - a0)*2/dx;
251 /* Set dc to the derivative at the next point */
254 if (dc_new != dc_new || dc_new < -tab_max || dc_new > tab_max)
264 table_f[(i-1)] = -0.5*dc;
266 /* Currently the last value only contains half the force: double it */
269 if (table_v != NULL && table_fdv0 != NULL)
271 /* Copy to FDV0 table too. Allocation occurs in forcerec.c,
272 * init_ewald_f_table().
274 for (i = 0; i < ntab-1; i++)
276 table_fdv0[4*i] = table_f[i];
277 table_fdv0[4*i+1] = table_f[i+1]-table_f[i];
278 table_fdv0[4*i+2] = table_v[i];
279 table_fdv0[4*i+3] = 0.0;
281 table_fdv0[4*(ntab-1)] = table_f[(ntab-1)];
282 table_fdv0[4*(ntab-1)+1] = -table_f[(ntab-1)];
283 table_fdv0[4*(ntab-1)+2] = table_v[(ntab-1)];
284 table_fdv0[4*(ntab-1)+3] = 0.0;
288 /* The scale (1/spacing) for third order spline interpolation
289 * of the Ewald mesh contribution which needs to be subtracted
290 * from the non-bonded interactions.
292 real ewald_spline3_table_scale(real ewaldcoeff, real rc)
294 double erf_x_d3 = 1.0522; /* max of (erf(x)/x)''' */
298 /* Force tolerance: single precision accuracy */
299 ftol = GMX_FLOAT_EPS;
300 sc_f = sqrt(erf_x_d3/(6*4*ftol*ewaldcoeff))*ewaldcoeff;
302 /* Energy tolerance: 10x more accurate than the cut-off jump */
303 etol = 0.1*gmx_erfc(ewaldcoeff*rc);
304 etol = max(etol, GMX_REAL_EPS);
305 sc_e = pow(erf_x_d3/(6*12*sqrt(3)*etol), 1.0/3.0)*ewaldcoeff;
307 return max(sc_f, sc_e);
310 /* Calculate the potential and force for an r value
311 * in exactly the same way it is done in the inner loop.
312 * VFtab is a pointer to the table data, offset is
313 * the point where we should begin and stride is
314 * 4 if we have a buckingham table, 3 otherwise.
315 * If you want to evaluate table no N, set offset to 4*N.
317 * We use normal precision here, since that is what we
318 * will use in the inner loops.
320 static void evaluate_table(real VFtab[], int offset, int stride,
321 real tabscale, real r, real *y, real *yp)
325 real Y, F, Geps, Heps2, Fp;
334 Geps = eps*VFtab[n+2];
335 Heps2 = eps2*VFtab[n+3];
338 *yp = (Fp+Geps+2.0*Heps2)*tabscale;
341 static void copy2table(int n, int offset, int stride,
342 double x[], double Vtab[], double Ftab[], real scalefactor,
345 /* Use double prec. for the intermediary variables
346 * and temporary x/vtab/vtab2 data to avoid unnecessary
353 for (i = 0; (i < n); i++)
359 G = 3*(Vtab[i+1] - Vtab[i]) + (Ftab[i+1] + 2*Ftab[i])*h;
360 H = -2*(Vtab[i+1] - Vtab[i]) - (Ftab[i+1] + Ftab[i])*h;
364 /* Fill the last entry with a linear potential,
365 * this is mainly for rounding issues with angle and dihedral potentials.
371 nn0 = offset + i*stride;
372 dest[nn0] = scalefactor*Vtab[i];
373 dest[nn0+1] = scalefactor*F;
374 dest[nn0+2] = scalefactor*G;
375 dest[nn0+3] = scalefactor*H;
379 static void init_table(FILE *fp, int n, int nx0,
380 double tabscale, t_tabledata *td, gmx_bool bAlloc)
386 td->tabscale = tabscale;
393 for (i = 0; (i < td->nx); i++)
395 td->x[i] = i/tabscale;
399 static void spline_forces(int nx, double h, double v[], gmx_bool bS3, gmx_bool bE3,
403 double v3, b_s, b_e, b;
406 /* Formulas can be found in:
407 * H.J.C. Berendsen, Simulating the Physical World, Cambridge 2007
410 if (nx < 4 && (bS3 || bE3))
412 gmx_fatal(FARGS, "Can not generate splines with third derivative boundary conditions with less than 4 (%d) points", nx);
415 /* To make life easy we initially set the spacing to 1
416 * and correct for this at the end.
421 /* Fit V''' at the start */
422 v3 = v[3] - 3*v[2] + 3*v[1] - v[0];
425 fprintf(debug, "The left third derivative is %g\n", v3/(h*h*h));
427 b_s = 2*(v[1] - v[0]) + v3/6;
432 /* Fit V'' at the start */
435 v2 = -v[3] + 4*v[2] - 5*v[1] + 2*v[0];
436 /* v2 = v[2] - 2*v[1] + v[0]; */
439 fprintf(debug, "The left second derivative is %g\n", v2/(h*h));
441 b_s = 3*(v[1] - v[0]) - v2/2;
447 b_s = 3*(v[2] - v[0]) + f[0]*h;
452 /* Fit V''' at the end */
453 v3 = v[nx-1] - 3*v[nx-2] + 3*v[nx-3] - v[nx-4];
456 fprintf(debug, "The right third derivative is %g\n", v3/(h*h*h));
458 b_e = 2*(v[nx-1] - v[nx-2]) + v3/6;
463 /* V'=0 at the end */
464 b_e = 3*(v[nx-1] - v[nx-3]) + f[nx-1]*h;
469 beta = (bS3 ? 1 : 4);
471 /* For V'' fitting */
472 /* beta = (bS3 ? 2 : 4); */
475 for (i = start+1; i < end; i++)
479 b = 3*(v[i+1] - v[i-1]);
480 f[i] = (b - f[i-1])/beta;
482 gamma[end-1] = 1/beta;
483 beta = (bE3 ? 1 : 4) - gamma[end-1];
484 f[end-1] = (b_e - f[end-2])/beta;
486 for (i = end-2; i >= start; i--)
488 f[i] -= gamma[i+1]*f[i+1];
492 /* Correct for the minus sign and the spacing */
493 for (i = start; i < end; i++)
499 static void set_forces(FILE *fp, int angle,
500 int nx, double h, double v[], double f[],
508 "Force generation for dihedral tables is not (yet) implemented");
512 while (v[start] == 0)
518 while (v[end-1] == 0)
533 fprintf(fp, "Generating forces for table %d, boundary conditions: V''' at %g, %s at %g\n",
534 table+1, start*h, end == nx ? "V'''" : "V'=0", (end-1)*h);
536 spline_forces(end-start, h, v+start, TRUE, end == nx, f+start);
539 static void read_tables(FILE *fp, const char *fn,
540 int ntab, int angle, t_tabledata td[])
544 double **yy = NULL, start, end, dx0, dx1, ssd, vm, vp, f, numf;
545 int k, i, nx, nx0 = 0, ny, nny, ns;
546 gmx_bool bAllZero, bZeroV, bZeroF;
550 libfn = gmxlibfn(fn);
551 nx = read_xvg(libfn, &yy, &ny);
554 gmx_fatal(FARGS, "Trying to read file %s, but nr columns = %d, should be %d",
562 "The first distance in file %s is %f nm instead of %f nm",
563 libfn, yy[0][0], 0.0);
577 if (yy[0][0] != start || yy[0][nx-1] != end)
579 gmx_fatal(FARGS, "The angles in file %s should go from %f to %f instead of %f to %f\n",
580 libfn, start, end, yy[0][0], yy[0][nx-1]);
584 tabscale = (nx-1)/(yy[0][nx-1] - yy[0][0]);
588 fprintf(fp, "Read user tables from %s with %d data points.\n", libfn, nx);
591 fprintf(fp, "Tabscale = %g points/nm\n", tabscale);
596 for (k = 0; k < ntab; k++)
600 for (i = 0; (i < nx); i++)
604 dx0 = yy[0][i-1] - yy[0][i-2];
605 dx1 = yy[0][i] - yy[0][i-1];
606 /* Check for 1% deviation in spacing */
607 if (fabs(dx1 - dx0) >= 0.005*(fabs(dx0) + fabs(dx1)))
609 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]);
612 if (yy[1+k*2][i] != 0)
620 if (yy[1+k*2][i] > 0.01*GMX_REAL_MAX ||
621 yy[1+k*2][i] < -0.01*GMX_REAL_MAX)
623 gmx_fatal(FARGS, "Out of range potential value %g in file '%s'",
627 if (yy[1+k*2+1][i] != 0)
635 if (yy[1+k*2+1][i] > 0.01*GMX_REAL_MAX ||
636 yy[1+k*2+1][i] < -0.01*GMX_REAL_MAX)
638 gmx_fatal(FARGS, "Out of range force value %g in file '%s'",
644 if (!bZeroV && bZeroF)
646 set_forces(fp, angle, nx, 1/tabscale, yy[1+k*2], yy[1+k*2+1], k);
650 /* Check if the second column is close to minus the numerical
651 * derivative of the first column.
655 for (i = 1; (i < nx-1); i++)
660 if (vm != 0 && vp != 0 && f != 0)
662 /* Take the centered difference */
663 numf = -(vp - vm)*0.5*tabscale;
664 ssd += fabs(2*(f - numf)/(f + numf));
671 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));
674 fprintf(debug, "%s", buf);
680 fprintf(fp, "\nWARNING: %s\n", buf);
682 fprintf(stderr, "\nWARNING: %s\n", buf);
689 fprintf(fp, "\nNOTE: All elements in table %s are zero\n\n", libfn);
692 for (k = 0; (k < ntab); k++)
694 init_table(fp, nx, nx0, tabscale, &(td[k]), TRUE);
695 for (i = 0; (i < nx); i++)
697 td[k].x[i] = yy[0][i];
698 td[k].v[i] = yy[2*k+1][i];
699 td[k].f[i] = yy[2*k+2][i];
702 for (i = 0; (i < ny); i++)
710 static void done_tabledata(t_tabledata *td)
724 static void fill_table(t_tabledata *td, int tp, const t_forcerec *fr)
726 /* Fill the table according to the formulas in the manual.
727 * In principle, we only need the potential and the second
728 * derivative, but then we would have to do lots of calculations
729 * in the inner loop. By precalculating some terms (see manual)
730 * we get better eventual performance, despite a larger table.
732 * Since some of these higher-order terms are very small,
733 * we always use double precision to calculate them here, in order
734 * to avoid unnecessary loss of precision.
741 double r1, rc, r12, r13;
742 double r, r2, r6, rc6;
743 double expr, Vtab, Ftab;
744 /* Parameters for David's function */
745 double A = 0, B = 0, C = 0, A_3 = 0, B_4 = 0;
746 /* Parameters for the switching function */
747 double ksw, swi, swi1;
748 /* Temporary parameters */
749 gmx_bool bSwitch, bShift;
750 double ewc = fr->ewaldcoeff;
752 bSwitch = ((tp == etabLJ6Switch) || (tp == etabLJ12Switch) ||
753 (tp == etabCOULSwitch) ||
754 (tp == etabEwaldSwitch) || (tp == etabEwaldUserSwitch));
755 bShift = ((tp == etabLJ6Shift) || (tp == etabLJ12Shift) ||
760 if (tprops[tp].bCoulomb)
762 r1 = fr->rcoulomb_switch;
767 r1 = fr->rvdw_switch;
772 ksw = 1.0/(pow5(rc-r1));
784 else if (tp == etabLJ6Shift)
793 A = p * ((p+1)*r1-(p+4)*rc)/(pow(rc, p+2)*pow2(rc-r1));
794 B = -p * ((p+1)*r1-(p+3)*rc)/(pow(rc, p+2)*pow3(rc-r1));
795 C = 1.0/pow(rc, p)-A/3.0*pow3(rc-r1)-B/4.0*pow4(rc-r1);
796 if (tp == etabLJ6Shift)
807 fprintf(debug, "Setting up tables\n"); fflush(debug);
811 fp = xvgropen("switch.xvg", "switch", "r", "s");
814 for (i = td->nx0; (i < td->nx); i++)
819 if (gmx_within_tol(reppow, 12.0, 10*GMX_DOUBLE_EPS))
825 r12 = pow(r, -reppow);
831 /* swi is function, swi1 1st derivative and swi2 2nd derivative */
832 /* The switch function is 1 for r<r1, 0 for r>rc, and smooth for
833 * r1<=r<=rc. The 1st and 2nd derivatives are both zero at
835 * ksw is just the constant 1/(rc-r1)^5, to save some calculations...
849 swi = 1 - 10*pow3(r-r1)*ksw*pow2(rc-r1)
850 + 15*pow4(r-r1)*ksw*(rc-r1) - 6*pow5(r-r1)*ksw;
851 swi1 = -30*pow2(r-r1)*ksw*pow2(rc-r1)
852 + 60*pow3(r-r1)*ksw*(rc-r1) - 30*pow4(r-r1)*ksw;
855 else /* not really needed, but avoids compiler warnings... */
861 fprintf(fp, "%10g %10g %10g %10g\n", r, swi, swi1, swi2);
887 Ftab = reppow*Vtab/r;
895 Ftab = reppow*Vtab/r;
901 Vtab = -(r6-6.0*(rc-r)*rc6/rc-rc6);
902 Ftab = -(6.0*r6/r-6.0*rc6/rc);
913 Vtab = -(r6-6.0*(rc-r)*rc6/rc-rc6);
914 Ftab = -(6.0*r6/r-6.0*rc6/rc);
935 case etabEwaldSwitch:
936 Vtab = gmx_erfc(ewc*r)/r;
937 Ftab = gmx_erfc(ewc*r)/r2+exp(-(ewc*ewc*r2))*ewc*M_2_SQRTPI/r;
940 case etabEwaldUserSwitch:
941 /* Only calculate minus the reciprocal space contribution */
942 Vtab = -gmx_erf(ewc*r)/r;
943 Ftab = -gmx_erf(ewc*r)/r2+exp(-(ewc*ewc*r2))*ewc*M_2_SQRTPI/r;
947 Vtab = 1.0/r + fr->k_rf*r2 - fr->c_rf;
948 Ftab = 1.0/r2 - 2*fr->k_rf*r;
949 if (tp == etabRF_ZERO && r >= rc)
963 Vtab = 1.0/r-(rc-r)/(rc*rc)-1.0/rc;
964 Ftab = 1.0/r2-1.0/(rc*rc);
973 gmx_fatal(FARGS, "Table type %d not implemented yet. (%s,%d)",
974 tp, __FILE__, __LINE__);
978 /* Normal coulomb with cut-off correction for potential */
982 /* If in Shifting range add something to it */
987 Vtab += -A_3*r13 - B_4*r12*r12;
988 Ftab += A*r12 + B*r13;
999 if ((r > r1) && bSwitch)
1001 Ftab = Ftab*swi - Vtab*swi1;
1005 /* Convert to single precision when we store to mem */
1010 /* Continue the table linearly from nx0 to 0.
1011 * These values are only required for energy minimization with overlap or TPI.
1013 for (i = td->nx0-1; i >= 0; i--)
1015 td->v[i] = td->v[i+1] + td->f[i+1]*(td->x[i+1] - td->x[i]);
1016 td->f[i] = td->f[i+1];
1024 static void set_table_type(int tabsel[], const t_forcerec *fr, gmx_bool b14only)
1026 int eltype, vdwtype;
1028 /* Set the different table indices.
1035 switch (fr->eeltype)
1042 case eelPMEUSERSWITCH:
1051 eltype = fr->eeltype;
1057 tabsel[etiCOUL] = etabCOUL;
1060 tabsel[etiCOUL] = etabShift;
1063 if (fr->rcoulomb > fr->rcoulomb_switch)
1065 tabsel[etiCOUL] = etabShift;
1069 tabsel[etiCOUL] = etabCOUL;
1075 tabsel[etiCOUL] = etabEwald;
1078 tabsel[etiCOUL] = etabEwaldSwitch;
1081 tabsel[etiCOUL] = etabEwaldUser;
1083 case eelPMEUSERSWITCH:
1084 tabsel[etiCOUL] = etabEwaldUserSwitch;
1089 tabsel[etiCOUL] = etabRF;
1092 tabsel[etiCOUL] = etabRF_ZERO;
1095 tabsel[etiCOUL] = etabCOULSwitch;
1098 tabsel[etiCOUL] = etabUSER;
1101 tabsel[etiCOUL] = etabCOULEncad;
1104 gmx_fatal(FARGS, "Invalid eeltype %d", eltype);
1107 /* Van der Waals time */
1108 if (fr->bBHAM && !b14only)
1110 tabsel[etiLJ6] = etabLJ6;
1111 tabsel[etiLJ12] = etabEXPMIN;
1115 if (b14only && fr->vdwtype != evdwUSER)
1121 vdwtype = fr->vdwtype;
1127 tabsel[etiLJ6] = etabLJ6Switch;
1128 tabsel[etiLJ12] = etabLJ12Switch;
1131 tabsel[etiLJ6] = etabLJ6Shift;
1132 tabsel[etiLJ12] = etabLJ12Shift;
1135 tabsel[etiLJ6] = etabUSER;
1136 tabsel[etiLJ12] = etabUSER;
1139 tabsel[etiLJ6] = etabLJ6;
1140 tabsel[etiLJ12] = etabLJ12;
1142 case evdwENCADSHIFT:
1143 tabsel[etiLJ6] = etabLJ6Encad;
1144 tabsel[etiLJ12] = etabLJ12Encad;
1147 gmx_fatal(FARGS, "Invalid vdwtype %d in %s line %d", vdwtype,
1148 __FILE__, __LINE__);
1153 t_forcetable make_tables(FILE *out, const output_env_t oenv,
1154 const t_forcerec *fr,
1155 gmx_bool bVerbose, const char *fn,
1156 real rtab, int flags)
1158 const char *fns[3] = { "ctab.xvg", "dtab.xvg", "rtab.xvg" };
1159 const char *fns14[3] = { "ctab14.xvg", "dtab14.xvg", "rtab14.xvg" };
1162 gmx_bool b14only, bReadTab, bGenTab;
1164 int i, j, k, nx, nx0, tabsel[etiNR];
1169 b14only = (flags & GMX_MAKETABLES_14ONLY);
1171 if (flags & GMX_MAKETABLES_FORCEUSER)
1173 tabsel[etiCOUL] = etabUSER;
1174 tabsel[etiLJ6] = etabUSER;
1175 tabsel[etiLJ12] = etabUSER;
1179 set_table_type(tabsel, fr, b14only);
1185 table.scale_exp = 0;
1189 table.interaction = GMX_TABLE_INTERACTION_ELEC_VDWREP_VDWDISP;
1190 table.format = GMX_TABLE_FORMAT_CUBICSPLINE_YFGH;
1191 table.formatsize = 4;
1192 table.ninteractions = 3;
1193 table.stride = table.formatsize*table.ninteractions;
1195 /* Check whether we have to read or generate */
1198 for (i = 0; (i < etiNR); i++)
1200 if (ETAB_USER(tabsel[i]))
1204 if (tabsel[i] != etabUSER)
1211 read_tables(out, fn, etiNR, 0, td);
1212 if (rtab == 0 || (flags & GMX_MAKETABLES_14ONLY))
1214 rtab = td[0].x[td[0].nx-1];
1220 if (td[0].x[td[0].nx-1] < rtab)
1222 gmx_fatal(FARGS, "Tables in file %s not long enough for cut-off:\n"
1223 "\tshould be at least %f nm\n", fn, rtab);
1225 nx = table.n = (int)(rtab*td[0].tabscale + 0.5);
1227 table.scale = td[0].tabscale;
1235 table.scale = 2000.0;
1237 table.scale = 500.0;
1239 nx = table.n = rtab*table.scale;
1244 if (fr->bham_b_max != 0)
1246 table.scale_exp = table.scale/fr->bham_b_max;
1250 table.scale_exp = table.scale;
1254 /* Each table type (e.g. coul,lj6,lj12) requires four
1255 * numbers per nx+1 data points. For performance reasons we want
1256 * the table data to be aligned to 16-byte.
1258 snew_aligned(table.data, 12*(nx+1)*sizeof(real), 32);
1260 for (k = 0; (k < etiNR); k++)
1262 if (tabsel[k] != etabUSER)
1264 init_table(out, nx, nx0,
1265 (tabsel[k] == etabEXPMIN) ? table.scale_exp : table.scale,
1266 &(td[k]), !bReadTab);
1267 fill_table(&(td[k]), tabsel[k], fr);
1270 fprintf(out, "%s table with %d data points for %s%s.\n"
1271 "Tabscale = %g points/nm\n",
1272 ETAB_USER(tabsel[k]) ? "Modified" : "Generated",
1273 td[k].nx, b14only ? "1-4 " : "", tprops[tabsel[k]].name,
1278 /* Set scalefactor for c6/c12 tables. This is because we save flops in the non-table kernels
1279 * by including the derivative constants (6.0 or 12.0) in the parameters, since
1280 * we no longer calculate force in most steps. This means the c6/c12 parameters
1281 * have been scaled up, so we need to scale down the table interactions too.
1282 * It comes here since we need to scale user tables too.
1286 scalefactor = 1.0/6.0;
1288 else if (k == etiLJ12 && tabsel[k] != etabEXPMIN)
1290 scalefactor = 1.0/12.0;
1297 copy2table(table.n, k*4, 12, td[k].x, td[k].v, td[k].f, scalefactor, table.data);
1299 if (bDebugMode() && bVerbose)
1303 fp = xvgropen(fns14[k], fns14[k], "r", "V", oenv);
1307 fp = xvgropen(fns[k], fns[k], "r", "V", oenv);
1309 /* plot the output 5 times denser than the table data */
1310 for (i = 5*((nx0+1)/2); i < 5*table.n; i++)
1312 x0 = i*table.r/(5*(table.n-1));
1313 evaluate_table(table.data, 4*k, 12, table.scale, x0, &y0, &yp);
1314 fprintf(fp, "%15.10e %15.10e %15.10e\n", x0, y0, yp);
1318 done_tabledata(&(td[k]));
1325 t_forcetable make_gb_table(FILE *out, const output_env_t oenv,
1326 const t_forcerec *fr,
1330 const char *fns[3] = { "gbctab.xvg", "gbdtab.xvg", "gbrtab.xvg" };
1331 const char *fns14[3] = { "gbctab14.xvg", "gbdtab14.xvg", "gbrtab14.xvg" };
1334 gmx_bool bReadTab, bGenTab;
1336 int i, j, k, nx, nx0, tabsel[etiNR];
1337 double r, r2, Vtab, Ftab, expterm;
1341 double abs_error_r, abs_error_r2;
1342 double rel_error_r, rel_error_r2;
1343 double rel_error_r_old = 0, rel_error_r2_old = 0;
1344 double x0_r_error, x0_r2_error;
1347 /* Only set a Coulomb table for GB */
1354 /* Set the table dimensions for GB, not really necessary to
1355 * use etiNR (since we only have one table, but ...)
1358 table.interaction = GMX_TABLE_INTERACTION_ELEC;
1359 table.format = GMX_TABLE_FORMAT_CUBICSPLINE_YFGH;
1360 table.r = fr->gbtabr;
1361 table.scale = fr->gbtabscale;
1362 table.scale_exp = 0;
1363 table.n = table.scale*table.r;
1364 table.formatsize = 4;
1365 table.ninteractions = 1;
1366 table.stride = table.formatsize*table.ninteractions;
1368 nx = table.scale*table.r;
1370 /* Check whether we have to read or generate
1371 * We will always generate a table, so remove the read code
1372 * (Compare with original make_table function
1377 /* Each table type (e.g. coul,lj6,lj12) requires four
1378 * numbers per datapoint. For performance reasons we want
1379 * the table data to be aligned to 16-byte. This is accomplished
1380 * by allocating 16 bytes extra to a temporary pointer, and then
1381 * calculating an aligned pointer. This new pointer must not be
1382 * used in a free() call, but thankfully we're sloppy enough not
1386 snew_aligned(table.data, 4*nx, 32);
1388 init_table(out, nx, nx0, table.scale, &(td[0]), !bReadTab);
1390 /* Local implementation so we don't have to use the etabGB
1391 * enum above, which will cause problems later when
1392 * making the other tables (right now even though we are using
1393 * GB, the normal Coulomb tables will be created, but this
1394 * will cause a problem since fr->eeltype==etabGB which will not
1395 * be defined in fill_table and set_table_type
1398 for (i = nx0; i < nx; i++)
1404 expterm = exp(-0.25*r2);
1406 Vtab = 1/sqrt(r2+expterm);
1407 Ftab = (r-0.25*r*expterm)/((r2+expterm)*sqrt(r2+expterm));
1409 /* Convert to single precision when we store to mem */
1415 copy2table(table.n, 0, 4, td[0].x, td[0].v, td[0].f, 1.0, table.data);
1419 fp = xvgropen(fns[0], fns[0], "r", "V", oenv);
1420 /* plot the output 5 times denser than the table data */
1421 /* for(i=5*nx0;i<5*table.n;i++) */
1422 for (i = nx0; i < table.n; i++)
1424 /* x0=i*table.r/(5*table.n); */
1425 x0 = i*table.r/table.n;
1426 evaluate_table(table.data, 0, 4, table.scale, x0, &y0, &yp);
1427 fprintf(fp, "%15.10e %15.10e %15.10e\n", x0, y0, yp);
1434 for(i=100*nx0;i<99.81*table.n;i++)
1436 r = i*table.r/(100*table.n);
1438 expterm = exp(-0.25*r2);
1440 Vtab = 1/sqrt(r2+expterm);
1441 Ftab = (r-0.25*r*expterm)/((r2+expterm)*sqrt(r2+expterm));
1444 evaluate_table(table.data,0,4,table.scale,r,&y0,&yp);
1445 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);
1447 abs_error_r=fabs(y0-Vtab);
1448 abs_error_r2=fabs(yp-(-1)*Ftab);
1450 rel_error_r=abs_error_r/y0;
1451 rel_error_r2=fabs(abs_error_r2/yp);
1454 if(rel_error_r>rel_error_r_old)
1456 rel_error_r_old=rel_error_r;
1460 if(rel_error_r2>rel_error_r2_old)
1462 rel_error_r2_old=rel_error_r2;
1467 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);
1468 printf("gb: XO_R=%g, X0_R2=%g\n",x0_r_error, x0_r2_error);
1471 done_tabledata(&(td[0]));
1479 t_forcetable make_atf_table(FILE *out, const output_env_t oenv,
1480 const t_forcerec *fr,
1484 const char *fns[3] = { "tf_tab.xvg", "atfdtab.xvg", "atfrtab.xvg" };
1487 real x0, y0, yp, rtab;
1489 real rx, ry, rz, box_r;
1494 /* Set the table dimensions for ATF, not really necessary to
1495 * use etiNR (since we only have one table, but ...)
1499 if (fr->adress_type == eAdressSphere)
1501 /* take half box diagonal direction as tab range */
1502 rx = 0.5*box[0][0]+0.5*box[1][0]+0.5*box[2][0];
1503 ry = 0.5*box[0][1]+0.5*box[1][1]+0.5*box[2][1];
1504 rz = 0.5*box[0][2]+0.5*box[1][2]+0.5*box[2][2];
1505 box_r = sqrt(rx*rx+ry*ry+rz*rz);
1510 /* xsplit: take half box x direction as tab range */
1511 box_r = box[0][0]/2;
1516 table.scale_exp = 0;
1520 read_tables(out, fn, 1, 0, td);
1521 rtab = td[0].x[td[0].nx-1];
1523 if (fr->adress_type == eAdressXSplit && (rtab < box[0][0]/2))
1525 gmx_fatal(FARGS, "AdResS full box therm force table in file %s extends to %f:\n"
1526 "\tshould extend to at least half the length of the box in x-direction"
1527 "%f\n", fn, rtab, box[0][0]/2);
1531 gmx_fatal(FARGS, "AdResS full box therm force table in file %s extends to %f:\n"
1532 "\tshould extend to at least for spherical adress"
1533 "%f (=distance from center to furthermost point in box \n", fn, rtab, box_r);
1539 table.scale = td[0].tabscale;
1542 /* Each table type (e.g. coul,lj6,lj12) requires four
1543 * numbers per datapoint. For performance reasons we want
1544 * the table data to be aligned to 16-byte. This is accomplished
1545 * by allocating 16 bytes extra to a temporary pointer, and then
1546 * calculating an aligned pointer. This new pointer must not be
1547 * used in a free() call, but thankfully we're sloppy enough not
1551 snew_aligned(table.data, 4*nx, 32);
1553 copy2table(table.n, 0, 4, td[0].x, td[0].v, td[0].f, 1.0, table.data);
1557 fp = xvgropen(fns[0], fns[0], "r", "V", oenv);
1558 /* plot the output 5 times denser than the table data */
1559 /* for(i=5*nx0;i<5*table.n;i++) */
1561 for (i = 5*((nx0+1)/2); i < 5*table.n; i++)
1563 /* x0=i*table.r/(5*table.n); */
1564 x0 = i*table.r/(5*(table.n-1));
1565 evaluate_table(table.data, 0, 4, table.scale, x0, &y0, &yp);
1566 fprintf(fp, "%15.10e %15.10e %15.10e\n", x0, y0, yp);
1572 done_tabledata(&(td[0]));
1575 table.interaction = GMX_TABLE_INTERACTION_ELEC_VDWREP_VDWDISP;
1576 table.format = GMX_TABLE_FORMAT_CUBICSPLINE_YFGH;
1577 table.formatsize = 4;
1578 table.ninteractions = 3;
1579 table.stride = table.formatsize*table.ninteractions;
1585 bondedtable_t make_bonded_table(FILE *fplog, char *fn, int angle)
1600 read_tables(fplog, fn, 1, angle, &td);
1603 /* Convert the table from degrees to radians */
1604 for (i = 0; i < td.nx; i++)
1609 td.tabscale *= RAD2DEG;
1612 tab.scale = td.tabscale;
1613 snew(tab.data, tab.n*4);
1614 copy2table(tab.n, 0, 4, td.x, td.v, td.f, 1.0, tab.data);
1615 done_tabledata(&td);