73a874e620366bd5d281f48c68be691384686052
[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, 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 /*! \internal \file
36  * \brief
37  * Tests for the surface area calculation used by the `sasa` analysis module.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_trajectoryanalysis
41  */
42 #include "gmxpre.h"
43
44 #include "gromacs/trajectoryanalysis/modules/surfacearea.h"
45
46 #include <cstdlib>
47
48 #include <gtest/gtest.h>
49
50 #include "gromacs/math/utilities.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/pbcutil/pbc.h"
53 #include "gromacs/random/threefry.h"
54 #include "gromacs/random/uniformrealdistribution.h"
55 #include "gromacs/utility/arrayref.h"
56 #include "gromacs/utility/gmxassert.h"
57 #include "gromacs/utility/smalloc.h"
58
59 #include "testutils/refdata.h"
60 #include "testutils/testasserts.h"
61
62 namespace
63 {
64
65 /********************************************************************
66  * SurfaceAreaTest
67  */
68
69 class SurfaceAreaTest : public ::testing::Test
70 {
71     public:
72         SurfaceAreaTest()
73             : rng_(12345), area_(0.0), volume_(0.0),
74               atomArea_(nullptr), dotCount_(0), dots_(nullptr)
75         {
76             clear_mat(box_);
77         }
78         ~SurfaceAreaTest()
79         {
80             sfree(atomArea_);
81             sfree(dots_);
82         }
83
84         void addSphere(real x, real y, real z, real radius,
85                        bool bAddToIndex = true)
86         {
87             if (bAddToIndex)
88             {
89                 index_.push_back(x_.size());
90             }
91             x_.emplace_back(x, y, z);
92             radius_.push_back(radius);
93         }
94
95         void generateRandomPosition(rvec x, real *radius)
96         {
97             rvec fx;
98             gmx::UniformRealDistribution<real> dist;
99
100             fx[XX]  = dist(rng_);
101             fx[YY]  = dist(rng_);
102             fx[ZZ]  = dist(rng_);
103             mvmul(box_, fx, x);
104             *radius = 1.5*dist(rng_) + 0.5;
105         }
106
107         void generateRandomPositions(int count)
108         {
109             x_.reserve(count);
110             radius_.reserve(count);
111             index_.reserve(count);
112             for (int i = 0; i < count; ++i)
113             {
114                 rvec x;
115                 real radius;
116                 generateRandomPosition(x, &radius);
117                 addSphere(x[XX], x[YY], x[ZZ], radius);
118             }
119         }
120         void translatePoints(real x, real y, real z)
121         {
122             for (size_t i = 0; i < x_.size(); ++i)
123             {
124                 x_[i][XX] += x;
125                 x_[i][YY] += y;
126                 x_[i][ZZ] += z;
127             }
128         }
129
130         void calculate(int ndots, int flags, bool bPBC)
131         {
132             volume_   = 0.0;
133             sfree(atomArea_);
134             atomArea_ = nullptr;
135             dotCount_ = 0;
136             sfree(dots_);
137             dots_     = nullptr;
138             t_pbc       pbc;
139             if (bPBC)
140             {
141                 set_pbc(&pbc, epbcXYZ, box_);
142             }
143             ASSERT_NO_THROW_GMX(
144                     {
145                         gmx::SurfaceAreaCalculator calculator;
146                         calculator.setDotCount(ndots);
147                         calculator.setRadii(radius_);
148                         calculator.calculate(as_rvec_array(x_.data()), bPBC ? &pbc : nullptr,
149                                              index_.size(), index_.data(), flags,
150                                              &area_, &volume_, &atomArea_,
151                                              &dots_, &dotCount_);
152                     });
153         }
154         real resultArea() const { return area_; }
155         real resultVolume() const { return volume_; }
156         real atomArea(int index) const { return atomArea_[index]; }
157
158         void checkReference(gmx::test::TestReferenceChecker *checker, const char *id,
159                             bool checkDotCoordinates)
160         {
161             gmx::test::TestReferenceChecker compound(
162                     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(
230             gmx::test::defaultRealTolerance());
231     addSphere(1, 1, 1, 1);
232     ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA, false));
233     EXPECT_REAL_EQ_TOL(4*M_PI, resultArea(), tolerance);
234     EXPECT_REAL_EQ_TOL(4*M_PI, atomArea(0), tolerance);
235     EXPECT_REAL_EQ_TOL(4*M_PI/3, resultVolume(), tolerance);
236 }
237
238 TEST_F(SurfaceAreaTest, ComputesTwoPoints)
239 {
240     gmx::test::FloatingPointTolerance tolerance(
241             gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
242     addSphere(1, 1, 1, 1);
243     addSphere(2, 1, 1, 1);
244     ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
245     EXPECT_REAL_EQ_TOL(2*2*M_PI*1.5, resultArea(), tolerance);
246     EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(0), tolerance);
247     EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(1), tolerance);
248 }
249
250 TEST_F(SurfaceAreaTest, ComputesTwoPointsOfUnequalRadius)
251 {
252     gmx::test::FloatingPointTolerance tolerance(
253             gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
254     // Spheres of radius 1 and 2 with intersection at 1.5
255     const real dist = 0.5 + sqrt(3.25);
256     addSphere(1.0, 1.0, 1.0, 1);
257     addSphere(1.0 + dist, 1.0, 1.0, 2);
258     ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
259     EXPECT_REAL_EQ_TOL(2*M_PI*(1.5 + (dist - 0.5 + 2)*2), resultArea(), tolerance);
260     EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(0), tolerance);
261     EXPECT_REAL_EQ_TOL(2*M_PI*(dist - 0.5 + 2)*2, atomArea(1), tolerance);
262 }
263
264 TEST_F(SurfaceAreaTest, SurfacePoints12)
265 {
266     gmx::test::TestReferenceChecker checker(data_.rootChecker());
267     addSphere(0, 0, 0, 1);
268     ASSERT_NO_FATAL_FAILURE(calculate(12, FLAG_DOTS, false));
269     checkReference(&checker, "Surface", true);
270 }
271
272 TEST_F(SurfaceAreaTest, SurfacePoints32)
273 {
274     gmx::test::TestReferenceChecker checker(data_.rootChecker());
275     addSphere(0, 0, 0, 1);
276     ASSERT_NO_FATAL_FAILURE(calculate(32, FLAG_DOTS, false));
277     checkReference(&checker, "Surface", true);
278 }
279
280 TEST_F(SurfaceAreaTest, SurfacePoints42)
281 {
282     gmx::test::TestReferenceChecker checker(data_.rootChecker());
283     addSphere(0, 0, 0, 1);
284     ASSERT_NO_FATAL_FAILURE(calculate(42, FLAG_DOTS, false));
285     checkReference(&checker, "Surface", true);
286 }
287
288 TEST_F(SurfaceAreaTest, SurfacePoints122)
289 {
290     gmx::test::TestReferenceChecker checker(data_.rootChecker());
291     addSphere(0, 0, 0, 1);
292     ASSERT_NO_FATAL_FAILURE(calculate(122, FLAG_DOTS, false));
293     checkReference(&checker, "Surface", true);
294 }
295
296 TEST_F(SurfaceAreaTest, Computes100Points)
297 {
298     gmx::test::TestReferenceChecker checker(data_.rootChecker());
299     checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
300     box_[XX][XX] = 10.0;
301     box_[YY][YY] = 10.0;
302     box_[ZZ][ZZ] = 10.0;
303     generateRandomPositions(100);
304     ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA | FLAG_DOTS, false));
305     checkReference(&checker, "100Points", false);
306 }
307
308 TEST_F(SurfaceAreaTest, Computes100PointsWithRectangularPBC)
309 {
310     // TODO: It would be nice to check that this produces the same result as
311     // without PBC, without duplicating the reference files.
312     gmx::test::TestReferenceChecker checker(data_.rootChecker());
313     checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
314     box_[XX][XX] = 10.0;
315     box_[YY][YY] = 10.0;
316     box_[ZZ][ZZ] = 10.0;
317     generateRandomPositions(100);
318     box_[XX][XX] = 20.0;
319     box_[YY][YY] = 20.0;
320     box_[ZZ][ZZ] = 20.0;
321     const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
322     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
323     checkReference(&checker, "100Points", false);
324
325     translatePoints(15.0, 0, 0);
326     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
327     checkReference(&checker, "100Points", false);
328
329     translatePoints(-15.0, 15.0, 0);
330     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
331     checkReference(&checker, "100Points", false);
332
333     translatePoints(0, -15.0, 15.0);
334     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
335     checkReference(&checker, "100Points", false);
336 }
337
338 TEST_F(SurfaceAreaTest, Computes100PointsWithTriclinicPBC)
339 {
340     // TODO: It would be nice to check that this produces the same result as
341     // without PBC, without duplicating the reference files.
342     gmx::test::TestReferenceChecker checker(data_.rootChecker());
343     checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
344     box_[XX][XX] = 10.0;
345     box_[YY][YY] = 10.0;
346     box_[ZZ][ZZ] = 10.0;
347     generateRandomPositions(100);
348     box_[XX][XX] = 20.0;
349     box_[YY][XX] = 10.0;
350     box_[YY][YY] = 10.0*sqrt(3.0);
351     box_[ZZ][XX] = 10.0;
352     box_[ZZ][YY] = 10.0*sqrt(1.0/3.0);
353     box_[ZZ][ZZ] = 20.0*sqrt(2.0/3.0);
354
355     const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
356     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
357     checkReference(&checker, "100Points", false);
358
359     translatePoints(15.0, 0, 0);
360     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
361     checkReference(&checker, "100Points", false);
362
363     translatePoints(-15.0, box_[YY][YY] - 5.0, 0);
364     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
365     checkReference(&checker, "100Points", false);
366
367     translatePoints(0, -(box_[YY][YY] - 5.0), 15.0);
368     ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
369     checkReference(&checker, "100Points", false);
370 }
371
372 } // namespace