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.
41 #include "gromacs/math/utilities.h"
42 #include "gromacs/legacyheaders/typedefs.h"
43 #include "gromacs/legacyheaders/types/commrec.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/coulomb.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/legacyheaders/txtdump.h"
49 #include "gromacs/utility/futil.h"
50 #include "gromacs/legacyheaders/names.h"
51 #include "gromacs/legacyheaders/macros.h"
53 real calc_ewaldcoeff_q(real rc, real dtol)
55 real x = 5, low, high;
64 while (gmx_erfc(x*rc) > dtol);
66 n = i+60; /* search tolerance is 2^-60 */
69 for (i = 0; i < n; i++)
72 if (gmx_erfc(x*rc) > dtol)
84 static real ewald_function_lj(real x, real rc)
86 real xrc, xrc2, xrc4, factor;
91 factor = exp(-xrc2)*(1 + xrc2 + xrc4/2.0);
93 factor = expf(-xrc2)*(1 + xrc2 + xrc4/2.0);
99 real calc_ewaldcoeff_lj(real rc, real dtol)
101 real x = 5, low, high;
109 while (ewald_function_lj(x, rc) > dtol);
111 n = i + 60; /* search tolerance is 2^-60 */
114 for (i = 0; i < n; ++i)
116 x = (low + high) / 2.0;
117 if (ewald_function_lj(x, rc) > dtol)
129 void ewald_LRcorrection(int start, int end,
130 t_commrec *cr, int thread, t_forcerec *fr,
131 real *chargeA, real *chargeB,
132 real *C6A, real *C6B,
133 real *sigmaA, real *sigmaB,
134 real *sigma3A, real *sigma3B,
135 gmx_bool calc_excl_corr,
136 t_blocka *excl, rvec x[],
137 matrix box, rvec mu_tot[],
138 int ewald_geometry, real epsilon_surface,
139 rvec *f, tensor vir_q, tensor vir_lj,
140 real *Vcorr_q, real *Vcorr_lj,
141 real lambda_q, real lambda_lj,
142 real *dvdlambda_q, real *dvdlambda_lj)
144 int i, i1, i2, j, k, m, iv, jv, q;
146 double Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
149 real v, vc, qiA, qiB, dr2, rinv, enercorr;
150 real Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ewaldcoeff_q, ewcdr;
151 real ewc_lj = fr->ewaldcoeff_lj, ewc_lj2 = ewc_lj * ewc_lj;
152 real c6Ai = 0, c6Bi = 0, c6A = 0, c6B = 0, ewcdr2, ewcdr4, c6L = 0, rinv6;
153 rvec df, dx, mutot[2], dipcorrA, dipcorrB;
154 tensor dxdf_q, dxdf_lj;
155 real vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
156 real L1_q, L1_lj, dipole_coeff, qqA, qqB, qqL, vr0_q, vr0_lj = 0;
157 gmx_bool bFreeEnergy = (chargeB != NULL);
158 gmx_bool bMolPBC = fr->bMolPBC;
159 gmx_bool bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
161 /* This routine can be made faster by using tables instead of analytical interactions
162 * However, that requires a thorough verification that they are correct in all cases.
165 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
166 vr0_q = ewc_q*M_2_SQRTPI;
167 if (EVDW_PME(fr->vdwtype))
169 vr0_lj = -pow(ewc_lj, 6)/6.0;
180 L1_lj = 1.0-lambda_lj;
181 /* Note that we have to transform back to gromacs units, since
182 * mu_tot contains the dipole in debye units (for output).
184 for (i = 0; (i < DIM); i++)
186 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
187 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
192 switch (ewald_geometry)
195 if (epsilon_surface != 0)
198 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
199 for (i = 0; (i < DIM); i++)
201 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
202 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
207 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
208 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
209 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
212 gmx_incons("Unsupported Ewald geometry");
217 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n",
218 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
219 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n",
220 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
223 if (EVDW_PME(fr->vdwtype))
227 if ((calc_excl_corr || dipole_coeff != 0) && !bFreeEnergy)
229 for (i = start; (i < end); i++)
231 /* Initiate local variables (for this i-particle) to 0 */
232 qiA = chargeA[i]*one_4pi_eps;
233 if (EVDW_PME(fr->vdwtype))
244 i2 = excl->index[i+1];
246 /* Loop over excluded neighbours */
247 for (j = i1; (j < i2); j++)
251 * First we must test whether k <> i, and then,
252 * because the exclusions are all listed twice i->k
253 * and k->i we must select just one of the two. As
254 * a minor optimization we only compute forces when
255 * the charges are non-zero.
259 qqA = qiA*chargeA[k];
260 if (EVDW_PME(fr->vdwtype))
265 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
268 if (qqA != 0.0 || c6A != 0.0)
273 rvec_sub(x[i], x[k], dx);
276 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
277 for (m = DIM-1; (m >= 0); m--)
279 if (dx[m] > 0.5*box[m][m])
281 rvec_dec(dx, box[m]);
283 else if (dx[m] < -0.5*box[m][m])
285 rvec_inc(dx, box[m]);
290 /* Distance between two excluded particles
291 * may be zero in the case of shells
295 rinv = gmx_invsqrt(dr2);
303 vc = qqA*gmx_erf(ewcdr)*rinv;
306 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
307 #define R_ERF_R_INACC 0.006
309 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
310 #define R_ERF_R_INACC 0.1
312 /* fscal is the scalar force pre-multiplied by rinv,
313 * to normalise the relative position vector dx */
314 if (ewcdr > R_ERF_R_INACC)
316 fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*exp(-ewcdr*ewcdr));
320 /* Use a fourth order series expansion for small ewcdr */
321 fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
324 /* The force vector is obtained by multiplication with
325 * the relative position vector
327 svmul(fscal, dx, df);
330 for (iv = 0; (iv < DIM); iv++)
332 for (jv = 0; (jv < DIM); jv++)
334 dxdf_q[iv][jv] += dx[iv]*df[jv];
341 rinv6 = rinv2*rinv2*rinv2;
342 ewcdr2 = ewc_lj2*dr2;
343 ewcdr4 = ewcdr2*ewcdr2;
344 /* We get the excluded long-range contribution from -C6*(1-g(r))
345 * g(r) is also defined in the manual under LJ-PME
347 vc = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
349 /* The force is the derivative of the potential vc.
350 * fscal is the scalar force pre-multiplied by rinv,
351 * to normalise the relative position vector dx */
352 fscal = 6.0*vc*rinv2 + c6A*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
354 /* The force vector is obtained by multiplication with
355 * the relative position vector
357 svmul(fscal, dx, df);
360 for (iv = 0; (iv < DIM); iv++)
362 for (jv = 0; (jv < DIM); jv++)
364 dxdf_lj[iv][jv] += dx[iv]*df[jv];
371 Vexcl_q += qqA*vr0_q;
372 Vexcl_lj += c6A*vr0_lj;
378 /* Dipole correction on force */
379 if (dipole_coeff != 0)
381 for (j = 0; (j < DIM); j++)
383 f[i][j] -= dipcorrA[j]*chargeA[i];
388 else if (calc_excl_corr || dipole_coeff != 0)
390 for (i = start; (i < end); i++)
392 /* Initiate local variables (for this i-particle) to 0 */
393 qiA = chargeA[i]*one_4pi_eps;
394 qiB = chargeB[i]*one_4pi_eps;
395 if (EVDW_PME(fr->vdwtype))
408 i2 = excl->index[i+1];
410 /* Loop over excluded neighbours */
411 for (j = i1; (j < i2); j++)
416 qqA = qiA*chargeA[k];
417 qqB = qiB*chargeB[k];
418 if (EVDW_PME(fr->vdwtype))
424 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
425 c6B *= pow(0.5*(sigmaB[i]+sigmaB[k]), 6)*sigma3B[k];
428 if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
433 qqL = L1_q*qqA + lambda_q*qqB;
434 if (EVDW_PME(fr->vdwtype))
436 c6L = L1_lj*c6A + lambda_lj*c6B;
438 rvec_sub(x[i], x[k], dx);
441 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
442 for (m = DIM-1; (m >= 0); m--)
444 if (dx[m] > 0.5*box[m][m])
446 rvec_dec(dx, box[m]);
448 else if (dx[m] < -0.5*box[m][m])
450 rvec_inc(dx, box[m]);
457 rinv = gmx_invsqrt(dr2);
459 if (qqA != 0.0 || qqB != 0.0)
464 v = gmx_erf(ewc_q*dr)*rinv;
467 /* fscal is the scalar force pre-multiplied by rinv,
468 * to normalise the relative position vector dx */
469 fscal = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*exp(-ewc_q*ewc_q*dr2));
470 dvdl_excl_q += (qqB - qqA)*v;
472 /* The force vector is obtained by multiplication with
473 * the relative position vector
475 svmul(fscal, dx, df);
478 for (iv = 0; (iv < DIM); iv++)
480 for (jv = 0; (jv < DIM); jv++)
482 dxdf_q[iv][jv] += dx[iv]*df[jv];
487 if ((c6A != 0.0 || c6B != 0.0) && EVDW_PME(fr->vdwtype))
489 rinv6 = rinv2*rinv2*rinv2;
490 ewcdr2 = ewc_lj2*dr2;
491 ewcdr4 = ewcdr2*ewcdr2;
492 v = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
495 /* fscal is the scalar force pre-multiplied by rinv,
496 * to normalise the relative position vector dx */
497 fscal = 6.0*vc*rinv2 + c6L*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
498 dvdl_excl_lj += (c6B - c6A)*v;
500 /* The force vector is obtained by multiplication with
501 * the relative position vector
503 svmul(fscal, dx, df);
506 for (iv = 0; (iv < DIM); iv++)
508 for (jv = 0; (jv < DIM); jv++)
510 dxdf_lj[iv][jv] += dx[iv]*df[jv];
517 Vexcl_q += qqL*vr0_q;
518 dvdl_excl_q += (qqB - qqA)*vr0_q;
519 Vexcl_lj += c6L*vr0_lj;
520 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
526 /* Dipole correction on force */
527 if (dipole_coeff != 0)
529 for (j = 0; (j < DIM); j++)
531 f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
532 + lambda_q*dipcorrB[j]*chargeB[i];
537 for (iv = 0; (iv < DIM); iv++)
539 for (jv = 0; (jv < DIM); jv++)
541 vir_q[iv][jv] += 0.5*dxdf_q[iv][jv];
542 vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
551 /* Global corrections only on master process */
552 if (MASTER(cr) && thread == 0)
554 for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
558 /* Self-energy correction */
559 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
560 if (EVDW_PME(fr->vdwtype))
562 Vself_lj[q] = fr->c6sum[q]*0.5*vr0_lj;
566 /* Apply surface dipole correction:
567 * correction = dipole_coeff * (dipole)^2
569 if (dipole_coeff != 0)
571 if (ewald_geometry == eewg3D)
573 Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
575 else if (ewald_geometry == eewg3DC)
577 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
584 *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
585 if (EVDW_PME(fr->vdwtype))
587 *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
592 *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
593 + lambda_q*(Vdipole[1] - Vself_q[1])
595 *dvdlambda_q += Vdipole[1] - Vself_q[1]
596 - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
597 if (EVDW_PME(fr->vdwtype))
599 *Vcorr_lj = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
600 *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
606 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
607 fprintf(debug, "start=%d,natoms=%d\n", start, end-start);
608 fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
609 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]);
610 fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
611 fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
612 if (MASTER(cr) && thread == 0)
614 if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
616 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
617 L1_q*Vdipole[0]+lambda_q*Vdipole[1]);
623 real ewald_charge_correction(t_commrec *cr, t_forcerec *fr, real lambda,
625 real *dvdlambda, tensor vir)
628 real vol, fac, qs2A, qs2B, vc, enercorr;
633 /* Apply charge correction */
634 vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
636 fac = M_PI*ONE_4PI_EPS0/(fr->epsilon_r*2.0*vol*vol*sqr(fr->ewaldcoeff_q));
638 qs2A = fr->qsum[0]*fr->qsum[0];
639 qs2B = fr->qsum[1]*fr->qsum[1];
641 vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
645 *dvdlambda += -vol*(qs2B - qs2A)*fac;
647 for (d = 0; d < DIM; d++)
654 fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);