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,2016,2017,2018,2019, 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/ewald/ewald_utils.h"
44 #include "gromacs/math/functions.h"
45 #include "gromacs/math/units.h"
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/mdtypes/forcerec.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/gmxassert.h"
55 #include "pme_internal.h"
57 /* There's nothing special to do here if just masses are perturbed,
58 * but if either charge or type is perturbed then the implementation
59 * requires that B states are defined for both charge and type, and
60 * does not optimize for the cases where only one changes.
62 * The parameter vectors for B states are left undefined in atoms2md()
63 * when either FEP is inactive, or when there are no mass/charge/type
64 * perturbations. The parameter vectors for LJ-PME are likewise
65 * undefined when LJ-PME is not active. This works because
66 * bHaveChargeOrTypePerturbed handles the control flow. */
67 void ewald_LRcorrection(int numAtomsLocal,
69 int numThreads, int thread,
72 const real *chargeA, const real *chargeB,
73 const real *C6A, const real *C6B,
74 const real *sigmaA, const real *sigmaB,
75 const real *sigma3A, const real *sigma3B,
76 gmx_bool bHaveChargeOrTypePerturbed,
77 gmx_bool calc_excl_corr,
80 matrix box, rvec mu_tot[],
81 int ewald_geometry, real epsilon_surface,
82 rvec *f, tensor vir_q, tensor vir_lj,
83 real *Vcorr_q, real *Vcorr_lj,
84 real lambda_q, real lambda_lj,
85 real *dvdlambda_q, real *dvdlambda_lj)
87 int numAtomsToBeCorrected;
90 /* We need to correct all exclusion pairs (cutoff-scheme = group) */
91 numAtomsToBeCorrected = excl->nr;
93 GMX_RELEASE_ASSERT(numAtomsToBeCorrected >= numAtomsLocal, "We might need to do self-corrections");
97 /* We need to correct only self interactions */
98 numAtomsToBeCorrected = numAtomsLocal;
100 int start = (numAtomsToBeCorrected* thread )/numThreads;
101 int end = (numAtomsToBeCorrected*(thread + 1))/numThreads;
103 int i, i1, i2, j, k, m, iv, jv, q;
105 double Vexcl_q, dvdl_excl_q, dvdl_excl_lj; /* Necessary for precision */
108 real v, vc, qiA, qiB, dr2, rinv;
109 real Vself_q[2], Vself_lj[2], Vdipole[2], rinv2, ewc_q = fr->ic->ewaldcoeff_q, ewcdr;
110 real ewc_lj = fr->ic->ewaldcoeff_lj, ewc_lj2 = ewc_lj * ewc_lj;
111 real c6Ai = 0, c6Bi = 0, c6A = 0, c6B = 0, ewcdr2, ewcdr4, c6L = 0, rinv6;
112 rvec df, dx, mutot[2], dipcorrA, dipcorrB;
113 tensor dxdf_q = {{0}}, dxdf_lj = {{0}};
114 real L1_q, L1_lj, dipole_coeff, qqA, qqB, qqL, vr0_q, vr0_lj = 0;
115 real chargecorr[2] = { 0, 0 };
116 gmx_bool bMolPBC = fr->bMolPBC;
117 gmx_bool bDoingLBRule = (fr->ljpme_combination_rule == eljpmeLB);
118 gmx_bool bNeedLongRangeCorrection;
120 GMX_ASSERT(ir, "Invalid inputrec pointer");
122 EwaldBoxZScaler boxScaler(*ir);
123 boxScaler.scaleBox(box, scaledBox);
125 /* This routine can be made faster by using tables instead of analytical interactions
126 * However, that requires a thorough verification that they are correct in all cases.
129 bool vdwPme = EVDW_PME(fr->ic->vdwtype);
131 one_4pi_eps = ONE_4PI_EPS0/fr->ic->epsilon_r;
132 vr0_q = ewc_q*M_2_SQRTPI;
135 vr0_lj = -gmx::power6(ewc_lj)/6.0;
146 L1_lj = 1.0-lambda_lj;
147 /* Note that we have to transform back to gromacs units, since
148 * mu_tot contains the dipole in debye units (for output).
150 for (i = 0; (i < DIM); i++)
152 mutot[0][i] = mu_tot[0][i]*DEBYE2ENM;
153 mutot[1][i] = mu_tot[1][i]*DEBYE2ENM;
159 real boxVolume = scaledBox[XX][XX]*scaledBox[YY][YY]*scaledBox[ZZ][ZZ];
160 switch (ewald_geometry)
163 if (epsilon_surface != 0)
166 2*M_PI*ONE_4PI_EPS0/((2*epsilon_surface + fr->ic->epsilon_r)*boxVolume);
167 for (i = 0; (i < DIM); i++)
169 dipcorrA[i] = 2*dipole_coeff*mutot[0][i];
170 dipcorrB[i] = 2*dipole_coeff*mutot[1][i];
175 dipole_coeff = 2*M_PI*one_4pi_eps/boxVolume;
176 dipcorrA[ZZ] = 2*dipole_coeff*mutot[0][ZZ];
177 dipcorrB[ZZ] = 2*dipole_coeff*mutot[1][ZZ];
178 for (int q = 0; q < (bHaveChargeOrTypePerturbed ? 2 : 1); q++)
180 /* Avoid charge corrections with near-zero net charge */
181 if (fabs(fr->qsum[q]) > 1e-4)
183 chargecorr[q] = 2*dipole_coeff*fr->qsum[q];
188 gmx_incons("Unsupported Ewald geometry");
192 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n",
193 dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
194 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n",
195 mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
197 bNeedLongRangeCorrection = (calc_excl_corr || dipole_coeff != 0);
198 if (bNeedLongRangeCorrection && !bHaveChargeOrTypePerturbed)
200 for (i = start; (i < end); i++)
202 /* Initiate local variables (for this i-particle) to 0 */
203 qiA = chargeA[i]*one_4pi_eps;
215 i2 = excl->index[i+1];
217 /* Loop over excluded neighbours */
218 for (j = i1; (j < i2); j++)
222 * First we must test whether k <> i, and then,
223 * because the exclusions are all listed twice i->k
224 * and k->i we must select just one of the two. As
225 * a minor optimization we only compute forces when
226 * the charges are non-zero.
230 qqA = qiA*chargeA[k];
236 c6A *= gmx::power6(0.5*(sigmaA[i]+sigmaA[k]))*sigma3A[k];
239 if (qqA != 0.0 || c6A != 0.0)
241 rvec_sub(x[i], x[k], dx);
244 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
245 for (m = DIM-1; (m >= 0); m--)
247 if (dx[m] > 0.5*box[m][m])
249 rvec_dec(dx, box[m]);
251 else if (dx[m] < -0.5*box[m][m])
253 rvec_inc(dx, box[m]);
258 /* Distance between two excluded particles
259 * may be zero in the case of shells
263 rinv = gmx::invsqrt(dr2);
271 vc = qqA*std::erf(ewcdr)*rinv;
274 /* Relative accuracy at R_ERF_R_INACC of 3e-10 */
275 #define R_ERF_R_INACC 0.006
277 /* Relative accuracy at R_ERF_R_INACC of 2e-5 */
278 #define R_ERF_R_INACC 0.1
280 /* fscal is the scalar force pre-multiplied by rinv,
281 * to normalise the relative position vector dx */
282 if (ewcdr > R_ERF_R_INACC)
284 fscal = rinv2*(vc - qqA*ewc_q*M_2_SQRTPI*std::exp(-ewcdr*ewcdr));
288 /* Use a fourth order series expansion for small ewcdr */
289 fscal = ewc_q*ewc_q*qqA*vr0_q*(2.0/3.0 - 0.4*ewcdr*ewcdr);
292 /* The force vector is obtained by multiplication with
293 * the relative position vector
295 svmul(fscal, dx, df);
298 for (iv = 0; (iv < DIM); iv++)
300 for (jv = 0; (jv < DIM); jv++)
302 dxdf_q[iv][jv] += dx[iv]*df[jv];
311 rinv6 = rinv2*rinv2*rinv2;
312 ewcdr2 = ewc_lj2*dr2;
313 ewcdr4 = ewcdr2*ewcdr2;
314 /* We get the excluded long-range contribution from -C6*(1-g(r))
315 * g(r) is also defined in the manual under LJ-PME
317 vc = -c6A*rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
319 /* The force is the derivative of the potential vc.
320 * fscal is the scalar force pre-multiplied by rinv,
321 * to normalise the relative position vector dx */
322 fscal = 6.0*vc*rinv2 + c6A*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
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_lj[iv][jv] += dx[iv]*df[jv];
341 Vexcl_q += qqA*vr0_q;
342 Vexcl_lj += c6A*vr0_lj;
348 /* Dipole correction on force */
349 if (dipole_coeff != 0 && i < numAtomsLocal)
351 for (j = 0; (j < DIM); j++)
353 f[i][j] -= dipcorrA[j]*chargeA[i];
355 if (chargecorr[0] != 0)
357 f[i][ZZ] += chargecorr[0]*chargeA[i]*x[i][ZZ];
362 else if (bNeedLongRangeCorrection)
364 for (i = start; (i < end); i++)
366 /* Initiate local variables (for this i-particle) to 0 */
367 qiA = chargeA[i]*one_4pi_eps;
368 qiB = chargeB[i]*one_4pi_eps;
382 i2 = excl->index[i+1];
384 /* Loop over excluded neighbours */
385 for (j = i1; (j < i2); j++)
390 qqA = qiA*chargeA[k];
391 qqB = qiB*chargeB[k];
398 c6A *= gmx::power6(0.5*(sigmaA[i]+sigmaA[k]))*sigma3A[k];
399 c6B *= gmx::power6(0.5*(sigmaB[i]+sigmaB[k]))*sigma3B[k];
402 if (qqA != 0.0 || qqB != 0.0 || c6A != 0.0 || c6B != 0.0)
406 qqL = L1_q*qqA + lambda_q*qqB;
409 c6L = L1_lj*c6A + lambda_lj*c6B;
411 rvec_sub(x[i], x[k], dx);
414 /* Cheap pbc_dx, assume excluded pairs are at short distance. */
415 for (m = DIM-1; (m >= 0); m--)
417 if (dx[m] > 0.5*box[m][m])
419 rvec_dec(dx, box[m]);
421 else if (dx[m] < -0.5*box[m][m])
423 rvec_inc(dx, box[m]);
430 rinv = gmx::invsqrt(dr2);
432 if (qqA != 0.0 || qqB != 0.0)
437 v = std::erf(ewc_q*dr)*rinv;
440 /* fscal is the scalar force pre-multiplied by rinv,
441 * to normalise the relative position vector dx */
442 fscal = rinv2*(vc-qqL*ewc_q*M_2_SQRTPI*std::exp(-ewc_q*ewc_q*dr2));
443 dvdl_excl_q += (qqB - qqA)*v;
445 /* The force vector is obtained by multiplication with
446 * the relative position vector
448 svmul(fscal, dx, df);
451 for (iv = 0; (iv < DIM); iv++)
453 for (jv = 0; (jv < DIM); jv++)
455 dxdf_q[iv][jv] += dx[iv]*df[jv];
460 if ((c6A != 0.0 || c6B != 0.0) && vdwPme)
462 rinv6 = rinv2*rinv2*rinv2;
463 ewcdr2 = ewc_lj2*dr2;
464 ewcdr4 = ewcdr2*ewcdr2;
465 v = -rinv6*(1.0 - exp(-ewcdr2)*(1 + ewcdr2 + 0.5*ewcdr4));
468 /* fscal is the scalar force pre-multiplied by rinv,
469 * to normalise the relative position vector dx */
470 fscal = 6.0*vc*rinv2 + c6L*rinv6*exp(-ewcdr2)*ewc_lj2*ewcdr4;
471 dvdl_excl_lj += (c6B - c6A)*v;
473 /* The force vector is obtained by multiplication with
474 * the relative position vector
476 svmul(fscal, dx, df);
479 for (iv = 0; (iv < DIM); iv++)
481 for (jv = 0; (jv < DIM); jv++)
483 dxdf_lj[iv][jv] += dx[iv]*df[jv];
490 Vexcl_q += qqL*vr0_q;
491 dvdl_excl_q += (qqB - qqA)*vr0_q;
492 Vexcl_lj += c6L*vr0_lj;
493 dvdl_excl_lj += (c6B - c6A)*vr0_lj;
499 /* Dipole correction on force */
500 if (dipole_coeff != 0 && i < numAtomsLocal)
502 for (j = 0; (j < DIM); j++)
504 f[i][j] -= L1_q*dipcorrA[j]*chargeA[i]
505 + lambda_q*dipcorrB[j]*chargeB[i];
507 if (chargecorr[0] != 0 || chargecorr[1] != 0)
509 f[i][ZZ] += (L1_q*chargecorr[0]*chargeA[i]
510 + lambda_q*chargecorr[1])*x[i][ZZ];
515 for (iv = 0; (iv < DIM); iv++)
517 for (jv = 0; (jv < DIM); jv++)
519 vir_q[iv][jv] += 0.5*dxdf_q[iv][jv];
520 vir_lj[iv][jv] += 0.5*dxdf_lj[iv][jv];
529 /* Global corrections only on master process */
530 if (MASTER(cr) && thread == 0)
532 for (q = 0; q < (bHaveChargeOrTypePerturbed ? 2 : 1); q++)
536 /* Self-energy correction */
537 Vself_q[q] = ewc_q*one_4pi_eps*fr->q2sum[q]*M_1_SQRTPI;
540 Vself_lj[q] = fr->c6sum[q]*0.5*vr0_lj;
544 /* Apply surface and charged surface dipole correction:
545 * correction = dipole_coeff * ( (dipole)^2
546 * - qsum*sum_i q_i z_i^2 - qsum^2 * box_z^2 / 12 )
548 if (dipole_coeff != 0)
550 if (ewald_geometry == eewg3D)
552 Vdipole[q] = dipole_coeff*iprod(mutot[q], mutot[q]);
554 else if (ewald_geometry == eewg3DC)
556 Vdipole[q] = dipole_coeff*mutot[q][ZZ]*mutot[q][ZZ];
558 if (chargecorr[q] != 0)
560 /* Here we use a non thread-parallelized loop,
561 * because this is the only loop over atoms for
562 * energies and they need reduction (unlike forces).
563 * We could implement a reduction over threads,
564 * but this case is rarely used.
566 const real *qPtr = (q == 0 ? chargeA : chargeB);
568 for (int i = 0; i < numAtomsLocal; i++)
570 sumQZ2 += qPtr[i]*x[i][ZZ]*x[i][ZZ];
572 Vdipole[q] -= dipole_coeff*fr->qsum[q]*(sumQZ2 + fr->qsum[q]*box[ZZ][ZZ]*box[ZZ][ZZ]/12);
578 if (!bHaveChargeOrTypePerturbed)
580 *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
583 *Vcorr_lj = -Vself_lj[0] - Vexcl_lj;
588 *Vcorr_q = L1_q*(Vdipole[0] - Vself_q[0])
589 + lambda_q*(Vdipole[1] - Vself_q[1])
591 *dvdlambda_q += Vdipole[1] - Vself_q[1]
592 - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
595 *Vcorr_lj = -(L1_lj*Vself_lj[0] + lambda_lj*Vself_lj[1]) - Vexcl_lj;
596 *dvdlambda_lj += -Vself_lj[1] + Vself_lj[0] - dvdl_excl_lj;
602 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
603 fprintf(debug, "q2sum = %g, Vself_q=%g c6sum = %g, Vself_lj=%g\n",
604 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]);
605 fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
606 fprintf(debug, "Lennard-Jones Long Range correction: Vexcl=%g\n", Vexcl_lj);
607 if (MASTER(cr) && thread == 0)
609 if (epsilon_surface > 0 || ewald_geometry == eewg3DC)
611 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
612 L1_q*Vdipole[0]+lambda_q*Vdipole[1]);