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.
39 #include "nb_kernel_allvsall.h"
45 #include "gromacs/legacyheaders/nrnb.h"
46 #include "gromacs/legacyheaders/types/simple.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/smalloc.h"
54 int ** exclusion_mask;
59 calc_maxoffset(int i, int natoms)
63 if ((natoms % 2) == 1)
65 /* Odd number of atoms, easy */
68 else if ((natoms % 4) == 0)
70 /* Multiple of four is hard */
79 maxoffset = natoms/2-1;
90 maxoffset = natoms/2-1;
103 maxoffset = natoms/2-1;
112 setup_exclusions_and_indices(gmx_allvsall_data_t * aadata,
122 /* This routine can appear to be a bit complex, but it is mostly book-keeping.
123 * To enable the fast all-vs-all kernel we need to be able to stream through all coordinates
124 * whether they should interact or not.
126 * To avoid looping over the exclusions, we create a simple mask that is 1 if the interaction
127 * should be present, otherwise 0. Since exclusions typically only occur when i & j are close,
128 * we create a jindex array with three elements per i atom: the starting point, the point to
129 * which we need to check exclusions, and the end point.
130 * This way we only have to allocate a short exclusion mask per i atom.
133 /* Allocate memory for our modified jindex array */
134 snew(aadata->jindex, 3*natoms);
136 /* Pointer to lists with exclusion masks */
137 snew(aadata->exclusion_mask, natoms);
139 for (i = 0; i < natoms; i++)
142 aadata->jindex[3*i] = i+1;
143 max_offset = calc_maxoffset(i, natoms);
146 nj0 = excl->index[i];
147 nj1 = excl->index[i+1];
149 /* first check the max range */
150 max_excl_offset = -1;
152 for (j = nj0; j < nj1; j++)
158 if (k+natoms <= max_offset)
163 max_excl_offset = (k > max_excl_offset) ? k : max_excl_offset;
166 max_excl_offset = (max_offset < max_excl_offset) ? max_offset : max_excl_offset;
168 aadata->jindex[3*i+1] = i+1+max_excl_offset;
171 snew(aadata->exclusion_mask[i], max_excl_offset);
172 /* Include everything by default */
173 for (j = 0; j < max_excl_offset; j++)
175 /* Use all-ones to mark interactions that should be present, compatible with SSE */
176 aadata->exclusion_mask[i][j] = 0xFFFFFFFF;
179 /* Go through exclusions again */
180 for (j = nj0; j < nj1; j++)
186 if (k+natoms <= max_offset)
191 if (k > 0 && k <= max_excl_offset)
193 /* Excluded, kill it! */
194 aadata->exclusion_mask[i][k-1] = 0;
199 aadata->jindex[3*i+2] = i+1+max_offset;
204 setup_aadata(gmx_allvsall_data_t ** p_aadata,
212 gmx_allvsall_data_t *aadata;
218 /* Generate vdw params */
219 snew(aadata->pvdwparam, ntype);
221 for (i = 0; i < ntype; i++)
223 snew(aadata->pvdwparam[i], 2*natoms);
224 p = aadata->pvdwparam[i];
226 /* Lets keep it simple and use multiple steps - first create temp. c6/c12 arrays */
227 for (j = 0; j < natoms; j++)
229 idx = i*ntype+type[j];
230 p[2*j] = pvdwparam[2*idx];
231 p[2*j+1] = pvdwparam[2*idx+1];
235 setup_exclusions_and_indices(aadata, excl, natoms);
241 nb_kernel_allvsall(t_nblist gmx_unused * nlist,
246 nb_kernel_data_t * kernel_data,
249 gmx_allvsall_data_t *aadata;
266 real rsq, rinv, rinvsq, rinvsix;
268 real c6, c12, Vvdw6, Vvdw12, Vvdwtot;
278 charge = mdatoms->chargeA;
279 type = mdatoms->typeA;
281 natoms = mdatoms->nr;
283 ni1 = mdatoms->homenr;
284 aadata = fr->AllvsAll_work;
285 excl = kernel_data->exclusions;
287 Vc = kernel_data->energygrp_elec;
288 Vvdw = kernel_data->energygrp_vdw;
292 setup_aadata(&aadata, excl, natoms, type, fr->ntype, fr->nbfp);
293 fr->AllvsAll_work = aadata;
296 for (i = ni0; i < ni1; i++)
298 /* We assume shifts are NOT used for all-vs-all interactions */
300 /* Load i atom data */
304 iq = facel*charge[i];
306 pvdw = aadata->pvdwparam[type[i]];
308 /* Zero the potential energy for this list */
312 /* Clear i atom forces */
317 /* Load limits for loop over neighbors */
318 nj0 = aadata->jindex[3*i];
319 nj1 = aadata->jindex[3*i+1];
320 nj2 = aadata->jindex[3*i+2];
322 mask = aadata->exclusion_mask[i];
324 /* Prologue part, including exclusion mask */
325 for (j = nj0; j < nj1; j++, mask++)
331 /* load j atom coordinates */
336 /* Calculate distance */
340 rsq = dx*dx+dy*dy+dz*dz;
342 /* Calculate 1/r and 1/r2 */
343 rinv = gmx_invsqrt(rsq);
346 /* Load parameters for j atom */
351 /* Coulomb interaction */
355 /* Lennard-Jones interaction */
356 rinvsix = rinvsq*rinvsq*rinvsq;
358 Vvdw12 = c12*rinvsix*rinvsix;
359 Vvdwtot = Vvdwtot+Vvdw12-Vvdw6;
360 fscal = (vcoul+12.0*Vvdw12-6.0*Vvdw6)*rinvsq;
362 /* Calculate temporary vectorial force */
367 /* Increment i atom force */
372 /* Decrement j atom force */
373 f[3*k] = f[3*k] - tx;
374 f[3*k+1] = f[3*k+1] - ty;
375 f[3*k+2] = f[3*k+2] - tz;
377 /* Inner loop uses 38 flops/iteration */
380 /* Main part, no exclusions */
381 for (j = nj1; j < nj2; j++)
385 /* load j atom coordinates */
390 /* Calculate distance */
394 rsq = dx*dx+dy*dy+dz*dz;
396 /* Calculate 1/r and 1/r2 */
397 rinv = gmx_invsqrt(rsq);
400 /* Load parameters for j atom */
405 /* Coulomb interaction */
409 /* Lennard-Jones interaction */
410 rinvsix = rinvsq*rinvsq*rinvsq;
412 Vvdw12 = c12*rinvsix*rinvsix;
413 Vvdwtot = Vvdwtot+Vvdw12-Vvdw6;
414 fscal = (vcoul+12.0*Vvdw12-6.0*Vvdw6)*rinvsq;
416 /* Calculate temporary vectorial force */
421 /* Increment i atom force */
426 /* Decrement j atom force */
427 f[3*k] = f[3*k] - tx;
428 f[3*k+1] = f[3*k+1] - ty;
429 f[3*k+2] = f[3*k+2] - tz;
431 /* Inner loop uses 38 flops/iteration */
438 /* Add potential energies to the group for this list */
441 Vc[ggid] = Vc[ggid] + vctot;
442 Vvdw[ggid] = Vvdw[ggid] + Vvdwtot;
444 /* Outer loop uses 6 flops/iteration */
447 /* 12 flops per outer iteration
448 * 19 flops per inner iteration
450 inc_nrnb(nrnb, eNR_NBKERNEL_ELEC_VDW_VF, (ni1-ni0)*12 + ((ni1-ni0)*natoms/2)*19);