2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017,2018, 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 = nullptr;
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", nullptr));
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_, nullptr, nullptr, INDEX_UNKNOWN, false);
149 done_blocka(&blocka_);
150 gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
154 TEST_F(IndexBlockTest, CreatesAtomBlock)
156 const int group[] = { 0, 1, 3, 4, 6 };
158 gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, false);
160 done_blocka(&blocka_);
161 gmx_ana_index_make_block(&blocka_, nullptr, &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_.initUniformResidues(1);
181 topManager_.initUniformMolecules(3);
183 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
188 TEST_F(IndexBlockTest, CreatesResidueBlockWithCompletion)
190 const int group[] = { 3, 4, 7, 8, 13 };
191 topManager_.initAtoms(18);
192 topManager_.initUniformResidues(3);
194 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
199 TEST_F(IndexBlockTest, CreatesMoleculeBlockWithCompletion)
201 const int group[] = { 3, 4, 7, 8, 13 };
202 topManager_.initAtoms(18);
203 topManager_.initUniformResidues(1);
204 topManager_.initUniformMolecules(3);
206 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
211 TEST_F(IndexBlockTest, CreatesSingleBlock)
213 const int group[] = { 0, 1, 3, 4, 6 };
215 gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, false);
217 done_blocka(&blocka_);
218 gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, true);
222 /********************************************************************
223 * gmx_ana_index_has_full_ablocks() tests
226 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
228 const int maxGroup[] = {
229 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
231 const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
232 topManager_.initAtoms(18);
233 topManager_.initUniformResidues(3);
235 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
238 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
241 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
243 const int maxGroup[] = {
244 15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
246 const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
247 topManager_.initAtoms(18);
248 topManager_.initUniformResidues(3);
250 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
253 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
256 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
258 const int maxGroup[] = {
259 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
261 const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
262 const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
263 const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
265 topManager_.initAtoms(18);
266 topManager_.initUniformResidues(3);
268 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
271 setGroup(testGroup1);
272 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
274 setGroup(testGroup2);
275 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
277 setGroup(testGroup3);
278 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
281 /********************************************************************
282 * gmx_ana_index_has_complete_elems() tests
285 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
287 const int group[] = { 0, 1, 2 };
289 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, nullptr));
290 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, nullptr));
291 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, nullptr));
294 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
296 const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
297 const int group2[] = { 3, 4, 5, 6, 7, 8 };
299 topManager_.initAtoms(15);
300 topManager_.initUniformResidues(3);
301 gmx_mtop_t *top = topManager_.topology();
304 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
307 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
310 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
312 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
313 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
314 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
315 const int group4[] = { 3, 4, 5, 6, 8, 12, 13, 14 };
317 topManager_.initAtoms(18);
318 topManager_.initUniformResidues(3);
319 gmx_mtop_t *top = topManager_.topology();
322 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
325 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
328 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
331 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
334 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
336 const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
338 topManager_.initAtoms(15);
339 topManager_.initUniformResidues(1);
340 topManager_.initUniformMolecules(3);
341 gmx_mtop_t *top = topManager_.topology();
344 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
347 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
349 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
350 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
351 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
353 topManager_.initAtoms(18);
354 topManager_.initUniformResidues(1);
355 topManager_.initUniformMolecules(3);
356 gmx_mtop_t *top = topManager_.topology();
359 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
362 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
365 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
368 /********************************************************************
372 class IndexMapTest : public ::testing::Test
378 void testInit(int atomCount, const int atoms[], e_index_t type);
379 void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
381 void testOrgIdGroup(e_index_t type, const char *name);
383 void testInit(const int (&atoms)[count], e_index_t type)
385 testInit(count, atoms, type);
388 void testUpdate(const int (&atoms)[count], bool bMaskOnly,
391 testUpdate(count, atoms, bMaskOnly, name);
394 void checkMapping(int atomCount, const int atoms[], const char *name);
396 gmx::test::TestReferenceData data_;
397 gmx::test::TestReferenceChecker checker_;
398 gmx::test::TopologyManager topManager_;
399 gmx_ana_indexmap_t map_;
402 gmx_ana_index_t initGroup_;
405 IndexMapTest::IndexMapTest()
406 : checker_(data_.rootChecker())
408 gmx_ana_indexmap_clear(&map_);
409 gmx_ana_index_clear(&initGroup_);
412 IndexMapTest::~IndexMapTest()
414 gmx_ana_indexmap_deinit(&map_);
417 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
419 initGroup_.isize = atomCount;
420 initGroup_.index = const_cast<int *>(atoms);
421 gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
422 EXPECT_EQ(type, map_.type);
423 checkMapping(atomCount, atoms, "Initialized");
426 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
431 g.index = const_cast<int *>(atoms);
432 gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
439 checkMapping(initGroup_.isize, initGroup_.index, name);
443 checkMapping(atomCount, atoms, name);
447 void IndexMapTest::testOrgIdGroup(e_index_t type, const char *name)
449 gmx::test::TestReferenceChecker compound(
450 checker_.checkCompound("OrgIdGroups", name));
452 = gmx_ana_indexmap_init_orgid_group(&map_, topManager_.topology(), type);
453 compound.checkInteger(count, "GroupCount");
454 compound.checkSequenceArray(map_.mapb.nr, map_.orgid, "OrgId");
455 for (int i = 0; i < map_.mapb.nr; ++i)
457 EXPECT_EQ(map_.orgid[i], map_.mapid[i]);
461 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
464 gmx::test::TestReferenceChecker compound(
465 checker_.checkCompound("IndexMapping", name));
466 compound.checkSequenceArray(atomCount, atoms, "Input");
467 compound.checkInteger(map_.mapb.nr, "Count");
468 for (int i = 0; i < map_.mapb.nr; ++i)
470 gmx::test::TestReferenceChecker blockCompound(
471 compound.checkCompound("Block", nullptr));
472 blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
473 &atoms[map_.mapb.index[i+1]],
475 blockCompound.checkInteger(map_.refid[i], "RefId");
476 blockCompound.checkInteger(map_.mapid[i], "MapId");
477 int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
478 EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
482 /********************************************************************
483 * gmx_ana_indexmap_t tests
486 TEST_F(IndexMapTest, InitializesAtomBlock)
488 const int maxGroup[] = { 1, 2, 4, 5 };
489 testInit(maxGroup, INDEX_ATOM);
492 TEST_F(IndexMapTest, InitializesOrgIdGroupAtom)
494 const int maxGroup[] = { 2, 5, 7 };
495 testInit(maxGroup, INDEX_ATOM);
496 testOrgIdGroup(INDEX_ATOM, "Atoms");
499 TEST_F(IndexMapTest, InitializesOrgIdGroupSingle)
501 const int maxGroup[] = { 3, 4, 7, 8, 13 };
502 topManager_.initAtoms(18);
503 topManager_.initUniformResidues(3);
504 testInit(maxGroup, INDEX_RES);
505 testOrgIdGroup(INDEX_ATOM, "Single");
508 TEST_F(IndexMapTest, InitializesOrgIdGroupResidue)
510 const int maxGroup[] = { 3, 4, 7, 8, 13 };
511 topManager_.initAtoms(18);
512 topManager_.initUniformResidues(3);
513 testInit(maxGroup, INDEX_ATOM);
514 testOrgIdGroup(INDEX_RES, "Residues");
517 TEST_F(IndexMapTest, InitializesOrgIdGroupMolecule)
519 const int maxGroup[] = { 1, 2, 3, 4, 7, 8, 13 };
520 topManager_.initAtoms(18);
521 topManager_.initUniformResidues(3);
522 topManager_.initUniformMolecules(6);
523 testInit(maxGroup, INDEX_RES);
524 testOrgIdGroup(INDEX_MOL, "Molecules");
527 TEST_F(IndexMapTest, InitializesOrgIdGroupAll)
529 const int maxGroup[] = { 3, 4, 7, 8, 13 };
530 testInit(maxGroup, INDEX_ATOM);
531 testOrgIdGroup(INDEX_ALL, "All");
534 TEST_F(IndexMapTest, InitializesMoleculeBlock)
536 const int maxGroup[] = { 3, 4, 7, 8, 13 };
537 topManager_.initAtoms(18);
538 topManager_.initUniformResidues(1);
539 topManager_.initUniformMolecules(3);
540 testInit(maxGroup, INDEX_MOL);
543 TEST_F(IndexMapTest, MapsSingleBlock)
545 const int maxGroup[] = { 0, 1, 2, 3 };
546 const int evalGroup[] = { 0, 2 };
547 testInit(maxGroup, INDEX_ALL);
548 testUpdate(evalGroup, false, nullptr);
551 TEST_F(IndexMapTest, MapsResidueBlocks)
553 const int maxGroup[] = {
554 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
556 const int evalGroup[] = { 3, 4, 7, 8, 13 };
557 topManager_.initAtoms(18);
558 topManager_.initUniformResidues(3);
559 testInit(maxGroup, INDEX_RES);
560 testUpdate(evalGroup, false, nullptr);
563 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
565 const int maxGroup[] = {
566 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
568 const int evalGroup[] = { 3, 4, 7, 8, 13 };
569 topManager_.initAtoms(18);
570 topManager_.initUniformResidues(3);
571 testInit(maxGroup, INDEX_RES);
572 testUpdate(evalGroup, true, nullptr);
575 TEST_F(IndexMapTest, HandlesMultipleRequests)
577 const int maxGroup[] = {
578 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
580 const int evalGroup[] = { 3, 4, 7, 8, 13 };
581 topManager_.initAtoms(18);
582 topManager_.initUniformResidues(3);
583 testInit(maxGroup, INDEX_RES);
584 testUpdate(evalGroup, false, "EvaluatedNoMask");
585 testUpdate(evalGroup, true, "EvaluatedMask");
586 testUpdate(maxGroup, true, "Initialized");
587 testUpdate(evalGroup, true, "EvaluatedMask");
588 testUpdate(evalGroup, false, "EvaluatedNoMask");
589 testUpdate(maxGroup, false, "Initialized");