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,2017,2018 by the GROMACS development team.
7 * Copyright (c) 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 * \brief This file contains function definitions necessary for
41 * computing energies and forces for the plain-Ewald long-ranged part,
42 * and the correction for overall system charge for all Ewald-family
45 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
46 * \author Mark Abraham <mark.j.abraham@gmail.com>
47 * \ingroup module_ewald
59 #include "gromacs/ewald/ewald_utils.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/gmxcomplex.h"
62 #include "gromacs/math/units.h"
63 #include "gromacs/math/utilities.h"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/math/vectypes.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/forcerec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/interaction_const.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/smalloc.h"
74 using cvec = std::array<t_complex, DIM>;
76 gmx_ewald_tab_t::gmx_ewald_tab_t(const t_inputrec& ir, FILE* fp)
80 fprintf(fp, "Will do ordinary reciprocal space Ewald sum.\n");
86 kmax = std::max(nx, std::max(ny, nz));
89 gmx_ewald_tab_t::~gmx_ewald_tab_t() = default;
91 //! Calculates wave vectors.
92 static void calc_lll(const rvec box, rvec lll)
94 lll[XX] = 2.0 * M_PI / box[XX];
95 lll[YY] = 2.0 * M_PI / box[YY];
96 lll[ZZ] = 2.0 * M_PI / box[ZZ];
99 //! Make tables for the structure factor parts
100 static void tabulateStructureFactors(int natom, const rvec x[], int kmax, cvec** eir, const rvec lll)
106 printf("Go away! kmax = %d\n", kmax);
110 for (i = 0; (i < natom); i++)
112 for (m = 0; (m < 3); m++)
118 for (m = 0; (m < 3); m++)
120 eir[1][i][m].re = std::cos(x[i][m] * lll[m]);
121 eir[1][i][m].im = std::sin(x[i][m] * lll[m]);
123 for (j = 2; (j < kmax); j++)
125 for (m = 0; (m < 3); m++)
127 eir[j][i][m] = cmul(eir[j - 1][i][m], eir[1][i][m]);
133 real do_ewald(const t_inputrec& ir,
136 const real chargeA[],
137 const real chargeB[],
147 real factor = -1.0 / (4 * ewaldcoeff * ewaldcoeff);
149 real energy_AB[2], energy;
151 int lowiy, lowiz, ix, iy, iz, n, q;
152 real tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
153 gmx_bool bFreeEnergy;
160 gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
164 /* Scale box with Ewald wall factor */
166 EwaldBoxZScaler boxScaler(ir);
167 boxScaler.scaleBox(box, scaledBox);
170 for (int i = 0; (i < DIM); i++)
172 boxDiag[i] = scaledBox[i][i];
176 real scaleRecip = 4.0 * M_PI / (boxDiag[XX] * boxDiag[YY] * boxDiag[ZZ]) * ONE_4PI_EPS0 / ir.epsilon_r;
179 for (n = 0; n < et->kmax; n++)
181 snew(eir[n], natoms);
183 et->tab_xy.resize(natoms);
184 et->tab_qxyz.resize(natoms);
186 bFreeEnergy = (ir.efep != FreeEnergyPerturbationType::No);
190 calc_lll(boxDiag, lll);
191 tabulateStructureFactors(natoms, x, et->kmax, eir, lll);
193 for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
203 scale = 1.0 - lambda;
213 for (ix = 0; ix < et->nx; ix++)
216 for (iy = lowiy; iy < et->ny; iy++)
221 for (n = 0; n < natoms; n++)
223 et->tab_xy[n] = cmul(eir[ix][n][XX], eir[iy][n][YY]);
228 for (n = 0; n < natoms; n++)
230 et->tab_xy[n] = cmul(eir[ix][n][XX], conjugate(eir[-iy][n][YY]));
233 for (iz = lowiz; iz < et->nz; iz++)
236 m2 = mx * mx + my * my + mz * mz;
237 ak = std::exp(m2 * factor) / m2;
238 akv = 2.0 * ak * (1.0 / m2 - factor);
241 for (n = 0; n < natoms; n++)
243 et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n], eir[iz][n][ZZ]));
248 for (n = 0; n < natoms; n++)
251 rcmul(charge[n], cmul(et->tab_xy[n], conjugate(eir[-iz][n][ZZ])));
256 for (n = 0; n < natoms; n++)
258 cs += et->tab_qxyz[n].re;
259 ss += et->tab_qxyz[n].im;
261 energy_AB[q] += ak * (cs * cs + ss * ss);
262 tmp = scale * akv * (cs * cs + ss * ss);
263 lrvir[XX][XX] -= tmp * mx * mx;
264 lrvir[XX][YY] -= tmp * mx * my;
265 lrvir[XX][ZZ] -= tmp * mx * mz;
266 lrvir[YY][YY] -= tmp * my * my;
267 lrvir[YY][ZZ] -= tmp * my * mz;
268 lrvir[ZZ][ZZ] -= tmp * mz * mz;
269 for (n = 0; n < natoms; n++)
271 /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
272 tmp = scale * ak * (cs * et->tab_qxyz[n].im - ss * et->tab_qxyz[n].re);
273 f[n][XX] += tmp * mx * 2 * scaleRecip;
274 f[n][YY] += tmp * my * 2 * scaleRecip;
275 f[n][ZZ] += tmp * mz * 2 * scaleRecip;
291 energy = energy_AB[0];
295 energy = (1.0 - lambda) * energy_AB[0] + lambda * energy_AB[1];
296 *dvdlambda += scaleRecip * (energy_AB[1] - energy_AB[0]);
299 lrvir[XX][XX] = -0.5 * scaleRecip * (lrvir[XX][XX] + energy);
300 lrvir[XX][YY] = -0.5 * scaleRecip * (lrvir[XX][YY]);
301 lrvir[XX][ZZ] = -0.5 * scaleRecip * (lrvir[XX][ZZ]);
302 lrvir[YY][YY] = -0.5 * scaleRecip * (lrvir[YY][YY] + energy);
303 lrvir[YY][ZZ] = -0.5 * scaleRecip * (lrvir[YY][ZZ]);
304 lrvir[ZZ][ZZ] = -0.5 * scaleRecip * (lrvir[ZZ][ZZ] + energy);
306 lrvir[YY][XX] = lrvir[XX][YY];
307 lrvir[ZZ][XX] = lrvir[XX][ZZ];
308 lrvir[ZZ][YY] = lrvir[YY][ZZ];
310 energy *= scaleRecip;
315 real ewald_charge_correction(const t_commrec* cr,
316 const t_forcerec* fr,
323 real vol, fac, qs2A, qs2B, vc, enercorr;
328 /* Apply charge correction */
329 vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
331 fac = M_PI * ONE_4PI_EPS0
332 / (fr->ic->epsilon_r * 2.0 * vol * vol * gmx::square(fr->ic->ewaldcoeff_q));
334 qs2A = fr->qsum[0] * fr->qsum[0];
335 qs2B = fr->qsum[1] * fr->qsum[1];
337 vc = (qs2A * (1 - lambda) + qs2B * lambda) * fac;
339 enercorr = -vol * vc;
341 *dvdlambda += -vol * (qs2B - qs2A) * fac;
343 for (d = 0; d < DIM; d++)
350 fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);