2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014, 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 selection neighborhood searching.
40 * Increase coverage of these tests for different corner cases: other PBC cases
41 * than full 3D, large cutoffs (larger than half the box size), etc.
42 * At least some of these probably don't work correctly.
44 * \author Teemu Murtola <teemu.murtola@gmail.com>
45 * \ingroup module_selection
47 #include <gtest/gtest.h>
55 #include "gromacs/legacyheaders/pbc.h"
56 #include "gromacs/legacyheaders/vec.h"
58 #include "gromacs/selection/nbsearch.h"
59 #include "gromacs/random/random.h"
60 #include "gromacs/utility/smalloc.h"
62 #include "testutils/testasserts.h"
67 /********************************************************************
68 * NeighborhoodSearchTestData
71 class NeighborhoodSearchTestData
76 TestPosition() : refMinDist(0.0), refNearestPoint(-1)
80 explicit TestPosition(const rvec x)
81 : refMinDist(0.0), refNearestPoint(-1)
83 copy_rvec(x, this->x);
89 std::set<int> refPairs;
91 typedef std::vector<TestPosition> TestPositionList;
93 NeighborhoodSearchTestData(int seed, real cutoff);
94 ~NeighborhoodSearchTestData();
96 gmx::AnalysisNeighborhoodPositions refPositions() const
98 return gmx::AnalysisNeighborhoodPositions(refPos_, refPosCount_);
100 gmx::AnalysisNeighborhoodPositions testPositions() const
102 if (testPos_ == NULL)
104 snew(testPos_, testPositions_.size());
105 for (size_t i = 0; i < testPositions_.size(); ++i)
107 copy_rvec(testPositions_[i].x, testPos_[i]);
110 return gmx::AnalysisNeighborhoodPositions(testPos_,
111 testPositions_.size());
113 gmx::AnalysisNeighborhoodPositions testPosition(int index) const
115 return testPositions().selectSingleFromArray(index);
118 void addTestPosition(const rvec x)
120 GMX_RELEASE_ASSERT(testPos_ == NULL,
121 "Cannot add positions after testPositions() call");
122 testPositions_.push_back(TestPosition(x));
124 void generateRandomPosition(rvec x);
125 void generateRandomRefPositions(int count);
126 void generateRandomTestPositions(int count);
127 void computeReferences(t_pbc *pbc);
135 TestPositionList testPositions_;
138 mutable rvec *testPos_;
141 NeighborhoodSearchTestData::NeighborhoodSearchTestData(int seed, real cutoff)
142 : rng_(NULL), cutoff_(cutoff), refPosCount_(0), refPos_(NULL), testPos_(NULL)
144 // TODO: Handle errors.
145 rng_ = gmx_rng_init(seed);
147 set_pbc(&pbc_, epbcNONE, box_);
150 NeighborhoodSearchTestData::~NeighborhoodSearchTestData()
154 gmx_rng_destroy(rng_);
160 void NeighborhoodSearchTestData::generateRandomPosition(rvec x)
163 fx[XX] = gmx_rng_uniform_real(rng_);
164 fx[YY] = gmx_rng_uniform_real(rng_);
165 fx[ZZ] = gmx_rng_uniform_real(rng_);
167 // Add a small displacement to allow positions outside the box
168 x[XX] += 0.2 * gmx_rng_uniform_real(rng_) - 0.1;
169 x[YY] += 0.2 * gmx_rng_uniform_real(rng_) - 0.1;
170 x[ZZ] += 0.2 * gmx_rng_uniform_real(rng_) - 0.1;
173 void NeighborhoodSearchTestData::generateRandomRefPositions(int count)
175 refPosCount_ = count;
176 snew(refPos_, refPosCount_);
177 for (int i = 0; i < refPosCount_; ++i)
179 generateRandomPosition(refPos_[i]);
183 void NeighborhoodSearchTestData::generateRandomTestPositions(int count)
185 testPositions_.reserve(count);
186 for (int i = 0; i < count; ++i)
189 generateRandomPosition(x);
194 void NeighborhoodSearchTestData::computeReferences(t_pbc *pbc)
196 real cutoff = cutoff_;
199 cutoff = std::numeric_limits<real>::max();
201 TestPositionList::iterator i;
202 for (i = testPositions_.begin(); i != testPositions_.end(); ++i)
204 i->refMinDist = cutoff;
205 i->refNearestPoint = -1;
207 for (int j = 0; j < refPosCount_; ++j)
212 pbc_dx(pbc, i->x, refPos_[j], dx);
216 rvec_sub(i->x, refPos_[j], dx);
218 const real dist = norm(dx);
219 if (dist < i->refMinDist)
221 i->refMinDist = dist;
222 i->refNearestPoint = j;
226 i->refPairs.insert(j);
232 /********************************************************************
233 * NeighborhoodSearchTest
236 class NeighborhoodSearchTest : public ::testing::Test
239 void testIsWithin(gmx::AnalysisNeighborhoodSearch *search,
240 const NeighborhoodSearchTestData &data);
241 void testMinimumDistance(gmx::AnalysisNeighborhoodSearch *search,
242 const NeighborhoodSearchTestData &data);
243 void testNearestPoint(gmx::AnalysisNeighborhoodSearch *search,
244 const NeighborhoodSearchTestData &data);
245 void testPairSearch(gmx::AnalysisNeighborhoodSearch *search,
246 const NeighborhoodSearchTestData &data);
248 gmx::AnalysisNeighborhood nb_;
251 void NeighborhoodSearchTest::testIsWithin(
252 gmx::AnalysisNeighborhoodSearch *search,
253 const NeighborhoodSearchTestData &data)
255 NeighborhoodSearchTestData::TestPositionList::const_iterator i;
256 for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
258 const bool bWithin = (i->refMinDist <= data.cutoff_);
259 EXPECT_EQ(bWithin, search->isWithin(i->x))
260 << "Distance is " << i->refMinDist;
264 void NeighborhoodSearchTest::testMinimumDistance(
265 gmx::AnalysisNeighborhoodSearch *search,
266 const NeighborhoodSearchTestData &data)
268 NeighborhoodSearchTestData::TestPositionList::const_iterator i;
269 for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
271 const real refDist = i->refMinDist;
272 EXPECT_REAL_EQ_TOL(refDist, search->minimumDistance(i->x),
273 gmx::test::ulpTolerance(20));
277 void NeighborhoodSearchTest::testNearestPoint(
278 gmx::AnalysisNeighborhoodSearch *search,
279 const NeighborhoodSearchTestData &data)
281 NeighborhoodSearchTestData::TestPositionList::const_iterator i;
282 for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
284 const gmx::AnalysisNeighborhoodPair pair = search->nearestPoint(i->x);
287 EXPECT_EQ(i->refNearestPoint, pair.refIndex());
288 EXPECT_EQ(0, pair.testIndex());
292 EXPECT_EQ(i->refNearestPoint, -1);
297 void NeighborhoodSearchTest::testPairSearch(
298 gmx::AnalysisNeighborhoodSearch *search,
299 const NeighborhoodSearchTestData &data)
301 NeighborhoodSearchTestData::TestPositionList::const_iterator i;
302 for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
304 std::set<int> checkSet = i->refPairs;
305 gmx::AnalysisNeighborhoodPairSearch pairSearch =
306 search->startPairSearch(i->x);
307 gmx::AnalysisNeighborhoodPair pair;
308 while (pairSearch.findNextPair(&pair))
310 EXPECT_EQ(0, pair.testIndex());
311 if (checkSet.erase(pair.refIndex()) == 0)
313 // TODO: Check whether the same pair was returned more than
314 // once and give a better error message if so.
316 << "Expected: Position " << pair.refIndex()
317 << " is within cutoff.\n"
318 << " Actual: It is not.";
321 EXPECT_TRUE(checkSet.empty()) << "Some positions were not returned by the pair search.";
325 /********************************************************************
326 * Test data generation
329 class TrivialTestData
332 static const NeighborhoodSearchTestData &get()
334 static TrivialTestData singleton;
335 return singleton.data_;
338 TrivialTestData() : data_(12345, 1.0)
340 data_.box_[XX][XX] = 5.0;
341 data_.box_[YY][YY] = 5.0;
342 data_.box_[ZZ][ZZ] = 5.0;
343 data_.generateRandomRefPositions(10);
344 data_.generateRandomTestPositions(5);
345 set_pbc(&data_.pbc_, epbcXYZ, data_.box_);
346 data_.computeReferences(&data_.pbc_);
350 NeighborhoodSearchTestData data_;
353 class RandomBoxFullPBCData
356 static const NeighborhoodSearchTestData &get()
358 static RandomBoxFullPBCData singleton;
359 return singleton.data_;
362 RandomBoxFullPBCData() : data_(12345, 1.0)
364 data_.box_[XX][XX] = 10.0;
365 data_.box_[YY][YY] = 5.0;
366 data_.box_[ZZ][ZZ] = 7.0;
367 // TODO: Consider whether manually picking some positions would give better
369 data_.generateRandomRefPositions(1000);
370 data_.generateRandomTestPositions(100);
371 set_pbc(&data_.pbc_, epbcXYZ, data_.box_);
372 data_.computeReferences(&data_.pbc_);
376 NeighborhoodSearchTestData data_;
379 class RandomTriclinicFullPBCData
382 static const NeighborhoodSearchTestData &get()
384 static RandomTriclinicFullPBCData singleton;
385 return singleton.data_;
388 RandomTriclinicFullPBCData() : data_(12345, 1.0)
390 data_.box_[XX][XX] = 5.0;
391 data_.box_[YY][XX] = 2.5;
392 data_.box_[YY][YY] = 2.5*sqrt(3.0);
393 data_.box_[ZZ][XX] = 2.5;
394 data_.box_[ZZ][YY] = 2.5*sqrt(1.0/3.0);
395 data_.box_[ZZ][ZZ] = 5.0*sqrt(2.0/3.0);
396 // TODO: Consider whether manually picking some positions would give better
398 data_.generateRandomRefPositions(1000);
399 data_.generateRandomTestPositions(100);
400 set_pbc(&data_.pbc_, epbcXYZ, data_.box_);
401 data_.computeReferences(&data_.pbc_);
405 NeighborhoodSearchTestData data_;
408 class RandomBox2DPBCData
411 static const NeighborhoodSearchTestData &get()
413 static RandomBox2DPBCData singleton;
414 return singleton.data_;
417 RandomBox2DPBCData() : data_(12345, 1.0)
419 data_.box_[XX][XX] = 10.0;
420 data_.box_[YY][YY] = 7.0;
421 data_.box_[ZZ][ZZ] = 5.0;
422 // TODO: Consider whether manually picking some positions would give better
424 data_.generateRandomRefPositions(1000);
425 data_.generateRandomTestPositions(100);
426 set_pbc(&data_.pbc_, epbcXY, data_.box_);
427 data_.computeReferences(&data_.pbc_);
431 NeighborhoodSearchTestData data_;
434 /********************************************************************
438 TEST_F(NeighborhoodSearchTest, SimpleSearch)
440 const NeighborhoodSearchTestData &data = RandomBoxFullPBCData::get();
442 nb_.setCutoff(data.cutoff_);
443 nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Simple);
444 gmx::AnalysisNeighborhoodSearch search =
445 nb_.initSearch(&data.pbc_, data.refPositions());
446 ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Simple, search.mode());
448 testIsWithin(&search, data);
449 testMinimumDistance(&search, data);
450 testNearestPoint(&search, data);
451 testPairSearch(&search, data);
454 TEST_F(NeighborhoodSearchTest, GridSearchBox)
456 const NeighborhoodSearchTestData &data = RandomBoxFullPBCData::get();
458 nb_.setCutoff(data.cutoff_);
459 nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
460 gmx::AnalysisNeighborhoodSearch search =
461 nb_.initSearch(&data.pbc_, data.refPositions());
462 ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
464 testIsWithin(&search, data);
465 testMinimumDistance(&search, data);
466 testNearestPoint(&search, data);
467 testPairSearch(&search, data);
470 TEST_F(NeighborhoodSearchTest, GridSearchTriclinic)
472 const NeighborhoodSearchTestData &data = RandomTriclinicFullPBCData::get();
474 nb_.setCutoff(data.cutoff_);
475 nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
476 gmx::AnalysisNeighborhoodSearch search =
477 nb_.initSearch(&data.pbc_, data.refPositions());
478 ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
480 testPairSearch(&search, data);
483 TEST_F(NeighborhoodSearchTest, GridSearch2DPBC)
485 const NeighborhoodSearchTestData &data = RandomBox2DPBCData::get();
487 nb_.setCutoff(data.cutoff_);
488 nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
489 gmx::AnalysisNeighborhoodSearch search =
490 nb_.initSearch(&data.pbc_, data.refPositions());
491 // Currently, grid searching not supported with 2D PBC.
492 //ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
494 testIsWithin(&search, data);
495 testMinimumDistance(&search, data);
496 testNearestPoint(&search, data);
497 testPairSearch(&search, data);
500 TEST_F(NeighborhoodSearchTest, HandlesConcurrentSearches)
502 const NeighborhoodSearchTestData &data = TrivialTestData::get();
504 nb_.setCutoff(data.cutoff_);
505 gmx::AnalysisNeighborhoodSearch search1 =
506 nb_.initSearch(&data.pbc_, data.refPositions());
507 gmx::AnalysisNeighborhoodSearch search2 =
508 nb_.initSearch(&data.pbc_, data.refPositions());
510 gmx::AnalysisNeighborhoodPairSearch pairSearch1 =
511 search1.startPairSearch(data.testPosition(0));
512 gmx::AnalysisNeighborhoodPairSearch pairSearch2 =
513 search1.startPairSearch(data.testPosition(1));
515 testPairSearch(&search2, data);
517 gmx::AnalysisNeighborhoodPair pair;
518 pairSearch1.findNextPair(&pair);
519 EXPECT_EQ(0, pair.testIndex());
520 EXPECT_TRUE(data.testPositions_[0].refPairs.count(pair.refIndex()) == 1);
522 pairSearch2.findNextPair(&pair);
523 EXPECT_EQ(1, pair.testIndex());
524 EXPECT_TRUE(data.testPositions_[1].refPairs.count(pair.refIndex()) == 1);
527 TEST_F(NeighborhoodSearchTest, HandlesSkippingPairs)
529 const NeighborhoodSearchTestData &data = TrivialTestData::get();
531 nb_.setCutoff(data.cutoff_);
532 gmx::AnalysisNeighborhoodSearch search =
533 nb_.initSearch(&data.pbc_, data.refPositions());
534 gmx::AnalysisNeighborhoodPairSearch pairSearch =
535 search.startPairSearch(data.testPositions());
536 gmx::AnalysisNeighborhoodPair pair;
537 // TODO: This test needs to be adjusted if the grid search gets optimized
538 // to loop over the test positions in cell order (first, the ordering
539 // assumption here breaks, and second, it then needs to be tested
540 // separately for simple and grid searches).
541 int currentIndex = 0;
542 while (pairSearch.findNextPair(&pair))
544 while (currentIndex < pair.testIndex())
548 EXPECT_EQ(currentIndex, pair.testIndex());
549 EXPECT_TRUE(data.testPositions_[currentIndex].refPairs.count(pair.refIndex()) == 1);
550 pairSearch.skipRemainingPairsForTestPosition();