2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018,2019, 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
82 ~SurfaceAreaTest() override
88 void addSphere(real x, real y, real z, real radius, bool bAddToIndex = true)
92 index_.push_back(x_.size());
94 x_.emplace_back(x, y, z);
95 radius_.push_back(radius);
98 void generateRandomPosition(rvec x, real* radius)
101 gmx::UniformRealDistribution<real> dist;
107 *radius = 1.5 * dist(rng_) + 0.5;
110 void generateRandomPositions(int count)
113 radius_.reserve(count);
114 index_.reserve(count);
115 for (int i = 0; i < count; ++i)
119 generateRandomPosition(x, &radius);
120 addSphere(x[XX], x[YY], x[ZZ], radius);
123 void translatePoints(real x, real y, real z)
125 for (size_t i = 0; i < x_.size(); ++i)
133 void calculate(int ndots, int flags, bool bPBC)
144 set_pbc(&pbc, epbcXYZ, box_);
146 ASSERT_NO_THROW_GMX({
147 gmx::SurfaceAreaCalculator calculator;
148 calculator.setDotCount(ndots);
149 calculator.setRadii(radius_);
150 calculator.calculate(as_rvec_array(x_.data()), bPBC ? &pbc : nullptr, index_.size(),
151 index_.data(), flags, &area_, &volume_, &atomArea_, &dots_, &dotCount_);
154 real resultArea() const { return area_; }
155 real resultVolume() const { return volume_; }
156 real atomArea(int index) const { return atomArea_[index]; }
158 void checkReference(gmx::test::TestReferenceChecker* checker, const char* id, bool checkDotCoordinates)
160 gmx::test::TestReferenceChecker compound(checker->checkCompound("SASA", id));
161 compound.checkReal(area_, "Area");
164 compound.checkReal(volume_, "Volume");
166 if (atomArea_ != nullptr)
168 compound.checkSequenceArray(index_.size(), atomArea_, "AtomArea");
170 if (dots_ != nullptr)
172 if (checkDotCoordinates)
174 // The algorithm may produce the dots in different order in
175 // single and double precision due to some internal
177 std::qsort(dots_, dotCount_, sizeof(rvec), &dotComparer);
178 compound.checkSequenceArray(3 * dotCount_, dots_, "Dots");
182 compound.checkInteger(dotCount_, "DotCount");
187 gmx::test::TestReferenceData data_;
191 static int dotComparer(const void* a, const void* b)
193 for (int d = DIM - 1; d >= 0; --d)
195 const real ad = reinterpret_cast<const real*>(a)[d];
196 const real bd = reinterpret_cast<const real*>(b)[d];
197 // A fudge factor is needed to get an ordering that is the same
198 // in single and double precision, since the points are not
199 // exactly on the same Z plane even though in exact arithmetic
200 // they probably would be.
205 else if (ad > bd + 0.001)
213 gmx::DefaultRandomEngine rng_;
214 std::vector<gmx::RVec> x_;
215 std::vector<real> radius_;
216 std::vector<int> index_;
225 TEST_F(SurfaceAreaTest, ComputesSinglePoint)
227 gmx::test::FloatingPointTolerance tolerance(gmx::test::defaultRealTolerance());
228 addSphere(1, 1, 1, 1);
229 ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA, false));
230 EXPECT_REAL_EQ_TOL(4 * M_PI, resultArea(), tolerance);
231 EXPECT_REAL_EQ_TOL(4 * M_PI, atomArea(0), tolerance);
232 EXPECT_REAL_EQ_TOL(4 * M_PI / 3, resultVolume(), tolerance);
235 TEST_F(SurfaceAreaTest, ComputesTwoPoints)
237 gmx::test::FloatingPointTolerance tolerance(gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
238 addSphere(1, 1, 1, 1);
239 addSphere(2, 1, 1, 1);
240 ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
241 EXPECT_REAL_EQ_TOL(2 * 2 * M_PI * 1.5, resultArea(), tolerance);
242 EXPECT_REAL_EQ_TOL(2 * M_PI * 1.5, atomArea(0), tolerance);
243 EXPECT_REAL_EQ_TOL(2 * M_PI * 1.5, atomArea(1), tolerance);
246 TEST_F(SurfaceAreaTest, ComputesTwoPointsOfUnequalRadius)
248 gmx::test::FloatingPointTolerance tolerance(gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
249 // Spheres of radius 1 and 2 with intersection at 1.5
250 const real dist = 0.5 + sqrt(3.25);
251 addSphere(1.0, 1.0, 1.0, 1);
252 addSphere(1.0 + dist, 1.0, 1.0, 2);
253 ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
254 EXPECT_REAL_EQ_TOL(2 * M_PI * (1.5 + (dist - 0.5 + 2) * 2), resultArea(), tolerance);
255 EXPECT_REAL_EQ_TOL(2 * M_PI * 1.5, atomArea(0), tolerance);
256 EXPECT_REAL_EQ_TOL(2 * M_PI * (dist - 0.5 + 2) * 2, atomArea(1), tolerance);
259 TEST_F(SurfaceAreaTest, SurfacePoints12)
261 gmx::test::TestReferenceChecker checker(data_.rootChecker());
262 addSphere(0, 0, 0, 1);
263 ASSERT_NO_FATAL_FAILURE(calculate(12, FLAG_DOTS, false));
264 checkReference(&checker, "Surface", true);
267 TEST_F(SurfaceAreaTest, SurfacePoints32)
269 gmx::test::TestReferenceChecker checker(data_.rootChecker());
270 addSphere(0, 0, 0, 1);
271 ASSERT_NO_FATAL_FAILURE(calculate(32, FLAG_DOTS, false));
272 checkReference(&checker, "Surface", true);
275 TEST_F(SurfaceAreaTest, SurfacePoints42)
277 gmx::test::TestReferenceChecker checker(data_.rootChecker());
278 addSphere(0, 0, 0, 1);
279 ASSERT_NO_FATAL_FAILURE(calculate(42, FLAG_DOTS, false));
280 checkReference(&checker, "Surface", true);
283 TEST_F(SurfaceAreaTest, SurfacePoints122)
285 gmx::test::TestReferenceChecker checker(data_.rootChecker());
286 addSphere(0, 0, 0, 1);
287 ASSERT_NO_FATAL_FAILURE(calculate(122, FLAG_DOTS, false));
288 checkReference(&checker, "Surface", true);
291 TEST_F(SurfaceAreaTest, Computes100Points)
293 gmx::test::TestReferenceChecker checker(data_.rootChecker());
294 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
298 generateRandomPositions(100);
299 ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA | FLAG_DOTS, false));
300 checkReference(&checker, "100Points", false);
303 TEST_F(SurfaceAreaTest, Computes100PointsWithRectangularPBC)
305 // TODO: It would be nice to check that this produces the same result as
306 // without PBC, without duplicating the reference files.
307 gmx::test::TestReferenceChecker checker(data_.rootChecker());
308 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
312 generateRandomPositions(100);
316 const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
317 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
318 checkReference(&checker, "100Points", false);
320 translatePoints(15.0, 0, 0);
321 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
322 checkReference(&checker, "100Points", false);
324 translatePoints(-15.0, 15.0, 0);
325 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
326 checkReference(&checker, "100Points", false);
328 translatePoints(0, -15.0, 15.0);
329 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
330 checkReference(&checker, "100Points", false);
333 TEST_F(SurfaceAreaTest, Computes100PointsWithTriclinicPBC)
335 // TODO: It would be nice to check that this produces the same result as
336 // without PBC, without duplicating the reference files.
337 gmx::test::TestReferenceChecker checker(data_.rootChecker());
338 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
342 generateRandomPositions(100);
345 box_[YY][YY] = 10.0 * sqrt(3.0);
347 box_[ZZ][YY] = 10.0 * sqrt(1.0 / 3.0);
348 box_[ZZ][ZZ] = 20.0 * sqrt(2.0 / 3.0);
350 const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
351 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
352 checkReference(&checker, "100Points", false);
354 translatePoints(15.0, 0, 0);
355 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
356 checkReference(&checker, "100Points", false);
358 translatePoints(-15.0, box_[YY][YY] - 5.0, 0);
359 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
360 checkReference(&checker, "100Points", false);
362 translatePoints(0, -(box_[YY][YY] - 5.0), 15.0);
363 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
364 checkReference(&checker, "100Points", false);