Avoid allocating SYCL buffer on each call to PME solve
[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(bool                           havePbcXY2Walls,
135               real                           wallEwaldZfac,
136               real                           epsilonR,
137               FreeEnergyPerturbationType     freeEnergyPerturbationType,
138               gmx::ArrayRef<const gmx::RVec> coords,
139               gmx::ArrayRef<gmx::RVec>       forces,
140               gmx::ArrayRef<const real>      chargeA,
141               gmx::ArrayRef<const real>      chargeB,
142               const matrix                   box,
143               const t_commrec*               commrec,
144               int                            natoms,
145               matrix                         lrvir,
146               real                           ewaldcoeff,
147               real                           lambda,
148               real*                          dvdlambda,
149               gmx_ewald_tab_t*               et)
150 {
151     real   factor = -1.0 / (4 * ewaldcoeff * ewaldcoeff);
152     real   energy_AB[2], energy;
153     rvec   lll;
154     int    lowiy, lowiz, ix, iy, iz, n, q;
155     real   tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
156     cvec** eir;
157     bool   bFreeEnergy;
158
159     if (commrec != nullptr)
160     {
161         if (PAR(commrec))
162         {
163             gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
164         }
165     }
166
167     /* Scale box with Ewald wall factor */
168     matrix          scaledBox;
169     EwaldBoxZScaler boxScaler(havePbcXY2Walls, wallEwaldZfac);
170     boxScaler.scaleBox(box, scaledBox);
171
172     rvec boxDiag;
173     for (int i = 0; (i < DIM); i++)
174     {
175         boxDiag[i] = scaledBox[i][i];
176     }
177
178     /* 1/(Vol*e0) */
179     real scaleRecip = 4.0 * M_PI / (boxDiag[XX] * boxDiag[YY] * boxDiag[ZZ]) * gmx::c_one4PiEps0 / epsilonR;
180
181     snew(eir, et->kmax);
182     for (n = 0; n < et->kmax; n++)
183     {
184         snew(eir[n], natoms);
185     }
186     et->tab_xy.resize(natoms);
187     et->tab_qxyz.resize(natoms);
188
189     bFreeEnergy = (freeEnergyPerturbationType != FreeEnergyPerturbationType::No);
190
191     clear_mat(lrvir);
192
193     calc_lll(boxDiag, lll);
194     tabulateStructureFactors(natoms, coords, et->kmax, eir, lll);
195
196     gmx::ArrayRef<const real> charge;
197     for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
198     {
199         if (!bFreeEnergy)
200         {
201             charge = chargeA;
202             scale  = 1.0;
203         }
204         else if (q == 0)
205         {
206             charge = chargeA;
207             scale  = 1.0 - lambda;
208         }
209         else
210         {
211             charge = chargeB;
212             scale  = lambda;
213         }
214         lowiy        = 0;
215         lowiz        = 1;
216         energy_AB[q] = 0;
217         for (ix = 0; ix < et->nx; ix++)
218         {
219             mx = ix * lll[XX];
220             for (iy = lowiy; iy < et->ny; iy++)
221             {
222                 my = iy * lll[YY];
223                 if (iy >= 0)
224                 {
225                     for (n = 0; n < natoms; n++)
226                     {
227                         et->tab_xy[n] = cmul(eir[ix][n][XX], eir[iy][n][YY]);
228                     }
229                 }
230                 else
231                 {
232                     for (n = 0; n < natoms; n++)
233                     {
234                         et->tab_xy[n] = cmul(eir[ix][n][XX], conjugate(eir[-iy][n][YY]));
235                     }
236                 }
237                 for (iz = lowiz; iz < et->nz; iz++)
238                 {
239                     mz  = iz * lll[ZZ];
240                     m2  = mx * mx + my * my + mz * mz;
241                     ak  = std::exp(m2 * factor) / m2;
242                     akv = 2.0 * ak * (1.0 / m2 - factor);
243                     if (iz >= 0)
244                     {
245                         for (n = 0; n < natoms; n++)
246                         {
247                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n], eir[iz][n][ZZ]));
248                         }
249                     }
250                     else
251                     {
252                         for (n = 0; n < natoms; n++)
253                         {
254                             et->tab_qxyz[n] =
255                                     rcmul(charge[n], cmul(et->tab_xy[n], conjugate(eir[-iz][n][ZZ])));
256                         }
257                     }
258
259                     cs = ss = 0;
260                     for (n = 0; n < natoms; n++)
261                     {
262                         cs += et->tab_qxyz[n].re;
263                         ss += et->tab_qxyz[n].im;
264                     }
265                     energy_AB[q] += ak * (cs * cs + ss * ss);
266                     tmp = scale * akv * (cs * cs + ss * ss);
267                     lrvir[XX][XX] -= tmp * mx * mx;
268                     lrvir[XX][YY] -= tmp * mx * my;
269                     lrvir[XX][ZZ] -= tmp * mx * mz;
270                     lrvir[YY][YY] -= tmp * my * my;
271                     lrvir[YY][ZZ] -= tmp * my * mz;
272                     lrvir[ZZ][ZZ] -= tmp * mz * mz;
273                     for (n = 0; n < natoms; n++)
274                     {
275                         /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
276                         tmp = scale * ak * (cs * et->tab_qxyz[n].im - ss * et->tab_qxyz[n].re);
277                         forces[n][XX] += tmp * mx * 2 * scaleRecip;
278                         forces[n][YY] += tmp * my * 2 * scaleRecip;
279                         forces[n][ZZ] += tmp * mz * 2 * scaleRecip;
280                     }
281                     lowiz = 1 - et->nz;
282                 }
283                 lowiy = 1 - et->ny;
284             }
285         }
286     }
287
288     if (!bFreeEnergy)
289     {
290         energy = energy_AB[0];
291     }
292     else
293     {
294         energy = (1.0 - lambda) * energy_AB[0] + lambda * energy_AB[1];
295         *dvdlambda += scaleRecip * (energy_AB[1] - energy_AB[0]);
296     }
297
298     lrvir[XX][XX] = -0.5 * scaleRecip * (lrvir[XX][XX] + energy);
299     lrvir[XX][YY] = -0.5 * scaleRecip * (lrvir[XX][YY]);
300     lrvir[XX][ZZ] = -0.5 * scaleRecip * (lrvir[XX][ZZ]);
301     lrvir[YY][YY] = -0.5 * scaleRecip * (lrvir[YY][YY] + energy);
302     lrvir[YY][ZZ] = -0.5 * scaleRecip * (lrvir[YY][ZZ]);
303     lrvir[ZZ][ZZ] = -0.5 * scaleRecip * (lrvir[ZZ][ZZ] + energy);
304
305     lrvir[YY][XX] = lrvir[XX][YY];
306     lrvir[ZZ][XX] = lrvir[XX][ZZ];
307     lrvir[ZZ][YY] = lrvir[YY][ZZ];
308
309     energy *= scaleRecip;
310
311     return energy;
312 }
313
314 real ewald_charge_correction(const t_commrec*            commrec,
315                              const real                  epsilonR,
316                              const real                  ewaldcoeffQ,
317                              gmx::ArrayRef<const double> qsum,
318                              const real                  lambda,
319                              const matrix                box,
320                              real*                       dvdlambda,
321                              tensor                      vir)
322
323 {
324     real enercorr = 0;
325
326     if (MASTER(commrec))
327     {
328         /* Apply charge correction */
329         real vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
330
331         real fac = M_PI * gmx::c_one4PiEps0 / (epsilonR * 2.0 * vol * vol * gmx::square(ewaldcoeffQ));
332
333         real qs2A = qsum[0] * qsum[0];
334         real qs2B = qsum[1] * qsum[1];
335
336         real vc = (qs2A * (1 - lambda) + qs2B * lambda) * fac;
337
338         enercorr = -vol * vc;
339
340         *dvdlambda += -vol * (qs2B - qs2A) * fac;
341
342         for (int d = 0; d < DIM; d++)
343         {
344             vir[d][d] += vc;
345         }
346     }
347
348     return enercorr;
349 }