50fc7d163daef0b7a4a998279e0aae299f7e9687
[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,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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief
39  * Implements utility functions used by all nbnxm CPU kernels.
40  *
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_nbnxm
43  */
44
45 #include "gmxpre.h"
46
47 #include "kernel_common.h"
48
49 #include "gromacs/pbcutil/ishift.h"
50 #include "gromacs/utility/gmxassert.h"
51
52 //! Clears all elements of buffer
53 static void clearBufferAll(gmx::ArrayRef<real> buffer)
54 {
55     for (real& elem : buffer)
56     {
57         elem = 0;
58     }
59 }
60
61 /*! \brief Clears elements of size and stride \p numComponentsPerElement
62  *
63  * Only elements with flags in \p nbat set for index \p outputIndex
64  * are cleared.
65  */
66 template<int numComponentsPerElement>
67 static void clearBufferFlagged(const nbnxn_atomdata_t& nbat, int outputIndex, gmx::ArrayRef<real> buffer)
68 {
69     const nbnxn_buffer_flags_t& flags = nbat.buffer_flags;
70     gmx_bitmask_t               our_flag;
71     bitmask_init_bit(&our_flag, outputIndex);
72
73     constexpr size_t numComponentsPerBlock = NBNXN_BUFFERFLAG_SIZE * numComponentsPerElement;
74
75     for (int b = 0; b < flags.nflag; b++)
76     {
77         if (!bitmask_is_disjoint(flags.flag[b], our_flag))
78         {
79             clearBufferAll(buffer.subArray(b * numComponentsPerBlock, numComponentsPerBlock));
80         }
81     }
82 }
83
84 void clearForceBuffer(nbnxn_atomdata_t* nbat, int outputIndex)
85 {
86     if (nbat->bUseBufferFlags)
87     {
88         GMX_ASSERT(nbat->fstride == DIM, "Only fstride=3 is currently handled here");
89
90         clearBufferFlagged<DIM>(*nbat, outputIndex, nbat->out[outputIndex].f);
91     }
92     else
93     {
94         clearBufferAll(nbat->out[outputIndex].f);
95     }
96 }
97
98 void clear_fshift(real* fshift)
99 {
100     int i;
101
102     for (i = 0; i < SHIFTS * DIM; i++)
103     {
104         fshift[i] = 0;
105     }
106 }
107
108 void reduce_energies_over_lists(const nbnxn_atomdata_t* nbat, int nlist, real* Vvdw, real* Vc)
109 {
110     const int nenergrp = nbat->params().nenergrp;
111
112     for (int nb = 0; nb < nlist; nb++)
113     {
114         for (int i = 0; i < nenergrp; i++)
115         {
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];
120
121             /* Reduce the off-diagonal terms */
122             for (int j = i + 1; j < nenergrp; j++)
123             {
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];
129             }
130         }
131     }
132 }