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