Update instructions in containers.rst
[gromacs.git] / src / gromacs / topology / mtop_util.h
blob7d0527e654cfa4f844f39acb30ed06111e7b1a91
1 /*
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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef GMX_TOPOLOGY_MTOP_UTIL_H
39 #define GMX_TOPOLOGY_MTOP_UTIL_H
41 #include <cstddef>
43 #include <array>
44 #include <vector>
46 #include <boost/stl_interfaces/iterator_interface.hpp>
48 #include "gromacs/topology/topology.h"
49 #include "gromacs/utility/basedefinitions.h"
51 struct gmx_localtop_t;
52 struct t_atom;
53 struct t_atoms;
54 struct t_block;
55 struct t_symtab;
57 // TODO All of the functions taking a const gmx_mtop * are deprecated
58 // and should be replaced by versions taking const gmx_mtop & when
59 // their callers are refactored similarly.
61 /* Counts the number of atoms of each type. State should be 0 for
62 * state A and 1 for state B types. typecount should have at
63 * least mtop->ffparams.atnr elements.
65 void gmx_mtop_count_atomtypes(const gmx_mtop_t* mtop, int state, int typecount[]);
67 /*!\brief Returns the total number of molecules in mtop
69 * \param[in] mtop The global topology
71 int gmx_mtop_num_molecules(const gmx_mtop_t& mtop);
73 /* Returns the total number of residues in mtop. */
74 int gmx_mtop_nres(const gmx_mtop_t* mtop);
76 class AtomIterator;
78 //! Proxy object returned from AtomIterator
79 class AtomProxy
81 public:
82 //! Default constructor.
83 AtomProxy(const AtomIterator* it) : it_(it) {}
84 //! Access current global atom number.
85 int globalAtomNumber() const;
86 //! Access current t_atom struct.
87 const t_atom& atom() const;
88 //! Access current name of the atom.
89 const char* atomName() const;
90 //! Access current name of the residue the atom is in.
91 const char* residueName() const;
92 //! Access current residue number.
93 int residueNumber() const;
94 //! Access current molecule type.
95 const gmx_moltype_t& moleculeType() const;
96 //! Access the position of the current atom in the molecule.
97 int atomNumberInMol() const;
99 private:
100 const AtomIterator* it_;
103 /*! \brief
104 * Object that allows looping over all atoms in an mtop.
106 class AtomIterator :
107 public boost::stl_interfaces::proxy_iterator_interface<AtomIterator, std::forward_iterator_tag, t_atom, AtomProxy>
109 using Base =
110 boost::stl_interfaces::proxy_iterator_interface<AtomIterator, std::forward_iterator_tag, t_atom, AtomProxy>;
112 public:
113 //! Construct from topology and optionalally a global atom number.
114 explicit AtomIterator(const gmx_mtop_t& mtop, int globalAtomNumber = 0);
116 //! Prefix increment.
117 AtomIterator& operator++();
118 using Base:: operator++;
120 //! Equality comparison.
121 bool operator==(const AtomIterator& o) const;
123 //! Dereference operator. Returns proxy.
124 AtomProxy operator*() const { return { this }; }
126 private:
127 //! Global topology.
128 const gmx_mtop_t* mtop_;
129 //! Current molecule block.
130 size_t mblock_;
131 //! The atoms of the current molecule.
132 const t_atoms* atoms_;
133 //! The current molecule.
134 int currentMolecule_;
135 //! Current highest number for residues.
136 int highestResidueNumber_;
137 //! Current local atom number.
138 int localAtomNumber_;
139 //! Global current atom number.
140 int globalAtomNumber_;
142 friend class AtomProxy;
145 //! Range over all atoms of topology.
146 class AtomRange
148 public:
149 //! Default constructor.
150 explicit AtomRange(const gmx_mtop_t& mtop) : begin_(mtop), end_(mtop, mtop.natoms) {}
151 //! Iterator to begin of range.
152 AtomIterator& begin() { return begin_; }
153 //! Iterator to end of range.
154 AtomIterator& end() { return end_; }
156 private:
157 AtomIterator begin_, end_;
160 /* Abstract type for atom loop over atoms in all molecule blocks */
161 typedef struct gmx_mtop_atomloop_block* gmx_mtop_atomloop_block_t;
163 /* Initialize an atom loop over atoms in all molecule blocks the system.
165 gmx_mtop_atomloop_block_t gmx_mtop_atomloop_block_init(const gmx_mtop_t* mtop);
167 /* Loop to the next atom.
168 * When not at the end:
169 * returns TRUE
170 * sets the pointer atom to the t_atom struct of that atom
171 * and return the number of molecules corresponding to this atom.
172 * When at the end, destroys aloop and returns FALSE.
173 * Use as:
174 * gmx_mtop_atomloop_block_t aloop;
175 * aloop = gmx_mtop_atomloop_block_init(mtop)
176 * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
177 * ...
180 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop, const t_atom** atom, int* nmol);
183 /* Abstract type for ilist loop over all ilists */
184 typedef struct gmx_mtop_ilistloop* gmx_mtop_ilistloop_t;
186 /* Initialize an ilist loop over all molecule types in the system. */
187 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t* mtop);
189 /* Initialize an ilist loop over all molecule types in the system. */
190 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t& mtop);
192 /* Loop to the next molecule,
193 * When not at the end:
194 * returns a valid pointer to the next array ilist_mol[F_NRE],
195 * writes the number of molecules for this ilist in *nmol.
196 * When at the end, destroys iloop and returns nullptr.
198 const InteractionLists* gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop, int* nmol);
200 /* Returns the total number of interactions in the system of type ftype */
201 int gmx_mtop_ftype_count(const gmx_mtop_t* mtop, int ftype);
203 /* Returns the total number of interactions in the system of type ftype */
204 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype);
206 /* Returns the total number of interactions in the system with all interaction flags that are set in \p if_flags set */
207 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, int unsigned if_flags);
209 /* Returns the count of atoms for each particle type */
210 std::array<int, eptNR> gmx_mtop_particletype_count(const gmx_mtop_t& mtop);
212 /* Returns a single t_atoms struct for the whole system */
213 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t* mtop);
216 /*! \brief
217 * Populate a 'local' topology for the whole system.
219 * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
220 * be sorted to the end.
222 * \param[in] mtop The global topology used to populate the local one.
223 * \param[in,out] top New local topology populated from global \p mtop.
224 * \param[in] freeEnergyInteractionsAtEnd If free energy interactions will be sorted.
226 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd);
229 /*!\brief Creates and returns a struct with begin/end atom indices of all molecules
231 * \param[in] mtop The global topology
232 * \returns A RangePartitioning object with numBlocks() equal to the number
233 * of molecules and atom indices such that molecule m contains atoms a with:
234 * index[m] <= a < index[m+1].
236 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop);
238 /*! \brief
239 * Returns the index range from residue begin to end for each residue in a molecule block.
241 * Note that residues will always have consecutive atoms numbers internally.
243 * \param[in] moltype Molecule Type to parse for start and end.
244 * \returns Vector of ranges for all residues.
246 std::vector<gmx::Range<int>> atomRangeOfEachResidue(const gmx_moltype_t& moltype);
248 /* Converts a gmx_mtop_t struct to t_topology.
250 * If the lifetime of the returned topology should be longer than that
251 * of mtop, your need to pass freeMtop==true.
252 * If freeMTop == true, memory related to mtop will be freed so that done_top()
253 * on the result value will free all memory.
254 * If freeMTop == false, mtop and the return value will share some of their
255 * memory, and there is currently no way to consistently free all the memory.
257 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t* mtop, bool freeMTop);
259 /*! \brief Get vector of atoms indices from topology
261 * This function returns the indices of all particles with type
262 * eptAtom, that is shells, vsites etc. are left out.
263 * \param[in] mtop Molecular topology
264 * \returns Vector that will be filled with the atom indices
266 std::vector<int> get_atom_index(const gmx_mtop_t* mtop);
268 /*! \brief Converts a t_atoms struct to an mtop struct
270 * All pointers contained in \p atoms will be copied into \p mtop.
271 * Note that this will produce one moleculetype encompassing the whole system.
273 * \param[in] symtab The symbol table
274 * \param[in] name Pointer to the name for the topology
275 * \param[in] atoms The atoms to convert
276 * \param[out] mtop The molecular topology output containing atoms.
278 void convertAtomsToMtop(t_symtab* symtab, char** name, t_atoms* atoms, gmx_mtop_t* mtop);
280 //! Checks and returns whether non-bonded interactions are perturbed for free-energy calculations
281 bool haveFepPerturbedNBInteractions(const gmx_mtop_t& mtop);
283 //! Checks whether masses are perturbed for free-energy calculations
284 bool haveFepPerturbedMasses(const gmx_mtop_t& mtop);
286 //! Checks whether constraints are perturbed for free-energy calculations
287 bool havePerturbedConstraints(const gmx_mtop_t& mtop);
289 #endif