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