Use ListOfLists in gmx_mtop_t and gmx_localtop_t
[alexxy/gromacs.git] / src / gromacs / topology / tests / exclusionblocks.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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  * Implements test of exclusionblock routines
38  *
39  * \author Paul Bauer <paul.bauer.q@gmail.com>
40  * \ingroup module_topology
41  */
42 #include "gmxpre.h"
43
44 #include "gromacs/topology/exclusionblocks.h"
45
46 #include <gtest/gtest.h>
47
48 #include "gromacs/topology/block.h"
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/listoflists.h"
51 #include "gromacs/utility/smalloc.h"
52
53 #include "testutils/cmdlinetest.h"
54
55 namespace gmx
56 {
57 namespace testing
58 {
59 namespace
60 {
61
62 //! Add a new group to t_blocka
63 void addGroupToBlocka(t_blocka* b, gmx::ArrayRef<const int> indices)
64 {
65     srenew(b->index, b->nr + 2);
66     srenew(b->a, b->nra + indices.size());
67     for (index i = 0; i < indices.ssize(); i++)
68     {
69         b->a[b->nra++] = indices[i];
70     }
71     b->nr++;
72     b->index[b->nr] = b->nra;
73 }
74
75 //! Fill ExclusionBlock with data.
76 int fillExclusionBlock(gmx::ArrayRef<ExclusionBlock> b)
77 {
78     std::vector<std::vector<int>> indices = { { 0, 4, 7 }, { 1, 5, 8, 10 }, { 2, 6, 9, 11, 12 } };
79     int                           nra     = 0;
80     for (index i = 0; i < b.ssize(); i++)
81     {
82         b[i].atomNumber.clear();
83         for (const auto& j : indices[i])
84         {
85             b[i].atomNumber.push_back(j);
86         }
87         nra += b[i].nra();
88     }
89     return nra;
90 }
91
92 //! Fill the t_blocka with some datastructures
93 void makeTestBlockAData(t_blocka* ba)
94 {
95     init_blocka(ba);
96
97     std::vector<int> indices = { 12, 11, 9, 6, 2 };
98     addGroupToBlocka(ba, indices);
99     indices = { 10, 8, 5, 1 };
100     addGroupToBlocka(ba, indices);
101     indices = { 7, 4, 0 };
102     addGroupToBlocka(ba, indices);
103 }
104
105 //! Return ListOfLists filled with some datastructures
106 ListOfLists<int> makeTestListOfLists()
107 {
108     ListOfLists<int> list;
109
110     std::vector<int> indices = { 12, 11, 9, 6, 2 };
111     list.pushBack(indices);
112     indices = { 10, 8, 5, 1 };
113     list.pushBack(indices);
114     indices = { 7, 4, 0 };
115     list.pushBack(indices);
116
117     return list;
118 }
119
120 class ExclusionBlockTest : public ::testing::Test
121 {
122 public:
123     ExclusionBlockTest()
124     {
125         const int natom = 3;
126         makeTestBlockAData(&ba_);
127         list_ = makeTestListOfLists();
128         b_.resize(natom);
129     }
130     ~ExclusionBlockTest() override { done_blocka(&ba_); }
131
132     void compareBlocks()
133     {
134         for (index i = 0; i < ssize(b_); i++)
135         {
136             int index = ba_.index[i];
137             for (int j = 0; j < b_[i].nra(); j++)
138             {
139                 int pos = index + j;
140                 EXPECT_EQ(b_[i].atomNumber[j], ba_.a[pos])
141                         << "Block mismatch at " << i << " , " << j << ".";
142             }
143         }
144     }
145
146     void compareBlocksAndList()
147     {
148         GMX_RELEASE_ASSERT(ssize(b_) == list_.ssize(), "The list counts should match");
149         for (index i = 0; i < ssize(b_); i++)
150         {
151             gmx::ArrayRef<const int> jList = list_[i];
152             EXPECT_EQ(b_[i].nra(), jList.ssize()) << "Block size mismatch at " << i << ".";
153             for (int j = 0; j < b_[i].nra(); j++)
154             {
155                 EXPECT_EQ(b_[i].atomNumber[j], jList[j]) << "Block mismatch at " << i << " , " << j << ".";
156             }
157         }
158     }
159
160 protected:
161     t_blocka                    ba_;
162     ListOfLists<int>            list_;
163     std::vector<ExclusionBlock> b_;
164 };
165
166 TEST_F(ExclusionBlockTest, ConvertBlockAToExclusionBlocks)
167 {
168     blockaToExclusionBlocks(&ba_, b_);
169     compareBlocks();
170 }
171
172 TEST_F(ExclusionBlockTest, ConvertExclusionBlockToBlocka)
173 {
174     int nra = fillExclusionBlock(b_);
175     srenew(ba_.a, nra + 1);
176     srenew(ba_.index, b_.size() + 1);
177     exclusionBlocksToBlocka(b_, &ba_);
178     compareBlocks();
179 }
180
181 TEST_F(ExclusionBlockTest, MergeExclusions)
182 {
183     mergeExclusions(&list_, b_);
184     compareBlocksAndList();
185 }
186
187 } // namespace
188
189 } // namespace testing
190
191 } // namespace gmx