1ab116afdedea512b464366c89920b17df4d889a
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_util.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,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 #ifndef GMX_MDLIB_NBNXN_UTIL_H
37 #define GMX_MDLIB_NBNXN_UTIL_H
38
39 #include <cmath>
40
41 #include "gromacs/mdlib/nb_verlet.h"
42 #include "gromacs/mdlib/nbnxn_consts.h"
43 #include "gromacs/mdlib/nbnxn_pairlist.h"
44 #include "gromacs/simd/simd.h"
45 #include "gromacs/utility/fatalerror.h"
46
47
48 /* Returns the base-2 log of n.
49  * Generates a fatal error when n is not an integer power of 2.
50  */
51 static inline int get_2log(int n)
52 {
53     int log2;
54
55     log2 = 0;
56     while ((1 << log2) < n)
57     {
58         log2++;
59     }
60     if ((1 << log2) != n)
61     {
62         gmx_fatal(FARGS, "nbnxn na_c (%d) is not a power of 2", n);
63     }
64
65     return log2;
66 }
67
68 /* Returns the nbnxn i-cluster size in atoms for the nbnxn kernel type */
69 static inline int nbnxn_kernel_to_cluster_i_size(int nb_kernel_type)
70 {
71     switch (nb_kernel_type)
72     {
73         case nbnxnk4x4_PlainC:
74         case nbnxnk4xN_SIMD_4xN:
75         case nbnxnk4xN_SIMD_2xNN:
76             return c_nbnxnCpuIClusterSize;
77         case nbnxnk8x8x8_GPU:
78         case nbnxnk8x8x8_PlainC:
79             /* The cluster size for super/sub lists is only set here.
80              * Any value should work for the pair-search and atomdata code.
81              * The kernels, of course, might require a particular value.
82              */
83             return c_nbnxnGpuClusterSize;
84         default:
85             gmx_incons("unknown kernel type");
86     }
87
88     return 0;
89 }
90
91 /* Returns the nbnxn i-cluster size in atoms for the nbnxn kernel type */
92 static inline int nbnxn_kernel_to_cluster_j_size(int nb_kernel_type)
93 {
94     int nbnxn_simd_width = 0;
95     int cj_size          = 0;
96
97 #if GMX_SIMD
98     nbnxn_simd_width = GMX_SIMD_REAL_WIDTH;
99 #endif
100
101     switch (nb_kernel_type)
102     {
103         case nbnxnk4x4_PlainC:
104             cj_size = c_nbnxnCpuIClusterSize;
105             break;
106         case nbnxnk4xN_SIMD_4xN:
107             cj_size = nbnxn_simd_width;
108             break;
109         case nbnxnk4xN_SIMD_2xNN:
110             cj_size = nbnxn_simd_width/2;
111             break;
112         case nbnxnk8x8x8_GPU:
113         case nbnxnk8x8x8_PlainC:
114             cj_size = nbnxn_kernel_to_cluster_i_size(nb_kernel_type);
115             break;
116         default:
117             gmx_incons("unknown kernel type");
118     }
119
120     return cj_size;
121 }
122
123
124 #endif