0733d9412e2577869b9b1afd1304387bc09c334f
[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 /*
46    Kernel launch parameters:
47     - #blocks   = #pair lists, blockId = pair list Id
48     - #threads  = CL_SIZE^2
49     - shmem     = CL_SIZE^2 * sizeof(float)
50
51     Each thread calculates an i force-component taking one pair of i-j atoms.
52  */
53 #ifdef PRUNE_NBL
54 #ifdef CALC_ENERGIES
55 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_prune_legacy)
56 #else
57 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _prune_legacy)
58 #endif
59 #else
60 #ifdef CALC_ENERGIES
61 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _ener_legacy)
62 #else
63 __global__ void NB_KERNEL_FUNC_NAME(k_nbnxn, _legacy)
64 #endif
65 #endif
66             (const cu_atomdata_t atdat,
67              const cu_nbparam_t nbparam,
68              const cu_plist_t plist,
69              bool bCalcFshift)
70 {
71     /* convenience variables */
72     const nbnxn_sci_t *pl_sci   = plist.sci;
73 #ifndef PRUNE_NBL
74     const
75 #endif
76     nbnxn_cj4_t *pl_cj4         = plist.cj4;
77     const nbnxn_excl_t *excl    = plist.excl;
78     const int *atom_types       = atdat.atom_types;
79     int ntypes                  = atdat.ntypes;
80     const float4 *xq            = atdat.xq;
81     float3 *f                   = atdat.f;
82     const float3 *shift_vec     = atdat.shift_vec;
83     float rcoulomb_sq           = nbparam.rcoulomb_sq;
84 #ifdef VDW_CUTOFF_CHECK
85     float rvdw_sq               = nbparam.rvdw_sq;
86     float vdw_in_range;
87 #endif
88 #ifdef EL_RF
89     float two_k_rf              = nbparam.two_k_rf;
90 #endif
91 #ifdef EL_EWALD
92     float coulomb_tab_scale     = nbparam.coulomb_tab_scale;
93 #endif
94 #ifdef PRUNE_NBL
95     float rlist_sq              = nbparam.rlist_sq;
96 #endif
97
98 #ifdef CALC_ENERGIES
99     float lj_shift    = nbparam.sh_invrc6;
100 #ifdef EL_EWALD
101     float beta        = nbparam.ewald_beta;
102     float ewald_shift = nbparam.sh_ewald;
103 #else
104     float c_rf        = nbparam.c_rf;
105 #endif
106     float *e_lj       = atdat.e_lj;
107     float *e_el       = atdat.e_el;
108 #endif
109
110     /* thread/block/warp id-s */
111     unsigned int tidxi  = threadIdx.x;
112     unsigned int tidxj  = threadIdx.y;
113     unsigned int tidx   = threadIdx.y * blockDim.x + threadIdx.x;
114     unsigned int bidx   = blockIdx.x;
115     unsigned int widx   = tidx / WARP_SIZE; /* warp index */
116
117     int sci, ci, cj, ci_offset,
118         ai, aj,
119         cij4_start, cij4_end,
120         typei, typej,
121         i, cii, jm, j4, nsubi, wexcl_idx;
122     float qi, qj_f,
123           r2, inv_r, inv_r2, inv_r6,
124           c6, c12,
125 #ifdef CALC_ENERGIES
126           E_lj, E_el, E_lj_p,
127 #endif
128           F_invr;
129     unsigned int wexcl, int_bit, imask, imask_j;
130 #ifdef PRUNE_NBL
131     unsigned int imask_prune;
132 #endif
133     float4 xqbuf;
134     float3 xi, xj, rv, f_ij, fcj_buf, fshift_buf;
135     float3 fci_buf[NCL_PER_SUPERCL];    /* i force buffer */
136     nbnxn_sci_t nb_sci;
137
138     /* shmem buffer for i x+q pre-loading */
139     extern __shared__  float4 xqib[];
140     /* shmem j force buffer */
141     float *f_buf = (float *)(xqib + NCL_PER_SUPERCL * CL_SIZE);
142
143     nb_sci      = pl_sci[bidx];         /* my i super-cluster's index = current bidx */
144     sci         = nb_sci.sci;           /* super-cluster */
145     cij4_start  = nb_sci.cj4_ind_start; /* first ...*/
146     cij4_end    = nb_sci.cj4_ind_end;   /* and last index of j clusters */
147
148     /* Store the i-atom x and q in shared memory */
149     /* Note: the thread indexing here is inverted with respect to the
150        inner-loop as this results in slightly higher performance */
151     ci = sci * NCL_PER_SUPERCL + tidxi;
152     ai = ci * CL_SIZE + tidxj;
153     xqib[tidxi * CL_SIZE + tidxj] = xq[ai] + shift_vec[nb_sci.shift];
154     __syncthreads();
155
156     for(ci_offset = 0; ci_offset < NCL_PER_SUPERCL; ci_offset++)
157     {
158         fci_buf[ci_offset] = make_float3(0.0f);
159     }
160
161 #ifdef CALC_ENERGIES
162     E_lj = 0.0f;
163     E_el = 0.0f;
164
165 #if defined EL_EWALD || defined EL_RF
166     if (nb_sci.shift == CENTRAL && pl_cj4[cij4_start].cj[0] == sci*NCL_PER_SUPERCL)
167     {
168         /* we have the diagonal: add the charge self interaction energy term */
169         for (i = 0; i < NCL_PER_SUPERCL; i++)
170         {
171             qi    = xqib[i * CL_SIZE + tidxi].w;
172             E_el += qi*qi;
173         }
174         /* divide the self term equally over the j-threads */
175         E_el /= CL_SIZE;
176 #ifdef EL_RF
177         E_el *= -nbparam.epsfac*0.5f*c_rf;
178 #else
179         E_el *= -nbparam.epsfac*beta*M_FLOAT_1_SQRTPI; /* last factor 1/sqrt(pi) */
180 #endif
181     }
182 #endif
183 #endif
184
185     /* skip central shifts when summing shift forces */
186     if (nb_sci.shift == CENTRAL)
187     {
188         bCalcFshift = false;
189     }
190
191     fshift_buf = make_float3(0.0f);
192
193     /* loop over the j clusters = seen by any of the atoms in the current super-cluster */
194     for (j4 = cij4_start; j4 < cij4_end; j4++)
195     {
196         wexcl_idx   = pl_cj4[j4].imei[widx].excl_ind;
197         imask       = pl_cj4[j4].imei[widx].imask;
198         wexcl       = excl[wexcl_idx].pair[(tidx) & (WARP_SIZE - 1)];
199
200 #ifndef PRUNE_NBL
201         if (imask)
202 #endif
203         {
204 #ifdef PRUNE_NBL
205             imask_prune = 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                 imask_j = (imask >> (jm * CL_SIZE)) & supercl_interaction_mask;
215                 if (imask_j)
216                 {
217                     nsubi = __popc(imask_j);
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(imask_j) - 1;
237                         imask_j &= ~(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_prune &= ~(1U << (jm * NCL_PER_SUPERCL + i));
259                         }
260 #endif
261
262                         int_bit = ((wexcl >> (jm * NCL_PER_SUPERCL + i)) & 1);
263
264                         /* cutoff & exclusion check */
265 #if defined EL_EWALD || 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(tex_nbfp, 2 * (ntypes * typei + typej));
278                             c12     = tex1Dfetch(tex_nbfp, 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 || 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
318                             F_invr  += qi * qj_f * (int_bit*inv_r2 - interpolate_coulomb_force_r(r2 * inv_r, coulomb_tab_scale)) * inv_r;
319 #endif
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
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_prune;
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 }