3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * GROwing Monsters And Cloning Shrimps
43 #include "gmxcomplex.h"
46 #include "gmx_fatal.h"
57 t_complex *tab_xy, *tab_qxyz;
62 /* TODO: fix thread-safety */
64 /* the other routines are in complex.h */
65 static t_complex conjmul(t_complex a, t_complex b)
69 c.re = a.re*b.re + a.im*b.im;
70 c.im = a.im*b.re - a.re*b.im;
78 static void tabulate_eir(int natom, rvec x[], int kmax, cvec **eir, rvec lll)
84 printf("Go away! kmax = %d\n", kmax);
88 for (i = 0; (i < natom); i++)
90 for (m = 0; (m < 3); m++)
96 for (m = 0; (m < 3); m++)
98 eir[1][i][m].re = cos(x[i][m]*lll[m]);
99 eir[1][i][m].im = sin(x[i][m]*lll[m]);
101 for (j = 2; (j < kmax); j++)
103 for (m = 0; (m < 3); m++)
105 eir[j][i][m] = cmul(eir[j-1][i][m], eir[1][i][m]);
111 void init_ewald_tab(ewald_tab_t *et, const t_inputrec *ir, FILE *fp)
118 fprintf(fp, "Will do ordinary reciprocal space Ewald sum.\n");
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));
126 (*et)->tab_xy = NULL;
127 (*et)->tab_qxyz = NULL;
132 real do_ewald(t_inputrec *ir,
134 real chargeA[], real chargeB[],
136 t_commrec *cr, int natoms,
137 matrix lrvir, real ewaldcoeff,
138 real lambda, real *dvdlambda,
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;
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;
153 gmx_fatal(FARGS, "No parallel Ewald. Use PME instead.\n");
158 if (!et->eir) /* allocate if we need to */
160 snew(et->eir, et->kmax);
161 for (n = 0; n < et->kmax; n++)
163 snew(et->eir[n], natoms);
165 snew(et->tab_xy, natoms);
166 snew(et->tab_qxyz, natoms);
169 bFreeEnergy = (ir->efep != efepNO);
174 /* make tables for the structure factor parts */
175 tabulate_eir(natoms, x, et->kmax, et->eir, lll);
177 for (q = 0; q < (bFreeEnergy ? 2 : 1); q++)
187 scale = 1.0 - lambda;
197 for (ix = 0; ix < et->nx; ix++)
200 for (iy = lowiy; iy < et->ny; iy++)
205 for (n = 0; n < natoms; n++)
207 et->tab_xy[n] = cmul(et->eir[ix][n][XX], et->eir[iy][n][YY]);
212 for (n = 0; n < natoms; n++)
214 et->tab_xy[n] = conjmul(et->eir[ix][n][XX], et->eir[-iy][n][YY]);
217 for (iz = lowiz; iz < et->nz; iz++)
220 m2 = mx*mx+my*my+mz*mz;
221 ak = exp(m2*factor)/m2;
222 akv = 2.0*ak*(1.0/m2-factor);
225 for (n = 0; n < natoms; n++)
227 et->tab_qxyz[n] = rcmul(charge[n], cmul(et->tab_xy[n],
228 et->eir[iz][n][ZZ]));
233 for (n = 0; n < natoms; n++)
235 et->tab_qxyz[n] = rcmul(charge[n], conjmul(et->tab_xy[n],
236 et->eir[-iz][n][ZZ]));
241 for (n = 0; n < natoms; n++)
243 cs += et->tab_qxyz[n].re;
244 ss += et->tab_qxyz[n].im;
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++)
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;
276 energy = energy_AB[0];
280 energy = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
281 *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
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);
291 lrvir[YY][XX] = lrvir[XX][YY];
292 lrvir[ZZ][XX] = lrvir[XX][ZZ];
293 lrvir[ZZ][YY] = lrvir[YY][ZZ];
295 energy *= scaleRecip;