Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_generic_cg.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  * 
4  *                This source code is part of
5  * 
6  *                 G   R   O   M   A   C   S
7  * 
8  *          GROningen MAchine for Chemical Simulations
9  * 
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  * 
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  * 
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  * 
31  * For more info, check our website at http://www.gromacs.org
32  * 
33  * And Hey:
34  * GROningen Mixture of Alchemy and Childrens' Stories
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <math.h>
41
42 #include "types/simple.h"
43 #include "vec.h"
44 #include "typedefs.h"
45 #include "nb_generic_cg.h"
46 #include "nonbonded.h"
47 #include "nb_kernel.h"
48 #include "nrnb.h"
49
50  void
51 gmx_nb_generic_cg_kernel(t_nblist *                nlist,
52                          rvec *                    xx,
53                          rvec *                    ff,
54                          t_forcerec *              fr,
55                          t_mdatoms *               mdatoms,
56                          nb_kernel_data_t *        kernel_data,
57                          t_nrnb *                  nrnb)
58 {
59      int           nri,ntype,table_nelements,ielec,ivdw;
60      real          facel,gbtabscale;
61      int           n,is3,i3,k,nj0,nj1,j3,ggid,nnn,n0;
62      int           ai0,ai1,ai,aj0,aj1,aj;
63      real          shX,shY,shZ;
64      real          fscal,tx,ty,tz;
65      real          rinvsq;
66      real          iq;
67      real          qq,vcoul,krsq,vctot;
68      int           nti,nvdwparam;
69      int           tj;
70      real          rt,r,eps,eps2,Y,F,Geps,Heps2,VV,FF,Fp,fijD,fijR;
71      real          rinvsix;
72      real          Vvdwtot;
73      real          Vvdw_rep,Vvdw_disp;
74      real          ix,iy,iz,fix,fiy,fiz;
75      real          jx,jy,jz;
76      real          dx,dy,dz,rsq,rinv;
77      real          c6,c12,cexp1,cexp2,br;
78      real *        charge;
79      real *        shiftvec;
80      real *        vdwparam;
81      int *         shift;
82      int *         type;
83      t_excl *      excl;
84      real *        fshift;
85      real *        Vc;
86      real *        Vvdw;
87      real          tabscale;
88      real *        VFtab;
89      real *        x;
90      real *        f;
91
92      x                   = xx[0];
93      f                   = ff[0];
94      ielec               = nlist->ielec;
95      ivdw                = nlist->ivdw;
96
97      fshift              = fr->fshift[0];
98      Vc                  = kernel_data->energygrp_elec;
99      Vvdw                = kernel_data->energygrp_vdw;
100      tabscale            = kernel_data->table_elec_vdw->scale;
101      VFtab               = kernel_data->table_elec_vdw->data;
102
103      /* avoid compiler warnings for cases that cannot happen */
104      nnn                 = 0;
105      vcoul               = 0.0;
106      eps                 = 0.0;
107      eps2                = 0.0;
108      
109      /* 3 VdW parameters for buckingham, otherwise 2 */
110      nvdwparam           = (nlist->ivdw==2) ? 3 : 2;
111      table_nelements     = (ielec==3) ? 4 : 0;
112      table_nelements    += (ivdw==3) ? 8 : 0;
113      
114      charge              = mdatoms->chargeA;
115      type                = mdatoms->typeA;
116      facel               = fr->epsfac;
117      shiftvec            = fr->shift_vec[0];
118      vdwparam            = fr->nbfp;
119      ntype               = fr->ntype;
120      
121      for(n=0; (n<nlist->nri); n++)
122      {
123          is3              = 3*nlist->shift[n];     
124          shX              = shiftvec[is3];  
125          shY              = shiftvec[is3+1];
126          shZ              = shiftvec[is3+2];
127          nj0              = nlist->jindex[n];      
128          nj1              = nlist->jindex[n+1];    
129          ai0              = nlist->iinr[n];
130          ai1              = nlist->iinr_end[n];
131          vctot            = 0;              
132          Vvdwtot          = 0;              
133          fix              = 0;
134          fiy              = 0;
135          fiz              = 0;
136          
137          for(k=nj0; (k<nj1); k++)
138          {
139              aj0              = nlist->jjnr[k];
140              aj1              = nlist->jjnr_end[k];
141              excl             = &nlist->excl[k*MAX_CGCGSIZE];
142
143              for(ai=ai0; (ai<ai1); ai++)
144              {
145                  i3               = ai*3;
146                  ix               = shX + x[i3+0];
147                  iy               = shY + x[i3+1];
148                  iz               = shZ + x[i3+2];
149                  iq               = facel*charge[ai];
150                  nti              = nvdwparam*ntype*type[ai];
151
152                  /* Note that this code currently calculates
153                   * all LJ and Coulomb interactions,
154                   * even if the LJ parameters or charges are zero.
155                   * If required, this can be optimized.
156                   */
157
158                  for(aj=aj0; (aj<aj1); aj++)
159                  {
160                      /* Check if this interaction is excluded */
161                      if (excl[aj-aj0] & (1<<(ai-ai0)))
162                      {
163                          continue;
164                      }
165
166                      j3               = aj*3;
167                      jx               = x[j3+0];
168                      jy               = x[j3+1];
169                      jz               = x[j3+2];
170                      dx               = ix - jx;      
171                      dy               = iy - jy;      
172                      dz               = iz - jz;      
173                      rsq              = dx*dx+dy*dy+dz*dz;
174                      rinv             = gmx_invsqrt(rsq);
175                      rinvsq           = rinv*rinv;  
176                      fscal            = 0;
177
178                      if (ielec==3 || ivdw==3)
179                      {
180                          r                = rsq*rinv;
181                          rt               = r*tabscale;     
182                          n0               = rt;             
183                          eps              = rt-n0;          
184                          eps2             = eps*eps;        
185                          nnn              = table_nelements*n0;                                         
186                      }
187                      
188                      /* Coulomb interaction. ielec==0 means no interaction */
189                      if (ielec > 0)
190                      {
191                          qq               = iq*charge[aj]; 
192                          
193                          switch(ielec)
194                          {
195                          case 1:
196                              /* Vanilla cutoff coulomb */
197                              vcoul            = qq*rinv;      
198                              fscal            = vcoul*rinvsq; 
199                              break;
200                              
201                          case 2:
202                              /* Reaction-field */
203                              krsq             = fr->k_rf*rsq;      
204                              vcoul            = qq*(rinv+krsq-fr->c_rf);
205                              fscal            = qq*(rinv-2.0*krsq)*rinvsq;
206                              break;
207                              
208                          case 3:
209                              /* Tabulated coulomb */
210                              Y                = VFtab[nnn];     
211                              F                = VFtab[nnn+1];   
212                              Geps             = eps*VFtab[nnn+2];
213                              Heps2            = eps2*VFtab[nnn+3];
214                              nnn             += 4;
215                              Fp               = F+Geps+Heps2;   
216                              VV               = Y+eps*Fp;       
217                              FF               = Fp+Geps+2.0*Heps2;
218                              vcoul            = qq*VV;          
219                              fscal            = -qq*FF*tabscale*rinv;
220                              break;
221                              
222                          case 4:
223                              /* GB */
224                            gmx_fatal(FARGS,"Death & horror! GB generic interaction not implemented.\n");
225                            break;
226                            
227                          default:
228                              gmx_fatal(FARGS,"Death & horror! No generic coulomb interaction for ielec=%d.\n",ielec);
229                              break;
230                          }
231                          vctot            = vctot+vcoul;    
232                      } /* End of coulomb interactions */
233                      
234                      
235                      /* VdW interaction. ivdw==0 means no interaction */
236                      if (ivdw > 0)
237                      {
238                          tj               = nti+nvdwparam*type[aj];
239                          
240                          switch(ivdw)
241                          {
242                          case 1:
243                              /* Vanilla Lennard-Jones cutoff */
244                              c6               = vdwparam[tj];   
245                              c12              = vdwparam[tj+1]; 
246                              
247                              rinvsix          = rinvsq*rinvsq*rinvsq;
248                              Vvdw_disp        = c6*rinvsix;     
249                              Vvdw_rep         = c12*rinvsix*rinvsix;
250                              fscal           += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
251                              Vvdwtot          = Vvdwtot+Vvdw_rep-Vvdw_disp;
252                              break;
253                              
254                          case 2:
255                              /* Buckingham */
256                              c6               = vdwparam[tj];   
257                              cexp1            = vdwparam[tj+1]; 
258                              cexp2            = vdwparam[tj+2]; 
259                              
260                              rinvsix          = rinvsq*rinvsq*rinvsq;
261                              Vvdw_disp        = c6*rinvsix;     
262                              br               = cexp2*rsq*rinv;
263                              Vvdw_rep         = cexp1*exp(-br); 
264                              fscal           += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
265                              Vvdwtot          = Vvdwtot+Vvdw_rep-Vvdw_disp;
266                              break;
267                              
268                          case 3:
269                              /* Tabulated VdW */
270                              c6               = vdwparam[tj];   
271                              c12              = vdwparam[tj+1]; 
272                              
273                              Y                = VFtab[nnn];     
274                              F                = VFtab[nnn+1];   
275                              Geps             = eps*VFtab[nnn+2];
276                              Heps2            = eps2*VFtab[nnn+3];
277                              Fp               = F+Geps+Heps2;   
278                              VV               = Y+eps*Fp;       
279                              FF               = Fp+Geps+2.0*Heps2;
280                              Vvdw_disp        = c6*VV;          
281                              fijD             = c6*FF;          
282                              nnn             += 4;          
283                              Y                = VFtab[nnn];     
284                              F                = VFtab[nnn+1];   
285                              Geps             = eps*VFtab[nnn+2];
286                              Heps2            = eps2*VFtab[nnn+3];
287                              Fp               = F+Geps+Heps2;   
288                              VV               = Y+eps*Fp;       
289                              FF               = Fp+Geps+2.0*Heps2;
290                              Vvdw_rep         = c12*VV;         
291                              fijR             = c12*FF;         
292                              fscal           += -(fijD+fijR)*tabscale*rinv;
293                              Vvdwtot          = Vvdwtot + Vvdw_disp + Vvdw_rep;                                         
294                              break;
295                              
296                          default:
297                              gmx_fatal(FARGS,"Death & horror! No generic VdW interaction for ivdw=%d.\n",ivdw);
298                              break;
299                          }
300                      } /* end VdW interactions */
301                      
302                      
303                      tx               = fscal*dx;     
304                      ty               = fscal*dy;     
305                      tz               = fscal*dz;     
306                      f[i3+0]         += tx;
307                      f[i3+1]         += ty;
308                      f[i3+2]         += tz;
309                      f[j3+0]         -= tx;
310                      f[j3+1]         -= ty;
311                      f[j3+2]         -= tz;
312                      fix             += tx;
313                      fiy             += ty;
314                      fiz             += tz;
315                  }
316              }
317          }
318
319          fshift[is3]     += fix;
320          fshift[is3+1]   += fiy;
321          fshift[is3+2]   += fiz;
322          ggid             = nlist->gid[n];         
323          Vc[ggid]        += vctot;
324          Vvdw[ggid]      += Vvdwtot;
325      }
326     /* Estimate flops, average for generic cg kernel:
327      * 12  flops per outer iteration
328      * 100 flops per inner iteration
329      */
330     inc_nrnb(nrnb,eNR_NBKERNEL_FREE_ENERGY,nlist->nri*12 + nlist->jindex[n]*100);
331 }
332