Remove unused function generate_excls and make clean_excls static
[gromacs.git] / src / gromacs / topology / mtop_util.h
blob7f3aa256960aeee53f43be66520fdbee5457e489
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,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
40 #include <cstddef>
42 #include <vector>
44 #include "gromacs/topology/topology.h"
45 #include "gromacs/utility/basedefinitions.h"
47 struct gmx_localtop_t;
48 struct t_atom;
49 struct t_atoms;
50 struct t_block;
51 struct t_symtab;
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.
62 void
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.
69 void
70 gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[]);
72 /*!\brief Returns the total number of molecules in mtop
74 * \param[in] mtop The global topology
76 int gmx_mtop_num_molecules(const gmx_mtop_t &mtop);
78 /* Returns the total number of residues in mtop. */
79 int gmx_mtop_nres(const gmx_mtop_t *mtop);
81 class AtomIterator;
83 //! Proxy object returned from AtomIterator
84 class AtomProxy
86 public:
87 //! Default constructor.
88 AtomProxy(const AtomIterator* it) : it_(it) {}
89 //! Access current global atom number.
90 int globalAtomNumber() const;
91 //! Access current t_atom struct.
92 const t_atom &atom() const;
93 //! Access current name of the atom.
94 const char *atomName() const;
95 //! Access current name of the residue the atom is in.
96 const char *residueName() const;
97 //! Access current residue number.
98 int residueNumber() const;
99 //! Access current molecule type.
100 const gmx_moltype_t &moleculeType() const;
101 //! Access the position of the current atom in the molecule.
102 int atomNumberInMol() const;
103 private:
104 const AtomIterator* it_;
107 //! Wrapper around proxy object to implement operator->
108 template <typename T>
109 class ProxyPtr
111 public:
112 //! Construct with proxy object.
113 ProxyPtr(T t) : t_(t) {}
114 //! Member of pointer operator.
115 T* operator->() { return &t_; }
116 private:
117 T t_;
120 /*! \brief
121 * Object that allows looping over all atoms in an mtop.
123 class AtomIterator
125 public:
126 //! Construct from topology and optionalally a global atom number.
127 explicit AtomIterator(const gmx_mtop_t &mtop, int globalAtomNumber = 0);
129 //! Prefix increment.
130 AtomIterator &operator++();
131 //! Postfix increment.
132 AtomIterator operator++(int);
134 //! Equality comparison.
135 bool operator==(const AtomIterator &o) const;
136 //! Non-equal comparison.
137 bool operator!=(const AtomIterator &o) const;
139 //! Dereference operator. Returns proxy.
140 AtomProxy operator*() const { return {this}; }
141 //! Member of pointer operator.
142 ProxyPtr<AtomProxy> operator->() const { return {this}; }
144 private:
145 //! Global topology.
146 const gmx_mtop_t *mtop_;
147 //! Current molecule block.
148 size_t mblock_;
149 //! The atoms of the current molecule.
150 const t_atoms *atoms_;
151 //! The current molecule.
152 int currentMolecule_;
153 //! Current highest number for residues.
154 int highestResidueNumber_;
155 //! Current local atom number.
156 int localAtomNumber_;
157 //! Global current atom number.
158 int globalAtomNumber_;
160 friend class AtomProxy;
163 //! Range over all atoms of topology.
164 class AtomRange
166 public:
167 //! Default constructor.
168 explicit AtomRange(const gmx_mtop_t &mtop) :
169 begin_(mtop), end_(mtop, mtop.natoms) {}
170 //! Iterator to begin of range.
171 AtomIterator &begin() { return begin_; }
172 //! Iterator to end of range.
173 AtomIterator &end() { return end_; }
174 private:
175 AtomIterator begin_, end_;
178 /* Abstract type for atom loop over atoms in all molecule blocks */
179 typedef struct gmx_mtop_atomloop_block *gmx_mtop_atomloop_block_t;
181 /* Initialize an atom loop over atoms in all molecule blocks the system.
183 gmx_mtop_atomloop_block_t
184 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop);
186 /* Loop to the next atom.
187 * When not at the end:
188 * returns TRUE
189 * sets the pointer atom to the t_atom struct of that atom
190 * and return the number of molecules corresponding to this atom.
191 * When at the end, destroys aloop and returns FALSE.
192 * Use as:
193 * gmx_mtop_atomloop_block_t aloop;
194 * aloop = gmx_mtop_atomloop_block_init(mtop)
195 * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
196 * ...
199 gmx_bool
200 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
201 const t_atom **atom, int *nmol);
204 /* Abstract type for ilist loop over all ilists */
205 typedef struct gmx_mtop_ilistloop *gmx_mtop_ilistloop_t;
207 /* Initialize an ilist loop over all molecule types in the system. */
208 gmx_mtop_ilistloop_t
209 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop);
211 /* Initialize an ilist loop over all molecule types in the system. */
212 gmx_mtop_ilistloop_t
213 gmx_mtop_ilistloop_init(const gmx_mtop_t &mtop);
215 /* Loop to the next molecule,
216 * When not at the end:
217 * returns a valid pointer to the next array ilist_mol[F_NRE],
218 * writes the number of molecules for this ilist in *nmol.
219 * When at the end, destroys iloop and returns nullptr.
221 const InteractionLists *
222 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
223 int *nmol);
225 /* Abstract type for ilist loop over all ilists of all molecules */
226 typedef struct gmx_mtop_ilistloop_all *gmx_mtop_ilistloop_all_t;
228 /* Initialize an ilist loop over all molecule types in the system.
229 * Only use this when you really need to loop over all molecules,
230 * i.e. when you use groups which might differ per molecule,
231 * otherwise use gmx_mtop_ilistloop.
233 gmx_mtop_ilistloop_all_t
234 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop);
236 /* Loop to the next molecule,
237 * When not at the end:
238 * returns a valid pointer to the next array ilist_mol[F_NRE],
239 * writes the atom offset which should be added to iatoms in atnr_offset.
240 * When at the end, destroys iloop and returns nullptr.
242 const InteractionLists *
243 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
244 int *atnr_offset);
247 /* Returns the total number of interactions in the system of type ftype */
249 gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype);
251 /* Returns the total number of interactions in the system of type ftype */
253 gmx_mtop_ftype_count(const gmx_mtop_t &mtop, int ftype);
255 /* Returns a single t_atoms struct for the whole system */
256 t_atoms
257 gmx_mtop_global_atoms(const gmx_mtop_t *mtop);
260 /*! \brief
261 * Populate a 'local' topology for the whole system.
263 * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
264 * be sorted to the end.
266 * \param[in] mtop The global topology used to populate the local one.
267 * \param[in,out] top New local topology populated from global \p mtop.
268 * \param[in] freeEnergyInteractionsAtEnd If free energy interactions will be sorted.
270 void
271 gmx_mtop_generate_local_top(const gmx_mtop_t &mtop,
272 gmx_localtop_t *top,
273 bool freeEnergyInteractionsAtEnd);
276 /*!\brief Creates and returns a struct with begin/end atom indices of all molecules
278 * \param[in] mtop The global topology
279 * \returns A RangePartitioning object with numBlocks() equal to the number
280 * of molecules and atom indices such that molecule m contains atoms a with:
281 * index[m] <= a < index[m+1].
283 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t &mtop);
286 /* Converts a gmx_mtop_t struct to t_topology.
288 * If the lifetime of the returned topology should be longer than that
289 * of mtop, your need to pass freeMtop==true.
290 * If freeMTop == true, memory related to mtop will be freed so that done_top()
291 * on the result value will free all memory.
292 * If freeMTop == false, mtop and the return value will share some of their
293 * memory, and there is currently no way to consistently free all the memory.
295 t_topology
296 gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop, bool freeMTop);
298 /*! \brief Get vector of atoms indices from topology
300 * This function returns the indices of all particles with type
301 * eptAtom, that is shells, vsites etc. are left out.
302 * \param[in] mtop Molecular topology
303 * \returns Vector that will be filled with the atom indices
305 std::vector<int> get_atom_index(const gmx_mtop_t *mtop);
307 /*! \brief Converts a t_atoms struct to an mtop struct
309 * All pointers contained in \p atoms will be copied into \p mtop.
310 * Note that this will produce one moleculetype encompassing the whole system.
312 * \param[in] symtab The symbol table
313 * \param[in] name Pointer to the name for the topology
314 * \param[in] atoms The atoms to convert
315 * \param[out] mtop The molecular topology output containing atoms.
317 void
318 convertAtomsToMtop(t_symtab *symtab,
319 char **name,
320 t_atoms *atoms,
321 gmx_mtop_t *mtop);
323 #endif