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"
56 t_complex *tab_xy, *tab_qxyz;
61 /* TODO: fix thread-safety */
63 /* the other routines are in complex.h */
64 static t_complex conjmul(t_complex a,t_complex b)
68 c.re = a.re*b.re + a.im*b.im;
69 c.im = a.im*b.re - a.re*b.im;
77 static void tabulate_eir(int natom,rvec x[],int kmax,cvec **eir,rvec lll)
82 printf("Go away! kmax = %d\n",kmax);
86 for(i=0; (i<natom); i++) {
87 for(m=0; (m<3); m++) {
92 for(m=0; (m<3); m++) {
93 eir[1][i][m].re = cos(x[i][m]*lll[m]);
94 eir[1][i][m].im = sin(x[i][m]*lll[m]);
96 for(j=2; (j<kmax); j++)
98 eir[j][i][m] = cmul(eir[j-1][i][m],eir[1][i][m]);
102 void init_ewald_tab(ewald_tab_t *et, const t_commrec *cr, const t_inputrec *ir,
109 fprintf(fp,"Will do ordinary reciprocal space Ewald sum.\n");
111 (*et)->nx = ir->nkx+1;
112 (*et)->ny = ir->nky+1;
113 (*et)->nz = ir->nkz+1;
114 (*et)->kmax = max((*et)->nx,max((*et)->ny,(*et)->nz));
116 (*et)->tab_xy = NULL;
117 (*et)->tab_qxyz = NULL;
122 real do_ewald(FILE *log, gmx_bool bVerbose,
125 real chargeA[], real chargeB[],
127 t_commrec *cr, int natoms,
128 matrix lrvir, real ewaldcoeff,
129 real lambda, real *dvdlambda,
132 real factor=-1.0/(4*ewaldcoeff*ewaldcoeff);
133 real scaleRecip =4.0*M_PI/(box[XX]*box[YY]*box[ZZ])*ONE_4PI_EPS0/ir->epsilon_r; /* 1/(Vol*e0) */
134 real *charge,energy_AB[2],energy;
136 int lowiy,lowiz,ix,iy,iz,n,q;
137 real tmp,cs,ss,ak,akv,mx,my,mz,m2,scale;
138 gmx_bool bFreeEnergy;
144 gmx_fatal(FARGS,"No parallel Ewald. Use PME instead.\n");
149 if (!et->eir) /* allocate if we need to */
151 snew(et->eir,et->kmax);
152 for(n=0;n<et->kmax;n++)
153 snew(et->eir[n],natoms);
154 snew(et->tab_xy,natoms);
155 snew(et->tab_qxyz,natoms);
158 bFreeEnergy = (ir->efep != efepNO);
163 /* make tables for the structure factor parts */
164 tabulate_eir(natoms,x,et->kmax,et->eir,lll);
166 for(q=0; q<(bFreeEnergy ? 2 : 1); q++) {
172 scale = 1.0 - lambda;
180 for(ix=0;ix<et->nx;ix++) {
182 for(iy=lowiy;iy<et->ny;iy++) {
185 for(n=0;n<natoms;n++)
186 et->tab_xy[n]=cmul(et->eir[ix][n][XX],et->eir[iy][n][YY]);
188 for(n=0;n<natoms;n++)
189 et->tab_xy[n]=conjmul(et->eir[ix][n][XX],et->eir[-iy][n][YY]);
190 for(iz=lowiz;iz<et->nz;iz++) {
192 m2=mx*mx+my*my+mz*mz;
193 ak=exp(m2*factor)/m2;
194 akv=2.0*ak*(1.0/m2-factor);
196 for(n=0;n<natoms;n++)
197 et->tab_qxyz[n]=rcmul(charge[n],cmul(et->tab_xy[n],
198 et->eir[iz][n][ZZ]));
200 for(n=0;n<natoms;n++)
201 et->tab_qxyz[n]=rcmul(charge[n],conjmul(et->tab_xy[n],
202 et->eir[-iz][n][ZZ]));
205 for(n=0;n<natoms;n++) {
206 cs+=et->tab_qxyz[n].re;
207 ss+=et->tab_qxyz[n].im;
209 energy_AB[q]+=ak*(cs*cs+ss*ss);
210 tmp=scale*akv*(cs*cs+ss*ss);
211 lrvir[XX][XX]-=tmp*mx*mx;
212 lrvir[XX][YY]-=tmp*mx*my;
213 lrvir[XX][ZZ]-=tmp*mx*mz;
214 lrvir[YY][YY]-=tmp*my*my;
215 lrvir[YY][ZZ]-=tmp*my*mz;
216 lrvir[ZZ][ZZ]-=tmp*mz*mz;
217 for(n=0;n<natoms;n++) {
218 /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
219 tmp=scale*ak*(cs*et->tab_qxyz[n].im-ss*et->tab_qxyz[n].re);
220 f[n][XX]+=tmp*mx*2*scaleRecip;
221 f[n][YY]+=tmp*my*2*scaleRecip;
222 f[n][ZZ]+=tmp*mz*2*scaleRecip;
237 energy = energy_AB[0];
239 energy = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
240 *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
243 lrvir[XX][XX]=-0.5*scaleRecip*(lrvir[XX][XX]+energy);
244 lrvir[XX][YY]=-0.5*scaleRecip*(lrvir[XX][YY]);
245 lrvir[XX][ZZ]=-0.5*scaleRecip*(lrvir[XX][ZZ]);
246 lrvir[YY][YY]=-0.5*scaleRecip*(lrvir[YY][YY]+energy);
247 lrvir[YY][ZZ]=-0.5*scaleRecip*(lrvir[YY][ZZ]);
248 lrvir[ZZ][ZZ]=-0.5*scaleRecip*(lrvir[ZZ][ZZ]+energy);
250 lrvir[YY][XX]=lrvir[XX][YY];
251 lrvir[ZZ][XX]=lrvir[XX][ZZ];
252 lrvir[ZZ][YY]=lrvir[YY][ZZ];