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