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 <gtest/gtest.h>
49 #include "gromacs/selection/indexutil.h"
50 #include "gromacs/topology/block.h"
52 #include "testutils/refdata.h"
59 /********************************************************************
63 class IndexBlockTest : public ::testing::Test
69 void setGroup(int count, const int atoms[]);
71 void setGroup(const int (&atoms)[count])
73 setGroup(count, atoms);
78 gmx::test::TestReferenceData data_;
79 gmx::test::TopologyManager topManager_;
84 IndexBlockTest::IndexBlockTest()
88 blocka_.nalloc_index = 0;
92 gmx_ana_index_clear(&g_);
95 IndexBlockTest::~IndexBlockTest()
97 done_blocka(&blocka_);
100 void IndexBlockTest::setGroup(int count, const int atoms[])
103 g_.index = const_cast<int *>(atoms);
106 void IndexBlockTest::checkBlocka()
108 gmx::test::TestReferenceChecker compound(
109 data_.rootChecker().checkCompound("BlockAtoms", "Block"));
110 compound.checkSequenceArray(g_.isize, g_.index, "Input");
111 compound.checkInteger(blocka_.nr, "Count");
112 for (int i = 0; i < blocka_.nr; ++i)
114 gmx::test::TestReferenceChecker blockCompound(
115 compound.checkCompound("Block", NULL));
116 blockCompound.checkSequence(&blocka_.a[blocka_.index[i]],
117 &blocka_.a[blocka_.index[i+1]],
122 /********************************************************************
123 * gmx_ana_index_make_block() tests
126 TEST_F(IndexBlockTest, CreatesUnknownBlock)
128 gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
130 done_blocka(&blocka_);
131 gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
135 TEST_F(IndexBlockTest, CreatesAtomBlock)
137 const int group[] = { 0, 1, 3, 4, 6 };
139 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, false);
141 done_blocka(&blocka_);
142 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, true);
146 TEST_F(IndexBlockTest, CreatesResidueBlock)
148 const int group[] = { 0, 1, 3, 6, 7 };
149 topManager_.initAtoms(9);
150 topManager_.initUniformResidues(3);
152 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
157 TEST_F(IndexBlockTest, CreatesMoleculeBlock)
159 const int group[] = { 3, 4, 7, 8, 13 };
160 topManager_.initAtoms(18);
161 topManager_.initUniformMolecules(3);
163 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
168 TEST_F(IndexBlockTest, CreatesResidueBlockWithCompletion)
170 const int group[] = { 3, 4, 7, 8, 13 };
171 topManager_.initAtoms(18);
172 topManager_.initUniformResidues(3);
174 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
179 TEST_F(IndexBlockTest, CreatesMoleculeBlockWithCompletion)
181 const int group[] = { 3, 4, 7, 8, 13 };
182 topManager_.initAtoms(18);
183 topManager_.initUniformMolecules(3);
185 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
190 TEST_F(IndexBlockTest, CreatesSingleBlock)
192 const int group[] = { 0, 1, 3, 4, 6 };
194 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, false);
196 done_blocka(&blocka_);
197 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, true);
201 /********************************************************************
202 * gmx_ana_index_has_full_ablocks() tests
205 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
207 const int maxGroup[] = {
208 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
210 const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
211 topManager_.initAtoms(18);
212 topManager_.initUniformResidues(3);
214 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
217 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
220 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
222 const int maxGroup[] = {
223 15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
225 const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
226 topManager_.initAtoms(18);
227 topManager_.initUniformResidues(3);
229 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
232 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
235 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
237 const int maxGroup[] = {
238 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
240 const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
241 const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
242 const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
244 topManager_.initAtoms(18);
245 topManager_.initUniformResidues(3);
247 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
250 setGroup(testGroup1);
251 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
253 setGroup(testGroup2);
254 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
256 setGroup(testGroup3);
257 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
260 /********************************************************************
261 * gmx_ana_index_has_complete_elems() tests
264 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
266 const int group[] = { 0, 1, 2 };
268 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, NULL));
269 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, NULL));
270 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, NULL));
273 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
275 const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
276 const int group2[] = { 3, 4, 5, 6, 7, 8 };
278 topManager_.initAtoms(15);
279 topManager_.initUniformResidues(3);
280 t_topology *top = topManager_.topology();
283 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
286 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
289 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
291 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
292 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
293 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
295 topManager_.initAtoms(18);
296 topManager_.initUniformResidues(3);
297 t_topology *top = topManager_.topology();
300 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
303 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
306 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
309 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
311 const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
313 topManager_.initAtoms(15);
314 topManager_.initUniformMolecules(3);
315 t_topology *top = topManager_.topology();
318 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
321 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
323 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
324 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
325 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
327 topManager_.initAtoms(18);
328 topManager_.initUniformMolecules(3);
329 t_topology *top = topManager_.topology();
332 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
335 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
338 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
341 /********************************************************************
345 class IndexMapTest : public ::testing::Test
351 void testInit(int atomCount, const int atoms[], e_index_t type);
352 void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
355 void testInit(const int (&atoms)[count], e_index_t type)
357 testInit(count, atoms, type);
360 void testUpdate(const int (&atoms)[count], bool bMaskOnly,
363 testUpdate(count, atoms, bMaskOnly, name);
366 void checkMapping(int atomCount, const int atoms[], const char *name);
368 gmx::test::TestReferenceData data_;
369 gmx::test::TestReferenceChecker checker_;
370 gmx::test::TopologyManager topManager_;
371 gmx_ana_indexmap_t map_;
374 gmx_ana_index_t initGroup_;
377 IndexMapTest::IndexMapTest()
378 : checker_(data_.rootChecker())
380 gmx_ana_indexmap_clear(&map_);
381 gmx_ana_index_clear(&initGroup_);
384 IndexMapTest::~IndexMapTest()
386 gmx_ana_indexmap_deinit(&map_);
389 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
391 initGroup_.isize = atomCount;
392 initGroup_.index = const_cast<int *>(atoms);
393 gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
394 EXPECT_EQ(type, map_.type);
395 checkMapping(atomCount, atoms, "Initialized");
398 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
403 g.index = const_cast<int *>(atoms);
404 gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
411 checkMapping(initGroup_.isize, initGroup_.index, name);
415 checkMapping(atomCount, atoms, name);
419 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
422 gmx::test::TestReferenceChecker compound(
423 checker_.checkCompound("IndexMapping", name));
424 compound.checkSequenceArray(atomCount, atoms, "Input");
425 compound.checkInteger(map_.mapb.nr, "Count");
426 for (int i = 0; i < map_.mapb.nr; ++i)
428 gmx::test::TestReferenceChecker blockCompound(
429 compound.checkCompound("Block", NULL));
430 blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
431 &atoms[map_.mapb.index[i+1]],
433 blockCompound.checkInteger(map_.refid[i], "RefId");
434 blockCompound.checkInteger(map_.mapid[i], "MapId");
435 int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
436 EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
440 /********************************************************************
441 * gmx_ana_indexmap_t tests
444 TEST_F(IndexMapTest, InitializesAtomBlock)
446 const int maxGroup[] = { 1, 2, 4, 5 };
447 testInit(maxGroup, INDEX_ATOM);
450 TEST_F(IndexMapTest, InitializesMoleculeBlock)
452 const int maxGroup[] = { 3, 4, 7, 8, 13 };
453 topManager_.initAtoms(18);
454 topManager_.initUniformMolecules(3);
455 testInit(maxGroup, INDEX_MOL);
458 TEST_F(IndexMapTest, MapsSingleBlock)
460 const int maxGroup[] = { 0, 1, 2, 3 };
461 const int evalGroup[] = { 0, 2 };
462 testInit(maxGroup, INDEX_ALL);
463 testUpdate(evalGroup, false, NULL);
466 TEST_F(IndexMapTest, MapsResidueBlocks)
468 const int maxGroup[] = {
469 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
471 const int evalGroup[] = { 3, 4, 7, 8, 13 };
472 topManager_.initAtoms(18);
473 topManager_.initUniformResidues(3);
474 testInit(maxGroup, INDEX_RES);
475 testUpdate(evalGroup, false, NULL);
478 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
480 const int maxGroup[] = {
481 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
483 const int evalGroup[] = { 3, 4, 7, 8, 13 };
484 topManager_.initAtoms(18);
485 topManager_.initUniformResidues(3);
486 testInit(maxGroup, INDEX_RES);
487 testUpdate(evalGroup, true, NULL);
490 TEST_F(IndexMapTest, HandlesMultipleRequests)
492 const int maxGroup[] = {
493 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
495 const int evalGroup[] = { 3, 4, 7, 8, 13 };
496 topManager_.initAtoms(18);
497 topManager_.initUniformResidues(3);
498 testInit(maxGroup, INDEX_RES);
499 testUpdate(evalGroup, false, "EvaluatedNoMask");
500 testUpdate(evalGroup, true, "EvaluatedMask");
501 testUpdate(maxGroup, true, "Initialized");
502 testUpdate(evalGroup, true, "EvaluatedMask");
503 testUpdate(evalGroup, false, "EvaluatedNoMask");
504 testUpdate(maxGroup, false, "Initialized");