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.
43 #include "types/simple.h"
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/smalloc.h"
48 #include "nb_kernel_allvsallgb.h"
55 int ** exclusion_mask;
60 calc_maxoffset(int i, int natoms)
64 if ((natoms % 2) == 1)
66 /* Odd number of atoms, easy */
69 else if ((natoms % 4) == 0)
71 /* Multiple of four is hard */
80 maxoffset = natoms/2-1;
91 maxoffset = natoms/2-1;
100 maxoffset = natoms/2;
104 maxoffset = natoms/2-1;
113 setup_exclusions_and_indices(gmx_allvsall_data_t * aadata,
124 /* This routine can appear to be a bit complex, but it is mostly book-keeping.
125 * To enable the fast all-vs-all kernel we need to be able to stream through all coordinates
126 * whether they should interact or not.
128 * To avoid looping over the exclusions, we create a simple mask that is 1 if the interaction
129 * should be present, otherwise 0. Since exclusions typically only occur when i & j are close,
130 * we create a jindex array with three elements per i atom: the starting point, the point to
131 * which we need to check exclusions, and the end point.
132 * This way we only have to allocate a short exclusion mask per i atom.
135 /* Allocate memory for our modified jindex array */
136 snew(aadata->jindex, 3*natoms);
138 /* Pointer to lists with exclusion masks */
139 snew(aadata->exclusion_mask, natoms);
141 for (i = 0; i < natoms; i++)
144 aadata->jindex[3*i] = i+1;
145 max_offset = calc_maxoffset(i, natoms);
148 nj0 = excl->index[i];
149 nj1 = excl->index[i+1];
151 /* first check the max range */
152 max_excl_offset = -1;
154 for (j = nj0; j < nj1; j++)
160 if (k+natoms <= max_offset)
165 max_excl_offset = (k > max_excl_offset) ? k : max_excl_offset;
168 max_excl_offset = (max_offset < max_excl_offset) ? max_offset : max_excl_offset;
170 aadata->jindex[3*i+1] = i+1+max_excl_offset;
172 snew(aadata->exclusion_mask[i], max_excl_offset);
173 /* Include everything by default */
174 for (j = 0; j < max_excl_offset; j++)
176 /* Use all-ones to mark interactions that should be present, compatible with SSE */
177 aadata->exclusion_mask[i][j] = 0xFFFFFFFF;
180 /* Go through exclusions again */
181 for (j = nj0; j < nj1; j++)
187 if (k+natoms <= max_offset)
192 if (k > 0 && k <= max_excl_offset)
194 /* Excluded, kill it! */
195 aadata->exclusion_mask[i][k-1] = 0;
200 aadata->jindex[3*i+2] = i+1+max_offset;
206 setup_aadata(gmx_allvsall_data_t ** p_aadata,
214 gmx_allvsall_data_t *aadata;
220 /* Generate vdw params */
221 snew(aadata->pvdwparam, ntype);
223 for (i = 0; i < ntype; i++)
225 snew(aadata->pvdwparam[i], 2*natoms);
226 p = aadata->pvdwparam[i];
228 /* Lets keep it simple and use multiple steps - first create temp. c6/c12 arrays */
229 for (j = 0; j < natoms; j++)
231 idx = i*ntype+type[j];
232 p[2*j] = pvdwparam[2*idx];
233 p[2*j+1] = pvdwparam[2*idx+1];
237 setup_exclusions_and_indices(aadata, excl, natoms);
243 nb_kernel_allvsallgb(t_nblist gmx_unused * nlist,
248 nb_kernel_data_t * kernel_data,
251 gmx_allvsall_data_t *aadata;
266 real vgbtot, dvdasum;
274 real rsq, rinv, rinvsq, rinvsix;
276 real c6, c12, Vvdw6, Vvdw12, Vvdwtot;
277 real fscal, dvdatmp, fijC, vgb;
278 real Y, F, Fp, Geps, Heps2, VV, FF, eps, eps2, r, rt;
279 real dvdaj, gbscale, isaprod, isai, isaj, gbtabscale;
289 charge = mdatoms->chargeA;
290 type = mdatoms->typeA;
291 gbfactor = ((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
293 GBtab = fr->gbtab.data;
294 gbtabscale = fr->gbtab.scale;
295 invsqrta = fr->invsqrta;
297 vpol = kernel_data->energygrp_polarization;
299 natoms = mdatoms->nr;
301 ni1 = mdatoms->homenr;
303 aadata = fr->AllvsAll_work;
304 excl = kernel_data->exclusions;
306 Vc = kernel_data->energygrp_elec;
307 Vvdw = kernel_data->energygrp_vdw;
311 setup_aadata(&aadata, excl, natoms, type, fr->ntype, fr->nbfp);
312 fr->AllvsAll_work = aadata;
315 for (i = ni0; i < ni1; i++)
317 /* We assume shifts are NOT used for all-vs-all interactions */
319 /* Load i atom data */
323 iq = facel*charge[i];
327 pvdw = aadata->pvdwparam[type[i]];
329 /* Zero the potential energy for this list */
335 /* Clear i atom forces */
340 /* Load limits for loop over neighbors */
341 nj0 = aadata->jindex[3*i];
342 nj1 = aadata->jindex[3*i+1];
343 nj2 = aadata->jindex[3*i+2];
345 mask = aadata->exclusion_mask[i];
347 /* Prologue part, including exclusion mask */
348 for (j = nj0; j < nj1; j++, mask++)
354 /* load j atom coordinates */
359 /* Calculate distance */
363 rsq = dx*dx+dy*dy+dz*dz;
365 /* Calculate 1/r and 1/r2 */
366 rinv = gmx_invsqrt(rsq);
368 /* Load parameters for j atom */
374 qq = isaprod*(-qq)*gbfactor;
375 gbscale = isaprod*gbtabscale;
380 /* Tabulated Generalized-Born interaction */
384 /* Calculate table index */
392 Geps = eps*GBtab[nnn+2];
393 Heps2 = eps2*GBtab[nnn+3];
396 FF = Fp+Geps+2.0*Heps2;
398 fijC = qq*FF*gbscale;
399 dvdatmp = -0.5*(vgb+fijC*r);
400 dvdasum = dvdasum + dvdatmp;
401 dvda[k] = dvdaj+dvdatmp*isaj*isaj;
402 vctot = vctot + vcoul;
403 vgbtot = vgbtot + vgb;
405 /* Lennard-Jones interaction */
406 rinvsix = rinvsq*rinvsq*rinvsq;
408 Vvdw12 = c12*rinvsix*rinvsix;
409 Vvdwtot = Vvdwtot+Vvdw12-Vvdw6;
410 fscal = (12.0*Vvdw12-6.0*Vvdw6)*rinvsq-(fijC-fscal)*rinv;
412 /* Calculate temporary vectorial force */
417 /* Increment i atom force */
422 /* Decrement j atom force */
423 f[3*k] = f[3*k] - tx;
424 f[3*k+1] = f[3*k+1] - ty;
425 f[3*k+2] = f[3*k+2] - tz;
427 /* Inner loop uses 38 flops/iteration */
430 /* Main part, no exclusions */
431 for (j = nj1; j < nj2; j++)
435 /* load j atom coordinates */
440 /* Calculate distance */
444 rsq = dx*dx+dy*dy+dz*dz;
446 /* Calculate 1/r and 1/r2 */
447 rinv = gmx_invsqrt(rsq);
449 /* Load parameters for j atom */
455 qq = isaprod*(-qq)*gbfactor;
456 gbscale = isaprod*gbtabscale;
461 /* Tabulated Generalized-Born interaction */
465 /* Calculate table index */
473 Geps = eps*GBtab[nnn+2];
474 Heps2 = eps2*GBtab[nnn+3];
477 FF = Fp+Geps+2.0*Heps2;
479 fijC = qq*FF*gbscale;
480 dvdatmp = -0.5*(vgb+fijC*r);
481 dvdasum = dvdasum + dvdatmp;
482 dvda[k] = dvdaj+dvdatmp*isaj*isaj;
483 vctot = vctot + vcoul;
484 vgbtot = vgbtot + vgb;
486 /* Lennard-Jones interaction */
487 rinvsix = rinvsq*rinvsq*rinvsq;
489 Vvdw12 = c12*rinvsix*rinvsix;
490 Vvdwtot = Vvdwtot+Vvdw12-Vvdw6;
491 fscal = (12.0*Vvdw12-6.0*Vvdw6)*rinvsq-(fijC-fscal)*rinv;
493 /* Calculate temporary vectorial force */
498 /* Increment i atom force */
503 /* Decrement j atom force */
504 f[3*k] = f[3*k] - tx;
505 f[3*k+1] = f[3*k+1] - ty;
506 f[3*k+2] = f[3*k+2] - tz;
508 /* Inner loop uses 38 flops/iteration */
515 /* Add potential energies to the group for this list */
518 Vc[ggid] = Vc[ggid] + vctot;
519 Vvdw[ggid] = Vvdw[ggid] + Vvdwtot;
520 vpol[ggid] = vpol[ggid] + vgbtot;
521 dvda[i] = dvda[i] + dvdasum*isai*isai;
523 /* Outer loop uses 6 flops/iteration */
526 /* 12 flops per outer iteration
527 * 19 flops per inner iteration
529 inc_nrnb(nrnb, eNR_NBKERNEL_ELEC_VDW_VF, (ni1-ni0)*12 + ((ni1-ni0)*natoms/2)*19);