2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017, 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.
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.
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.
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.
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.
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.
37 * Tests for the surface area calculation used by the `sasa` analysis module.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_trajectoryanalysis
44 #include "gromacs/trajectoryanalysis/modules/surfacearea.h"
48 #include <gtest/gtest.h>
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"
59 #include "testutils/refdata.h"
60 #include "testutils/testasserts.h"
65 /********************************************************************
69 class SurfaceAreaTest : public ::testing::Test
73 : rng_(12345), area_(0.0), volume_(0.0),
74 atomArea_(nullptr), dotCount_(0), dots_(nullptr)
84 void addSphere(real x, real y, real z, real radius,
85 bool bAddToIndex = true)
89 index_.push_back(x_.size());
91 x_.emplace_back(x, y, z);
92 radius_.push_back(radius);
95 void generateRandomPosition(rvec x, real *radius)
98 gmx::UniformRealDistribution<real> dist;
104 *radius = 1.5*dist(rng_) + 0.5;
107 void addDummySpheres(int count)
109 for (int i = 0; i < count; ++i)
113 generateRandomPosition(x, &radius);
114 addSphere(x[XX], x[YY], x[ZZ], radius, false);
118 void generateRandomPositions(int count)
121 radius_.reserve(count);
122 index_.reserve(count);
123 for (int i = 0; i < count; ++i)
127 generateRandomPosition(x, &radius);
128 addSphere(x[XX], x[YY], x[ZZ], radius);
131 void translatePoints(real x, real y, real z)
133 for (size_t i = 0; i < x_.size(); ++i)
141 void calculate(int ndots, int flags, bool bPBC)
152 set_pbc(&pbc, epbcXYZ, box_);
156 gmx::SurfaceAreaCalculator calculator;
157 calculator.setDotCount(ndots);
158 calculator.setRadii(radius_);
159 calculator.calculate(as_rvec_array(x_.data()), bPBC ? &pbc : nullptr,
160 index_.size(), index_.data(), flags,
161 &area_, &volume_, &atomArea_,
165 real resultArea() const { return area_; }
166 real resultVolume() const { return volume_; }
167 real atomArea(int index) const { return atomArea_[index]; }
169 void checkReference(gmx::test::TestReferenceChecker *checker, const char *id,
170 bool checkDotCoordinates)
172 gmx::test::TestReferenceChecker compound(
173 checker->checkCompound("SASA", id));
174 compound.checkReal(area_, "Area");
177 compound.checkReal(volume_, "Volume");
179 if (atomArea_ != nullptr)
181 compound.checkSequenceArray(index_.size(), atomArea_, "AtomArea");
183 if (dots_ != nullptr)
185 if (checkDotCoordinates)
187 // The algorithm may produce the dots in different order in
188 // single and double precision due to some internal
190 std::qsort(dots_, dotCount_, sizeof(rvec), &dotComparer);
191 compound.checkSequenceArray(3*dotCount_, dots_, "Dots");
195 compound.checkInteger(dotCount_, "DotCount");
200 gmx::test::TestReferenceData data_;
204 static int dotComparer(const void *a, const void *b)
206 for (int d = DIM - 1; d >= 0; --d)
208 const real ad = reinterpret_cast<const real *>(a)[d];
209 const real bd = reinterpret_cast<const real *>(b)[d];
210 // A fudge factor is needed to get an ordering that is the same
211 // in single and double precision, since the points are not
212 // exactly on the same Z plane even though in exact arithmetic
213 // they probably would be.
218 else if (ad > bd + 0.001)
226 gmx::DefaultRandomEngine rng_;
227 std::vector<gmx::RVec> x_;
228 std::vector<real> radius_;
229 std::vector<int> index_;
238 TEST_F(SurfaceAreaTest, ComputesSinglePoint)
240 gmx::test::FloatingPointTolerance tolerance(
241 gmx::test::defaultRealTolerance());
242 addSphere(1, 1, 1, 1);
243 ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA, false));
244 EXPECT_REAL_EQ_TOL(4*M_PI, resultArea(), tolerance);
245 EXPECT_REAL_EQ_TOL(4*M_PI, atomArea(0), tolerance);
246 EXPECT_REAL_EQ_TOL(4*M_PI/3, resultVolume(), tolerance);
249 TEST_F(SurfaceAreaTest, ComputesTwoPoints)
251 gmx::test::FloatingPointTolerance tolerance(
252 gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
253 addSphere(1, 1, 1, 1);
254 addSphere(2, 1, 1, 1);
255 ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
256 EXPECT_REAL_EQ_TOL(2*2*M_PI*1.5, resultArea(), tolerance);
257 EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(0), tolerance);
258 EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(1), tolerance);
261 TEST_F(SurfaceAreaTest, ComputesTwoPointsOfUnequalRadius)
263 gmx::test::FloatingPointTolerance tolerance(
264 gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
265 // Spheres of radius 1 and 2 with intersection at 1.5
266 const real dist = 0.5 + sqrt(3.25);
267 addSphere(1.0, 1.0, 1.0, 1);
268 addSphere(1.0 + dist, 1.0, 1.0, 2);
269 ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
270 EXPECT_REAL_EQ_TOL(2*M_PI*(1.5 + (dist - 0.5 + 2)*2), resultArea(), tolerance);
271 EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(0), tolerance);
272 EXPECT_REAL_EQ_TOL(2*M_PI*(dist - 0.5 + 2)*2, atomArea(1), tolerance);
275 TEST_F(SurfaceAreaTest, SurfacePoints12)
277 gmx::test::TestReferenceChecker checker(data_.rootChecker());
278 addSphere(0, 0, 0, 1);
279 ASSERT_NO_FATAL_FAILURE(calculate(12, FLAG_DOTS, false));
280 checkReference(&checker, "Surface", true);
283 TEST_F(SurfaceAreaTest, SurfacePoints32)
285 gmx::test::TestReferenceChecker checker(data_.rootChecker());
286 addSphere(0, 0, 0, 1);
287 ASSERT_NO_FATAL_FAILURE(calculate(32, FLAG_DOTS, false));
288 checkReference(&checker, "Surface", true);
291 TEST_F(SurfaceAreaTest, SurfacePoints42)
293 gmx::test::TestReferenceChecker checker(data_.rootChecker());
294 addSphere(0, 0, 0, 1);
295 ASSERT_NO_FATAL_FAILURE(calculate(42, FLAG_DOTS, false));
296 checkReference(&checker, "Surface", true);
299 TEST_F(SurfaceAreaTest, SurfacePoints122)
301 gmx::test::TestReferenceChecker checker(data_.rootChecker());
302 addSphere(0, 0, 0, 1);
303 ASSERT_NO_FATAL_FAILURE(calculate(122, FLAG_DOTS, false));
304 checkReference(&checker, "Surface", true);
307 TEST_F(SurfaceAreaTest, Computes100Points)
309 gmx::test::TestReferenceChecker checker(data_.rootChecker());
310 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
314 generateRandomPositions(100);
315 ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA | FLAG_DOTS, false));
316 checkReference(&checker, "100Points", false);
319 TEST_F(SurfaceAreaTest, Computes100PointsWithRectangularPBC)
321 // TODO: It would be nice to check that this produces the same result as
322 // without PBC, without duplicating the reference files.
323 gmx::test::TestReferenceChecker checker(data_.rootChecker());
324 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
328 generateRandomPositions(100);
332 const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
333 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
334 checkReference(&checker, "100Points", false);
336 translatePoints(15.0, 0, 0);
337 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
338 checkReference(&checker, "100Points", false);
340 translatePoints(-15.0, 15.0, 0);
341 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
342 checkReference(&checker, "100Points", false);
344 translatePoints(0, -15.0, 15.0);
345 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
346 checkReference(&checker, "100Points", false);
349 TEST_F(SurfaceAreaTest, Computes100PointsWithTriclinicPBC)
351 // TODO: It would be nice to check that this produces the same result as
352 // without PBC, without duplicating the reference files.
353 gmx::test::TestReferenceChecker checker(data_.rootChecker());
354 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
358 generateRandomPositions(100);
361 box_[YY][YY] = 10.0*sqrt(3.0);
363 box_[ZZ][YY] = 10.0*sqrt(1.0/3.0);
364 box_[ZZ][ZZ] = 20.0*sqrt(2.0/3.0);
366 const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
367 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
368 checkReference(&checker, "100Points", false);
370 translatePoints(15.0, 0, 0);
371 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
372 checkReference(&checker, "100Points", false);
374 translatePoints(-15.0, box_[YY][YY] - 5.0, 0);
375 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
376 checkReference(&checker, "100Points", false);
378 translatePoints(0, -(box_[YY][YY] - 5.0), 15.0);
379 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
380 checkReference(&checker, "100Points", false);