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_allvsallgb.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,
123 /* This routine can appear to be a bit complex, but it is mostly book-keeping.
124 * To enable the fast all-vs-all kernel we need to be able to stream through all coordinates
125 * whether they should interact or not.
127 * To avoid looping over the exclusions, we create a simple mask that is 1 if the interaction
128 * should be present, otherwise 0. Since exclusions typically only occur when i & j are close,
129 * we create a jindex array with three elements per i atom: the starting point, the point to
130 * which we need to check exclusions, and the end point.
131 * This way we only have to allocate a short exclusion mask per i atom.
134 /* Allocate memory for our modified jindex array */
135 snew(aadata->jindex, 3*natoms);
137 /* Pointer to lists with exclusion masks */
138 snew(aadata->exclusion_mask, natoms);
140 for (i = 0; i < natoms; i++)
143 aadata->jindex[3*i] = i+1;
144 max_offset = calc_maxoffset(i, natoms);
147 nj0 = excl->index[i];
148 nj1 = excl->index[i+1];
150 /* first check the max range */
151 max_excl_offset = -1;
153 for (j = nj0; j < nj1; j++)
159 if (k+natoms <= max_offset)
164 max_excl_offset = (k > max_excl_offset) ? k : max_excl_offset;
167 max_excl_offset = (max_offset < max_excl_offset) ? max_offset : max_excl_offset;
169 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;
205 setup_aadata(gmx_allvsall_data_t ** p_aadata,
213 gmx_allvsall_data_t *aadata;
219 /* Generate vdw params */
220 snew(aadata->pvdwparam, ntype);
222 for (i = 0; i < ntype; i++)
224 snew(aadata->pvdwparam[i], 2*natoms);
225 p = aadata->pvdwparam[i];
227 /* Lets keep it simple and use multiple steps - first create temp. c6/c12 arrays */
228 for (j = 0; j < natoms; j++)
230 idx = i*ntype+type[j];
231 p[2*j] = pvdwparam[2*idx];
232 p[2*j+1] = pvdwparam[2*idx+1];
236 setup_exclusions_and_indices(aadata, excl, natoms);
242 nb_kernel_allvsallgb(t_nblist gmx_unused * nlist,
247 nb_kernel_data_t * kernel_data,
250 gmx_allvsall_data_t *aadata;
265 real vgbtot, dvdasum;
273 real rsq, rinv, rinvsq, rinvsix;
275 real c6, c12, Vvdw6, Vvdw12, Vvdwtot;
276 real fscal, dvdatmp, fijC, vgb;
277 real Y, F, Fp, Geps, Heps2, VV, FF, eps, eps2, r, rt;
278 real dvdaj, gbscale, isaprod, isai, isaj, gbtabscale;
288 charge = mdatoms->chargeA;
289 type = mdatoms->typeA;
290 gbfactor = ((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
292 GBtab = fr->gbtab.data;
293 gbtabscale = fr->gbtab.scale;
294 invsqrta = fr->invsqrta;
296 vpol = kernel_data->energygrp_polarization;
298 natoms = mdatoms->nr;
300 ni1 = mdatoms->homenr;
302 aadata = fr->AllvsAll_work;
303 excl = kernel_data->exclusions;
305 Vc = kernel_data->energygrp_elec;
306 Vvdw = kernel_data->energygrp_vdw;
310 setup_aadata(&aadata, excl, natoms, type, fr->ntype, fr->nbfp);
311 fr->AllvsAll_work = aadata;
314 for (i = ni0; i < ni1; i++)
316 /* We assume shifts are NOT used for all-vs-all interactions */
318 /* Load i atom data */
322 iq = facel*charge[i];
326 pvdw = aadata->pvdwparam[type[i]];
328 /* Zero the potential energy for this list */
334 /* Clear i atom forces */
339 /* Load limits for loop over neighbors */
340 nj0 = aadata->jindex[3*i];
341 nj1 = aadata->jindex[3*i+1];
342 nj2 = aadata->jindex[3*i+2];
344 mask = aadata->exclusion_mask[i];
346 /* Prologue part, including exclusion mask */
347 for (j = nj0; j < nj1; j++, mask++)
353 /* load j atom coordinates */
358 /* Calculate distance */
362 rsq = dx*dx+dy*dy+dz*dz;
364 /* Calculate 1/r and 1/r2 */
365 rinv = gmx_invsqrt(rsq);
367 /* Load parameters for j atom */
373 qq = isaprod*(-qq)*gbfactor;
374 gbscale = isaprod*gbtabscale;
379 /* Tabulated Generalized-Born interaction */
383 /* Calculate table index */
391 Geps = eps*GBtab[nnn+2];
392 Heps2 = eps2*GBtab[nnn+3];
395 FF = Fp+Geps+2.0*Heps2;
397 fijC = qq*FF*gbscale;
398 dvdatmp = -0.5*(vgb+fijC*r);
399 dvdasum = dvdasum + dvdatmp;
400 dvda[k] = dvdaj+dvdatmp*isaj*isaj;
401 vctot = vctot + vcoul;
402 vgbtot = vgbtot + vgb;
404 /* Lennard-Jones interaction */
405 rinvsix = rinvsq*rinvsq*rinvsq;
407 Vvdw12 = c12*rinvsix*rinvsix;
408 Vvdwtot = Vvdwtot+Vvdw12-Vvdw6;
409 fscal = (12.0*Vvdw12-6.0*Vvdw6)*rinvsq-(fijC-fscal)*rinv;
411 /* Calculate temporary vectorial force */
416 /* Increment i atom force */
421 /* Decrement j atom force */
422 f[3*k] = f[3*k] - tx;
423 f[3*k+1] = f[3*k+1] - ty;
424 f[3*k+2] = f[3*k+2] - tz;
426 /* Inner loop uses 38 flops/iteration */
429 /* Main part, no exclusions */
430 for (j = nj1; j < nj2; j++)
434 /* load j atom coordinates */
439 /* Calculate distance */
443 rsq = dx*dx+dy*dy+dz*dz;
445 /* Calculate 1/r and 1/r2 */
446 rinv = gmx_invsqrt(rsq);
448 /* Load parameters for j atom */
454 qq = isaprod*(-qq)*gbfactor;
455 gbscale = isaprod*gbtabscale;
460 /* Tabulated Generalized-Born interaction */
464 /* Calculate table index */
472 Geps = eps*GBtab[nnn+2];
473 Heps2 = eps2*GBtab[nnn+3];
476 FF = Fp+Geps+2.0*Heps2;
478 fijC = qq*FF*gbscale;
479 dvdatmp = -0.5*(vgb+fijC*r);
480 dvdasum = dvdasum + dvdatmp;
481 dvda[k] = dvdaj+dvdatmp*isaj*isaj;
482 vctot = vctot + vcoul;
483 vgbtot = vgbtot + vgb;
485 /* Lennard-Jones interaction */
486 rinvsix = rinvsq*rinvsq*rinvsq;
488 Vvdw12 = c12*rinvsix*rinvsix;
489 Vvdwtot = Vvdwtot+Vvdw12-Vvdw6;
490 fscal = (12.0*Vvdw12-6.0*Vvdw6)*rinvsq-(fijC-fscal)*rinv;
492 /* Calculate temporary vectorial force */
497 /* Increment i atom force */
502 /* Decrement j atom force */
503 f[3*k] = f[3*k] - tx;
504 f[3*k+1] = f[3*k+1] - ty;
505 f[3*k+2] = f[3*k+2] - tz;
507 /* Inner loop uses 38 flops/iteration */
514 /* Add potential energies to the group for this list */
517 Vc[ggid] = Vc[ggid] + vctot;
518 Vvdw[ggid] = Vvdw[ggid] + Vvdwtot;
519 vpol[ggid] = vpol[ggid] + vgbtot;
520 dvda[i] = dvda[i] + dvdasum*isai*isai;
522 /* Outer loop uses 6 flops/iteration */
525 /* 12 flops per outer iteration
526 * 19 flops per inner iteration
528 inc_nrnb(nrnb, eNR_NBKERNEL_ELEC_VDW_VF, (ni1-ni0)*12 + ((ni1-ni0)*natoms/2)*19);