bd057371920f7d529b8362dedde69760e22ed4fc
[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, 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/md_enums.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/smalloc.h"
72
73 struct gmx_ewald_tab_t
74 {
75     int        nx, ny, nz, kmax;
76     cvec**     eir;
77     t_complex *tab_xy, *tab_qxyz;
78 };
79
80 void init_ewald_tab(struct gmx_ewald_tab_t** et, const t_inputrec* ir, FILE* fp)
81 {
82     snew(*et, 1);
83     if (fp)
84     {
85         fprintf(fp, "Will do ordinary reciprocal space Ewald sum.\n");
86     }
87
88     (*et)->nx       = ir->nkx + 1;
89     (*et)->ny       = ir->nky + 1;
90     (*et)->nz       = ir->nkz + 1;
91     (*et)->kmax     = std::max((*et)->nx, std::max((*et)->ny, (*et)->nz));
92     (*et)->eir      = nullptr;
93     (*et)->tab_xy   = nullptr;
94     (*et)->tab_qxyz = nullptr;
95 }
96
97 //! Calculates wave vectors.
98 static void calc_lll(const rvec box, rvec lll)
99 {
100     lll[XX] = 2.0 * M_PI / box[XX];
101     lll[YY] = 2.0 * M_PI / box[YY];
102     lll[ZZ] = 2.0 * M_PI / box[ZZ];
103 }
104
105 //! Make tables for the structure factor parts
106 static void tabulateStructureFactors(int natom, const rvec x[], int kmax, cvec** eir, const rvec lll)
107 {
108     int i, j, m;
109
110     if (kmax < 1)
111     {
112         printf("Go away! kmax = %d\n", kmax);
113         exit(1);
114     }
115
116     for (i = 0; (i < natom); i++)
117     {
118         for (m = 0; (m < 3); m++)
119         {
120             eir[0][i][m].re = 1;
121             eir[0][i][m].im = 0;
122         }
123
124         for (m = 0; (m < 3); m++)
125         {
126             eir[1][i][m].re = std::cos(x[i][m] * lll[m]);
127             eir[1][i][m].im = std::sin(x[i][m] * lll[m]);
128         }
129         for (j = 2; (j < kmax); j++)
130         {
131             for (m = 0; (m < 3); m++)
132             {
133                 eir[j][i][m] = cmul(eir[j - 1][i][m], eir[1][i][m]);
134             }
135         }
136     }
137 }
138
139 real do_ewald(const t_inputrec* ir,
140               const rvec        x[],
141               rvec              f[],
142               const real        chargeA[],
143               const real        chargeB[],
144               const matrix      box,
145               const t_commrec*  cr,
146               int               natoms,
147               matrix            lrvir,
148               real              ewaldcoeff,
149               real              lambda,
150               real*             dvdlambda,
151               gmx_ewald_tab_t*  et)
152 {
153     real        factor = -1.0 / (4 * ewaldcoeff * ewaldcoeff);
154     const real* charge;
155     real        energy_AB[2], energy;
156     rvec        lll;
157     int         lowiy, lowiz, ix, iy, iz, n, q;
158     real        tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
159     gmx_bool    bFreeEnergy;
160
161     if (cr != nullptr)
162     {
163         if (PAR(cr))
164         {
165             gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
166         }
167     }
168
169     /* Scale box with Ewald wall factor */
170     matrix          scaledBox;
171     EwaldBoxZScaler boxScaler(*ir);
172     boxScaler.scaleBox(box, scaledBox);
173
174     rvec boxDiag;
175     for (int i = 0; (i < DIM); i++)
176     {
177         boxDiag[i] = scaledBox[i][i];
178     }
179
180     /* 1/(Vol*e0) */
181     real scaleRecip = 4.0 * M_PI / (boxDiag[XX] * boxDiag[YY] * boxDiag[ZZ]) * ONE_4PI_EPS0 / ir->epsilon_r;
182
183     if (!et->eir) /* allocate if we need to */
184     {
185         snew(et->eir, et->kmax);
186         for (n = 0; n < et->kmax; n++)
187         {
188             snew(et->eir[n], natoms);
189         }
190         snew(et->tab_xy, natoms);
191         snew(et->tab_qxyz, natoms);
192     }
193
194     bFreeEnergy = (ir->efep != efepNO);
195
196     clear_mat(lrvir);
197
198     calc_lll(boxDiag, lll);
199     tabulateStructureFactors(natoms, x, et->kmax, et->eir, lll);
200
201     for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
202     {
203         if (!bFreeEnergy)
204         {
205             charge = chargeA;
206             scale  = 1.0;
207         }
208         else if (q == 0)
209         {
210             charge = chargeA;
211             scale  = 1.0 - lambda;
212         }
213         else
214         {
215             charge = chargeB;
216             scale  = lambda;
217         }
218         lowiy        = 0;
219         lowiz        = 1;
220         energy_AB[q] = 0;
221         for (ix = 0; ix < et->nx; ix++)
222         {
223             mx = ix * lll[XX];
224             for (iy = lowiy; iy < et->ny; iy++)
225             {
226                 my = iy * lll[YY];
227                 if (iy >= 0)
228                 {
229                     for (n = 0; n < natoms; n++)
230                     {
231                         et->tab_xy[n] = cmul(et->eir[ix][n][XX], et->eir[iy][n][YY]);
232                     }
233                 }
234                 else
235                 {
236                     for (n = 0; n < natoms; n++)
237                     {
238                         et->tab_xy[n] = cmul(et->eir[ix][n][XX], conjugate(et->eir[-iy][n][YY]));
239                     }
240                 }
241                 for (iz = lowiz; iz < et->nz; iz++)
242                 {
243                     mz  = iz * lll[ZZ];
244                     m2  = mx * mx + my * my + mz * mz;
245                     ak  = std::exp(m2 * factor) / m2;
246                     akv = 2.0 * ak * (1.0 / m2 - factor);
247                     if (iz >= 0)
248                     {
249                         for (n = 0; n < natoms; n++)
250                         {
251                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n], et->eir[iz][n][ZZ]));
252                         }
253                     }
254                     else
255                     {
256                         for (n = 0; n < natoms; n++)
257                         {
258                             et->tab_qxyz[n] = rcmul(
259                                     charge[n], cmul(et->tab_xy[n], conjugate(et->eir[-iz][n][ZZ])));
260                         }
261                     }
262
263                     cs = ss = 0;
264                     for (n = 0; n < natoms; n++)
265                     {
266                         cs += et->tab_qxyz[n].re;
267                         ss += et->tab_qxyz[n].im;
268                     }
269                     energy_AB[q] += ak * (cs * cs + ss * ss);
270                     tmp = scale * akv * (cs * cs + ss * ss);
271                     lrvir[XX][XX] -= tmp * mx * mx;
272                     lrvir[XX][YY] -= tmp * mx * my;
273                     lrvir[XX][ZZ] -= tmp * mx * mz;
274                     lrvir[YY][YY] -= tmp * my * my;
275                     lrvir[YY][ZZ] -= tmp * my * mz;
276                     lrvir[ZZ][ZZ] -= tmp * mz * mz;
277                     for (n = 0; n < natoms; n++)
278                     {
279                         /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
280                         tmp = scale * ak * (cs * et->tab_qxyz[n].im - ss * et->tab_qxyz[n].re);
281                         f[n][XX] += tmp * mx * 2 * scaleRecip;
282                         f[n][YY] += tmp * my * 2 * scaleRecip;
283                         f[n][ZZ] += tmp * mz * 2 * scaleRecip;
284 #if 0
285                         f[n][XX] += tmp*mx;
286                         f[n][YY] += tmp*my;
287                         f[n][ZZ] += tmp*mz;
288 #endif
289                     }
290                     lowiz = 1 - et->nz;
291                 }
292                 lowiy = 1 - et->ny;
293             }
294         }
295     }
296
297     if (!bFreeEnergy)
298     {
299         energy = energy_AB[0];
300     }
301     else
302     {
303         energy = (1.0 - lambda) * energy_AB[0] + lambda * energy_AB[1];
304         *dvdlambda += scaleRecip * (energy_AB[1] - energy_AB[0]);
305     }
306
307     lrvir[XX][XX] = -0.5 * scaleRecip * (lrvir[XX][XX] + energy);
308     lrvir[XX][YY] = -0.5 * scaleRecip * (lrvir[XX][YY]);
309     lrvir[XX][ZZ] = -0.5 * scaleRecip * (lrvir[XX][ZZ]);
310     lrvir[YY][YY] = -0.5 * scaleRecip * (lrvir[YY][YY] + energy);
311     lrvir[YY][ZZ] = -0.5 * scaleRecip * (lrvir[YY][ZZ]);
312     lrvir[ZZ][ZZ] = -0.5 * scaleRecip * (lrvir[ZZ][ZZ] + energy);
313
314     lrvir[YY][XX] = lrvir[XX][YY];
315     lrvir[ZZ][XX] = lrvir[XX][ZZ];
316     lrvir[ZZ][YY] = lrvir[YY][ZZ];
317
318     energy *= scaleRecip;
319
320     return energy;
321 }
322
323 real ewald_charge_correction(const t_commrec*  cr,
324                              const t_forcerec* fr,
325                              const real        lambda,
326                              const matrix      box,
327                              real*             dvdlambda,
328                              tensor            vir)
329
330 {
331     real vol, fac, qs2A, qs2B, vc, enercorr;
332     int  d;
333
334     if (MASTER(cr))
335     {
336         /* Apply charge correction */
337         vol = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
338
339         fac = M_PI * ONE_4PI_EPS0
340               / (fr->ic->epsilon_r * 2.0 * vol * vol * gmx::square(fr->ic->ewaldcoeff_q));
341
342         qs2A = fr->qsum[0] * fr->qsum[0];
343         qs2B = fr->qsum[1] * fr->qsum[1];
344
345         vc = (qs2A * (1 - lambda) + qs2B * lambda) * fac;
346
347         enercorr = -vol * vc;
348
349         *dvdlambda += -vol * (qs2B - qs2A) * fac;
350
351         for (d = 0; d < DIM; d++)
352         {
353             vir[d][d] += vc;
354         }
355
356         if (debug)
357         {
358             fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);
359         }
360     }
361     else
362     {
363         enercorr = 0;
364     }
365
366     return enercorr;
367 }