Fix tune_pme
[gromacs.git] / src / gromacs / trajectoryanalysis / tests / surfacearea.cpp
blobf1473e479ffa4dd018080c7007a3654d14db3ac4
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
35 /*! \internal \file
36 * \brief
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
42 #include "gmxpre.h"
44 #include "gromacs/trajectoryanalysis/modules/surfacearea.h"
46 #include <cstdlib>
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"
62 namespace
65 /********************************************************************
66 * SurfaceAreaTest
69 class SurfaceAreaTest : public ::testing::Test
71 public:
72 SurfaceAreaTest()
73 : box_(), rng_(12345), area_(0.0), volume_(0.0),
74 atomArea_(nullptr), dotCount_(0), dots_(nullptr)
77 ~SurfaceAreaTest() override
79 sfree(atomArea_);
80 sfree(dots_);
83 void addSphere(real x, real y, real z, real radius,
84 bool bAddToIndex = true)
86 if (bAddToIndex)
88 index_.push_back(x_.size());
90 x_.emplace_back(x, y, z);
91 radius_.push_back(radius);
94 void generateRandomPosition(rvec x, real *radius)
96 rvec fx;
97 gmx::UniformRealDistribution<real> dist;
99 fx[XX] = dist(rng_);
100 fx[YY] = dist(rng_);
101 fx[ZZ] = dist(rng_);
102 mvmul(box_, fx, x);
103 *radius = 1.5*dist(rng_) + 0.5;
106 void generateRandomPositions(int count)
108 x_.reserve(count);
109 radius_.reserve(count);
110 index_.reserve(count);
111 for (int i = 0; i < count; ++i)
113 rvec x;
114 real radius;
115 generateRandomPosition(x, &radius);
116 addSphere(x[XX], x[YY], x[ZZ], radius);
119 void translatePoints(real x, real y, real z)
121 for (size_t i = 0; i < x_.size(); ++i)
123 x_[i][XX] += x;
124 x_[i][YY] += y;
125 x_[i][ZZ] += z;
129 void calculate(int ndots, int flags, bool bPBC)
131 volume_ = 0.0;
132 sfree(atomArea_);
133 atomArea_ = nullptr;
134 dotCount_ = 0;
135 sfree(dots_);
136 dots_ = nullptr;
137 t_pbc pbc;
138 if (bPBC)
140 set_pbc(&pbc, epbcXYZ, box_);
142 ASSERT_NO_THROW_GMX(
144 gmx::SurfaceAreaCalculator calculator;
145 calculator.setDotCount(ndots);
146 calculator.setRadii(radius_);
147 calculator.calculate(as_rvec_array(x_.data()), bPBC ? &pbc : nullptr,
148 index_.size(), index_.data(), flags,
149 &area_, &volume_, &atomArea_,
150 &dots_, &dotCount_);
153 real resultArea() const { return area_; }
154 real resultVolume() const { return volume_; }
155 real atomArea(int index) const { return atomArea_[index]; }
157 void checkReference(gmx::test::TestReferenceChecker *checker, const char *id,
158 bool checkDotCoordinates)
160 gmx::test::TestReferenceChecker compound(
161 checker->checkCompound("SASA", id));
162 compound.checkReal(area_, "Area");
163 if (volume_ > 0.0)
165 compound.checkReal(volume_, "Volume");
167 if (atomArea_ != nullptr)
169 compound.checkSequenceArray(index_.size(), atomArea_, "AtomArea");
171 if (dots_ != nullptr)
173 if (checkDotCoordinates)
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");
181 else
183 compound.checkInteger(dotCount_, "DotCount");
188 gmx::test::TestReferenceData data_;
189 matrix box_;
191 private:
192 static int dotComparer(const void *a, const void *b)
194 for (int d = DIM - 1; d >= 0; --d)
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)
204 return -1;
206 else if (ad > bd + 0.001)
208 return 1;
211 return 0;
214 gmx::DefaultRandomEngine rng_;
215 std::vector<gmx::RVec> x_;
216 std::vector<real> radius_;
217 std::vector<int> index_;
219 real area_;
220 real volume_;
221 real *atomArea_;
222 int dotCount_;
223 real *dots_;
226 TEST_F(SurfaceAreaTest, ComputesSinglePoint)
228 gmx::test::FloatingPointTolerance tolerance(
229 gmx::test::defaultRealTolerance());
230 addSphere(1, 1, 1, 1);
231 ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA, false));
232 EXPECT_REAL_EQ_TOL(4*M_PI, resultArea(), tolerance);
233 EXPECT_REAL_EQ_TOL(4*M_PI, atomArea(0), tolerance);
234 EXPECT_REAL_EQ_TOL(4*M_PI/3, resultVolume(), tolerance);
237 TEST_F(SurfaceAreaTest, ComputesTwoPoints)
239 gmx::test::FloatingPointTolerance tolerance(
240 gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
241 addSphere(1, 1, 1, 1);
242 addSphere(2, 1, 1, 1);
243 ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
244 EXPECT_REAL_EQ_TOL(2*2*M_PI*1.5, resultArea(), tolerance);
245 EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(0), tolerance);
246 EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(1), tolerance);
249 TEST_F(SurfaceAreaTest, ComputesTwoPointsOfUnequalRadius)
251 gmx::test::FloatingPointTolerance tolerance(
252 gmx::test::relativeToleranceAsFloatingPoint(1.0, 0.005));
253 // Spheres of radius 1 and 2 with intersection at 1.5
254 const real dist = 0.5 + sqrt(3.25);
255 addSphere(1.0, 1.0, 1.0, 1);
256 addSphere(1.0 + dist, 1.0, 1.0, 2);
257 ASSERT_NO_FATAL_FAILURE(calculate(1000, FLAG_ATOM_AREA, false));
258 EXPECT_REAL_EQ_TOL(2*M_PI*(1.5 + (dist - 0.5 + 2)*2), resultArea(), tolerance);
259 EXPECT_REAL_EQ_TOL(2*M_PI*1.5, atomArea(0), tolerance);
260 EXPECT_REAL_EQ_TOL(2*M_PI*(dist - 0.5 + 2)*2, atomArea(1), tolerance);
263 TEST_F(SurfaceAreaTest, SurfacePoints12)
265 gmx::test::TestReferenceChecker checker(data_.rootChecker());
266 addSphere(0, 0, 0, 1);
267 ASSERT_NO_FATAL_FAILURE(calculate(12, FLAG_DOTS, false));
268 checkReference(&checker, "Surface", true);
271 TEST_F(SurfaceAreaTest, SurfacePoints32)
273 gmx::test::TestReferenceChecker checker(data_.rootChecker());
274 addSphere(0, 0, 0, 1);
275 ASSERT_NO_FATAL_FAILURE(calculate(32, FLAG_DOTS, false));
276 checkReference(&checker, "Surface", true);
279 TEST_F(SurfaceAreaTest, SurfacePoints42)
281 gmx::test::TestReferenceChecker checker(data_.rootChecker());
282 addSphere(0, 0, 0, 1);
283 ASSERT_NO_FATAL_FAILURE(calculate(42, FLAG_DOTS, false));
284 checkReference(&checker, "Surface", true);
287 TEST_F(SurfaceAreaTest, SurfacePoints122)
289 gmx::test::TestReferenceChecker checker(data_.rootChecker());
290 addSphere(0, 0, 0, 1);
291 ASSERT_NO_FATAL_FAILURE(calculate(122, FLAG_DOTS, false));
292 checkReference(&checker, "Surface", true);
295 TEST_F(SurfaceAreaTest, Computes100Points)
297 gmx::test::TestReferenceChecker checker(data_.rootChecker());
298 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
299 box_[XX][XX] = 10.0;
300 box_[YY][YY] = 10.0;
301 box_[ZZ][ZZ] = 10.0;
302 generateRandomPositions(100);
303 ASSERT_NO_FATAL_FAILURE(calculate(24, FLAG_VOLUME | FLAG_ATOM_AREA | FLAG_DOTS, false));
304 checkReference(&checker, "100Points", false);
307 TEST_F(SurfaceAreaTest, Computes100PointsWithRectangularPBC)
309 // TODO: It would be nice to check that this produces the same result as
310 // without PBC, without duplicating the reference files.
311 gmx::test::TestReferenceChecker checker(data_.rootChecker());
312 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
313 box_[XX][XX] = 10.0;
314 box_[YY][YY] = 10.0;
315 box_[ZZ][ZZ] = 10.0;
316 generateRandomPositions(100);
317 box_[XX][XX] = 20.0;
318 box_[YY][YY] = 20.0;
319 box_[ZZ][ZZ] = 20.0;
320 const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
321 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
322 checkReference(&checker, "100Points", false);
324 translatePoints(15.0, 0, 0);
325 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
326 checkReference(&checker, "100Points", false);
328 translatePoints(-15.0, 15.0, 0);
329 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
330 checkReference(&checker, "100Points", false);
332 translatePoints(0, -15.0, 15.0);
333 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
334 checkReference(&checker, "100Points", false);
337 TEST_F(SurfaceAreaTest, Computes100PointsWithTriclinicPBC)
339 // TODO: It would be nice to check that this produces the same result as
340 // without PBC, without duplicating the reference files.
341 gmx::test::TestReferenceChecker checker(data_.rootChecker());
342 checker.setDefaultTolerance(gmx::test::absoluteTolerance(0.001));
343 box_[XX][XX] = 10.0;
344 box_[YY][YY] = 10.0;
345 box_[ZZ][ZZ] = 10.0;
346 generateRandomPositions(100);
347 box_[XX][XX] = 20.0;
348 box_[YY][XX] = 10.0;
349 box_[YY][YY] = 10.0*sqrt(3.0);
350 box_[ZZ][XX] = 10.0;
351 box_[ZZ][YY] = 10.0*sqrt(1.0/3.0);
352 box_[ZZ][ZZ] = 20.0*sqrt(2.0/3.0);
354 const int flags = FLAG_ATOM_AREA | FLAG_VOLUME | FLAG_DOTS;
355 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
356 checkReference(&checker, "100Points", false);
358 translatePoints(15.0, 0, 0);
359 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
360 checkReference(&checker, "100Points", false);
362 translatePoints(-15.0, box_[YY][YY] - 5.0, 0);
363 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
364 checkReference(&checker, "100Points", false);
366 translatePoints(0, -(box_[YY][YY] - 5.0), 15.0);
367 ASSERT_NO_FATAL_FAILURE(calculate(24, flags, true));
368 checkReference(&checker, "100Points", false);
371 } // namespace