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