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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 #include "long_range_correction.h"
44 #include "gromacs/ewald/ewald_utils.h"
45 #include "gromacs/math/functions.h"
46 #include "gromacs/math/units.h"
47 #include "gromacs/math/utilities.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/mdtypes/forcerec.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/gmxassert.h"
56 #include "pme_internal.h"
58 /* There's nothing special to do here if just masses are perturbed,
59 * but if either charge or type is perturbed then the implementation
60 * requires that B states are defined for both charge and type, and
61 * does not optimize for the cases where only one changes.
63 * The parameter vectors for B states are left undefined in atoms2md()
64 * when either FEP is inactive, or when there are no mass/charge/type
65 * perturbations. The parameter vectors for LJ-PME are likewise
66 * undefined when LJ-PME is not active. This works because
67 * bHaveChargeOrTypePerturbed handles the control flow. */
68 void ewald_LRcorrection(const int numAtomsLocal,
76 gmx_bool bHaveChargePerturbed,
85 /* We need to correct only self interactions */
86 const int start = (numAtomsLocal * thread) / numThreads;
87 const int end = (numAtomsLocal * (thread + 1)) / numThreads;
90 double Vexcl_q, dvdl_excl_q; /* Necessary for precision */
92 real Vself_q[2], Vdipole[2];
93 rvec mutot[2], dipcorrA, dipcorrB;
94 real L1_q, dipole_coeff;
95 real chargecorr[2] = { 0, 0 };
97 /* Scale the Ewald unit cell when dimension z is not periodic */
99 EwaldBoxZScaler boxScaler(ir);
100 boxScaler.scaleBox(box, scaledBox);
102 one_4pi_eps = ONE_4PI_EPS0 / fr.ic->epsilon_r;
107 L1_q = 1.0 - lambda_q;
108 /* Note that we have to transform back to gromacs units, since
109 * mu_tot contains the dipole in debye units (for output).
111 for (i = 0; (i < DIM); i++)
113 mutot[0][i] = mu_tot[0][i] * DEBYE2ENM;
114 mutot[1][i] = mu_tot[1][i] * DEBYE2ENM;
120 real boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
121 switch (ir.ewald_geometry)
124 if (ir.epsilon_surface != 0)
126 dipole_coeff = 2 * M_PI * ONE_4PI_EPS0
127 / ((2 * ir.epsilon_surface + fr.ic->epsilon_r) * boxVolume);
128 for (i = 0; (i < DIM); i++)
130 dipcorrA[i] = 2 * dipole_coeff * mutot[0][i];
131 dipcorrB[i] = 2 * dipole_coeff * mutot[1][i];
136 dipole_coeff = 2 * M_PI * one_4pi_eps / boxVolume;
137 dipcorrA[ZZ] = 2 * dipole_coeff * mutot[0][ZZ];
138 dipcorrB[ZZ] = 2 * dipole_coeff * mutot[1][ZZ];
139 for (int q = 0; q < (bHaveChargePerturbed ? 2 : 1); q++)
141 /* Avoid charge corrections with near-zero net charge */
142 if (fabs(fr.qsum[q]) > 1e-4)
144 chargecorr[q] = 2 * dipole_coeff * fr.qsum[q];
148 default: gmx_incons("Unsupported Ewald geometry");
152 fprintf(debug, "dipcorr = %8.3f %8.3f %8.3f\n", dipcorrA[XX], dipcorrA[YY], dipcorrA[ZZ]);
153 fprintf(debug, "mutot = %8.3f %8.3f %8.3f\n", mutot[0][XX], mutot[0][YY], mutot[0][ZZ]);
155 const bool bNeedLongRangeCorrection = (dipole_coeff != 0);
156 if (bNeedLongRangeCorrection && !bHaveChargePerturbed)
158 for (i = start; (i < end); i++)
160 for (j = 0; (j < DIM); j++)
162 f[i][j] -= dipcorrA[j] * chargeA[i];
164 if (chargecorr[0] != 0)
166 f[i][ZZ] += chargecorr[0] * chargeA[i] * x[i][ZZ];
170 else if (bNeedLongRangeCorrection)
172 for (i = start; (i < end); i++)
174 for (j = 0; (j < DIM); j++)
176 f[i][j] -= L1_q * dipcorrA[j] * chargeA[i] + lambda_q * dipcorrB[j] * chargeB[i];
178 if (chargecorr[0] != 0 || chargecorr[1] != 0)
180 f[i][ZZ] += (L1_q * chargecorr[0] * chargeA[i] + lambda_q * chargecorr[1]) * x[i][ZZ];
188 /* Global corrections only on master process */
189 if (MASTER(cr) && thread == 0)
191 for (q = 0; q < (bHaveChargePerturbed ? 2 : 1); q++)
193 /* Apply surface and charged surface dipole correction:
194 * correction = dipole_coeff * ( (dipole)^2
195 * - qsum*sum_i q_i z_i^2 - qsum^2 * box_z^2 / 12 )
197 if (dipole_coeff != 0)
199 if (ir.ewald_geometry == eewg3D)
201 Vdipole[q] = dipole_coeff * iprod(mutot[q], mutot[q]);
203 else if (ir.ewald_geometry == eewg3DC)
205 Vdipole[q] = dipole_coeff * mutot[q][ZZ] * mutot[q][ZZ];
207 if (chargecorr[q] != 0)
209 /* Here we use a non thread-parallelized loop,
210 * because this is the only loop over atoms for
211 * energies and they need reduction (unlike forces).
212 * We could implement a reduction over threads,
213 * but this case is rarely used.
215 const real* qPtr = (q == 0 ? chargeA : chargeB);
217 for (int i = 0; i < numAtomsLocal; i++)
219 sumQZ2 += qPtr[i] * x[i][ZZ] * x[i][ZZ];
221 Vdipole[q] -= dipole_coeff * fr.qsum[q]
222 * (sumQZ2 + fr.qsum[q] * box[ZZ][ZZ] * box[ZZ][ZZ] / 12);
228 if (!bHaveChargePerturbed)
230 *Vcorr_q = Vdipole[0] - Vself_q[0] - Vexcl_q;
234 *Vcorr_q = L1_q * (Vdipole[0] - Vself_q[0]) + lambda_q * (Vdipole[1] - Vself_q[1]) - Vexcl_q;
235 *dvdlambda_q += Vdipole[1] - Vself_q[1] - (Vdipole[0] - Vself_q[0]) - dvdl_excl_q;
240 fprintf(debug, "Long Range corrections for Ewald interactions:\n");
241 fprintf(debug, "q2sum = %g, Vself_q=%g\n", L1_q * fr.q2sum[0] + lambda_q * fr.q2sum[1],
242 L1_q * Vself_q[0] + lambda_q * Vself_q[1]);
243 fprintf(debug, "Electrostatic Long Range correction: Vexcl=%g\n", Vexcl_q);
244 if (MASTER(cr) && thread == 0)
246 if (ir.epsilon_surface > 0 || ir.ewald_geometry == eewg3DC)
248 fprintf(debug, "Total dipole correction: Vdipole=%g\n",
249 L1_q * Vdipole[0] + lambda_q * Vdipole[1]);