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