Split lines with many copyright years
[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 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief
38  * Tests the index group handling in the selection engine.
39  *
40  * \todo
41  * Tests for other functions, at least the set operations.
42  *
43  * \author Teemu Murtola <teemu.murtola@gmail.com>
44  * \ingroup module_selection
45  */
46 #include "gmxpre.h"
47
48 #include "gromacs/selection/indexutil.h"
49
50 #include <gmock/gmock.h>
51 #include <gtest/gtest.h>
52
53 #include "gromacs/topology/block.h"
54 #include "gromacs/utility/arrayref.h"
55 #include "gromacs/utility/smalloc.h"
56
57 #include "testutils/refdata.h"
58 #include "testutils/testasserts.h"
59
60 #include "toputils.h"
61
62 namespace
63 {
64
65 //! Helper for creating groups from an array.
66 gmx_ana_index_t initGroup(gmx::ArrayRef<int> index)
67 {
68     gmx_ana_index_t g = { static_cast<int>(index.size()), index.data(), 0 };
69     return g;
70 }
71
72 TEST(IndexGroupTest, RemovesDuplicates)
73 {
74     int             index[]    = { 1, 1, 2, 3, 4, 4 };
75     int             expected[] = { 1, 2, 3, 4 };
76     gmx_ana_index_t g          = initGroup(index);
77     gmx_ana_index_t e          = initGroup(expected);
78     gmx_ana_index_remove_duplicates(&g);
79     EXPECT_TRUE(gmx_ana_index_equals(&g, &e));
80 }
81
82 //! Text fixture for index block operations
83 class IndexBlockTest : public ::testing::Test
84 {
85 public:
86     IndexBlockTest();
87     ~IndexBlockTest() override;
88
89     //@{
90     //! Set the input group for the index block operation
91     void setGroup(int count, const int atoms[]);
92     template<int count>
93     void setGroup(const int (&atoms)[count])
94     {
95         setGroup(count, atoms);
96     }
97     //@}
98     /*! \brief Check the input group and output with refdata, with
99      * an optional \c id to name the refdata block */
100     void checkBlocka(const char* id = "Block");
101     //! Make a simple topology to check with
102     void makeSimpleTopology();
103     //! Make a complex topology to check with
104     void makeComplexTopology();
105
106     //! Managers reference data for the tests
107     gmx::test::TestReferenceData data_;
108     //! Manages setting up a topology and matching data structures
109     gmx::test::TopologyManager topManager_;
110     //! The input group to test with
111     gmx_ana_index_t g_;
112     //! The output block to test on
113     t_blocka blocka_;
114 };
115
116 IndexBlockTest::IndexBlockTest()
117 {
118     blocka_.nr           = 0;
119     blocka_.index        = nullptr;
120     blocka_.nalloc_index = 0;
121     blocka_.nra          = 0;
122     blocka_.a            = nullptr;
123     blocka_.nalloc_a     = 0;
124     gmx_ana_index_clear(&g_);
125 }
126
127 IndexBlockTest::~IndexBlockTest()
128 {
129     done_blocka(&blocka_);
130 }
131
132 void IndexBlockTest::setGroup(int count, const int atoms[])
133 {
134     g_.isize = count;
135     g_.index = const_cast<int*>(atoms);
136 }
137
138 void IndexBlockTest::checkBlocka(const char* id)
139 {
140     gmx::test::TestReferenceChecker compound(data_.rootChecker().checkCompound("BlockAtoms", id));
141     compound.checkSequenceArray(g_.isize, g_.index, "Input");
142     compound.checkInteger(blocka_.nr, "Count");
143     for (int i = 0; i < blocka_.nr; ++i)
144     {
145         gmx::test::TestReferenceChecker blockCompound(compound.checkCompound("Block", nullptr));
146         blockCompound.checkSequence(&blocka_.a[blocka_.index[i]], &blocka_.a[blocka_.index[i + 1]],
147                                     "Atoms");
148     }
149 }
150
151 void IndexBlockTest::makeSimpleTopology()
152 {
153     topManager_.initTopology(1, 1);
154     {
155         int              moleculeTypeIndex   = 0;
156         std::vector<int> numAtomsInResidues  = { 3, 3, 3 };
157         int              moleculeBlockIndex  = 0;
158         int              numMoleculesInBlock = 1;
159         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
160         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
161     }
162     topManager_.finalizeTopology();
163 }
164
165 void IndexBlockTest::makeComplexTopology()
166 {
167     topManager_.initTopology(3, 4);
168     {
169         int              moleculeTypeIndex   = 0;
170         std::vector<int> numAtomsInResidues  = { 2, 2, 1 };
171         int              moleculeBlockIndex  = 0;
172         int              numMoleculesInBlock = 1;
173         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
174         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
175     }
176     {
177         int              moleculeTypeIndex   = 1;
178         std::vector<int> numAtomsInResidues  = { 1 };
179         int              moleculeBlockIndex  = 1;
180         int              numMoleculesInBlock = 3;
181         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
182         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
183     }
184     {
185         int              moleculeTypeIndex   = 2;
186         std::vector<int> numAtomsInResidues  = { 3 };
187         int              moleculeBlockIndex  = 2;
188         int              numMoleculesInBlock = 1;
189         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
190         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
191     }
192     {
193         int moleculeTypeIndex   = 0; // Re-using a moltype definition
194         int moleculeBlockIndex  = 3;
195         int numMoleculesInBlock = 1;
196         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
197     }
198     topManager_.finalizeTopology();
199 }
200
201 /********************************************************************
202  * gmx_ana_index_make_block() tests
203  */
204
205 TEST_F(IndexBlockTest, CreatesUnknownBlock)
206 {
207     gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
208     checkBlocka();
209     done_blocka(&blocka_);
210     gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
211     checkBlocka();
212 }
213
214 TEST_F(IndexBlockTest, CreatesAtomBlock)
215 {
216     const int group[] = { 0, 1, 3, 4, 6 };
217     setGroup(group);
218     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, false);
219     checkBlocka();
220     done_blocka(&blocka_);
221     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, true);
222     checkBlocka();
223 }
224
225 TEST_F(IndexBlockTest, CreatesResidueBlocksForSimpleTopology)
226 {
227     makeSimpleTopology();
228
229     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8 };
230     setGroup(group);
231     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
232     checkBlocka("FromAllAtoms");
233     done_blocka(&blocka_);
234     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
235     checkBlocka("FromAllAtoms");
236 }
237
238 TEST_F(IndexBlockTest, CreatesResidueBlocksForComplexTopology)
239 {
240     makeComplexTopology();
241
242     SCOPED_TRACE("Group with all atoms without completion");
243     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15 };
244     setGroup(group);
245     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
246     checkBlocka("FromAllAtoms");
247     done_blocka(&blocka_);
248     // SCOPED_TRACE("Group with all atoms with completion");
249     // gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
250     // checkBlocka("FromAllAtoms");
251     // done_blocka(&blocka_);
252
253     SCOPED_TRACE("Group with some atoms without completion");
254     const int subgroup[] = { 0, 3, 4, 5, 6, 7, 8, 12, 13, 15 };
255     setGroup(subgroup);
256     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
257     checkBlocka("FromSomeAtomsWithoutCompletion");
258     // done_blocka(&blocka_);
259     // SCOPED_TRACE("Group with some atoms with completion");
260     // gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
261     // checkBlocka("FromSomeAtomsWithCompletion");
262 }
263
264 TEST_F(IndexBlockTest, CreatesMoleculeBlocksForSimpleTopology)
265 {
266     makeSimpleTopology();
267
268     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8 };
269     setGroup(group);
270     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
271     checkBlocka("FromAllAtoms");
272     done_blocka(&blocka_);
273     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
274     checkBlocka("FromAllAtoms");
275 }
276
277 TEST_F(IndexBlockTest, CreatesMoleculeBlocksForComplexTopology)
278 {
279     makeComplexTopology();
280
281     SCOPED_TRACE("Group with all atoms without completion");
282     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15 };
283     setGroup(group);
284     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
285     checkBlocka("FromAllAtoms");
286     done_blocka(&blocka_);
287     // SCOPED_TRACE("Group with all atoms with completion");
288     // gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
289     // checkBlocka("FromAllAtoms");
290     // done_blocka(&blocka_);
291
292     SCOPED_TRACE("Group with some atoms without completion");
293     const int subgroup[] = { 1, 5, 6, 7, 11, 12 };
294     setGroup(subgroup);
295     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
296     checkBlocka("FromSomeAtomsWithoutCompletion");
297     // done_blocka(&blocka_);
298     // SCOPED_TRACE("Group with some atoms with completion");
299     // gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
300     // checkBlocka("FromSomeAtomsWithCompletion");
301 }
302
303 TEST_F(IndexBlockTest, CreatesSingleBlock)
304 {
305     const int group[] = { 0, 1, 3, 4, 6 };
306     setGroup(group);
307     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, false);
308     checkBlocka();
309     done_blocka(&blocka_);
310     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, true);
311     checkBlocka();
312 }
313
314 /********************************************************************
315  * gmx_ana_index_has_full_ablocks() tests
316  */
317
318 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
319 {
320     const int maxGroup[]  = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
321     const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
322     topManager_.initAtoms(18);
323     topManager_.initUniformResidues(3);
324     setGroup(maxGroup);
325     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
326     setGroup(testGroup);
327     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
328 }
329
330 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
331 {
332     const int maxGroup[]  = { 15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6 };
333     const int testGroup[] = {
334         2, 1, 0, 5, 4, 3, 8, 7, 6,
335     };
336     topManager_.initAtoms(18);
337     topManager_.initUniformResidues(3);
338     setGroup(maxGroup);
339     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
340     setGroup(testGroup);
341     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
342 }
343
344 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
345 {
346     const int maxGroup[]   = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
347     const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
348     const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
349     const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
350
351     topManager_.initAtoms(18);
352     topManager_.initUniformResidues(3);
353     setGroup(maxGroup);
354     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
355
356     setGroup(testGroup1);
357     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
358
359     setGroup(testGroup2);
360     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
361
362     setGroup(testGroup3);
363     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
364 }
365
366 /********************************************************************
367  * gmx_ana_index_has_complete_elems() tests
368  */
369
370 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
371 {
372     const int group[] = { 0, 1, 2 };
373     setGroup(group);
374     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, nullptr));
375     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, nullptr));
376     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, nullptr));
377 }
378
379 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
380 {
381     const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
382     const int group2[] = { 3, 4, 5, 6, 7, 8 };
383
384     topManager_.initAtoms(15);
385     topManager_.initUniformResidues(3);
386     gmx_mtop_t* top = topManager_.topology();
387
388     setGroup(group1);
389     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
390
391     setGroup(group2);
392     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
393 }
394
395 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
396 {
397     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
398     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
399     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
400     const int group4[] = { 3, 4, 5, 6, 8, 12, 13, 14 };
401
402     topManager_.initAtoms(18);
403     topManager_.initUniformResidues(3);
404     gmx_mtop_t* top = topManager_.topology();
405
406     setGroup(group1);
407     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
408
409     setGroup(group2);
410     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
411
412     setGroup(group3);
413     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
414
415     setGroup(group4);
416     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
417 }
418
419 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
420 {
421     const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
422
423     topManager_.initAtoms(15);
424     topManager_.initUniformResidues(1);
425     topManager_.initUniformMolecules(3);
426     gmx_mtop_t* top = topManager_.topology();
427
428     setGroup(group);
429     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
430 }
431
432 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
433 {
434     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
435     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
436     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
437
438     topManager_.initAtoms(18);
439     topManager_.initUniformResidues(1);
440     topManager_.initUniformMolecules(3);
441     gmx_mtop_t* top = topManager_.topology();
442
443     setGroup(group1);
444     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
445
446     setGroup(group2);
447     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
448
449     setGroup(group3);
450     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
451 }
452
453 /********************************************************************
454  * IndexMapTest
455  */
456
457 class IndexMapTest : public ::testing::Test
458 {
459 public:
460     IndexMapTest();
461     ~IndexMapTest() override;
462
463     void testInit(int atomCount, const int atoms[], e_index_t type);
464     void testUpdate(int atomCount, const int atoms[], bool bMaskOnly, const char* name);
465     void testOrgIdGroup(e_index_t type, const char* name);
466     template<int count>
467     void testInit(const int (&atoms)[count], e_index_t type)
468     {
469         testInit(count, atoms, type);
470     }
471     template<int count>
472     void testUpdate(const int (&atoms)[count], bool bMaskOnly, const char* name)
473     {
474         testUpdate(count, atoms, bMaskOnly, name);
475     }
476
477     void checkMapping(int atomCount, const int atoms[], const char* name);
478
479     gmx::test::TestReferenceData    data_;
480     gmx::test::TestReferenceChecker checker_;
481     gmx::test::TopologyManager      topManager_;
482     gmx_ana_indexmap_t              map_;
483
484 private:
485     gmx_ana_index_t initGroup_;
486 };
487
488 IndexMapTest::IndexMapTest() : checker_(data_.rootChecker())
489 {
490     gmx_ana_indexmap_clear(&map_);
491     gmx_ana_index_clear(&initGroup_);
492 }
493
494 IndexMapTest::~IndexMapTest()
495 {
496     gmx_ana_indexmap_deinit(&map_);
497 }
498
499 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
500 {
501     initGroup_.isize = atomCount;
502     initGroup_.index = const_cast<int*>(atoms);
503     gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
504     EXPECT_EQ(type, map_.type);
505     checkMapping(atomCount, atoms, "Initialized");
506 }
507
508 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly, const char* name)
509 {
510     gmx_ana_index_t g;
511     g.isize = atomCount;
512     g.index = const_cast<int*>(atoms);
513     gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
514     if (name == nullptr)
515     {
516         name = "Updated";
517     }
518     if (bMaskOnly)
519     {
520         checkMapping(initGroup_.isize, initGroup_.index, name);
521     }
522     else
523     {
524         checkMapping(atomCount, atoms, name);
525     }
526 }
527
528 void IndexMapTest::testOrgIdGroup(e_index_t type, const char* name)
529 {
530     gmx::test::TestReferenceChecker compound(checker_.checkCompound("OrgIdGroups", name));
531     const int count = gmx_ana_indexmap_init_orgid_group(&map_, topManager_.topology(), type);
532     compound.checkInteger(count, "GroupCount");
533     compound.checkSequenceArray(map_.mapb.nr, map_.orgid, "OrgId");
534     for (int i = 0; i < map_.mapb.nr; ++i)
535     {
536         EXPECT_EQ(map_.orgid[i], map_.mapid[i]);
537     }
538 }
539
540 void IndexMapTest::checkMapping(int atomCount, const int atoms[], const char* name)
541 {
542     gmx::test::TestReferenceChecker compound(checker_.checkCompound("IndexMapping", name));
543     compound.checkSequenceArray(atomCount, atoms, "Input");
544     compound.checkInteger(map_.mapb.nr, "Count");
545     for (int i = 0; i < map_.mapb.nr; ++i)
546     {
547         gmx::test::TestReferenceChecker blockCompound(compound.checkCompound("Block", nullptr));
548         blockCompound.checkSequence(&atoms[map_.mapb.index[i]], &atoms[map_.mapb.index[i + 1]],
549                                     "Atoms");
550         blockCompound.checkInteger(map_.refid[i], "RefId");
551         blockCompound.checkInteger(map_.mapid[i], "MapId");
552         int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
553         EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
554     }
555 }
556
557 /********************************************************************
558  * gmx_ana_indexmap_t tests
559  */
560
561 TEST_F(IndexMapTest, InitializesAtomBlock)
562 {
563     const int maxGroup[] = { 1, 2, 4, 5 };
564     testInit(maxGroup, INDEX_ATOM);
565 }
566
567 TEST_F(IndexMapTest, InitializesOrgIdGroupAtom)
568 {
569     const int maxGroup[] = { 2, 5, 7 };
570     testInit(maxGroup, INDEX_ATOM);
571     testOrgIdGroup(INDEX_ATOM, "Atoms");
572 }
573
574 TEST_F(IndexMapTest, InitializesOrgIdGroupSingle)
575 {
576     const int maxGroup[] = { 3, 4, 7, 8, 13 };
577     topManager_.initAtoms(18);
578     topManager_.initUniformResidues(3);
579     testInit(maxGroup, INDEX_RES);
580     testOrgIdGroup(INDEX_ATOM, "Single");
581 }
582
583 TEST_F(IndexMapTest, InitializesOrgIdGroupResidue)
584 {
585     const int maxGroup[] = { 3, 4, 7, 8, 13 };
586     topManager_.initAtoms(18);
587     topManager_.initUniformResidues(3);
588     testInit(maxGroup, INDEX_ATOM);
589     testOrgIdGroup(INDEX_RES, "Residues");
590 }
591
592 TEST_F(IndexMapTest, InitializesOrgIdGroupMolecule)
593 {
594     const int maxGroup[] = { 1, 2, 3, 4, 7, 8, 13 };
595     topManager_.initAtoms(18);
596     topManager_.initUniformResidues(3);
597     topManager_.initUniformMolecules(6);
598     testInit(maxGroup, INDEX_RES);
599     testOrgIdGroup(INDEX_MOL, "Molecules");
600 }
601
602 TEST_F(IndexMapTest, InitializesOrgIdGroupAll)
603 {
604     const int maxGroup[] = { 3, 4, 7, 8, 13 };
605     testInit(maxGroup, INDEX_ATOM);
606     testOrgIdGroup(INDEX_ALL, "All");
607 }
608
609 TEST_F(IndexMapTest, InitializesMoleculeBlock)
610 {
611     const int maxGroup[] = { 3, 4, 7, 8, 13 };
612     topManager_.initAtoms(18);
613     topManager_.initUniformResidues(1);
614     topManager_.initUniformMolecules(3);
615     testInit(maxGroup, INDEX_MOL);
616 }
617
618 TEST_F(IndexMapTest, MapsSingleBlock)
619 {
620     const int maxGroup[]  = { 0, 1, 2, 3 };
621     const int evalGroup[] = { 0, 2 };
622     testInit(maxGroup, INDEX_ALL);
623     testUpdate(evalGroup, false, nullptr);
624 }
625
626 TEST_F(IndexMapTest, MapsResidueBlocks)
627 {
628     const int maxGroup[]  = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
629     const int evalGroup[] = { 3, 4, 7, 8, 13 };
630     topManager_.initAtoms(18);
631     topManager_.initUniformResidues(3);
632     testInit(maxGroup, INDEX_RES);
633     testUpdate(evalGroup, false, nullptr);
634 }
635
636 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
637 {
638     const int maxGroup[]  = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
639     const int evalGroup[] = { 3, 4, 7, 8, 13 };
640     topManager_.initAtoms(18);
641     topManager_.initUniformResidues(3);
642     testInit(maxGroup, INDEX_RES);
643     testUpdate(evalGroup, true, nullptr);
644 }
645
646 TEST_F(IndexMapTest, HandlesMultipleRequests)
647 {
648     const int maxGroup[]  = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17 };
649     const int evalGroup[] = { 3, 4, 7, 8, 13 };
650     topManager_.initAtoms(18);
651     topManager_.initUniformResidues(3);
652     testInit(maxGroup, INDEX_RES);
653     testUpdate(evalGroup, false, "EvaluatedNoMask");
654     testUpdate(evalGroup, true, "EvaluatedMask");
655     testUpdate(maxGroup, true, "Initialized");
656     testUpdate(evalGroup, true, "EvaluatedMask");
657     testUpdate(evalGroup, false, "EvaluatedNoMask");
658     testUpdate(maxGroup, false, "Initialized");
659 }
660
661 /***********************************************************************
662  * IndexGroupsAndNames tests
663  */
664
665 class IndexGroupsAndNamesTest : public ::testing::Test
666 {
667 public:
668     IndexGroupsAndNamesTest()
669     {
670         init_blocka(&blockA_);
671         addGroupToBlocka_(indicesGroupA_);
672         addGroupToBlocka_(indicesGroupB_);
673         addGroupToBlocka_(indicesGroupSecondA_);
674         addGroupToBlocka_(indicesGroupC_);
675
676         const char* const namesAsConstCharArray[4] = { groupNames[0].c_str(), groupNames[1].c_str(),
677                                                        groupNames[2].c_str(), groupNames[3].c_str() };
678         indexGroupAndNames_ = std::make_unique<gmx::IndexGroupsAndNames>(blockA_, namesAsConstCharArray);
679     }
680     ~IndexGroupsAndNamesTest() override { done_blocka(&blockA_); }
681
682 protected:
683     std::unique_ptr<gmx::IndexGroupsAndNames> indexGroupAndNames_;
684     const std::vector<std::string>            groupNames           = { "A", "B - Name", "A", "C" };
685     const std::vector<gmx::index>             indicesGroupA_       = {};
686     const std::vector<gmx::index>             indicesGroupB_       = { 1, 2 };
687     const std::vector<gmx::index>             indicesGroupSecondA_ = { 5 };
688     const std::vector<gmx::index>             indicesGroupC_       = { 10 };
689
690 private:
691     //! Add a new group to t_blocka
692     void addGroupToBlocka_(gmx::ArrayRef<const gmx::index> index)
693     {
694         srenew(blockA_.index, blockA_.nr + 2);
695         srenew(blockA_.a, blockA_.nra + index.size());
696         for (int i = 0; (i < index.ssize()); i++)
697         {
698             blockA_.a[blockA_.nra++] = index[i];
699         }
700         blockA_.nr++;
701         blockA_.index[blockA_.nr] = blockA_.nra;
702     }
703
704     t_blocka blockA_;
705 };
706
707 TEST_F(IndexGroupsAndNamesTest, containsNames)
708 {
709     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("a"));
710     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("A"));
711     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("B - Name"));
712     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("b - Name"));
713     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("B - naMe"));
714     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("C"));
715     EXPECT_FALSE(indexGroupAndNames_->containsGroupName("D"));
716 }
717
718 TEST_F(IndexGroupsAndNamesTest, throwsWhenNameMissing)
719 {
720     EXPECT_ANY_THROW(indexGroupAndNames_->indices("D"));
721 }
722
723 TEST_F(IndexGroupsAndNamesTest, groupIndicesCorrect)
724 {
725     using ::testing::Eq;
726     using ::testing::Pointwise;
727     EXPECT_THAT(indicesGroupA_, Pointwise(Eq(), indexGroupAndNames_->indices("A")));
728     EXPECT_THAT(indicesGroupB_, Pointwise(Eq(), indexGroupAndNames_->indices("B - name")));
729     EXPECT_THAT(indicesGroupC_, Pointwise(Eq(), indexGroupAndNames_->indices("C")));
730 }
731
732
733 } // namespace