657030789acfe7621eae89d48afd9690ce7c9d3e
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_kernel_legacy.cuh
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2012, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14  *
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35
36 #include "maths.h"
37 /* Note that floating-point constants in CUDA code should be suffixed
38  * with f (e.g. 0.5f), to stop the compiler producing intermediate
39  * code that is in double precision.
40  */
41
42 /* Analytical Ewald is not implemented for the legacy kernels (as it is anyway
43    slower than the tabulated kernel on Fermi). */
44 #ifdef EL_EWALD_ANA
45 #error Trying to generate Analytical Ewald legacy kernels which is not implemented in the legacy kernels!
46 #endif
47
48 /*
49    Kernel launch parameters:
50     - #blocks   = #pair lists, blockId = pair list Id
51     - #threads  = CL_SIZE^2
52     - shmem     = CL_SIZE^2 * sizeof(float)
53
54     Each thread calculates an i force-component taking one pair of i-j atoms.
55  */
56 #if __CUDA_ARCH__ >= 350
57 __launch_bounds__(64,16)
58 #endif
59 #ifdef PRUNE_NBL
60 #ifdef CALC_ENERGIES
61 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_prune_legacy)
62 #else
63 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _prune_legacy)
64 #endif
65 #else
66 #ifdef CALC_ENERGIES
67 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_legacy)
68 #else
69 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _legacy)
70 #endif
71 #endif
72             (const cu_atomdata_t atdat,
73              const cu_nbparam_t nbparam,
74              const cu_plist_t plist,
75              bool bCalcFshift)
76 {
77     /* convenience variables */
78     const nbnxn_sci_t *pl_sci   = plist.sci;
79 #ifndef PRUNE_NBL
80     const
81 #endif
82     nbnxn_cj4_t *pl_cj4         = plist.cj4;
83     const nbnxn_excl_t *excl    = plist.excl;
84     const int *atom_types       = atdat.atom_types;
85     int ntypes                  = atdat.ntypes;
86     const float4 *xq            = atdat.xq;
87     float3 *f                   = atdat.f;
88     const float3 *shift_vec     = atdat.shift_vec;
89     float rcoulomb_sq           = nbparam.rcoulomb_sq;
90 #ifdef VDW_CUTOFF_CHECK
91     float rvdw_sq               = nbparam.rvdw_sq;
92     float vdw_in_range;
93 #endif
94 #ifdef EL_RF
95     float two_k_rf              = nbparam.two_k_rf;
96 #endif
97 #ifdef EL_EWALD_TAB
98     float coulomb_tab_scale     = nbparam.coulomb_tab_scale;
99 #endif
100 #ifdef PRUNE_NBL
101     float rlist_sq              = nbparam.rlist_sq;
102 #endif
103
104 #ifdef CALC_ENERGIES
105     float lj_shift    = nbparam.sh_invrc6;
106 #ifdef EL_EWALD_TAB
107     float beta        = nbparam.ewald_beta;
108     float ewald_shift = nbparam.sh_ewald;
109 #else
110     float c_rf        = nbparam.c_rf;
111 #endif
112     float *e_lj       = atdat.e_lj;
113     float *e_el       = atdat.e_el;
114 #endif
115
116     /* thread/block/warp id-s */
117     unsigned int tidxi  = threadIdx.x;
118     unsigned int tidxj  = threadIdx.y;
119     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
120     unsigned int bidx   = blockIdx.x;
121     unsigned int widx   = tidx / WARP_SIZE; /* warp index */
122
123     int sci, ci, cj, ci_offset,
124         ai, aj,
125         cij4_start, cij4_end,
126         typei, typej,
127         i, cii, jm, j4, nsubi, wexcl_idx;
128     float qi, qj_f,
129           r2, inv_r, inv_r2, inv_r6,
130           c6, c12,
131           int_bit,
132 #ifdef CALC_ENERGIES
133           E_lj, E_el, E_lj_p,
134 #endif
135           F_invr;
136     unsigned int wexcl, imask, mask_ji;
137     float4 xqbuf;
138     float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
139     float3 fci_buf[NCL_PER_SUPERCL];    /* i force buffer */
140     nbnxn_sci_t nb_sci;
141
142     /* shmem buffer for i x+q pre-loading */
143     extern __shared__  float4 xqib[];
144     /* shmem j force buffer */
145     float *f_buf = (float *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
146
147     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
148     sci         = nb_sci.sci;           /* super-cluster */
149     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
150     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
151
152     /* Store the i-atom x and q in shared memory */
153     /* Note: the thread indexing here is inverted with respect to the
154        inner-loop as this results in slightly higher performance */
155     ci = sci * NCL_PER_SUPERCL + tidxi;
156     ai = ci * CL_SIZE + tidxj;
157     xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
158     __syncthreads();
159
160     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
161     {
162         fci_buf[ci_offset] = make_float3(0.0f);
163     }
164
165 #ifdef CALC_ENERGIES
166     E_lj = 0.0f;
167     E_el = 0.0f;
168
169 #if defined EL_EWALD_TAB || defined EL_RF
170     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
171     {
172         /* we have the diagonal: add the charge self interaction energy term */
173         for (i = 0; i < NCL_PER_SUPERCL; i++)
174         {
175             qi    = xqib[i * CL_SIZE + tidxi].w;
176             E_el += qi*qi;
177         }
178         /* divide the self term equally over the j-threads */
179         E_el /= CL_SIZE;
180 #ifdef EL_RF
181         E_el *= -nbparam.epsfac*0.5f*c_rf;
182 #else
183         E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
184 #endif
185     }
186 #endif
187 #endif
188
189     /* skip central shifts when summing shift forces */
190     if (nb_sci.shift == CENTRAL)
191     {
192         bCalcFshift = false;
193     }
194
195     fshift_buf = make_float3(0.0f);
196
197     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
198     for (j4 = cij4_start; j4 < cij4_end; j4++)
199     {
200         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
201         imask       = pl_cj4[j4].imei[widx].imask;
202         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
203
204 #ifndef PRUNE_NBL
205         if (imask)
206 #endif
207         {
208             /* nvcc >v4.1 doesn't like this loop, it refuses to unroll it */
209 #if CUDA_VERSION >= 4010
210             #pragma unroll 4
211 #endif
212             for (jm = 0; jm < NBNXN_GPU_JGROUP_SIZE; jm++)
213             {
214                 mask_ji = (imask >> (jm * CL_SIZE)) & supercl_interaction_mask;
215                 if (mask_ji)
216                 {
217                     nsubi = __popc(mask_ji);
218
219                     cj      = pl_cj4[j4].cj[jm];
220                     aj      = cj * CL_SIZE + tidxj;
221
222                     /* load j atom data */
223                     xqbuf   = xq[aj];
224                     xj      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
225                     qj_f    = nbparam.epsfac * xqbuf.w;
226                     typej   = atom_types[aj];
227
228                     fcj_buf = make_float3(0.0f);
229
230                     /* loop over the i-clusters in sci */
231                     /* #pragma unroll 8
232                        -- nvcc doesn't like my code, it refuses to unroll it
233                        which is a pity because here unrolling could help.  */
234                     for (cii = 0; cii < nsubi; cii++)
235                     {
236                         i = __ffs(mask_ji) - 1;
237                         mask_ji &= ~(1U << i);
238
239                         ci_offset   = i;    /* i force buffer offset */
240
241                         ci      = sci * NCL_PER_SUPERCL + i; /* i cluster index */
242                         ai      = ci * CL_SIZE + tidxi;      /* i atom index */
243
244                         /* all threads load an atom from i cluster ci into shmem! */
245                         xqbuf   = xqib[i * CL_SIZE + tidxi];
246                         xi      = make_float3(xqbuf.x, xqbuf.y, xqbuf.z);
247
248                         /* distance between i and j atoms */
249                         rv      = xi - xj;
250                         r2      = norm2(rv);
251
252 #ifdef PRUNE_NBL
253                         /* If _none_ of the atoms pairs are in cutoff range,
254                                the bit corresponding to the current
255                                cluster-pair in imask gets set to 0. */
256                         if (!__any(r2 < rlist_sq))
257                         {
258                             imask &= ~(1U << (jm * NCL_PER_SUPERCL + i));
259                         }
260 #endif
261
262                         int_bit = ((wexcl >> (jm * NCL_PER_SUPERCL + i)) & 1) ? 1.0f : 0.0f;
263
264                         /* cutoff & exclusion check */
265 #if defined EL_EWALD_TAB || defined EL_RF
266                         if (r2 < rcoulomb_sq *
267                             (nb_sci.shift != CENTRAL || ci != cj || tidxj > tidxi))
268 #else
269                         if (r2 < rcoulomb_sq * int_bit)
270 #endif
271                         {
272                             /* load the rest of the i-atom parameters */
273                             qi      = xqbuf.w;
274                             typei   = atom_types[ai];
275
276                             /* LJ 6*C6 and 12*C12 */
277                             c6      = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej));
278                             c12     = tex1Dfetch(nbfp_texref, 2 * (ntypes * typei + typej) + 1);
279
280                             /* avoid NaN for excluded pairs at r=0 */
281                             r2      += (1.0f - int_bit) * NBNXN_AVOID_SING_R2_INC;
282
283                             inv_r   = rsqrt(r2);
284                             inv_r2  = inv_r * inv_r;
285                             inv_r6  = inv_r2 * inv_r2 * inv_r2;
286 #if defined EL_EWALD_TAB || defined EL_RF
287                             /* We could mask inv_r2, but with Ewald
288                              * masking both inv_r6 and F_invr is faster */
289                             inv_r6  *= int_bit;
290 #endif
291
292                             F_invr  = inv_r6 * (c12 * inv_r6 - c6) * inv_r2;
293
294 #ifdef CALC_ENERGIES
295                             E_lj_p  = int_bit * (c12 * (inv_r6 * inv_r6 - lj_shift * lj_shift) * 0.08333333f - c6 * (inv_r6 - lj_shift) * 0.16666667f);
296 #endif
297
298 #ifdef VDW_CUTOFF_CHECK
299                             /* this enables twin-range cut-offs (rvdw < rcoulomb <= rlist) */
300                             vdw_in_range = (r2 < rvdw_sq) ? 1.0f : 0.0f;
301                             F_invr  *= vdw_in_range;
302 #ifdef CALC_ENERGIES
303                             E_lj_p  *= vdw_in_range;
304 #endif
305 #endif
306 #ifdef CALC_ENERGIES
307                             E_lj    += E_lj_p;
308 #endif
309
310
311 #ifdef EL_CUTOFF
312                             F_invr  += qi * qj_f * inv_r2 * inv_r;
313 #endif
314 #ifdef EL_RF
315                             F_invr  += qi * qj_f * (int_bit*inv_r2 * inv_r - two_k_rf);
316 #endif
317 #ifdef EL_EWALD_TAB
318                             F_invr  += qi * qj_f * (int_bit*inv_r2 - interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)) * inv_r;
319 #endif /* EL_EWALD_TAB */
320
321 #ifdef CALC_ENERGIES
322 #ifdef EL_CUTOFF
323                             E_el    += qi * qj_f * (inv_r - c_rf);
324 #endif
325 #ifdef EL_RF
326                             E_el    += qi * qj_f * (int_bit*inv_r + 0.5f * two_k_rf * r2 - c_rf);
327 #endif
328 #ifdef EL_EWALD_TAB
329                             /* 1.0f - erff is faster than erfcf */
330                             E_el    += qi * qj_f * (inv_r * (int_bit - erff(r2 * inv_r * beta)) - int_bit * ewald_shift);
331 #endif
332 #endif
333                             f_ij    = rv * F_invr;
334
335                             /* accumulate j forces in registers */
336                             fcj_buf -= f_ij;
337
338                             /* accumulate i forces in registers */
339                             fci_buf[ci_offset] += f_ij;
340                         }
341                     }
342
343                     /* store j forces in shmem */
344                     f_buf[                  tidx] = fcj_buf.x;
345                     f_buf[    FBUF_STRIDE + tidx] = fcj_buf.y;
346                     f_buf[2 * FBUF_STRIDE + tidx] = fcj_buf.z;
347
348                     /* reduce j forces */
349                     reduce_force_j_generic(f_buf, f, tidxi, tidxj, aj);
350                 }
351             }
352 #ifdef PRUNE_NBL
353             /* Update the imask with the new one which does not contain the
354                out of range clusters anymore. */
355             pl_cj4[j4].imei[widx].imask = imask;
356 #endif
357         }
358     }
359
360     /* reduce i forces */
361     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
362     {
363         ai  = (sci * NCL_PER_SUPERCL + ci_offset) * CL_SIZE + tidxi;
364         f_buf[                  tidx] = fci_buf[ci_offset].x;
365         f_buf[    FBUF_STRIDE + tidx] = fci_buf[ci_offset].y;
366         f_buf[2 * FBUF_STRIDE + tidx] = fci_buf[ci_offset].z;
367         __syncthreads();
368         reduce_force_i(f_buf, f,
369                        &fshift_buf, bCalcFshift,
370                        tidxi, tidxj, ai);
371         __syncthreads();
372     }
373
374     /* add up local shift forces into global mem */
375     if (bCalcFshift && tidxj == 0)
376     {
377         atomicAdd(&atdat.fshift[nb_sci.shift].x, fshift_buf.x);
378         atomicAdd(&atdat.fshift[nb_sci.shift].y, fshift_buf.y);
379         atomicAdd(&atdat.fshift[nb_sci.shift].z, fshift_buf.z);
380     }
381
382 #ifdef CALC_ENERGIES
383     /* flush the energies to shmem and reduce them */
384     f_buf[              tidx] = E_lj;
385     f_buf[FBUF_STRIDE + tidx] = E_el;
386     reduce_energy_pow2(f_buf + (tidx & WARP_SIZE), e_lj, e_el, tidx & ~WARP_SIZE);
387 #endif
388 }