e27d793e5fa5bf9b402afe4ac8785515545d58c6
[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 <stdio.h>
53
54 #include <cmath>
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, 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               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],
252                                                                     et->eir[iz][n][ZZ]));
253                         }
254                     }
255                     else
256                     {
257                         for (n = 0; n < natoms; n++)
258                         {
259                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n],
260                                                                     conjugate(et->eir[-iz][n][ZZ])));
261                         }
262                     }
263
264                     cs = ss = 0;
265                     for (n = 0; n < natoms; n++)
266                     {
267                         cs += et->tab_qxyz[n].re;
268                         ss += et->tab_qxyz[n].im;
269                     }
270                     energy_AB[q]  += ak*(cs*cs+ss*ss);
271                     tmp            = scale*akv*(cs*cs+ss*ss);
272                     lrvir[XX][XX] -= tmp*mx*mx;
273                     lrvir[XX][YY] -= tmp*mx*my;
274                     lrvir[XX][ZZ] -= tmp*mx*mz;
275                     lrvir[YY][YY] -= tmp*my*my;
276                     lrvir[YY][ZZ] -= tmp*my*mz;
277                     lrvir[ZZ][ZZ] -= tmp*mz*mz;
278                     for (n = 0; n < natoms; n++)
279                     {
280                         /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
281                         tmp       = scale*ak*(cs*et->tab_qxyz[n].im-ss*et->tab_qxyz[n].re);
282                         f[n][XX] += tmp*mx*2*scaleRecip;
283                         f[n][YY] += tmp*my*2*scaleRecip;
284                         f[n][ZZ] += tmp*mz*2*scaleRecip;
285 #if 0
286                         f[n][XX] += tmp*mx;
287                         f[n][YY] += tmp*my;
288                         f[n][ZZ] += tmp*mz;
289 #endif
290                     }
291                     lowiz = 1-et->nz;
292                 }
293                 lowiy = 1-et->ny;
294             }
295         }
296     }
297
298     if (!bFreeEnergy)
299     {
300         energy = energy_AB[0];
301     }
302     else
303     {
304         energy      = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
305         *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
306     }
307
308     lrvir[XX][XX] = -0.5*scaleRecip*(lrvir[XX][XX]+energy);
309     lrvir[XX][YY] = -0.5*scaleRecip*(lrvir[XX][YY]);
310     lrvir[XX][ZZ] = -0.5*scaleRecip*(lrvir[XX][ZZ]);
311     lrvir[YY][YY] = -0.5*scaleRecip*(lrvir[YY][YY]+energy);
312     lrvir[YY][ZZ] = -0.5*scaleRecip*(lrvir[YY][ZZ]);
313     lrvir[ZZ][ZZ] = -0.5*scaleRecip*(lrvir[ZZ][ZZ]+energy);
314
315     lrvir[YY][XX] = lrvir[XX][YY];
316     lrvir[ZZ][XX] = lrvir[XX][ZZ];
317     lrvir[ZZ][YY] = lrvir[YY][ZZ];
318
319     energy *= scaleRecip;
320
321     return energy;
322 }
323
324 real ewald_charge_correction(const t_commrec *cr, t_forcerec *fr, real lambda,
325                              matrix box,
326                              real *dvdlambda, tensor vir)
327
328 {
329     real vol, fac, qs2A, qs2B, vc, enercorr;
330     int  d;
331
332     if (MASTER(cr))
333     {
334         /* Apply charge correction */
335         vol = box[XX][XX]*box[YY][YY]*box[ZZ][ZZ];
336
337         fac = M_PI*ONE_4PI_EPS0/(fr->ic->epsilon_r*2.0*vol*vol*gmx::square(fr->ic->ewaldcoeff_q));
338
339         qs2A = fr->qsum[0]*fr->qsum[0];
340         qs2B = fr->qsum[1]*fr->qsum[1];
341
342         vc = (qs2A*(1 - lambda) + qs2B*lambda)*fac;
343
344         enercorr = -vol*vc;
345
346         *dvdlambda += -vol*(qs2B - qs2A)*fac;
347
348         for (d = 0; d < DIM; d++)
349         {
350             vir[d][d] += vc;
351         }
352
353         if (debug)
354         {
355             fprintf(debug, "Total charge correction: Vcharge=%g\n", enercorr);
356         }
357     }
358     else
359     {
360         enercorr = 0;
361     }
362
363     return enercorr;
364 }