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