Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / nbnxm / kernel_common.cpp
blob55a470b08c849884b380f89f80c5af2fe8f3ddbc
1 /*
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.
36 /*! \internal \file
38 * \brief
39 * Implements utility functions used by all nbnxm CPU kernels.
41 * \author Berk Hess <hess@kth.se>
44 #include "gmxpre.h"
46 #include "kernel_common.h"
48 #include "gromacs/pbcutil/ishift.h"
49 #include "gromacs/utility/gmxassert.h"
51 //! Clears all elements of buffer
52 static void
53 clearBufferAll(gmx::ArrayRef<real> buffer)
55 for (real &elem : buffer)
57 elem = 0;
61 /*! \brief Clears elements of size and stride \p numComponentsPerElement
63 * Only elements with flags in \p nbat set for index \p outputIndex
64 * are cleared.
66 template <int numComponentsPerElement>
67 static void
68 clearBufferFlagged(const nbnxn_atomdata_t &nbat,
69 int outputIndex,
70 gmx::ArrayRef<real> buffer)
72 const nbnxn_buffer_flags_t &flags = nbat.buffer_flags;
73 gmx_bitmask_t our_flag;
74 bitmask_init_bit(&our_flag, outputIndex);
76 constexpr size_t numComponentsPerBlock = NBNXN_BUFFERFLAG_SIZE*numComponentsPerElement;
78 for (int b = 0; b < flags.nflag; b++)
80 if (!bitmask_is_disjoint(flags.flag[b], our_flag))
82 clearBufferAll(buffer.subArray(b*numComponentsPerBlock, numComponentsPerBlock));
87 void
88 clearForceBuffer(nbnxn_atomdata_t *nbat,
89 int outputIndex)
91 if (nbat->bUseBufferFlags)
93 GMX_ASSERT(nbat->fstride == DIM, "Only fstride=3 is currently handled here");
95 clearBufferFlagged<DIM>(*nbat, outputIndex, nbat->out[outputIndex].f);
97 else
99 clearBufferAll(nbat->out[outputIndex].f);
103 void
104 clear_fshift(real *fshift)
106 int i;
108 for (i = 0; i < SHIFTS*DIM; i++)
110 fshift[i] = 0;
114 void
115 reduce_energies_over_lists(const nbnxn_atomdata_t *nbat,
116 int nlist,
117 real *Vvdw,
118 real *Vc)
120 const int nenergrp = nbat->params().nenergrp;
122 for (int nb = 0; nb < nlist; nb++)
124 for (int i = 0; i < nenergrp; i++)
126 /* Reduce the diagonal terms */
127 int ind = i*nenergrp + i;
128 Vvdw[ind] += nbat->out[nb].Vvdw[ind];
129 Vc[ind] += nbat->out[nb].Vc[ind];
131 /* Reduce the off-diagonal terms */
132 for (int j = i + 1; j < nenergrp; j++)
134 /* The output should contain only one off-diagonal part */
135 int ind = i*nenergrp + j;
136 int indr = j*nenergrp + i;
137 Vvdw[ind] += nbat->out[nb].Vvdw[ind] + nbat->out[nb].Vvdw[indr];
138 Vc[ind] += nbat->out[nb].Vc[ind] + nbat->out[nb].Vc[indr];