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