2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
46 #include "gmxcomplex.h"
49 #include "gmx_fatal.h"
59 t_complex *tab_xy, *tab_qxyz;
64 /* TODO: fix thread-safety */
66 /* the other routines are in complex.h */
67 static t_complex conjmul(t_complex a,t_complex b)
71 c.re = a.re*b.re + a.im*b.im;
72 c.im = a.im*b.re - a.re*b.im;
80 static void tabulate_eir(int natom,rvec x[],int kmax,cvec **eir,rvec lll)
85 printf("Go away! kmax = %d\n",kmax);
89 for(i=0; (i<natom); i++) {
90 for(m=0; (m<3); m++) {
95 for(m=0; (m<3); m++) {
96 eir[1][i][m].re = cos(x[i][m]*lll[m]);
97 eir[1][i][m].im = sin(x[i][m]*lll[m]);
99 for(j=2; (j<kmax); j++)
101 eir[j][i][m] = cmul(eir[j-1][i][m],eir[1][i][m]);
105 void init_ewald_tab(ewald_tab_t *et, const t_commrec *cr, const t_inputrec *ir,
112 fprintf(fp,"Will do ordinary reciprocal space Ewald sum.\n");
114 (*et)->nx = ir->nkx+1;
115 (*et)->ny = ir->nky+1;
116 (*et)->nz = ir->nkz+1;
117 (*et)->kmax = max((*et)->nx,max((*et)->ny,(*et)->nz));
119 (*et)->tab_xy = NULL;
120 (*et)->tab_qxyz = NULL;
125 real do_ewald(FILE *log, gmx_bool bVerbose,
128 real chargeA[], real chargeB[],
130 t_commrec *cr, int natoms,
131 matrix lrvir, real ewaldcoeff,
132 real lambda, real *dvdlambda,
135 real factor=-1.0/(4*ewaldcoeff*ewaldcoeff);
136 real scaleRecip =4.0*M_PI/(box[XX]*box[YY]*box[ZZ])*ONE_4PI_EPS0/ir->epsilon_r; /* 1/(Vol*e0) */
137 real *charge,energy_AB[2],energy;
139 int lowiy,lowiz,ix,iy,iz,n,q;
140 real tmp,cs,ss,ak,akv,mx,my,mz,m2,scale;
141 gmx_bool bFreeEnergy;
147 gmx_fatal(FARGS,"No parallel Ewald. Use PME instead.\n");
152 if (!et->eir) /* allocate if we need to */
154 snew(et->eir,et->kmax);
155 for(n=0;n<et->kmax;n++)
156 snew(et->eir[n],natoms);
157 snew(et->tab_xy,natoms);
158 snew(et->tab_qxyz,natoms);
161 bFreeEnergy = (ir->efep != efepNO);
166 /* make tables for the structure factor parts */
167 tabulate_eir(natoms,x,et->kmax,et->eir,lll);
169 for(q=0; q<(bFreeEnergy ? 2 : 1); q++) {
175 scale = 1.0 - lambda;
183 for(ix=0;ix<et->nx;ix++) {
185 for(iy=lowiy;iy<et->ny;iy++) {
188 for(n=0;n<natoms;n++)
189 et->tab_xy[n]=cmul(et->eir[ix][n][XX],et->eir[iy][n][YY]);
191 for(n=0;n<natoms;n++)
192 et->tab_xy[n]=conjmul(et->eir[ix][n][XX],et->eir[-iy][n][YY]);
193 for(iz=lowiz;iz<et->nz;iz++) {
195 m2=mx*mx+my*my+mz*mz;
196 ak=exp(m2*factor)/m2;
197 akv=2.0*ak*(1.0/m2-factor);
199 for(n=0;n<natoms;n++)
200 et->tab_qxyz[n]=rcmul(charge[n],cmul(et->tab_xy[n],
201 et->eir[iz][n][ZZ]));
203 for(n=0;n<natoms;n++)
204 et->tab_qxyz[n]=rcmul(charge[n],conjmul(et->tab_xy[n],
205 et->eir[-iz][n][ZZ]));
208 for(n=0;n<natoms;n++) {
209 cs+=et->tab_qxyz[n].re;
210 ss+=et->tab_qxyz[n].im;
212 energy_AB[q]+=ak*(cs*cs+ss*ss);
213 tmp=scale*akv*(cs*cs+ss*ss);
214 lrvir[XX][XX]-=tmp*mx*mx;
215 lrvir[XX][YY]-=tmp*mx*my;
216 lrvir[XX][ZZ]-=tmp*mx*mz;
217 lrvir[YY][YY]-=tmp*my*my;
218 lrvir[YY][ZZ]-=tmp*my*mz;
219 lrvir[ZZ][ZZ]-=tmp*mz*mz;
220 for(n=0;n<natoms;n++) {
221 /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
222 tmp=scale*ak*(cs*et->tab_qxyz[n].im-ss*et->tab_qxyz[n].re);
223 f[n][XX]+=tmp*mx*2*scaleRecip;
224 f[n][YY]+=tmp*my*2*scaleRecip;
225 f[n][ZZ]+=tmp*mz*2*scaleRecip;
240 energy = energy_AB[0];
242 energy = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
243 *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
246 lrvir[XX][XX]=-0.5*scaleRecip*(lrvir[XX][XX]+energy);
247 lrvir[XX][YY]=-0.5*scaleRecip*(lrvir[XX][YY]);
248 lrvir[XX][ZZ]=-0.5*scaleRecip*(lrvir[XX][ZZ]);
249 lrvir[YY][YY]=-0.5*scaleRecip*(lrvir[YY][YY]+energy);
250 lrvir[YY][ZZ]=-0.5*scaleRecip*(lrvir[YY][ZZ]);
251 lrvir[ZZ][ZZ]=-0.5*scaleRecip*(lrvir[ZZ][ZZ]+energy);
253 lrvir[YY][XX]=lrvir[XX][YY];
254 lrvir[ZZ][XX]=lrvir[XX][ZZ];
255 lrvir[ZZ][YY]=lrvir[YY][ZZ];