Fixed typo in CMakeLists for mdrun-gpu
[gromacs/qmmm-gamess-us.git] / src / mdlib / ewald.c
bloba00885fe6d4a8c5fe7cdc2e8f7ec38521976ab5a
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.
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
32 * And Hey:
33 * GROwing Monsters And Cloning Shrimps
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
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 "fftgrid.h"
47 #include "gmx_fatal.h"
48 #include "physics.h"
49 #include "coulomb.h"
51 #define TOL 2e-5
53 struct ewald_tab
55 int nx,ny,nz,kmax;
56 cvec **eir;
57 t_complex *tab_xy, *tab_qxyz;
58 };
62 /* TODO: fix thread-safety */
64 /* the other routines are in complex.h */
65 static t_complex conjmul(t_complex a,t_complex b)
67 t_complex c;
69 c.re = a.re*b.re + a.im*b.im;
70 c.im = a.im*b.re - a.re*b.im;
72 return c;
78 static void tabulate_eir(int natom,rvec x[],int kmax,cvec **eir,rvec lll)
80 int i,j,m;
82 if (kmax < 1) {
83 printf("Go away! kmax = %d\n",kmax);
84 exit(1);
87 for(i=0; (i<natom); i++) {
88 for(m=0; (m<3); m++) {
89 eir[0][i][m].re = 1;
90 eir[0][i][m].im = 0;
93 for(m=0; (m<3); m++) {
94 eir[1][i][m].re = cos(x[i][m]*lll[m]);
95 eir[1][i][m].im = sin(x[i][m]*lll[m]);
97 for(j=2; (j<kmax); j++)
98 for(m=0; (m<3); m++)
99 eir[j][i][m] = cmul(eir[j-1][i][m],eir[1][i][m]);
103 void init_ewald_tab(ewald_tab_t *et, const t_commrec *cr, const t_inputrec *ir,
104 FILE *fp)
106 int n;
108 snew(*et, 1);
109 if (fp)
110 fprintf(fp,"Will do ordinary reciprocal space Ewald sum.\n");
112 (*et)->nx = ir->nkx+1;
113 (*et)->ny = ir->nky+1;
114 (*et)->nz = ir->nkz+1;
115 (*et)->kmax = max((*et)->nx,max((*et)->ny,(*et)->nz));
116 (*et)->eir = NULL;
117 (*et)->tab_xy = NULL;
118 (*et)->tab_qxyz = NULL;
123 real do_ewald(FILE *log, bool bVerbose,
124 t_inputrec *ir,
125 rvec x[], rvec f[],
126 real chargeA[], real chargeB[],
127 rvec box,
128 t_commrec *cr, int natoms,
129 matrix lrvir, real ewaldcoeff,
130 real lambda, real *dvdlambda,
131 ewald_tab_t et)
133 real factor=-1.0/(4*ewaldcoeff*ewaldcoeff);
134 real scaleRecip =4.0*M_PI/(box[XX]*box[YY]*box[ZZ])*ONE_4PI_EPS0/ir->epsilon_r; // 1/(Vol*e0) //
135 real *charge,energy_AB[2],energy;
136 rvec lll;
137 int lowiy,lowiz,ix,iy,iz,n,q;
138 real tmp,cs,ss,ak,akv,mx,my,mz,m2,scale;
139 bool bFreeEnergy;
141 if (cr != NULL)
143 if (cr->nnodes > 1 || cr->nthreads>1)
144 gmx_fatal(FARGS,"No parallel Ewald. Use PME instead.\n");
148 if (!et->eir) /* allocate if we need to */
150 snew(et->eir,et->kmax);
151 for(n=0;n<et->kmax;n++)
152 snew(et->eir[n],natoms);
153 snew(et->tab_xy,natoms);
154 snew(et->tab_qxyz,natoms);
157 bFreeEnergy = (ir->efep != efepNO);
159 clear_mat(lrvir);
161 calc_lll(box,lll);
162 /* make tables for the structure factor parts */
163 tabulate_eir(natoms,x,et->kmax,et->eir,lll);
165 for(q=0; q<(bFreeEnergy ? 2 : 1); q++) {
166 if (!bFreeEnergy) {
167 charge = chargeA;
168 scale = 1.0;
169 } else if (q==0) {
170 charge = chargeA;
171 scale = 1.0 - lambda;
172 } else {
173 charge = chargeB;
174 scale = lambda;
176 lowiy=0;
177 lowiz=1;
178 energy_AB[q]=0;
179 for(ix=0;ix<et->nx;ix++) {
180 mx=ix*lll[XX];
181 for(iy=lowiy;iy<et->ny;iy++) {
182 my=iy*lll[YY];
183 if(iy>=0)
184 for(n=0;n<natoms;n++)
185 et->tab_xy[n]=cmul(et->eir[ix][n][XX],et->eir[iy][n][YY]);
186 else
187 for(n=0;n<natoms;n++)
188 et->tab_xy[n]=conjmul(et->eir[ix][n][XX],et->eir[-iy][n][YY]);
189 for(iz=lowiz;iz<et->nz;iz++) {
190 mz=iz*lll[ZZ];
191 m2=mx*mx+my*my+mz*mz;
192 ak=exp(m2*factor)/m2;
193 akv=2.0*ak*(1.0/m2-factor);
194 if(iz>=0)
195 for(n=0;n<natoms;n++)
196 et->tab_qxyz[n]=rcmul(charge[n],cmul(et->tab_xy[n],
197 et->eir[iz][n][ZZ]));
198 else
199 for(n=0;n<natoms;n++)
200 et->tab_qxyz[n]=rcmul(charge[n],conjmul(et->tab_xy[n],
201 et->eir[-iz][n][ZZ]));
203 cs=ss=0;
204 for(n=0;n<natoms;n++) {
205 cs+=et->tab_qxyz[n].re;
206 ss+=et->tab_qxyz[n].im;
208 energy_AB[q]+=ak*(cs*cs+ss*ss);
209 tmp=scale*akv*(cs*cs+ss*ss);
210 lrvir[XX][XX]-=tmp*mx*mx;
211 lrvir[XX][YY]-=tmp*mx*my;
212 lrvir[XX][ZZ]-=tmp*mx*mz;
213 lrvir[YY][YY]-=tmp*my*my;
214 lrvir[YY][ZZ]-=tmp*my*mz;
215 lrvir[ZZ][ZZ]-=tmp*mz*mz;
216 for(n=0;n<natoms;n++) {
217 /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
218 tmp=scale*ak*(cs*et->tab_qxyz[n].im-ss*et->tab_qxyz[n].re);
219 f[n][XX]+=tmp*mx*2*scaleRecip;
220 f[n][YY]+=tmp*my*2*scaleRecip;
221 f[n][ZZ]+=tmp*mz*2*scaleRecip;
222 #if 0
223 f[n][XX]+=tmp*mx;
224 f[n][YY]+=tmp*my;
225 f[n][ZZ]+=tmp*mz;
226 #endif
228 lowiz=1-et->nz;
230 lowiy=1-et->ny;
235 if (!bFreeEnergy) {
236 energy = energy_AB[0];
237 } else {
238 energy = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
239 *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
242 lrvir[XX][XX]=-0.5*scaleRecip*(lrvir[XX][XX]+energy);
243 lrvir[XX][YY]=-0.5*scaleRecip*(lrvir[XX][YY]);
244 lrvir[XX][ZZ]=-0.5*scaleRecip*(lrvir[XX][ZZ]);
245 lrvir[YY][YY]=-0.5*scaleRecip*(lrvir[YY][YY]+energy);
246 lrvir[YY][ZZ]=-0.5*scaleRecip*(lrvir[YY][ZZ]);
247 lrvir[ZZ][ZZ]=-0.5*scaleRecip*(lrvir[ZZ][ZZ]+energy);
249 lrvir[YY][XX]=lrvir[XX][YY];
250 lrvir[ZZ][XX]=lrvir[XX][ZZ];
251 lrvir[ZZ][YY]=lrvir[YY][ZZ];
253 energy*=scaleRecip;
255 return energy;