Better concurrency support for analysis nbsearch.
[alexxy/gromacs.git] / src / gromacs / selection / tests / nbsearch.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Tests selection neighborhood searching.
38  *
39  * \todo
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.
43  *
44  * \author Teemu Murtola <teemu.murtola@gmail.com>
45  * \ingroup module_selection
46  */
47 #include <gtest/gtest.h>
48
49 #include <cmath>
50
51 #include <limits>
52 #include <set>
53
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"
58
59 #include "gromacs/selection/nbsearch.h"
60
61 #include "testutils/testasserts.h"
62
63 namespace
64 {
65
66 /********************************************************************
67  * NeighborhoodSearchTestData
68  */
69
70 class NeighborhoodSearchTestData
71 {
72     public:
73         struct TestPosition
74         {
75             TestPosition() : refMinDist(0.0), refNearestPoint(-1)
76             {
77                 clear_rvec(x);
78             }
79             explicit TestPosition(const rvec x)
80                 : refMinDist(0.0), refNearestPoint(-1)
81             {
82                 copy_rvec(x, this->x);
83             }
84
85             rvec                x;
86             real                refMinDist;
87             int                 refNearestPoint;
88             std::set<int>       refPairs;
89         };
90         typedef std::vector<TestPosition> TestPositionList;
91
92         NeighborhoodSearchTestData(int seed, real cutoff);
93         ~NeighborhoodSearchTestData();
94
95         void addTestPosition(const rvec x)
96         {
97             testPositions_.push_back(TestPosition(x));
98         }
99         void generateRandomPosition(rvec x);
100         void generateRandomRefPositions(int count);
101         void generateRandomTestPositions(int count);
102         void computeReferences(t_pbc *pbc);
103
104         gmx_rng_t                        rng_;
105         real                             cutoff_;
106         matrix                           box_;
107         t_pbc                            pbc_;
108         int                              refPosCount_;
109         rvec                            *refPos_;
110         TestPositionList                 testPositions_;
111 };
112
113 NeighborhoodSearchTestData::NeighborhoodSearchTestData(int seed, real cutoff)
114     : rng_(NULL), cutoff_(cutoff), refPosCount_(0), refPos_(NULL)
115 {
116     // TODO: Handle errors.
117     rng_ = gmx_rng_init(seed);
118     clear_mat(box_);
119     set_pbc(&pbc_, epbcNONE, box_);
120 }
121
122 NeighborhoodSearchTestData::~NeighborhoodSearchTestData()
123 {
124     if (rng_ != NULL)
125     {
126         gmx_rng_destroy(rng_);
127     }
128     sfree(refPos_);
129 }
130
131 void NeighborhoodSearchTestData::generateRandomPosition(rvec x)
132 {
133     rvec fx;
134     fx[XX] = gmx_rng_uniform_real(rng_);
135     fx[YY] = gmx_rng_uniform_real(rng_);
136     fx[ZZ] = gmx_rng_uniform_real(rng_);
137     mvmul(box_, fx, x);
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;
142 }
143
144 void NeighborhoodSearchTestData::generateRandomRefPositions(int count)
145 {
146     refPosCount_ = count;
147     snew(refPos_, refPosCount_);
148     for (int i = 0; i < refPosCount_; ++i)
149     {
150         generateRandomPosition(refPos_[i]);
151     }
152 }
153
154 void NeighborhoodSearchTestData::generateRandomTestPositions(int count)
155 {
156     testPositions_.reserve(count);
157     for (int i = 0; i < count; ++i)
158     {
159         rvec x;
160         generateRandomPosition(x);
161         addTestPosition(x);
162     }
163 }
164
165 void NeighborhoodSearchTestData::computeReferences(t_pbc *pbc)
166 {
167     real cutoff = cutoff_;
168     if (cutoff <= 0)
169     {
170         cutoff = std::numeric_limits<real>::max();
171     }
172     TestPositionList::iterator i;
173     for (i = testPositions_.begin(); i != testPositions_.end(); ++i)
174     {
175         i->refMinDist      = cutoff;
176         i->refNearestPoint = -1;
177         i->refPairs.clear();
178         for (int j = 0; j < refPosCount_; ++j)
179         {
180             rvec dx;
181             if (pbc != NULL)
182             {
183                 pbc_dx(pbc, i->x, refPos_[j], dx);
184             }
185             else
186             {
187                 rvec_sub(i->x, refPos_[j], dx);
188             }
189             const real dist = norm(dx);
190             if (dist < i->refMinDist)
191             {
192                 i->refMinDist      = dist;
193                 i->refNearestPoint = j;
194             }
195             if (dist <= cutoff)
196             {
197                 i->refPairs.insert(j);
198             }
199         }
200     }
201 }
202
203 /********************************************************************
204  * NeighborhoodSearchTest
205  */
206
207 class NeighborhoodSearchTest : public ::testing::Test
208 {
209     public:
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);
218
219         gmx::AnalysisNeighborhood        nb_;
220 };
221
222 void NeighborhoodSearchTest::testIsWithin(
223         gmx::AnalysisNeighborhoodSearch  *search,
224         const NeighborhoodSearchTestData &data)
225 {
226     NeighborhoodSearchTestData::TestPositionList::const_iterator i;
227     for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
228     {
229         const bool bWithin = (i->refMinDist <= data.cutoff_);
230         EXPECT_EQ(bWithin, search->isWithin(i->x))
231         << "Distance is " << i->refMinDist;
232     }
233 }
234
235 void NeighborhoodSearchTest::testMinimumDistance(
236         gmx::AnalysisNeighborhoodSearch  *search,
237         const NeighborhoodSearchTestData &data)
238 {
239     NeighborhoodSearchTestData::TestPositionList::const_iterator i;
240     for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
241     {
242         const real refDist = i->refMinDist;
243         EXPECT_NEAR_REL(refDist, search->minimumDistance(i->x), 20*GMX_REAL_EPS);
244     }
245 }
246
247 void NeighborhoodSearchTest::testNearestPoint(
248         gmx::AnalysisNeighborhoodSearch  *search,
249         const NeighborhoodSearchTestData &data)
250 {
251     NeighborhoodSearchTestData::TestPositionList::const_iterator i;
252     for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
253     {
254         const gmx::AnalysisNeighborhoodPair pair = search->nearestPoint(i->x);
255         if (pair.isValid())
256         {
257             EXPECT_EQ(i->refNearestPoint, pair.refIndex());
258             EXPECT_EQ(0, pair.testIndex());
259         }
260         else
261         {
262             EXPECT_EQ(i->refNearestPoint, -1);
263         }
264     }
265 }
266
267 void NeighborhoodSearchTest::testPairSearch(
268         gmx::AnalysisNeighborhoodSearch  *search,
269         const NeighborhoodSearchTestData &data)
270 {
271     NeighborhoodSearchTestData::TestPositionList::const_iterator i;
272     for (i = data.testPositions_.begin(); i != data.testPositions_.end(); ++i)
273     {
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))
279         {
280             EXPECT_EQ(0, pair.testIndex());
281             if (checkSet.erase(pair.refIndex()) == 0)
282             {
283                 // TODO: Check whether the same pair was returned more than
284                 // once and give a better error message if so.
285                 ADD_FAILURE()
286                 << "Expected: Position " << pair.refIndex()
287                 << " is within cutoff.\n"
288                 << "  Actual: It is not.";
289             }
290         }
291         EXPECT_TRUE(checkSet.empty()) << "Some positions were not returned by the pair search.";
292     }
293 }
294
295 /********************************************************************
296  * Test data generation
297  */
298
299 class TrivialTestData
300 {
301     public:
302         static const NeighborhoodSearchTestData &get()
303         {
304             static TrivialTestData singleton;
305             return singleton.data_;
306         }
307
308         TrivialTestData() : data_(12345, 1.0)
309         {
310             data_.box_[XX][XX] = 5.0;
311             data_.box_[YY][YY] = 5.0;
312             data_.box_[ZZ][ZZ] = 5.0;
313             data_.generateRandomRefPositions(10);
314             data_.generateRandomTestPositions(5);
315             set_pbc(&data_.pbc_, epbcXYZ, data_.box_);
316             data_.computeReferences(&data_.pbc_);
317         }
318
319     private:
320         NeighborhoodSearchTestData data_;
321 };
322
323 class RandomBoxFullPBCData
324 {
325     public:
326         static const NeighborhoodSearchTestData &get()
327         {
328             static RandomBoxFullPBCData singleton;
329             return singleton.data_;
330         }
331
332         RandomBoxFullPBCData() : data_(12345, 1.0)
333         {
334             data_.box_[XX][XX] = 10.0;
335             data_.box_[YY][YY] = 5.0;
336             data_.box_[ZZ][ZZ] = 7.0;
337             // TODO: Consider whether manually picking some positions would give better
338             // test coverage.
339             data_.generateRandomRefPositions(1000);
340             data_.generateRandomTestPositions(100);
341             set_pbc(&data_.pbc_, epbcXYZ, data_.box_);
342             data_.computeReferences(&data_.pbc_);
343         }
344
345     private:
346         NeighborhoodSearchTestData data_;
347 };
348
349 class RandomTriclinicFullPBCData
350 {
351     public:
352         static const NeighborhoodSearchTestData &get()
353         {
354             static RandomTriclinicFullPBCData singleton;
355             return singleton.data_;
356         }
357
358         RandomTriclinicFullPBCData() : data_(12345, 1.0)
359         {
360             data_.box_[XX][XX] = 5.0;
361             data_.box_[YY][XX] = 2.5;
362             data_.box_[YY][YY] = 2.5*sqrt(3.0);
363             data_.box_[ZZ][XX] = 2.5;
364             data_.box_[ZZ][YY] = 2.5*sqrt(1.0/3.0);
365             data_.box_[ZZ][ZZ] = 5.0*sqrt(2.0/3.0);
366             // TODO: Consider whether manually picking some positions would give better
367             // test coverage.
368             data_.generateRandomRefPositions(1000);
369             data_.generateRandomTestPositions(100);
370             set_pbc(&data_.pbc_, epbcXYZ, data_.box_);
371             data_.computeReferences(&data_.pbc_);
372         }
373
374     private:
375         NeighborhoodSearchTestData data_;
376 };
377
378 class RandomBox2DPBCData
379 {
380     public:
381         static const NeighborhoodSearchTestData &get()
382         {
383             static RandomBox2DPBCData singleton;
384             return singleton.data_;
385         }
386
387         RandomBox2DPBCData() : data_(12345, 1.0)
388         {
389             data_.box_[XX][XX] = 10.0;
390             data_.box_[YY][YY] = 7.0;
391             data_.box_[ZZ][ZZ] = 5.0;
392             // TODO: Consider whether manually picking some positions would give better
393             // test coverage.
394             data_.generateRandomRefPositions(1000);
395             data_.generateRandomTestPositions(100);
396             set_pbc(&data_.pbc_, epbcXY, data_.box_);
397             data_.computeReferences(&data_.pbc_);
398         }
399
400     private:
401         NeighborhoodSearchTestData data_;
402 };
403
404 /********************************************************************
405  * Actual tests
406  */
407
408 TEST_F(NeighborhoodSearchTest, SimpleSearch)
409 {
410     const NeighborhoodSearchTestData &data = RandomBoxFullPBCData::get();
411
412     nb_.setCutoff(data.cutoff_);
413     nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Simple);
414     gmx::AnalysisNeighborhoodSearch search =
415         nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
416     ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Simple, search.mode());
417
418     testIsWithin(&search, data);
419     testMinimumDistance(&search, data);
420     testNearestPoint(&search, data);
421     testPairSearch(&search, data);
422 }
423
424 TEST_F(NeighborhoodSearchTest, GridSearchBox)
425 {
426     const NeighborhoodSearchTestData &data = RandomBoxFullPBCData::get();
427
428     nb_.setCutoff(data.cutoff_);
429     nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
430     gmx::AnalysisNeighborhoodSearch search =
431         nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
432     ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
433
434     testIsWithin(&search, data);
435     testMinimumDistance(&search, data);
436     testNearestPoint(&search, data);
437     testPairSearch(&search, data);
438 }
439
440 TEST_F(NeighborhoodSearchTest, GridSearchTriclinic)
441 {
442     const NeighborhoodSearchTestData &data = RandomTriclinicFullPBCData::get();
443
444     nb_.setCutoff(data.cutoff_);
445     nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
446     gmx::AnalysisNeighborhoodSearch search =
447         nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
448     ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
449
450     testPairSearch(&search, data);
451 }
452
453 TEST_F(NeighborhoodSearchTest, GridSearch2DPBC)
454 {
455     const NeighborhoodSearchTestData &data = RandomBox2DPBCData::get();
456
457     nb_.setCutoff(data.cutoff_);
458     nb_.setMode(gmx::AnalysisNeighborhood::eSearchMode_Grid);
459     gmx::AnalysisNeighborhoodSearch search =
460         nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
461     // Currently, grid searching not supported with 2D PBC.
462     //ASSERT_EQ(gmx::AnalysisNeighborhood::eSearchMode_Grid, search.mode());
463
464     testIsWithin(&search, data);
465     testMinimumDistance(&search, data);
466     testNearestPoint(&search, data);
467     testPairSearch(&search, data);
468 }
469
470 TEST_F(NeighborhoodSearchTest, HandlesConcurrentSearches)
471 {
472     const NeighborhoodSearchTestData &data = TrivialTestData::get();
473
474     nb_.setCutoff(data.cutoff_);
475     gmx::AnalysisNeighborhoodSearch search1 =
476         nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
477     gmx::AnalysisNeighborhoodSearch search2 =
478         nb_.initSearch(&data.pbc_, data.refPosCount_, data.refPos_);
479
480     testPairSearch(&search1, data);
481     testPairSearch(&search2, data);
482 }
483
484 } // namespace