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
45 #include <gtest/gtest.h>
47 #include "gromacs/selection/indexutil.h"
48 #include "gromacs/topology/block.h"
50 #include "testutils/refdata.h"
57 /********************************************************************
61 class IndexBlockTest : public ::testing::Test
67 void setGroup(int count, const int atoms[]);
69 void setGroup(const int (&atoms)[count])
71 setGroup(count, atoms);
76 gmx::test::TestReferenceData data_;
77 gmx::test::TopologyManager topManager_;
82 IndexBlockTest::IndexBlockTest()
86 blocka_.nalloc_index = 0;
90 gmx_ana_index_clear(&g_);
93 IndexBlockTest::~IndexBlockTest()
95 done_blocka(&blocka_);
98 void IndexBlockTest::setGroup(int count, const int atoms[])
101 g_.index = const_cast<int *>(atoms);
104 void IndexBlockTest::checkBlocka()
106 gmx::test::TestReferenceChecker compound(
107 data_.rootChecker().checkCompound("BlockAtoms", "Block"));
108 compound.checkSequenceArray(g_.isize, g_.index, "Input");
109 compound.checkInteger(blocka_.nr, "Count");
110 for (int i = 0; i < blocka_.nr; ++i)
112 gmx::test::TestReferenceChecker blockCompound(
113 compound.checkCompound("Block", NULL));
114 blockCompound.checkSequence(&blocka_.a[blocka_.index[i]],
115 &blocka_.a[blocka_.index[i+1]],
120 /********************************************************************
121 * gmx_ana_index_make_block() tests
124 TEST_F(IndexBlockTest, CreatesUnknownBlock)
126 gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
128 done_blocka(&blocka_);
129 gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
133 TEST_F(IndexBlockTest, CreatesAtomBlock)
135 const int group[] = { 0, 1, 3, 4, 6 };
137 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, false);
139 done_blocka(&blocka_);
140 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, true);
144 TEST_F(IndexBlockTest, CreatesResidueBlock)
146 const int group[] = { 0, 1, 3, 6, 7 };
147 topManager_.initAtoms(9);
148 topManager_.initUniformResidues(3);
150 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
155 TEST_F(IndexBlockTest, CreatesMoleculeBlock)
157 const int group[] = { 3, 4, 7, 8, 13 };
158 topManager_.initAtoms(18);
159 topManager_.initUniformMolecules(3);
161 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
166 TEST_F(IndexBlockTest, CreatesResidueBlockWithCompletion)
168 const int group[] = { 3, 4, 7, 8, 13 };
169 topManager_.initAtoms(18);
170 topManager_.initUniformResidues(3);
172 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
177 TEST_F(IndexBlockTest, CreatesMoleculeBlockWithCompletion)
179 const int group[] = { 3, 4, 7, 8, 13 };
180 topManager_.initAtoms(18);
181 topManager_.initUniformMolecules(3);
183 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
188 TEST_F(IndexBlockTest, CreatesSingleBlock)
190 const int group[] = { 0, 1, 3, 4, 6 };
192 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, false);
194 done_blocka(&blocka_);
195 gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, true);
199 /********************************************************************
200 * gmx_ana_index_has_full_ablocks() tests
203 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
205 const int maxGroup[] = {
206 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
208 const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
209 topManager_.initAtoms(18);
210 topManager_.initUniformResidues(3);
212 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
215 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
218 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
220 const int maxGroup[] = {
221 15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
223 const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
224 topManager_.initAtoms(18);
225 topManager_.initUniformResidues(3);
227 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
230 EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
233 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
235 const int maxGroup[] = {
236 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
238 const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
239 const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
240 const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
242 topManager_.initAtoms(18);
243 topManager_.initUniformResidues(3);
245 gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
248 setGroup(testGroup1);
249 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
251 setGroup(testGroup2);
252 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
254 setGroup(testGroup3);
255 EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
258 /********************************************************************
259 * gmx_ana_index_has_complete_elems() tests
262 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
264 const int group[] = { 0, 1, 2 };
266 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, NULL));
267 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, NULL));
268 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, NULL));
271 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
273 const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
274 const int group2[] = { 3, 4, 5, 6, 7, 8 };
276 topManager_.initAtoms(15);
277 topManager_.initUniformResidues(3);
278 t_topology *top = topManager_.topology();
281 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
284 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
287 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
289 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
290 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
291 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
293 topManager_.initAtoms(18);
294 topManager_.initUniformResidues(3);
295 t_topology *top = topManager_.topology();
298 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
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 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
309 const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
311 topManager_.initAtoms(15);
312 topManager_.initUniformMolecules(3);
313 t_topology *top = topManager_.topology();
316 EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
319 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
321 const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
322 const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
323 const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
325 topManager_.initAtoms(18);
326 topManager_.initUniformMolecules(3);
327 t_topology *top = topManager_.topology();
330 EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
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 /********************************************************************
343 class IndexMapTest : public ::testing::Test
349 void testInit(int atomCount, const int atoms[], e_index_t type);
350 void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
353 void testInit(const int (&atoms)[count], e_index_t type)
355 testInit(count, atoms, type);
358 void testUpdate(const int (&atoms)[count], bool bMaskOnly,
361 testUpdate(count, atoms, bMaskOnly, name);
364 void checkMapping(int atomCount, const int atoms[], const char *name);
366 gmx::test::TestReferenceData data_;
367 gmx::test::TestReferenceChecker checker_;
368 gmx::test::TopologyManager topManager_;
369 gmx_ana_indexmap_t map_;
372 gmx_ana_index_t initGroup_;
375 IndexMapTest::IndexMapTest()
376 : checker_(data_.rootChecker())
378 gmx_ana_indexmap_clear(&map_);
379 gmx_ana_index_clear(&initGroup_);
382 IndexMapTest::~IndexMapTest()
384 gmx_ana_indexmap_deinit(&map_);
387 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
389 initGroup_.isize = atomCount;
390 initGroup_.index = const_cast<int *>(atoms);
391 gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
392 EXPECT_EQ(type, map_.type);
393 checkMapping(atomCount, atoms, "Initialized");
396 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
401 g.index = const_cast<int *>(atoms);
402 gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
409 checkMapping(initGroup_.isize, initGroup_.index, name);
413 checkMapping(atomCount, atoms, name);
417 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
420 gmx::test::TestReferenceChecker compound(
421 checker_.checkCompound("IndexMapping", name));
422 compound.checkSequenceArray(atomCount, atoms, "Input");
423 compound.checkInteger(map_.mapb.nr, "Count");
424 for (int i = 0; i < map_.mapb.nr; ++i)
426 gmx::test::TestReferenceChecker blockCompound(
427 compound.checkCompound("Block", NULL));
428 blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
429 &atoms[map_.mapb.index[i+1]],
431 blockCompound.checkInteger(map_.refid[i], "RefId");
432 blockCompound.checkInteger(map_.mapid[i], "MapId");
433 int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
434 EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
438 /********************************************************************
439 * gmx_ana_indexmap_t tests
442 TEST_F(IndexMapTest, InitializesAtomBlock)
444 const int maxGroup[] = { 1, 2, 4, 5 };
445 testInit(maxGroup, INDEX_ATOM);
448 TEST_F(IndexMapTest, InitializesMoleculeBlock)
450 const int maxGroup[] = { 3, 4, 7, 8, 13 };
451 topManager_.initAtoms(18);
452 topManager_.initUniformMolecules(3);
453 testInit(maxGroup, INDEX_MOL);
456 TEST_F(IndexMapTest, MapsSingleBlock)
458 const int maxGroup[] = { 0, 1, 2, 3 };
459 const int evalGroup[] = { 0, 2 };
460 testInit(maxGroup, INDEX_ALL);
461 testUpdate(evalGroup, false, NULL);
464 TEST_F(IndexMapTest, MapsResidueBlocks)
466 const int maxGroup[] = {
467 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
469 const int evalGroup[] = { 3, 4, 7, 8, 13 };
470 topManager_.initAtoms(18);
471 topManager_.initUniformResidues(3);
472 testInit(maxGroup, INDEX_RES);
473 testUpdate(evalGroup, false, NULL);
476 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
478 const int maxGroup[] = {
479 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
481 const int evalGroup[] = { 3, 4, 7, 8, 13 };
482 topManager_.initAtoms(18);
483 topManager_.initUniformResidues(3);
484 testInit(maxGroup, INDEX_RES);
485 testUpdate(evalGroup, true, NULL);
488 TEST_F(IndexMapTest, HandlesMultipleRequests)
490 const int maxGroup[] = {
491 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
493 const int evalGroup[] = { 3, 4, 7, 8, 13 };
494 topManager_.initAtoms(18);
495 topManager_.initUniformResidues(3);
496 testInit(maxGroup, INDEX_RES);
497 testUpdate(evalGroup, false, "EvaluatedNoMask");
498 testUpdate(evalGroup, true, "EvaluatedMask");
499 testUpdate(maxGroup, true, "Initialized");
500 testUpdate(evalGroup, true, "EvaluatedMask");
501 testUpdate(evalGroup, false, "EvaluatedNoMask");
502 testUpdate(maxGroup, false, "Initialized");