Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / topology / tests / exclusionblocks.cpp
blob1bc41d82dea999cc68dcddde52cb03e66b11921d
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
35 /*! \internal \file
36 * \brief
37 * Implements test of exclusionblock routines
39 * \author Paul Bauer <paul.bauer.q@gmail.com>
40 * \ingroup module_topology
42 #include "gmxpre.h"
44 #include "gromacs/topology/exclusionblocks.h"
46 #include <gtest/gtest.h>
48 #include "gromacs/topology/block.h"
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/smalloc.h"
52 #include "testutils/cmdlinetest.h"
54 namespace gmx
56 namespace testing
58 namespace
61 //! Add a new group to t_blocka
62 void addGroupToBlocka(t_blocka *b, gmx::ArrayRef<const int> index)
64 srenew(b->index, b->nr+2);
65 srenew(b->a, b->nra+index.size());
66 for (int i = 0; (i < index.ssize()); i++)
68 b->a[b->nra++] = index[i];
70 b->nr++;
71 b->index[b->nr] = b->nra;
74 //! Fill ExclusionBlock with data.
75 int fillExclusionBlock(gmx::ArrayRef<ExclusionBlock> b)
77 std::vector < std::vector < int>> index = {{0, 4, 7}, {1, 5, 8, 10}, {2, 6, 9, 11, 12}};
78 int nra = 0;
79 for (int i = 0; i < b.ssize(); i++)
81 b[i].atomNumber.clear();
82 for (const auto &j : index[i])
84 b[i].atomNumber.push_back(j);
86 nra += b[i].nra();
88 return nra;
91 //! Fill the t_blocka with some datastructures
92 void makeTestBlockAData(t_blocka *ba)
94 init_blocka(ba);
96 std::vector<int> index = {12, 11, 9, 6, 2};
97 addGroupToBlocka(ba, index);
98 index = {10, 8, 5, 1};
99 addGroupToBlocka(ba, index);
100 index = {7, 4, 0};
101 addGroupToBlocka(ba, index);
104 class ExclusionBlockTest : public ::testing::Test
106 public:
107 ExclusionBlockTest()
109 const int natom = 3;
110 makeTestBlockAData(&ba_);
111 b_.resize(natom);
113 ~ExclusionBlockTest() override
115 done_blocka(&ba_);
118 void compareBlocks()
120 for (unsigned i = 0; i < b_.size(); i++)
122 int index = ba_.index[i];
123 for (int j = 0; j < b_[i].nra(); j++)
125 int pos = index + j;
126 EXPECT_EQ(b_[i].atomNumber[j], ba_.a[pos])<< "Block mismatch at " << i << " , " << j << ".";
130 protected:
131 t_blocka ba_;
132 std::vector<ExclusionBlock> b_;
135 TEST_F(ExclusionBlockTest, ConvertBlockAToExclusionBlocks)
137 blockaToExclusionBlocks(&ba_, b_);
138 compareBlocks();
141 TEST_F(ExclusionBlockTest, ConvertExclusionBlockToBlocka)
143 int nra = fillExclusionBlock(b_);
144 srenew(ba_.a, nra+1);
145 srenew(ba_.index, b_.size()+1);
146 exclusionBlocksToBlocka(b_, &ba_);
147 compareBlocks();
150 TEST_F(ExclusionBlockTest, MergeExclusions)
152 mergeExclusions(&ba_, b_);
153 compareBlocks();
156 } // namespace
158 } // namespace testing
160 } // namespace gmx