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