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