Apply clang-format-11
[alexxy/gromacs.git] / src / gromacs / trajectoryanalysis / tests / surfacearea.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5  * Copyright (c) 2019,2020,2021, 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.
9  *
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.
14  *
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.
19  *
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.
24  *
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.
32  *
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.
35  */
36 /*! \internal \file
37  * \brief
38  * Tests for the surface area calculation used by the `sasa` analysis module.
39  *
40  * \author Teemu Murtola <teemu.murtola@gmail.com>
41  * \ingroup module_trajectoryanalysis
42  */
43 #include "gmxpre.h"
44
45 #include "gromacs/trajectoryanalysis/modules/surfacearea.h"
46
47 #include <cstdlib>
48
49 #include <gtest/gtest.h>
50
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/random/threefry.h"
56 #include "gromacs/random/uniformrealdistribution.h"
57 #include "gromacs/utility/arrayref.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/smalloc.h"
60
61 #include "testutils/refdata.h"
62 #include "testutils/testasserts.h"
63
64 namespace
65 {
66
67 /********************************************************************
68  * SurfaceAreaTest
69  */
70
71 class SurfaceAreaTest : public ::testing::Test
72 {
73 public:
74     SurfaceAreaTest() :
75         box_(), rng_(12345), area_(0.0), volume_(0.0), atomArea_(nullptr), dotCount_(0), dots_(nullptr)
76     {
77     }
78     ~SurfaceAreaTest() override
79     {
80         sfree(atomArea_);
81         sfree(dots_);
82     }
83
84     void addSphere(real x, real y, real z, real radius, bool bAddToIndex = true)
85     {
86         if (bAddToIndex)
87         {
88             index_.push_back(x_.size());
89         }
90         x_.emplace_back(x, y, z);
91         radius_.push_back(radius);
92     }
93
94     void generateRandomPosition(rvec x, real* radius)
95     {
96         rvec                               fx;
97         gmx::UniformRealDistribution<real> dist;
98
99         fx[XX] = dist(rng_);
100         fx[YY] = dist(rng_);
101         fx[ZZ] = dist(rng_);
102         mvmul(box_, fx, x);
103         *radius = 1.5 * dist(rng_) + 0.5;
104     }
105
106     void generateRandomPositions(int count)
107     {
108         x_.reserve(count);
109         radius_.reserve(count);
110         index_.reserve(count);
111         for (int i = 0; i < count; ++i)
112         {
113             rvec x;
114             real radius = 0;
115             generateRandomPosition(x, &radius);
116             addSphere(x[XX], x[YY], x[ZZ], radius);
117         }
118     }
119     void translatePoints(real x, real y, real z)
120     {
121         for (size_t i = 0; i < x_.size(); ++i)
122         {
123             x_[i][XX] += x;
124             x_[i][YY] += y;
125             x_[i][ZZ] += z;
126         }
127     }
128
129     void calculate(int ndots, int flags, bool bPBC)
130     {
131         volume_ = 0.0;
132         sfree(atomArea_);
133         atomArea_ = nullptr;
134         dotCount_ = 0;
135         sfree(dots_);
136         dots_ = nullptr;
137         t_pbc pbc;
138         if (bPBC)
139         {
140             set_pbc(&pbc, PbcType::Xyz, box_);
141         }
142         gmx::SurfaceAreaCalculator calculator;
143         calculator.setDotCount(ndots);
144         calculator.setRadii(radius_);
145         calculator.calculate(as_rvec_array(x_.data()),
146                              bPBC ? &pbc : nullptr,
147                              index_.size(),
148                              index_.data(),
149                              flags,
150                              &area_,
151                              &volume_,
152                              &atomArea_,
153                              &dots_,
154                              &dotCount_);
155     }
156     real resultArea() const { return area_; }
157     real resultVolume() const { return volume_; }
158     real atomArea(int index) const { return atomArea_[index]; }
159
160     void checkReference(gmx::test::TestReferenceChecker* checker, const char* id, bool checkDotCoordinates)
161     {
162         gmx::test::TestReferenceChecker compound(checker->checkCompound("SASA", id));
163         compound.checkReal(area_, "Area");
164         if (volume_ > 0.0)
165         {
166             compound.checkReal(volume_, "Volume");
167         }
168         if (atomArea_ != nullptr)
169         {
170             compound.checkSequenceArray(index_.size(), atomArea_, "AtomArea");
171         }
172         if (dots_ != nullptr)
173         {
174             if (checkDotCoordinates)
175             {
176                 // The algorithm may produce the dots in different order in
177                 // single and double precision due to some internal
178                 // sorting...
179                 std::qsort(dots_, dotCount_, sizeof(rvec), &dotComparer);
180                 compound.checkSequenceArray(3 * dotCount_, dots_, "Dots");
181             }
182             else
183             {
184                 compound.checkInteger(dotCount_, "DotCount");
185             }
186         }
187     }
188
189     gmx::test::TestReferenceData data_;
190     matrix                       box_;
191
192 private:
193     static int dotComparer(const void* a, const void* b)
194     {
195         for (int d = DIM - 1; d >= 0; --d)
196         {
197             const real ad = reinterpret_cast<const real*>(a)[d];
198             const real bd = reinterpret_cast<const real*>(b)[d];
199             // A fudge factor is needed to get an ordering that is the same
200             // in single and double precision, since the points are not
201             // exactly on the same Z plane even though in exact arithmetic
202             // they probably would be.
203             if (ad < bd - 0.001)
204             {
205                 return -1;
206             }
207             else if (ad > bd + 0.001)
208             {
209                 return 1;
210             }
211         }
212         return 0;
213     }
214
215     gmx::DefaultRandomEngine rng_;
216     std::vector<gmx::RVec>   x_;
217     std::vector<real>        radius_;
218     std::vector<int>         index_;
219
220     real  area_;
221     real  volume_;
222     real* atomArea_;
223     int   dotCount_;
224     real* dots_;
225 };
226
227 TEST_F(SurfaceAreaTest, ComputesSinglePoint)
228 {
229     gmx::test::FloatingPointTolerance tolerance(gmx::test::defaultRealTolerance());
230     addSphere(1, 1, 1, 1);
231     ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA, false));
232     EXPECT_REAL_EQ_TOL(4 * M_PI, resultArea(), tolerance);
233     EXPECT_REAL_EQ_TOL(4 * M_PI, atomArea(0), tolerance);
234     EXPECT_REAL_EQ_TOL(4 * M_PI / 3, resultVolume(), tolerance);
235 }
236
237 TEST_F(SurfaceAreaTest, ComputesTwoPoints)
238 {
239     gmx::test::FloatingPointTolerance tolerance(gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
240     addSphere(1, 1, 1, 1);
241     addSphere(2, 1, 1, 1);
242     ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
243     EXPECT_REAL_EQ_TOL(2 * 2 * M_PI * 1.5, resultArea(), tolerance);
244     EXPECT_REAL_EQ_TOL(2 * M_PI * 1.5, atomArea(0), tolerance);
245     EXPECT_REAL_EQ_TOL(2 * M_PI * 1.5, atomArea(1), tolerance);
246 }
247
248 TEST_F(SurfaceAreaTest, ComputesTwoPointsOfUnequalRadius)
249 {
250     gmx::test::FloatingPointTolerance tolerance(gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
251     // Spheres of radius 1 and 2 with intersection at 1.5
252     const real dist = 0.5 + sqrt(3.25);
253     addSphere(1.0, 1.0, 1.0, 1);
254     addSphere(1.0 + dist, 1.0, 1.0, 2);
255     ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
256     EXPECT_REAL_EQ_TOL(2 * M_PI * (1.5 + (dist - 0.5 + 2) * 2), resultArea(), tolerance);
257     EXPECT_REAL_EQ_TOL(2 * M_PI * 1.5, atomArea(0), tolerance);
258     EXPECT_REAL_EQ_TOL(2 * M_PI * (dist - 0.5 + 2) * 2, atomArea(1), tolerance);
259 }
260
261 TEST_F(SurfaceAreaTest, SurfacePoints12)
262 {
263     gmx::test::TestReferenceChecker checker(data_.rootChecker());
264     addSphere(0, 0, 0, 1);
265     ASSERT_NO_FATAL_FAILURE(calculate(12, FLAG_DOTS, false));
266     checkReference(&checker, "Surface", true);
267 }
268
269 TEST_F(SurfaceAreaTest, SurfacePoints32)
270 {
271     gmx::test::TestReferenceChecker checker(data_.rootChecker());
272     addSphere(0, 0, 0, 1);
273     ASSERT_NO_FATAL_FAILURE(calculate(32, FLAG_DOTS, false));
274     checkReference(&checker, "Surface", true);
275 }
276
277 TEST_F(SurfaceAreaTest, SurfacePoints42)
278 {
279     gmx::test::TestReferenceChecker checker(data_.rootChecker());
280     addSphere(0, 0, 0, 1);
281     ASSERT_NO_FATAL_FAILURE(calculate(42, FLAG_DOTS, false));
282     checkReference(&checker, "Surface", true);
283 }
284
285 TEST_F(SurfaceAreaTest, SurfacePoints122)
286 {
287     gmx::test::TestReferenceChecker checker(data_.rootChecker());
288     addSphere(0, 0, 0, 1);
289     ASSERT_NO_FATAL_FAILURE(calculate(122, FLAG_DOTS, false));
290     checkReference(&checker, "Surface", true);
291 }
292
293 TEST_F(SurfaceAreaTest, Computes100Points)
294 {
295     gmx::test::TestReferenceChecker checker(data_.rootChecker());
296     checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
297     box_[XX][XX] = 10.0;
298     box_[YY][YY] = 10.0;
299     box_[ZZ][ZZ] = 10.0;
300     generateRandomPositions(100);
301     ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA | FLAG_DOTS, false));
302     checkReference(&checker, "100Points", false);
303 }
304
305 TEST_F(SurfaceAreaTest, Computes100PointsWithRectangularPBC)
306 {
307     // TODO: It would be nice to check that this produces the same result as
308     // without PBC, without duplicating the reference files.
309     gmx::test::TestReferenceChecker checker(data_.rootChecker());
310     checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
311     box_[XX][XX] = 10.0;
312     box_[YY][YY] = 10.0;
313     box_[ZZ][ZZ] = 10.0;
314     generateRandomPositions(100);
315     box_[XX][XX]    = 20.0;
316     box_[YY][YY]    = 20.0;
317     box_[ZZ][ZZ]    = 20.0;
318     const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
319     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
320     checkReference(&checker, "100Points", false);
321
322     translatePoints(15.0, 0, 0);
323     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
324     checkReference(&checker, "100Points", false);
325
326     translatePoints(-15.0, 15.0, 0);
327     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
328     checkReference(&checker, "100Points", false);
329
330     translatePoints(0, -15.0, 15.0);
331     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
332     checkReference(&checker, "100Points", false);
333 }
334
335 TEST_F(SurfaceAreaTest, Computes100PointsWithTriclinicPBC)
336 {
337     // TODO: It would be nice to check that this produces the same result as
338     // without PBC, without duplicating the reference files.
339     gmx::test::TestReferenceChecker checker(data_.rootChecker());
340     checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
341     box_[XX][XX] = 10.0;
342     box_[YY][YY] = 10.0;
343     box_[ZZ][ZZ] = 10.0;
344     generateRandomPositions(100);
345     box_[XX][XX] = 20.0;
346     box_[YY][XX] = 10.0;
347     box_[YY][YY] = 10.0 * sqrt(3.0);
348     box_[ZZ][XX] = 10.0;
349     box_[ZZ][YY] = 10.0 * sqrt(1.0 / 3.0);
350     box_[ZZ][ZZ] = 20.0 * sqrt(2.0 / 3.0);
351
352     const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
353     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
354     checkReference(&checker, "100Points", false);
355
356     translatePoints(15.0, 0, 0);
357     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
358     checkReference(&checker, "100Points", false);
359
360     translatePoints(-15.0, box_[YY][YY] - 5.0, 0);
361     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
362     checkReference(&checker, "100Points", false);
363
364     translatePoints(0, -(box_[YY][YY] - 5.0), 15.0);
365     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
366     checkReference(&checker, "100Points", false);
367 }
368
369 } // namespace