Merge release-5-0 into master
[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, 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 <gtest/gtest.h>
46
47 #include "gromacs/selection/indexutil.h"
48 #include "gromacs/topology/block.h"
49
50 #include "testutils/refdata.h"
51
52 #include "toputils.h"
53
54 namespace
55 {
56
57 /********************************************************************
58  * IndexBlockTest
59  */
60
61 class IndexBlockTest : public ::testing::Test
62 {
63     public:
64         IndexBlockTest();
65         ~IndexBlockTest();
66
67         void setGroup(int count, const int atoms[]);
68         template <int count>
69         void setGroup(const int (&atoms)[count])
70         {
71             setGroup(count, atoms);
72         }
73
74         void checkBlocka();
75
76         gmx::test::TestReferenceData    data_;
77         gmx::test::TopologyManager      topManager_;
78         gmx_ana_index_t                 g_;
79         t_blocka                        blocka_;
80 };
81
82 IndexBlockTest::IndexBlockTest()
83 {
84     blocka_.nr           = 0;
85     blocka_.index        = NULL;
86     blocka_.nalloc_index = 0;
87     blocka_.nra          = 0;
88     blocka_.a            = NULL;
89     blocka_.nalloc_a     = 0;
90     gmx_ana_index_clear(&g_);
91 }
92
93 IndexBlockTest::~IndexBlockTest()
94 {
95     done_blocka(&blocka_);
96 }
97
98 void IndexBlockTest::setGroup(int count, const int atoms[])
99 {
100     g_.isize = count;
101     g_.index = const_cast<int *>(atoms);
102 }
103
104 void IndexBlockTest::checkBlocka()
105 {
106     gmx::test::TestReferenceChecker compound(
107             data_.rootChecker().checkCompound("BlockAtoms", "Block"));
108     compound.checkSequenceArray(g_.isize, g_.index, "Input");
109     compound.checkInteger(blocka_.nr, "Count");
110     for (int i = 0; i < blocka_.nr; ++i)
111     {
112         gmx::test::TestReferenceChecker blockCompound(
113                 compound.checkCompound("Block", NULL));
114         blockCompound.checkSequence(&blocka_.a[blocka_.index[i]],
115                                     &blocka_.a[blocka_.index[i+1]],
116                                     "Atoms");
117     }
118 }
119
120 /********************************************************************
121  * gmx_ana_index_make_block() tests
122  */
123
124 TEST_F(IndexBlockTest, CreatesUnknownBlock)
125 {
126     gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
127     checkBlocka();
128     done_blocka(&blocka_);
129     gmx_ana_index_make_block(&blocka_, NULL, NULL, INDEX_UNKNOWN, false);
130     checkBlocka();
131 }
132
133 TEST_F(IndexBlockTest, CreatesAtomBlock)
134 {
135     const int group[] = { 0, 1, 3, 4, 6 };
136     setGroup(group);
137     gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, false);
138     checkBlocka();
139     done_blocka(&blocka_);
140     gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ATOM, true);
141     checkBlocka();
142 }
143
144 TEST_F(IndexBlockTest, CreatesResidueBlock)
145 {
146     const int group[] = { 0, 1, 3, 6, 7 };
147     topManager_.initAtoms(9);
148     topManager_.initUniformResidues(3);
149     setGroup(group);
150     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
151                              INDEX_RES, false);
152     checkBlocka();
153 }
154
155 TEST_F(IndexBlockTest, CreatesMoleculeBlock)
156 {
157     const int group[] = { 3, 4, 7, 8, 13 };
158     topManager_.initAtoms(18);
159     topManager_.initUniformMolecules(3);
160     setGroup(group);
161     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
162                              INDEX_MOL, false);
163     checkBlocka();
164 }
165
166 TEST_F(IndexBlockTest, CreatesResidueBlockWithCompletion)
167 {
168     const int group[] = { 3, 4, 7, 8, 13 };
169     topManager_.initAtoms(18);
170     topManager_.initUniformResidues(3);
171     setGroup(group);
172     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
173                              INDEX_RES, true);
174     checkBlocka();
175 }
176
177 TEST_F(IndexBlockTest, CreatesMoleculeBlockWithCompletion)
178 {
179     const int group[] = { 3, 4, 7, 8, 13 };
180     topManager_.initAtoms(18);
181     topManager_.initUniformMolecules(3);
182     setGroup(group);
183     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
184                              INDEX_MOL, true);
185     checkBlocka();
186 }
187
188 TEST_F(IndexBlockTest, CreatesSingleBlock)
189 {
190     const int group[] = { 0, 1, 3, 4, 6 };
191     setGroup(group);
192     gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, false);
193     checkBlocka();
194     done_blocka(&blocka_);
195     gmx_ana_index_make_block(&blocka_, NULL, &g_, INDEX_ALL, true);
196     checkBlocka();
197 }
198
199 /********************************************************************
200  * gmx_ana_index_has_full_ablocks() tests
201  */
202
203 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksPositive)
204 {
205     const int maxGroup[] = {
206         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
207     };
208     const int testGroup[] = { 3, 4, 5, 6, 7, 8, 12, 13, 14 };
209     topManager_.initAtoms(18);
210     topManager_.initUniformResidues(3);
211     setGroup(maxGroup);
212     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
213                              INDEX_RES, false);
214     setGroup(testGroup);
215     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
216 }
217
218 TEST_F(IndexBlockTest, ChecksOutOfOrderGroupForFullBlocksPositive)
219 {
220     const int maxGroup[] = {
221         15, 16, 17, 2, 1, 0, 12, 13, 14, 5, 4, 3, 9, 10, 11, 8, 7, 6
222     };
223     const int testGroup[] = { 2, 1, 0, 5, 4, 3, 8, 7, 6, };
224     topManager_.initAtoms(18);
225     topManager_.initUniformResidues(3);
226     setGroup(maxGroup);
227     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
228                              INDEX_RES, false);
229     setGroup(testGroup);
230     EXPECT_TRUE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
231 }
232
233 TEST_F(IndexBlockTest, ChecksGroupForFullBlocksNegative)
234 {
235     const int maxGroup[] = {
236         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
237     };
238     const int testGroup1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
239     const int testGroup2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
240     const int testGroup3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
241
242     topManager_.initAtoms(18);
243     topManager_.initUniformResidues(3);
244     setGroup(maxGroup);
245     gmx_ana_index_make_block(&blocka_, topManager_.topology(), &g_,
246                              INDEX_RES, false);
247
248     setGroup(testGroup1);
249     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
250
251     setGroup(testGroup2);
252     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
253
254     setGroup(testGroup3);
255     EXPECT_FALSE(gmx_ana_index_has_full_ablocks(&g_, &blocka_));
256 }
257
258 /********************************************************************
259  * gmx_ana_index_has_complete_elems() tests
260  */
261
262 TEST_F(IndexBlockTest, ChecksGroupForCompleteElementsTrivial)
263 {
264     const int group[] = { 0, 1, 2 };
265     setGroup(group);
266     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_ATOM, NULL));
267     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_ALL, NULL));
268     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_UNKNOWN, NULL));
269 }
270
271 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesPositive)
272 {
273     const int group1[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
274     const int group2[] = { 3, 4, 5, 6, 7, 8 };
275
276     topManager_.initAtoms(15);
277     topManager_.initUniformResidues(3);
278     t_topology *top = topManager_.topology();
279
280     setGroup(group1);
281     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
282
283     setGroup(group2);
284     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
285 }
286
287 TEST_F(IndexBlockTest, ChecksGroupForCompleteResiduesNegative)
288 {
289     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
290     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
291     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
292
293     topManager_.initAtoms(18);
294     topManager_.initUniformResidues(3);
295     t_topology *top = topManager_.topology();
296
297     setGroup(group1);
298     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
299
300     setGroup(group2);
301     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
302
303     setGroup(group3);
304     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_RES, top));
305 }
306
307 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesPositive)
308 {
309     const int group[] = { 0, 1, 2, 6, 7, 8, 12, 13, 14 };
310
311     topManager_.initAtoms(15);
312     topManager_.initUniformMolecules(3);
313     t_topology *top = topManager_.topology();
314
315     setGroup(group);
316     EXPECT_TRUE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
317 }
318
319 TEST_F(IndexBlockTest, ChecksGroupForCompleteMoleculesNegative)
320 {
321     const int group1[] = { 3, 4, 5, 6, 7, 8, 12, 13 };
322     const int group2[] = { 3, 4, 5, 6, 7, 12, 13, 14 };
323     const int group3[] = { 4, 5, 6, 7, 8, 12, 13, 14 };
324
325     topManager_.initAtoms(18);
326     topManager_.initUniformMolecules(3);
327     t_topology *top = topManager_.topology();
328
329     setGroup(group1);
330     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
331
332     setGroup(group2);
333     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
334
335     setGroup(group3);
336     EXPECT_FALSE(gmx_ana_index_has_complete_elems(&g_, INDEX_MOL, top));
337 }
338
339 /********************************************************************
340  * IndexMapTest
341  */
342
343 class IndexMapTest : public ::testing::Test
344 {
345     public:
346         IndexMapTest();
347         ~IndexMapTest();
348
349         void testInit(int atomCount, const int atoms[], e_index_t type);
350         void testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
351                         const char *name);
352         template <int count>
353         void testInit(const int (&atoms)[count], e_index_t type)
354         {
355             testInit(count, atoms, type);
356         }
357         template <int count>
358         void testUpdate(const int (&atoms)[count], bool bMaskOnly,
359                         const char *name)
360         {
361             testUpdate(count, atoms, bMaskOnly, name);
362         }
363
364         void checkMapping(int atomCount, const int atoms[], const char *name);
365
366         gmx::test::TestReferenceData    data_;
367         gmx::test::TestReferenceChecker checker_;
368         gmx::test::TopologyManager      topManager_;
369         gmx_ana_indexmap_t              map_;
370
371     private:
372         gmx_ana_index_t                 initGroup_;
373 };
374
375 IndexMapTest::IndexMapTest()
376     : checker_(data_.rootChecker())
377 {
378     gmx_ana_indexmap_clear(&map_);
379     gmx_ana_index_clear(&initGroup_);
380 }
381
382 IndexMapTest::~IndexMapTest()
383 {
384     gmx_ana_indexmap_deinit(&map_);
385 }
386
387 void IndexMapTest::testInit(int atomCount, const int atoms[], e_index_t type)
388 {
389     initGroup_.isize = atomCount;
390     initGroup_.index = const_cast<int *>(atoms);
391     gmx_ana_indexmap_init(&map_, &initGroup_, topManager_.topology(), type);
392     EXPECT_EQ(type, map_.type);
393     checkMapping(atomCount, atoms, "Initialized");
394 }
395
396 void IndexMapTest::testUpdate(int atomCount, const int atoms[], bool bMaskOnly,
397                               const char *name)
398 {
399     gmx_ana_index_t g;
400     g.isize = atomCount;
401     g.index = const_cast<int *>(atoms);
402     gmx_ana_indexmap_update(&map_, &g, bMaskOnly);
403     if (name == NULL)
404     {
405         name = "Updated";
406     }
407     if (bMaskOnly)
408     {
409         checkMapping(initGroup_.isize, initGroup_.index, name);
410     }
411     else
412     {
413         checkMapping(atomCount, atoms, name);
414     }
415 }
416
417 void IndexMapTest::checkMapping(int atomCount, const int atoms[],
418                                 const char *name)
419 {
420     gmx::test::TestReferenceChecker compound(
421             checker_.checkCompound("IndexMapping", name));
422     compound.checkSequenceArray(atomCount, atoms, "Input");
423     compound.checkInteger(map_.mapb.nr, "Count");
424     for (int i = 0; i < map_.mapb.nr; ++i)
425     {
426         gmx::test::TestReferenceChecker blockCompound(
427                 compound.checkCompound("Block", NULL));
428         blockCompound.checkSequence(&atoms[map_.mapb.index[i]],
429                                     &atoms[map_.mapb.index[i+1]],
430                                     "Atoms");
431         blockCompound.checkInteger(map_.refid[i], "RefId");
432         blockCompound.checkInteger(map_.mapid[i], "MapId");
433         int originalIdIndex = (map_.refid[i] != -1 ? map_.refid[i] : i);
434         EXPECT_EQ(map_.orgid[originalIdIndex], map_.mapid[i]);
435     }
436 }
437
438 /********************************************************************
439  * gmx_ana_indexmap_t tests
440  */
441
442 TEST_F(IndexMapTest, InitializesAtomBlock)
443 {
444     const int maxGroup[] = { 1, 2, 4, 5 };
445     testInit(maxGroup, INDEX_ATOM);
446 }
447
448 TEST_F(IndexMapTest, InitializesMoleculeBlock)
449 {
450     const int maxGroup[] = { 3, 4, 7, 8, 13 };
451     topManager_.initAtoms(18);
452     topManager_.initUniformMolecules(3);
453     testInit(maxGroup, INDEX_MOL);
454 }
455
456 TEST_F(IndexMapTest, MapsSingleBlock)
457 {
458     const int maxGroup[]  = { 0, 1, 2, 3 };
459     const int evalGroup[] = { 0, 2 };
460     testInit(maxGroup, INDEX_ALL);
461     testUpdate(evalGroup, false, NULL);
462 }
463
464 TEST_F(IndexMapTest, MapsResidueBlocks)
465 {
466     const int maxGroup[] = {
467         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
468     };
469     const int evalGroup[] = { 3, 4, 7, 8, 13 };
470     topManager_.initAtoms(18);
471     topManager_.initUniformResidues(3);
472     testInit(maxGroup, INDEX_RES);
473     testUpdate(evalGroup, false, NULL);
474 }
475
476 TEST_F(IndexMapTest, MapsResidueBlocksWithMask)
477 {
478     const int maxGroup[] = {
479         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
480     };
481     const int evalGroup[] = { 3, 4, 7, 8, 13 };
482     topManager_.initAtoms(18);
483     topManager_.initUniformResidues(3);
484     testInit(maxGroup, INDEX_RES);
485     testUpdate(evalGroup, true, NULL);
486 }
487
488 TEST_F(IndexMapTest, HandlesMultipleRequests)
489 {
490     const int maxGroup[] = {
491         0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17
492     };
493     const int evalGroup[] = { 3, 4, 7, 8, 13 };
494     topManager_.initAtoms(18);
495     topManager_.initUniformResidues(3);
496     testInit(maxGroup, INDEX_RES);
497     testUpdate(evalGroup, false, "EvaluatedNoMask");
498     testUpdate(evalGroup, true, "EvaluatedMask");
499     testUpdate(maxGroup, true, "Initialized");
500     testUpdate(evalGroup, true, "EvaluatedMask");
501     testUpdate(evalGroup, false, "EvaluatedNoMask");
502     testUpdate(maxGroup, false, "Initialized");
503 }
504
505 } // namespace