Merge branch 'release-4-6' into master
[alexxy/gromacs.git] / src / gromacs / mdlib / ewald.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * GROwing Monsters And Cloning Shrimps
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <stdio.h>
40 #include <math.h>
41 #include "typedefs.h"
42 #include "vec.h"
43 #include "gmxcomplex.h"
44 #include "smalloc.h"
45 #include "futil.h"
46 #include "gmx_fatal.h"
47 #include "physics.h"
48 #include "coulomb.h"
49 #include "macros.h"
50
51 #define TOL 2e-5
52
53 struct ewald_tab
54 {
55     int        nx, ny, nz, kmax;
56     cvec     **eir;
57     t_complex *tab_xy, *tab_qxyz;
58 };
59
60
61
62 /* TODO: fix thread-safety */
63
64 /* the other routines are in complex.h */
65 static t_complex conjmul(t_complex a, t_complex b)
66 {
67     t_complex c;
68
69     c.re = a.re*b.re + a.im*b.im;
70     c.im = a.im*b.re - a.re*b.im;
71
72     return c;
73 }
74
75
76
77
78 static void tabulate_eir(int natom, rvec x[], int kmax, cvec **eir, rvec lll)
79 {
80     int  i, j, m;
81
82     if (kmax < 1)
83     {
84         printf("Go away! kmax = %d\n", kmax);
85         exit(1);
86     }
87
88     for (i = 0; (i < natom); i++)
89     {
90         for (m = 0; (m < 3); m++)
91         {
92             eir[0][i][m].re = 1;
93             eir[0][i][m].im = 0;
94         }
95
96         for (m = 0; (m < 3); m++)
97         {
98             eir[1][i][m].re = cos(x[i][m]*lll[m]);
99             eir[1][i][m].im = sin(x[i][m]*lll[m]);
100         }
101         for (j = 2; (j < kmax); j++)
102         {
103             for (m = 0; (m < 3); m++)
104             {
105                 eir[j][i][m] = cmul(eir[j-1][i][m], eir[1][i][m]);
106             }
107         }
108     }
109 }
110
111 void init_ewald_tab(ewald_tab_t *et, const t_commrec *cr, const t_inputrec *ir,
112                     FILE *fp)
113 {
114     int n;
115
116     snew(*et, 1);
117     if (fp)
118     {
119         fprintf(fp, "Will do ordinary reciprocal space Ewald sum.\n");
120     }
121
122     (*et)->nx       = ir->nkx+1;
123     (*et)->ny       = ir->nky+1;
124     (*et)->nz       = ir->nkz+1;
125     (*et)->kmax     = max((*et)->nx, max((*et)->ny, (*et)->nz));
126     (*et)->eir      = NULL;
127     (*et)->tab_xy   = NULL;
128     (*et)->tab_qxyz = NULL;
129 }
130
131
132
133 real do_ewald(FILE *log,       gmx_bool bVerbose,
134               t_inputrec *ir,
135               rvec x[],        rvec f[],
136               real chargeA[],  real chargeB[],
137               rvec box,
138               t_commrec *cr,   int natoms,
139               matrix lrvir,    real ewaldcoeff,
140               real lambda,     real *dvdlambda,
141               ewald_tab_t et)
142 {
143     real     factor     = -1.0/(4*ewaldcoeff*ewaldcoeff);
144     real     scaleRecip = 4.0*M_PI/(box[XX]*box[YY]*box[ZZ])*ONE_4PI_EPS0/ir->epsilon_r; /* 1/(Vol*e0) */
145     real    *charge, energy_AB[2], energy;
146     rvec     lll;
147     int      lowiy, lowiz, ix, iy, iz, n, q;
148     real     tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
149     gmx_bool bFreeEnergy;
150
151     if (cr != NULL)
152     {
153         if (PAR(cr))
154         {
155             gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
156         }
157     }
158
159
160     if (!et->eir) /* allocate if we need to */
161     {
162         snew(et->eir, et->kmax);
163         for (n = 0; n < et->kmax; n++)
164         {
165             snew(et->eir[n], natoms);
166         }
167         snew(et->tab_xy, natoms);
168         snew(et->tab_qxyz, natoms);
169     }
170
171     bFreeEnergy = (ir->efep != efepNO);
172
173     clear_mat(lrvir);
174
175     calc_lll(box, lll);
176     /* make tables for the structure factor parts */
177     tabulate_eir(natoms, x, et->kmax, et->eir, lll);
178
179     for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
180     {
181         if (!bFreeEnergy)
182         {
183             charge = chargeA;
184             scale  = 1.0;
185         }
186         else if (q == 0)
187         {
188             charge = chargeA;
189             scale  = 1.0 - lambda;
190         }
191         else
192         {
193             charge = chargeB;
194             scale  = lambda;
195         }
196         lowiy        = 0;
197         lowiz        = 1;
198         energy_AB[q] = 0;
199         for (ix = 0; ix < et->nx; ix++)
200         {
201             mx = ix*lll[XX];
202             for (iy = lowiy; iy < et->ny; iy++)
203             {
204                 my = iy*lll[YY];
205                 if (iy >= 0)
206                 {
207                     for (n = 0; n < natoms; n++)
208                     {
209                         et->tab_xy[n] = cmul(et->eir[ix][n][XX], et->eir[iy][n][YY]);
210                     }
211                 }
212                 else
213                 {
214                     for (n = 0; n < natoms; n++)
215                     {
216                         et->tab_xy[n] = conjmul(et->eir[ix][n][XX], et->eir[-iy][n][YY]);
217                     }
218                 }
219                 for (iz = lowiz; iz < et->nz; iz++)
220                 {
221                     mz  = iz*lll[ZZ];
222                     m2  = mx*mx+my*my+mz*mz;
223                     ak  = exp(m2*factor)/m2;
224                     akv = 2.0*ak*(1.0/m2-factor);
225                     if (iz >= 0)
226                     {
227                         for (n = 0; n < natoms; n++)
228                         {
229                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n],
230                                                                     et->eir[iz][n][ZZ]));
231                         }
232                     }
233                     else
234                     {
235                         for (n = 0; n < natoms; n++)
236                         {
237                             et->tab_qxyz[n] = rcmul(charge[n], conjmul(et->tab_xy[n],
238                                                                        et->eir[-iz][n][ZZ]));
239                         }
240                     }
241
242                     cs = ss = 0;
243                     for (n = 0; n < natoms; n++)
244                     {
245                         cs += et->tab_qxyz[n].re;
246                         ss += et->tab_qxyz[n].im;
247                     }
248                     energy_AB[q]  += ak*(cs*cs+ss*ss);
249                     tmp            = scale*akv*(cs*cs+ss*ss);
250                     lrvir[XX][XX] -= tmp*mx*mx;
251                     lrvir[XX][YY] -= tmp*mx*my;
252                     lrvir[XX][ZZ] -= tmp*mx*mz;
253                     lrvir[YY][YY] -= tmp*my*my;
254                     lrvir[YY][ZZ] -= tmp*my*mz;
255                     lrvir[ZZ][ZZ] -= tmp*mz*mz;
256                     for (n = 0; n < natoms; n++)
257                     {
258                         /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
259                         tmp       = scale*ak*(cs*et->tab_qxyz[n].im-ss*et->tab_qxyz[n].re);
260                         f[n][XX] += tmp*mx*2*scaleRecip;
261                         f[n][YY] += tmp*my*2*scaleRecip;
262                         f[n][ZZ] += tmp*mz*2*scaleRecip;
263 #if 0
264                         f[n][XX] += tmp*mx;
265                         f[n][YY] += tmp*my;
266                         f[n][ZZ] += tmp*mz;
267 #endif
268                     }
269                     lowiz = 1-et->nz;
270                 }
271                 lowiy = 1-et->ny;
272             }
273         }
274     }
275
276     if (!bFreeEnergy)
277     {
278         energy = energy_AB[0];
279     }
280     else
281     {
282         energy      = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
283         *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
284     }
285
286     lrvir[XX][XX] = -0.5*scaleRecip*(lrvir[XX][XX]+energy);
287     lrvir[XX][YY] = -0.5*scaleRecip*(lrvir[XX][YY]);
288     lrvir[XX][ZZ] = -0.5*scaleRecip*(lrvir[XX][ZZ]);
289     lrvir[YY][YY] = -0.5*scaleRecip*(lrvir[YY][YY]+energy);
290     lrvir[YY][ZZ] = -0.5*scaleRecip*(lrvir[YY][ZZ]);
291     lrvir[ZZ][ZZ] = -0.5*scaleRecip*(lrvir[ZZ][ZZ]+energy);
292
293     lrvir[YY][XX] = lrvir[XX][YY];
294     lrvir[ZZ][XX] = lrvir[XX][ZZ];
295     lrvir[ZZ][YY] = lrvir[YY][ZZ];
296
297     energy *= scaleRecip;
298
299     return energy;
300 }