Continue removing -nb gpu_cpu
[gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_kernel_common.cpp
blobbd199329983d6931ed8666a39962db4e68e62b02
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 #include "gmxpre.h"
37 #include "nbnxn_kernel_common.h"
39 #include "gromacs/pbcutil/ishift.h"
41 static void
42 clear_f_all(const nbnxn_atomdata_t *nbat, real *f)
44 int i;
46 for (i = 0; i < nbat->natoms*nbat->fstride; i++)
48 f[i] = 0;
52 static void
53 clear_f_flagged(const nbnxn_atomdata_t *nbat, int output_index, real *f)
55 const nbnxn_buffer_flags_t *flags;
56 gmx_bitmask_t our_flag;
57 int b, a0, a1, i;
59 flags = &nbat->buffer_flags;
61 bitmask_init_bit(&our_flag, output_index);
63 for (b = 0; b < flags->nflag; b++)
65 if (!bitmask_is_disjoint(flags->flag[b], our_flag))
67 a0 = b*NBNXN_BUFFERFLAG_SIZE;
68 a1 = a0 + NBNXN_BUFFERFLAG_SIZE;
69 for (i = a0*nbat->fstride; i < a1*nbat->fstride; i++)
71 f[i] = 0;
77 void
78 clear_f(const nbnxn_atomdata_t *nbat, int output_index, real *f)
80 if (nbat->bUseBufferFlags)
82 clear_f_flagged(nbat, output_index, f);
84 else
86 clear_f_all(nbat, f);
90 void
91 clear_fshift(real *fshift)
93 int i;
95 for (i = 0; i < SHIFTS*DIM; i++)
97 fshift[i] = 0;
101 void
102 reduce_energies_over_lists(const nbnxn_atomdata_t *nbat,
103 int nlist,
104 real *Vvdw,
105 real *Vc)
107 int nb;
108 int i, j, ind, indr;
110 for (nb = 0; nb < nlist; nb++)
112 for (i = 0; i < nbat->nenergrp; i++)
114 /* Reduce the diagonal terms */
115 ind = i*nbat->nenergrp + i;
116 Vvdw[ind] += nbat->out[nb].Vvdw[ind];
117 Vc[ind] += nbat->out[nb].Vc[ind];
119 /* Reduce the off-diagonal terms */
120 for (j = i+1; j < nbat->nenergrp; j++)
122 /* The output should contain only one off-diagonal part */
123 ind = i*nbat->nenergrp + j;
124 indr = j*nbat->nenergrp + i;
125 Vvdw[ind] += nbat->out[nb].Vvdw[ind] + nbat->out[nb].Vvdw[indr];
126 Vc[ind] += nbat->out[nb].Vc[ind] + nbat->out[nb].Vc[indr];