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.
39 #include "long-range-correction.h"
43 #include "gromacs/legacyheaders/macros.h"
44 #include "gromacs/legacyheaders/names.h"
45 #include "gromacs/legacyheaders/types/commrec.h"
46 #include "gromacs/legacyheaders/types/forcerec.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/math/utilities.h"
49 #include "gromacs/math/vec.h"
51 /* There's nothing special to do here if just masses are perturbed,
52 * but if either charge or type is perturbed then the implementation
53 * requires that B states are defined for both charge and type, and
54 * does not optimize for the cases where only one changes.
56 * The parameter vectors for B states are left undefined in atoms2md()
57 * when either FEP is inactive, or when there are no mass/charge/type
58 * perturbations. The parameter vectors for LJ-PME are likewise
59 * undefined when LJ-PME is not active. This works because
60 * bHaveChargeOrTypePerturbed handles the control flow. */
61 void ewald_LRcorrection(int start, int end,
62 t_commrec *cr, int thread, t_forcerec *fr,
63 real *chargeA, real *chargeB,
65 real *sigmaA, real *sigmaB,
66 real *sigma3A, real *sigma3B,
67 gmx_bool bHaveChargeOrTypePerturbed,
68 gmx_bool calc_excl_corr,
69 t_blocka *excl, rvec x[],
70 matrix box, rvec mu_tot[],
71 int ewald_geometry, real epsilon_surface,
72 rvec *f, tensor vir_q, tensor vir_lj,
73 real *Vcorr_q, real *Vcorr_lj,
74 real lambda_q, real lambda_lj,
75 real *dvdlambda_q, real *dvdlambda_lj)
77 int i, i1, i2, j, k, m, iv, jv, q;
79 double Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
82 real v, vc, qiA, qiB, dr2, rinv, enercorr;
83 real Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ewaldcoeff_q, ewcdr;
84 real ewc_lj = fr->ewaldcoeff_lj, ewc_lj2 = ewc_lj * ewc_lj;
85 real c6Ai = 0, c6Bi = 0, c6A = 0, c6B = 0, ewcdr2, ewcdr4, c6L = 0, rinv6;
86 rvec df, dx, mutot[2], dipcorrA, dipcorrB;
87 tensor dxdf_q, dxdf_lj;
88 real vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
89 real L1_q, L1_lj, dipole_coeff, qqA, qqB, qqL, vr0_q, vr0_lj = 0;
90 gmx_bool bMolPBC = fr->bMolPBC;
91 gmx_bool bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
93 /* This routine can be made faster by using tables instead of analytical interactions
94 * However, that requires a thorough verification that they are correct in all cases.
97 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
98 vr0_q = ewc_q*M_2_SQRTPI;
99 if (EVDW_PME(fr->vdwtype))
101 vr0_lj = -pow(ewc_lj, 6)/6.0;
112 L1_lj = 1.0-lambda_lj;
113 /* Note that we have to transform back to gromacs units, since
114 * mu_tot contains the dipole in debye units (for output).
116 for (i = 0; (i < DIM); i++)
118 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
119 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
124 switch (ewald_geometry)
127 if (epsilon_surface != 0)
130 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
131 for (i = 0; (i < DIM); i++)
133 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
134 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
139 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
140 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
141 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
144 gmx_incons("Unsupported Ewald geometry");
149 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n",
150 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
151 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n",
152 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
155 if (EVDW_PME(fr->vdwtype))
159 if ((calc_excl_corr || dipole_coeff != 0) && !bHaveChargeOrTypePerturbed)
161 for (i = start; (i < end); i++)
163 /* Initiate local variables (for this i-particle) to 0 */
164 qiA = chargeA[i]*one_4pi_eps;
165 if (EVDW_PME(fr->vdwtype))
176 i2 = excl->index[i+1];
178 /* Loop over excluded neighbours */
179 for (j = i1; (j < i2); j++)
183 * First we must test whether k <> i, and then,
184 * because the exclusions are all listed twice i->k
185 * and k->i we must select just one of the two. As
186 * a minor optimization we only compute forces when
187 * the charges are non-zero.
191 qqA = qiA*chargeA[k];
192 if (EVDW_PME(fr->vdwtype))
197 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
200 if (qqA != 0.0 || c6A != 0.0)
205 rvec_sub(x[i], x[k], dx);
208 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
209 for (m = DIM-1; (m >= 0); m--)
211 if (dx[m] > 0.5*box[m][m])
213 rvec_dec(dx, box[m]);
215 else if (dx[m] < -0.5*box[m][m])
217 rvec_inc(dx, box[m]);
222 /* Distance between two excluded particles
223 * may be zero in the case of shells
227 rinv = gmx_invsqrt(dr2);
235 vc = qqA*gmx_erf(ewcdr)*rinv;
238 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
239 #define R_ERF_R_INACC 0.006
241 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
242 #define R_ERF_R_INACC 0.1
244 /* fscal is the scalar force pre-multiplied by rinv,
245 * to normalise the relative position vector dx */
246 if (ewcdr > R_ERF_R_INACC)
248 fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*exp(-ewcdr*ewcdr));
252 /* Use a fourth order series expansion for small ewcdr */
253 fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
256 /* The force vector is obtained by multiplication with
257 * the relative position vector
259 svmul(fscal, dx, df);
262 for (iv = 0; (iv < DIM); iv++)
264 for (jv = 0; (jv < DIM); jv++)
266 dxdf_q[iv][jv] += dx[iv]*df[jv];
273 rinv6 = rinv2*rinv2*rinv2;
274 ewcdr2 = ewc_lj2*dr2;
275 ewcdr4 = ewcdr2*ewcdr2;
276 /* We get the excluded long-range contribution from -C6*(1-g(r))
277 * g(r) is also defined in the manual under LJ-PME
279 vc = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
281 /* The force is the derivative of the potential vc.
282 * fscal is the scalar force pre-multiplied by rinv,
283 * to normalise the relative position vector dx */
284 fscal = 6.0*vc*rinv2 + c6A*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
286 /* The force vector is obtained by multiplication with
287 * the relative position vector
289 svmul(fscal, dx, df);
292 for (iv = 0; (iv < DIM); iv++)
294 for (jv = 0; (jv < DIM); jv++)
296 dxdf_lj[iv][jv] += dx[iv]*df[jv];
303 Vexcl_q += qqA*vr0_q;
304 Vexcl_lj += c6A*vr0_lj;
310 /* Dipole correction on force */
311 if (dipole_coeff != 0)
313 for (j = 0; (j < DIM); j++)
315 f[i][j] -= dipcorrA[j]*chargeA[i];
320 else if (calc_excl_corr || dipole_coeff != 0)
322 for (i = start; (i < end); i++)
324 /* Initiate local variables (for this i-particle) to 0 */
325 qiA = chargeA[i]*one_4pi_eps;
326 qiB = chargeB[i]*one_4pi_eps;
327 if (EVDW_PME(fr->vdwtype))
340 i2 = excl->index[i+1];
342 /* Loop over excluded neighbours */
343 for (j = i1; (j < i2); j++)
348 qqA = qiA*chargeA[k];
349 qqB = qiB*chargeB[k];
350 if (EVDW_PME(fr->vdwtype))
356 c6A *= pow(0.5*(sigmaA[i]+sigmaA[k]), 6)*sigma3A[k];
357 c6B *= pow(0.5*(sigmaB[i]+sigmaB[k]), 6)*sigma3B[k];
360 if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
365 qqL = L1_q*qqA + lambda_q*qqB;
366 if (EVDW_PME(fr->vdwtype))
368 c6L = L1_lj*c6A + lambda_lj*c6B;
370 rvec_sub(x[i], x[k], dx);
373 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
374 for (m = DIM-1; (m >= 0); m--)
376 if (dx[m] > 0.5*box[m][m])
378 rvec_dec(dx, box[m]);
380 else if (dx[m] < -0.5*box[m][m])
382 rvec_inc(dx, box[m]);
389 rinv = gmx_invsqrt(dr2);
391 if (qqA != 0.0 || qqB != 0.0)
396 v = gmx_erf(ewc_q*dr)*rinv;
399 /* fscal is the scalar force pre-multiplied by rinv,
400 * to normalise the relative position vector dx */
401 fscal = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*exp(-ewc_q*ewc_q*dr2));
402 dvdl_excl_q += (qqB - qqA)*v;
404 /* The force vector is obtained by multiplication with
405 * the relative position vector
407 svmul(fscal, dx, df);
410 for (iv = 0; (iv < DIM); iv++)
412 for (jv = 0; (jv < DIM); jv++)
414 dxdf_q[iv][jv] += dx[iv]*df[jv];
419 if ((c6A != 0.0 || c6B != 0.0) && EVDW_PME(fr->vdwtype))
421 rinv6 = rinv2*rinv2*rinv2;
422 ewcdr2 = ewc_lj2*dr2;
423 ewcdr4 = ewcdr2*ewcdr2;
424 v = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
427 /* fscal is the scalar force pre-multiplied by rinv,
428 * to normalise the relative position vector dx */
429 fscal = 6.0*vc*rinv2 + c6L*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
430 dvdl_excl_lj += (c6B - c6A)*v;
432 /* The force vector is obtained by multiplication with
433 * the relative position vector
435 svmul(fscal, dx, df);
438 for (iv = 0; (iv < DIM); iv++)
440 for (jv = 0; (jv < DIM); jv++)
442 dxdf_lj[iv][jv] += dx[iv]*df[jv];
449 Vexcl_q += qqL*vr0_q;
450 dvdl_excl_q += (qqB - qqA)*vr0_q;
451 Vexcl_lj += c6L*vr0_lj;
452 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
458 /* Dipole correction on force */
459 if (dipole_coeff != 0)
461 for (j = 0; (j < DIM); j++)
463 f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
464 + lambda_q*dipcorrB[j]*chargeB[i];
469 for (iv = 0; (iv < DIM); iv++)
471 for (jv = 0; (jv < DIM); jv++)
473 vir_q[iv][jv] += 0.5*dxdf_q[iv][jv];
474 vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
483 /* Global corrections only on master process */
484 if (MASTER(cr) && thread == 0)
486 for (q = 0; q < (bHaveChargeOrTypePerturbed ? 2 : 1); q++)
490 /* Self-energy correction */
491 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
492 if (EVDW_PME(fr->vdwtype))
494 Vself_lj[q] = fr->c6sum[q]*0.5*vr0_lj;
498 /* Apply surface dipole correction:
499 * correction = dipole_coeff * (dipole)^2
501 if (dipole_coeff != 0)
503 if (ewald_geometry == eewg3D)
505 Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
507 else if (ewald_geometry == eewg3DC)
509 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
514 if (!bHaveChargeOrTypePerturbed)
516 *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
517 if (EVDW_PME(fr->vdwtype))
519 *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
524 *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
525 + lambda_q*(Vdipole[1] - Vself_q[1])
527 *dvdlambda_q += Vdipole[1] - Vself_q[1]
528 - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
529 if (EVDW_PME(fr->vdwtype))
531 *Vcorr_lj = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
532 *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
538 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
539 fprintf(debug, "start=%d,natoms=%d\n", start, end-start);
540 fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
541 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]);
542 fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
543 fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
544 if (MASTER(cr) && thread == 0)
546 if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
548 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
549 L1_q*Vdipole[0]+lambda_q*Vdipole[1]);