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