Split lines with many copyright years
[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 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #ifndef GMX_NBNXM_NBNXM_GEOMETRY_H
38 #define GMX_NBNXM_NBNXM_GEOMETRY_H
39
40 #include "gromacs/math/vectypes.h"
41 #include "gromacs/nbnxm/nbnxm.h"
42 #include "gromacs/simd/simd.h"
43 #include "gromacs/utility/fatalerror.h"
44
45 #include "pairlist.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 namespace Nbnxm
69 {
70
71 /* The nbnxn i-cluster size in atoms for each nbnxn kernel type */
72 static constexpr gmx::EnumerationArray<KernelType, int> IClusterSizePerKernelType = {
73     { 0, c_nbnxnCpuIClusterSize, c_nbnxnCpuIClusterSize, c_nbnxnCpuIClusterSize,
74       c_nbnxnGpuClusterSize, c_nbnxnGpuClusterSize }
75 };
76
77 /* The nbnxn j-cluster size in atoms for each nbnxn kernel type */
78 static constexpr gmx::EnumerationArray<KernelType, int> JClusterSizePerKernelType = {
79     { 0, c_nbnxnCpuIClusterSize,
80 #if GMX_SIMD
81       GMX_SIMD_REAL_WIDTH, GMX_SIMD_REAL_WIDTH / 2,
82 #else
83       0, 0,
84 #endif
85       c_nbnxnGpuClusterSize, c_nbnxnGpuClusterSize }
86 };
87
88 /* Returns whether the pair-list corresponding to nb_kernel_type is simple */
89 static inline bool kernelTypeUsesSimplePairlist(const KernelType kernelType)
90 {
91     return (kernelType == KernelType::Cpu4x4_PlainC || kernelType == KernelType::Cpu4xN_Simd_4xN
92             || kernelType == KernelType::Cpu4xN_Simd_2xNN);
93 }
94
95 static inline bool kernelTypeIsSimd(const KernelType kernelType)
96 {
97     return (kernelType == KernelType::Cpu4xN_Simd_4xN || kernelType == KernelType::Cpu4xN_Simd_2xNN);
98 }
99
100 } // namespace Nbnxm
101
102 /* Returns the effective list radius of the pair-list
103  *
104  * Due to the cluster size the effective pair-list is longer than
105  * that of a simple atom pair-list. This function gives the extra distance.
106  *
107  * NOTE: If the i- and j-cluster sizes are identical and you know
108  *       the physical dimensions of the clusters, use the next function
109  *       for more accurate results
110  */
111 real nbnxn_get_rlist_effective_inc(int jClusterSize, real atomDensity);
112
113 /* Returns the effective list radius of the pair-list
114  *
115  * Due to the cluster size the effective pair-list is longer than
116  * that of a simple atom pair-list. This function gives the extra distance.
117  */
118 real nbnxn_get_rlist_effective_inc(int clusterSize, const gmx::RVec& averageClusterBoundingBox);
119
120 #endif