2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2012,2013,2014,2015,2016,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_TOPOLOGY_MTOP_UTIL_H
38 #define GMX_TOPOLOGY_MTOP_UTIL_H
44 #include "gromacs/topology/topology.h"
45 #include "gromacs/utility/basedefinitions.h"
47 struct gmx_localtop_t;
52 enum struct GmxQmmmMode;
54 // TODO All of the functions taking a const gmx_mtop * are deprecated
55 // and should be replaced by versions taking const gmx_mtop & when
56 // their callers are refactored similarly.
58 /* Should be called after generating or reading mtop,
59 * to set some compute intesive variables to avoid
60 * N^2 operations later on.
63 gmx_mtop_finalize(gmx_mtop_t *mtop);
65 /* Counts the number of atoms of each type. State should be 0 for
66 * state A and 1 for state B types. typecount should have at
67 * least mtop->ffparams.atnr elements.
70 gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[]);
72 /* Returns the total number of charge groups in mtop */
74 ncg_mtop(const gmx_mtop_t *mtop);
76 /*!\brief Returns the total number of molecules in mtop
78 * \param[in] mtop The global topology
80 int gmx_mtop_num_molecules(const gmx_mtop_t &mtop);
82 /* Returns the total number of residues in mtop. */
83 int gmx_mtop_nres(const gmx_mtop_t *mtop);
85 /* Removes the charge groups, i.e. makes single atom charge groups, in mtop */
86 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop);
90 //! Proxy object returned from AtomIterator
94 //! Default constructor.
95 AtomProxy(const AtomIterator* it) : it_(it) {}
96 //! Access current global atom number.
97 int globalAtomNumber() const;
98 //! Access current t_atom struct.
99 const t_atom &atom() const;
100 //! Access current name of the atom.
101 const char *atomName() const;
102 //! Access current name of the residue the atom is in.
103 const char *residueName() const;
104 //! Access current residue number.
105 int residueNumber() const;
106 //! Access current molecule type.
107 const gmx_moltype_t &moleculeType() const;
108 //! Access the position of the current atom in the molecule.
109 int atomNumberInMol() const;
111 const AtomIterator* it_;
114 //! Wrapper around proxy object to implement operator->
115 template <typename T>
119 //! Construct with proxy object.
120 ProxyPtr(T t) : t_(t) {}
121 //! Member of pointer operator.
122 T* operator->() { return &t_; }
128 * Object that allows looping over all atoms in an mtop.
133 //! Construct from topology and optionalally a global atom number.
134 explicit AtomIterator(const gmx_mtop_t &mtop, int globalAtomNumber = 0);
136 //! Prefix increment.
137 AtomIterator &operator++();
138 //! Postfix increment.
139 AtomIterator operator++(int);
141 //! Equality comparison.
142 bool operator==(const AtomIterator &o) const;
143 //! Non-equal comparison.
144 bool operator!=(const AtomIterator &o) const;
146 //! Dereference operator. Returns proxy.
147 AtomProxy operator*() const { return {this}; }
148 //! Member of pointer operator.
149 ProxyPtr<AtomProxy> operator->() const { return {this}; };
153 const gmx_mtop_t *mtop_;
154 //! Current molecule block.
156 //! The atoms of the current molecule.
157 const t_atoms *atoms_;
158 //! The current molecule.
159 int currentMolecule_;
160 //! Current highest number for residues.
161 int highestResidueNumber_;
162 //! Current local atom number.
163 int localAtomNumber_;
164 //! Global current atom number.
165 int globalAtomNumber_;
167 friend class AtomProxy;
170 //! Range over all atoms of topology.
174 //! Default constructor.
175 explicit AtomRange(const gmx_mtop_t &mtop) :
176 begin_(mtop), end_(mtop, mtop.natoms) {}
177 //! Iterator to begin of range.
178 AtomIterator &begin() { return begin_; }
179 //! Iterator to end of range.
180 AtomIterator &end() { return end_; }
182 AtomIterator begin_, end_;
185 /* Abstract type for atom loop over atoms in all molecule blocks */
186 typedef struct gmx_mtop_atomloop_block *gmx_mtop_atomloop_block_t;
188 /* Initialize an atom loop over atoms in all molecule blocks the system.
190 gmx_mtop_atomloop_block_t
191 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop);
193 /* Loop to the next atom.
194 * When not at the end:
196 * sets the pointer atom to the t_atom struct of that atom
197 * and return the number of molecules corresponding to this atom.
198 * When at the end, destroys aloop and returns FALSE.
200 * gmx_mtop_atomloop_block_t aloop;
201 * aloop = gmx_mtop_atomloop_block_init(mtop)
202 * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
207 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
208 const t_atom **atom, int *nmol);
211 /* Abstract type for ilist loop over all ilists */
212 typedef struct gmx_mtop_ilistloop *gmx_mtop_ilistloop_t;
214 /* Initialize an ilist loop over all molecule types in the system. */
216 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop);
218 /* Initialize an ilist loop over all molecule types in the system. */
220 gmx_mtop_ilistloop_init(const gmx_mtop_t &mtop);
222 /* Loop to the next molecule,
223 * When not at the end:
224 * returns a valid pointer to the next array ilist_mol[F_NRE],
225 * writes the number of molecules for this ilist in *nmol.
226 * When at the end, destroys iloop and returns nullptr.
228 const InteractionLists *
229 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
232 /* Abstract type for ilist loop over all ilists of all molecules */
233 typedef struct gmx_mtop_ilistloop_all *gmx_mtop_ilistloop_all_t;
235 /* Initialize an ilist loop over all molecule types in the system.
236 * Only use this when you really need to loop over all molecules,
237 * i.e. when you use groups which might differ per molecule,
238 * otherwise use gmx_mtop_ilistloop.
240 gmx_mtop_ilistloop_all_t
241 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop);
243 /* Loop to the next molecule,
244 * When not at the end:
245 * returns a valid pointer to the next array ilist_mol[F_NRE],
246 * writes the atom offset which should be added to iatoms in atnr_offset.
247 * When at the end, destroys iloop and returns nullptr.
249 const InteractionLists *
250 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
254 /* Returns the total number of interactions in the system of type ftype */
256 gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype);
258 /* Returns the total number of interactions in the system of type ftype */
260 gmx_mtop_ftype_count(const gmx_mtop_t &mtop, int ftype);
262 /* Returns a charge group index for the whole system */
264 gmx_mtop_global_cgs(const gmx_mtop_t *mtop);
267 /* Returns a single t_atoms struct for the whole system */
269 gmx_mtop_global_atoms(const gmx_mtop_t *mtop);
273 * Populate a 'local' topology for the whole system.
275 * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
276 * be sorted to the end.
278 * \param[in] mtop The global topology used to populate the local one.
279 * \param[in,out] top New local topology populated from global \p mtop.
280 * \param[in] freeEnergyInteractionsAtEnd If free energy interactions will be sorted.
283 gmx_mtop_generate_local_top(const gmx_mtop_t &mtop,
285 bool freeEnergyInteractionsAtEnd);
288 /*!\brief Creates and returns a struct with begin/end atom indices of all molecules
290 * \param[in] mtop The global topology
291 * \returns A RangePartitioning object with numBlocks() equal to the number
292 * of molecules and atom indices such that molecule m contains atoms a with:
293 * index[m] <= a < index[m+1].
295 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t &mtop);
298 /* Converts a gmx_mtop_t struct to t_topology.
300 * If the lifetime of the returned topology should be longer than that
301 * of mtop, your need to pass freeMtop==true.
302 * If freeMTop == true, memory related to mtop will be freed so that done_top()
303 * on the result value will free all memory.
304 * If freeMTop == false, mtop and the return value will share some of their
305 * memory, and there is currently no way to consistently free all the memory.
308 gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop, bool freeMTop);
310 /*! \brief Get vector of atoms indices from topology
312 * This function returns the indices of all particles with type
313 * eptAtom, that is shells, vsites etc. are left out.
314 * \param[in] mtop Molecular topology
315 * \returns Vector that will be filled with the atom indices
317 std::vector<int> get_atom_index(const gmx_mtop_t *mtop);
319 /*! \brief Converts a t_atoms struct to an mtop struct
321 * All pointers contained in \p atoms will be copied into \p mtop.
322 * Note that this will produce one moleculetype encompassing the whole system.
324 * \param[in] symtab The symbol table
325 * \param[in] name Pointer to the name for the topology
326 * \param[in] atoms The atoms to convert
327 * \param[out] mtop The molecular topology output containing atoms.
330 convertAtomsToMtop(t_symtab *symtab,