Create fileio module
[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 "gromacs/fileio/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_inputrec *ir, FILE *fp)
112 {
113     int n;
114
115     snew(*et, 1);
116     if (fp)
117     {
118         fprintf(fp, "Will do ordinary reciprocal space Ewald sum.\n");
119     }
120
121     (*et)->nx       = ir->nkx+1;
122     (*et)->ny       = ir->nky+1;
123     (*et)->nz       = ir->nkz+1;
124     (*et)->kmax     = max((*et)->nx, max((*et)->ny, (*et)->nz));
125     (*et)->eir      = NULL;
126     (*et)->tab_xy   = NULL;
127     (*et)->tab_qxyz = NULL;
128 }
129
130
131
132 real do_ewald(t_inputrec *ir,
133               rvec x[],        rvec f[],
134               real chargeA[],  real chargeB[],
135               rvec box,
136               t_commrec *cr,   int natoms,
137               matrix lrvir,    real ewaldcoeff,
138               real lambda,     real *dvdlambda,
139               ewald_tab_t et)
140 {
141     real     factor     = -1.0/(4*ewaldcoeff*ewaldcoeff);
142     real     scaleRecip = 4.0*M_PI/(box[XX]*box[YY]*box[ZZ])*ONE_4PI_EPS0/ir->epsilon_r; /* 1/(Vol*e0) */
143     real    *charge, energy_AB[2], energy;
144     rvec     lll;
145     int      lowiy, lowiz, ix, iy, iz, n, q;
146     real     tmp, cs, ss, ak, akv, mx, my, mz, m2, scale;
147     gmx_bool bFreeEnergy;
148
149     if (cr != NULL)
150     {
151         if (PAR(cr))
152         {
153             gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
154         }
155     }
156
157
158     if (!et->eir) /* allocate if we need to */
159     {
160         snew(et->eir, et->kmax);
161         for (n = 0; n < et->kmax; n++)
162         {
163             snew(et->eir[n], natoms);
164         }
165         snew(et->tab_xy, natoms);
166         snew(et->tab_qxyz, natoms);
167     }
168
169     bFreeEnergy = (ir->efep != efepNO);
170
171     clear_mat(lrvir);
172
173     calc_lll(box, lll);
174     /* make tables for the structure factor parts */
175     tabulate_eir(natoms, x, et->kmax, et->eir, lll);
176
177     for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
178     {
179         if (!bFreeEnergy)
180         {
181             charge = chargeA;
182             scale  = 1.0;
183         }
184         else if (q == 0)
185         {
186             charge = chargeA;
187             scale  = 1.0 - lambda;
188         }
189         else
190         {
191             charge = chargeB;
192             scale  = lambda;
193         }
194         lowiy        = 0;
195         lowiz        = 1;
196         energy_AB[q] = 0;
197         for (ix = 0; ix < et->nx; ix++)
198         {
199             mx = ix*lll[XX];
200             for (iy = lowiy; iy < et->ny; iy++)
201             {
202                 my = iy*lll[YY];
203                 if (iy >= 0)
204                 {
205                     for (n = 0; n < natoms; n++)
206                     {
207                         et->tab_xy[n] = cmul(et->eir[ix][n][XX], et->eir[iy][n][YY]);
208                     }
209                 }
210                 else
211                 {
212                     for (n = 0; n < natoms; n++)
213                     {
214                         et->tab_xy[n] = conjmul(et->eir[ix][n][XX], et->eir[-iy][n][YY]);
215                     }
216                 }
217                 for (iz = lowiz; iz < et->nz; iz++)
218                 {
219                     mz  = iz*lll[ZZ];
220                     m2  = mx*mx+my*my+mz*mz;
221                     ak  = exp(m2*factor)/m2;
222                     akv = 2.0*ak*(1.0/m2-factor);
223                     if (iz >= 0)
224                     {
225                         for (n = 0; n < natoms; n++)
226                         {
227                             et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n],
228                                                                     et->eir[iz][n][ZZ]));
229                         }
230                     }
231                     else
232                     {
233                         for (n = 0; n < natoms; n++)
234                         {
235                             et->tab_qxyz[n] = rcmul(charge[n], conjmul(et->tab_xy[n],
236                                                                        et->eir[-iz][n][ZZ]));
237                         }
238                     }
239
240                     cs = ss = 0;
241                     for (n = 0; n < natoms; n++)
242                     {
243                         cs += et->tab_qxyz[n].re;
244                         ss += et->tab_qxyz[n].im;
245                     }
246                     energy_AB[q]  += ak*(cs*cs+ss*ss);
247                     tmp            = scale*akv*(cs*cs+ss*ss);
248                     lrvir[XX][XX] -= tmp*mx*mx;
249                     lrvir[XX][YY] -= tmp*mx*my;
250                     lrvir[XX][ZZ] -= tmp*mx*mz;
251                     lrvir[YY][YY] -= tmp*my*my;
252                     lrvir[YY][ZZ] -= tmp*my*mz;
253                     lrvir[ZZ][ZZ] -= tmp*mz*mz;
254                     for (n = 0; n < natoms; n++)
255                     {
256                         /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
257                         tmp       = scale*ak*(cs*et->tab_qxyz[n].im-ss*et->tab_qxyz[n].re);
258                         f[n][XX] += tmp*mx*2*scaleRecip;
259                         f[n][YY] += tmp*my*2*scaleRecip;
260                         f[n][ZZ] += tmp*mz*2*scaleRecip;
261 #if 0
262                         f[n][XX] += tmp*mx;
263                         f[n][YY] += tmp*my;
264                         f[n][ZZ] += tmp*mz;
265 #endif
266                     }
267                     lowiz = 1-et->nz;
268                 }
269                 lowiy = 1-et->ny;
270             }
271         }
272     }
273
274     if (!bFreeEnergy)
275     {
276         energy = energy_AB[0];
277     }
278     else
279     {
280         energy      = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
281         *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
282     }
283
284     lrvir[XX][XX] = -0.5*scaleRecip*(lrvir[XX][XX]+energy);
285     lrvir[XX][YY] = -0.5*scaleRecip*(lrvir[XX][YY]);
286     lrvir[XX][ZZ] = -0.5*scaleRecip*(lrvir[XX][ZZ]);
287     lrvir[YY][YY] = -0.5*scaleRecip*(lrvir[YY][YY]+energy);
288     lrvir[YY][ZZ] = -0.5*scaleRecip*(lrvir[YY][ZZ]);
289     lrvir[ZZ][ZZ] = -0.5*scaleRecip*(lrvir[ZZ][ZZ]+energy);
290
291     lrvir[YY][XX] = lrvir[XX][YY];
292     lrvir[ZZ][XX] = lrvir[XX][ZZ];
293     lrvir[ZZ][YY] = lrvir[YY][ZZ];
294
295     energy *= scaleRecip;
296
297     return energy;
298 }