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,2021, 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/interaction_const.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/gmxassert.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(const int numAtomsLocal,
68 const t_commrec* commrec,
72 gmx::ArrayRef<const double> qsum,
73 EwaldGeometry ewaldGeometry,
74 const real epsilonSurface,
77 gmx::ArrayRef<const real> chargeA,
78 gmx::ArrayRef<const real> chargeB,
79 bool bHaveChargePerturbed,
80 gmx::ArrayRef<const gmx::RVec> coords,
82 gmx::ArrayRef<const gmx::RVec> mu_tot,
83 gmx::ArrayRef<gmx::RVec> forces,
88 /* We need to correct only self interactions */
89 const int start = (numAtomsLocal * thread) / numThreads;
90 const int end = (numAtomsLocal * (thread + 1)) / numThreads;
92 std::array<rvec, 2> mutot;
93 gmx::RVec dipcorrA = { 0, 0, 0 };
94 gmx::RVec dipcorrB = { 0, 0, 0 };
95 std::array<real, 2> chargecorr = { 0 };
96 std::array<real, 2> Vdipole = { 0 };
98 /* Scale the Ewald unit cell when dimension z is not periodic */
100 EwaldBoxZScaler boxScaler(havePbcXY2Walls, wallEwaldZfac);
101 boxScaler.scaleBox(box, scaledBox);
103 const real one_4pi_eps = gmx::c_one4PiEps0 / epsilonR;
104 const real L1_q = 1.0 - lambda_q;
105 /* Note that we have to transform back to gromacs units, since
106 * mu_tot contains the dipole in debye units (for output).
108 for (int i = 0; (i < DIM); i++)
110 mutot[0][i] = mu_tot[0][i] * gmx::c_debye2Enm;
111 mutot[1][i] = mu_tot[1][i] * gmx::c_debye2Enm;
113 real dipole_coeff = 0;
115 real boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
116 switch (ewaldGeometry)
118 case EwaldGeometry::ThreeD:
119 if (epsilonSurface != 0)
122 2 * M_PI * gmx::c_one4PiEps0 / ((2 * epsilonSurface + epsilonR) * boxVolume);
123 for (int i = 0; (i < DIM); i++)
125 dipcorrA[i] = 2 * dipole_coeff * mutot[0][i];
126 dipcorrB[i] = 2 * dipole_coeff * mutot[1][i];
130 case EwaldGeometry::ThreeDC:
131 dipole_coeff = 2 * M_PI * one_4pi_eps / boxVolume;
132 dipcorrA[ZZ] = 2 * dipole_coeff * mutot[0][ZZ];
133 dipcorrB[ZZ] = 2 * dipole_coeff * mutot[1][ZZ];
134 for (int q = 0; q < (bHaveChargePerturbed ? 2 : 1); q++)
136 /* Avoid charge corrections with near-zero net charge */
137 if (fabs(qsum[q]) > 1e-4)
139 chargecorr[q] = 2 * dipole_coeff * qsum[q];
143 default: gmx_incons("Unsupported Ewald geometry");
145 const bool bNeedLongRangeCorrection = (dipole_coeff != 0);
146 if (bNeedLongRangeCorrection && !bHaveChargePerturbed)
148 for (int i = start; (i < end); i++)
150 for (int j = 0; (j < DIM); j++)
152 forces[i][j] -= dipcorrA[j] * chargeA[i];
154 if (chargecorr[0] != 0)
156 forces[i][ZZ] += chargecorr[0] * chargeA[i] * coords[i][ZZ];
160 else if (bNeedLongRangeCorrection)
162 for (int i = start; (i < end); i++)
164 for (int j = 0; (j < DIM); j++)
166 forces[i][j] -= L1_q * dipcorrA[j] * chargeA[i] + lambda_q * dipcorrB[j] * chargeB[i];
168 if (chargecorr[0] != 0 || chargecorr[1] != 0)
171 (L1_q * chargecorr[0] * chargeA[i] + lambda_q * chargecorr[1]) * coords[i][ZZ];
176 /* Global corrections only on master process */
177 if (MASTER(commrec) && thread == 0)
179 for (int q = 0; q < (bHaveChargePerturbed ? 2 : 1); q++)
181 /* Apply surface and charged surface dipole correction:
182 * correction = dipole_coeff * ( (dipole)^2
183 * - qsum*sum_i q_i z_i^2 - qsum^2 * box_z^2 / 12 )
185 if (dipole_coeff != 0)
187 if (ewaldGeometry == EwaldGeometry::ThreeD)
189 Vdipole[q] = dipole_coeff * iprod(mutot[q], mutot[q]);
191 else if (ewaldGeometry == EwaldGeometry::ThreeDC)
193 Vdipole[q] = dipole_coeff * mutot[q][ZZ] * mutot[q][ZZ];
195 if (chargecorr[q] != 0)
197 /* Here we use a non thread-parallelized loop,
198 * because this is the only loop over atoms for
199 * energies and they need reduction (unlike forces).
200 * We could implement a reduction over threads,
201 * but this case is rarely used.
203 gmx::ArrayRef<const real> charge = (q == 0 ? chargeA : chargeB);
205 for (int i = 0; i < numAtomsLocal; i++)
207 sumQZ2 += charge[i] * coords[i][ZZ] * coords[i][ZZ];
209 Vdipole[q] -= dipole_coeff * qsum[q]
210 * (sumQZ2 + qsum[q] * box[ZZ][ZZ] * box[ZZ][ZZ] / 12);
216 if (!bHaveChargePerturbed)
218 *Vcorr_q = Vdipole[0];
222 *Vcorr_q = L1_q * Vdipole[0] + lambda_q * Vdipole[1];
223 *dvdlambda_q += Vdipole[1] - Vdipole[0];