2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
43 #include "gromacs/math/utilities.h"
50 #include "gromacs/fileio/futil.h"
54 real calc_ewaldcoeff_q(real rc, real dtol)
56 real x = 5, low, high;
65 while (gmx_erfc(x*rc) > dtol);
67 n = i+60; /* search tolerance is 2^-60 */
70 for (i = 0; i < n; i++)
73 if (gmx_erfc(x*rc) > dtol)
85 static real ewald_function_lj(real x, real rc)
87 real xrc, xrc2, xrc4, factor;
92 factor = exp(-xrc2)*(1 + xrc2 + xrc4/2.0);
94 factor = expf(-xrc2)*(1 + xrc2 + xrc4/2.0);
100 real calc_ewaldcoeff_lj(real rc, real dtol)
102 real x = 5, low, high;
110 while (ewald_function_lj(x, rc) > dtol);
112 n = i + 60; /* search tolerance is 2^-60 */
115 for (i = 0; i < n; ++i)
117 x = (low + high) / 2.0;
118 if (ewald_function_lj(x, rc) > dtol)
130 void ewald_LRcorrection(int start, int end,
131 t_commrec *cr, int thread, t_forcerec *fr,
132 real *chargeA, real *chargeB,
133 real *C6A, real *C6B,
134 real *sigmaA, real *sigmaB,
135 real *sigma3A, real *sigma3B,
136 gmx_bool calc_excl_corr,
137 t_blocka *excl, rvec x[],
138 matrix box, rvec mu_tot[],
139 int ewald_geometry, real epsilon_surface,
140 rvec *f, tensor vir_q, tensor vir_lj,
141 real *Vcorr_q, real *Vcorr_lj,
142 real lambda_q, real lambda_lj,
143 real *dvdlambda_q, real *dvdlambda_lj)
145 int i, i1, i2, j, k, m, iv, jv, q;
147 double Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
150 real v, vc, qiA, qiB, dr, dr2, rinv, fscal, enercorr;
151 real Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ewaldcoeff_q, ewcdr;
152 real ewc_lj = fr->ewaldcoeff_lj;
153 real c6Ai = 0, c6Bi = 0, c6A = 0, c6B = 0, ewcdr2, ewcdr4, ewcdr5, c6L = 0, rinv6;
154 rvec df, dx, mutot[2], dipcorrA, dipcorrB;
155 tensor dxdf_q, dxdf_lj;
156 real vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
157 real L1_q, L1_lj, dipole_coeff, qqA, qqB, qqL, vr0_q, vr0_lj = 0;
158 gmx_bool bFreeEnergy = (chargeB != NULL);
159 gmx_bool bMolPBC = fr->bMolPBC;
160 gmx_bool bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
162 /* This routine can be made faster by using tables instead of analytical interactions
163 * However, that requires a thorough verification that they are correct in all cases.
166 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
167 vr0_q = ewc_q*M_2_SQRTPI;
168 if (EVDW_PME(fr->vdwtype))
170 vr0_lj = -pow(ewc_lj, 6)/6.0;
181 L1_lj = 1.0-lambda_lj;
182 /* Note that we have to transform back to gromacs units, since
183 * mu_tot contains the dipole in debye units (for output).
185 for (i = 0; (i < DIM); i++)
187 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
188 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
193 switch (ewald_geometry)
196 if (epsilon_surface != 0)
199 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
200 for (i = 0; (i < DIM); i++)
202 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
203 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
208 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
209 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
210 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
213 gmx_incons("Unsupported Ewald geometry");
218 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n",
219 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
220 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n",
221 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
224 if (EVDW_PME(fr->vdwtype))
228 if ((calc_excl_corr || dipole_coeff != 0 || EVDW_PME(fr->vdwtype)) && !bFreeEnergy)
230 for (i = start; (i < end); i++)
232 /* Initiate local variables (for this i-particle) to 0 */
233 qiA = chargeA[i]*one_4pi_eps;
234 if (EVDW_PME(fr->vdwtype))
245 i2 = excl->index[i+1];
247 /* Loop over excluded neighbours */
248 for (j = i1; (j < i2); j++)
252 * First we must test whether k <> i, and then,
253 * because the exclusions are all listed twice i->k
254 * and k->i we must select just one of the two. As
255 * a minor optimization we only compute forces when
256 * the charges are non-zero.
260 qqA = qiA*chargeA[k];
261 if (EVDW_PME(fr->vdwtype))
266 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
269 if (qqA != 0.0 || c6A != 0.0)
271 rvec_sub(x[i], x[k], dx);
275 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
276 for (m = DIM-1; (m >= 0); m--)
278 if (dx[m] > 0.5*box[m][m])
280 rvec_dec(dx, box[m]);
282 else if (dx[m] < -0.5*box[m][m])
284 rvec_inc(dx, box[m]);
289 /* Distance between two excluded particles
290 * may be zero in the case of shells
294 rinv = gmx_invsqrt(dr2);
300 vc = qqA*gmx_erf(ewcdr)*rinv;
303 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
304 #define R_ERF_R_INACC 0.006
306 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
307 #define R_ERF_R_INACC 0.1
309 if (ewcdr > R_ERF_R_INACC)
311 fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*exp(-ewcdr*ewcdr));
315 /* Use a fourth order series expansion for small ewcdr */
316 fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
319 /* The force vector is obtained by multiplication with
320 * the distance vector
322 svmul(fscal, dx, df);
325 for (iv = 0; (iv < DIM); iv++)
327 for (jv = 0; (jv < DIM); jv++)
329 dxdf_q[iv][jv] += dx[iv]*df[jv];
336 rinv6 = rinv2*rinv2*rinv2;
338 ewcdr2 = ewcdr*ewcdr;
339 ewcdr4 = ewcdr2*ewcdr2;
340 ewcdr5 = ewcdr4*ewcdr;
341 /* We get the excluded long-range contribution from -C6*(1-g(r))
342 * g(r) is also defined in the manual under LJ-PME
344 vc = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
346 /* The force is the derivative of the potential vc */
347 fscal = 6.0*vc*rinv + c6A*rinv6*exp(-ewcdr2)*ewc_lj*ewcdr5;
349 /* The force vector is obtained by
350 * multiplication with the distance
353 svmul(fscal, dx, df);
356 for (iv = 0; (iv < DIM); iv++)
358 for (jv = 0; (jv < DIM); jv++)
360 dxdf_lj[iv][jv] += dx[iv]*df[jv];
367 Vexcl_q += qqA*vr0_q;
368 Vexcl_lj += c6A*vr0_lj;
374 /* Dipole correction on force */
375 if (dipole_coeff != 0)
377 for (j = 0; (j < DIM); j++)
379 f[i][j] -= dipcorrA[j]*chargeA[i];
384 else if (calc_excl_corr || dipole_coeff != 0 || EVDW_PME(fr->vdwtype))
386 for (i = start; (i < end); i++)
388 /* Initiate local variables (for this i-particle) to 0 */
389 qiA = chargeA[i]*one_4pi_eps;
390 qiB = chargeB[i]*one_4pi_eps;
391 if (EVDW_PME(fr->vdwtype))
404 i2 = excl->index[i+1];
406 /* Loop over excluded neighbours */
407 for (j = i1; (j < i2); j++)
412 qqA = qiA*chargeA[k];
413 qqB = qiB*chargeB[k];
414 if (EVDW_PME(fr->vdwtype))
420 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
421 c6B *= pow(0.5*(sigmaB[i]+sigmaB[k]), 6)*sigma3B[k];
424 if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
427 qqL = L1_q*qqA + lambda_q*qqB;
428 if (EVDW_PME(fr->vdwtype))
430 c6L = L1_lj*c6A + lambda_lj*c6B;
432 rvec_sub(x[i], x[k], dx);
435 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
436 for (m = DIM-1; (m >= 0); m--)
438 if (dx[m] > 0.5*box[m][m])
440 rvec_dec(dx, box[m]);
442 else if (dx[m] < -0.5*box[m][m])
444 rvec_inc(dx, box[m]);
451 rinv = gmx_invsqrt(dr2);
454 if (qqA != 0.0 || qqB != 0.0)
456 v = gmx_erf(ewc_q*dr)*rinv;
459 fscal = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*exp(-ewc_q*ewc_q*dr2));
460 dvdl_excl_q += (qqB - qqA)*v;
462 svmul(fscal, dx, df);
465 for (iv = 0; (iv < DIM); iv++)
467 for (jv = 0; (jv < DIM); jv++)
469 dxdf_q[iv][jv] += dx[iv]*df[jv];
474 if ((c6A != 0.0 || c6B != 0.0) && EVDW_PME(fr->vdwtype))
476 rinv6 = rinv2*rinv2*rinv2;
478 ewcdr2 = ewcdr*ewcdr;
479 ewcdr4 = ewcdr2*ewcdr2;
480 ewcdr5 = ewcdr4*ewcdr;
481 v = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
484 fscal = 6.0*vc*rinv + c6L*rinv6*exp(-ewcdr2)*ewc_lj*ewcdr5;
485 dvdl_excl_lj += (c6B - c6A)*v;
486 svmul(fscal, dx, df);
489 for (iv = 0; (iv < DIM); iv++)
491 for (jv = 0; (jv < DIM); jv++)
493 dxdf_lj[iv][jv] += dx[iv]*df[jv];
500 Vexcl_q += qqL*vr0_q;
501 dvdl_excl_q += (qqB - qqA)*vr0_q;
502 Vexcl_lj += c6L*vr0_lj;
503 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
509 /* Dipole correction on force */
510 if (dipole_coeff != 0)
512 for (j = 0; (j < DIM); j++)
514 f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
515 + lambda_q*dipcorrB[j]*chargeB[i];
520 for (iv = 0; (iv < DIM); iv++)
522 for (jv = 0; (jv < DIM); jv++)
524 vir_q[iv][jv] += 0.5*dxdf_q[iv][jv];
525 vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
534 /* Global corrections only on master process */
535 if (MASTER(cr) && thread == 0)
537 for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
541 /* Self-energy correction */
542 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
543 if (EVDW_PME(fr->vdwtype))
545 Vself_lj[q] = -pow(ewc_lj, 6)*fr->c6sum[q];
549 /* Apply surface dipole correction:
550 * correction = dipole_coeff * (dipole)^2
552 if (dipole_coeff != 0)
554 if (ewald_geometry == eewg3D)
556 Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
558 else if (ewald_geometry == eewg3DC)
560 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
567 *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
568 if (EVDW_PME(fr->vdwtype))
570 *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
575 *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
576 + lambda_q*(Vdipole[1] - Vself_q[1])
578 *dvdlambda_q += Vdipole[1] - Vself_q[1]
579 - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
580 if (EVDW_PME(fr->vdwtype))
582 *Vcorr_lj = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
583 *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
589 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
590 fprintf(debug, "start=%d,natoms=%d\n", start, end-start);
591 fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
592 L1_q*fr->q2sum[0]+lambda_q*fr->q2sum[1], L1_q*Vself_q[0]+lambda_q*Vself_q[1], L1_lj*fr->c6sum[0]+lambda_lj*fr->c6sum[1], L1_lj*Vself_lj[0]+lambda_lj*Vself_lj[1]);
593 fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
594 fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
595 if (MASTER(cr) && thread == 0)
597 if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
599 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
600 L1_q*Vdipole[0]+lambda_q*Vdipole[1]);
606 real ewald_charge_correction(t_commrec *cr, t_forcerec *fr, real lambda,
608 real *dvdlambda, tensor vir)
611 real vol, fac, qs2A, qs2B, vc, enercorr;
616 /* Apply charge correction */
617 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
619 fac = M_PI*ONE_4PI_EPS0/(fr->epsilon_r*2.0*vol*vol*sqr(fr->ewaldcoeff_q));
621 qs2A = fr->qsum[0]*fr->qsum[0];
622 qs2B = fr->qsum[1]*fr->qsum[1];
624 vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
628 *dvdlambda += -vol*(qs2B - qs2A)*fac;
630 for (d = 0; d < DIM; d++)
637 fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);