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