2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013, by the GROMACS development team, led by
5 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6 * others, as listed in the AUTHORS file in the top-level source
7 * 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>
54 #include "gromacs/legacyheaders/gmx_random.h"
55 #include "gromacs/legacyheaders/pbc.h"
56 #include "gromacs/legacyheaders/smalloc.h"
57 #include "gromacs/legacyheaders/vec.h"
59 #include "gromacs/selection/nbsearch.h"
61 #include "testutils/testasserts.h"
66 /********************************************************************
67 * NeighborhoodSearchTestData
70 class NeighborhoodSearchTestData
75 TestPosition() : refMinDist(0.0), refNearestPoint(-1)
79 explicit TestPosition(const rvec x)
80 : refMinDist(0.0), refNearestPoint(-1)
82 copy_rvec(x, this->x);
88 std::set<int> refPairs;
90 typedef std::vector<TestPosition> TestPositionList;
92 NeighborhoodSearchTestData(int seed, real cutoff);
93 ~NeighborhoodSearchTestData();
95 void addTestPosition(const rvec x)
97 testPositions_.push_back(TestPosition(x));
99 void generateRandomPosition(rvec x);
100 void generateRandomRefPositions(int count);
101 void generateRandomTestPositions(int count);
102 void computeReferences(t_pbc *pbc);
110 TestPositionList testPositions_;
113 NeighborhoodSearchTestData::NeighborhoodSearchTestData(int seed, real cutoff)
114 : rng_(NULL), cutoff_(cutoff), refPosCount_(0), refPos_(NULL)
116 // TODO: Handle errors.
117 rng_ = gmx_rng_init(seed);
119 set_pbc(&pbc_, epbcNONE, box_);
122 NeighborhoodSearchTestData::~NeighborhoodSearchTestData()
126 gmx_rng_destroy(rng_);
131 void NeighborhoodSearchTestData::generateRandomPosition(rvec x)
134 fx[XX] = gmx_rng_uniform_real(rng_);
135 fx[YY] = gmx_rng_uniform_real(rng_);
136 fx[ZZ] = gmx_rng_uniform_real(rng_);
138 // Add a small displacement to allow positions outside the box
139 x[XX] += 0.2 * gmx_rng_uniform_real(rng_) - 0.1;
140 x[YY] += 0.2 * gmx_rng_uniform_real(rng_) - 0.1;
141 x[ZZ] += 0.2 * gmx_rng_uniform_real(rng_) - 0.1;
144 void NeighborhoodSearchTestData::generateRandomRefPositions(int count)
146 refPosCount_ = count;
147 snew(refPos_, refPosCount_);
148 for (int i = 0; i < refPosCount_; ++i)
150 generateRandomPosition(refPos_[i]);
154 void NeighborhoodSearchTestData::generateRandomTestPositions(int count)
156 testPositions_.reserve(count);
157 for (int i = 0; i < count; ++i)
160 generateRandomPosition(x);
165 void NeighborhoodSearchTestData::computeReferences(t_pbc *pbc)
167 real cutoff = cutoff_;
170 cutoff = std::numeric_limits<real>::max();
172 TestPositionList::iterator i;
173 for (i = testPositions_.begin(); i != testPositions_.end(); ++i)
175 i->refMinDist = cutoff;
176 i->refNearestPoint = -1;
178 for (int j = 0; j < refPosCount_; ++j)
183 pbc_dx(pbc, i->x, refPos_[j], dx);
187 rvec_sub(i->x, refPos_[j], dx);
189 const real dist = norm(dx);
190 if (dist < i->refMinDist)
192 i->refMinDist = dist;
193 i->refNearestPoint = j;
197 i->refPairs.insert(j);
203 /********************************************************************
204 * NeighborhoodSearchTest
207 class NeighborhoodSearchTest : public ::testing::Test
210 void testIsWithin(gmx::AnalysisNeighborhoodSearch *search,
211 const NeighborhoodSearchTestData &data);
212 void testMinimumDistance(gmx::AnalysisNeighborhoodSearch *search,
213 const NeighborhoodSearchTestData &data);
214 void testNearestPoint(gmx::AnalysisNeighborhoodSearch *search,
215 const NeighborhoodSearchTestData &data);
216 void testPairSearch(gmx::AnalysisNeighborhoodSearch *search,
217 const NeighborhoodSearchTestData &data);
219 gmx::AnalysisNeighborhood nb_;
222 void NeighborhoodSearchTest::testIsWithin(
223 gmx::AnalysisNeighborhoodSearch *search,
224 const NeighborhoodSearchTestData &data)
226 NeighborhoodSearchTestData::TestPositionList::const_iterator i;
227 for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
229 const bool bWithin = (i->refMinDist <= data.cutoff_);
230 EXPECT_EQ(bWithin, search->isWithin(i->x))
231 << "Distance is " << i->refMinDist;
235 void NeighborhoodSearchTest::testMinimumDistance(
236 gmx::AnalysisNeighborhoodSearch *search,
237 const NeighborhoodSearchTestData &data)
239 NeighborhoodSearchTestData::TestPositionList::const_iterator i;
240 for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
242 const real refDist = i->refMinDist;
243 EXPECT_NEAR_REL(refDist, search->minimumDistance(i->x), 20*GMX_REAL_EPS);
247 void NeighborhoodSearchTest::testNearestPoint(
248 gmx::AnalysisNeighborhoodSearch *search,
249 const NeighborhoodSearchTestData &data)
251 NeighborhoodSearchTestData::TestPositionList::const_iterator i;
252 for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
254 const gmx::AnalysisNeighborhoodPair pair = search->nearestPoint(i->x);
257 EXPECT_EQ(i->refNearestPoint, pair.refIndex());
258 EXPECT_EQ(0, pair.testIndex());
262 EXPECT_EQ(i->refNearestPoint, -1);
267 void NeighborhoodSearchTest::testPairSearch(
268 gmx::AnalysisNeighborhoodSearch *search,
269 const NeighborhoodSearchTestData &data)
271 NeighborhoodSearchTestData::TestPositionList::const_iterator i;
272 for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
274 std::set<int> checkSet = i->refPairs;
275 gmx::AnalysisNeighborhoodPairSearch pairSearch =
276 search->startPairSearch(i->x);
277 gmx::AnalysisNeighborhoodPair pair;
278 while (pairSearch.findNextPair(&pair))
280 EXPECT_EQ(0, pair.testIndex());
281 if (checkSet.erase(pair.refIndex()) == 0)
283 // TODO: Check whether the same pair was returned more than
284 // once and give a better error message if so.
286 << "Expected: Position " << pair.refIndex()
287 << " is within cutoff.\n"
288 << " Actual: It is not.";
291 EXPECT_TRUE(checkSet.empty()) << "Some positions were not returned by the pair search.";
295 /********************************************************************
296 * Test data generation
299 class RandomBoxFullPBCData
302 static const NeighborhoodSearchTestData &get()
304 static RandomBoxFullPBCData singleton;
305 return singleton.data_;
308 RandomBoxFullPBCData() : data_(12345, 1.0)
310 data_.box_[XX][XX] = 10.0;
311 data_.box_[YY][YY] = 5.0;
312 data_.box_[ZZ][ZZ] = 7.0;
313 // TODO: Consider whether manually picking some positions would give better
315 data_.generateRandomRefPositions(1000);
316 data_.generateRandomTestPositions(100);
317 set_pbc(&data_.pbc_, epbcXYZ, data_.box_);
318 data_.computeReferences(&data_.pbc_);
322 NeighborhoodSearchTestData data_;
325 class RandomTriclinicFullPBCData
328 static const NeighborhoodSearchTestData &get()
330 static RandomTriclinicFullPBCData singleton;
331 return singleton.data_;
334 RandomTriclinicFullPBCData() : data_(12345, 1.0)
336 data_.box_[XX][XX] = 5.0;
337 data_.box_[YY][XX] = 2.5;
338 data_.box_[YY][YY] = 2.5*sqrt(3.0);
339 data_.box_[ZZ][XX] = 2.5;
340 data_.box_[ZZ][YY] = 2.5*sqrt(1.0/3.0);
341 data_.box_[ZZ][ZZ] = 5.0*sqrt(2.0/3.0);
342 // TODO: Consider whether manually picking some positions would give better
344 data_.generateRandomRefPositions(1000);
345 data_.generateRandomTestPositions(100);
346 set_pbc(&data_.pbc_, epbcXYZ, data_.box_);
347 data_.computeReferences(&data_.pbc_);
351 NeighborhoodSearchTestData data_;
354 class RandomBox2DPBCData
357 static const NeighborhoodSearchTestData &get()
359 static RandomBox2DPBCData singleton;
360 return singleton.data_;
363 RandomBox2DPBCData() : data_(12345, 1.0)
365 data_.box_[XX][XX] = 10.0;
366 data_.box_[YY][YY] = 7.0;
367 data_.box_[ZZ][ZZ] = 5.0;
368 // TODO: Consider whether manually picking some positions would give better
370 data_.generateRandomRefPositions(1000);
371 data_.generateRandomTestPositions(100);
372 set_pbc(&data_.pbc_, epbcXY, data_.box_);
373 data_.computeReferences(&data_.pbc_);
377 NeighborhoodSearchTestData data_;
380 /********************************************************************
384 TEST_F(NeighborhoodSearchTest, SimpleSearch)
386 const NeighborhoodSearchTestData &data = RandomBoxFullPBCData::get();
388 nb_.setCutoff(data.cutoff_);
389 nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Simple);
390 gmx::AnalysisNeighborhoodSearch search =
391 nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
392 ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Simple, search.mode());
394 testIsWithin(&search, data);
395 testMinimumDistance(&search, data);
396 testNearestPoint(&search, data);
397 testPairSearch(&search, data);
400 TEST_F(NeighborhoodSearchTest, GridSearchBox)
402 const NeighborhoodSearchTestData &data = RandomBoxFullPBCData::get();
404 nb_.setCutoff(data.cutoff_);
405 nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
406 gmx::AnalysisNeighborhoodSearch search =
407 nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
408 ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
410 testIsWithin(&search, data);
411 testMinimumDistance(&search, data);
412 testNearestPoint(&search, data);
413 testPairSearch(&search, data);
416 TEST_F(NeighborhoodSearchTest, GridSearchTriclinic)
418 const NeighborhoodSearchTestData &data = RandomTriclinicFullPBCData::get();
420 nb_.setCutoff(data.cutoff_);
421 nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
422 gmx::AnalysisNeighborhoodSearch search =
423 nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
424 ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
426 testPairSearch(&search, data);
429 TEST_F(NeighborhoodSearchTest, GridSearch2DPBC)
431 const NeighborhoodSearchTestData &data = RandomBox2DPBCData::get();
433 nb_.setCutoff(data.cutoff_);
434 nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
435 gmx::AnalysisNeighborhoodSearch search =
436 nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
437 // Currently, grid searching not supported with 2D PBC.
438 //ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
440 testIsWithin(&search, data);
441 testMinimumDistance(&search, data);
442 testNearestPoint(&search, data);
443 testPairSearch(&search, data);