Simplify vsite PBC handling
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_generic_cg.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) 2013,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_cg.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/utility/fatalerror.h"
49
50 void
51 gmx_nb_generic_cg_kernel(t_nblist *                nlist,
52                          rvec *                    xx,
53                          rvec *                    ff,
54                          t_forcerec *              fr,
55                          t_mdatoms *               mdatoms,
56                          nb_kernel_data_t *        kernel_data,
57                          t_nrnb *                  nrnb)
58 {
59     int           ntype, table_nelements, ielec, ivdw;
60     real          facel;
61     int           n, is3, i3, k, nj0, nj1, j3, ggid, nnn, n0;
62     int           ai0, ai1, ai, aj0, aj1, aj;
63     real          shX, shY, shZ;
64     real          fscal, tx, ty, tz;
65     real          rinvsq;
66     real          iq;
67     real          qq, vcoul, krsq, 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 *         type;
82     t_excl *      excl;
83     real *        fshift;
84     real *        Vc;
85     real *        Vvdw;
86     real          tabscale;
87     real *        VFtab;
88     real *        x;
89     real *        f;
90     gmx_bool      do_tab;
91
92     x                   = xx[0];
93     f                   = ff[0];
94     ielec               = nlist->ielec;
95     ivdw                = nlist->ivdw;
96
97     fshift              = fr->fshift[0];
98     Vc                  = kernel_data->energygrp_elec;
99     Vvdw                = kernel_data->energygrp_vdw;
100
101     do_tab = (ielec == GMX_NBKERNEL_ELEC_CUBICSPLINETABLE ||
102               ivdw == GMX_NBKERNEL_VDW_CUBICSPLINETABLE);
103     if (do_tab)
104     {
105         tabscale         = kernel_data->table_elec_vdw->scale;
106         VFtab            = kernel_data->table_elec_vdw->data;
107     }
108     else
109     {
110         tabscale        = 0;
111         VFtab           = nullptr;
112     }
113
114     /* avoid compiler warnings for cases that cannot happen */
115     nnn                 = 0;
116     eps                 = 0.0;
117     eps2                = 0.0;
118
119     /* 3 VdW parameters for buckingham, otherwise 2 */
120     nvdwparam           = (nlist->ivdw == 2) ? 3 : 2;
121     table_nelements     = (ielec == 3) ? 4 : 0;
122     table_nelements    += (ivdw == 3) ? 8 : 0;
123
124     charge              = mdatoms->chargeA;
125     type                = mdatoms->typeA;
126     facel               = fr->ic->epsfac;
127     shiftvec            = fr->shift_vec[0];
128     vdwparam            = fr->nbfp;
129     ntype               = fr->ntype;
130
131     for (n = 0; (n < nlist->nri); n++)
132     {
133         is3              = 3*nlist->shift[n];
134         shX              = shiftvec[is3];
135         shY              = shiftvec[is3+1];
136         shZ              = shiftvec[is3+2];
137         nj0              = nlist->jindex[n];
138         nj1              = nlist->jindex[n+1];
139         ai0              = nlist->iinr[n];
140         ai1              = nlist->iinr_end[n];
141         vctot            = 0;
142         Vvdwtot          = 0;
143         fix              = 0;
144         fiy              = 0;
145         fiz              = 0;
146
147         for (k = nj0; (k < nj1); k++)
148         {
149             aj0              = nlist->jjnr[k];
150             aj1              = nlist->jjnr_end[k];
151             excl             = &nlist->excl[k*MAX_CGCGSIZE];
152
153             for (ai = ai0; (ai < ai1); ai++)
154             {
155                 i3               = ai*3;
156                 ix               = shX + x[i3+0];
157                 iy               = shY + x[i3+1];
158                 iz               = shZ + x[i3+2];
159                 iq               = facel*charge[ai];
160                 nti              = nvdwparam*ntype*type[ai];
161
162                 /* Note that this code currently calculates
163                  * all LJ and Coulomb interactions,
164                  * even if the LJ parameters or charges are zero.
165                  * If required, this can be optimized.
166                  */
167
168                 for (aj = aj0; (aj < aj1); aj++)
169                 {
170                     /* Check if this interaction is excluded */
171                     if (excl[aj-aj0] & (1<<(ai-ai0)))
172                     {
173                         continue;
174                     }
175
176                     j3               = aj*3;
177                     jx               = x[j3+0];
178                     jy               = x[j3+1];
179                     jz               = x[j3+2];
180                     dx               = ix - jx;
181                     dy               = iy - jy;
182                     dz               = iz - jz;
183                     rsq              = dx*dx+dy*dy+dz*dz;
184                     rinv             = gmx::invsqrt(rsq);
185                     rinvsq           = rinv*rinv;
186                     fscal            = 0;
187
188                     if (ielec == 3 || ivdw == 3)
189                     {
190                         r                = rsq*rinv;
191                         rt               = r*tabscale;
192                         n0               = rt;
193                         eps              = rt-n0;
194                         eps2             = eps*eps;
195                         nnn              = table_nelements*n0;
196                     }
197
198                     /* Coulomb interaction. ielec==0 means no interaction */
199                     if (ielec > 0)
200                     {
201                         qq               = iq*charge[aj];
202
203                         switch (ielec)
204                         {
205                             case 1:
206                                 /* Vanilla cutoff coulomb */
207                                 vcoul            = qq*rinv;
208                                 fscal            = vcoul*rinvsq;
209                                 break;
210
211                             case 2:
212                                 /* Reaction-field */
213                                 krsq             = fr->ic->k_rf*rsq;
214                                 vcoul            = qq*(rinv + krsq - fr->ic->c_rf);
215                                 fscal            = qq*(rinv-2.0*krsq)*rinvsq;
216                                 break;
217
218                             case 3:
219                                 /* Tabulated coulomb */
220                                 Y                = VFtab[nnn];
221                                 F                = VFtab[nnn+1];
222                                 Geps             = eps*VFtab[nnn+2];
223                                 Heps2            = eps2*VFtab[nnn+3];
224                                 nnn             += 4;
225                                 Fp               = F+Geps+Heps2;
226                                 VV               = Y+eps*Fp;
227                                 FF               = Fp+Geps+2.0*Heps2;
228                                 vcoul            = qq*VV;
229                                 fscal            = -qq*FF*tabscale*rinv;
230                                 break;
231
232                             case 4:
233                                 /* GB */
234                                 gmx_fatal(FARGS, "Death & horror! GB generic interaction not implemented.\n");
235                             default:
236                                 gmx_fatal(FARGS, "Death & horror! No generic coulomb interaction for ielec=%d.\n", ielec);
237                         }
238                         vctot            = vctot+vcoul;
239                     }  /* End of coulomb interactions */
240
241
242                     /* VdW interaction. ivdw==0 means no interaction */
243                     if (ivdw > 0)
244                     {
245                         tj               = nti+nvdwparam*type[aj];
246
247                         switch (ivdw)
248                         {
249                             case 1:
250                                 /* Vanilla Lennard-Jones cutoff */
251                                 c6               = vdwparam[tj];
252                                 c12              = vdwparam[tj+1];
253
254                                 rinvsix          = rinvsq*rinvsq*rinvsq;
255                                 Vvdw_disp        = c6*rinvsix;
256                                 Vvdw_rep         = c12*rinvsix*rinvsix;
257                                 fscal           += (12.0*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
258                                 Vvdwtot          = Vvdwtot+Vvdw_rep-Vvdw_disp;
259                                 break;
260
261                             case 2:
262                                 /* Buckingham */
263                                 c6               = vdwparam[tj];
264                                 cexp1            = vdwparam[tj+1];
265                                 cexp2            = vdwparam[tj+2];
266
267                                 rinvsix          = rinvsq*rinvsq*rinvsq;
268                                 Vvdw_disp        = c6*rinvsix;
269                                 br               = cexp2*rsq*rinv;
270                                 Vvdw_rep         = cexp1*std::exp(-br);
271                                 fscal           += (br*Vvdw_rep-6.0*Vvdw_disp)*rinvsq;
272                                 Vvdwtot          = Vvdwtot+Vvdw_rep-Vvdw_disp;
273                                 break;
274
275                             case 3:
276                                 /* Tabulated VdW */
277                                 c6               = vdwparam[tj];
278                                 c12              = vdwparam[tj+1];
279
280                                 Y                = VFtab[nnn];
281                                 F                = VFtab[nnn+1];
282                                 Geps             = eps*VFtab[nnn+2];
283                                 Heps2            = eps2*VFtab[nnn+3];
284                                 Fp               = F+Geps+Heps2;
285                                 VV               = Y+eps*Fp;
286                                 FF               = Fp+Geps+2.0*Heps2;
287                                 Vvdw_disp        = c6*VV;
288                                 fijD             = c6*FF;
289                                 nnn             += 4;
290                                 Y                = VFtab[nnn];
291                                 F                = VFtab[nnn+1];
292                                 Geps             = eps*VFtab[nnn+2];
293                                 Heps2            = eps2*VFtab[nnn+3];
294                                 Fp               = F+Geps+Heps2;
295                                 VV               = Y+eps*Fp;
296                                 FF               = Fp+Geps+2.0*Heps2;
297                                 Vvdw_rep         = c12*VV;
298                                 fijR             = c12*FF;
299                                 fscal           += -(fijD+fijR)*tabscale*rinv;
300                                 Vvdwtot          = Vvdwtot + Vvdw_disp + Vvdw_rep;
301                                 break;
302
303                             default:
304                                 gmx_fatal(FARGS, "Death & horror! No generic VdW interaction for ivdw=%d.\n", ivdw);
305                         }
306                     }  /* end VdW interactions */
307
308
309                     tx               = fscal*dx;
310                     ty               = fscal*dy;
311                     tz               = fscal*dz;
312                     f[i3+0]         += tx;
313                     f[i3+1]         += ty;
314                     f[i3+2]         += tz;
315                     f[j3+0]         -= tx;
316                     f[j3+1]         -= ty;
317                     f[j3+2]         -= tz;
318                     fix             += tx;
319                     fiy             += ty;
320                     fiz             += tz;
321                 }
322             }
323         }
324
325         fshift[is3]     += fix;
326         fshift[is3+1]   += fiy;
327         fshift[is3+2]   += fiz;
328         ggid             = nlist->gid[n];
329         Vc[ggid]        += vctot;
330         Vvdw[ggid]      += Vvdwtot;
331     }
332     /* Estimate flops, average for generic cg kernel:
333      * 12  flops per outer iteration
334      * 100 flops per inner iteration
335      */
336     inc_nrnb(nrnb, eNR_NBKERNEL_GENERIC_CG, nlist->nri*12 + nlist->jindex[n]*100);
337 }