ebc2791074f31d4496ccff54529e45df24a7b260
[alexxy/gromacs.git] / src / gromacs / awh / tests / biasgrid.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,2018,2019,2020, 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 #include "gmxpre.h"
36
37 #include "gromacs/awh/biasgrid.h"
38
39 #include <cmath>
40
41 #include <memory>
42 #include <vector>
43
44 #include <gmock/gmock.h>
45 #include <gtest/gtest.h>
46
47 #include "gromacs/mdtypes/awh_params.h"
48 #include "gromacs/utility/gmxassert.h"
49
50 #include "testutils/testasserts.h"
51
52 namespace gmx
53 {
54
55 namespace test
56 {
57
58 TEST(biasGridTest, neighborhood)
59 {
60     constexpr double pointsPerScope = BiasGrid::c_scopeCutoff * BiasGrid::c_numPointsPerSigma;
61     GMX_RELEASE_ASSERT(std::abs(pointsPerScope - std::round(pointsPerScope)) > 1e-4,
62                        "If the scope is close to an integer number of points, this test can be "
63                        "unstable due to rounding issues");
64
65     const int scopeInPoints = static_cast<int>(pointsPerScope);
66
67     /* We test a 2-dimensional grid with dim0 periodic, dim1 not periodic.
68      * This test should cover most relevant cases: multi-dimension grids
69      * and non-periodic and periodic dimensions.
70      * The only thing not tested here is a periodic dimension with a gap.
71      */
72
73     const int                 numDim = 2;
74     std::vector<AwhDimParams> awhDimParams(numDim);
75
76     awhDimParams[0].origin = -5;
77     awhDimParams[0].end    = 5;
78     awhDimParams[0].period = 10;
79
80     awhDimParams[1].origin = 0.5;
81     awhDimParams[1].end    = 2.0;
82     awhDimParams[1].period = 0;
83
84     const real conversionFactor = 1;
85     const real beta             = 3.0;
86
87     /* Set up dimParams to get about 15 points along each dimension */
88     std::vector<DimParams> dimParams;
89     dimParams.emplace_back(conversionFactor, 1 / (beta * 0.7 * 0.7), beta);
90     dimParams.emplace_back(conversionFactor, 1 / (beta * 0.1 * 0.1), beta);
91
92     BiasGrid grid(dimParams, awhDimParams.data());
93
94     const int numPoints = grid.numPoints();
95
96     int numPointsDim[numDim];
97     for (int d = 0; d < numDim; d++)
98     {
99         numPointsDim[d] = grid.axis(d).numPoints();
100
101         GMX_RELEASE_ASSERT(numPointsDim[d] >= 2 * scopeInPoints + 1,
102                            "The test code currently assume that the grid is larger of equal to the "
103                            "scope in each dimension");
104     }
105     int halfNumPoints0 = numPointsDim[0] / 2;
106
107     bool haveOutOfGridNeighbors  = false;
108     bool haveDuplicateNeighbors  = false;
109     bool haveIncorrectNeighbors  = false;
110     bool haveCorrectNumNeighbors = true;
111
112     /* Set up a grid for checking for duplicate neighbors */
113     std::vector<bool> isInNeighborhood(grid.numPoints(), false);
114
115     /* Checking for all points is overkill, we check every 7th */
116     for (size_t i = 0; i < grid.numPoints(); i += 7)
117     {
118         const GridPoint& point = grid.point(i);
119
120         /* NOTE: This code relies on major-minor index ordering in Grid */
121         int pointIndex0 = i / numPointsDim[1];
122         int pointIndex1 = i - pointIndex0 * numPointsDim[1];
123
124         /* To check if we have the correct neighbors, we check the expected
125          * number of neighbors, if all neighbors are within the grid bounds
126          * and if they are within scope.
127          */
128         int    distanceFromEdge1 = std::min(pointIndex1, numPointsDim[1] - 1 - pointIndex1);
129         size_t numNeighbors      = (2 * scopeInPoints + 1)
130                               * (scopeInPoints + std::min(scopeInPoints, distanceFromEdge1) + 1);
131         if (point.neighbor.size() != numNeighbors)
132         {
133             haveCorrectNumNeighbors = false;
134         }
135
136         for (auto& j : point.neighbor)
137         {
138             if (j >= 0 && j < numPoints)
139             {
140                 if (isInNeighborhood[j])
141                 {
142                     haveDuplicateNeighbors = true;
143                 }
144                 isInNeighborhood[j] = true;
145
146                 int neighborIndex0 = j / numPointsDim[1];
147                 int neighborIndex1 = j - neighborIndex0 * numPointsDim[1];
148
149                 int distance0 = neighborIndex0 - pointIndex0;
150                 int distance1 = neighborIndex1 - pointIndex1;
151                 /* Adjust distance for periodicity of dimension 0 */
152                 if (distance0 < -halfNumPoints0)
153                 {
154                     distance0 += numPointsDim[0];
155                 }
156                 else if (distance0 > halfNumPoints0)
157                 {
158                     distance0 -= numPointsDim[0];
159                 }
160                 /* Check if the distance is within scope */
161                 if (distance0 < -scopeInPoints || distance0 > scopeInPoints
162                     || distance1 < -scopeInPoints || distance1 > scopeInPoints)
163                 {
164                     haveIncorrectNeighbors = true;
165                 }
166             }
167             else
168             {
169                 haveOutOfGridNeighbors = true;
170             }
171         }
172
173         /* Clear the marked points in the checking grid */
174         for (auto& neighbor : point.neighbor)
175         {
176             if (neighbor >= 0 && neighbor < numPoints)
177             {
178                 isInNeighborhood[neighbor] = false;
179             }
180         }
181     }
182
183     EXPECT_FALSE(haveOutOfGridNeighbors);
184     EXPECT_FALSE(haveDuplicateNeighbors);
185     EXPECT_FALSE(haveIncorrectNeighbors);
186     EXPECT_TRUE(haveCorrectNumNeighbors);
187 }
188
189 } // namespace test
190 } // namespace gmx