Index groups for MdModules from grompp
[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(
140             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(
146                 compound.checkCompound("Block", nullptr));
147         blockCompound.checkSequence(&blocka_.a[blocka_.index[i]],
148                                     &blocka_.a[blocka_.index[i+1]],
149                                     "Atoms");
150     }
151 }
152
153 void IndexBlockTest::makeSimpleTopology()
154 {
155     topManager_.initTopology(1, 1);
156     {
157         int              moleculeTypeIndex   = 0;
158         std::vector<int> numAtomsInResidues  = {3, 3, 3};
159         int              moleculeBlockIndex  = 0;
160         int              numMoleculesInBlock = 1;
161         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
162         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
163     }
164     topManager_.finalizeTopology();
165 }
166
167 void IndexBlockTest::makeComplexTopology()
168 {
169     topManager_.initTopology(3, 4);
170     {
171         int              moleculeTypeIndex   = 0;
172         std::vector<int> numAtomsInResidues  = {2, 2, 1};
173         int              moleculeBlockIndex  = 0;
174         int              numMoleculesInBlock = 1;
175         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
176         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
177     }
178     {
179         int              moleculeTypeIndex   = 1;
180         std::vector<int> numAtomsInResidues  = {1};
181         int              moleculeBlockIndex  = 1;
182         int              numMoleculesInBlock = 3;
183         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
184         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
185     }
186     {
187         int              moleculeTypeIndex   = 2;
188         std::vector<int> numAtomsInResidues  = {3};
189         int              moleculeBlockIndex  = 2;
190         int              numMoleculesInBlock = 1;
191         topManager_.setMoleculeType(moleculeTypeIndex, numAtomsInResidues);
192         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
193     }
194     {
195         int moleculeTypeIndex   = 0; // Re-using a moltype definition
196         int moleculeBlockIndex  = 3;
197         int numMoleculesInBlock = 1;
198         topManager_.setMoleculeBlock(moleculeBlockIndex, moleculeTypeIndex, numMoleculesInBlock);
199     }
200     topManager_.finalizeTopology();
201 }
202
203 /********************************************************************
204  * gmx_ana_index_make_block() tests
205  */
206
207 TEST_F(IndexBlockTest, CreatesUnknownBlock)
208 {
209     gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
210     checkBlocka();
211     done_blocka(&blocka_);
212     gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
213     checkBlocka();
214 }
215
216 TEST_F(IndexBlockTest, CreatesAtomBlock)
217 {
218     const int group[] = { 0, 1, 3, 4, 6 };
219     setGroup(group);
220     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, false);
221     checkBlocka();
222     done_blocka(&blocka_);
223     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, true);
224     checkBlocka();
225 }
226
227 TEST_F(IndexBlockTest, CreatesResidueBlocksForSimpleTopology)
228 {
229     makeSimpleTopology();
230
231     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8 };
232     setGroup(group);
233     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
234     checkBlocka("FromAllAtoms");
235     done_blocka(&blocka_);
236     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
237     checkBlocka("FromAllAtoms");
238 }
239
240 TEST_F(IndexBlockTest, CreatesResidueBlocksForComplexTopology)
241 {
242     makeComplexTopology();
243
244     SCOPED_TRACE("Group with all atoms without completion");
245     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15 };
246     setGroup(group);
247     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
248     checkBlocka("FromAllAtoms");
249     done_blocka(&blocka_);
250     //SCOPED_TRACE("Group with all atoms with completion");
251     //gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
252     //checkBlocka("FromAllAtoms");
253     //done_blocka(&blocka_);
254
255     SCOPED_TRACE("Group with some atoms without completion");
256     const int subgroup[] = { 0, 3, 4, 5, 6, 7, 8, 12, 13, 15 };
257     setGroup(subgroup);
258     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, false);
259     checkBlocka("FromSomeAtomsWithoutCompletion");
260     //done_blocka(&blocka_);
261     //SCOPED_TRACE("Group with some atoms with completion");
262     //gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_RES, true);
263     //checkBlocka("FromSomeAtomsWithCompletion");
264 }
265
266 TEST_F(IndexBlockTest, CreatesMoleculeBlocksForSimpleTopology)
267 {
268     makeSimpleTopology();
269
270     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8 };
271     setGroup(group);
272     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
273     checkBlocka("FromAllAtoms");
274     done_blocka(&blocka_);
275     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
276     checkBlocka("FromAllAtoms");
277 }
278
279 TEST_F(IndexBlockTest, CreatesMoleculeBlocksForComplexTopology)
280 {
281     makeComplexTopology();
282
283     SCOPED_TRACE("Group with all atoms without completion");
284     const int group[] = { 0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15 };
285     setGroup(group);
286     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
287     checkBlocka("FromAllAtoms");
288     done_blocka(&blocka_);
289     //SCOPED_TRACE("Group with all atoms with completion");
290     //gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
291     //checkBlocka("FromAllAtoms");
292     //done_blocka(&blocka_);
293
294     SCOPED_TRACE("Group with some atoms without completion");
295     const int subgroup[] = { 1, 5, 6, 7, 11, 12 };
296     setGroup(subgroup);
297     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, false);
298     checkBlocka("FromSomeAtomsWithoutCompletion");
299     //done_blocka(&blocka_);
300     //SCOPED_TRACE("Group with some atoms with completion");
301     //gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_, INDEX_MOL, true);
302     //checkBlocka("FromSomeAtomsWithCompletion");
303 }
304
305 TEST_F(IndexBlockTest, CreatesSingleBlock)
306 {
307     const int group[] = { 0, 1, 3, 4, 6 };
308     setGroup(group);
309     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, false);
310     checkBlocka();
311     done_blocka(&blocka_);
312     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, true);
313     checkBlocka();
314 }
315
316 /********************************************************************
317  * gmx_ana_index_has_full_ablocks() tests
318  */
319
320 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
321 {
322     const int maxGroup[] = {
323         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
324     };
325     const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
326     topManager_.initAtoms(18);
327     topManager_.initUniformResidues(3);
328     setGroup(maxGroup);
329     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
330                              INDEX_RES, false);
331     setGroup(testGroup);
332     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
333 }
334
335 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
336 {
337     const int maxGroup[] = {
338         15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
339     };
340     const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
341     topManager_.initAtoms(18);
342     topManager_.initUniformResidues(3);
343     setGroup(maxGroup);
344     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
345                              INDEX_RES, false);
346     setGroup(testGroup);
347     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
348 }
349
350 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
351 {
352     const int maxGroup[] = {
353         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
354     };
355     const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
356     const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
357     const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
358
359     topManager_.initAtoms(18);
360     topManager_.initUniformResidues(3);
361     setGroup(maxGroup);
362     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
363                              INDEX_RES, false);
364
365     setGroup(testGroup1);
366     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
367
368     setGroup(testGroup2);
369     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
370
371     setGroup(testGroup3);
372     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
373 }
374
375 /********************************************************************
376  * gmx_ana_index_has_complete_elems() tests
377  */
378
379 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
380 {
381     const int group[] = { 0, 1, 2 };
382     setGroup(group);
383     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, nullptr));
384     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, nullptr));
385     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, nullptr));
386 }
387
388 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
389 {
390     const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
391     const int group2[] = { 3, 4, 5, 6, 7, 8 };
392
393     topManager_.initAtoms(15);
394     topManager_.initUniformResidues(3);
395     gmx_mtop_t *top = topManager_.topology();
396
397     setGroup(group1);
398     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
399
400     setGroup(group2);
401     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
402 }
403
404 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
405 {
406     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
407     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
408     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
409     const int group4[] = { 3, 4, 5, 6, 8, 12, 13, 14 };
410
411     topManager_.initAtoms(18);
412     topManager_.initUniformResidues(3);
413     gmx_mtop_t *top = topManager_.topology();
414
415     setGroup(group1);
416     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
417
418     setGroup(group2);
419     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
420
421     setGroup(group3);
422     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
423
424     setGroup(group4);
425     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
426 }
427
428 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
429 {
430     const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
431
432     topManager_.initAtoms(15);
433     topManager_.initUniformResidues(1);
434     topManager_.initUniformMolecules(3);
435     gmx_mtop_t *top = topManager_.topology();
436
437     setGroup(group);
438     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
439 }
440
441 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
442 {
443     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
444     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
445     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
446
447     topManager_.initAtoms(18);
448     topManager_.initUniformResidues(1);
449     topManager_.initUniformMolecules(3);
450     gmx_mtop_t *top = topManager_.topology();
451
452     setGroup(group1);
453     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
454
455     setGroup(group2);
456     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
457
458     setGroup(group3);
459     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
460 }
461
462 /********************************************************************
463  * IndexMapTest
464  */
465
466 class IndexMapTest : public ::testing::Test
467 {
468     public:
469         IndexMapTest();
470         ~IndexMapTest() override;
471
472         void testInit(int atomCount, const int atoms[], e_index_t type);
473         void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
474                         const char *name);
475         void testOrgIdGroup(e_index_t type, const char *name);
476         template <int count>
477         void testInit(const int (&atoms)[count], e_index_t type)
478         {
479             testInit(count, atoms, type);
480         }
481         template <int count>
482         void testUpdate(const int (&atoms)[count], bool bMaskOnly,
483                         const char *name)
484         {
485             testUpdate(count, atoms, bMaskOnly, name);
486         }
487
488         void checkMapping(int atomCount, const int atoms[], const char *name);
489
490         gmx::test::TestReferenceData    data_;
491         gmx::test::TestReferenceChecker checker_;
492         gmx::test::TopologyManager      topManager_;
493         gmx_ana_indexmap_t              map_;
494
495     private:
496         gmx_ana_index_t                 initGroup_;
497 };
498
499 IndexMapTest::IndexMapTest()
500     : checker_(data_.rootChecker())
501 {
502     gmx_ana_indexmap_clear(&map_);
503     gmx_ana_index_clear(&initGroup_);
504 }
505
506 IndexMapTest::~IndexMapTest()
507 {
508     gmx_ana_indexmap_deinit(&map_);
509 }
510
511 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
512 {
513     initGroup_.isize = atomCount;
514     initGroup_.index = const_cast<int *>(atoms);
515     gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
516     EXPECT_EQ(type, map_.type);
517     checkMapping(atomCount, atoms, "Initialized");
518 }
519
520 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
521                               const char *name)
522 {
523     gmx_ana_index_t g;
524     g.isize = atomCount;
525     g.index = const_cast<int *>(atoms);
526     gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
527     if (name == nullptr)
528     {
529         name = "Updated";
530     }
531     if (bMaskOnly)
532     {
533         checkMapping(initGroup_.isize, initGroup_.index, name);
534     }
535     else
536     {
537         checkMapping(atomCount, atoms, name);
538     }
539 }
540
541 void IndexMapTest::testOrgIdGroup(e_index_t type, const char *name)
542 {
543     gmx::test::TestReferenceChecker compound(
544             checker_.checkCompound("OrgIdGroups", name));
545     const int count
546         = gmx_ana_indexmap_init_orgid_group(&map_, topManager_.topology(), type);
547     compound.checkInteger(count, "GroupCount");
548     compound.checkSequenceArray(map_.mapb.nr, map_.orgid, "OrgId");
549     for (int i = 0; i < map_.mapb.nr; ++i)
550     {
551         EXPECT_EQ(map_.orgid[i], map_.mapid[i]);
552     }
553 }
554
555 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
556                                 const char *name)
557 {
558     gmx::test::TestReferenceChecker compound(
559             checker_.checkCompound("IndexMapping", name));
560     compound.checkSequenceArray(atomCount, atoms, "Input");
561     compound.checkInteger(map_.mapb.nr, "Count");
562     for (int i = 0; i < map_.mapb.nr; ++i)
563     {
564         gmx::test::TestReferenceChecker blockCompound(
565                 compound.checkCompound("Block", nullptr));
566         blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
567                                     &atoms[map_.mapb.index[i+1]],
568                                     "Atoms");
569         blockCompound.checkInteger(map_.refid[i], "RefId");
570         blockCompound.checkInteger(map_.mapid[i], "MapId");
571         int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
572         EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
573     }
574 }
575
576 /********************************************************************
577  * gmx_ana_indexmap_t tests
578  */
579
580 TEST_F(IndexMapTest, InitializesAtomBlock)
581 {
582     const int maxGroup[] = { 1, 2, 4, 5 };
583     testInit(maxGroup, INDEX_ATOM);
584 }
585
586 TEST_F(IndexMapTest, InitializesOrgIdGroupAtom)
587 {
588     const int maxGroup[] = { 2, 5, 7 };
589     testInit(maxGroup, INDEX_ATOM);
590     testOrgIdGroup(INDEX_ATOM, "Atoms");
591 }
592
593 TEST_F(IndexMapTest, InitializesOrgIdGroupSingle)
594 {
595     const int maxGroup[] = { 3, 4, 7, 8, 13 };
596     topManager_.initAtoms(18);
597     topManager_.initUniformResidues(3);
598     testInit(maxGroup, INDEX_RES);
599     testOrgIdGroup(INDEX_ATOM, "Single");
600 }
601
602 TEST_F(IndexMapTest, InitializesOrgIdGroupResidue)
603 {
604     const int maxGroup[] = { 3, 4, 7, 8, 13 };
605     topManager_.initAtoms(18);
606     topManager_.initUniformResidues(3);
607     testInit(maxGroup, INDEX_ATOM);
608     testOrgIdGroup(INDEX_RES, "Residues");
609 }
610
611 TEST_F(IndexMapTest, InitializesOrgIdGroupMolecule)
612 {
613     const int maxGroup[] = { 1, 2, 3, 4, 7, 8, 13 };
614     topManager_.initAtoms(18);
615     topManager_.initUniformResidues(3);
616     topManager_.initUniformMolecules(6);
617     testInit(maxGroup, INDEX_RES);
618     testOrgIdGroup(INDEX_MOL, "Molecules");
619 }
620
621 TEST_F(IndexMapTest, InitializesOrgIdGroupAll)
622 {
623     const int maxGroup[] = { 3, 4, 7, 8, 13 };
624     testInit(maxGroup, INDEX_ATOM);
625     testOrgIdGroup(INDEX_ALL, "All");
626 }
627
628 TEST_F(IndexMapTest, InitializesMoleculeBlock)
629 {
630     const int maxGroup[] = { 3, 4, 7, 8, 13 };
631     topManager_.initAtoms(18);
632     topManager_.initUniformResidues(1);
633     topManager_.initUniformMolecules(3);
634     testInit(maxGroup, INDEX_MOL);
635 }
636
637 TEST_F(IndexMapTest, MapsSingleBlock)
638 {
639     const int maxGroup[]  = { 0, 1, 2, 3 };
640     const int evalGroup[] = { 0, 2 };
641     testInit(maxGroup, INDEX_ALL);
642     testUpdate(evalGroup, false, nullptr);
643 }
644
645 TEST_F(IndexMapTest, MapsResidueBlocks)
646 {
647     const int maxGroup[] = {
648         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
649     };
650     const int evalGroup[] = { 3, 4, 7, 8, 13 };
651     topManager_.initAtoms(18);
652     topManager_.initUniformResidues(3);
653     testInit(maxGroup, INDEX_RES);
654     testUpdate(evalGroup, false, nullptr);
655 }
656
657 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
658 {
659     const int maxGroup[] = {
660         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
661     };
662     const int evalGroup[] = { 3, 4, 7, 8, 13 };
663     topManager_.initAtoms(18);
664     topManager_.initUniformResidues(3);
665     testInit(maxGroup, INDEX_RES);
666     testUpdate(evalGroup, true, nullptr);
667 }
668
669 TEST_F(IndexMapTest, HandlesMultipleRequests)
670 {
671     const int maxGroup[] = {
672         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
673     };
674     const int evalGroup[] = { 3, 4, 7, 8, 13 };
675     topManager_.initAtoms(18);
676     topManager_.initUniformResidues(3);
677     testInit(maxGroup, INDEX_RES);
678     testUpdate(evalGroup, false, "EvaluatedNoMask");
679     testUpdate(evalGroup, true, "EvaluatedMask");
680     testUpdate(maxGroup, true, "Initialized");
681     testUpdate(evalGroup, true, "EvaluatedMask");
682     testUpdate(evalGroup, false, "EvaluatedNoMask");
683     testUpdate(maxGroup, false, "Initialized");
684 }
685
686 /***********************************************************************
687  * IndexGroupsAndNames tests
688  */
689
690 class IndexGroupsAndNamesTest : public ::testing::Test
691 {
692     public:
693         IndexGroupsAndNamesTest()
694         {
695             init_blocka(&blockA_);
696             addGroupToBlocka_(indicesGroupA_);
697             addGroupToBlocka_(indicesGroupB_);
698             addGroupToBlocka_(indicesGroupSecondA_);
699             addGroupToBlocka_(indicesGroupC_);
700
701             const char *const namesAsConstCharArray[4] = {
702                 groupNames[0].c_str(), groupNames[1].c_str(),
703                 groupNames[2].c_str(), groupNames[3].c_str()
704             };
705             indexGroupAndNames_ = std::make_unique<gmx::IndexGroupsAndNames>(blockA_, namesAsConstCharArray);
706         }
707         ~IndexGroupsAndNamesTest() override
708         {
709             done_blocka(&blockA_);
710         }
711     protected:
712         std::unique_ptr<gmx::IndexGroupsAndNames> indexGroupAndNames_;
713         const std::vector<std::string>            groupNames     = { "A", "B - Name", "A", "C" };
714         const std::vector<gmx::index>             indicesGroupA_ = { };
715         const std::vector<gmx::index>             indicesGroupB_       = { 1, 2 };
716         const std::vector<gmx::index>             indicesGroupSecondA_ = { 5 };
717         const std::vector<gmx::index>             indicesGroupC_       = { 10 };
718
719     private:
720         //! Add a new group to t_blocka
721         void addGroupToBlocka_(gmx::ArrayRef<const gmx::index> index)
722         {
723             srenew(blockA_.index, blockA_.nr + 2);
724             srenew(blockA_.a, blockA_.nra + index.size());
725             for (int i = 0; (i < index.ssize()); i++)
726             {
727                 blockA_.a[blockA_.nra++] = index[i];
728             }
729             blockA_.nr++;
730             blockA_.index[blockA_.nr] = blockA_.nra;
731         }
732
733         t_blocka               blockA_;
734 };
735
736 TEST_F(IndexGroupsAndNamesTest, containsNames)
737 {
738     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("a"));
739     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("A"));
740     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("B - Name"));
741     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("b - Name"));
742     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("B - naMe"));
743     EXPECT_TRUE(indexGroupAndNames_->containsGroupName("C"));
744     EXPECT_FALSE(indexGroupAndNames_->containsGroupName("D"));
745 }
746
747 TEST_F(IndexGroupsAndNamesTest, throwsWhenNameMissing)
748 {
749     EXPECT_ANY_THROW(indexGroupAndNames_->indices("D"));
750 }
751
752 TEST_F(IndexGroupsAndNamesTest, groupIndicesCorrect)
753 {
754     using ::testing::Pointwise;
755     using ::testing::Eq;
756     EXPECT_THAT(indicesGroupA_, Pointwise(Eq(), indexGroupAndNames_->indices("A")));
757     EXPECT_THAT(indicesGroupB_, Pointwise(Eq(), indexGroupAndNames_->indices("B - name")));
758     EXPECT_THAT(indicesGroupC_, Pointwise(Eq(), indexGroupAndNames_->indices("C")));
759 }
760
761
762 } // namespace