54a7d3af7137cb9b41cb565b44e22faa071765a4
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_generic.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include <math.h>
42
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/typedefs.h"
46 #include "nb_generic.h"
47 #include "gromacs/legacyheaders/nrnb.h"
48
49 #include "gromacs/utility/fatalerror.h"
50
51 #include "gromacs/legacyheaders/nonbonded.h"
52 #include "nb_kernel.h"
53
54 void
55 gmx_nb_generic_kernel(t_nblist *                nlist,
56                       rvec *                    xx,
57                       rvec *                    ff,
58                       t_forcerec *              fr,
59                       t_mdatoms *               mdatoms,
60                       nb_kernel_data_t *        kernel_data,
61                       t_nrnb *                  nrnb)
62 {
63     int           nri, ntype, table_nelements, ielec, ivdw;
64     real          facel, gbtabscale;
65     int           n, ii, is3, ii3, k, nj0, nj1, jnr, j3, ggid, nnn, n0;
66     real          shX, shY, shZ;
67     real          fscal, felec, fvdw, velec, vvdw, tx, ty, tz;
68     real          rinvsq;
69     real          iq;
70     real          qq, vctot;
71     int           nti, nvdwparam;
72     int           tj;
73     real          rt, r, eps, eps2, Y, F, Geps, Heps2, VV, FF, Fp, fijD, fijR;
74     real          rinvsix;
75     real          vvdwtot;
76     real          vvdw_rep, vvdw_disp;
77     real          ix, iy, iz, fix, fiy, fiz;
78     real          jx, jy, jz;
79     real          dx, dy, dz, rsq, rinv;
80     real          c6, c12, c6grid, cexp1, cexp2, br;
81     real *        charge;
82     real *        shiftvec;
83     real *        vdwparam, *vdwgridparam;
84     int *         shift;
85     int *         type;
86     real *        fshift;
87     real *        velecgrp;
88     real *        vvdwgrp;
89     real          tabscale;
90     real *        VFtab;
91     real *        x;
92     real *        f;
93     int           ewitab;
94     real          ewtabscale, eweps, sh_ewald, ewrt, ewtabhalfspace;
95     real *        ewtab;
96     real          rcoulomb2, rvdw, rvdw2, sh_dispersion, sh_repulsion;
97     real          rcutoff, rcutoff2;
98     real          rswitch_elec, rswitch_vdw, d, d2, sw, dsw, rinvcorr;
99     real          elec_swV3, elec_swV4, elec_swV5, elec_swF2, elec_swF3, elec_swF4;
100     real          vdw_swV3, vdw_swV4, vdw_swV5, vdw_swF2, vdw_swF3, vdw_swF4;
101     real          ewclj, ewclj2, ewclj6, ewcljrsq, poly, exponent, sh_lj_ewald;
102     gmx_bool      bExactElecCutoff, bExactVdwCutoff, bExactCutoff;
103
104     x                   = xx[0];
105     f                   = ff[0];
106     ielec               = nlist->ielec;
107     ivdw                = nlist->ivdw;
108
109     fshift              = fr->fshift[0];
110     velecgrp            = kernel_data->energygrp_elec;
111     vvdwgrp             = kernel_data->energygrp_vdw;
112     tabscale            = kernel_data->table_elec_vdw->scale;
113     VFtab               = kernel_data->table_elec_vdw->data;
114
115     sh_ewald            = fr->ic->sh_ewald;
116     ewtab               = fr->ic->tabq_coul_FDV0;
117     ewtabscale          = fr->ic->tabq_scale;
118     ewtabhalfspace      = 0.5/ewtabscale;
119
120     rcoulomb2           = fr->rcoulomb*fr->rcoulomb;
121     rvdw                = fr->rvdw;
122     rvdw2               = rvdw*rvdw;
123     sh_dispersion       = fr->ic->dispersion_shift.cpot;
124     sh_repulsion        = fr->ic->repulsion_shift.cpot;
125     sh_lj_ewald         = fr->ic->sh_lj_ewald;
126
127     ewclj               = fr->ewaldcoeff_lj;
128     ewclj2              = ewclj*ewclj;
129     ewclj6              = ewclj2*ewclj2*ewclj2;
130
131     if (fr->coulomb_modifier == eintmodPOTSWITCH)
132     {
133         d               = fr->rcoulomb-fr->rcoulomb_switch;
134         elec_swV3       = -10.0/(d*d*d);
135         elec_swV4       =  15.0/(d*d*d*d);
136         elec_swV5       =  -6.0/(d*d*d*d*d);
137         elec_swF2       = -30.0/(d*d*d);
138         elec_swF3       =  60.0/(d*d*d*d);
139         elec_swF4       = -30.0/(d*d*d*d*d);
140     }
141     else
142     {
143         /* Avoid warnings from stupid compilers (looking at you, Clang!) */
144         elec_swV3 = elec_swV4 = elec_swV5 = elec_swF2 = elec_swF3 = elec_swF4 = 0.0;
145     }
146     if (fr->vdw_modifier == eintmodPOTSWITCH)
147     {
148         d               = fr->rvdw-fr->rvdw_switch;
149         vdw_swV3        = -10.0/(d*d*d);
150         vdw_swV4        =  15.0/(d*d*d*d);
151         vdw_swV5        =  -6.0/(d*d*d*d*d);
152         vdw_swF2        = -30.0/(d*d*d);
153         vdw_swF3        =  60.0/(d*d*d*d);
154         vdw_swF4        = -30.0/(d*d*d*d*d);
155     }
156     else
157     {
158         /* Avoid warnings from stupid compilers (looking at you, Clang!) */
159         vdw_swV3 = vdw_swV4 = vdw_swV5 = vdw_swF2 = vdw_swF3 = vdw_swF4 = 0.0;
160     }
161
162     bExactElecCutoff    = (fr->coulomb_modifier != eintmodNONE) || fr->eeltype == eelRF_ZERO;
163     bExactVdwCutoff     = (fr->vdw_modifier != eintmodNONE);
164     bExactCutoff        = bExactElecCutoff && bExactVdwCutoff;
165
166     if (bExactCutoff)
167     {
168         rcutoff  = ( fr->rcoulomb > fr->rvdw ) ? fr->rcoulomb : fr->rvdw;
169         rcutoff2 = rcutoff*rcutoff;
170     }
171     else
172     {
173         /* Fix warnings for stupid compilers */
174         rcutoff = rcutoff2 = 1e30;
175     }
176
177     /* avoid compiler warnings for cases that cannot happen */
178     nnn                 = 0;
179     eps                 = 0.0;
180     eps2                = 0.0;
181
182     /* 3 VdW parameters for Buckingham, otherwise 2 */
183     nvdwparam           = (ivdw == GMX_NBKERNEL_VDW_BUCKINGHAM) ? 3 : 2;
184     table_nelements     = 12;
185
186     charge              = mdatoms->chargeA;
187     type                = mdatoms->typeA;
188     facel               = fr->epsfac;
189     shiftvec            = fr->shift_vec[0];
190     vdwparam            = fr->nbfp;
191     ntype               = fr->ntype;
192     vdwgridparam        = fr->ljpme_c6grid;
193
194     for (n = 0; (n < nlist->nri); n++)
195     {
196         is3              = 3*nlist->shift[n];
197         shX              = shiftvec[is3];
198         shY              = shiftvec[is3+1];
199         shZ              = shiftvec[is3+2];
200         nj0              = nlist->jindex[n];
201         nj1              = nlist->jindex[n+1];
202         ii               = nlist->iinr[n];
203         ii3              = 3*ii;
204         ix               = shX + x[ii3+0];
205         iy               = shY + x[ii3+1];
206         iz               = shZ + x[ii3+2];
207         iq               = facel*charge[ii];
208         nti              = nvdwparam*ntype*type[ii];
209         vctot            = 0;
210         vvdwtot          = 0;
211         fix              = 0;
212         fiy              = 0;
213         fiz              = 0;
214
215         for (k = nj0; (k < nj1); k++)
216         {
217             jnr              = nlist->jjnr[k];
218             j3               = 3*jnr;
219             jx               = x[j3+0];
220             jy               = x[j3+1];
221             jz               = x[j3+2];
222             dx               = ix - jx;
223             dy               = iy - jy;
224             dz               = iz - jz;
225             rsq              = dx*dx+dy*dy+dz*dz;
226             rinv             = gmx_invsqrt(rsq);
227             rinvsq           = rinv*rinv;
228             felec            = 0;
229             fvdw             = 0;
230             velec            = 0;
231             vvdw             = 0;
232
233             if (bExactCutoff && rsq >= rcutoff2)
234             {
235                 continue;
236             }
237
238             if (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE || ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE)
239             {
240                 r                = rsq*rinv;
241                 rt               = r*tabscale;
242                 n0               = rt;
243                 eps              = rt-n0;
244                 eps2             = eps*eps;
245                 nnn              = table_nelements*n0;
246             }
247
248             /* Coulomb interaction. ielec==0 means no interaction */
249             if (ielec != GMX_NBKERNEL_ELEC_NONE)
250             {
251                 qq               = iq*charge[jnr];
252
253                 switch (ielec)
254                 {
255                     case GMX_NBKERNEL_ELEC_NONE:
256                         break;
257
258                     case GMX_NBKERNEL_ELEC_COULOMB:
259                         /* Vanilla cutoff coulomb */
260                         velec            = qq*rinv;
261                         felec            = velec*rinvsq;
262                         /* The shift for the Coulomb potential is stored in
263                          * the RF parameter c_rf, which is 0 without shift
264                          */
265                         velec           -= qq*fr->ic->c_rf;
266                         break;
267
268                     case GMX_NBKERNEL_ELEC_REACTIONFIELD:
269                         /* Reaction-field */
270                         velec            = qq*(rinv+fr->k_rf*rsq-fr->c_rf);
271                         felec            = qq*(rinv*rinvsq-2.0*fr->k_rf);
272                         break;
273
274                     case GMX_NBKERNEL_ELEC_CUBICSPLINETABLE:
275                         /* Tabulated coulomb */
276                         Y                = VFtab[nnn];
277                         F                = VFtab[nnn+1];
278                         Geps             = eps*VFtab[nnn+2];
279                         Heps2            = eps2*VFtab[nnn+3];
280                         Fp               = F+Geps+Heps2;
281                         VV               = Y+eps*Fp;
282                         FF               = Fp+Geps+2.0*Heps2;
283                         velec            = qq*VV;
284                         felec            = -qq*FF*tabscale*rinv;
285                         break;
286
287                     case GMX_NBKERNEL_ELEC_GENERALIZEDBORN:
288                         /* GB */
289                         gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
290                         break;
291
292                     case GMX_NBKERNEL_ELEC_EWALD:
293                         ewrt             = rsq*rinv*ewtabscale;
294                         ewitab           = ewrt;
295                         eweps            = ewrt-ewitab;
296                         ewitab           = 4*ewitab;
297                         felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
298                         rinvcorr         = (fr->coulomb_modifier == eintmodPOTSHIFT) ? rinv-fr->ic->sh_ewald : rinv;
299                         velec            = qq*(rinvcorr-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
300                         felec            = qq*rinv*(rinvsq-felec);
301                         break;
302
303                     default:
304                         gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
305                         break;
306                 }
307                 if (fr->coulomb_modifier == eintmodPOTSWITCH)
308                 {
309                     d                = rsq*rinv-fr->rcoulomb_switch;
310                     d                = (d > 0.0) ? d : 0.0;
311                     d2               = d*d;
312                     sw               = 1.0+d2*d*(elec_swV3+d*(elec_swV4+d*elec_swV5));
313                     dsw              = d2*(elec_swF2+d*(elec_swF3+d*elec_swF4));
314                     /* Apply switch function. Note that felec=f/r since it will be multiplied
315                      * by the i-j displacement vector. This means felec'=f'/r=-(v*sw)'/r=
316                      * -(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=felec*sw-v*dsw/r
317                      */
318                     felec            = felec*sw - rinv*velec*dsw;
319                     /* Once we have used velec to update felec we can modify velec too */
320                     velec           *= sw;
321                 }
322                 if (bExactElecCutoff)
323                 {
324                     felec            = (rsq < rcoulomb2) ? felec : 0.0;
325                     velec            = (rsq < rcoulomb2) ? velec : 0.0;
326                 }
327                 vctot           += velec;
328             } /* End of coulomb interactions */
329
330
331             /* VdW interaction. ivdw==0 means no interaction */
332             if (ivdw != GMX_NBKERNEL_VDW_NONE)
333             {
334                 tj               = nti+nvdwparam*type[jnr];
335
336                 switch (ivdw)
337                 {
338                     case GMX_NBKERNEL_VDW_NONE:
339                         break;
340
341                     case GMX_NBKERNEL_VDW_LENNARDJONES:
342                         /* Vanilla Lennard-Jones cutoff */
343                         c6               = vdwparam[tj];
344                         c12              = vdwparam[tj+1];
345                         rinvsix          = rinvsq*rinvsq*rinvsq;
346                         vvdw_disp        = c6*rinvsix;
347                         vvdw_rep         = c12*rinvsix*rinvsix;
348                         fvdw             = (vvdw_rep-vvdw_disp)*rinvsq;
349                         if (fr->vdw_modifier == eintmodPOTSHIFT)
350                         {
351                             vvdw             = (vvdw_rep + c12*sh_repulsion)/12.0 - (vvdw_disp + c6*sh_dispersion)/6.0;
352                         }
353                         else
354                         {
355                             vvdw             = vvdw_rep/12.0-vvdw_disp/6.0;
356                         }
357                         break;
358
359                     case GMX_NBKERNEL_VDW_BUCKINGHAM:
360                         /* Buckingham */
361                         c6               = vdwparam[tj];
362                         cexp1            = vdwparam[tj+1];
363                         cexp2            = vdwparam[tj+2];
364
365                         rinvsix          = rinvsq*rinvsq*rinvsq;
366                         vvdw_disp        = c6*rinvsix;
367                         br               = cexp2*rsq*rinv;
368                         vvdw_rep         = cexp1*exp(-br);
369                         fvdw             = (br*vvdw_rep-vvdw_disp)*rinvsq;
370                         if (fr->vdw_modifier == eintmodPOTSHIFT)
371                         {
372                             vvdw             = (vvdw_rep-cexp1*exp(-cexp2*rvdw))-(vvdw_disp + c6*sh_dispersion)/6.0;
373                         }
374                         else
375                         {
376                             vvdw             = vvdw_rep-vvdw_disp/6.0;
377                         }
378                         break;
379
380                     case GMX_NBKERNEL_VDW_CUBICSPLINETABLE:
381                         /* Tabulated VdW */
382                         c6               = vdwparam[tj];
383                         c12              = vdwparam[tj+1];
384                         Y                = VFtab[nnn+4];
385                         F                = VFtab[nnn+5];
386                         Geps             = eps*VFtab[nnn+6];
387                         Heps2            = eps2*VFtab[nnn+7];
388                         Fp               = F+Geps+Heps2;
389                         VV               = Y+eps*Fp;
390                         FF               = Fp+Geps+2.0*Heps2;
391                         vvdw_disp        = c6*VV;
392                         fijD             = c6*FF;
393                         Y                = VFtab[nnn+8];
394                         F                = VFtab[nnn+9];
395                         Geps             = eps*VFtab[nnn+10];
396                         Heps2            = eps2*VFtab[nnn+11];
397                         Fp               = F+Geps+Heps2;
398                         VV               = Y+eps*Fp;
399                         FF               = Fp+Geps+2.0*Heps2;
400                         vvdw_rep         = c12*VV;
401                         fijR             = c12*FF;
402                         fvdw             = -(fijD+fijR)*tabscale*rinv;
403                         vvdw             = vvdw_disp + vvdw_rep;
404                         break;
405
406
407                     case GMX_NBKERNEL_VDW_LJEWALD:
408                         /* LJ-PME */
409                         rinvsix          = rinvsq*rinvsq*rinvsq;
410                         ewcljrsq         = ewclj2*rsq;
411                         exponent         = exp(-ewcljrsq);
412                         poly             = exponent*(1.0 + ewcljrsq + ewcljrsq*ewcljrsq*0.5);
413                         c6               = vdwparam[tj];
414                         c12              = vdwparam[tj+1];
415                         c6grid           = vdwgridparam[tj];
416                         vvdw_disp        = (c6-c6grid*(1.0-poly))*rinvsix;
417                         vvdw_rep         = c12*rinvsix*rinvsix;
418                         fvdw             = (vvdw_rep - vvdw_disp - c6grid*(1.0/6.0)*exponent*ewclj6)*rinvsq;
419                         if (fr->vdw_modifier == eintmodPOTSHIFT)
420                         {
421                             vvdw             = (vvdw_rep + c12*sh_repulsion)/12.0 - (vvdw_disp + c6*sh_dispersion - c6grid*sh_lj_ewald)/6.0;
422                         }
423                         else
424                         {
425                             vvdw             = vvdw_rep/12.0-vvdw_disp/6.0;
426                         }
427                         break;
428
429                     default:
430                         gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
431                         break;
432                 }
433                 if (fr->vdw_modifier == eintmodPOTSWITCH)
434                 {
435                     d                = rsq*rinv-fr->rvdw_switch;
436                     d                = (d > 0.0) ? d : 0.0;
437                     d2               = d*d;
438                     sw               = 1.0+d2*d*(vdw_swV3+d*(vdw_swV4+d*vdw_swV5));
439                     dsw              = d2*(vdw_swF2+d*(vdw_swF3+d*vdw_swF4));
440                     /* See coulomb interaction for the force-switch formula */
441                     fvdw             = fvdw*sw - rinv*vvdw*dsw;
442                     vvdw            *= sw;
443                 }
444                 if (bExactVdwCutoff)
445                 {
446                     fvdw             = (rsq < rvdw2) ? fvdw : 0.0;
447                     vvdw             = (rsq < rvdw2) ? vvdw : 0.0;
448                 }
449                 vvdwtot         += vvdw;
450             } /* end VdW interactions */
451
452             fscal            = felec+fvdw;
453
454             tx               = fscal*dx;
455             ty               = fscal*dy;
456             tz               = fscal*dz;
457             fix              = fix + tx;
458             fiy              = fiy + ty;
459             fiz              = fiz + tz;
460             f[j3+0]          = f[j3+0] - tx;
461             f[j3+1]          = f[j3+1] - ty;
462             f[j3+2]          = f[j3+2] - tz;
463         }
464
465         f[ii3+0]         = f[ii3+0] + fix;
466         f[ii3+1]         = f[ii3+1] + fiy;
467         f[ii3+2]         = f[ii3+2] + fiz;
468         fshift[is3]      = fshift[is3]+fix;
469         fshift[is3+1]    = fshift[is3+1]+fiy;
470         fshift[is3+2]    = fshift[is3+2]+fiz;
471         ggid             = nlist->gid[n];
472         velecgrp[ggid]  += vctot;
473         vvdwgrp[ggid]   += vvdwtot;
474     }
475     /* Estimate flops, average for generic kernel:
476      * 12 flops per outer iteration
477      * 50 flops per inner iteration
478      */
479     inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC, nlist->nri*12 + nlist->jindex[n]*50);
480 }