696a47c9d38e7b77b1ddd8f6b6c22fcf342e6465
[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, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 /*! \internal \file
38  *
39  * \brief This file contains function definitions necessary for
40  * computing energies and forces for the plain-Ewald long-ranged part,
41  * and the correction for overall system charge for all Ewald-family
42  * methods.
43  *
44  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
45  * \author Mark Abraham <mark.j.abraham@gmail.com>
46  * \ingroup module_ewald
47  */
48 #include "gmxpre.h"
49
50 #include "ewald.h"
51
52 #include <cmath>
53 #include <cstdio>
54 #include <cstdlib>
55
56 #include <algorithm>
57
58 #include "gromacs/ewald/ewald-utils.h"
59 #include "gromacs/math/functions.h"
60 #include "gromacs/math/gmxcomplex.h"
61 #include "gromacs/math/units.h"
62 #include "gromacs/math/utilities.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/math/vectypes.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/forcerec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/smalloc.h"
71
72 struct gmx_ewald_tab_t
73 {
74     int        nx, ny, nz, kmax;
75     cvec     **eir;
76     t_complex *tab_xy, *tab_qxyz;
77 };
78
79 void init_ewald_tab(struct gmx_ewald_tab_t **et, const t_inputrec *ir, FILE *fp)
80 {
81     snew(*et, 1);
82     if (fp)
83     {
84         fprintf(fp, "Will do ordinary reciprocal space Ewald sum.\n");
85     }
86
87     (*et)->nx       = ir->nkx+1;
88     (*et)->ny       = ir->nky+1;
89     (*et)->nz       = ir->nkz+1;
90     (*et)->kmax     = std::max((*et)->nx, std::max((*et)->ny, (*et)->nz));
91     (*et)->eir      = nullptr;
92     (*et)->tab_xy   = nullptr;
93     (*et)->tab_qxyz = nullptr;
94 }
95
96 //! Calculates wave vectors.
97 static void calc_lll(const rvec box, rvec lll)
98 {
99     lll[XX] = 2.0*M_PI/box[XX];
100     lll[YY] = 2.0*M_PI/box[YY];
101     lll[ZZ] = 2.0*M_PI/box[ZZ];
102 }
103
104 //! Make tables for the structure factor parts
105 static void tabulateStructureFactors(int natom, const rvec x[], int kmax, cvec **eir, const rvec lll)
106 {
107     int  i, j, m;
108
109     if (kmax < 1)
110     {
111         printf("Go away! kmax = %d\n", kmax);
112         exit(1);
113     }
114
115     for (i = 0; (i < natom); i++)
116     {
117         for (m = 0; (m < 3); m++)
118         {
119             eir[0][i][m].re = 1;
120             eir[0][i][m].im = 0;
121         }
122
123         for (m = 0; (m < 3); m++)
124         {
125             eir[1][i][m].re = std::cos(x[i][m]*lll[m]);
126             eir[1][i][m].im = std::sin(x[i][m]*lll[m]);
127         }
128         for (j = 2; (j < kmax); j++)
129         {
130             for (m = 0; (m < 3); m++)
131             {
132                 eir[j][i][m] = cmul(eir[j-1][i][m], eir[1][i][m]);
133             }
134         }
135     }
136 }
137
138 real do_ewald(const t_inputrec *ir,
139               const rvec        x[],
140               rvec              f[],
141               const real        chargeA[],
142               const real        chargeB[],
143               matrix            box,
144               const t_commrec  *cr,
145               int               natoms,
146               matrix            lrvir,
147               real              ewaldcoeff,
148               real              lambda,
149               real             *dvdlambda,
150               gmx_ewald_tab_t  *et)
151 {
152     real        factor     = -1.0/(4*ewaldcoeff*ewaldcoeff);
153     const real *charge;
154     real        energy_AB[2], energy;
155     rvec        lll;
156     int         lowiy, lowiz, ix, iy, iz, n, q;
157     real        tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
158     gmx_bool    bFreeEnergy;
159
160     if (cr != nullptr)
161     {
162         if (PAR(cr))
163         {
164             gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
165         }
166     }
167
168     /* Scale box with Ewald wall factor */
169     matrix          scaledBox;
170     EwaldBoxZScaler boxScaler(*ir);
171     boxScaler.scaleBox(box, scaledBox);
172
173     rvec boxDiag;
174     for (int i = 0; (i < DIM); i++)
175     {
176         boxDiag[i] = scaledBox[i][i];
177     }
178
179     /* 1/(Vol*e0) */
180     real scaleRecip = 4.0*M_PI/(boxDiag[XX]*boxDiag[YY]*boxDiag[ZZ])*ONE_4PI_EPS0/ir->epsilon_r;
181
182     if (!et->eir) /* allocate if we need to */
183     {
184         snew(et->eir, et->kmax);
185         for (n = 0; n < et->kmax; n++)
186         {
187             snew(et->eir[n], natoms);
188         }
189         snew(et->tab_xy, natoms);
190         snew(et->tab_qxyz, natoms);
191     }
192
193     bFreeEnergy = (ir->efep != efepNO);
194
195     clear_mat(lrvir);
196
197     calc_lll(boxDiag, lll);
198     tabulateStructureFactors(natoms, x, et->kmax, et->eir, lll);
199
200     for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
201     {
202         if (!bFreeEnergy)
203         {
204             charge = chargeA;
205             scale  = 1.0;
206         }
207         else if (q == 0)
208         {
209             charge = chargeA;
210             scale  = 1.0 - lambda;
211         }
212         else
213         {
214             charge = chargeB;
215             scale  = lambda;
216         }
217         lowiy        = 0;
218         lowiz        = 1;
219         energy_AB[q] = 0;
220         for (ix = 0; ix < et->nx; ix++)
221         {
222             mx = ix*lll[XX];
223             for (iy = lowiy; iy < et->ny; iy++)
224             {
225                 my = iy*lll[YY];
226                 if (iy >= 0)
227                 {
228                     for (n = 0; n < natoms; n++)
229                     {
230                         et->tab_xy[n] = cmul(et->eir[ix][n][XX], et->eir[iy][n][YY]);
231                     }
232                 }
233                 else
234                 {
235                     for (n = 0; n < natoms; n++)
236                     {
237                         et->tab_xy[n] = cmul(et->eir[ix][n][XX], conjugate(et->eir[-iy][n][YY]));
238                     }
239                 }
240                 for (iz = lowiz; iz < et->nz; iz++)
241                 {
242                     mz  = iz*lll[ZZ];
243                     m2  = mx*mx+my*my+mz*mz;
244                     ak  = std::exp(m2*factor)/m2;
245                     akv = 2.0*ak*(1.0/m2-factor);
246                     if (iz >= 0)
247                     {
248                         for (n = 0; n < natoms; n++)
249                         {
250                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n],
251                                                                     et->eir[iz][n][ZZ]));
252                         }
253                     }
254                     else
255                     {
256                         for (n = 0; n < natoms; n++)
257                         {
258                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n],
259                                                                     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, t_forcerec *fr, real lambda,
324                              matrix box,
325                              real *dvdlambda, tensor vir)
326
327 {
328     real vol, fac, qs2A, qs2B, vc, enercorr;
329     int  d;
330
331     if (MASTER(cr))
332     {
333         /* Apply charge correction */
334         vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
335
336         fac = M_PI*ONE_4PI_EPS0/(fr->ic->epsilon_r*2.0*vol*vol*gmx::square(fr->ic->ewaldcoeff_q));
337
338         qs2A = fr->qsum[0]*fr->qsum[0];
339         qs2B = fr->qsum[1]*fr->qsum[1];
340
341         vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
342
343         enercorr = -vol*vc;
344
345         *dvdlambda += -vol*(qs2B - qs2A)*fac;
346
347         for (d = 0; d < DIM; d++)
348         {
349             vir[d][d] += vc;
350         }
351
352         if (debug)
353         {
354             fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);
355         }
356     }
357     else
358     {
359         enercorr = 0;
360     }
361
362     return enercorr;
363 }