Merge branch release-2021 into merge-2021-into-master
[alexxy/gromacs.git] / src / gromacs / topology / mtop_util.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2021, 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #ifndef GMX_TOPOLOGY_MTOP_UTIL_H
39 #define GMX_TOPOLOGY_MTOP_UTIL_H
40
41 #include <cstddef>
42
43 #include <array>
44 #include <vector>
45
46 #include <boost/stl_interfaces/iterator_interface.hpp>
47
48 #include "gromacs/topology/topology.h"
49 #include "gromacs/utility/enumerationhelpers.h"
50
51 struct gmx_localtop_t;
52 struct t_atom;
53 struct t_atoms;
54 struct t_block;
55 struct t_symtab;
56
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.
60
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.
64  */
65 void gmx_mtop_count_atomtypes(const gmx_mtop_t& mtop, int state, int typecount[]);
66
67 /*!\brief Returns the total number of molecules in mtop
68  *
69  * \param[in] mtop  The global topology
70  */
71 int gmx_mtop_num_molecules(const gmx_mtop_t& mtop);
72
73 /* Returns the total number of residues in mtop. */
74 int gmx_mtop_nres(const gmx_mtop_t& mtop);
75
76 class AtomIterator;
77
78 //! Proxy object returned from AtomIterator
79 class AtomProxy
80 {
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;
98
99 private:
100     const AtomIterator* it_;
101 };
102
103 /*! \brief
104  * Object that allows looping over all atoms in an mtop.
105  */
106 class AtomIterator :
107     public boost::stl_interfaces::proxy_iterator_interface<AtomIterator, std::forward_iterator_tag, t_atom, AtomProxy>
108 {
109     using Base =
110             boost::stl_interfaces::proxy_iterator_interface<AtomIterator, std::forward_iterator_tag, t_atom, AtomProxy>;
111
112 public:
113     //! Construct from topology and optionalally a global atom number.
114     explicit AtomIterator(const gmx_mtop_t& mtop, int globalAtomNumber = 0);
115
116     //! Prefix increment.
117     AtomIterator& operator++();
118     using Base::  operator++;
119
120     //! Equality comparison.
121     bool operator==(const AtomIterator& o) const;
122
123     //! Dereference operator. Returns proxy.
124     AtomProxy operator*() const { return { this }; }
125
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_;
141
142     friend class AtomProxy;
143 };
144
145 //! Range over all atoms of topology.
146 class AtomRange
147 {
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_; }
155
156 private:
157     AtomIterator begin_, end_;
158 };
159
160 class IListIterator;
161
162 //! Proxy object returned from IListIterator
163 class IListProxy
164 {
165 public:
166     //! Default constructor.
167     IListProxy(const IListIterator* it) : it_(it) {}
168     //! Access current global atom number.
169     const InteractionLists& list() const;
170     //! Access current molecule.
171     int nmol() const;
172
173 private:
174     const IListIterator* it_;
175 };
176
177 /*! \brief
178  * Object that allows looping over all atoms in an mtop.
179  */
180 class IListIterator :
181     public boost::stl_interfaces::proxy_iterator_interface<IListIterator, std::forward_iterator_tag, InteractionLists, IListProxy>
182 {
183     using Base =
184             boost::stl_interfaces::proxy_iterator_interface<IListIterator, std::forward_iterator_tag, InteractionLists, IListProxy>;
185
186 public:
187     //! Construct from topology.
188     explicit IListIterator(const gmx_mtop_t& mtop, size_t mblock = 0);
189
190     //! Prefix increment.
191     IListIterator& operator++();
192     using Base::   operator++;
193
194     //! Equality comparison.
195     bool operator==(const IListIterator& o) const;
196
197     //! Dereference operator. Returns proxy.
198     IListProxy operator*() const { return { this }; }
199
200 private:
201     //! Global topology.
202     const gmx_mtop_t* mtop_;
203     //! Index of molecule block corresponding to the current location.
204     size_t mblock_;
205
206     friend class IListProxy;
207 };
208
209
210 /*! \brief
211  * Range over all interaction lists of topology.
212  *
213  * Includes the intermolecular interactions as the final element in the
214  * range if present.
215  */
216 class IListRange
217 {
218 public:
219     //! Default constructor.
220     explicit IListRange(const gmx_mtop_t& mtop);
221     //! Iterator to begin of range.
222     IListIterator& begin() { return begin_; }
223     //! Iterator to end of range.
224     IListIterator& end() { return end_; }
225
226 private:
227     IListIterator begin_, end_;
228 };
229
230 /* Abstract type for atom loop over atoms in all molecule blocks */
231 typedef struct gmx_mtop_atomloop_block* gmx_mtop_atomloop_block_t;
232
233 /* Initialize an atom loop over atoms in all molecule blocks the system.
234  */
235 gmx_mtop_atomloop_block_t gmx_mtop_atomloop_block_init(const gmx_mtop_t& mtop);
236
237 /* Loop to the next atom.
238  * When not at the end:
239  *   returns TRUE
240  *   sets the pointer atom to the t_atom struct of that atom
241  *   and return the number of molecules corresponding to this atom.
242  * When at the end, destroys aloop and returns FALSE.
243  * Use as:
244  * gmx_mtop_atomloop_block_t aloop;
245  * aloop = gmx_mtop_atomloop_block_init(mtop)
246  * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
247  *     ...
248  * }
249  */
250 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop, const t_atom** atom, int* nmol);
251
252
253 /* Returns the total number of interactions in the system of type ftype */
254 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype);
255
256 /* Returns the total number of interactions in the system with all interaction flags that are set in \p if_flags set */
257 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, int unsigned if_flags);
258
259 /* Returns the count of atoms for each particle type */
260 gmx::EnumerationArray<ParticleType, int> gmx_mtop_particletype_count(const gmx_mtop_t& mtop);
261
262 /* Returns a single t_atoms struct for the whole system */
263 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t& mtop);
264
265
266 /*! \brief
267  * Populate a 'local' topology for the whole system.
268  *
269  * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
270  * be sorted to the end.
271  *
272  * \param[in]     mtop                        The global topology used to populate the local one.
273  * \param[in,out] top                         New local topology populated from global \p mtop.
274  * \param[in]     freeEnergyInteractionsAtEnd If free energy interactions will be sorted.
275  */
276 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd);
277
278
279 /*!\brief Creates and returns a struct with begin/end atom indices of all molecules
280  *
281  * \param[in] mtop  The global topology
282  * \returns A RangePartitioning object with numBlocks() equal to the number
283  * of molecules and atom indices such that molecule m contains atoms a with:
284  * index[m] <= a < index[m+1].
285  */
286 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop);
287
288 /*! \brief
289  * Returns the index range from residue begin to end for each residue in a molecule block.
290  *
291  * Note that residues will always have consecutive atoms numbers internally.
292  *
293  * \param[in] moltype  Molecule Type to parse for start and end.
294  * \returns Vector of ranges for all residues.
295  */
296 std::vector<gmx::Range<int>> atomRangeOfEachResidue(const gmx_moltype_t& moltype);
297
298 /* Converts a gmx_mtop_t struct to t_topology.
299  *
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.
306  */
307 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t* mtop, bool freeMTop);
308
309 /*! \brief Get vector of atoms indices from topology
310  *
311  * This function returns the indices of all particles with type
312  * eptAtom, that is shells, vsites etc. are left out.
313  * \param[in]  mtop Molecular topology
314  * \returns Vector that will be filled with the atom indices
315  */
316 std::vector<int> get_atom_index(const gmx_mtop_t& mtop);
317
318 /*! \brief Converts a t_atoms struct to an mtop struct
319  *
320  * All pointers contained in \p atoms will be copied into \p mtop.
321  * Note that this will produce one moleculetype encompassing the whole system.
322  *
323  * \param[in]  symtab  The symbol table
324  * \param[in]  name    Pointer to the name for the topology
325  * \param[in]  atoms   The atoms to convert
326  * \param[out] mtop    The molecular topology output containing atoms.
327  */
328 void convertAtomsToMtop(t_symtab* symtab, char** name, t_atoms* atoms, gmx_mtop_t* mtop);
329
330 //! Checks and returns whether non-bonded interactions are perturbed for free-energy calculations
331 bool haveFepPerturbedNBInteractions(const gmx_mtop_t& mtop);
332
333 //! Checks whether masses are perturbed for free-energy calculations
334 bool haveFepPerturbedMasses(const gmx_mtop_t& mtop);
335
336 //! Checks whether masses are perturbed for free-energy calculations in SETTLE interactions
337 bool haveFepPerturbedMassesInSettles(const gmx_mtop_t& mtop);
338
339 //! Checks whether constraints are perturbed for free-energy calculations
340 bool havePerturbedConstraints(const gmx_mtop_t& mtop);
341
342 #endif