9d7330226a07cb119c9cc6c85f39b6cfbe18c580
[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/basedefinitions.h"
50 #include "gromacs/utility/enumerationhelpers.h"
51
52 struct gmx_localtop_t;
53 struct t_atom;
54 struct t_atoms;
55 struct t_block;
56 struct t_symtab;
57
58 // TODO All of the functions taking a const gmx_mtop * are deprecated
59 // and should be replaced by versions taking const gmx_mtop & when
60 // their callers are refactored similarly.
61
62 /* Counts the number of atoms of each type. State should be 0 for
63  * state A and 1 for state B types.  typecount should have at
64  * least mtop->ffparams.atnr elements.
65  */
66 void gmx_mtop_count_atomtypes(const gmx_mtop_t& mtop, int state, int typecount[]);
67
68 /*!\brief Returns the total number of molecules in mtop
69  *
70  * \param[in] mtop  The global topology
71  */
72 int gmx_mtop_num_molecules(const gmx_mtop_t& mtop);
73
74 /* Returns the total number of residues in mtop. */
75 int gmx_mtop_nres(const gmx_mtop_t& mtop);
76
77 class AtomIterator;
78
79 //! Proxy object returned from AtomIterator
80 class AtomProxy
81 {
82 public:
83     //! Default constructor.
84     AtomProxy(const AtomIterator* it) : it_(it) {}
85     //! Access current global atom number.
86     int globalAtomNumber() const;
87     //! Access current t_atom struct.
88     const t_atom& atom() const;
89     //! Access current name of the atom.
90     const char* atomName() const;
91     //! Access current name of the residue the atom is in.
92     const char* residueName() const;
93     //! Access current residue number.
94     int residueNumber() const;
95     //! Access current molecule type.
96     const gmx_moltype_t& moleculeType() const;
97     //! Access the position of the current atom in the molecule.
98     int atomNumberInMol() const;
99
100 private:
101     const AtomIterator* it_;
102 };
103
104 /*! \brief
105  * Object that allows looping over all atoms in an mtop.
106  */
107 class AtomIterator :
108     public boost::stl_interfaces::proxy_iterator_interface<AtomIterator, std::forward_iterator_tag, t_atom, AtomProxy>
109 {
110     using Base =
111             boost::stl_interfaces::proxy_iterator_interface<AtomIterator, std::forward_iterator_tag, t_atom, AtomProxy>;
112
113 public:
114     //! Construct from topology and optionalally a global atom number.
115     explicit AtomIterator(const gmx_mtop_t& mtop, int globalAtomNumber = 0);
116
117     //! Prefix increment.
118     AtomIterator& operator++();
119     using Base::  operator++;
120
121     //! Equality comparison.
122     bool operator==(const AtomIterator& o) const;
123
124     //! Dereference operator. Returns proxy.
125     AtomProxy operator*() const { return { this }; }
126
127 private:
128     //! Global topology.
129     const gmx_mtop_t* mtop_;
130     //! Current molecule block.
131     size_t mblock_;
132     //! The atoms of the current molecule.
133     const t_atoms* atoms_;
134     //! The current molecule.
135     int currentMolecule_;
136     //! Current highest number for residues.
137     int highestResidueNumber_;
138     //! Current local atom number.
139     int localAtomNumber_;
140     //! Global current atom number.
141     int globalAtomNumber_;
142
143     friend class AtomProxy;
144 };
145
146 //! Range over all atoms of topology.
147 class AtomRange
148 {
149 public:
150     //! Default constructor.
151     explicit AtomRange(const gmx_mtop_t& mtop) : begin_(mtop), end_(mtop, mtop.natoms) {}
152     //! Iterator to begin of range.
153     AtomIterator& begin() { return begin_; }
154     //! Iterator to end of range.
155     AtomIterator& end() { return end_; }
156
157 private:
158     AtomIterator begin_, end_;
159 };
160
161 /* Abstract type for atom loop over atoms in all molecule blocks */
162 typedef struct gmx_mtop_atomloop_block* gmx_mtop_atomloop_block_t;
163
164 /* Initialize an atom loop over atoms in all molecule blocks the system.
165  */
166 gmx_mtop_atomloop_block_t gmx_mtop_atomloop_block_init(const gmx_mtop_t& mtop);
167
168 /* Loop to the next atom.
169  * When not at the end:
170  *   returns TRUE
171  *   sets the pointer atom to the t_atom struct of that atom
172  *   and return the number of molecules corresponding to this atom.
173  * When at the end, destroys aloop and returns FALSE.
174  * Use as:
175  * gmx_mtop_atomloop_block_t aloop;
176  * aloop = gmx_mtop_atomloop_block_init(mtop)
177  * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
178  *     ...
179  * }
180  */
181 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop, const t_atom** atom, int* nmol);
182
183
184 /* Abstract type for ilist loop over all ilists */
185 typedef struct gmx_mtop_ilistloop* gmx_mtop_ilistloop_t;
186
187 /* Initialize an ilist loop over all molecule types in the system. */
188 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t& mtop);
189
190 /* Loop to the next molecule,
191  * When not at the end:
192  *   returns a valid pointer to the next array ilist_mol[F_NRE],
193  *   writes the number of molecules for this ilist in *nmol.
194  * When at the end, destroys iloop and returns nullptr.
195  */
196 const InteractionLists* gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop, int* nmol);
197
198 /* Returns the total number of interactions in the system of type ftype */
199 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype);
200
201 /* Returns the total number of interactions in the system with all interaction flags that are set in \p if_flags set */
202 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, int unsigned if_flags);
203
204 /* Returns the count of atoms for each particle type */
205 gmx::EnumerationArray<ParticleType, int> gmx_mtop_particletype_count(const gmx_mtop_t& mtop);
206
207 /* Returns a single t_atoms struct for the whole system */
208 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t& mtop);
209
210
211 /*! \brief
212  * Populate a 'local' topology for the whole system.
213  *
214  * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
215  * be sorted to the end.
216  *
217  * \param[in]     mtop                        The global topology used to populate the local one.
218  * \param[in,out] top                         New local topology populated from global \p mtop.
219  * \param[in]     freeEnergyInteractionsAtEnd If free energy interactions will be sorted.
220  */
221 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd);
222
223
224 /*!\brief Creates and returns a struct with begin/end atom indices of all molecules
225  *
226  * \param[in] mtop  The global topology
227  * \returns A RangePartitioning object with numBlocks() equal to the number
228  * of molecules and atom indices such that molecule m contains atoms a with:
229  * index[m] <= a < index[m+1].
230  */
231 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop);
232
233 /*! \brief
234  * Returns the index range from residue begin to end for each residue in a molecule block.
235  *
236  * Note that residues will always have consecutive atoms numbers internally.
237  *
238  * \param[in] moltype  Molecule Type to parse for start and end.
239  * \returns Vector of ranges for all residues.
240  */
241 std::vector<gmx::Range<int>> atomRangeOfEachResidue(const gmx_moltype_t& moltype);
242
243 /* Converts a gmx_mtop_t struct to t_topology.
244  *
245  * If the lifetime of the returned topology should be longer than that
246  * of mtop, your need to pass freeMtop==true.
247  * If freeMTop == true, memory related to mtop will be freed so that done_top()
248  * on the result value will free all memory.
249  * If freeMTop == false, mtop and the return value will share some of their
250  * memory, and there is currently no way to consistently free all the memory.
251  */
252 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t* mtop, bool freeMTop);
253
254 /*! \brief Get vector of atoms indices from topology
255  *
256  * This function returns the indices of all particles with type
257  * eptAtom, that is shells, vsites etc. are left out.
258  * \param[in]  mtop Molecular topology
259  * \returns Vector that will be filled with the atom indices
260  */
261 std::vector<int> get_atom_index(const gmx_mtop_t& mtop);
262
263 /*! \brief Converts a t_atoms struct to an mtop struct
264  *
265  * All pointers contained in \p atoms will be copied into \p mtop.
266  * Note that this will produce one moleculetype encompassing the whole system.
267  *
268  * \param[in]  symtab  The symbol table
269  * \param[in]  name    Pointer to the name for the topology
270  * \param[in]  atoms   The atoms to convert
271  * \param[out] mtop    The molecular topology output containing atoms.
272  */
273 void convertAtomsToMtop(t_symtab* symtab, char** name, t_atoms* atoms, gmx_mtop_t* mtop);
274
275 //! Checks and returns whether non-bonded interactions are perturbed for free-energy calculations
276 bool haveFepPerturbedNBInteractions(const gmx_mtop_t& mtop);
277
278 //! Checks whether masses are perturbed for free-energy calculations
279 bool haveFepPerturbedMasses(const gmx_mtop_t& mtop);
280
281 //! Checks whether constraints are perturbed for free-energy calculations
282 bool havePerturbedConstraints(const gmx_mtop_t& mtop);
283
284 #endif