Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_generic.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  * GROningen Mixture of Alchemy and Childrens' Stories
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40
41 #include "types/simple.h"
42 #include "vec.h"
43 #include "typedefs.h"
44 #include "nb_generic.h"
45 #include "nrnb.h"
46
47 #include "nonbonded.h"
48 #include "nb_kernel.h"
49
50
51 void
52 gmx_nb_generic_kernel(t_nblist *                nlist,
53                       rvec *                    xx,
54                       rvec *                    ff,
55                       t_forcerec *              fr,
56                       t_mdatoms *               mdatoms,
57                       nb_kernel_data_t *        kernel_data,
58                       t_nrnb *                  nrnb)
59 {
60     int           nri, ntype, table_nelements, ielec, ivdw;
61     real          facel, gbtabscale;
62     int           n, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid, nnn, n0;
63     real          shX, shY, shZ;
64     real          fscal, felec, fvdw, velec, vvdw, tx, ty, tz;
65     real          rinvsq;
66     real          iq;
67     real          qq, 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     real *        fshift;
84     real *        velecgrp;
85     real *        vvdwgrp;
86     real          tabscale;
87     real *        VFtab;
88     real *        x;
89     real *        f;
90     int           ewitab;
91     real          ewtabscale, eweps, sh_ewald, ewrt, ewtabhalfspace;
92     real *        ewtab;
93     real          rcoulomb2, rvdw, rvdw2, sh_invrc6;
94     real          rcutoff, rcutoff2;
95     real          rswitch_elec, rswitch_vdw, d, d2, sw, dsw, rinvcorr;
96     real          elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
97     real          vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
98     gmx_bool      bExactElecCutoff, bExactVdwCutoff, bExactCutoff;
99
100     x                   = xx[0];
101     f                   = ff[0];
102     ielec               = nlist->ielec;
103     ivdw                = nlist->ivdw;
104
105     fshift              = fr->fshift[0];
106     velecgrp            = kernel_data->energygrp_elec;
107     vvdwgrp             = kernel_data->energygrp_vdw;
108     tabscale            = kernel_data->table_elec_vdw->scale;
109     VFtab               = kernel_data->table_elec_vdw->data;
110
111     sh_ewald            = fr->ic->sh_ewald;
112     ewtab               = fr->ic->tabq_coul_FDV0;
113     ewtabscale          = fr->ic->tabq_scale;
114     ewtabhalfspace      = 0.5/ewtabscale;
115
116     rcoulomb2           = fr->rcoulomb*fr->rcoulomb;
117     rvdw                = fr->rvdw;
118     rvdw2               = rvdw*rvdw;
119     sh_invrc6           = fr->ic->sh_invrc6;
120
121     if (fr->coulomb_modifier == eintmodPOTSWITCH)
122     {
123         d               = fr->rcoulomb-fr->rcoulomb_switch;
124         elec_swV3       = -10.0/(d*d*d);
125         elec_swV4       =  15.0/(d*d*d*d);
126         elec_swV5       =  -6.0/(d*d*d*d*d);
127         elec_swF2       = -30.0/(d*d*d);
128         elec_swF3       =  60.0/(d*d*d*d);
129         elec_swF4       = -30.0/(d*d*d*d*d);
130     }
131     else
132     {
133         /* Avoid warnings from stupid compilers (looking at you, Clang!) */
134         elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
135     }
136     if (fr->vdw_modifier == eintmodPOTSWITCH)
137     {
138         d               = fr->rvdw-fr->rvdw_switch;
139         vdw_swV3        = -10.0/(d*d*d);
140         vdw_swV4        =  15.0/(d*d*d*d);
141         vdw_swV5        =  -6.0/(d*d*d*d*d);
142         vdw_swF2        = -30.0/(d*d*d);
143         vdw_swF3        =  60.0/(d*d*d*d);
144         vdw_swF4        = -30.0/(d*d*d*d*d);
145     }
146     else
147     {
148         /* Avoid warnings from stupid compilers (looking at you, Clang!) */
149         vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
150     }
151
152     bExactElecCutoff    = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
153     bExactVdwCutoff     = (fr->vdw_modifier != eintmodNONE);
154     bExactCutoff        = bExactElecCutoff || bExactVdwCutoff;
155
156     if (bExactCutoff)
157     {
158         rcutoff  = ( fr->rcoulomb > fr->rvdw ) ? fr->rcoulomb : fr->rvdw;
159         rcutoff2 = rcutoff*rcutoff;
160     }
161     else
162     {
163         /* Fix warnings for stupid compilers */
164         rcutoff = rcutoff2 = 1e30;
165     }
166
167     /* avoid compiler warnings for cases that cannot happen */
168     nnn                 = 0;
169     eps                 = 0.0;
170     eps2                = 0.0;
171
172     /* 3 VdW parameters for buckingham, otherwise 2 */
173     nvdwparam           = (ivdw == GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
174     table_nelements     = 12;
175
176     charge              = mdatoms->chargeA;
177     type                = mdatoms->typeA;
178     facel               = fr->epsfac;
179     shiftvec            = fr->shift_vec[0];
180     vdwparam            = fr->nbfp;
181     ntype               = fr->ntype;
182
183     for (n = 0; (n < nlist->nri); n++)
184     {
185         is3              = 3*nlist->shift[n];
186         shX              = shiftvec[is3];
187         shY              = shiftvec[is3+1];
188         shZ              = shiftvec[is3+2];
189         nj0              = nlist->jindex[n];
190         nj1              = nlist->jindex[n+1];
191         ii               = nlist->iinr[n];
192         ii3              = 3*ii;
193         ix               = shX + x[ii3+0];
194         iy               = shY + x[ii3+1];
195         iz               = shZ + x[ii3+2];
196         iq               = facel*charge[ii];
197         nti              = nvdwparam*ntype*type[ii];
198         vctot            = 0;
199         vvdwtot          = 0;
200         fix              = 0;
201         fiy              = 0;
202         fiz              = 0;
203
204         for (k = nj0; (k < nj1); k++)
205         {
206             jnr              = nlist->jjnr[k];
207             j3               = 3*jnr;
208             jx               = x[j3+0];
209             jy               = x[j3+1];
210             jz               = x[j3+2];
211             dx               = ix - jx;
212             dy               = iy - jy;
213             dz               = iz - jz;
214             rsq              = dx*dx+dy*dy+dz*dz;
215             rinv             = gmx_invsqrt(rsq);
216             rinvsq           = rinv*rinv;
217             felec            = 0;
218             fvdw             = 0;
219             velec            = 0;
220             vvdw             = 0;
221
222             if (bExactCutoff && rsq > rcutoff2)
223             {
224                 continue;
225             }
226
227             if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
228             {
229                 r                = rsq*rinv;
230                 rt               = r*tabscale;
231                 n0               = rt;
232                 eps              = rt-n0;
233                 eps2             = eps*eps;
234                 nnn              = table_nelements*n0;
235             }
236
237             /* Coulomb interaction. ielec==0 means no interaction */
238             if (ielec != GMX_NBKERNEL_ELEC_NONE)
239             {
240                 qq               = iq*charge[jnr];
241
242                 switch (ielec)
243                 {
244                     case GMX_NBKERNEL_ELEC_NONE:
245                         break;
246
247                     case GMX_NBKERNEL_ELEC_COULOMB:
248                         /* Vanilla cutoff coulomb */
249                         velec            = qq*rinv;
250                         felec            = velec*rinvsq;
251                         break;
252
253                     case GMX_NBKERNEL_ELEC_REACTIONFIELD:
254                         /* Reaction-field */
255                         velec            = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
256                         felec            = qq*(rinv*rinvsq-2.0*fr->k_rf);
257                         break;
258
259                     case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
260                         /* Tabulated coulomb */
261                         Y                = VFtab[nnn];
262                         F                = VFtab[nnn+1];
263                         Geps             = eps*VFtab[nnn+2];
264                         Heps2            = eps2*VFtab[nnn+3];
265                         Fp               = F+Geps+Heps2;
266                         VV               = Y+eps*Fp;
267                         FF               = Fp+Geps+2.0*Heps2;
268                         velec            = qq*VV;
269                         felec            = -qq*FF*tabscale*rinv;
270                         break;
271
272                     case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
273                         /* GB */
274                         gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
275                         break;
276
277                     case GMX_NBKERNEL_ELEC_EWALD:
278                         ewrt             = rsq*rinv*ewtabscale;
279                         ewitab           = ewrt;
280                         eweps            = ewrt-ewitab;
281                         ewitab           = 4*ewitab;
282                         felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
283                         rinvcorr         = (fr->coulomb_modifier == eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
284                         velec            = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
285                         felec            = qq*rinv*(rinvsq-felec);
286                         break;
287
288                     default:
289                         gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
290                         break;
291                 }
292                 if (fr->coulomb_modifier == eintmodPOTSWITCH)
293                 {
294                     d                = rsq*rinv-fr->rcoulomb_switch;
295                     d                = (d > 0.0) ? d : 0.0;
296                     d2               = d*d;
297                     sw               = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
298                     dsw              = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
299                     /* Apply switch function. Note that felec=f/r since it will be multiplied
300                      * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
301                      * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
302                      */
303                     felec            = felec*sw - rinv*velec*dsw;
304                     /* Once we have used velec to update felec we can modify velec too */
305                     velec           *= sw;
306                 }
307                 if (bExactElecCutoff)
308                 {
309                     felec            = (rsq <= rcoulomb2) ? felec : 0.0;
310                     velec            = (rsq <= rcoulomb2) ? velec : 0.0;
311                 }
312                 vctot           += velec;
313             } /* End of coulomb interactions */
314
315
316             /* VdW interaction. ivdw==0 means no interaction */
317             if (ivdw != GMX_NBKERNEL_VDW_NONE)
318             {
319                 tj               = nti+nvdwparam*type[jnr];
320
321                 switch (ivdw)
322                 {
323                     case GMX_NBKERNEL_VDW_NONE:
324                         break;
325
326                     case GMX_NBKERNEL_VDW_LENNARDJONES:
327                         /* Vanilla Lennard-Jones cutoff */
328                         c6               = vdwparam[tj];
329                         c12              = vdwparam[tj+1];
330                         rinvsix          = rinvsq*rinvsq*rinvsq;
331                         vvdw_disp        = c6*rinvsix;
332                         vvdw_rep         = c12*rinvsix*rinvsix;
333                         fvdw             = (vvdw_rep-vvdw_disp)*rinvsq;
334                         if (fr->vdw_modifier == eintmodPOTSHIFT)
335                         {
336                             vvdw             = (vvdw_rep-c12*sh_invrc6*sh_invrc6)*(1.0/12.0)-(vvdw_disp-c6*sh_invrc6)*(1.0/6.0);
337                         }
338                         else
339                         {
340                             vvdw             = vvdw_rep/12.0-vvdw_disp/6.0;
341                         }
342                         break;
343
344                     case GMX_NBKERNEL_VDW_BUCKINGHAM:
345                         /* Buckingham */
346                         c6               = vdwparam[tj];
347                         cexp1            = vdwparam[tj+1];
348                         cexp2            = vdwparam[tj+2];
349
350                         rinvsix          = rinvsq*rinvsq*rinvsq;
351                         vvdw_disp        = c6*rinvsix;
352                         br               = cexp2*rsq*rinv;
353                         vvdw_rep         = cexp1*exp(-br);
354                         fvdw             = (br*vvdw_rep-vvdw_disp)*rinvsq;
355                         if (fr->vdw_modifier == eintmodPOTSHIFT)
356                         {
357                             vvdw             = (vvdw_rep-cexp1*exp(-cexp2*rvdw))-(vvdw_disp-c6*sh_invrc6)/6.0;
358                         }
359                         else
360                         {
361                             vvdw             = vvdw_rep-vvdw_disp/6.0;
362                         }
363                         break;
364
365                     case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
366                         /* Tabulated VdW */
367                         c6               = vdwparam[tj];
368                         c12              = vdwparam[tj+1];
369                         Y                = VFtab[nnn+4];
370                         F                = VFtab[nnn+5];
371                         Geps             = eps*VFtab[nnn+6];
372                         Heps2            = eps2*VFtab[nnn+7];
373                         Fp               = F+Geps+Heps2;
374                         VV               = Y+eps*Fp;
375                         FF               = Fp+Geps+2.0*Heps2;
376                         vvdw_disp        = c6*VV;
377                         fijD             = c6*FF;
378                         Y                = VFtab[nnn+8];
379                         F                = VFtab[nnn+9];
380                         Geps             = eps*VFtab[nnn+10];
381                         Heps2            = eps2*VFtab[nnn+11];
382                         Fp               = F+Geps+Heps2;
383                         VV               = Y+eps*Fp;
384                         FF               = Fp+Geps+2.0*Heps2;
385                         vvdw_rep         = c12*VV;
386                         fijR             = c12*FF;
387                         fvdw             = -(fijD+fijR)*tabscale*rinv;
388                         vvdw             = vvdw_disp + vvdw_rep;
389                         break;
390
391                     default:
392                         gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
393                         break;
394                 }
395                 if (fr->vdw_modifier == eintmodPOTSWITCH)
396                 {
397                     d                = rsq*rinv-fr->rvdw_switch;
398                     d                = (d > 0.0) ? d : 0.0;
399                     d2               = d*d;
400                     sw               = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
401                     dsw              = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
402                     /* See coulomb interaction for the force-switch formula */
403                     fvdw             = fvdw*sw - rinv*vvdw*dsw;
404                     vvdw            *= sw;
405                 }
406                 if (bExactVdwCutoff)
407                 {
408                     fvdw             = (rsq <= rvdw2) ? fvdw : 0.0;
409                     vvdw             = (rsq <= rvdw2) ? vvdw : 0.0;
410                 }
411                 vvdwtot         += vvdw;
412             } /* end VdW interactions */
413
414             fscal            = felec+fvdw;
415
416             tx               = fscal*dx;
417             ty               = fscal*dy;
418             tz               = fscal*dz;
419             fix              = fix + tx;
420             fiy              = fiy + ty;
421             fiz              = fiz + tz;
422             f[j3+0]          = f[j3+0] - tx;
423             f[j3+1]          = f[j3+1] - ty;
424             f[j3+2]          = f[j3+2] - tz;
425         }
426
427         f[ii3+0]         = f[ii3+0] + fix;
428         f[ii3+1]         = f[ii3+1] + fiy;
429         f[ii3+2]         = f[ii3+2] + fiz;
430         fshift[is3]      = fshift[is3]+fix;
431         fshift[is3+1]    = fshift[is3+1]+fiy;
432         fshift[is3+2]    = fshift[is3+2]+fiz;
433         ggid             = nlist->gid[n];
434         velecgrp[ggid]  += vctot;
435         vvdwgrp[ggid]   += vvdwtot;
436     }
437     /* Estimate flops, average for generic kernel:
438      * 12 flops per outer iteration
439      * 50 flops per inner iteration
440      */
441     inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC, nlist->nri*12 + nlist->jindex[n]*50);
442 }