Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / nbnxm / nbnxm_geometry.h
blobc43e22e7ba1ed522993b8ddcbea9435cc0efc367
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
36 #ifndef GMX_NBNXM_NBNXM_GEOMETRY_H
37 #define GMX_NBNXM_NBNXM_GEOMETRY_H
39 #include "gromacs/math/vectypes.h"
40 #include "gromacs/nbnxm/nbnxm.h"
41 #include "gromacs/nbnxm/pairlist.h"
42 #include "gromacs/simd/simd.h"
43 #include "gromacs/utility/fatalerror.h"
45 /* Returns the base-2 log of n.
46 * Generates a fatal error when n is not an integer power of 2.
48 static inline int get_2log(int n)
50 int log2;
52 log2 = 0;
53 while ((1 << log2) < n)
55 log2++;
57 if ((1 << log2) != n)
59 gmx_fatal(FARGS, "nbnxn na_c (%d) is not a power of 2", n);
62 return log2;
65 namespace Nbnxm
68 /* The nbnxn i-cluster size in atoms for each nbnxn kernel type */
69 static constexpr gmx::EnumerationArray<KernelType, int> IClusterSizePerKernelType =
70 { {
72 c_nbnxnCpuIClusterSize,
73 c_nbnxnCpuIClusterSize,
74 c_nbnxnCpuIClusterSize,
75 c_nbnxnGpuClusterSize,
76 c_nbnxnGpuClusterSize
77 } };
79 /* The nbnxn j-cluster size in atoms for each nbnxn kernel type */
80 static constexpr gmx::EnumerationArray<KernelType, int> JClusterSizePerKernelType =
81 { {
83 c_nbnxnCpuIClusterSize,
84 #if GMX_SIMD
85 GMX_SIMD_REAL_WIDTH,
86 GMX_SIMD_REAL_WIDTH/2,
87 #else
90 #endif
91 c_nbnxnGpuClusterSize,
92 c_nbnxnGpuClusterSize
93 } };
95 /* Returns whether the pair-list corresponding to nb_kernel_type is simple */
96 static inline bool kernelTypeUsesSimplePairlist(const KernelType kernelType)
98 return (kernelType == KernelType::Cpu4x4_PlainC ||
99 kernelType == KernelType::Cpu4xN_Simd_4xN ||
100 kernelType == KernelType::Cpu4xN_Simd_2xNN);
103 static inline bool kernelTypeIsSimd(const KernelType kernelType)
105 return (kernelType == KernelType::Cpu4xN_Simd_4xN ||
106 kernelType == KernelType::Cpu4xN_Simd_2xNN);
109 } // namespace Nbnxm
111 /* Returns the effective list radius of the pair-list
113 * Due to the cluster size the effective pair-list is longer than
114 * that of a simple atom pair-list. This function gives the extra distance.
116 * NOTE: If the i- and j-cluster sizes are identical and you know
117 * the physical dimensions of the clusters, use the next function
118 * for more accurate results
120 real nbnxn_get_rlist_effective_inc(int jClusterSize,
121 real atomDensity);
123 /* Returns the effective list radius of the pair-list
125 * Due to the cluster size the effective pair-list is longer than
126 * that of a simple atom pair-list. This function gives the extra distance.
128 real nbnxn_get_rlist_effective_inc(int clusterSize,
129 const gmx::RVec &averageClusterBoundingBox);
131 #endif