b4a3f1327fcabd45e3def92d72e9cb9a140e4749
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_kernels / nbnxn_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 #include "gmxpre.h"
36
37 #include "nbnxn_kernel_common.h"
38
39 #include "gromacs/pbcutil/ishift.h"
40 #include "gromacs/utility/gmxassert.h"
41
42 static void
43 clear_f_all(const nbnxn_atomdata_t *nbat, real *f)
44 {
45     for (int i = 0; i < nbat->numAtoms()*nbat->fstride; i++)
46     {
47         f[i] = 0;
48     }
49 }
50
51 static void
52 clear_f_flagged(const nbnxn_atomdata_t *nbat, int output_index, real *f)
53 {
54     GMX_ASSERT(nbat->fstride == DIM, "For performance we use compile time constant DIM instead of nbat->fstride");
55
56     const nbnxn_buffer_flags_t &flags = nbat->buffer_flags;
57     gmx_bitmask_t               our_flag;
58     bitmask_init_bit(&our_flag, output_index);
59
60     for (int b = 0; b < flags.nflag; b++)
61     {
62         if (!bitmask_is_disjoint(flags.flag[b], our_flag))
63         {
64             int a0 = b*NBNXN_BUFFERFLAG_SIZE;
65             int a1 = a0 + NBNXN_BUFFERFLAG_SIZE;
66             for (int i = a0*DIM; i < a1*DIM; i++)
67             {
68                 f[i] = 0;
69             }
70         }
71     }
72 }
73
74 void
75 clear_f(const nbnxn_atomdata_t *nbat, int output_index, real *f)
76 {
77     if (nbat->bUseBufferFlags)
78     {
79         clear_f_flagged(nbat, output_index, f);
80     }
81     else
82     {
83         clear_f_all(nbat, f);
84     }
85 }
86
87 void
88 clear_fshift(real *fshift)
89 {
90     int i;
91
92     for (i = 0; i < SHIFTS*DIM; i++)
93     {
94         fshift[i] = 0;
95     }
96 }
97
98 void
99 reduce_energies_over_lists(const nbnxn_atomdata_t     *nbat,
100                            int                         nlist,
101                            real                       *Vvdw,
102                            real                       *Vc)
103 {
104     const int nenergrp = nbat->params().nenergrp;
105
106     for (int nb = 0; nb < nlist; nb++)
107     {
108         for (int i = 0; i < nenergrp; i++)
109         {
110             /* Reduce the diagonal terms */
111             int ind    = i*nenergrp + i;
112             Vvdw[ind] += nbat->out[nb].Vvdw[ind];
113             Vc[ind]   += nbat->out[nb].Vc[ind];
114
115             /* Reduce the off-diagonal terms */
116             for (int j = i + 1; j < nenergrp; j++)
117             {
118                 /* The output should contain only one off-diagonal part */
119                 int ind    = i*nenergrp + j;
120                 int indr   = j*nenergrp + i;
121                 Vvdw[ind] += nbat->out[nb].Vvdw[ind] + nbat->out[nb].Vvdw[indr];
122                 Vc[ind]   += nbat->out[nb].Vc[ind]   + nbat->out[nb].Vc[indr];
123             }
124         }
125     }
126 }