c91519e34377a3b32a874d153c783025a762f19e
[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, 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 /********************************************************************
79  * IndexBlockTest
80  */
81
82 class IndexBlockTest : public ::testing::Test
83 {
84     public:
85         IndexBlockTest();
86         ~IndexBlockTest();
87
88         void setGroup(int count, const int atoms[]);
89         template <int count>
90         void setGroup(const int (&atoms)[count])
91         {
92             setGroup(count, atoms);
93         }
94
95         void checkBlocka();
96
97         gmx::test::TestReferenceData    data_;
98         gmx::test::TopologyManager      topManager_;
99         gmx_ana_index_t                 g_;
100         t_blocka                        blocka_;
101 };
102
103 IndexBlockTest::IndexBlockTest()
104 {
105     blocka_.nr           = 0;
106     blocka_.index        = nullptr;
107     blocka_.nalloc_index = 0;
108     blocka_.nra          = 0;
109     blocka_.a            = nullptr;
110     blocka_.nalloc_a     = 0;
111     gmx_ana_index_clear(&g_);
112 }
113
114 IndexBlockTest::~IndexBlockTest()
115 {
116     done_blocka(&blocka_);
117 }
118
119 void IndexBlockTest::setGroup(int count, const int atoms[])
120 {
121     g_.isize = count;
122     g_.index = const_cast<int *>(atoms);
123 }
124
125 void IndexBlockTest::checkBlocka()
126 {
127     gmx::test::TestReferenceChecker compound(
128             data_.rootChecker().checkCompound("BlockAtoms", "Block"));
129     compound.checkSequenceArray(g_.isize, g_.index, "Input");
130     compound.checkInteger(blocka_.nr, "Count");
131     for (int i = 0; i < blocka_.nr; ++i)
132     {
133         gmx::test::TestReferenceChecker blockCompound(
134                 compound.checkCompound("Block", nullptr));
135         blockCompound.checkSequence(&blocka_.a[blocka_.index[i]],
136                                     &blocka_.a[blocka_.index[i+1]],
137                                     "Atoms");
138     }
139 }
140
141 /********************************************************************
142  * gmx_ana_index_make_block() tests
143  */
144
145 TEST_F(IndexBlockTest, CreatesUnknownBlock)
146 {
147     gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
148     checkBlocka();
149     done_blocka(&blocka_);
150     gmx_ana_index_make_block(&blocka_, nullptr, nullptr, INDEX_UNKNOWN, false);
151     checkBlocka();
152 }
153
154 TEST_F(IndexBlockTest, CreatesAtomBlock)
155 {
156     const int group[] = { 0, 1, 3, 4, 6 };
157     setGroup(group);
158     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, false);
159     checkBlocka();
160     done_blocka(&blocka_);
161     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ATOM, true);
162     checkBlocka();
163 }
164
165 TEST_F(IndexBlockTest, CreatesResidueBlock)
166 {
167     const int group[] = { 0, 1, 3, 6, 7 };
168     topManager_.initAtoms(9);
169     topManager_.initUniformResidues(3);
170     setGroup(group);
171     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
172                              INDEX_RES, false);
173     checkBlocka();
174 }
175
176 TEST_F(IndexBlockTest, CreatesMoleculeBlock)
177 {
178     const int group[] = { 3, 4, 7, 8, 13 };
179     topManager_.initAtoms(18);
180     topManager_.initUniformResidues(1);
181     topManager_.initUniformMolecules(3);
182     setGroup(group);
183     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
184                              INDEX_MOL, false);
185     checkBlocka();
186 }
187
188 TEST_F(IndexBlockTest, CreatesResidueBlockWithCompletion)
189 {
190     const int group[] = { 3, 4, 7, 8, 13 };
191     topManager_.initAtoms(18);
192     topManager_.initUniformResidues(3);
193     setGroup(group);
194     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
195                              INDEX_RES, true);
196     checkBlocka();
197 }
198
199 TEST_F(IndexBlockTest, CreatesMoleculeBlockWithCompletion)
200 {
201     const int group[] = { 3, 4, 7, 8, 13 };
202     topManager_.initAtoms(18);
203     topManager_.initUniformResidues(1);
204     topManager_.initUniformMolecules(3);
205     setGroup(group);
206     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
207                              INDEX_MOL, true);
208     checkBlocka();
209 }
210
211 TEST_F(IndexBlockTest, CreatesSingleBlock)
212 {
213     const int group[] = { 0, 1, 3, 4, 6 };
214     setGroup(group);
215     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, false);
216     checkBlocka();
217     done_blocka(&blocka_);
218     gmx_ana_index_make_block(&blocka_, nullptr, &g_, INDEX_ALL, true);
219     checkBlocka();
220 }
221
222 /********************************************************************
223  * gmx_ana_index_has_full_ablocks() tests
224  */
225
226 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
227 {
228     const int maxGroup[] = {
229         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
230     };
231     const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
232     topManager_.initAtoms(18);
233     topManager_.initUniformResidues(3);
234     setGroup(maxGroup);
235     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
236                              INDEX_RES, false);
237     setGroup(testGroup);
238     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
239 }
240
241 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
242 {
243     const int maxGroup[] = {
244         15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
245     };
246     const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
247     topManager_.initAtoms(18);
248     topManager_.initUniformResidues(3);
249     setGroup(maxGroup);
250     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
251                              INDEX_RES, false);
252     setGroup(testGroup);
253     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
254 }
255
256 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
257 {
258     const int maxGroup[] = {
259         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
260     };
261     const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
262     const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
263     const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
264
265     topManager_.initAtoms(18);
266     topManager_.initUniformResidues(3);
267     setGroup(maxGroup);
268     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
269                              INDEX_RES, false);
270
271     setGroup(testGroup1);
272     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
273
274     setGroup(testGroup2);
275     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
276
277     setGroup(testGroup3);
278     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
279 }
280
281 /********************************************************************
282  * gmx_ana_index_has_complete_elems() tests
283  */
284
285 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
286 {
287     const int group[] = { 0, 1, 2 };
288     setGroup(group);
289     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, nullptr));
290     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, nullptr));
291     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, nullptr));
292 }
293
294 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
295 {
296     const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
297     const int group2[] = { 3, 4, 5, 6, 7, 8 };
298
299     topManager_.initAtoms(15);
300     topManager_.initUniformResidues(3);
301     gmx_mtop_t *top = topManager_.topology();
302
303     setGroup(group1);
304     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
305
306     setGroup(group2);
307     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
308 }
309
310 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
311 {
312     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
313     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
314     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
315     const int group4[] = { 3, 4, 5, 6, 8, 12, 13, 14 };
316
317     topManager_.initAtoms(18);
318     topManager_.initUniformResidues(3);
319     gmx_mtop_t *top = topManager_.topology();
320
321     setGroup(group1);
322     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
323
324     setGroup(group2);
325     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
326
327     setGroup(group3);
328     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
329
330     setGroup(group4);
331     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
332 }
333
334 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
335 {
336     const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
337
338     topManager_.initAtoms(15);
339     topManager_.initUniformResidues(1);
340     topManager_.initUniformMolecules(3);
341     gmx_mtop_t *top = topManager_.topology();
342
343     setGroup(group);
344     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
345 }
346
347 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
348 {
349     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
350     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
351     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
352
353     topManager_.initAtoms(18);
354     topManager_.initUniformResidues(1);
355     topManager_.initUniformMolecules(3);
356     gmx_mtop_t *top = topManager_.topology();
357
358     setGroup(group1);
359     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
360
361     setGroup(group2);
362     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
363
364     setGroup(group3);
365     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
366 }
367
368 /********************************************************************
369  * IndexMapTest
370  */
371
372 class IndexMapTest : public ::testing::Test
373 {
374     public:
375         IndexMapTest();
376         ~IndexMapTest();
377
378         void testInit(int atomCount, const int atoms[], e_index_t type);
379         void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
380                         const char *name);
381         void testOrgIdGroup(e_index_t type, const char *name);
382         template <int count>
383         void testInit(const int (&atoms)[count], e_index_t type)
384         {
385             testInit(count, atoms, type);
386         }
387         template <int count>
388         void testUpdate(const int (&atoms)[count], bool bMaskOnly,
389                         const char *name)
390         {
391             testUpdate(count, atoms, bMaskOnly, name);
392         }
393
394         void checkMapping(int atomCount, const int atoms[], const char *name);
395
396         gmx::test::TestReferenceData    data_;
397         gmx::test::TestReferenceChecker checker_;
398         gmx::test::TopologyManager      topManager_;
399         gmx_ana_indexmap_t              map_;
400
401     private:
402         gmx_ana_index_t                 initGroup_;
403 };
404
405 IndexMapTest::IndexMapTest()
406     : checker_(data_.rootChecker())
407 {
408     gmx_ana_indexmap_clear(&map_);
409     gmx_ana_index_clear(&initGroup_);
410 }
411
412 IndexMapTest::~IndexMapTest()
413 {
414     gmx_ana_indexmap_deinit(&map_);
415 }
416
417 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
418 {
419     initGroup_.isize = atomCount;
420     initGroup_.index = const_cast<int *>(atoms);
421     gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
422     EXPECT_EQ(type, map_.type);
423     checkMapping(atomCount, atoms, "Initialized");
424 }
425
426 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
427                               const char *name)
428 {
429     gmx_ana_index_t g;
430     g.isize = atomCount;
431     g.index = const_cast<int *>(atoms);
432     gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
433     if (name == nullptr)
434     {
435         name = "Updated";
436     }
437     if (bMaskOnly)
438     {
439         checkMapping(initGroup_.isize, initGroup_.index, name);
440     }
441     else
442     {
443         checkMapping(atomCount, atoms, name);
444     }
445 }
446
447 void IndexMapTest::testOrgIdGroup(e_index_t type, const char *name)
448 {
449     gmx::test::TestReferenceChecker compound(
450             checker_.checkCompound("OrgIdGroups", name));
451     const int count
452         = gmx_ana_indexmap_init_orgid_group(&map_, topManager_.topology(), type);
453     compound.checkInteger(count, "GroupCount");
454     compound.checkSequenceArray(map_.mapb.nr, map_.orgid, "OrgId");
455     for (int i = 0; i < map_.mapb.nr; ++i)
456     {
457         EXPECT_EQ(map_.orgid[i], map_.mapid[i]);
458     }
459 }
460
461 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
462                                 const char *name)
463 {
464     gmx::test::TestReferenceChecker compound(
465             checker_.checkCompound("IndexMapping", name));
466     compound.checkSequenceArray(atomCount, atoms, "Input");
467     compound.checkInteger(map_.mapb.nr, "Count");
468     for (int i = 0; i < map_.mapb.nr; ++i)
469     {
470         gmx::test::TestReferenceChecker blockCompound(
471                 compound.checkCompound("Block", nullptr));
472         blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
473                                     &atoms[map_.mapb.index[i+1]],
474                                     "Atoms");
475         blockCompound.checkInteger(map_.refid[i], "RefId");
476         blockCompound.checkInteger(map_.mapid[i], "MapId");
477         int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
478         EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
479     }
480 }
481
482 /********************************************************************
483  * gmx_ana_indexmap_t tests
484  */
485
486 TEST_F(IndexMapTest, InitializesAtomBlock)
487 {
488     const int maxGroup[] = { 1, 2, 4, 5 };
489     testInit(maxGroup, INDEX_ATOM);
490 }
491
492 TEST_F(IndexMapTest, InitializesOrgIdGroupAtom)
493 {
494     const int maxGroup[] = { 2, 5, 7 };
495     testInit(maxGroup, INDEX_ATOM);
496     testOrgIdGroup(INDEX_ATOM, "Atoms");
497 }
498
499 TEST_F(IndexMapTest, InitializesOrgIdGroupSingle)
500 {
501     const int maxGroup[] = { 3, 4, 7, 8, 13 };
502     topManager_.initAtoms(18);
503     topManager_.initUniformResidues(3);
504     testInit(maxGroup, INDEX_RES);
505     testOrgIdGroup(INDEX_ATOM, "Single");
506 }
507
508 TEST_F(IndexMapTest, InitializesOrgIdGroupResidue)
509 {
510     const int maxGroup[] = { 3, 4, 7, 8, 13 };
511     topManager_.initAtoms(18);
512     topManager_.initUniformResidues(3);
513     testInit(maxGroup, INDEX_ATOM);
514     testOrgIdGroup(INDEX_RES, "Residues");
515 }
516
517 TEST_F(IndexMapTest, InitializesOrgIdGroupMolecule)
518 {
519     const int maxGroup[] = { 1, 2, 3, 4, 7, 8, 13 };
520     topManager_.initAtoms(18);
521     topManager_.initUniformResidues(3);
522     topManager_.initUniformMolecules(6);
523     testInit(maxGroup, INDEX_RES);
524     testOrgIdGroup(INDEX_MOL, "Molecules");
525 }
526
527 TEST_F(IndexMapTest, InitializesOrgIdGroupAll)
528 {
529     const int maxGroup[] = { 3, 4, 7, 8, 13 };
530     testInit(maxGroup, INDEX_ATOM);
531     testOrgIdGroup(INDEX_ALL, "All");
532 }
533
534 TEST_F(IndexMapTest, InitializesMoleculeBlock)
535 {
536     const int maxGroup[] = { 3, 4, 7, 8, 13 };
537     topManager_.initAtoms(18);
538     topManager_.initUniformResidues(1);
539     topManager_.initUniformMolecules(3);
540     testInit(maxGroup, INDEX_MOL);
541 }
542
543 TEST_F(IndexMapTest, MapsSingleBlock)
544 {
545     const int maxGroup[]  = { 0, 1, 2, 3 };
546     const int evalGroup[] = { 0, 2 };
547     testInit(maxGroup, INDEX_ALL);
548     testUpdate(evalGroup, false, nullptr);
549 }
550
551 TEST_F(IndexMapTest, MapsResidueBlocks)
552 {
553     const int maxGroup[] = {
554         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
555     };
556     const int evalGroup[] = { 3, 4, 7, 8, 13 };
557     topManager_.initAtoms(18);
558     topManager_.initUniformResidues(3);
559     testInit(maxGroup, INDEX_RES);
560     testUpdate(evalGroup, false, nullptr);
561 }
562
563 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
564 {
565     const int maxGroup[] = {
566         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
567     };
568     const int evalGroup[] = { 3, 4, 7, 8, 13 };
569     topManager_.initAtoms(18);
570     topManager_.initUniformResidues(3);
571     testInit(maxGroup, INDEX_RES);
572     testUpdate(evalGroup, true, nullptr);
573 }
574
575 TEST_F(IndexMapTest, HandlesMultipleRequests)
576 {
577     const int maxGroup[] = {
578         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
579     };
580     const int evalGroup[] = { 3, 4, 7, 8, 13 };
581     topManager_.initAtoms(18);
582     topManager_.initUniformResidues(3);
583     testInit(maxGroup, INDEX_RES);
584     testUpdate(evalGroup, false, "EvaluatedNoMask");
585     testUpdate(evalGroup, true, "EvaluatedMask");
586     testUpdate(maxGroup, true, "Initialized");
587     testUpdate(evalGroup, true, "EvaluatedMask");
588     testUpdate(evalGroup, false, "EvaluatedNoMask");
589     testUpdate(maxGroup, false, "Initialized");
590 }
591
592 } // namespace