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