Add function to get residue start and end
[gromacs.git] / src / gromacs / topology / tests / mtop.cpp
blob3bae46fab56015706fad9a0c5708089f07daa258
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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 mtop routines
39 * \author Roland Schulz <roland.schulz@intel.com>
40 * \ingroup module_topology
42 #include "gmxpre.h"
44 #include <gtest/gtest.h>
46 #include "gromacs/topology/mtop_util.h"
47 #include "gromacs/utility/basedefinitions.h"
49 namespace gmx
51 namespace
54 /*! \brief Initializes a basic topology with 9 atoms with settle*/
55 void createBasicTop(gmx_mtop_t* mtop)
57 gmx_moltype_t moltype;
58 moltype.atoms.nr = NRAL(F_SETTLE);
59 std::vector<int>& iatoms = moltype.ilist[F_SETTLE].iatoms;
60 const int settleType = 0;
61 iatoms.push_back(settleType);
62 iatoms.push_back(0);
63 iatoms.push_back(1);
64 iatoms.push_back(2);
65 mtop->moltype.push_back(moltype);
67 mtop->molblock.resize(1);
68 mtop->molblock[0].type = 0;
69 mtop->molblock[0].nmol = 3;
70 mtop->natoms = moltype.atoms.nr * mtop->molblock[0].nmol;
71 gmx_mtop_finalize(mtop);
74 /*! \brief
75 * Creates dummy topology with two differently sized residues.
77 * Residue begin and end are set to allow checking routines
78 * that make use of the boundaries.
80 * \returns The residue ranges.
82 std::vector<gmx::Range<int>> createTwoResidueTopology(gmx_mtop_t* mtop)
84 auto& moltype = mtop->moltype.emplace_back();
85 int residueOneSize = 5;
86 int residueTwoSize = 4;
87 moltype.atoms.nr = residueOneSize + residueTwoSize;
88 snew(moltype.atoms.atom, residueOneSize + residueTwoSize);
89 for (int i = 0; i < residueOneSize; i++)
91 moltype.atoms.atom[i].resind = 0;
93 for (int i = residueOneSize; i < residueOneSize + residueTwoSize; i++)
95 moltype.atoms.atom[i].resind = 1;
98 mtop->molblock.resize(1);
99 mtop->molblock[0].type = 0;
100 mtop->molblock[0].nmol = 1;
101 mtop->natoms = moltype.atoms.nr * mtop->molblock[0].nmol;
102 gmx_mtop_finalize(mtop);
103 std::vector<gmx::Range<int>> residueRange;
104 residueRange.emplace_back(0, residueOneSize);
105 residueRange.emplace_back(residueOneSize, residueOneSize + residueTwoSize);
106 return residueRange;
109 TEST(MtopTest, RangeBasedLoop)
111 gmx_mtop_t mtop;
112 createBasicTop(&mtop);
113 int count = 0;
114 for (const AtomProxy atomP : AtomRange(mtop))
116 EXPECT_EQ(atomP.globalAtomNumber(), count);
117 ++count;
119 EXPECT_EQ(count, 9);
122 TEST(MtopTest, Operators)
124 gmx_mtop_t mtop;
125 createBasicTop(&mtop);
126 AtomIterator it(mtop);
127 AtomIterator otherIt(mtop);
128 EXPECT_EQ((*it).globalAtomNumber(), 0);
129 EXPECT_EQ(it->globalAtomNumber(), 0);
130 EXPECT_TRUE(it == otherIt);
131 EXPECT_FALSE(it != otherIt);
132 ++it;
133 EXPECT_EQ(it->globalAtomNumber(), 1);
134 it++;
135 EXPECT_EQ(it->globalAtomNumber(), 2);
136 EXPECT_TRUE(it != otherIt);
137 EXPECT_FALSE(it == otherIt);
140 TEST(MtopTest, CanFindResidueStartAndEndAtoms)
142 gmx_mtop_t mtop;
143 auto expectedResidueRange = createTwoResidueTopology(&mtop);
145 auto atomRanges = atomRangeOfEachResidue(mtop.molblock[0], mtop);
146 ASSERT_EQ(atomRanges.size(), expectedResidueRange.size());
147 for (gmx::index i = 0; i < gmx::ssize(atomRanges); i++)
149 EXPECT_EQ(atomRanges[i].begin(), expectedResidueRange[i].begin());
150 EXPECT_EQ(atomRanges[i].end(), expectedResidueRange[i].end());
151 ASSERT_EQ(atomRanges[i].size(), expectedResidueRange[i].size());
153 int rangeIndex = 0;
154 for (const auto& range : atomRanges)
156 auto referenceRangeIt = expectedResidueRange[rangeIndex].begin();
157 for (int i : range)
159 EXPECT_EQ(i, *referenceRangeIt);
160 referenceRangeIt++;
162 rangeIndex++;
166 } // namespace
168 } // namespace gmx