2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017,2019, 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.
39 * Implements utility functions used by all nbnxm CPU kernels.
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_nbnxm
47 #include "kernel_common.h"
49 #include "gromacs/pbcutil/ishift.h"
50 #include "gromacs/utility/gmxassert.h"
52 //! Clears all elements of buffer
53 static void clearBufferAll(gmx::ArrayRef<real> buffer)
55 for (real& elem : buffer)
61 /*! \brief Clears elements of size and stride \p numComponentsPerElement
63 * Only elements with flags in \p nbat set for index \p outputIndex
66 template<int numComponentsPerElement>
67 static void clearBufferFlagged(const nbnxn_atomdata_t& nbat, int outputIndex, gmx::ArrayRef<real> buffer)
69 const nbnxn_buffer_flags_t& flags = nbat.buffer_flags;
70 gmx_bitmask_t our_flag;
71 bitmask_init_bit(&our_flag, outputIndex);
73 constexpr size_t numComponentsPerBlock = NBNXN_BUFFERFLAG_SIZE * numComponentsPerElement;
75 for (int b = 0; b < flags.nflag; b++)
77 if (!bitmask_is_disjoint(flags.flag[b], our_flag))
79 clearBufferAll(buffer.subArray(b * numComponentsPerBlock, numComponentsPerBlock));
84 void clearForceBuffer(nbnxn_atomdata_t* nbat, int outputIndex)
86 if (nbat->bUseBufferFlags)
88 GMX_ASSERT(nbat->fstride == DIM, "Only fstride=3 is currently handled here");
90 clearBufferFlagged<DIM>(*nbat, outputIndex, nbat->out[outputIndex].f);
94 clearBufferAll(nbat->out[outputIndex].f);
98 void clear_fshift(real* fshift)
102 for (i = 0; i < SHIFTS * DIM; i++)
108 void reduce_energies_over_lists(const nbnxn_atomdata_t* nbat, int nlist, real* Vvdw, real* Vc)
110 const int nenergrp = nbat->params().nenergrp;
112 for (int nb = 0; nb < nlist; nb++)
114 for (int i = 0; i < nenergrp; i++)
116 /* Reduce the diagonal terms */
117 int ind = i * nenergrp + i;
118 Vvdw[ind] += nbat->out[nb].Vvdw[ind];
119 Vc[ind] += nbat->out[nb].Vc[ind];
121 /* Reduce the off-diagonal terms */
122 for (int j = i + 1; j < nenergrp; j++)
124 /* The output should contain only one off-diagonal part */
125 int ind = i * nenergrp + j;
126 int indr = j * nenergrp + i;
127 Vvdw[ind] += nbat->out[nb].Vvdw[ind] + nbat->out[nb].Vvdw[indr];
128 Vc[ind] += nbat->out[nb].Vc[ind] + nbat->out[nb].Vc[indr];