Redefine the default boolean type to gmx_bool.
[alexxy/gromacs.git] / src / 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 "futil.h"
46 #include "gmx_fatal.h"
47 #include "physics.h"
48 #include "coulomb.h"
49
50 #define TOL 2e-5
51
52 struct ewald_tab
53 {
54     int nx,ny,nz,kmax;
55     cvec **eir;
56     t_complex *tab_xy, *tab_qxyz;
57 }; 
58
59
60
61 /* TODO: fix thread-safety */
62
63 /* the other routines are in complex.h */
64 static t_complex conjmul(t_complex a,t_complex b)
65 {
66   t_complex c;
67   
68   c.re = a.re*b.re + a.im*b.im;
69   c.im = a.im*b.re - a.re*b.im;
70   
71   return c;
72 }
73
74   
75           
76
77 static void tabulate_eir(int natom,rvec x[],int kmax,cvec **eir,rvec lll)
78 {
79   int  i,j,m;
80   
81   if (kmax < 1) {
82       printf("Go away! kmax = %d\n",kmax);
83       exit(1);
84   }
85   
86   for(i=0; (i<natom); i++) {
87     for(m=0; (m<3); m++) {
88       eir[0][i][m].re = 1;
89       eir[0][i][m].im = 0;
90     }
91     
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]); 
95     }
96     for(j=2; (j<kmax); j++) 
97       for(m=0; (m<3); m++)
98         eir[j][i][m] = cmul(eir[j-1][i][m],eir[1][i][m]);
99   }
100 }
101
102 void init_ewald_tab(ewald_tab_t *et, const t_commrec *cr, const t_inputrec *ir,
103                     FILE *fp)
104 {
105     int n;
106
107     snew(*et, 1);
108     if (fp)
109         fprintf(fp,"Will do ordinary reciprocal space Ewald sum.\n");
110
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));
115     (*et)->eir = NULL;
116     (*et)->tab_xy = NULL;
117     (*et)->tab_qxyz = NULL;
118 }
119
120
121
122 real do_ewald(FILE *log,       gmx_bool bVerbose,
123               t_inputrec *ir,
124               rvec x[],        rvec f[],
125               real chargeA[],  real chargeB[],
126               rvec box,
127               t_commrec *cr,   int natoms,
128               matrix lrvir,    real ewaldcoeff,
129               real lambda,     real *dvdlambda,
130               ewald_tab_t et)
131 {
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;
135   rvec lll;
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;
139
140     if (cr != NULL) 
141     {
142         if (PAR(cr))
143         {
144             gmx_fatal(FARGS,"No parallel Ewald. Use PME instead.\n");
145         }
146     }
147
148
149   if (!et->eir) /* allocate if we need to */
150   {
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);
156   }
157
158   bFreeEnergy = (ir->efep != efepNO);
159
160   clear_mat(lrvir);
161   
162   calc_lll(box,lll);
163   /* make tables for the structure factor parts */
164   tabulate_eir(natoms,x,et->kmax,et->eir,lll);
165
166   for(q=0; q<(bFreeEnergy ? 2 : 1); q++) {
167     if (!bFreeEnergy) {
168       charge = chargeA;
169       scale = 1.0;
170     } else if (q==0) {
171       charge = chargeA;
172       scale = 1.0 - lambda;
173     } else {
174       charge = chargeB;
175       scale = lambda;
176     }
177     lowiy=0;
178     lowiz=1;
179     energy_AB[q]=0;
180     for(ix=0;ix<et->nx;ix++) {
181       mx=ix*lll[XX];
182       for(iy=lowiy;iy<et->ny;iy++) {
183         my=iy*lll[YY];
184         if(iy>=0) 
185           for(n=0;n<natoms;n++) 
186             et->tab_xy[n]=cmul(et->eir[ix][n][XX],et->eir[iy][n][YY]);
187         else 
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++) {
191           mz=iz*lll[ZZ];               
192           m2=mx*mx+my*my+mz*mz;
193           ak=exp(m2*factor)/m2;
194           akv=2.0*ak*(1.0/m2-factor);  
195           if(iz>=0) 
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]));
199           else 
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]));
203           
204           cs=ss=0;
205           for(n=0;n<natoms;n++) {
206             cs+=et->tab_qxyz[n].re;
207             ss+=et->tab_qxyz[n].im;
208           }
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;
223 #if 0
224             f[n][XX]+=tmp*mx;
225             f[n][YY]+=tmp*my;
226             f[n][ZZ]+=tmp*mz;
227 #endif
228           }
229           lowiz=1-et->nz;
230         }
231         lowiy=1-et->ny;
232       }
233     }
234   }
235
236   if (!bFreeEnergy) {
237     energy = energy_AB[0];
238   } else {
239     energy = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
240     *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
241   }
242
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);
249   
250   lrvir[YY][XX]=lrvir[XX][YY];
251   lrvir[ZZ][XX]=lrvir[XX][ZZ];
252   lrvir[ZZ][YY]=lrvir[YY][ZZ];
253   
254   energy*=scaleRecip;
255   
256   return energy;
257 }