Replace nbnxn_buffer_flags_t with vector
[alexxy/gromacs.git] / src / gromacs / nbnxm / kernel_common.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
5  * Copyright (c) 2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  *
39  * \brief
40  * Implements utility functions used by all nbnxm CPU kernels.
41  *
42  * \author Berk Hess <hess@kth.se>
43  * \ingroup module_nbnxm
44  */
45
46 #include "gmxpre.h"
47
48 #include "kernel_common.h"
49
50 #include "gromacs/pbcutil/ishift.h"
51 #include "gromacs/utility/gmxassert.h"
52
53 //! Clears all elements of buffer
54 static void clearBufferAll(gmx::ArrayRef<real> buffer)
55 {
56     for (real& elem : buffer)
57     {
58         elem = 0;
59     }
60 }
61
62 /*! \brief Clears elements of size and stride \p numComponentsPerElement
63  *
64  * Only elements with flags in \p nbat set for index \p outputIndex
65  * are cleared.
66  */
67 template<int numComponentsPerElement>
68 static void clearBufferFlagged(const nbnxn_atomdata_t& nbat, int outputIndex, gmx::ArrayRef<real> buffer)
69 {
70     gmx::ArrayRef<const gmx_bitmask_t> flags = nbat.buffer_flags;
71     gmx_bitmask_t                      our_flag;
72     bitmask_init_bit(&our_flag, outputIndex);
73
74     constexpr size_t numComponentsPerBlock = NBNXN_BUFFERFLAG_SIZE * numComponentsPerElement;
75
76     for (size_t b = 0; b < flags.size(); b++)
77     {
78         if (!bitmask_is_disjoint(flags[b], our_flag))
79         {
80             clearBufferAll(buffer.subArray(b * numComponentsPerBlock, numComponentsPerBlock));
81         }
82     }
83 }
84
85 void clearForceBuffer(nbnxn_atomdata_t* nbat, int outputIndex)
86 {
87     if (nbat->bUseBufferFlags)
88     {
89         GMX_ASSERT(nbat->fstride == DIM, "Only fstride=3 is currently handled here");
90
91         clearBufferFlagged<DIM>(*nbat, outputIndex, nbat->out[outputIndex].f);
92     }
93     else
94     {
95         clearBufferAll(nbat->out[outputIndex].f);
96     }
97 }
98
99 void clear_fshift(real* fshift)
100 {
101     int i;
102
103     for (i = 0; i < SHIFTS * DIM; i++)
104     {
105         fshift[i] = 0;
106     }
107 }
108
109 void reduce_energies_over_lists(const nbnxn_atomdata_t* nbat, int nlist, real* Vvdw, real* Vc)
110 {
111     const int nenergrp = nbat->params().nenergrp;
112
113     for (int nb = 0; nb < nlist; nb++)
114     {
115         for (int i = 0; i < nenergrp; i++)
116         {
117             /* Reduce the diagonal terms */
118             int ind = i * nenergrp + i;
119             Vvdw[ind] += nbat->out[nb].Vvdw[ind];
120             Vc[ind] += nbat->out[nb].Vc[ind];
121
122             /* Reduce the off-diagonal terms */
123             for (int j = i + 1; j < nenergrp; j++)
124             {
125                 /* The output should contain only one off-diagonal part */
126                 int ind  = i * nenergrp + j;
127                 int indr = j * nenergrp + i;
128                 Vvdw[ind] += nbat->out[nb].Vvdw[ind] + nbat->out[nb].Vvdw[indr];
129                 Vc[ind] += nbat->out[nb].Vc[ind] + nbat->out[nb].Vc[indr];
130             }
131         }
132     }
133 }