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