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 the index group handling in the selection engine.
40 * Tests for other functions, at least the set operations.
42 * \author Teemu Murtola <teemu.murtola@gmail.com>
43 * \ingroup module_selection
47 #include "gromacs/selection/indexutil.h"
49 #include <gtest/gtest.h>
51 #include "gromacs/topology/block.h"
53 #include "testutils/refdata.h"
60 /********************************************************************
64 class IndexBlockTest : public ::testing::Test
70 void setGroup(int count, const int atoms[]);
72 void setGroup(const int (&atoms)[count])
74 setGroup(count, atoms);
79 gmx::test::TestReferenceData data_;
80 gmx::test::TopologyManager topManager_;
85 IndexBlockTest::IndexBlockTest()
89 blocka_.nalloc_index = 0;
93 gmx_ana_index_clear(&g_);
96 IndexBlockTest::~IndexBlockTest()
98 done_blocka(&blocka_);
101 void IndexBlockTest::setGroup(int count, const int atoms[])
104 g_.index = const_cast<int *>(atoms);
107 void IndexBlockTest::checkBlocka()
109 gmx::test::TestReferenceChecker compound(
110 data_.rootChecker().checkCompound("BlockAtoms", "Block"));
111 compound.checkSequenceArray(g_.isize, g_.index, "Input");
112 compound.checkInteger(blocka_.nr, "Count");
113 for (int i = 0; i < blocka_.nr; ++i)
115 gmx::test::TestReferenceChecker blockCompound(
116 compound.checkCompound("Block", NULL));
117 blockCompound.checkSequence(&blocka_.a[blocka_.index[i]],
118 &blocka_.a[blocka_.index[i+1]],
123 /********************************************************************
124 * gmx_ana_index_make_block() tests
127 TEST_F(IndexBlockTest, CreatesUnknownBlock)
129 gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
131 done_blocka(&blocka_);
132 gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
136 TEST_F(IndexBlockTest, CreatesAtomBlock)
138 const int group[] = { 0, 1, 3, 4, 6 };
140 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, false);
142 done_blocka(&blocka_);
143 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, true);
147 TEST_F(IndexBlockTest, CreatesResidueBlock)
149 const int group[] = { 0, 1, 3, 6, 7 };
150 topManager_.initAtoms(9);
151 topManager_.initUniformResidues(3);
153 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
158 TEST_F(IndexBlockTest, CreatesMoleculeBlock)
160 const int group[] = { 3, 4, 7, 8, 13 };
161 topManager_.initAtoms(18);
162 topManager_.initUniformMolecules(3);
164 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
169 TEST_F(IndexBlockTest, CreatesResidueBlockWithCompletion)
171 const int group[] = { 3, 4, 7, 8, 13 };
172 topManager_.initAtoms(18);
173 topManager_.initUniformResidues(3);
175 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
180 TEST_F(IndexBlockTest, CreatesMoleculeBlockWithCompletion)
182 const int group[] = { 3, 4, 7, 8, 13 };
183 topManager_.initAtoms(18);
184 topManager_.initUniformMolecules(3);
186 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
191 TEST_F(IndexBlockTest, CreatesSingleBlock)
193 const int group[] = { 0, 1, 3, 4, 6 };
195 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, false);
197 done_blocka(&blocka_);
198 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, true);
202 /********************************************************************
203 * gmx_ana_index_has_full_ablocks() tests
206 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
208 const int maxGroup[] = {
209 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
211 const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
212 topManager_.initAtoms(18);
213 topManager_.initUniformResidues(3);
215 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
218 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
221 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
223 const int maxGroup[] = {
224 15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
226 const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
227 topManager_.initAtoms(18);
228 topManager_.initUniformResidues(3);
230 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
233 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
236 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
238 const int maxGroup[] = {
239 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
241 const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
242 const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
243 const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
245 topManager_.initAtoms(18);
246 topManager_.initUniformResidues(3);
248 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
251 setGroup(testGroup1);
252 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
254 setGroup(testGroup2);
255 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
257 setGroup(testGroup3);
258 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
261 /********************************************************************
262 * gmx_ana_index_has_complete_elems() tests
265 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
267 const int group[] = { 0, 1, 2 };
269 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, NULL));
270 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, NULL));
271 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, NULL));
274 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
276 const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
277 const int group2[] = { 3, 4, 5, 6, 7, 8 };
279 topManager_.initAtoms(15);
280 topManager_.initUniformResidues(3);
281 t_topology *top = topManager_.topology();
284 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
287 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
290 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
292 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
293 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
294 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
296 topManager_.initAtoms(18);
297 topManager_.initUniformResidues(3);
298 t_topology *top = topManager_.topology();
301 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
304 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
307 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
310 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
312 const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
314 topManager_.initAtoms(15);
315 topManager_.initUniformMolecules(3);
316 t_topology *top = topManager_.topology();
319 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
322 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
324 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
325 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
326 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
328 topManager_.initAtoms(18);
329 topManager_.initUniformMolecules(3);
330 t_topology *top = topManager_.topology();
333 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
336 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
339 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
342 /********************************************************************
346 class IndexMapTest : public ::testing::Test
352 void testInit(int atomCount, const int atoms[], e_index_t type);
353 void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
356 void testInit(const int (&atoms)[count], e_index_t type)
358 testInit(count, atoms, type);
361 void testUpdate(const int (&atoms)[count], bool bMaskOnly,
364 testUpdate(count, atoms, bMaskOnly, name);
367 void checkMapping(int atomCount, const int atoms[], const char *name);
369 gmx::test::TestReferenceData data_;
370 gmx::test::TestReferenceChecker checker_;
371 gmx::test::TopologyManager topManager_;
372 gmx_ana_indexmap_t map_;
375 gmx_ana_index_t initGroup_;
378 IndexMapTest::IndexMapTest()
379 : checker_(data_.rootChecker())
381 gmx_ana_indexmap_clear(&map_);
382 gmx_ana_index_clear(&initGroup_);
385 IndexMapTest::~IndexMapTest()
387 gmx_ana_indexmap_deinit(&map_);
390 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
392 initGroup_.isize = atomCount;
393 initGroup_.index = const_cast<int *>(atoms);
394 gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
395 EXPECT_EQ(type, map_.type);
396 checkMapping(atomCount, atoms, "Initialized");
399 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
404 g.index = const_cast<int *>(atoms);
405 gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
412 checkMapping(initGroup_.isize, initGroup_.index, name);
416 checkMapping(atomCount, atoms, name);
420 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
423 gmx::test::TestReferenceChecker compound(
424 checker_.checkCompound("IndexMapping", name));
425 compound.checkSequenceArray(atomCount, atoms, "Input");
426 compound.checkInteger(map_.mapb.nr, "Count");
427 for (int i = 0; i < map_.mapb.nr; ++i)
429 gmx::test::TestReferenceChecker blockCompound(
430 compound.checkCompound("Block", NULL));
431 blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
432 &atoms[map_.mapb.index[i+1]],
434 blockCompound.checkInteger(map_.refid[i], "RefId");
435 blockCompound.checkInteger(map_.mapid[i], "MapId");
436 int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
437 EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
441 /********************************************************************
442 * gmx_ana_indexmap_t tests
445 TEST_F(IndexMapTest, InitializesAtomBlock)
447 const int maxGroup[] = { 1, 2, 4, 5 };
448 testInit(maxGroup, INDEX_ATOM);
451 TEST_F(IndexMapTest, InitializesMoleculeBlock)
453 const int maxGroup[] = { 3, 4, 7, 8, 13 };
454 topManager_.initAtoms(18);
455 topManager_.initUniformMolecules(3);
456 testInit(maxGroup, INDEX_MOL);
459 TEST_F(IndexMapTest, MapsSingleBlock)
461 const int maxGroup[] = { 0, 1, 2, 3 };
462 const int evalGroup[] = { 0, 2 };
463 testInit(maxGroup, INDEX_ALL);
464 testUpdate(evalGroup, false, NULL);
467 TEST_F(IndexMapTest, MapsResidueBlocks)
469 const int maxGroup[] = {
470 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
472 const int evalGroup[] = { 3, 4, 7, 8, 13 };
473 topManager_.initAtoms(18);
474 topManager_.initUniformResidues(3);
475 testInit(maxGroup, INDEX_RES);
476 testUpdate(evalGroup, false, NULL);
479 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
481 const int maxGroup[] = {
482 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
484 const int evalGroup[] = { 3, 4, 7, 8, 13 };
485 topManager_.initAtoms(18);
486 topManager_.initUniformResidues(3);
487 testInit(maxGroup, INDEX_RES);
488 testUpdate(evalGroup, true, NULL);
491 TEST_F(IndexMapTest, HandlesMultipleRequests)
493 const int maxGroup[] = {
494 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
496 const int evalGroup[] = { 3, 4, 7, 8, 13 };
497 topManager_.initAtoms(18);
498 topManager_.initUniformResidues(3);
499 testInit(maxGroup, INDEX_RES);
500 testUpdate(evalGroup, false, "EvaluatedNoMask");
501 testUpdate(evalGroup, true, "EvaluatedMask");
502 testUpdate(maxGroup, true, "Initialized");
503 testUpdate(evalGroup, true, "EvaluatedMask");
504 testUpdate(evalGroup, false, "EvaluatedNoMask");
505 testUpdate(maxGroup, false, "Initialized");