7bce7c26b9f9e48bbe69f1b94767ffef2fcd8cad
[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 <gmock/gmock.h>
50 #include <gtest/gtest.h>
51
52 #include "gromacs/topology/block.h"
53 #include "gromacs/utility/arrayref.h"
54 #include "gromacs/utility/smalloc.h"
55
56 #include "testutils/refdata.h"
57 #include "testutils/testasserts.h"
58
59 #include "toputils.h"
60
61 namespace
62 {
63
64 //! Helper for creating groups from an array.
65 gmx_ana_index_t initGroup(gmx::ArrayRef<int> index)
66 {
67     gmx_ana_index_t g = { static_cast<int>(index.size()), index.data(), 0 };
68     return g;
69 }
70
71 TEST(IndexGroupTest, RemovesDuplicates)
72 {
73     int             index[]    = { 1, 1, 2, 3, 4, 4 };
74     int             expected[] = { 1, 2, 3, 4 };
75     gmx_ana_index_t g          = initGroup(index);
76     gmx_ana_index_t e          = initGroup(expected);
77     gmx_ana_index_remove_duplicates(&g);
78     EXPECT_TRUE(gmx_ana_index_equals(&g, &e));
79 }
80
81 //! Text fixture for index block operations
82 class IndexBlockTest : public ::testing::Test
83 {
84 public:
85     IndexBlockTest();
86     ~IndexBlockTest() override;
87
88     //@{
89     //! Set the input group for the index block operation
90     void setGroup(int count, const int atoms[]);
91     template<int count>
92     void setGroup(const int (&atoms)[count])
93     {
94         setGroup(count, atoms);
95     }
96     //@}
97     /*! \brief Check the input group and output with refdata, with
98      * an optional \c id to name the refdata block */
99     void checkBlocka(const char* id = "Block");
100     //! Make a simple topology to check with
101     void makeSimpleTopology();
102     //! Make a complex topology to check with
103     void makeComplexTopology();
104
105     //! Managers reference data for the tests
106     gmx::test::TestReferenceData data_;
107     //! Manages setting up a topology and matching data structures
108     gmx::test::TopologyManager topManager_;
109     //! The input group to test with
110     gmx_ana_index_t g_;
111     //! The output block to test on
112     t_blocka blocka_;
113 };
114
115 IndexBlockTest::IndexBlockTest()
116 {
117     blocka_.nr           = 0;
118     blocka_.index        = nullptr;
119     blocka_.nalloc_index = 0;
120     blocka_.nra          = 0;
121     blocka_.a            = nullptr;
122     blocka_.nalloc_a     = 0;
123     gmx_ana_index_clear(&g_);
124 }
125
126 IndexBlockTest::~IndexBlockTest()
127 {
128     done_blocka(&blocka_);
129 }
130
131 void IndexBlockTest::setGroup(int count, const int atoms[])
132 {
133     g_.isize = count;
134     g_.index = const_cast<int*>(atoms);
135 }
136
137 void IndexBlockTest::checkBlocka(const char* id)
138 {
139     gmx::test::TestReferenceChecker compound(data_.rootChecker().checkCompound("BlockAtoms", id));
140     compound.checkSequenceArray(g_.isize, g_.index, "Input");
141     compound.checkInteger(blocka_.nr, "Count");
142     for (int i = 0; i < blocka_.nr; ++i)
143     {
144         gmx::test::TestReferenceChecker blockCompound(compound.checkCompound("Block", nullptr));
145         blockCompound.checkSequence(&blocka_.a[blocka_.index[i]], &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[]  = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
320     const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
321     topManager_.initAtoms(18);
322     topManager_.initUniformResidues(3);
323     setGroup(maxGroup);
324     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
325     setGroup(testGroup);
326     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
327 }
328
329 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
330 {
331     const int maxGroup[]  = { 15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6 };
332     const int testGroup[] = {
333         2, 1, 0, 5, 4, 3, 8, 7, 6,
334     };
335     topManager_.initAtoms(18);
336     topManager_.initUniformResidues(3);
337     setGroup(maxGroup);
338     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
339     setGroup(testGroup);
340     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
341 }
342
343 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
344 {
345     const int maxGroup[]   = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
346     const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
347     const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
348     const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
349
350     topManager_.initAtoms(18);
351     topManager_.initUniformResidues(3);
352     setGroup(maxGroup);
353     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
354
355     setGroup(testGroup1);
356     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
357
358     setGroup(testGroup2);
359     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
360
361     setGroup(testGroup3);
362     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
363 }
364
365 /********************************************************************
366  * gmx_ana_index_has_complete_elems() tests
367  */
368
369 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
370 {
371     const int group[] = { 0, 1, 2 };
372     setGroup(group);
373     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, nullptr));
374     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, nullptr));
375     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, nullptr));
376 }
377
378 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
379 {
380     const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
381     const int group2[] = { 3, 4, 5, 6, 7, 8 };
382
383     topManager_.initAtoms(15);
384     topManager_.initUniformResidues(3);
385     gmx_mtop_t* top = topManager_.topology();
386
387     setGroup(group1);
388     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
389
390     setGroup(group2);
391     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
392 }
393
394 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
395 {
396     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
397     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
398     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
399     const int group4[] = { 3, 4, 5, 6, 8, 12, 13, 14 };
400
401     topManager_.initAtoms(18);
402     topManager_.initUniformResidues(3);
403     gmx_mtop_t* top = topManager_.topology();
404
405     setGroup(group1);
406     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
407
408     setGroup(group2);
409     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
410
411     setGroup(group3);
412     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
413
414     setGroup(group4);
415     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
416 }
417
418 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
419 {
420     const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
421
422     topManager_.initAtoms(15);
423     topManager_.initUniformResidues(1);
424     topManager_.initUniformMolecules(3);
425     gmx_mtop_t* top = topManager_.topology();
426
427     setGroup(group);
428     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
429 }
430
431 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
432 {
433     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
434     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
435     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
436
437     topManager_.initAtoms(18);
438     topManager_.initUniformResidues(1);
439     topManager_.initUniformMolecules(3);
440     gmx_mtop_t* top = topManager_.topology();
441
442     setGroup(group1);
443     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
444
445     setGroup(group2);
446     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
447
448     setGroup(group3);
449     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
450 }
451
452 /********************************************************************
453  * IndexMapTest
454  */
455
456 class IndexMapTest : public ::testing::Test
457 {
458 public:
459     IndexMapTest();
460     ~IndexMapTest() override;
461
462     void testInit(int atomCount, const int atoms[], e_index_t type);
463     void testUpdate(int atomCount, const int atoms[], bool bMaskOnly, const char* name);
464     void testOrgIdGroup(e_index_t type, const char* name);
465     template<int count>
466     void testInit(const int (&atoms)[count], e_index_t type)
467     {
468         testInit(count, atoms, type);
469     }
470     template<int count>
471     void testUpdate(const int (&atoms)[count], bool bMaskOnly, const char* name)
472     {
473         testUpdate(count, atoms, bMaskOnly, name);
474     }
475
476     void checkMapping(int atomCount, const int atoms[], const char* name);
477
478     gmx::test::TestReferenceData    data_;
479     gmx::test::TestReferenceChecker checker_;
480     gmx::test::TopologyManager      topManager_;
481     gmx_ana_indexmap_t              map_;
482
483 private:
484     gmx_ana_index_t initGroup_;
485 };
486
487 IndexMapTest::IndexMapTest() : checker_(data_.rootChecker())
488 {
489     gmx_ana_indexmap_clear(&map_);
490     gmx_ana_index_clear(&initGroup_);
491 }
492
493 IndexMapTest::~IndexMapTest()
494 {
495     gmx_ana_indexmap_deinit(&map_);
496 }
497
498 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
499 {
500     initGroup_.isize = atomCount;
501     initGroup_.index = const_cast<int*>(atoms);
502     gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
503     EXPECT_EQ(type, map_.type);
504     checkMapping(atomCount, atoms, "Initialized");
505 }
506
507 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly, const char* name)
508 {
509     gmx_ana_index_t g;
510     g.isize = atomCount;
511     g.index = const_cast<int*>(atoms);
512     gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
513     if (name == nullptr)
514     {
515         name = "Updated";
516     }
517     if (bMaskOnly)
518     {
519         checkMapping(initGroup_.isize, initGroup_.index, name);
520     }
521     else
522     {
523         checkMapping(atomCount, atoms, name);
524     }
525 }
526
527 void IndexMapTest::testOrgIdGroup(e_index_t type, const char* name)
528 {
529     gmx::test::TestReferenceChecker compound(checker_.checkCompound("OrgIdGroups", name));
530     const int count = gmx_ana_indexmap_init_orgid_group(&map_, topManager_.topology(), type);
531     compound.checkInteger(count, "GroupCount");
532     compound.checkSequenceArray(map_.mapb.nr, map_.orgid, "OrgId");
533     for (int i = 0; i < map_.mapb.nr; ++i)
534     {
535         EXPECT_EQ(map_.orgid[i], map_.mapid[i]);
536     }
537 }
538
539 void IndexMapTest::checkMapping(int atomCount, const int atoms[], const char* name)
540 {
541     gmx::test::TestReferenceChecker compound(checker_.checkCompound("IndexMapping", name));
542     compound.checkSequenceArray(atomCount, atoms, "Input");
543     compound.checkInteger(map_.mapb.nr, "Count");
544     for (int i = 0; i < map_.mapb.nr; ++i)
545     {
546         gmx::test::TestReferenceChecker blockCompound(compound.checkCompound("Block", nullptr));
547         blockCompound.checkSequence(&atoms[map_.mapb.index[i]], &atoms[map_.mapb.index[i + 1]],
548                                     "Atoms");
549         blockCompound.checkInteger(map_.refid[i], "RefId");
550         blockCompound.checkInteger(map_.mapid[i], "MapId");
551         int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
552         EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
553     }
554 }
555
556 /********************************************************************
557  * gmx_ana_indexmap_t tests
558  */
559
560 TEST_F(IndexMapTest, InitializesAtomBlock)
561 {
562     const int maxGroup[] = { 1, 2, 4, 5 };
563     testInit(maxGroup, INDEX_ATOM);
564 }
565
566 TEST_F(IndexMapTest, InitializesOrgIdGroupAtom)
567 {
568     const int maxGroup[] = { 2, 5, 7 };
569     testInit(maxGroup, INDEX_ATOM);
570     testOrgIdGroup(INDEX_ATOM, "Atoms");
571 }
572
573 TEST_F(IndexMapTest, InitializesOrgIdGroupSingle)
574 {
575     const int maxGroup[] = { 3, 4, 7, 8, 13 };
576     topManager_.initAtoms(18);
577     topManager_.initUniformResidues(3);
578     testInit(maxGroup, INDEX_RES);
579     testOrgIdGroup(INDEX_ATOM, "Single");
580 }
581
582 TEST_F(IndexMapTest, InitializesOrgIdGroupResidue)
583 {
584     const int maxGroup[] = { 3, 4, 7, 8, 13 };
585     topManager_.initAtoms(18);
586     topManager_.initUniformResidues(3);
587     testInit(maxGroup, INDEX_ATOM);
588     testOrgIdGroup(INDEX_RES, "Residues");
589 }
590
591 TEST_F(IndexMapTest, InitializesOrgIdGroupMolecule)
592 {
593     const int maxGroup[] = { 1, 2, 3, 4, 7, 8, 13 };
594     topManager_.initAtoms(18);
595     topManager_.initUniformResidues(3);
596     topManager_.initUniformMolecules(6);
597     testInit(maxGroup, INDEX_RES);
598     testOrgIdGroup(INDEX_MOL, "Molecules");
599 }
600
601 TEST_F(IndexMapTest, InitializesOrgIdGroupAll)
602 {
603     const int maxGroup[] = { 3, 4, 7, 8, 13 };
604     testInit(maxGroup, INDEX_ATOM);
605     testOrgIdGroup(INDEX_ALL, "All");
606 }
607
608 TEST_F(IndexMapTest, InitializesMoleculeBlock)
609 {
610     const int maxGroup[] = { 3, 4, 7, 8, 13 };
611     topManager_.initAtoms(18);
612     topManager_.initUniformResidues(1);
613     topManager_.initUniformMolecules(3);
614     testInit(maxGroup, INDEX_MOL);
615 }
616
617 TEST_F(IndexMapTest, MapsSingleBlock)
618 {
619     const int maxGroup[]  = { 0, 1, 2, 3 };
620     const int evalGroup[] = { 0, 2 };
621     testInit(maxGroup, INDEX_ALL);
622     testUpdate(evalGroup, false, nullptr);
623 }
624
625 TEST_F(IndexMapTest, MapsResidueBlocks)
626 {
627     const int maxGroup[]  = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
628     const int evalGroup[] = { 3, 4, 7, 8, 13 };
629     topManager_.initAtoms(18);
630     topManager_.initUniformResidues(3);
631     testInit(maxGroup, INDEX_RES);
632     testUpdate(evalGroup, false, nullptr);
633 }
634
635 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
636 {
637     const int maxGroup[]  = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
638     const int evalGroup[] = { 3, 4, 7, 8, 13 };
639     topManager_.initAtoms(18);
640     topManager_.initUniformResidues(3);
641     testInit(maxGroup, INDEX_RES);
642     testUpdate(evalGroup, true, nullptr);
643 }
644
645 TEST_F(IndexMapTest, HandlesMultipleRequests)
646 {
647     const int maxGroup[]  = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
648     const int evalGroup[] = { 3, 4, 7, 8, 13 };
649     topManager_.initAtoms(18);
650     topManager_.initUniformResidues(3);
651     testInit(maxGroup, INDEX_RES);
652     testUpdate(evalGroup, false, "EvaluatedNoMask");
653     testUpdate(evalGroup, true, "EvaluatedMask");
654     testUpdate(maxGroup, true, "Initialized");
655     testUpdate(evalGroup, true, "EvaluatedMask");
656     testUpdate(evalGroup, false, "EvaluatedNoMask");
657     testUpdate(maxGroup, false, "Initialized");
658 }
659
660 /***********************************************************************
661  * IndexGroupsAndNames tests
662  */
663
664 class IndexGroupsAndNamesTest : public ::testing::Test
665 {
666 public:
667     IndexGroupsAndNamesTest()
668     {
669         init_blocka(&blockA_);
670         addGroupToBlocka_(indicesGroupA_);
671         addGroupToBlocka_(indicesGroupB_);
672         addGroupToBlocka_(indicesGroupSecondA_);
673         addGroupToBlocka_(indicesGroupC_);
674
675         const char* const namesAsConstCharArray[4] = { groupNames[0].c_str(), groupNames[1].c_str(),
676                                                        groupNames[2].c_str(), groupNames[3].c_str() };
677         indexGroupAndNames_ = std::make_unique<gmx::IndexGroupsAndNames>(blockA_, namesAsConstCharArray);
678     }
679     ~IndexGroupsAndNamesTest() override { done_blocka(&blockA_); }
680
681 protected:
682     std::unique_ptr<gmx::IndexGroupsAndNames> indexGroupAndNames_;
683     const std::vector<std::string>            groupNames           = { "A", "B - Name", "A", "C" };
684     const std::vector<gmx::index>             indicesGroupA_       = {};
685     const std::vector<gmx::index>             indicesGroupB_       = { 1, 2 };
686     const std::vector<gmx::index>             indicesGroupSecondA_ = { 5 };
687     const std::vector<gmx::index>             indicesGroupC_       = { 10 };
688
689 private:
690     //! Add a new group to t_blocka
691     void addGroupToBlocka_(gmx::ArrayRef<const gmx::index> index)
692     {
693         srenew(blockA_.index, blockA_.nr + 2);
694         srenew(blockA_.a, blockA_.nra + index.size());
695         for (int i = 0; (i < index.ssize()); i++)
696         {
697             blockA_.a[blockA_.nra++] = index[i];
698         }
699         blockA_.nr++;
700         blockA_.index[blockA_.nr] = blockA_.nra;
701     }
702
703     t_blocka blockA_;
704 };
705
706 TEST_F(IndexGroupsAndNamesTest, containsNames)
707 {
708     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("a"));
709     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("A"));
710     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("B - Name"));
711     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("b - Name"));
712     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("B - naMe"));
713     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("C"));
714     EXPECT_FALSE(indexGroupAndNames_->containsGroupName("D"));
715 }
716
717 TEST_F(IndexGroupsAndNamesTest, throwsWhenNameMissing)
718 {
719     EXPECT_ANY_THROW(indexGroupAndNames_->indices("D"));
720 }
721
722 TEST_F(IndexGroupsAndNamesTest, groupIndicesCorrect)
723 {
724     using ::testing::Eq;
725     using ::testing::Pointwise;
726     EXPECT_THAT(indicesGroupA_, Pointwise(Eq(), indexGroupAndNames_->indices("A")));
727     EXPECT_THAT(indicesGroupB_, Pointwise(Eq(), indexGroupAndNames_->indices("B - name")));
728     EXPECT_THAT(indicesGroupC_, Pointwise(Eq(), indexGroupAndNames_->indices("C")));
729 }
730
731
732 } // namespace