Remove forcerec and inputrec from ewald_LRcorrection
[alexxy/gromacs.git] / src / gromacs / ewald / ewald.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 /*! \internal \file
39  *
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
43  * methods.
44  *
45  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
46  * \author Mark Abraham <mark.j.abraham@gmail.com>
47  * \ingroup module_ewald
48  */
49 #include "gmxpre.h"
50
51 #include "ewald.h"
52
53 #include <cmath>
54 #include <cstdio>
55 #include <cstdlib>
56
57 #include <algorithm>
58
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/arrayref.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/smalloc.h"
74
75 using cvec = std::array<t_complex, DIM>;
76
77 gmx_ewald_tab_t::gmx_ewald_tab_t(const t_inputrec& ir, FILE* fp)
78 {
79     if (fp)
80     {
81         fprintf(fp, "Will do ordinary reciprocal space Ewald sum.\n");
82     }
83
84     nx   = ir.nkx + 1;
85     ny   = ir.nky + 1;
86     nz   = ir.nkz + 1;
87     kmax = std::max(nx, std::max(ny, nz));
88 }
89
90 gmx_ewald_tab_t::~gmx_ewald_tab_t() = default;
91
92 //! Calculates wave vectors.
93 static void calc_lll(const rvec box, rvec lll)
94 {
95     lll[XX] = 2.0 * M_PI / box[XX];
96     lll[YY] = 2.0 * M_PI / box[YY];
97     lll[ZZ] = 2.0 * M_PI / box[ZZ];
98 }
99
100 //! Make tables for the structure factor parts
101 static void tabulateStructureFactors(int natom, gmx::ArrayRef<const gmx::RVec> x, int kmax, cvec** eir, const rvec lll)
102 {
103     int i, j, m;
104
105     if (kmax < 1)
106     {
107         printf("Go away! kmax = %d\n", kmax);
108         exit(1);
109     }
110
111     for (i = 0; (i < natom); i++)
112     {
113         for (m = 0; (m < 3); m++)
114         {
115             eir[0][i][m].re = 1;
116             eir[0][i][m].im = 0;
117         }
118
119         for (m = 0; (m < 3); m++)
120         {
121             eir[1][i][m].re = std::cos(x[i][m] * lll[m]);
122             eir[1][i][m].im = std::sin(x[i][m] * lll[m]);
123         }
124         for (j = 2; (j < kmax); j++)
125         {
126             for (m = 0; (m < 3); m++)
127             {
128                 eir[j][i][m] = cmul(eir[j - 1][i][m], eir[1][i][m]);
129             }
130         }
131     }
132 }
133
134 real do_ewald(const t_inputrec&              ir,
135               gmx::ArrayRef<const gmx::RVec> x,
136               gmx::ArrayRef<gmx::RVec>       f,
137               gmx::ArrayRef<const real>      chargeA,
138               gmx::ArrayRef<const real>      chargeB,
139               const matrix                   box,
140               const t_commrec*               cr,
141               int                            natoms,
142               matrix                         lrvir,
143               real                           ewaldcoeff,
144               real                           lambda,
145               real*                          dvdlambda,
146               gmx_ewald_tab_t*               et)
147 {
148     real   factor = -1.0 / (4 * ewaldcoeff * ewaldcoeff);
149     real   energy_AB[2], energy;
150     rvec   lll;
151     int    lowiy, lowiz, ix, iy, iz, n, q;
152     real   tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
153     cvec** eir;
154     bool   bFreeEnergy;
155
156     if (cr != nullptr)
157     {
158         if (PAR(cr))
159         {
160             gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
161         }
162     }
163
164     /* Scale box with Ewald wall factor */
165     matrix          scaledBox;
166     EwaldBoxZScaler boxScaler(inputrecPbcXY2Walls(&ir), ir.wall_ewald_zfac);
167     boxScaler.scaleBox(box, scaledBox);
168
169     rvec boxDiag;
170     for (int i = 0; (i < DIM); i++)
171     {
172         boxDiag[i] = scaledBox[i][i];
173     }
174
175     /* 1/(Vol*e0) */
176     real scaleRecip =
177             4.0 * M_PI / (boxDiag[XX] * boxDiag[YY] * boxDiag[ZZ]) * gmx::c_one4PiEps0 / ir.epsilon_r;
178
179     snew(eir, et->kmax);
180     for (n = 0; n < et->kmax; n++)
181     {
182         snew(eir[n], natoms);
183     }
184     et->tab_xy.resize(natoms);
185     et->tab_qxyz.resize(natoms);
186
187     bFreeEnergy = (ir.efep != FreeEnergyPerturbationType::No);
188
189     clear_mat(lrvir);
190
191     calc_lll(boxDiag, lll);
192     tabulateStructureFactors(natoms, x, et->kmax, eir, lll);
193
194     gmx::ArrayRef<const real> charge;
195     for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
196     {
197         if (!bFreeEnergy)
198         {
199             charge = chargeA;
200             scale  = 1.0;
201         }
202         else if (q == 0)
203         {
204             charge = chargeA;
205             scale  = 1.0 - lambda;
206         }
207         else
208         {
209             charge = chargeB;
210             scale  = lambda;
211         }
212         lowiy        = 0;
213         lowiz        = 1;
214         energy_AB[q] = 0;
215         for (ix = 0; ix < et->nx; ix++)
216         {
217             mx = ix * lll[XX];
218             for (iy = lowiy; iy < et->ny; iy++)
219             {
220                 my = iy * lll[YY];
221                 if (iy >= 0)
222                 {
223                     for (n = 0; n < natoms; n++)
224                     {
225                         et->tab_xy[n] = cmul(eir[ix][n][XX], eir[iy][n][YY]);
226                     }
227                 }
228                 else
229                 {
230                     for (n = 0; n < natoms; n++)
231                     {
232                         et->tab_xy[n] = cmul(eir[ix][n][XX], conjugate(eir[-iy][n][YY]));
233                     }
234                 }
235                 for (iz = lowiz; iz < et->nz; iz++)
236                 {
237                     mz  = iz * lll[ZZ];
238                     m2  = mx * mx + my * my + mz * mz;
239                     ak  = std::exp(m2 * factor) / m2;
240                     akv = 2.0 * ak * (1.0 / m2 - factor);
241                     if (iz >= 0)
242                     {
243                         for (n = 0; n < natoms; n++)
244                         {
245                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n], eir[iz][n][ZZ]));
246                         }
247                     }
248                     else
249                     {
250                         for (n = 0; n < natoms; n++)
251                         {
252                             et->tab_qxyz[n] =
253                                     rcmul(charge[n], cmul(et->tab_xy[n], conjugate(eir[-iz][n][ZZ])));
254                         }
255                     }
256
257                     cs = ss = 0;
258                     for (n = 0; n < natoms; n++)
259                     {
260                         cs += et->tab_qxyz[n].re;
261                         ss += et->tab_qxyz[n].im;
262                     }
263                     energy_AB[q] += ak * (cs * cs + ss * ss);
264                     tmp = scale * akv * (cs * cs + ss * ss);
265                     lrvir[XX][XX] -= tmp * mx * mx;
266                     lrvir[XX][YY] -= tmp * mx * my;
267                     lrvir[XX][ZZ] -= tmp * mx * mz;
268                     lrvir[YY][YY] -= tmp * my * my;
269                     lrvir[YY][ZZ] -= tmp * my * mz;
270                     lrvir[ZZ][ZZ] -= tmp * mz * mz;
271                     for (n = 0; n < natoms; n++)
272                     {
273                         /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
274                         tmp = scale * ak * (cs * et->tab_qxyz[n].im - ss * et->tab_qxyz[n].re);
275                         f[n][XX] += tmp * mx * 2 * scaleRecip;
276                         f[n][YY] += tmp * my * 2 * scaleRecip;
277                         f[n][ZZ] += tmp * mz * 2 * scaleRecip;
278 #if 0
279                         f[n][XX] += tmp*mx;
280                         f[n][YY] += tmp*my;
281                         f[n][ZZ] += tmp*mz;
282 #endif
283                     }
284                     lowiz = 1 - et->nz;
285                 }
286                 lowiy = 1 - et->ny;
287             }
288         }
289     }
290
291     if (!bFreeEnergy)
292     {
293         energy = energy_AB[0];
294     }
295     else
296     {
297         energy = (1.0 - lambda) * energy_AB[0] + lambda * energy_AB[1];
298         *dvdlambda += scaleRecip * (energy_AB[1] - energy_AB[0]);
299     }
300
301     lrvir[XX][XX] = -0.5 * scaleRecip * (lrvir[XX][XX] + energy);
302     lrvir[XX][YY] = -0.5 * scaleRecip * (lrvir[XX][YY]);
303     lrvir[XX][ZZ] = -0.5 * scaleRecip * (lrvir[XX][ZZ]);
304     lrvir[YY][YY] = -0.5 * scaleRecip * (lrvir[YY][YY] + energy);
305     lrvir[YY][ZZ] = -0.5 * scaleRecip * (lrvir[YY][ZZ]);
306     lrvir[ZZ][ZZ] = -0.5 * scaleRecip * (lrvir[ZZ][ZZ] + energy);
307
308     lrvir[YY][XX] = lrvir[XX][YY];
309     lrvir[ZZ][XX] = lrvir[XX][ZZ];
310     lrvir[ZZ][YY] = lrvir[YY][ZZ];
311
312     energy *= scaleRecip;
313
314     return energy;
315 }
316
317 real ewald_charge_correction(const t_commrec*            commrec,
318                              const real                  epsilonR,
319                              const real                  ewaldcoeffQ,
320                              gmx::ArrayRef<const double> qsum,
321                              const real                  lambda,
322                              const matrix                box,
323                              real*                       dvdlambda,
324                              tensor                      vir)
325
326 {
327     real enercorr = 0;
328
329     if (MASTER(commrec))
330     {
331         /* Apply charge correction */
332         real vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
333
334         real fac = M_PI * gmx::c_one4PiEps0 / (epsilonR * 2.0 * vol * vol * gmx::square(ewaldcoeffQ));
335
336         real qs2A = qsum[0] * qsum[0];
337         real qs2B = qsum[1] * qsum[1];
338
339         real vc = (qs2A * (1 - lambda) + qs2B * lambda) * fac;
340
341         enercorr = -vol * vc;
342
343         *dvdlambda += -vol * (qs2B - qs2A) * fac;
344
345         for (int d = 0; d < DIM; d++)
346         {
347             vir[d][d] += vc;
348         }
349     }
350
351     return enercorr;
352 }