2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2009, The GROMACS Development Team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
41 #include "types/simple.h"
43 #include "gromacs/math/vec.h"
44 #include "gromacs/utility/smalloc.h"
46 #include "nb_kernel_allvsall.h"
53 int ** exclusion_mask;
58 calc_maxoffset(int i, int natoms)
62 if ((natoms % 2) == 1)
64 /* Odd number of atoms, easy */
67 else if ((natoms % 4) == 0)
69 /* Multiple of four is hard */
78 maxoffset = natoms/2-1;
89 maxoffset = natoms/2-1;
102 maxoffset = natoms/2-1;
111 setup_exclusions_and_indices(gmx_allvsall_data_t * aadata,
121 /* This routine can appear to be a bit complex, but it is mostly book-keeping.
122 * To enable the fast all-vs-all kernel we need to be able to stream through all coordinates
123 * whether they should interact or not.
125 * To avoid looping over the exclusions, we create a simple mask that is 1 if the interaction
126 * should be present, otherwise 0. Since exclusions typically only occur when i & j are close,
127 * we create a jindex array with three elements per i atom: the starting point, the point to
128 * which we need to check exclusions, and the end point.
129 * This way we only have to allocate a short exclusion mask per i atom.
132 /* Allocate memory for our modified jindex array */
133 snew(aadata->jindex, 3*natoms);
135 /* Pointer to lists with exclusion masks */
136 snew(aadata->exclusion_mask, natoms);
138 for (i = 0; i < natoms; i++)
141 aadata->jindex[3*i] = i+1;
142 max_offset = calc_maxoffset(i, natoms);
145 nj0 = excl->index[i];
146 nj1 = excl->index[i+1];
148 /* first check the max range */
149 max_excl_offset = -1;
151 for (j = nj0; j < nj1; j++)
157 if (k+natoms <= max_offset)
162 max_excl_offset = (k > max_excl_offset) ? k : max_excl_offset;
165 max_excl_offset = (max_offset < max_excl_offset) ? max_offset : max_excl_offset;
167 aadata->jindex[3*i+1] = i+1+max_excl_offset;
170 snew(aadata->exclusion_mask[i], max_excl_offset);
171 /* Include everything by default */
172 for (j = 0; j < max_excl_offset; j++)
174 /* Use all-ones to mark interactions that should be present, compatible with SSE */
175 aadata->exclusion_mask[i][j] = 0xFFFFFFFF;
178 /* Go through exclusions again */
179 for (j = nj0; j < nj1; j++)
185 if (k+natoms <= max_offset)
190 if (k > 0 && k <= max_excl_offset)
192 /* Excluded, kill it! */
193 aadata->exclusion_mask[i][k-1] = 0;
198 aadata->jindex[3*i+2] = i+1+max_offset;
203 setup_aadata(gmx_allvsall_data_t ** p_aadata,
211 gmx_allvsall_data_t *aadata;
217 /* Generate vdw params */
218 snew(aadata->pvdwparam, ntype);
220 for (i = 0; i < ntype; i++)
222 snew(aadata->pvdwparam[i], 2*natoms);
223 p = aadata->pvdwparam[i];
225 /* Lets keep it simple and use multiple steps - first create temp. c6/c12 arrays */
226 for (j = 0; j < natoms; j++)
228 idx = i*ntype+type[j];
229 p[2*j] = pvdwparam[2*idx];
230 p[2*j+1] = pvdwparam[2*idx+1];
234 setup_exclusions_and_indices(aadata, excl, natoms);
240 nb_kernel_allvsall(t_nblist gmx_unused * nlist,
245 nb_kernel_data_t * kernel_data,
248 gmx_allvsall_data_t *aadata;
265 real rsq, rinv, rinvsq, rinvsix;
267 real c6, c12, Vvdw6, Vvdw12, Vvdwtot;
277 charge = mdatoms->chargeA;
278 type = mdatoms->typeA;
280 natoms = mdatoms->nr;
282 ni1 = mdatoms->homenr;
283 aadata = fr->AllvsAll_work;
284 excl = kernel_data->exclusions;
286 Vc = kernel_data->energygrp_elec;
287 Vvdw = kernel_data->energygrp_vdw;
291 setup_aadata(&aadata, excl, natoms, type, fr->ntype, fr->nbfp);
292 fr->AllvsAll_work = aadata;
295 for (i = ni0; i < ni1; i++)
297 /* We assume shifts are NOT used for all-vs-all interactions */
299 /* Load i atom data */
303 iq = facel*charge[i];
305 pvdw = aadata->pvdwparam[type[i]];
307 /* Zero the potential energy for this list */
311 /* Clear i atom forces */
316 /* Load limits for loop over neighbors */
317 nj0 = aadata->jindex[3*i];
318 nj1 = aadata->jindex[3*i+1];
319 nj2 = aadata->jindex[3*i+2];
321 mask = aadata->exclusion_mask[i];
323 /* Prologue part, including exclusion mask */
324 for (j = nj0; j < nj1; j++, mask++)
330 /* load j atom coordinates */
335 /* Calculate distance */
339 rsq = dx*dx+dy*dy+dz*dz;
341 /* Calculate 1/r and 1/r2 */
342 rinv = gmx_invsqrt(rsq);
345 /* Load parameters for j atom */
350 /* Coulomb interaction */
354 /* Lennard-Jones interaction */
355 rinvsix = rinvsq*rinvsq*rinvsq;
357 Vvdw12 = c12*rinvsix*rinvsix;
358 Vvdwtot = Vvdwtot+Vvdw12-Vvdw6;
359 fscal = (vcoul+12.0*Vvdw12-6.0*Vvdw6)*rinvsq;
361 /* Calculate temporary vectorial force */
366 /* Increment i atom force */
371 /* Decrement j atom force */
372 f[3*k] = f[3*k] - tx;
373 f[3*k+1] = f[3*k+1] - ty;
374 f[3*k+2] = f[3*k+2] - tz;
376 /* Inner loop uses 38 flops/iteration */
379 /* Main part, no exclusions */
380 for (j = nj1; j < nj2; j++)
384 /* load j atom coordinates */
389 /* Calculate distance */
393 rsq = dx*dx+dy*dy+dz*dz;
395 /* Calculate 1/r and 1/r2 */
396 rinv = gmx_invsqrt(rsq);
399 /* Load parameters for j atom */
404 /* Coulomb interaction */
408 /* Lennard-Jones interaction */
409 rinvsix = rinvsq*rinvsq*rinvsq;
411 Vvdw12 = c12*rinvsix*rinvsix;
412 Vvdwtot = Vvdwtot+Vvdw12-Vvdw6;
413 fscal = (vcoul+12.0*Vvdw12-6.0*Vvdw6)*rinvsq;
415 /* Calculate temporary vectorial force */
420 /* Increment i atom force */
425 /* Decrement j atom force */
426 f[3*k] = f[3*k] - tx;
427 f[3*k+1] = f[3*k+1] - ty;
428 f[3*k+2] = f[3*k+2] - tz;
430 /* Inner loop uses 38 flops/iteration */
437 /* Add potential energies to the group for this list */
440 Vc[ggid] = Vc[ggid] + vctot;
441 Vvdw[ggid] = Vvdw[ggid] + Vvdwtot;
443 /* Outer loop uses 6 flops/iteration */
446 /* 12 flops per outer iteration
447 * 19 flops per inner iteration
449 inc_nrnb(nrnb, eNR_NBKERNEL_ELEC_VDW_VF, (ni1-ni0)*12 + ((ni1-ni0)*natoms/2)*19);