Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / nbnxm / nbnxm_geometry.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_NBNXM_NBNXM_GEOMETRY_H
37 #define GMX_NBNXM_NBNXM_GEOMETRY_H
38
39 #include "gromacs/math/vectypes.h"
40 #include "gromacs/nbnxm/nbnxm.h"
41 #include "gromacs/simd/simd.h"
42 #include "gromacs/utility/fatalerror.h"
43
44 #include "pairlist.h"
45
46
47 /* Returns the base-2 log of n.
48  * Generates a fatal error when n is not an integer power of 2.
49  */
50 static inline int get_2log(int n)
51 {
52     int log2;
53
54     log2 = 0;
55     while ((1 << log2) < n)
56     {
57         log2++;
58     }
59     if ((1 << log2) != n)
60     {
61         gmx_fatal(FARGS, "nbnxn na_c (%d) is not a power of 2", n);
62     }
63
64     return log2;
65 }
66
67 namespace Nbnxm
68 {
69
70 /* The nbnxn i-cluster size in atoms for each nbnxn kernel type */
71 static constexpr gmx::EnumerationArray<KernelType, int> IClusterSizePerKernelType = {
72     { 0, c_nbnxnCpuIClusterSize, c_nbnxnCpuIClusterSize, c_nbnxnCpuIClusterSize,
73       c_nbnxnGpuClusterSize, c_nbnxnGpuClusterSize }
74 };
75
76 /* The nbnxn j-cluster size in atoms for each nbnxn kernel type */
77 static constexpr gmx::EnumerationArray<KernelType, int> JClusterSizePerKernelType = {
78     { 0, c_nbnxnCpuIClusterSize,
79 #if GMX_SIMD
80       GMX_SIMD_REAL_WIDTH, GMX_SIMD_REAL_WIDTH / 2,
81 #else
82       0, 0,
83 #endif
84       c_nbnxnGpuClusterSize, c_nbnxnGpuClusterSize }
85 };
86
87 /* Returns whether the pair-list corresponding to nb_kernel_type is simple */
88 static inline bool kernelTypeUsesSimplePairlist(const KernelType kernelType)
89 {
90     return (kernelType == KernelType::Cpu4x4_PlainC || kernelType == KernelType::Cpu4xN_Simd_4xN
91             || kernelType == KernelType::Cpu4xN_Simd_2xNN);
92 }
93
94 static inline bool kernelTypeIsSimd(const KernelType kernelType)
95 {
96     return (kernelType == KernelType::Cpu4xN_Simd_4xN || kernelType == KernelType::Cpu4xN_Simd_2xNN);
97 }
98
99 } // namespace Nbnxm
100
101 /* Returns the effective list radius of the pair-list
102  *
103  * Due to the cluster size the effective pair-list is longer than
104  * that of a simple atom pair-list. This function gives the extra distance.
105  *
106  * NOTE: If the i- and j-cluster sizes are identical and you know
107  *       the physical dimensions of the clusters, use the next function
108  *       for more accurate results
109  */
110 real nbnxn_get_rlist_effective_inc(int jClusterSize, real atomDensity);
111
112 /* Returns the effective list radius of the pair-list
113  *
114  * Due to the cluster size the effective pair-list is longer than
115  * that of a simple atom pair-list. This function gives the extra distance.
116  */
117 real nbnxn_get_rlist_effective_inc(int clusterSize, const gmx::RVec& averageClusterBoundingBox);
118
119 #endif