4c55452e9caf49c95e83aaab8048d8b68407c6cd
[alexxy/gromacs.git] / src / gromacs / nbnxm / nbnxm_geometry.cpp
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 #include "gmxpre.h"
37
38 #include "nbnxm_geometry.h"
39
40 #include "gromacs/nbnxm/nbnxm.h"
41 #include "gromacs/nbnxm/pairlist.h"
42 #include "gromacs/simd/simd.h"
43 #include "gromacs/utility/fatalerror.h"
44 #include "gromacs/utility/real.h"
45
46 bool nbnxn_kernel_pairlist_simple(int nb_kernel_type)
47 {
48     if (nb_kernel_type == nbnxnkNotSet)
49     {
50         gmx_fatal(FARGS, "Non-bonded kernel type not set for Verlet-style pair-list.");
51     }
52
53     switch (nb_kernel_type)
54     {
55         case nbnxnk8x8x8_GPU:
56         case nbnxnk8x8x8_PlainC:
57             return false;
58
59         case nbnxnk4x4_PlainC:
60         case nbnxnk4xN_SIMD_4xN:
61         case nbnxnk4xN_SIMD_2xNN:
62             return true;
63
64         default:
65             gmx_incons("Invalid nonbonded kernel type passed!");
66             return false;
67     }
68 }
69
70 int nbnxn_kernel_to_cluster_i_size(int nb_kernel_type)
71 {
72     switch (nb_kernel_type)
73     {
74         case nbnxnk4x4_PlainC:
75         case nbnxnk4xN_SIMD_4xN:
76         case nbnxnk4xN_SIMD_2xNN:
77             return c_nbnxnCpuIClusterSize;
78         case nbnxnk8x8x8_GPU:
79         case nbnxnk8x8x8_PlainC:
80             /* The cluster size for super/sub lists is only set here.
81              * Any value should work for the pair-search and atomdata code.
82              * The kernels, of course, might require a particular value.
83              */
84             return c_nbnxnGpuClusterSize;
85         default:
86             gmx_incons("unknown kernel type");
87     }
88 }
89
90 int nbnxn_kernel_to_cluster_j_size(int nb_kernel_type)
91 {
92     int nbnxn_simd_width = 0;
93     int cj_size          = 0;
94
95 #if GMX_SIMD
96     nbnxn_simd_width = GMX_SIMD_REAL_WIDTH;
97 #endif
98
99     switch (nb_kernel_type)
100     {
101         case nbnxnk4x4_PlainC:
102             cj_size = c_nbnxnCpuIClusterSize;
103             break;
104         case nbnxnk4xN_SIMD_4xN:
105             cj_size = nbnxn_simd_width;
106             break;
107         case nbnxnk4xN_SIMD_2xNN:
108             cj_size = nbnxn_simd_width/2;
109             break;
110         case nbnxnk8x8x8_GPU:
111         case nbnxnk8x8x8_PlainC:
112             cj_size = nbnxn_kernel_to_cluster_i_size(nb_kernel_type);
113             break;
114         default:
115             gmx_incons("unknown kernel type");
116     }
117
118     return cj_size;
119 }
120
121 /* Clusters at the cut-off only increase rlist by 60% of their size */
122 static constexpr real c_nbnxnRlistIncreaseOutsideFactor = 0.6;
123
124 real nbnxn_get_rlist_effective_inc(const int  jClusterSize,
125                                    const real atomDensity)
126 {
127     /* We should get this from the setup, but currently it's the same for
128      * all setups, including GPUs.
129      */
130     const real iClusterSize    = c_nbnxnCpuIClusterSize;
131
132     const real iVolumeIncrease = (iClusterSize - 1)/atomDensity;
133     const real jVolumeIncrease = (jClusterSize - 1)/atomDensity;
134
135     return c_nbnxnRlistIncreaseOutsideFactor*std::cbrt(iVolumeIncrease +
136                                                        jVolumeIncrease);
137 }
138
139 real nbnxn_get_rlist_effective_inc(const int        clusterSize,
140                                    const gmx::RVec &averageClusterBoundingBox)
141 {
142     /* The average length of the diagonal of a sub cell */
143     const real diagonal    = std::sqrt(norm2(averageClusterBoundingBox));
144
145     const real volumeRatio = (clusterSize - 1.0_real)/clusterSize;
146
147     return c_nbnxnRlistIncreaseOutsideFactor*gmx::square(volumeRatio)*0.5_real*diagonal;
148 }