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,2015, 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/math/functions.h"
44 #include "gromacs/math/units.h"
45 #include "gromacs/math/utilities.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/mdtypes/commrec.h"
48 #include "gromacs/mdtypes/forcerec.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/utility/fatalerror.h"
52 /* There's nothing special to do here if just masses are perturbed,
53 * but if either charge or type is perturbed then the implementation
54 * requires that B states are defined for both charge and type, and
55 * does not optimize for the cases where only one changes.
57 * The parameter vectors for B states are left undefined in atoms2md()
58 * when either FEP is inactive, or when there are no mass/charge/type
59 * perturbations. The parameter vectors for LJ-PME are likewise
60 * undefined when LJ-PME is not active. This works because
61 * bHaveChargeOrTypePerturbed handles the control flow. */
62 void ewald_LRcorrection(int start, int end,
63 t_commrec *cr, int thread, t_forcerec *fr,
64 real *chargeA, real *chargeB,
66 real *sigmaA, real *sigmaB,
67 real *sigma3A, real *sigma3B,
68 gmx_bool bHaveChargeOrTypePerturbed,
69 gmx_bool calc_excl_corr,
70 t_blocka *excl, rvec x[],
71 matrix box, rvec mu_tot[],
72 int ewald_geometry, real epsilon_surface,
73 rvec *f, tensor vir_q, tensor vir_lj,
74 real *Vcorr_q, real *Vcorr_lj,
75 real lambda_q, real lambda_lj,
76 real *dvdlambda_q, real *dvdlambda_lj)
78 int i, i1, i2, j, k, m, iv, jv, q;
80 double Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
83 real v, vc, qiA, qiB, dr2, rinv;
84 real Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ewaldcoeff_q, ewcdr;
85 real ewc_lj = fr->ewaldcoeff_lj, ewc_lj2 = ewc_lj * ewc_lj;
86 real c6Ai = 0, c6Bi = 0, c6A = 0, c6B = 0, ewcdr2, ewcdr4, c6L = 0, rinv6;
87 rvec df, dx, mutot[2], dipcorrA, dipcorrB;
88 tensor dxdf_q = {{0}}, dxdf_lj = {{0}};
89 real vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
90 real L1_q, L1_lj, dipole_coeff, qqA, qqB, qqL, vr0_q, vr0_lj = 0;
91 gmx_bool bMolPBC = fr->bMolPBC;
92 gmx_bool bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
93 gmx_bool bNeedLongRangeCorrection;
95 /* This routine can be made faster by using tables instead of analytical interactions
96 * However, that requires a thorough verification that they are correct in all cases.
99 one_4pi_eps = ONE_4PI_EPS0/fr->epsilon_r;
100 vr0_q = ewc_q*M_2_SQRTPI;
101 if (EVDW_PME(fr->vdwtype))
103 vr0_lj = -gmx::power6(ewc_lj)/6.0;
114 L1_lj = 1.0-lambda_lj;
115 /* Note that we have to transform back to gromacs units, since
116 * mu_tot contains the dipole in debye units (for output).
118 for (i = 0; (i < DIM); i++)
120 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
121 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
126 switch (ewald_geometry)
129 if (epsilon_surface != 0)
132 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->epsilon_r)*vol);
133 for (i = 0; (i < DIM); i++)
135 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
136 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
141 dipole_coeff = 2*M_PI*one_4pi_eps/vol;
142 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
143 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
146 gmx_incons("Unsupported Ewald geometry");
151 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n",
152 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
153 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n",
154 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
156 bNeedLongRangeCorrection = (calc_excl_corr || dipole_coeff != 0);
157 if (bNeedLongRangeCorrection && !bHaveChargeOrTypePerturbed)
159 for (i = start; (i < end); i++)
161 /* Initiate local variables (for this i-particle) to 0 */
162 qiA = chargeA[i]*one_4pi_eps;
163 if (EVDW_PME(fr->vdwtype))
174 i2 = excl->index[i+1];
176 /* Loop over excluded neighbours */
177 for (j = i1; (j < i2); j++)
181 * First we must test whether k <> i, and then,
182 * because the exclusions are all listed twice i->k
183 * and k->i we must select just one of the two. As
184 * a minor optimization we only compute forces when
185 * the charges are non-zero.
189 qqA = qiA*chargeA[k];
190 if (EVDW_PME(fr->vdwtype))
195 c6A *= gmx::power6(0.5*(sigmaA[i]+sigmaA[k]))*sigma3A[k];
198 if (qqA != 0.0 || c6A != 0.0)
200 rvec_sub(x[i], x[k], dx);
203 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
204 for (m = DIM-1; (m >= 0); m--)
206 if (dx[m] > 0.5*box[m][m])
208 rvec_dec(dx, box[m]);
210 else if (dx[m] < -0.5*box[m][m])
212 rvec_inc(dx, box[m]);
217 /* Distance between two excluded particles
218 * may be zero in the case of shells
222 rinv = gmx::invsqrt(dr2);
230 vc = qqA*std::erf(ewcdr)*rinv;
233 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
234 #define R_ERF_R_INACC 0.006
236 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
237 #define R_ERF_R_INACC 0.1
239 /* fscal is the scalar force pre-multiplied by rinv,
240 * to normalise the relative position vector dx */
241 if (ewcdr > R_ERF_R_INACC)
243 fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*exp(-ewcdr*ewcdr));
247 /* Use a fourth order series expansion for small ewcdr */
248 fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
251 /* The force vector is obtained by multiplication with
252 * the relative position vector
254 svmul(fscal, dx, df);
257 for (iv = 0; (iv < DIM); iv++)
259 for (jv = 0; (jv < DIM); jv++)
261 dxdf_q[iv][jv] += dx[iv]*df[jv];
270 rinv6 = rinv2*rinv2*rinv2;
271 ewcdr2 = ewc_lj2*dr2;
272 ewcdr4 = ewcdr2*ewcdr2;
273 /* We get the excluded long-range contribution from -C6*(1-g(r))
274 * g(r) is also defined in the manual under LJ-PME
276 vc = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
278 /* The force is the derivative of the potential vc.
279 * fscal is the scalar force pre-multiplied by rinv,
280 * to normalise the relative position vector dx */
281 fscal = 6.0*vc*rinv2 + c6A*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
283 /* The force vector is obtained by multiplication with
284 * the relative position vector
286 svmul(fscal, dx, df);
289 for (iv = 0; (iv < DIM); iv++)
291 for (jv = 0; (jv < DIM); jv++)
293 dxdf_lj[iv][jv] += dx[iv]*df[jv];
300 Vexcl_q += qqA*vr0_q;
301 Vexcl_lj += c6A*vr0_lj;
307 /* Dipole correction on force */
308 if (dipole_coeff != 0)
310 for (j = 0; (j < DIM); j++)
312 f[i][j] -= dipcorrA[j]*chargeA[i];
317 else if (bNeedLongRangeCorrection)
319 for (i = start; (i < end); i++)
321 /* Initiate local variables (for this i-particle) to 0 */
322 qiA = chargeA[i]*one_4pi_eps;
323 qiB = chargeB[i]*one_4pi_eps;
324 if (EVDW_PME(fr->vdwtype))
337 i2 = excl->index[i+1];
339 /* Loop over excluded neighbours */
340 for (j = i1; (j < i2); j++)
345 qqA = qiA*chargeA[k];
346 qqB = qiB*chargeB[k];
347 if (EVDW_PME(fr->vdwtype))
353 c6A *= gmx::power6(0.5*(sigmaA[i]+sigmaA[k]))*sigma3A[k];
354 c6B *= gmx::power6(0.5*(sigmaB[i]+sigmaB[k]))*sigma3B[k];
357 if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
361 qqL = L1_q*qqA + lambda_q*qqB;
362 if (EVDW_PME(fr->vdwtype))
364 c6L = L1_lj*c6A + lambda_lj*c6B;
366 rvec_sub(x[i], x[k], dx);
369 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
370 for (m = DIM-1; (m >= 0); m--)
372 if (dx[m] > 0.5*box[m][m])
374 rvec_dec(dx, box[m]);
376 else if (dx[m] < -0.5*box[m][m])
378 rvec_inc(dx, box[m]);
385 rinv = gmx::invsqrt(dr2);
387 if (qqA != 0.0 || qqB != 0.0)
392 v = std::erf(ewc_q*dr)*rinv;
395 /* fscal is the scalar force pre-multiplied by rinv,
396 * to normalise the relative position vector dx */
397 fscal = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*exp(-ewc_q*ewc_q*dr2));
398 dvdl_excl_q += (qqB - qqA)*v;
400 /* The force vector is obtained by multiplication with
401 * the relative position vector
403 svmul(fscal, dx, df);
406 for (iv = 0; (iv < DIM); iv++)
408 for (jv = 0; (jv < DIM); jv++)
410 dxdf_q[iv][jv] += dx[iv]*df[jv];
415 if ((c6A != 0.0 || c6B != 0.0) && EVDW_PME(fr->vdwtype))
417 rinv6 = rinv2*rinv2*rinv2;
418 ewcdr2 = ewc_lj2*dr2;
419 ewcdr4 = ewcdr2*ewcdr2;
420 v = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
423 /* fscal is the scalar force pre-multiplied by rinv,
424 * to normalise the relative position vector dx */
425 fscal = 6.0*vc*rinv2 + c6L*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
426 dvdl_excl_lj += (c6B - c6A)*v;
428 /* The force vector is obtained by multiplication with
429 * the relative position vector
431 svmul(fscal, dx, df);
434 for (iv = 0; (iv < DIM); iv++)
436 for (jv = 0; (jv < DIM); jv++)
438 dxdf_lj[iv][jv] += dx[iv]*df[jv];
445 Vexcl_q += qqL*vr0_q;
446 dvdl_excl_q += (qqB - qqA)*vr0_q;
447 Vexcl_lj += c6L*vr0_lj;
448 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
454 /* Dipole correction on force */
455 if (dipole_coeff != 0)
457 for (j = 0; (j < DIM); j++)
459 f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
460 + lambda_q*dipcorrB[j]*chargeB[i];
465 for (iv = 0; (iv < DIM); iv++)
467 for (jv = 0; (jv < DIM); jv++)
469 vir_q[iv][jv] += 0.5*dxdf_q[iv][jv];
470 vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
479 /* Global corrections only on master process */
480 if (MASTER(cr) && thread == 0)
482 for (q = 0; q < (bHaveChargeOrTypePerturbed ? 2 : 1); q++)
486 /* Self-energy correction */
487 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
488 if (EVDW_PME(fr->vdwtype))
490 Vself_lj[q] = fr->c6sum[q]*0.5*vr0_lj;
494 /* Apply surface dipole correction:
495 * correction = dipole_coeff * (dipole)^2
497 if (dipole_coeff != 0)
499 if (ewald_geometry == eewg3D)
501 Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
503 else if (ewald_geometry == eewg3DC)
505 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
510 if (!bHaveChargeOrTypePerturbed)
512 *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
513 if (EVDW_PME(fr->vdwtype))
515 *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
520 *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
521 + lambda_q*(Vdipole[1] - Vself_q[1])
523 *dvdlambda_q += Vdipole[1] - Vself_q[1]
524 - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
525 if (EVDW_PME(fr->vdwtype))
527 *Vcorr_lj = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
528 *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
534 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
535 fprintf(debug, "start=%d,natoms=%d\n", start, end-start);
536 fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
537 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]);
538 fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
539 fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
540 if (MASTER(cr) && thread == 0)
542 if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
544 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
545 L1_q*Vdipole[0]+lambda_q*Vdipole[1]);