Merge remote-tracking branch 'origin/release-4-6' into HEAD
[alexxy/gromacs.git] / src / gromacs / 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 #include "macros.h"
50
51 #define TOL 2e-5
52
53 struct ewald_tab
54 {
55     int nx,ny,nz,kmax;
56     cvec **eir;
57     t_complex *tab_xy, *tab_qxyz;
58 }; 
59
60
61
62 /* TODO: fix thread-safety */
63
64 /* the other routines are in complex.h */
65 static t_complex conjmul(t_complex a,t_complex b)
66 {
67   t_complex c;
68   
69   c.re = a.re*b.re + a.im*b.im;
70   c.im = a.im*b.re - a.re*b.im;
71   
72   return c;
73 }
74
75   
76           
77
78 static void tabulate_eir(int natom,rvec x[],int kmax,cvec **eir,rvec lll)
79 {
80   int  i,j,m;
81   
82   if (kmax < 1) {
83       printf("Go away! kmax = %d\n",kmax);
84       exit(1);
85   }
86   
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;
91     }
92     
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]); 
96     }
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]);
100   }
101 }
102
103 void init_ewald_tab(ewald_tab_t *et, const t_commrec *cr, const t_inputrec *ir,
104                     FILE *fp)
105 {
106     int n;
107
108     snew(*et, 1);
109     if (fp)
110         fprintf(fp,"Will do ordinary reciprocal space Ewald sum.\n");
111
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;
119 }
120
121
122
123 real do_ewald(FILE *log,       gmx_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)
132 {
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   gmx_bool bFreeEnergy;
140
141     if (cr != NULL) 
142     {
143         if (PAR(cr))
144         {
145             gmx_fatal(FARGS,"No parallel Ewald. Use PME instead.\n");
146         }
147     }
148
149
150   if (!et->eir) /* allocate if we need to */
151   {
152       snew(et->eir,et->kmax);
153       for(n=0;n<et->kmax;n++)
154           snew(et->eir[n],natoms);
155       snew(et->tab_xy,natoms);
156       snew(et->tab_qxyz,natoms);
157   }
158
159   bFreeEnergy = (ir->efep != efepNO);
160
161   clear_mat(lrvir);
162   
163   calc_lll(box,lll);
164   /* make tables for the structure factor parts */
165   tabulate_eir(natoms,x,et->kmax,et->eir,lll);
166
167   for(q=0; q<(bFreeEnergy ? 2 : 1); q++) {
168     if (!bFreeEnergy) {
169       charge = chargeA;
170       scale = 1.0;
171     } else if (q==0) {
172       charge = chargeA;
173       scale = 1.0 - lambda;
174     } else {
175       charge = chargeB;
176       scale = lambda;
177     }
178     lowiy=0;
179     lowiz=1;
180     energy_AB[q]=0;
181     for(ix=0;ix<et->nx;ix++) {
182       mx=ix*lll[XX];
183       for(iy=lowiy;iy<et->ny;iy++) {
184         my=iy*lll[YY];
185         if(iy>=0) 
186           for(n=0;n<natoms;n++) 
187             et->tab_xy[n]=cmul(et->eir[ix][n][XX],et->eir[iy][n][YY]);
188         else 
189           for(n=0;n<natoms;n++) 
190             et->tab_xy[n]=conjmul(et->eir[ix][n][XX],et->eir[-iy][n][YY]); 
191         for(iz=lowiz;iz<et->nz;iz++) {
192           mz=iz*lll[ZZ];               
193           m2=mx*mx+my*my+mz*mz;
194           ak=exp(m2*factor)/m2;
195           akv=2.0*ak*(1.0/m2-factor);  
196           if(iz>=0) 
197             for(n=0;n<natoms;n++) 
198               et->tab_qxyz[n]=rcmul(charge[n],cmul(et->tab_xy[n],
199                                                    et->eir[iz][n][ZZ]));
200           else 
201             for(n=0;n<natoms;n++) 
202               et->tab_qxyz[n]=rcmul(charge[n],conjmul(et->tab_xy[n],
203                                                       et->eir[-iz][n][ZZ]));
204           
205           cs=ss=0;
206           for(n=0;n<natoms;n++) {
207             cs+=et->tab_qxyz[n].re;
208             ss+=et->tab_qxyz[n].im;
209           }
210           energy_AB[q]+=ak*(cs*cs+ss*ss);
211           tmp=scale*akv*(cs*cs+ss*ss);         
212           lrvir[XX][XX]-=tmp*mx*mx;
213           lrvir[XX][YY]-=tmp*mx*my;
214           lrvir[XX][ZZ]-=tmp*mx*mz;
215           lrvir[YY][YY]-=tmp*my*my;
216           lrvir[YY][ZZ]-=tmp*my*mz;
217           lrvir[ZZ][ZZ]-=tmp*mz*mz;
218           for(n=0;n<natoms;n++) {
219             /*tmp=scale*ak*(cs*tab_qxyz[n].im-ss*tab_qxyz[n].re);*/
220             tmp=scale*ak*(cs*et->tab_qxyz[n].im-ss*et->tab_qxyz[n].re);
221             f[n][XX]+=tmp*mx*2*scaleRecip;
222             f[n][YY]+=tmp*my*2*scaleRecip;
223             f[n][ZZ]+=tmp*mz*2*scaleRecip;
224 #if 0
225             f[n][XX]+=tmp*mx;
226             f[n][YY]+=tmp*my;
227             f[n][ZZ]+=tmp*mz;
228 #endif
229           }
230           lowiz=1-et->nz;
231         }
232         lowiy=1-et->ny;
233       }
234     }
235   }
236
237   if (!bFreeEnergy) {
238     energy = energy_AB[0];
239   } else {
240     energy = (1.0 - lambda)*energy_AB[0] + lambda*energy_AB[1];
241     *dvdlambda += scaleRecip*(energy_AB[1] - energy_AB[0]);
242   }
243
244   lrvir[XX][XX]=-0.5*scaleRecip*(lrvir[XX][XX]+energy);
245   lrvir[XX][YY]=-0.5*scaleRecip*(lrvir[XX][YY]);
246   lrvir[XX][ZZ]=-0.5*scaleRecip*(lrvir[XX][ZZ]);
247   lrvir[YY][YY]=-0.5*scaleRecip*(lrvir[YY][YY]+energy);
248   lrvir[YY][ZZ]=-0.5*scaleRecip*(lrvir[YY][ZZ]);
249   lrvir[ZZ][ZZ]=-0.5*scaleRecip*(lrvir[ZZ][ZZ]+energy);
250   
251   lrvir[YY][XX]=lrvir[XX][YY];
252   lrvir[ZZ][XX]=lrvir[XX][ZZ];
253   lrvir[ZZ][YY]=lrvir[YY][ZZ];
254   
255   energy*=scaleRecip;
256   
257   return energy;
258 }