Improved CUDA non-bonded kernel performance
[alexxy/gromacs.git] / src / mdlib / nbnxn_cuda / nbnxn_cuda_kernel.cuh
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-2012, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38
39 #include "maths.h"
40 /* Note that floating-point constants in CUDA code should be suffixed
41  * with f (e.g. 0.5f), to stop the compiler producing intermediate
42  * code that is in double precision.
43  */
44
45 #if __CUDA_ARCH__ >= 300
46 #define REDUCE_SHUFFLE
47 /* On Kepler pre-loading i-atom types to shmem gives a few %,
48    but on Fermi it does not */
49 #define IATYPE_SHMEM
50 #endif
51
52 #if defined EL_EWALD_ANA || defined EL_EWALD_TAB
53 /* Note: convenience macro, needs to be undef-ed at the end of the file. */
54 #define EL_EWALD_ANY
55 #endif
56
57 /*
58    Kernel launch parameters:
59     - #blocks   = #pair lists, blockId = pair list Id
60     - #threads  = CL_SIZE^2
61     - shmem     = CL_SIZE^2 * sizeof(float)
62
63     Each thread calculates an i force-component taking one pair of i-j atoms.
64  */
65 #if __CUDA_ARCH__ >= 350
66 __launch_bounds__(64,16)
67 #endif
68 #ifdef PRUNE_NBL
69 #ifdef CALC_ENERGIES
70 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_prune)
71 #else
72 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _prune)
73 #endif
74 #else
75 #ifdef CALC_ENERGIES
76 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener)
77 #else
78 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn)
79 #endif
80 #endif
81             (const cu_atomdata_t atdat,
82              const cu_nbparam_t nbparam,
83              const cu_plist_t plist,
84              bool bCalcFshift)
85 {
86     /* convenience variables */
87     const nbnxn_sci_t *pl_sci   = plist.sci;
88 #ifndef PRUNE_NBL
89     const
90 #endif
91     nbnxn_cj4_t *pl_cj4         = plist.cj4;
92     const nbnxn_excl_t *excl    = plist.excl;
93     const int *atom_types       = atdat.atom_types;
94     int ntypes                  = atdat.ntypes;
95     const float4 *xq            = atdat.xq;
96     float3 *f                   = atdat.f;
97     const float3 *shift_vec     = atdat.shift_vec;
98     float rcoulomb_sq           = nbparam.rcoulomb_sq;
99 #ifdef VDW_CUTOFF_CHECK
100     float rvdw_sq               = nbparam.rvdw_sq;
101     float vdw_in_range;
102 #endif
103 #ifdef EL_RF
104     float two_k_rf              = nbparam.two_k_rf;
105 #endif
106 #ifdef EL_EWALD_TAB
107     float coulomb_tab_scale     = nbparam.coulomb_tab_scale;
108 #endif
109 #ifdef EL_EWALD_ANA
110     float beta2                 = nbparam.ewald_beta*nbparam.ewald_beta;
111     float beta3                 = nbparam.ewald_beta*nbparam.ewald_beta*nbparam.ewald_beta;
112 #endif
113 #ifdef PRUNE_NBL
114     float rlist_sq              = nbparam.rlist_sq;
115 #endif
116
117 #ifdef CALC_ENERGIES
118     float lj_shift    = nbparam.sh_invrc6;
119 #ifdef EL_EWALD_ANY
120     float beta        = nbparam.ewald_beta;
121     float ewald_shift = nbparam.sh_ewald;
122 #else
123     float c_rf        = nbparam.c_rf;
124 #endif
125     float *e_lj       = atdat.e_lj;
126     float *e_el       = atdat.e_el;
127 #endif
128
129     /* thread/block/warp id-s */
130     unsigned int tidxi  = threadIdx.x;
131     unsigned int tidxj  = threadIdx.y;
132     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
133     unsigned int bidx   = blockIdx.x;
134     unsigned int widx   = tidx / WARP_SIZE; /* warp index */
135
136     int sci, ci, cj, ci_offset,
137         ai, aj,
138         cij4_start, cij4_end,
139         typei, typej,
140         i, jm, j4, wexcl_idx;
141     float qi, qj_f,
142           r2, inv_r, inv_r2, inv_r6,
143           c6, c12,
144           int_bit,
145 #ifdef CALC_ENERGIES
146           E_lj, E_el, E_lj_p,
147 #endif
148           F_invr;
149     unsigned int wexcl, imask, mask_ji;
150     float4 xqbuf;
151     float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
152     float3 fci_buf[NCL_PER_SUPERCL];    /* i force buffer */
153     nbnxn_sci_t nb_sci;
154
155     /* shmem buffer for i x+q pre-loading */
156     extern __shared__  float4 xqib[];
157     /* shmem buffer for cj, for both warps separately */
158     int *cjs     = (int *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
159 #ifdef IATYPE_SHMEM
160     /* shmem buffer for i atom-type pre-loading */
161     int *atib = (int *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
162 #endif
163
164 #ifndef REDUCE_SHUFFLE
165     /* shmem j force buffer */
166 #ifdef IATYPE_SHMEM
167     float *f_buf = (float *)(atib + NCL_PER_SUPERCL * CL_SIZE);
168 #else
169     float *f_buf = (float *)(cjs + 2 * NBNXN_GPU_JGROUP_SIZE);
170 #endif
171 #endif
172
173     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
174     sci         = nb_sci.sci;           /* super-cluster */
175     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
176     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
177
178     /* Pre-load i-atom x and q into shared memory */
179     ci = sci * NCL_PER_SUPERCL + tidxj;
180     ai = ci * CL_SIZE + tidxi;
181     xqib[tidxj * CL_SIZE + tidxi] = xq[ai] + shift_vec[nb_sci.shift];
182 #ifdef IATYPE_SHMEM
183     /* Pre-load the i-atom types into shared memory */
184     atib[tidxj * CL_SIZE + tidxi] = atom_types[ai];
185 #endif
186     __syncthreads();
187
188     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
189     {
190         fci_buf[ci_offset] = make_float3(0.0f);
191     }
192
193 #ifdef CALC_ENERGIES
194     E_lj = 0.0f;
195     E_el = 0.0f;
196
197 #if defined EL_EWALD_ANY || defined EL_RF
198     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
199     {
200         /* we have the diagonal: add the charge self interaction energy term */
201         for (i = 0; i < NCL_PER_SUPERCL; i++)
202         {
203             qi    = xqib[i * CL_SIZE + tidxi].w;
204             E_el += qi*qi;
205         }
206         /* divide the self term equally over the j-threads */
207         E_el /= CL_SIZE;
208 #ifdef EL_RF
209         E_el *= -nbparam.epsfac*0.5f*c_rf;
210 #else
211         E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
212 #endif
213     }
214 #endif
215 #endif
216
217     /* skip central shifts when summing shift forces */
218     if (nb_sci.shift == CENTRAL)
219     {
220         bCalcFshift = false;
221     }
222
223     fshift_buf = make_float3(0.0f);
224
225     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
226     for (j4 = cij4_start; j4 < cij4_end; j4++)
227     {
228         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
229         imask       = pl_cj4[j4].imei[widx].imask;
230         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
231
232 #ifndef PRUNE_NBL
233         if (imask)
234 #endif
235         {
236             /* Pre-load cj into shared memory on both warps separately */
237             if ((tidxj == 0 || tidxj == 4) && tidxi < NBNXN_GPU_JGROUP_SIZE)
238             {
239                 cjs[tidxi + tidxj * NBNXN_GPU_JGROUP_SIZE / 4] = pl_cj4[j4].cj[tidxi];
240             }
241
242             /* Unrolling this loop
243                - with pruning leads to register spilling;
244                - on Kepler is much slower;
245                - doesn't work on CUDA <v4.1
246                Tested with nvcc 3.2 - 5.0.7 */
247 #if !defined PRUNE_NBL && __CUDA_ARCH__ < 300 && CUDA_VERSION >= 4010
248 #pragma unroll 4
249 #endif
250             for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
251             {
252                 if (imask & (supercl_interaction_mask << (jm * NCL_PER_SUPERCL)))
253                 {
254                     mask_ji = (1U << (jm * NCL_PER_SUPERCL));
255
256                     cj      = cjs[jm + (tidxj & 4) * NBNXN_GPU_JGROUP_SIZE / 4];
257                     aj      = cj * CL_SIZE + tidxj;
258
259                     /* load j atom data */
260                     xqbuf   = xq[aj];
261                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
262                     qj_f    = nbparam.epsfac * xqbuf.w;
263                     typej   = atom_types[aj];
264
265                     fcj_buf = make_float3(0.0f);
266
267                     /* The PME and RF kernels don't unroll with CUDA <v4.1. */
268 #if !defined PRUNE_NBL && !(CUDA_VERSION < 4010 && (defined EL_EWALD_ANY || defined EL_RF))
269 #pragma unroll 8
270 #endif
271                     for(i = 0; i < NCL_PER_SUPERCL; i++)
272                     {
273                         if (imask & mask_ji)
274                         {
275                             ci_offset   = i;    /* i force buffer offset */
276
277                             ci      = sci * NCL_PER_SUPERCL + i; /* i cluster index */
278                             ai      = ci * CL_SIZE + tidxi;      /* i atom index */
279
280                             /* all threads load an atom from i cluster ci into shmem! */
281                             xqbuf   = xqib[i * CL_SIZE + tidxi];
282                             xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
283
284                             /* distance between i and j atoms */
285                             rv      = xi - xj;
286                             r2      = norm2(rv);
287
288 #ifdef PRUNE_NBL
289                             /* If _none_ of the atoms pairs are in cutoff range,
290                                the bit corresponding to the current
291                                cluster-pair in imask gets set to 0. */
292                             if (!__any(r2 < rlist_sq))
293                             {
294                                 imask &= ~mask_ji;
295                             }
296 #endif
297
298                             int_bit = (wexcl & mask_ji) ? 1.0f : 0.0f;
299
300                             /* cutoff & exclusion check */
301 #if defined EL_EWALD_ANY || defined EL_RF
302                             if (r2 < rcoulomb_sq *
303                                 (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
304 #else
305                             if (r2 < rcoulomb_sq * int_bit)
306 #endif
307                             {
308                                 /* load the rest of the i-atom parameters */
309                                 qi      = xqbuf.w;
310 #ifdef IATYPE_SHMEM
311                                 typei   = atib[i * CL_SIZE + tidxi];
312 #else
313                                 typei   = atom_types[ai];
314 #endif
315
316                                 /* LJ 6*C6 and 12*C12 */
317 #ifdef USE_TEXOBJ
318                                 c6      = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej));
319                                 c12     = tex1Dfetch<float>(nbparam.nbfp_texobj, 2 * (ntypes * typei + typej) + 1);
320 #else
321                                 c6      = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej));
322                                 c12     = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej) + 1);
323 #endif /* USE_TEXOBJ */
324
325
326                                 /* avoid NaN for excluded pairs at r=0 */
327                                 r2      += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
328
329                                 inv_r   = rsqrt(r2);
330                                 inv_r2  = inv_r * inv_r;
331                                 inv_r6  = inv_r2 * inv_r2 * inv_r2;
332 #if defined EL_EWALD_ANY || defined EL_RF
333                                 /* We could mask inv_r2, but with Ewald
334                                  * masking both inv_r6 and F_invr is faster */
335                                 inv_r6  *= int_bit;
336 #endif
337
338                                 F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
339
340 #ifdef CALC_ENERGIES
341                                 E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 - lj_shift * lj_shift) * 0.08333333f - c6 * (inv_r6 - lj_shift) * 0.16666667f);
342 #endif
343
344 #ifdef VDW_CUTOFF_CHECK
345                                 /* this enables twin-range cut-offs (rvdw < rcoulomb <= rlist) */
346                                 vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
347                                 F_invr  *= vdw_in_range;
348 #ifdef CALC_ENERGIES
349                                 E_lj_p  *= vdw_in_range;
350 #endif
351 #endif
352 #ifdef CALC_ENERGIES
353                                 E_lj    += E_lj_p;
354 #endif
355
356
357 #ifdef EL_CUTOFF
358                                 F_invr  += qi * qj_f * inv_r2 * inv_r;
359 #endif
360 #ifdef EL_RF
361                                 F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
362 #endif
363 #if defined EL_EWALD_ANA
364                                 F_invr  += qi * qj_f * (int_bit*inv_r2*inv_r + pmecorrF(beta2*r2)*beta3);
365 #elif defined EL_EWALD_TAB
366                                 F_invr  += qi * qj_f * (int_bit*inv_r2 -
367 #ifdef USE_TEXOBJ
368                                                         interpolate_coulomb_force_r(nbparam.coulomb_tab_texobj, r2 * inv_r, coulomb_tab_scale)
369 #else
370                                                         interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)
371 #endif /* USE_TEXOBJ */
372                                                         ) * inv_r;
373 #endif /* EL_EWALD_ANA/TAB */
374
375 #ifdef CALC_ENERGIES
376 #ifdef EL_CUTOFF
377                                 E_el    += qi * qj_f * (inv_r - c_rf);
378 #endif
379 #ifdef EL_RF
380                                 E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
381 #endif
382 #ifdef EL_EWALD_ANY
383                                 /* 1.0f - erff is faster than erfcf */
384                                 E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
385 #endif /* EL_EWALD_ANY */
386 #endif
387                                 f_ij    = rv * F_invr;
388
389                                 /* accumulate j forces in registers */
390                                 fcj_buf -= f_ij;
391
392                                 /* accumulate i forces in registers */
393                                 fci_buf[ci_offset] += f_ij;
394                             }
395                         }
396
397                         /* shift the mask bit by 1 */
398                         mask_ji += mask_ji;
399                     }
400
401                     /* reduce j forces */
402 #ifdef REDUCE_SHUFFLE
403                     reduce_force_j_warp_shfl(fcj_buf, f, tidxi, aj);
404 #else
405                     /* store j forces in shmem */
406                     f_buf[                  tidx] = fcj_buf.x;
407                     f_buf[    FBUF_STRIDE + tidx] = fcj_buf.y;
408                     f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
409
410                     reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
411 #endif
412                 }
413             }
414 #ifdef PRUNE_NBL
415             /* Update the imask with the new one which does not contain the
416                out of range clusters anymore. */
417             pl_cj4[j4].imei[widx].imask = imask;
418 #endif
419         }
420     }
421
422     /* reduce i forces */
423     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
424     {
425         ai  = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
426 #ifdef REDUCE_SHUFFLE
427         reduce_force_i_warp_shfl(fci_buf[ci_offset], f,
428                                  &fshift_buf, bCalcFshift,
429                                  tidxj, ai);
430 #else
431         f_buf[                  tidx] = fci_buf[ci_offset].x;
432         f_buf[    FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
433         f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
434         __syncthreads();
435         reduce_force_i(f_buf, f,
436                        &fshift_buf, bCalcFshift,
437                        tidxi, tidxj, ai);
438         __syncthreads();
439 #endif
440     }
441
442     /* add up local shift forces into global mem */
443 #ifdef REDUCE_SHUFFLE
444     if (bCalcFshift && (tidxj == 0 || tidxj == 4))
445 #else
446     if (bCalcFshift && tidxj == 0)
447 #endif
448     {
449         atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
450         atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
451         atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
452     }
453
454 #ifdef CALC_ENERGIES
455 #ifdef REDUCE_SHUFFLE
456     /* reduce the energies over warps and store into global memory */
457     reduce_energy_warp_shfl(E_lj, E_el, e_lj, e_el, tidx);
458 #else
459     /* flush the energies to shmem and reduce them */
460     f_buf[              tidx] = E_lj;
461     f_buf[FBUF_STRIDE + tidx] = E_el;
462     reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
463 #endif
464 #endif
465 }
466
467 #undef EL_EWALD_ANY