309d9a33e421e947944cf81a9e6697bc31d26853
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_gpu_ref.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #include "gmxpre.h"
36
37 #include "nbnxn_kernel_gpu_ref.h"
38
39 #include "config.h"
40
41 #include <cmath>
42
43 #include <algorithm>
44
45 #include "gromacs/math/functions.h"
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdlib/force.h"
49 #include "gromacs/mdlib/nb_verlet.h"
50 #include "gromacs/mdlib/nbnxn_consts.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/pbcutil/ishift.h"
53 #include "gromacs/utility/fatalerror.h"
54
55 #include "nbnxn_kernel_common.h"
56
57 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
58 static const int c_clSize          = c_nbnxnGpuClusterSize;
59
60 void
61 nbnxn_kernel_gpu_ref(const nbnxn_pairlist_t     *nbl,
62                      const nbnxn_atomdata_t     *nbat,
63                      const interaction_const_t  *iconst,
64                      rvec                       *shift_vec,
65                      int                         force_flags,
66                      int                         clearF,
67                      real  *                     f,
68                      real  *                     fshift,
69                      real  *                     Vc,
70                      real  *                     Vvdw)
71 {
72     const nbnxn_sci_t  *nbln;
73     const real         *x;
74     gmx_bool            bEner;
75     gmx_bool            bEwald;
76     const real         *Ftab = nullptr;
77     real                rcut2, rvdw2, rlist2;
78     int                 ntype;
79     real                facel;
80     int                 n;
81     int                 ish3;
82     int                 sci;
83     int                 cj4_ind0, cj4_ind1, cj4_ind;
84     int                 ci, cj;
85     int                 ic, jc, ia, ja, is, ifs, js, jfs, im, jm;
86     int                 n0;
87     int                 ggid;
88     real                shX, shY, shZ;
89     real                fscal, tx, ty, tz;
90     real                rinvsq;
91     real                iq;
92     real                qq, vcoul = 0, krsq, vctot;
93     int                 nti;
94     int                 tj;
95     real                rt, r, eps;
96     real                rinvsix;
97     real                Vvdwtot;
98     real                Vvdw_rep, Vvdw_disp;
99     real                ix, iy, iz, fix, fiy, fiz;
100     real                jx, jy, jz;
101     real                dx, dy, dz, rsq, rinv;
102     int                 int_bit;
103     real                fexcl;
104     real                c6, c12;
105     const real       *  shiftvec;
106     real       *        vdwparam;
107     int       *         type;
108     const nbnxn_excl_t *excl[2];
109
110     int                 npair_tot, npair;
111     int                 nhwu, nhwu_pruned;
112
113     if (nbl->na_ci != c_clSize)
114     {
115         gmx_fatal(FARGS, "The neighborlist cluster size in the GPU reference kernel is %d, expected it to be %d", nbl->na_ci, c_clSize);
116     }
117
118     if (clearF == enbvClearFYes)
119     {
120         clear_f(nbat, 0, f);
121     }
122
123     bEner = (force_flags & GMX_FORCE_ENERGY);
124
125     bEwald = EEL_FULL(iconst->eeltype);
126     if (bEwald)
127     {
128         Ftab = iconst->tabq_coul_F;
129     }
130
131     rcut2               = iconst->rcoulomb*iconst->rcoulomb;
132     rvdw2               = iconst->rvdw*iconst->rvdw;
133
134     rlist2              = nbl->rlist*nbl->rlist;
135
136     type                = nbat->type;
137     facel               = iconst->epsfac;
138     shiftvec            = shift_vec[0];
139     vdwparam            = nbat->nbfp;
140     ntype               = nbat->ntype;
141
142     x = nbat->x;
143
144     npair_tot   = 0;
145     nhwu        = 0;
146     nhwu_pruned = 0;
147
148     for (n = 0; n < nbl->nsci; n++)
149     {
150         nbln = &nbl->sci[n];
151
152         ish3             = 3*nbln->shift;
153         shX              = shiftvec[ish3];
154         shY              = shiftvec[ish3+1];
155         shZ              = shiftvec[ish3+2];
156         cj4_ind0         = nbln->cj4_ind_start;
157         cj4_ind1         = nbln->cj4_ind_end;
158         sci              = nbln->sci;
159         vctot            = 0;
160         Vvdwtot          = 0;
161
162         if (nbln->shift == CENTRAL &&
163             nbl->cj4[cj4_ind0].cj[0] == sci*c_numClPerSupercl)
164         {
165             /* we have the diagonal:
166              * add the charge self interaction energy term
167              */
168             for (im = 0; im < c_numClPerSupercl; im++)
169             {
170                 ci = sci*c_numClPerSupercl + im;
171                 for (ic = 0; ic < c_clSize; ic++)
172                 {
173                     ia     = ci*c_clSize + ic;
174                     iq     = x[ia*nbat->xstride+3];
175                     vctot += iq*iq;
176                 }
177             }
178             if (!bEwald)
179             {
180                 vctot *= -facel*0.5*iconst->c_rf;
181             }
182             else
183             {
184                 /* last factor 1/sqrt(pi) */
185                 vctot *= -facel*iconst->ewaldcoeff_q*M_1_SQRTPI;
186             }
187         }
188
189         for (cj4_ind = cj4_ind0; (cj4_ind < cj4_ind1); cj4_ind++)
190         {
191             excl[0]           = &nbl->excl[nbl->cj4[cj4_ind].imei[0].excl_ind];
192             excl[1]           = &nbl->excl[nbl->cj4[cj4_ind].imei[1].excl_ind];
193
194             for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
195             {
196                 cj               = nbl->cj4[cj4_ind].cj[jm];
197
198                 for (im = 0; im < c_numClPerSupercl; im++)
199                 {
200                     /* We're only using the first imask,
201                      * but here imei[1].imask is identical.
202                      */
203                     if ((nbl->cj4[cj4_ind].imei[0].imask >> (jm*c_numClPerSupercl + im)) & 1)
204                     {
205                         gmx_bool within_rlist;
206
207                         ci               = sci*c_numClPerSupercl + im;
208
209                         within_rlist     = FALSE;
210                         npair            = 0;
211                         for (ic = 0; ic < c_clSize; ic++)
212                         {
213                             ia               = ci*c_clSize + ic;
214
215                             is               = ia*nbat->xstride;
216                             ifs              = ia*nbat->fstride;
217                             ix               = shX + x[is+0];
218                             iy               = shY + x[is+1];
219                             iz               = shZ + x[is+2];
220                             iq               = facel*x[is+3];
221                             nti              = ntype*2*type[ia];
222
223                             fix              = 0;
224                             fiy              = 0;
225                             fiz              = 0;
226
227                             for (jc = 0; jc < c_clSize; jc++)
228                             {
229                                 ja               = cj*c_clSize + jc;
230
231                                 if (nbln->shift == CENTRAL &&
232                                     ci == cj && ja <= ia)
233                                 {
234                                     continue;
235                                 }
236
237                                 int_bit = ((excl[jc >> 2]->pair[(jc & 3)*c_clSize + ic] >> (jm*c_numClPerSupercl + im)) & 1);
238
239                                 js               = ja*nbat->xstride;
240                                 jfs              = ja*nbat->fstride;
241                                 jx               = x[js+0];
242                                 jy               = x[js+1];
243                                 jz               = x[js+2];
244                                 dx               = ix - jx;
245                                 dy               = iy - jy;
246                                 dz               = iz - jz;
247                                 rsq              = dx*dx + dy*dy + dz*dz;
248                                 if (rsq < rlist2)
249                                 {
250                                     within_rlist = TRUE;
251                                 }
252                                 if (rsq >= rcut2)
253                                 {
254                                     continue;
255                                 }
256
257                                 if (type[ia] != ntype-1 && type[ja] != ntype-1)
258                                 {
259                                     npair++;
260                                 }
261
262                                 // Ensure distance do not become so small that r^-12 overflows
263                                 rsq              = std::max(rsq, NBNXN_MIN_RSQ);
264
265                                 rinv             = gmx::invsqrt(rsq);
266                                 rinvsq           = rinv*rinv;
267
268                                 qq               = iq*x[js+3];
269                                 if (!bEwald)
270                                 {
271                                     /* Reaction-field */
272                                     krsq  = iconst->k_rf*rsq;
273                                     fscal = qq*(int_bit*rinv - 2*krsq)*rinvsq;
274                                     if (bEner)
275                                     {
276                                         vcoul = qq*(int_bit*rinv + krsq - iconst->c_rf);
277                                     }
278                                 }
279                                 else
280                                 {
281                                     r     = rsq*rinv;
282                                     rt    = r*iconst->tabq_scale;
283                                     n0    = rt;
284                                     eps   = rt - n0;
285
286                                     fexcl = (1 - eps)*Ftab[n0] + eps*Ftab[n0+1];
287
288                                     fscal = qq*(int_bit*rinvsq - fexcl)*rinv;
289
290                                     if (bEner)
291                                     {
292                                         vcoul = qq*((int_bit - std::erf(iconst->ewaldcoeff_q*r))*rinv - int_bit*iconst->sh_ewald);
293                                     }
294                                 }
295
296                                 if (rsq < rvdw2)
297                                 {
298                                     tj        = nti + 2*type[ja];
299
300                                     /* Vanilla Lennard-Jones cutoff */
301                                     c6        = vdwparam[tj];
302                                     c12       = vdwparam[tj+1];
303
304                                     rinvsix   = int_bit*rinvsq*rinvsq*rinvsq;
305                                     Vvdw_disp = c6*rinvsix;
306                                     Vvdw_rep  = c12*rinvsix*rinvsix;
307                                     fscal    += (Vvdw_rep - Vvdw_disp)*rinvsq;
308
309                                     if (bEner)
310                                     {
311                                         vctot   += vcoul;
312
313                                         Vvdwtot +=
314                                             (Vvdw_rep - int_bit*c12*iconst->sh_invrc6*iconst->sh_invrc6)/12 -
315                                             (Vvdw_disp - int_bit*c6*iconst->sh_invrc6)/6;
316                                     }
317                                 }
318
319                                 tx        = fscal*dx;
320                                 ty        = fscal*dy;
321                                 tz        = fscal*dz;
322                                 fix       = fix + tx;
323                                 fiy       = fiy + ty;
324                                 fiz       = fiz + tz;
325                                 f[jfs+0] -= tx;
326                                 f[jfs+1] -= ty;
327                                 f[jfs+2] -= tz;
328                             }
329
330                             f[ifs+0]        += fix;
331                             f[ifs+1]        += fiy;
332                             f[ifs+2]        += fiz;
333                             fshift[ish3]     = fshift[ish3]   + fix;
334                             fshift[ish3+1]   = fshift[ish3+1] + fiy;
335                             fshift[ish3+2]   = fshift[ish3+2] + fiz;
336
337                             /* Count in half work-units.
338                              * In CUDA one work-unit is 2 warps.
339                              */
340                             if ((ic+1) % (c_clSize/c_nbnxnGpuClusterpairSplit) == 0)
341                             {
342                                 npair_tot += npair;
343
344                                 nhwu++;
345                                 if (within_rlist)
346                                 {
347                                     nhwu_pruned++;
348                                 }
349
350                                 within_rlist = FALSE;
351                                 npair        = 0;
352                             }
353                         }
354                     }
355                 }
356             }
357         }
358
359         if (bEner)
360         {
361             ggid             = 0;
362             Vc[ggid]         = Vc[ggid]   + vctot;
363             Vvdw[ggid]       = Vvdw[ggid] + Vvdwtot;
364         }
365     }
366
367     if (debug)
368     {
369         fprintf(debug, "number of half %dx%d atom pairs: %d after pruning: %d fraction %4.2f\n",
370                 nbl->na_ci, nbl->na_ci,
371                 nhwu, nhwu_pruned, nhwu_pruned/(double)nhwu);
372         fprintf(debug, "generic kernel pair interactions:            %d\n",
373                 nhwu*nbl->na_ci/2*nbl->na_ci);
374         fprintf(debug, "generic kernel post-prune pair interactions: %d\n",
375                 nhwu_pruned*nbl->na_ci/2*nbl->na_ci);
376         fprintf(debug, "generic kernel non-zero pair interactions:   %d\n",
377                 npair_tot);
378         fprintf(debug, "ratio non-zero/post-prune pair interactions: %4.2f\n",
379                 npair_tot/(double)(nhwu_pruned*nbl->na_ci/2*nbl->na_ci));
380     }
381 }