Move nbnxn files to nbnxm directory
[alexxy/gromacs.git] / src / gromacs / nbnxm / kernels_reference / 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,2018,2019, 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 "kernel_gpu_ref.h"
38
39 #include <cmath>
40
41 #include <algorithm>
42
43 #include "gromacs/math/functions.h"
44 #include "gromacs/math/utilities.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/mdlib/force_flags.h"
47 #include "gromacs/mdtypes/md_enums.h"
48 #include "gromacs/nbnxm/kernel_common.h"
49 #include "gromacs/nbnxm/nbnxm.h"
50 #include "gromacs/pbcutil/ishift.h"
51 #include "gromacs/utility/fatalerror.h"
52
53 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
54 static const int c_clSize          = c_nbnxnGpuClusterSize;
55
56 void
57 nbnxn_kernel_gpu_ref(const NbnxnPairlistGpu     *nbl,
58                      const nbnxn_atomdata_t     *nbat,
59                      const interaction_const_t  *iconst,
60                      rvec                       *shift_vec,
61                      int                         force_flags,
62                      int                         clearF,
63                      gmx::ArrayRef<real>         f,
64                      real  *                     fshift,
65                      real  *                     Vc,
66                      real  *                     Vvdw)
67 {
68     gmx_bool            bEner;
69     gmx_bool            bEwald;
70     const real         *Ftab = nullptr;
71     real                rcut2, rvdw2, rlist2;
72     int                 ntype;
73     real                facel;
74     int                 ish3;
75     int                 sci;
76     int                 cj4_ind0, cj4_ind1, cj4_ind;
77     int                 ci, cj;
78     int                 ic, jc, ia, ja, is, ifs, js, jfs, im, jm;
79     int                 n0;
80     int                 ggid;
81     real                shX, shY, shZ;
82     real                fscal, tx, ty, tz;
83     real                rinvsq;
84     real                iq;
85     real                qq, vcoul = 0, krsq, vctot;
86     int                 nti;
87     int                 tj;
88     real                rt, r, eps;
89     real                rinvsix;
90     real                Vvdwtot;
91     real                Vvdw_rep, Vvdw_disp;
92     real                ix, iy, iz, fix, fiy, fiz;
93     real                jx, jy, jz;
94     real                dx, dy, dz, rsq, rinv;
95     int                 int_bit;
96     real                fexcl;
97     real                c6, c12;
98     const nbnxn_excl_t *excl[2];
99
100     int                 npair_tot, npair;
101     int                 nhwu, nhwu_pruned;
102
103     if (nbl->na_ci != c_clSize)
104     {
105         gmx_fatal(FARGS, "The neighborlist cluster size in the GPU reference kernel is %d, expected it to be %d", nbl->na_ci, c_clSize);
106     }
107
108     if (clearF == enbvClearFYes)
109     {
110         clear_f(nbat, 0, f.data());
111     }
112
113     bEner = ((force_flags & GMX_FORCE_ENERGY) != 0);
114
115     bEwald = EEL_FULL(iconst->eeltype);
116     if (bEwald)
117     {
118         Ftab = iconst->tabq_coul_F;
119     }
120
121     rcut2                = iconst->rcoulomb*iconst->rcoulomb;
122     rvdw2                = iconst->rvdw*iconst->rvdw;
123
124     rlist2               = nbl->rlist*nbl->rlist;
125
126     const int  *type     = nbat->params().type.data();
127     facel                = iconst->epsfac;
128     const real *shiftvec = shift_vec[0];
129     const real *vdwparam = nbat->params().nbfp.data();
130     ntype                = nbat->params().numTypes;
131
132     const real *x        = nbat->x().data();
133
134     npair_tot   = 0;
135     nhwu        = 0;
136     nhwu_pruned = 0;
137
138     for (const nbnxn_sci_t &nbln : nbl->sci)
139     {
140         ish3             = 3*nbln.shift;
141         shX              = shiftvec[ish3];
142         shY              = shiftvec[ish3+1];
143         shZ              = shiftvec[ish3+2];
144         cj4_ind0         = nbln.cj4_ind_start;
145         cj4_ind1         = nbln.cj4_ind_end;
146         sci              = nbln.sci;
147         vctot            = 0;
148         Vvdwtot          = 0;
149
150         if (nbln.shift == CENTRAL &&
151             nbl->cj4[cj4_ind0].cj[0] == sci*c_numClPerSupercl)
152         {
153             /* we have the diagonal:
154              * add the charge self interaction energy term
155              */
156             for (im = 0; im < c_numClPerSupercl; im++)
157             {
158                 ci = sci*c_numClPerSupercl + im;
159                 for (ic = 0; ic < c_clSize; ic++)
160                 {
161                     ia     = ci*c_clSize + ic;
162                     iq     = x[ia*nbat->xstride+3];
163                     vctot += iq*iq;
164                 }
165             }
166             if (!bEwald)
167             {
168                 vctot *= -facel*0.5*iconst->c_rf;
169             }
170             else
171             {
172                 /* last factor 1/sqrt(pi) */
173                 vctot *= -facel*iconst->ewaldcoeff_q*M_1_SQRTPI;
174             }
175         }
176
177         for (cj4_ind = cj4_ind0; (cj4_ind < cj4_ind1); cj4_ind++)
178         {
179             excl[0]           = &nbl->excl[nbl->cj4[cj4_ind].imei[0].excl_ind];
180             excl[1]           = &nbl->excl[nbl->cj4[cj4_ind].imei[1].excl_ind];
181
182             for (jm = 0; jm < c_nbnxnGpuJgroupSize; jm++)
183             {
184                 cj               = nbl->cj4[cj4_ind].cj[jm];
185
186                 for (im = 0; im < c_numClPerSupercl; im++)
187                 {
188                     /* We're only using the first imask,
189                      * but here imei[1].imask is identical.
190                      */
191                     if ((nbl->cj4[cj4_ind].imei[0].imask >> (jm*c_numClPerSupercl + im)) & 1)
192                     {
193                         gmx_bool within_rlist;
194
195                         ci               = sci*c_numClPerSupercl + im;
196
197                         within_rlist     = FALSE;
198                         npair            = 0;
199                         for (ic = 0; ic < c_clSize; ic++)
200                         {
201                             ia               = ci*c_clSize + ic;
202
203                             is               = ia*nbat->xstride;
204                             ifs              = ia*nbat->fstride;
205                             ix               = shX + x[is+0];
206                             iy               = shY + x[is+1];
207                             iz               = shZ + x[is+2];
208                             iq               = facel*x[is+3];
209                             nti              = ntype*2*type[ia];
210
211                             fix              = 0;
212                             fiy              = 0;
213                             fiz              = 0;
214
215                             for (jc = 0; jc < c_clSize; jc++)
216                             {
217                                 ja               = cj*c_clSize + jc;
218
219                                 if (nbln.shift == CENTRAL &&
220                                     ci == cj && ja <= ia)
221                                 {
222                                     continue;
223                                 }
224
225                                 constexpr int clusterPerSplit = c_nbnxnGpuClusterSize/c_nbnxnGpuClusterpairSplit;
226                                 int_bit = ((excl[jc/clusterPerSplit]->pair[(jc & (clusterPerSplit - 1))*c_clSize + ic]
227                                             >> (jm*c_numClPerSupercl + im)) & 1);
228
229                                 js               = ja*nbat->xstride;
230                                 jfs              = ja*nbat->fstride;
231                                 jx               = x[js+0];
232                                 jy               = x[js+1];
233                                 jz               = x[js+2];
234                                 dx               = ix - jx;
235                                 dy               = iy - jy;
236                                 dz               = iz - jz;
237                                 rsq              = dx*dx + dy*dy + dz*dz;
238                                 if (rsq < rlist2)
239                                 {
240                                     within_rlist = TRUE;
241                                 }
242                                 if (rsq >= rcut2)
243                                 {
244                                     continue;
245                                 }
246
247                                 if (type[ia] != ntype-1 && type[ja] != ntype-1)
248                                 {
249                                     npair++;
250                                 }
251
252                                 // Ensure distance do not become so small that r^-12 overflows
253                                 rsq              = std::max(rsq, NBNXN_MIN_RSQ);
254
255                                 rinv             = gmx::invsqrt(rsq);
256                                 rinvsq           = rinv*rinv;
257
258                                 qq               = iq*x[js+3];
259                                 if (!bEwald)
260                                 {
261                                     /* Reaction-field */
262                                     krsq  = iconst->k_rf*rsq;
263                                     fscal = qq*(int_bit*rinv - 2*krsq)*rinvsq;
264                                     if (bEner)
265                                     {
266                                         vcoul = qq*(int_bit*rinv + krsq - iconst->c_rf);
267                                     }
268                                 }
269                                 else
270                                 {
271                                     r     = rsq*rinv;
272                                     rt    = r*iconst->tabq_scale;
273                                     n0    = static_cast<int>(rt);
274                                     eps   = rt - n0;
275
276                                     fexcl = (1 - eps)*Ftab[n0] + eps*Ftab[n0+1];
277
278                                     fscal = qq*(int_bit*rinvsq - fexcl)*rinv;
279
280                                     if (bEner)
281                                     {
282                                         vcoul = qq*((int_bit - std::erf(iconst->ewaldcoeff_q*r))*rinv - int_bit*iconst->sh_ewald);
283                                     }
284                                 }
285
286                                 if (rsq < rvdw2)
287                                 {
288                                     tj        = nti + 2*type[ja];
289
290                                     /* Vanilla Lennard-Jones cutoff */
291                                     c6        = vdwparam[tj];
292                                     c12       = vdwparam[tj+1];
293
294                                     rinvsix   = int_bit*rinvsq*rinvsq*rinvsq;
295                                     Vvdw_disp = c6*rinvsix;
296                                     Vvdw_rep  = c12*rinvsix*rinvsix;
297                                     fscal    += (Vvdw_rep - Vvdw_disp)*rinvsq;
298
299                                     if (bEner)
300                                     {
301                                         vctot   += vcoul;
302
303                                         Vvdwtot +=
304                                             (Vvdw_rep - int_bit*c12*iconst->sh_invrc6*iconst->sh_invrc6)/12 -
305                                             (Vvdw_disp - int_bit*c6*iconst->sh_invrc6)/6;
306                                     }
307                                 }
308
309                                 tx        = fscal*dx;
310                                 ty        = fscal*dy;
311                                 tz        = fscal*dz;
312                                 fix       = fix + tx;
313                                 fiy       = fiy + ty;
314                                 fiz       = fiz + tz;
315                                 f[jfs+0] -= tx;
316                                 f[jfs+1] -= ty;
317                                 f[jfs+2] -= tz;
318                             }
319
320                             f[ifs+0]        += fix;
321                             f[ifs+1]        += fiy;
322                             f[ifs+2]        += fiz;
323                             fshift[ish3]     = fshift[ish3]   + fix;
324                             fshift[ish3+1]   = fshift[ish3+1] + fiy;
325                             fshift[ish3+2]   = fshift[ish3+2] + fiz;
326
327                             /* Count in half work-units.
328                              * In CUDA one work-unit is 2 warps.
329                              */
330                             if ((ic+1) % (c_clSize/c_nbnxnGpuClusterpairSplit) == 0)
331                             {
332                                 npair_tot += npair;
333
334                                 nhwu++;
335                                 if (within_rlist)
336                                 {
337                                     nhwu_pruned++;
338                                 }
339
340                                 within_rlist = FALSE;
341                                 npair        = 0;
342                             }
343                         }
344                     }
345                 }
346             }
347         }
348
349         if (bEner)
350         {
351             ggid             = 0;
352             Vc[ggid]         = Vc[ggid]   + vctot;
353             Vvdw[ggid]       = Vvdw[ggid] + Vvdwtot;
354         }
355     }
356
357     if (debug)
358     {
359         fprintf(debug, "number of half %dx%d atom pairs: %d after pruning: %d fraction %4.2f\n",
360                 nbl->na_ci, nbl->na_ci,
361                 nhwu, nhwu_pruned, nhwu_pruned/static_cast<double>(nhwu));
362         fprintf(debug, "generic kernel pair interactions:            %d\n",
363                 nhwu*nbl->na_ci/2*nbl->na_ci);
364         fprintf(debug, "generic kernel post-prune pair interactions: %d\n",
365                 nhwu_pruned*nbl->na_ci/2*nbl->na_ci);
366         fprintf(debug, "generic kernel non-zero pair interactions:   %d\n",
367                 npair_tot);
368         fprintf(debug, "ratio non-zero/post-prune pair interactions: %4.2f\n",
369                 npair_tot/static_cast<double>(nhwu_pruned*gmx::exactDiv(nbl->na_ci, 2)*nbl->na_ci));
370     }
371 }