allow compilation to optimize for CUDA compute cap. 3.5
[alexxy/gromacs.git] / src / mdlib / nbnxn_cuda / nbnxn_cuda_kernel_legacy.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 /* Analytical Ewald is not implemented for the legacy kernels (as it is anyway
46    slower than the tabulated kernel on Fermi). */
47 #ifdef EL_EWALD_ANA
48 #error Trying to generate Analytical Ewald legacy kernels which is not implemented in the legacy kernels!
49 #endif
50
51 /*
52    Kernel launch parameters:
53     - #blocks   = #pair lists, blockId = pair list Id
54     - #threads  = CL_SIZE^2
55     - shmem     = CL_SIZE^2 * sizeof(float)
56
57     Each thread calculates an i force-component taking one pair of i-j atoms.
58  */
59 #if __CUDA_ARCH__ >= 350
60 __launch_bounds__(64,16)
61 #endif
62 #ifdef PRUNE_NBL
63 #ifdef CALC_ENERGIES
64 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_prune_legacy)
65 #else
66 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _prune_legacy)
67 #endif
68 #else
69 #ifdef CALC_ENERGIES
70 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_legacy)
71 #else
72 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _legacy)
73 #endif
74 #endif
75             (const cu_atomdata_t atdat,
76              const cu_nbparam_t nbparam,
77              const cu_plist_t plist,
78              bool bCalcFshift)
79 {
80     /* convenience variables */
81     const nbnxn_sci_t *pl_sci   = plist.sci;
82 #ifndef PRUNE_NBL
83     const
84 #endif
85     nbnxn_cj4_t *pl_cj4         = plist.cj4;
86     const nbnxn_excl_t *excl    = plist.excl;
87     const int *atom_types       = atdat.atom_types;
88     int ntypes                  = atdat.ntypes;
89     const float4 *xq            = atdat.xq;
90     float3 *f                   = atdat.f;
91     const float3 *shift_vec     = atdat.shift_vec;
92     float rcoulomb_sq           = nbparam.rcoulomb_sq;
93 #ifdef VDW_CUTOFF_CHECK
94     float rvdw_sq               = nbparam.rvdw_sq;
95     float vdw_in_range;
96 #endif
97 #ifdef EL_RF
98     float two_k_rf              = nbparam.two_k_rf;
99 #endif
100 #ifdef EL_EWALD_TAB
101     float coulomb_tab_scale     = nbparam.coulomb_tab_scale;
102 #endif
103 #ifdef PRUNE_NBL
104     float rlist_sq              = nbparam.rlist_sq;
105 #endif
106
107 #ifdef CALC_ENERGIES
108     float lj_shift    = nbparam.sh_invrc6;
109 #ifdef EL_EWALD_TAB
110     float beta        = nbparam.ewald_beta;
111     float ewald_shift = nbparam.sh_ewald;
112 #else
113     float c_rf        = nbparam.c_rf;
114 #endif
115     float *e_lj       = atdat.e_lj;
116     float *e_el       = atdat.e_el;
117 #endif
118
119     /* thread/block/warp id-s */
120     unsigned int tidxi  = threadIdx.x;
121     unsigned int tidxj  = threadIdx.y;
122     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
123     unsigned int bidx   = blockIdx.x;
124     unsigned int widx   = tidx / WARP_SIZE; /* warp index */
125
126     int sci, ci, cj, ci_offset,
127         ai, aj,
128         cij4_start, cij4_end,
129         typei, typej,
130         i, cii, jm, j4, nsubi, wexcl_idx;
131     float qi, qj_f,
132           r2, inv_r, inv_r2, inv_r6,
133           c6, c12,
134           int_bit,
135 #ifdef CALC_ENERGIES
136           E_lj, E_el, E_lj_p,
137 #endif
138           F_invr;
139     unsigned int wexcl, imask, mask_ji;
140     float4 xqbuf;
141     float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
142     float3 fci_buf[NCL_PER_SUPERCL];    /* i force buffer */
143     nbnxn_sci_t nb_sci;
144
145     /* shmem buffer for i x+q pre-loading */
146     extern __shared__  float4 xqib[];
147     /* shmem j force buffer */
148     float *f_buf = (float *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
149
150     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
151     sci         = nb_sci.sci;           /* super-cluster */
152     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
153     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
154
155     /* Store the i-atom x and q in shared memory */
156     /* Note: the thread indexing here is inverted with respect to the
157        inner-loop as this results in slightly higher performance */
158     ci = sci * NCL_PER_SUPERCL + tidxi;
159     ai = ci * CL_SIZE + tidxj;
160     xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
161     __syncthreads();
162
163     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
164     {
165         fci_buf[ci_offset] = make_float3(0.0f);
166     }
167
168 #ifdef CALC_ENERGIES
169     E_lj = 0.0f;
170     E_el = 0.0f;
171
172 #if defined EL_EWALD_TAB || defined EL_RF
173     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
174     {
175         /* we have the diagonal: add the charge self interaction energy term */
176         for (i = 0; i < NCL_PER_SUPERCL; i++)
177         {
178             qi    = xqib[i * CL_SIZE + tidxi].w;
179             E_el += qi*qi;
180         }
181         /* divide the self term equally over the j-threads */
182         E_el /= CL_SIZE;
183 #ifdef EL_RF
184         E_el *= -nbparam.epsfac*0.5f*c_rf;
185 #else
186         E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
187 #endif
188     }
189 #endif
190 #endif
191
192     /* skip central shifts when summing shift forces */
193     if (nb_sci.shift == CENTRAL)
194     {
195         bCalcFshift = false;
196     }
197
198     fshift_buf = make_float3(0.0f);
199
200     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
201     for (j4 = cij4_start; j4 < cij4_end; j4++)
202     {
203         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
204         imask       = pl_cj4[j4].imei[widx].imask;
205         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
206
207 #ifndef PRUNE_NBL
208         if (imask)
209 #endif
210         {
211             /* nvcc >v4.1 doesn't like this loop, it refuses to unroll it */
212 #if CUDA_VERSION >= 4010
213             #pragma unroll 4
214 #endif
215             for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
216             {
217                 mask_ji = (imask >> (jm * CL_SIZE)) & supercl_interaction_mask;
218                 if (mask_ji)
219                 {
220                     nsubi = __popc(mask_ji);
221
222                     cj      = pl_cj4[j4].cj[jm];
223                     aj      = cj * CL_SIZE + tidxj;
224
225                     /* load j atom data */
226                     xqbuf   = xq[aj];
227                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
228                     qj_f    = nbparam.epsfac * xqbuf.w;
229                     typej   = atom_types[aj];
230
231                     fcj_buf = make_float3(0.0f);
232
233                     /* loop over the i-clusters in sci */
234                     /* #pragma unroll 8
235                        -- nvcc doesn't like my code, it refuses to unroll it
236                        which is a pity because here unrolling could help.  */
237                     for (cii = 0; cii < nsubi; cii++)
238                     {
239                         i = __ffs(mask_ji) - 1;
240                         mask_ji &= ~(1U << i);
241
242                         ci_offset   = i;    /* i force buffer offset */
243
244                         ci      = sci * NCL_PER_SUPERCL + i; /* i cluster index */
245                         ai      = ci * CL_SIZE + tidxi;      /* i atom index */
246
247                         /* all threads load an atom from i cluster ci into shmem! */
248                         xqbuf   = xqib[i * CL_SIZE + tidxi];
249                         xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
250
251                         /* distance between i and j atoms */
252                         rv      = xi - xj;
253                         r2      = norm2(rv);
254
255 #ifdef PRUNE_NBL
256                         /* If _none_ of the atoms pairs are in cutoff range,
257                                the bit corresponding to the current
258                                cluster-pair in imask gets set to 0. */
259                         if (!__any(r2 < rlist_sq))
260                         {
261                             imask &= ~(1U << (jm * NCL_PER_SUPERCL + i));
262                         }
263 #endif
264
265                         int_bit = ((wexcl >> (jm * NCL_PER_SUPERCL + i)) & 1) ? 1.0f : 0.0f;
266
267                         /* cutoff & exclusion check */
268 #if defined EL_EWALD_TAB || defined EL_RF
269                         if (r2 < rcoulomb_sq *
270                             (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
271 #else
272                         if (r2 < rcoulomb_sq * int_bit)
273 #endif
274                         {
275                             /* load the rest of the i-atom parameters */
276                             qi      = xqbuf.w;
277                             typei   = atom_types[ai];
278
279                             /* LJ 6*C6 and 12*C12 */
280                             c6      = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej));
281                             c12     = tex1Dfetch(tex_nbfp, 2 * (ntypes * typei + typej) + 1);
282
283                             /* avoid NaN for excluded pairs at r=0 */
284                             r2      += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
285
286                             inv_r   = rsqrt(r2);
287                             inv_r2  = inv_r * inv_r;
288                             inv_r6  = inv_r2 * inv_r2 * inv_r2;
289 #if defined EL_EWALD_TAB || defined EL_RF
290                             /* We could mask inv_r2, but with Ewald
291                              * masking both inv_r6 and F_invr is faster */
292                             inv_r6  *= int_bit;
293 #endif
294
295                             F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
296
297 #ifdef CALC_ENERGIES
298                             E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 - lj_shift * lj_shift) * 0.08333333f - c6 * (inv_r6 - lj_shift) * 0.16666667f);
299 #endif
300
301 #ifdef VDW_CUTOFF_CHECK
302                             /* this enables twin-range cut-offs (rvdw < rcoulomb <= rlist) */
303                             vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
304                             F_invr  *= vdw_in_range;
305 #ifdef CALC_ENERGIES
306                             E_lj_p  *= vdw_in_range;
307 #endif
308 #endif
309 #ifdef CALC_ENERGIES
310                             E_lj    += E_lj_p;
311 #endif
312
313
314 #ifdef EL_CUTOFF
315                             F_invr  += qi * qj_f * inv_r2 * inv_r;
316 #endif
317 #ifdef EL_RF
318                             F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
319 #endif
320 #ifdef EL_EWALD_TAB
321                             F_invr  += qi * qj_f * (int_bit*inv_r2 - interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)) * inv_r;
322 #endif /* EL_EWALD_TAB */
323
324 #ifdef CALC_ENERGIES
325 #ifdef EL_CUTOFF
326                             E_el    += qi * qj_f * (inv_r - c_rf);
327 #endif
328 #ifdef EL_RF
329                             E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
330 #endif
331 #ifdef EL_EWALD_TAB
332                             /* 1.0f - erff is faster than erfcf */
333                             E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
334 #endif
335 #endif
336                             f_ij    = rv * F_invr;
337
338                             /* accumulate j forces in registers */
339                             fcj_buf -= f_ij;
340
341                             /* accumulate i forces in registers */
342                             fci_buf[ci_offset] += f_ij;
343                         }
344                     }
345
346                     /* store j forces in shmem */
347                     f_buf[                  tidx] = fcj_buf.x;
348                     f_buf[    FBUF_STRIDE + tidx] = fcj_buf.y;
349                     f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
350
351                     /* reduce j forces */
352                     reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
353                 }
354             }
355 #ifdef PRUNE_NBL
356             /* Update the imask with the new one which does not contain the
357                out of range clusters anymore. */
358             pl_cj4[j4].imei[widx].imask = imask;
359 #endif
360         }
361     }
362
363     /* reduce i forces */
364     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
365     {
366         ai  = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
367         f_buf[                  tidx] = fci_buf[ci_offset].x;
368         f_buf[    FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
369         f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
370         __syncthreads();
371         reduce_force_i(f_buf, f,
372                        &fshift_buf, bCalcFshift,
373                        tidxi, tidxj, ai);
374         __syncthreads();
375     }
376
377     /* add up local shift forces into global mem */
378     if (bCalcFshift && tidxj == 0)
379     {
380         atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
381         atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
382         atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
383     }
384
385 #ifdef CALC_ENERGIES
386     /* flush the energies to shmem and reduce them */
387     f_buf[              tidx] = E_lj;
388     f_buf[FBUF_STRIDE + tidx] = E_el;
389     reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
390 #endif
391 }