6628cda62da473592dd0809d597d2d6da07a9412
[alexxy/gromacs.git] / src / gromacs / selection / tests / indexutil.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Tests the index group handling in the selection engine.
38  *
39  * \todo
40  * Tests for other functions, at least the set operations.
41  *
42  * \author Teemu Murtola <teemu.murtola@gmail.com>
43  * \ingroup module_selection
44  */
45 #include "gmxpre.h"
46
47 #include "gromacs/selection/indexutil.h"
48
49 #include <gtest/gtest.h>
50
51 #include "gromacs/topology/block.h"
52 #include "gromacs/utility/arrayref.h"
53
54 #include "testutils/refdata.h"
55
56 #include "toputils.h"
57
58 namespace
59 {
60
61 //! Helper for creating groups from an array.
62 gmx_ana_index_t initGroup(gmx::ArrayRef<int> index)
63 {
64     gmx_ana_index_t g = { static_cast<int>(index.size()), index.data(), 0 };
65     return g;
66 }
67
68 TEST(IndexGroupTest, RemovesDuplicates)
69 {
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));
76 }
77
78 //! Text fixture for index block operations
79 class IndexBlockTest : public ::testing::Test
80 {
81     public:
82         IndexBlockTest();
83         ~IndexBlockTest() override;
84
85         //@{
86         //! Set the input group for the index block operation
87         void setGroup(int count, const int atoms[]);
88         template <int count>
89         void setGroup(const int (&atoms)[count])
90         {
91             setGroup(count, atoms);
92         }
93         //@}
94         /*! \brief Check the input group and output with refdata, with
95          * an optional \c id to name the refdata block */
96         void checkBlocka(const char *id = "Block");
97         //! Make a simple topology to check with
98         void makeSimpleTopology();
99         //! Make a complex topology to check with
100         void makeComplexTopology();
101
102         //! Managers reference data for the tests
103         gmx::test::TestReferenceData    data_;
104         //! Manages setting up a topology and matching data structures
105         gmx::test::TopologyManager      topManager_;
106         //! The input group to test with
107         gmx_ana_index_t                 g_;
108         //! The output block to test on
109         t_blocka                        blocka_;
110 };
111
112 IndexBlockTest::IndexBlockTest()
113 {
114     blocka_.nr           = 0;
115     blocka_.index        = nullptr;
116     blocka_.nalloc_index = 0;
117     blocka_.nra          = 0;
118     blocka_.a            = nullptr;
119     blocka_.nalloc_a     = 0;
120     gmx_ana_index_clear(&g_);
121 }
122
123 IndexBlockTest::~IndexBlockTest()
124 {
125     done_blocka(&blocka_);
126 }
127
128 void IndexBlockTest::setGroup(int count, const int atoms[])
129 {
130     g_.isize = count;
131     g_.index = const_cast<int *>(atoms);
132 }
133
134 void IndexBlockTest::checkBlocka(const char *id)
135 {
136     gmx::test::TestReferenceChecker compound(
137             data_.rootChecker().checkCompound("BlockAtoms", id));
138     compound.checkSequenceArray(g_.isize, g_.index, "Input");
139     compound.checkInteger(blocka_.nr, "Count");
140     for (int i = 0; i < blocka_.nr; ++i)
141     {
142         gmx::test::TestReferenceChecker blockCompound(
143                 compound.checkCompound("Block", nullptr));
144         blockCompound.checkSequence(&blocka_.a[blocka_.index[i]],
145                                     &blocka_.a[blocka_.index[i+1]],
146                                     "Atoms");
147     }
148 }
149
150 void IndexBlockTest::makeSimpleTopology()
151 {
152     topManager_.initTopology(1, 1);
153     {
154         int              moleculeTypeIndex   = 0;
155         std::vector<int> numAtomsInResidues  = {3, 3, 3};
156         int              moleculeBlockIndex  = 0;
157         int              numMoleculesInBlock = 1;
158         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
159         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
160     }
161     topManager_.finalizeTopology();
162 }
163
164 void IndexBlockTest::makeComplexTopology()
165 {
166     topManager_.initTopology(3, 4);
167     {
168         int              moleculeTypeIndex   = 0;
169         std::vector<int> numAtomsInResidues  = {2, 2, 1};
170         int              moleculeBlockIndex  = 0;
171         int              numMoleculesInBlock = 1;
172         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
173         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
174     }
175     {
176         int              moleculeTypeIndex   = 1;
177         std::vector<int> numAtomsInResidues  = {1};
178         int              moleculeBlockIndex  = 1;
179         int              numMoleculesInBlock = 3;
180         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
181         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
182     }
183     {
184         int              moleculeTypeIndex   = 2;
185         std::vector<int> numAtomsInResidues  = {3};
186         int              moleculeBlockIndex  = 2;
187         int              numMoleculesInBlock = 1;
188         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
189         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
190     }
191     {
192         int moleculeTypeIndex   = 0; // Re-using a moltype definition
193         int moleculeBlockIndex  = 3;
194         int numMoleculesInBlock = 1;
195         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
196     }
197     topManager_.finalizeTopology();
198 }
199
200 /********************************************************************
201  * gmx_ana_index_make_block() tests
202  */
203
204 TEST_F(IndexBlockTest, CreatesUnknownBlock)
205 {
206     gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
207     checkBlocka();
208     done_blocka(&blocka_);
209     gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
210     checkBlocka();
211 }
212
213 TEST_F(IndexBlockTest, CreatesAtomBlock)
214 {
215     const int group[] = { 0, 1, 3, 4, 6 };
216     setGroup(group);
217     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, false);
218     checkBlocka();
219     done_blocka(&blocka_);
220     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, true);
221     checkBlocka();
222 }
223
224 TEST_F(IndexBlockTest, CreatesResidueBlocksForSimpleTopology)
225 {
226     makeSimpleTopology();
227
228     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8 };
229     setGroup(group);
230     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
231     checkBlocka("FromAllAtoms");
232     done_blocka(&blocka_);
233     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
234     checkBlocka("FromAllAtoms");
235 }
236
237 TEST_F(IndexBlockTest, CreatesResidueBlocksForComplexTopology)
238 {
239     makeComplexTopology();
240
241     SCOPED_TRACE("Group with all atoms without completion");
242     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15 };
243     setGroup(group);
244     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
245     checkBlocka("FromAllAtoms");
246     done_blocka(&blocka_);
247     //SCOPED_TRACE("Group with all atoms with completion");
248     //gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
249     //checkBlocka("FromAllAtoms");
250     //done_blocka(&blocka_);
251
252     SCOPED_TRACE("Group with some atoms without completion");
253     const int subgroup[] = { 0, 3, 4, 5, 6, 7, 8, 12, 13, 15 };
254     setGroup(subgroup);
255     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
256     checkBlocka("FromSomeAtomsWithoutCompletion");
257     //done_blocka(&blocka_);
258     //SCOPED_TRACE("Group with some atoms with completion");
259     //gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
260     //checkBlocka("FromSomeAtomsWithCompletion");
261 }
262
263 TEST_F(IndexBlockTest, CreatesMoleculeBlocksForSimpleTopology)
264 {
265     makeSimpleTopology();
266
267     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8 };
268     setGroup(group);
269     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
270     checkBlocka("FromAllAtoms");
271     done_blocka(&blocka_);
272     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
273     checkBlocka("FromAllAtoms");
274 }
275
276 TEST_F(IndexBlockTest, CreatesMoleculeBlocksForComplexTopology)
277 {
278     makeComplexTopology();
279
280     SCOPED_TRACE("Group with all atoms without completion");
281     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15 };
282     setGroup(group);
283     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
284     checkBlocka("FromAllAtoms");
285     done_blocka(&blocka_);
286     //SCOPED_TRACE("Group with all atoms with completion");
287     //gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
288     //checkBlocka("FromAllAtoms");
289     //done_blocka(&blocka_);
290
291     SCOPED_TRACE("Group with some atoms without completion");
292     const int subgroup[] = { 1, 5, 6, 7, 11, 12 };
293     setGroup(subgroup);
294     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
295     checkBlocka("FromSomeAtomsWithoutCompletion");
296     //done_blocka(&blocka_);
297     //SCOPED_TRACE("Group with some atoms with completion");
298     //gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
299     //checkBlocka("FromSomeAtomsWithCompletion");
300 }
301
302 TEST_F(IndexBlockTest, CreatesSingleBlock)
303 {
304     const int group[] = { 0, 1, 3, 4, 6 };
305     setGroup(group);
306     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, false);
307     checkBlocka();
308     done_blocka(&blocka_);
309     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, true);
310     checkBlocka();
311 }
312
313 /********************************************************************
314  * gmx_ana_index_has_full_ablocks() tests
315  */
316
317 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
318 {
319     const int maxGroup[] = {
320         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
321     };
322     const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
323     topManager_.initAtoms(18);
324     topManager_.initUniformResidues(3);
325     setGroup(maxGroup);
326     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
327                              INDEX_RES, false);
328     setGroup(testGroup);
329     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
330 }
331
332 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
333 {
334     const int maxGroup[] = {
335         15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
336     };
337     const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
338     topManager_.initAtoms(18);
339     topManager_.initUniformResidues(3);
340     setGroup(maxGroup);
341     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
342                              INDEX_RES, false);
343     setGroup(testGroup);
344     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
345 }
346
347 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
348 {
349     const int maxGroup[] = {
350         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
351     };
352     const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
353     const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
354     const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
355
356     topManager_.initAtoms(18);
357     topManager_.initUniformResidues(3);
358     setGroup(maxGroup);
359     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
360                              INDEX_RES, false);
361
362     setGroup(testGroup1);
363     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
364
365     setGroup(testGroup2);
366     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
367
368     setGroup(testGroup3);
369     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
370 }
371
372 /********************************************************************
373  * gmx_ana_index_has_complete_elems() tests
374  */
375
376 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
377 {
378     const int group[] = { 0, 1, 2 };
379     setGroup(group);
380     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, nullptr));
381     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, nullptr));
382     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, nullptr));
383 }
384
385 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
386 {
387     const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
388     const int group2[] = { 3, 4, 5, 6, 7, 8 };
389
390     topManager_.initAtoms(15);
391     topManager_.initUniformResidues(3);
392     gmx_mtop_t *top = topManager_.topology();
393
394     setGroup(group1);
395     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
396
397     setGroup(group2);
398     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
399 }
400
401 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
402 {
403     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
404     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
405     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
406     const int group4[] = { 3, 4, 5, 6, 8, 12, 13, 14 };
407
408     topManager_.initAtoms(18);
409     topManager_.initUniformResidues(3);
410     gmx_mtop_t *top = topManager_.topology();
411
412     setGroup(group1);
413     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
414
415     setGroup(group2);
416     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
417
418     setGroup(group3);
419     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
420
421     setGroup(group4);
422     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
423 }
424
425 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
426 {
427     const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
428
429     topManager_.initAtoms(15);
430     topManager_.initUniformResidues(1);
431     topManager_.initUniformMolecules(3);
432     gmx_mtop_t *top = topManager_.topology();
433
434     setGroup(group);
435     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
436 }
437
438 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
439 {
440     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
441     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
442     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
443
444     topManager_.initAtoms(18);
445     topManager_.initUniformResidues(1);
446     topManager_.initUniformMolecules(3);
447     gmx_mtop_t *top = topManager_.topology();
448
449     setGroup(group1);
450     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
451
452     setGroup(group2);
453     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
454
455     setGroup(group3);
456     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
457 }
458
459 /********************************************************************
460  * IndexMapTest
461  */
462
463 class IndexMapTest : public ::testing::Test
464 {
465     public:
466         IndexMapTest();
467         ~IndexMapTest() override;
468
469         void testInit(int atomCount, const int atoms[], e_index_t type);
470         void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
471                         const char *name);
472         void testOrgIdGroup(e_index_t type, const char *name);
473         template <int count>
474         void testInit(const int (&atoms)[count], e_index_t type)
475         {
476             testInit(count, atoms, type);
477         }
478         template <int count>
479         void testUpdate(const int (&atoms)[count], bool bMaskOnly,
480                         const char *name)
481         {
482             testUpdate(count, atoms, bMaskOnly, name);
483         }
484
485         void checkMapping(int atomCount, const int atoms[], const char *name);
486
487         gmx::test::TestReferenceData    data_;
488         gmx::test::TestReferenceChecker checker_;
489         gmx::test::TopologyManager      topManager_;
490         gmx_ana_indexmap_t              map_;
491
492     private:
493         gmx_ana_index_t                 initGroup_;
494 };
495
496 IndexMapTest::IndexMapTest()
497     : checker_(data_.rootChecker())
498 {
499     gmx_ana_indexmap_clear(&map_);
500     gmx_ana_index_clear(&initGroup_);
501 }
502
503 IndexMapTest::~IndexMapTest()
504 {
505     gmx_ana_indexmap_deinit(&map_);
506 }
507
508 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
509 {
510     initGroup_.isize = atomCount;
511     initGroup_.index = const_cast<int *>(atoms);
512     gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
513     EXPECT_EQ(type, map_.type);
514     checkMapping(atomCount, atoms, "Initialized");
515 }
516
517 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
518                               const char *name)
519 {
520     gmx_ana_index_t g;
521     g.isize = atomCount;
522     g.index = const_cast<int *>(atoms);
523     gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
524     if (name == nullptr)
525     {
526         name = "Updated";
527     }
528     if (bMaskOnly)
529     {
530         checkMapping(initGroup_.isize, initGroup_.index, name);
531     }
532     else
533     {
534         checkMapping(atomCount, atoms, name);
535     }
536 }
537
538 void IndexMapTest::testOrgIdGroup(e_index_t type, const char *name)
539 {
540     gmx::test::TestReferenceChecker compound(
541             checker_.checkCompound("OrgIdGroups", name));
542     const int count
543         = gmx_ana_indexmap_init_orgid_group(&map_, topManager_.topology(), type);
544     compound.checkInteger(count, "GroupCount");
545     compound.checkSequenceArray(map_.mapb.nr, map_.orgid, "OrgId");
546     for (int i = 0; i < map_.mapb.nr; ++i)
547     {
548         EXPECT_EQ(map_.orgid[i], map_.mapid[i]);
549     }
550 }
551
552 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
553                                 const char *name)
554 {
555     gmx::test::TestReferenceChecker compound(
556             checker_.checkCompound("IndexMapping", name));
557     compound.checkSequenceArray(atomCount, atoms, "Input");
558     compound.checkInteger(map_.mapb.nr, "Count");
559     for (int i = 0; i < map_.mapb.nr; ++i)
560     {
561         gmx::test::TestReferenceChecker blockCompound(
562                 compound.checkCompound("Block", nullptr));
563         blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
564                                     &atoms[map_.mapb.index[i+1]],
565                                     "Atoms");
566         blockCompound.checkInteger(map_.refid[i], "RefId");
567         blockCompound.checkInteger(map_.mapid[i], "MapId");
568         int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
569         EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
570     }
571 }
572
573 /********************************************************************
574  * gmx_ana_indexmap_t tests
575  */
576
577 TEST_F(IndexMapTest, InitializesAtomBlock)
578 {
579     const int maxGroup[] = { 1, 2, 4, 5 };
580     testInit(maxGroup, INDEX_ATOM);
581 }
582
583 TEST_F(IndexMapTest, InitializesOrgIdGroupAtom)
584 {
585     const int maxGroup[] = { 2, 5, 7 };
586     testInit(maxGroup, INDEX_ATOM);
587     testOrgIdGroup(INDEX_ATOM, "Atoms");
588 }
589
590 TEST_F(IndexMapTest, InitializesOrgIdGroupSingle)
591 {
592     const int maxGroup[] = { 3, 4, 7, 8, 13 };
593     topManager_.initAtoms(18);
594     topManager_.initUniformResidues(3);
595     testInit(maxGroup, INDEX_RES);
596     testOrgIdGroup(INDEX_ATOM, "Single");
597 }
598
599 TEST_F(IndexMapTest, InitializesOrgIdGroupResidue)
600 {
601     const int maxGroup[] = { 3, 4, 7, 8, 13 };
602     topManager_.initAtoms(18);
603     topManager_.initUniformResidues(3);
604     testInit(maxGroup, INDEX_ATOM);
605     testOrgIdGroup(INDEX_RES, "Residues");
606 }
607
608 TEST_F(IndexMapTest, InitializesOrgIdGroupMolecule)
609 {
610     const int maxGroup[] = { 1, 2, 3, 4, 7, 8, 13 };
611     topManager_.initAtoms(18);
612     topManager_.initUniformResidues(3);
613     topManager_.initUniformMolecules(6);
614     testInit(maxGroup, INDEX_RES);
615     testOrgIdGroup(INDEX_MOL, "Molecules");
616 }
617
618 TEST_F(IndexMapTest, InitializesOrgIdGroupAll)
619 {
620     const int maxGroup[] = { 3, 4, 7, 8, 13 };
621     testInit(maxGroup, INDEX_ATOM);
622     testOrgIdGroup(INDEX_ALL, "All");
623 }
624
625 TEST_F(IndexMapTest, InitializesMoleculeBlock)
626 {
627     const int maxGroup[] = { 3, 4, 7, 8, 13 };
628     topManager_.initAtoms(18);
629     topManager_.initUniformResidues(1);
630     topManager_.initUniformMolecules(3);
631     testInit(maxGroup, INDEX_MOL);
632 }
633
634 TEST_F(IndexMapTest, MapsSingleBlock)
635 {
636     const int maxGroup[]  = { 0, 1, 2, 3 };
637     const int evalGroup[] = { 0, 2 };
638     testInit(maxGroup, INDEX_ALL);
639     testUpdate(evalGroup, false, nullptr);
640 }
641
642 TEST_F(IndexMapTest, MapsResidueBlocks)
643 {
644     const int maxGroup[] = {
645         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
646     };
647     const int evalGroup[] = { 3, 4, 7, 8, 13 };
648     topManager_.initAtoms(18);
649     topManager_.initUniformResidues(3);
650     testInit(maxGroup, INDEX_RES);
651     testUpdate(evalGroup, false, nullptr);
652 }
653
654 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
655 {
656     const int maxGroup[] = {
657         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
658     };
659     const int evalGroup[] = { 3, 4, 7, 8, 13 };
660     topManager_.initAtoms(18);
661     topManager_.initUniformResidues(3);
662     testInit(maxGroup, INDEX_RES);
663     testUpdate(evalGroup, true, nullptr);
664 }
665
666 TEST_F(IndexMapTest, HandlesMultipleRequests)
667 {
668     const int maxGroup[] = {
669         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
670     };
671     const int evalGroup[] = { 3, 4, 7, 8, 13 };
672     topManager_.initAtoms(18);
673     topManager_.initUniformResidues(3);
674     testInit(maxGroup, INDEX_RES);
675     testUpdate(evalGroup, false, "EvaluatedNoMask");
676     testUpdate(evalGroup, true, "EvaluatedMask");
677     testUpdate(maxGroup, true, "Initialized");
678     testUpdate(evalGroup, true, "EvaluatedMask");
679     testUpdate(evalGroup, false, "EvaluatedNoMask");
680     testUpdate(maxGroup, false, "Initialized");
681 }
682
683 } // namespace