2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017,2018,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.
37 * Implements test helper routines from toputils.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
52 #include "gromacs/fileio/confio.h"
53 #include "gromacs/fileio/trxio.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/topology/atoms.h"
56 #include "gromacs/topology/mtop_util.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/arrayref.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/smalloc.h"
64 #include "testutils/testfilemanager.h"
71 TopologyManager::TopologyManager() : frame_(nullptr) {}
73 TopologyManager::~TopologyManager()
75 if (frame_ != nullptr)
84 for (char* atomtype : atomtypes_)
90 void TopologyManager::requestFrame()
92 GMX_RELEASE_ASSERT(mtop_ == nullptr, "Frame must be requested before initializing topology");
93 if (frame_ == nullptr)
99 void TopologyManager::requestVelocities()
101 GMX_RELEASE_ASSERT(frame_ != nullptr, "Velocities requested before requesting a frame");
103 if (frame_->natoms > 0)
105 snew(frame_->v, frame_->natoms);
109 void TopologyManager::requestForces()
111 GMX_RELEASE_ASSERT(frame_ != nullptr, "Forces requested before requesting a frame");
113 if (frame_->natoms > 0)
115 snew(frame_->f, frame_->natoms);
119 void TopologyManager::loadTopology(const char* filename)
123 rvec* xtop = nullptr;
126 GMX_RELEASE_ASSERT(mtop_ == nullptr, "Topology initialized more than once");
127 mtop_ = std::make_unique<gmx_mtop_t>();
128 readConfAndTopology(gmx::test::TestFileManager::getInputFilePath(filename).c_str(), &fullTopology,
129 mtop_.get(), &ePBC, frame_ != nullptr ? &xtop : nullptr, nullptr, box);
131 if (frame_ != nullptr)
133 GMX_ASSERT(xtop != nullptr, "Keep the static analyzer happy");
134 frame_->natoms = mtop_->natoms;
136 snew(frame_->x, frame_->natoms);
137 std::memcpy(frame_->x, xtop, sizeof(*frame_->x) * frame_->natoms);
139 copy_mat(box, frame_->box);
145 void TopologyManager::initAtoms(int count)
147 GMX_RELEASE_ASSERT(mtop_ == nullptr, "Topology initialized more than once");
148 mtop_ = std::make_unique<gmx_mtop_t>();
149 mtop_->moltype.resize(1);
150 init_t_atoms(&mtop_->moltype[0].atoms, count, FALSE);
151 mtop_->molblock.resize(1);
152 mtop_->molblock[0].type = 0;
153 mtop_->molblock[0].nmol = 1;
154 mtop_->natoms = count;
155 mtop_->maxres_renum = 0;
156 gmx_mtop_finalize(mtop_.get());
157 GMX_RELEASE_ASSERT(mtop_->maxres_renum == 0,
158 "maxres_renum in mtop can be modified by an env.var., that is not supported "
160 t_atoms& atoms = this->atoms();
161 for (int i = 0; i < count; ++i)
163 atoms.atom[i].m = (i % 3 == 0 ? 2.0 : 1.0);
165 atoms.haveMass = TRUE;
166 if (frame_ != nullptr)
168 frame_->natoms = count;
170 snew(frame_->x, count);
173 snew(frame_->v, count);
177 snew(frame_->f, count);
182 void TopologyManager::initAtomTypes(const ArrayRef<const char* const>& types)
184 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
185 atomtypes_.reserve(types.size());
186 for (const char* type : types)
188 atomtypes_.push_back(gmx_strdup(type));
190 t_atoms& atoms = this->atoms();
191 snew(atoms.atomtype, atoms.nr);
193 for (int i = 0; i < atoms.nr; ++i, ++j)
195 if (j == types.ssize())
199 atoms.atomtype[i] = &atomtypes_[j];
201 atoms.haveType = TRUE;
204 void TopologyManager::initTopology(const int numMoleculeTypes, const int numMoleculeBlocks)
206 GMX_RELEASE_ASSERT(mtop_ == nullptr, "Topology initialized more than once");
207 mtop_ = std::make_unique<gmx_mtop_t>();
208 mtop_->moltype.resize(numMoleculeTypes);
209 mtop_->molblock.resize(numMoleculeBlocks);
212 void TopologyManager::setMoleculeType(const int moleculeTypeIndex, const ArrayRef<const int> residueSizes)
214 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
216 // Make a molecule block that refers to a new molecule type
217 auto& newMoleculeType = mtop_->moltype[moleculeTypeIndex];
218 auto& atoms = newMoleculeType.atoms;
220 auto numAtomsInMolecule = std::accumulate(residueSizes.begin(), residueSizes.end(), 0);
221 init_t_atoms(&atoms, numAtomsInMolecule, FALSE);
222 atoms.nres = residueSizes.size();
224 // Fill the residues in the molecule type
225 int residueIndex = 0;
226 auto residueSize = std::begin(residueSizes);
227 for (int i = 0; i < atoms.nr && residueSize != std::end(residueSizes); ++residueSize, ++residueIndex)
229 for (int j = 0; i < atoms.nr && j < *residueSize; ++i, ++j)
231 atoms.atom[i].resind = residueIndex;
232 atoms.atom[i].m = (i % 3 == 0 ? 2.0 : 1.0);
235 atoms.nres = residueIndex;
236 atoms.haveMass = true;
239 void TopologyManager::setMoleculeBlock(const int moleculeBlockIndex,
240 const int moleculeTypeIndex,
241 const int numMoleculesToAdd)
243 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
245 auto& newMoleculeBlock = mtop_->molblock[moleculeBlockIndex];
246 newMoleculeBlock.type = moleculeTypeIndex;
247 newMoleculeBlock.nmol = numMoleculesToAdd;
249 mtop_->natoms += numMoleculesToAdd * mtop_->moltype[moleculeTypeIndex].atoms.nr;
252 void TopologyManager::finalizeTopology()
254 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
256 mtop_->maxres_renum = 0;
257 mtop_->haveMoleculeIndices = true;
258 gmx_mtop_finalize(mtop_.get());
261 void TopologyManager::initUniformResidues(int residueSize)
263 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
264 t_atoms& atoms = this->atoms();
265 int residueIndex = -1;
266 for (int i = 0; i < atoms.nr; ++i)
268 if (i % residueSize == 0)
272 atoms.atom[i].resind = residueIndex;
274 atoms.nres = residueIndex;
277 void TopologyManager::initUniformMolecules(int moleculeSize)
279 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
280 GMX_RELEASE_ASSERT(mtop_->molblock.size() == 1,
281 "initUniformMolecules only implemented for a single molblock");
282 gmx_molblock_t& molblock = mtop_->molblock[0];
283 t_atoms& atoms = mtop_->moltype[molblock.type].atoms;
284 GMX_RELEASE_ASSERT(atoms.nr % moleculeSize == 0,
285 "The number of atoms should be a multiple of moleculeSize");
286 molblock.nmol = atoms.nr / moleculeSize;
287 atoms.nr = moleculeSize;
288 const int nres = atoms.atom[atoms.nr].resind;
289 GMX_RELEASE_ASSERT(atoms.atom[atoms.nr - 1].resind != nres,
290 "The residues should break at molecule boundaries");
292 mtop_->haveMoleculeIndices = true;
293 gmx_mtop_finalize(mtop_.get());
296 void TopologyManager::initFrameIndices(const ArrayRef<const int>& index)
298 GMX_RELEASE_ASSERT(frame_ != nullptr, "Frame not initialized");
299 GMX_RELEASE_ASSERT(!frame_->bIndex, "Frame atom indices can only be set once");
301 frame_->bIndex = TRUE;
302 snew(frame_->index, index.size());
303 std::copy(index.begin(), index.end(), frame_->index);
305 frame_->natoms = index.size();
308 t_atoms& TopologyManager::atoms()
310 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
311 GMX_RELEASE_ASSERT(mtop_->natoms == mtop_->moltype[0].atoms.nr,
312 "Test setup assumes all atoms in a single molecule type");
313 return mtop_->moltype[0].atoms;