2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016, 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"
52 #include "gromacs/utility/arrayref.h"
54 #include "testutils/refdata.h"
61 //! Helper for creating groups from an array.
62 gmx_ana_index_t initGroup(gmx::ArrayRef<int> index)
64 gmx_ana_index_t g = { static_cast<int>(index.size()), index.data(), 0 };
68 TEST(IndexGroupTest, RemovesDuplicates)
70 int index[] = { 1, 1, 2, 3, 4, 4 };
71 int expected[] = { 1, 2, 3, 4 };
72 gmx_ana_index_t g = initGroup(index);
73 gmx_ana_index_t e = initGroup(expected);
74 gmx_ana_index_remove_duplicates(&g);
75 EXPECT_TRUE(gmx_ana_index_equals(&g, &e));
78 /********************************************************************
82 class IndexBlockTest : public ::testing::Test
88 void setGroup(int count, const int atoms[]);
90 void setGroup(const int (&atoms)[count])
92 setGroup(count, atoms);
97 gmx::test::TestReferenceData data_;
98 gmx::test::TopologyManager topManager_;
103 IndexBlockTest::IndexBlockTest()
106 blocka_.index = NULL;
107 blocka_.nalloc_index = 0;
110 blocka_.nalloc_a = 0;
111 gmx_ana_index_clear(&g_);
114 IndexBlockTest::~IndexBlockTest()
116 done_blocka(&blocka_);
119 void IndexBlockTest::setGroup(int count, const int atoms[])
122 g_.index = const_cast<int *>(atoms);
125 void IndexBlockTest::checkBlocka()
127 gmx::test::TestReferenceChecker compound(
128 data_.rootChecker().checkCompound("BlockAtoms", "Block"));
129 compound.checkSequenceArray(g_.isize, g_.index, "Input");
130 compound.checkInteger(blocka_.nr, "Count");
131 for (int i = 0; i < blocka_.nr; ++i)
133 gmx::test::TestReferenceChecker blockCompound(
134 compound.checkCompound("Block", NULL));
135 blockCompound.checkSequence(&blocka_.a[blocka_.index[i]],
136 &blocka_.a[blocka_.index[i+1]],
141 /********************************************************************
142 * gmx_ana_index_make_block() tests
145 TEST_F(IndexBlockTest, CreatesUnknownBlock)
147 gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
149 done_blocka(&blocka_);
150 gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
154 TEST_F(IndexBlockTest, CreatesAtomBlock)
156 const int group[] = { 0, 1, 3, 4, 6 };
158 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, false);
160 done_blocka(&blocka_);
161 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, true);
165 TEST_F(IndexBlockTest, CreatesResidueBlock)
167 const int group[] = { 0, 1, 3, 6, 7 };
168 topManager_.initAtoms(9);
169 topManager_.initUniformResidues(3);
171 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
176 TEST_F(IndexBlockTest, CreatesMoleculeBlock)
178 const int group[] = { 3, 4, 7, 8, 13 };
179 topManager_.initAtoms(18);
180 topManager_.initUniformMolecules(3);
182 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
187 TEST_F(IndexBlockTest, CreatesResidueBlockWithCompletion)
189 const int group[] = { 3, 4, 7, 8, 13 };
190 topManager_.initAtoms(18);
191 topManager_.initUniformResidues(3);
193 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
198 TEST_F(IndexBlockTest, CreatesMoleculeBlockWithCompletion)
200 const int group[] = { 3, 4, 7, 8, 13 };
201 topManager_.initAtoms(18);
202 topManager_.initUniformMolecules(3);
204 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
209 TEST_F(IndexBlockTest, CreatesSingleBlock)
211 const int group[] = { 0, 1, 3, 4, 6 };
213 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, false);
215 done_blocka(&blocka_);
216 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, true);
220 /********************************************************************
221 * gmx_ana_index_has_full_ablocks() tests
224 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
226 const int maxGroup[] = {
227 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
229 const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
230 topManager_.initAtoms(18);
231 topManager_.initUniformResidues(3);
233 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
236 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
239 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
241 const int maxGroup[] = {
242 15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
244 const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
245 topManager_.initAtoms(18);
246 topManager_.initUniformResidues(3);
248 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
251 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
254 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
256 const int maxGroup[] = {
257 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
259 const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
260 const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
261 const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
263 topManager_.initAtoms(18);
264 topManager_.initUniformResidues(3);
266 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
269 setGroup(testGroup1);
270 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
272 setGroup(testGroup2);
273 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
275 setGroup(testGroup3);
276 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
279 /********************************************************************
280 * gmx_ana_index_has_complete_elems() tests
283 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
285 const int group[] = { 0, 1, 2 };
287 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, NULL));
288 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, NULL));
289 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, NULL));
292 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
294 const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
295 const int group2[] = { 3, 4, 5, 6, 7, 8 };
297 topManager_.initAtoms(15);
298 topManager_.initUniformResidues(3);
299 gmx_mtop_t *top = topManager_.topology();
302 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
305 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
308 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
310 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
311 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
312 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
313 const int group4[] = { 3, 4, 5, 6, 8, 12, 13, 14 };
315 topManager_.initAtoms(18);
316 topManager_.initUniformResidues(3);
317 gmx_mtop_t *top = topManager_.topology();
320 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
323 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
326 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
329 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
332 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
334 const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
336 topManager_.initAtoms(15);
337 topManager_.initUniformMolecules(3);
338 gmx_mtop_t *top = topManager_.topology();
341 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
344 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
346 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
347 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
348 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
350 topManager_.initAtoms(18);
351 topManager_.initUniformMolecules(3);
352 gmx_mtop_t *top = topManager_.topology();
355 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
358 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
361 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
364 /********************************************************************
368 class IndexMapTest : public ::testing::Test
374 void testInit(int atomCount, const int atoms[], e_index_t type);
375 void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
377 void testOrgIdGroup(e_index_t type, const char *name);
379 void testInit(const int (&atoms)[count], e_index_t type)
381 testInit(count, atoms, type);
384 void testUpdate(const int (&atoms)[count], bool bMaskOnly,
387 testUpdate(count, atoms, bMaskOnly, name);
390 void checkMapping(int atomCount, const int atoms[], const char *name);
392 gmx::test::TestReferenceData data_;
393 gmx::test::TestReferenceChecker checker_;
394 gmx::test::TopologyManager topManager_;
395 gmx_ana_indexmap_t map_;
398 gmx_ana_index_t initGroup_;
401 IndexMapTest::IndexMapTest()
402 : checker_(data_.rootChecker())
404 gmx_ana_indexmap_clear(&map_);
405 gmx_ana_index_clear(&initGroup_);
408 IndexMapTest::~IndexMapTest()
410 gmx_ana_indexmap_deinit(&map_);
413 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
415 initGroup_.isize = atomCount;
416 initGroup_.index = const_cast<int *>(atoms);
417 gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
418 EXPECT_EQ(type, map_.type);
419 checkMapping(atomCount, atoms, "Initialized");
422 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
427 g.index = const_cast<int *>(atoms);
428 gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
435 checkMapping(initGroup_.isize, initGroup_.index, name);
439 checkMapping(atomCount, atoms, name);
443 void IndexMapTest::testOrgIdGroup(e_index_t type, const char *name)
445 gmx::test::TestReferenceChecker compound(
446 checker_.checkCompound("OrgIdGroups", name));
448 = gmx_ana_indexmap_init_orgid_group(&map_, topManager_.topology(), type);
449 compound.checkInteger(count, "GroupCount");
450 compound.checkSequenceArray(map_.mapb.nr, map_.orgid, "OrgId");
451 for (int i = 0; i < map_.mapb.nr; ++i)
453 EXPECT_EQ(map_.orgid[i], map_.mapid[i]);
457 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
460 gmx::test::TestReferenceChecker compound(
461 checker_.checkCompound("IndexMapping", name));
462 compound.checkSequenceArray(atomCount, atoms, "Input");
463 compound.checkInteger(map_.mapb.nr, "Count");
464 for (int i = 0; i < map_.mapb.nr; ++i)
466 gmx::test::TestReferenceChecker blockCompound(
467 compound.checkCompound("Block", NULL));
468 blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
469 &atoms[map_.mapb.index[i+1]],
471 blockCompound.checkInteger(map_.refid[i], "RefId");
472 blockCompound.checkInteger(map_.mapid[i], "MapId");
473 int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
474 EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
478 /********************************************************************
479 * gmx_ana_indexmap_t tests
482 TEST_F(IndexMapTest, InitializesAtomBlock)
484 const int maxGroup[] = { 1, 2, 4, 5 };
485 testInit(maxGroup, INDEX_ATOM);
488 TEST_F(IndexMapTest, InitializesOrgIdGroupAtom)
490 const int maxGroup[] = { 2, 5, 7 };
491 testInit(maxGroup, INDEX_ATOM);
492 testOrgIdGroup(INDEX_ATOM, "Atoms");
495 TEST_F(IndexMapTest, InitializesOrgIdGroupSingle)
497 const int maxGroup[] = { 3, 4, 7, 8, 13 };
498 topManager_.initAtoms(18);
499 topManager_.initUniformResidues(3);
500 testInit(maxGroup, INDEX_RES);
501 testOrgIdGroup(INDEX_ATOM, "Single");
504 TEST_F(IndexMapTest, InitializesOrgIdGroupResidue)
506 const int maxGroup[] = { 3, 4, 7, 8, 13 };
507 topManager_.initAtoms(18);
508 topManager_.initUniformResidues(3);
509 testInit(maxGroup, INDEX_ATOM);
510 testOrgIdGroup(INDEX_RES, "Residues");
513 TEST_F(IndexMapTest, InitializesOrgIdGroupMolecule)
515 const int maxGroup[] = { 1, 2, 3, 4, 7, 8, 13 };
516 topManager_.initAtoms(18);
517 topManager_.initUniformResidues(3);
518 topManager_.initUniformMolecules(6);
519 testInit(maxGroup, INDEX_RES);
520 testOrgIdGroup(INDEX_MOL, "Molecules");
523 TEST_F(IndexMapTest, InitializesOrgIdGroupAll)
525 const int maxGroup[] = { 3, 4, 7, 8, 13 };
526 testInit(maxGroup, INDEX_ATOM);
527 testOrgIdGroup(INDEX_ALL, "All");
530 TEST_F(IndexMapTest, InitializesMoleculeBlock)
532 const int maxGroup[] = { 3, 4, 7, 8, 13 };
533 topManager_.initAtoms(18);
534 topManager_.initUniformMolecules(3);
535 testInit(maxGroup, INDEX_MOL);
538 TEST_F(IndexMapTest, MapsSingleBlock)
540 const int maxGroup[] = { 0, 1, 2, 3 };
541 const int evalGroup[] = { 0, 2 };
542 testInit(maxGroup, INDEX_ALL);
543 testUpdate(evalGroup, false, NULL);
546 TEST_F(IndexMapTest, MapsResidueBlocks)
548 const int maxGroup[] = {
549 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
551 const int evalGroup[] = { 3, 4, 7, 8, 13 };
552 topManager_.initAtoms(18);
553 topManager_.initUniformResidues(3);
554 testInit(maxGroup, INDEX_RES);
555 testUpdate(evalGroup, false, NULL);
558 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
560 const int maxGroup[] = {
561 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
563 const int evalGroup[] = { 3, 4, 7, 8, 13 };
564 topManager_.initAtoms(18);
565 topManager_.initUniformResidues(3);
566 testInit(maxGroup, INDEX_RES);
567 testUpdate(evalGroup, true, NULL);
570 TEST_F(IndexMapTest, HandlesMultipleRequests)
572 const int maxGroup[] = {
573 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
575 const int evalGroup[] = { 3, 4, 7, 8, 13 };
576 topManager_.initAtoms(18);
577 topManager_.initUniformResidues(3);
578 testInit(maxGroup, INDEX_RES);
579 testUpdate(evalGroup, false, "EvaluatedNoMask");
580 testUpdate(evalGroup, true, "EvaluatedMask");
581 testUpdate(maxGroup, true, "Initialized");
582 testUpdate(evalGroup, true, "EvaluatedMask");
583 testUpdate(evalGroup, false, "EvaluatedNoMask");
584 testUpdate(maxGroup, false, "Initialized");