2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
5 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
38 * Implements test helper routines from toputils.h.
40 * \author Teemu Murtola <teemu.murtola@gmail.com>
41 * \ingroup module_selection
53 #include "gromacs/fileio/confio.h"
54 #include "gromacs/fileio/trxio.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/topology/atoms.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/trajectory/trajectoryframe.h"
60 #include "gromacs/utility/arrayref.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/gmxassert.h"
63 #include "gromacs/utility/smalloc.h"
65 #include "testutils/testfilemanager.h"
72 TopologyManager::TopologyManager() : frame_(nullptr) {}
74 TopologyManager::~TopologyManager()
76 if (frame_ != nullptr)
85 for (char* atomtype : atomtypes_)
91 void TopologyManager::requestFrame()
93 GMX_RELEASE_ASSERT(mtop_ == nullptr, "Frame must be requested before initializing topology");
94 if (frame_ == nullptr)
100 void TopologyManager::requestVelocities()
102 GMX_RELEASE_ASSERT(frame_ != nullptr, "Velocities requested before requesting a frame");
104 if (frame_->natoms > 0)
106 snew(frame_->v, frame_->natoms);
110 void TopologyManager::requestForces()
112 GMX_RELEASE_ASSERT(frame_ != nullptr, "Forces requested before requesting a frame");
114 if (frame_->natoms > 0)
116 snew(frame_->f, frame_->natoms);
120 void TopologyManager::loadTopology(const char* filename)
124 rvec* xtop = nullptr;
127 GMX_RELEASE_ASSERT(mtop_ == nullptr, "Topology initialized more than once");
128 mtop_ = std::make_unique<gmx_mtop_t>();
129 readConfAndTopology(gmx::test::TestFileManager::getInputFilePath(filename).c_str(),
133 frame_ != nullptr ? &xtop : nullptr,
137 if (frame_ != nullptr)
139 GMX_ASSERT(xtop != nullptr, "Keep the static analyzer happy");
140 frame_->natoms = mtop_->natoms;
142 snew(frame_->x, frame_->natoms);
143 std::memcpy(frame_->x, xtop, sizeof(*frame_->x) * frame_->natoms);
145 copy_mat(box, frame_->box);
151 void TopologyManager::initAtoms(int count)
153 GMX_RELEASE_ASSERT(mtop_ == nullptr, "Topology initialized more than once");
154 mtop_ = std::make_unique<gmx_mtop_t>();
155 mtop_->moltype.resize(1);
156 init_t_atoms(&mtop_->moltype[0].atoms, count, FALSE);
157 mtop_->molblock.resize(1);
158 mtop_->molblock[0].type = 0;
159 mtop_->molblock[0].nmol = 1;
160 mtop_->natoms = count;
162 GMX_RELEASE_ASSERT(mtop_->maxResiduesPerMoleculeToTriggerRenumber() == 0,
163 "maxres_renum in mtop can be modified by an env.var., that is not supported "
165 t_atoms& atoms = this->atoms();
166 for (int i = 0; i < count; ++i)
168 atoms.atom[i].m = (i % 3 == 0 ? 2.0 : 1.0);
170 atoms.haveMass = TRUE;
171 if (frame_ != nullptr)
173 frame_->natoms = count;
175 snew(frame_->x, count);
178 snew(frame_->v, count);
182 snew(frame_->f, count);
187 void TopologyManager::initAtomTypes(const ArrayRef<const char* const>& types)
189 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
190 atomtypes_.reserve(types.size());
191 for (const char* type : types)
193 atomtypes_.push_back(gmx_strdup(type));
195 t_atoms& atoms = this->atoms();
196 snew(atoms.atomtype, atoms.nr);
198 for (int i = 0; i < atoms.nr; ++i, ++j)
200 if (j == types.ssize())
204 atoms.atomtype[i] = &atomtypes_[j];
206 atoms.haveType = TRUE;
209 void TopologyManager::initTopology(const int numMoleculeTypes, const int numMoleculeBlocks)
211 GMX_RELEASE_ASSERT(mtop_ == nullptr, "Topology initialized more than once");
212 mtop_ = std::make_unique<gmx_mtop_t>();
213 mtop_->moltype.resize(numMoleculeTypes);
214 mtop_->molblock.resize(numMoleculeBlocks);
217 void TopologyManager::setMoleculeType(const int moleculeTypeIndex, const ArrayRef<const int> residueSizes)
219 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
221 // Make a molecule block that refers to a new molecule type
222 auto& newMoleculeType = mtop_->moltype[moleculeTypeIndex];
223 auto& atoms = newMoleculeType.atoms;
225 auto numAtomsInMolecule = std::accumulate(residueSizes.begin(), residueSizes.end(), 0);
226 init_t_atoms(&atoms, numAtomsInMolecule, FALSE);
227 atoms.nres = residueSizes.size();
229 // Fill the residues in the molecule type
230 int residueIndex = 0;
231 auto residueSize = std::begin(residueSizes);
232 for (int i = 0; i < atoms.nr && residueSize != std::end(residueSizes); ++residueSize, ++residueIndex)
234 for (int j = 0; i < atoms.nr && j < *residueSize; ++i, ++j)
236 atoms.atom[i].resind = residueIndex;
237 atoms.atom[i].m = (i % 3 == 0 ? 2.0 : 1.0);
240 atoms.nres = residueIndex;
241 atoms.haveMass = true;
244 void TopologyManager::setMoleculeBlock(const int moleculeBlockIndex,
245 const int moleculeTypeIndex,
246 const int numMoleculesToAdd)
248 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
250 auto& newMoleculeBlock = mtop_->molblock[moleculeBlockIndex];
251 newMoleculeBlock.type = moleculeTypeIndex;
252 newMoleculeBlock.nmol = numMoleculesToAdd;
254 mtop_->natoms += numMoleculesToAdd * mtop_->moltype[moleculeTypeIndex].atoms.nr;
257 void TopologyManager::finalizeTopology()
259 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
261 mtop_->haveMoleculeIndices = true;
265 void TopologyManager::initUniformResidues(int residueSize)
267 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
268 t_atoms& atoms = this->atoms();
269 int residueIndex = -1;
270 for (int i = 0; i < atoms.nr; ++i)
272 if (i % residueSize == 0)
276 atoms.atom[i].resind = residueIndex;
278 atoms.nres = residueIndex;
281 void TopologyManager::initUniformMolecules(int moleculeSize)
283 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
284 GMX_RELEASE_ASSERT(mtop_->molblock.size() == 1,
285 "initUniformMolecules only implemented for a single molblock");
286 gmx_molblock_t& molblock = mtop_->molblock[0];
287 t_atoms& atoms = mtop_->moltype[molblock.type].atoms;
288 GMX_RELEASE_ASSERT(atoms.nr % moleculeSize == 0,
289 "The number of atoms should be a multiple of moleculeSize");
290 molblock.nmol = atoms.nr / moleculeSize;
291 atoms.nr = moleculeSize;
292 const int nres = atoms.atom[atoms.nr].resind;
293 GMX_RELEASE_ASSERT(atoms.atom[atoms.nr - 1].resind != nres,
294 "The residues should break at molecule boundaries");
296 mtop_->haveMoleculeIndices = true;
300 void TopologyManager::initFrameIndices(const ArrayRef<const int>& index)
302 GMX_RELEASE_ASSERT(frame_ != nullptr, "Frame not initialized");
303 GMX_RELEASE_ASSERT(!frame_->bIndex, "Frame atom indices can only be set once");
305 frame_->bIndex = TRUE;
306 snew(frame_->index, index.size());
307 std::copy(index.begin(), index.end(), frame_->index);
309 frame_->natoms = index.size();
312 t_atoms& TopologyManager::atoms()
314 GMX_RELEASE_ASSERT(mtop_ != nullptr, "Topology not initialized");
315 GMX_RELEASE_ASSERT(mtop_->natoms == mtop_->moltype[0].atoms.nr,
316 "Test setup assumes all atoms in a single molecule type");
317 return mtop_->moltype[0].atoms;