e4b789df76f2e235f3741b21add66d24cd93d299
[alexxy/gromacs.git] / src / gromacs / mdlib / ewald.c
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, 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 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include <math.h>
42 #include <stdio.h>
43 #include <stdlib.h>
44
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "gromacs/legacyheaders/types/commrec.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/legacyheaders/coulomb.h"
52 #include "gromacs/legacyheaders/macros.h"
53
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/math/gmxcomplex.h"
56
57 #define TOL 2e-5
58
59 struct ewald_tab
60 {
61     int        nx, ny, nz, kmax;
62     cvec     **eir;
63     t_complex *tab_xy, *tab_qxyz;
64 };
65
66
67
68 /* TODO: fix thread-safety */
69
70 /* the other routines are in complex.h */
71 static t_complex conjmul(t_complex a, t_complex b)
72 {
73     t_complex c;
74
75     c.re = a.re*b.re + a.im*b.im;
76     c.im = a.im*b.re - a.re*b.im;
77
78     return c;
79 }
80
81
82
83
84 static void tabulate_eir(int natom, rvec x[], int kmax, cvec **eir, rvec lll)
85 {
86     int  i, j, m;
87
88     if (kmax < 1)
89     {
90         printf("Go away! kmax = %d\n", kmax);
91         exit(1);
92     }
93
94     for (i = 0; (i < natom); i++)
95     {
96         for (m = 0; (m < 3); m++)
97         {
98             eir[0][i][m].re = 1;
99             eir[0][i][m].im = 0;
100         }
101
102         for (m = 0; (m < 3); m++)
103         {
104             eir[1][i][m].re = cos(x[i][m]*lll[m]);
105             eir[1][i][m].im = sin(x[i][m]*lll[m]);
106         }
107         for (j = 2; (j < kmax); j++)
108         {
109             for (m = 0; (m < 3); m++)
110             {
111                 eir[j][i][m] = cmul(eir[j-1][i][m], eir[1][i][m]);
112             }
113         }
114     }
115 }
116
117 void init_ewald_tab(ewald_tab_t *et, const t_inputrec *ir, FILE *fp)
118 {
119     int n;
120
121     snew(*et, 1);
122     if (fp)
123     {
124         fprintf(fp, "Will do ordinary reciprocal space Ewald sum.\n");
125     }
126
127     (*et)->nx       = ir->nkx+1;
128     (*et)->ny       = ir->nky+1;
129     (*et)->nz       = ir->nkz+1;
130     (*et)->kmax     = max((*et)->nx, max((*et)->ny, (*et)->nz));
131     (*et)->eir      = NULL;
132     (*et)->tab_xy   = NULL;
133     (*et)->tab_qxyz = NULL;
134 }
135
136
137
138 real do_ewald(t_inputrec *ir,
139               rvec x[],        rvec f[],
140               real chargeA[],  real chargeB[],
141               rvec box,
142               t_commrec *cr,   int natoms,
143               matrix lrvir,    real ewaldcoeff,
144               real lambda,     real *dvdlambda,
145               ewald_tab_t et)
146 {
147     real     factor     = -1.0/(4*ewaldcoeff*ewaldcoeff);
148     real     scaleRecip = 4.0*M_PI/(box[XX]*box[YY]*box[ZZ])*ONE_4PI_EPS0/ir->epsilon_r; /* 1/(Vol*e0) */
149     real    *charge, 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     gmx_bool bFreeEnergy;
154
155     if (cr != NULL)
156     {
157         if (PAR(cr))
158         {
159             gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
160         }
161     }
162
163
164     if (!et->eir) /* allocate if we need to */
165     {
166         snew(et->eir, et->kmax);
167         for (n = 0; n < et->kmax; n++)
168         {
169             snew(et->eir[n], natoms);
170         }
171         snew(et->tab_xy, natoms);
172         snew(et->tab_qxyz, natoms);
173     }
174
175     bFreeEnergy = (ir->efep != efepNO);
176
177     clear_mat(lrvir);
178
179     calc_lll(box, lll);
180     /* make tables for the structure factor parts */
181     tabulate_eir(natoms, x, et->kmax, et->eir, lll);
182
183     for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
184     {
185         if (!bFreeEnergy)
186         {
187             charge = chargeA;
188             scale  = 1.0;
189         }
190         else if (q == 0)
191         {
192             charge = chargeA;
193             scale  = 1.0 - lambda;
194         }
195         else
196         {
197             charge = chargeB;
198             scale  = lambda;
199         }
200         lowiy        = 0;
201         lowiz        = 1;
202         energy_AB[q] = 0;
203         for (ix = 0; ix < et->nx; ix++)
204         {
205             mx = ix*lll[XX];
206             for (iy = lowiy; iy < et->ny; iy++)
207             {
208                 my = iy*lll[YY];
209                 if (iy >= 0)
210                 {
211                     for (n = 0; n < natoms; n++)
212                     {
213                         et->tab_xy[n] = cmul(et->eir[ix][n][XX], et->eir[iy][n][YY]);
214                     }
215                 }
216                 else
217                 {
218                     for (n = 0; n < natoms; n++)
219                     {
220                         et->tab_xy[n] = conjmul(et->eir[ix][n][XX], et->eir[-iy][n][YY]);
221                     }
222                 }
223                 for (iz = lowiz; iz < et->nz; iz++)
224                 {
225                     mz  = iz*lll[ZZ];
226                     m2  = mx*mx+my*my+mz*mz;
227                     ak  = exp(m2*factor)/m2;
228                     akv = 2.0*ak*(1.0/m2-factor);
229                     if (iz >= 0)
230                     {
231                         for (n = 0; n < natoms; n++)
232                         {
233                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n],
234                                                                     et->eir[iz][n][ZZ]));
235                         }
236                     }
237                     else
238                     {
239                         for (n = 0; n < natoms; n++)
240                         {
241                             et->tab_qxyz[n] = rcmul(charge[n], conjmul(et->tab_xy[n],
242                                                                        et->eir[-iz][n][ZZ]));
243                         }
244                     }
245
246                     cs = ss = 0;
247                     for (n = 0; n < natoms; n++)
248                     {
249                         cs += et->tab_qxyz[n].re;
250                         ss += et->tab_qxyz[n].im;
251                     }
252                     energy_AB[q]  += ak*(cs*cs+ss*ss);
253                     tmp            = scale*akv*(cs*cs+ss*ss);
254                     lrvir[XX][XX] -= tmp*mx*mx;
255                     lrvir[XX][YY] -= tmp*mx*my;
256                     lrvir[XX][ZZ] -= tmp*mx*mz;
257                     lrvir[YY][YY] -= tmp*my*my;
258                     lrvir[YY][ZZ] -= tmp*my*mz;
259                     lrvir[ZZ][ZZ] -= tmp*mz*mz;
260                     for (n = 0; n < natoms; n++)
261                     {
262                         /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
263                         tmp       = scale*ak*(cs*et->tab_qxyz[n].im-ss*et->tab_qxyz[n].re);
264                         f[n][XX] += tmp*mx*2*scaleRecip;
265                         f[n][YY] += tmp*my*2*scaleRecip;
266                         f[n][ZZ] += tmp*mz*2*scaleRecip;
267 #if 0
268                         f[n][XX] += tmp*mx;
269                         f[n][YY] += tmp*my;
270                         f[n][ZZ] += tmp*mz;
271 #endif
272                     }
273                     lowiz = 1-et->nz;
274                 }
275                 lowiy = 1-et->ny;
276             }
277         }
278     }
279
280     if (!bFreeEnergy)
281     {
282         energy = energy_AB[0];
283     }
284     else
285     {
286         energy      = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
287         *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
288     }
289
290     lrvir[XX][XX] = -0.5*scaleRecip*(lrvir[XX][XX]+energy);
291     lrvir[XX][YY] = -0.5*scaleRecip*(lrvir[XX][YY]);
292     lrvir[XX][ZZ] = -0.5*scaleRecip*(lrvir[XX][ZZ]);
293     lrvir[YY][YY] = -0.5*scaleRecip*(lrvir[YY][YY]+energy);
294     lrvir[YY][ZZ] = -0.5*scaleRecip*(lrvir[YY][ZZ]);
295     lrvir[ZZ][ZZ] = -0.5*scaleRecip*(lrvir[ZZ][ZZ]+energy);
296
297     lrvir[YY][XX] = lrvir[XX][YY];
298     lrvir[ZZ][XX] = lrvir[XX][ZZ];
299     lrvir[ZZ][YY] = lrvir[YY][ZZ];
300
301     energy *= scaleRecip;
302
303     return energy;
304 }