2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
39 * Declares the geometry-related functionality
41 * \author Berk Hess <hess@kth.se>
42 * \ingroup module_nbnxm
44 #ifndef GMX_NBNXM_NBNXM_GEOMETRY_H
45 #define GMX_NBNXM_NBNXM_GEOMETRY_H
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/nbnxm/nbnxm.h"
49 #include "gromacs/simd/simd.h"
50 #include "gromacs/utility/fatalerror.h"
55 /*! \brief Returns the base-2 log of n.
57 * Generates a fatal error when n is not an integer power of 2.
59 static inline int get_2log(int n
)
64 while ((1 << log2
) < n
)
70 gmx_fatal(FARGS
, "nbnxn na_c (%d) is not a power of 2", n
);
79 /*! \brief The nbnxn i-cluster size in atoms for each nbnxn kernel type */
80 static constexpr gmx::EnumerationArray
<KernelType
, int> IClusterSizePerKernelType
= {
81 { 0, c_nbnxnCpuIClusterSize
, c_nbnxnCpuIClusterSize
, c_nbnxnCpuIClusterSize
,
82 c_nbnxnGpuClusterSize
, c_nbnxnGpuClusterSize
}
85 /*! \brief The nbnxn j-cluster size in atoms for each nbnxn kernel type */
86 static constexpr gmx::EnumerationArray
<KernelType
, int> JClusterSizePerKernelType
= {
87 { 0, c_nbnxnCpuIClusterSize
,
89 GMX_SIMD_REAL_WIDTH
, GMX_SIMD_REAL_WIDTH
/ 2,
93 c_nbnxnGpuClusterSize
, c_nbnxnGpuClusterSize
/ 2 }
96 /*! \brief Returns whether the pair-list corresponding to nb_kernel_type is simple */
97 static inline bool kernelTypeUsesSimplePairlist(const KernelType kernelType
)
99 return (kernelType
== KernelType::Cpu4x4_PlainC
|| kernelType
== KernelType::Cpu4xN_Simd_4xN
100 || kernelType
== KernelType::Cpu4xN_Simd_2xNN
);
103 //! Returns whether a SIMD kernel is in use
104 static inline bool kernelTypeIsSimd(const KernelType kernelType
)
106 return (kernelType
== KernelType::Cpu4xN_Simd_4xN
|| kernelType
== KernelType::Cpu4xN_Simd_2xNN
);
111 /*! \brief 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
, real atomDensity
);
122 /*! \brief Returns the effective list radius of the pair-list
124 * Due to the cluster size the effective pair-list is longer than
125 * that of a simple atom pair-list. This function gives the extra distance.
127 real
nbnxn_get_rlist_effective_inc(int clusterSize
, const gmx::RVec
& averageClusterBoundingBox
);