0884c514dca465fe7d308946de8882bf562c5b93
[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,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef GMX_TOPOLOGY_MTOP_UTIL_H
38 #define GMX_TOPOLOGY_MTOP_UTIL_H
39
40 #include <cstddef>
41
42 #include <vector>
43
44 #include "gromacs/topology/topology.h"
45 #include "gromacs/utility/basedefinitions.h"
46
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;
53
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.
57
58 /* Should be called after generating or reading mtop,
59  * to set some compute intesive variables to avoid
60  * N^2 operations later on.
61  */
62 void gmx_mtop_finalize(gmx_mtop_t* mtop);
63
64 /* Counts the number of atoms of each type. State should be 0 for
65  * state A and 1 for state B types.  typecount should have at
66  * least mtop->ffparams.atnr elements.
67  */
68 void gmx_mtop_count_atomtypes(const gmx_mtop_t* mtop, int state, int typecount[]);
69
70 /*!\brief Returns the total number of molecules in mtop
71  *
72  * \param[in] mtop  The global topology
73  */
74 int gmx_mtop_num_molecules(const gmx_mtop_t& mtop);
75
76 /* Returns the total number of residues in mtop. */
77 int gmx_mtop_nres(const gmx_mtop_t* mtop);
78
79 class AtomIterator;
80
81 //! Proxy object returned from AtomIterator
82 class AtomProxy
83 {
84 public:
85     //! Default constructor.
86     AtomProxy(const AtomIterator* it) : it_(it) {}
87     //! Access current global atom number.
88     int globalAtomNumber() const;
89     //! Access current t_atom struct.
90     const t_atom& atom() const;
91     //! Access current name of the atom.
92     const char* atomName() const;
93     //! Access current name of the residue the atom is in.
94     const char* residueName() const;
95     //! Access current residue number.
96     int residueNumber() const;
97     //! Access current molecule type.
98     const gmx_moltype_t& moleculeType() const;
99     //! Access the position of the current atom in the molecule.
100     int atomNumberInMol() const;
101
102 private:
103     const AtomIterator* it_;
104 };
105
106 //! Wrapper around proxy object to implement operator->
107 template<typename T>
108 class ProxyPtr
109 {
110 public:
111     //! Construct with proxy object.
112     ProxyPtr(T t) : t_(t) {}
113     //! Member of pointer operator.
114     T* operator->() { return &t_; }
115
116 private:
117     T t_;
118 };
119
120 /*! \brief
121  * Object that allows looping over all atoms in an mtop.
122  */
123 class AtomIterator
124 {
125 public:
126     //! Construct from topology and optionalally a global atom number.
127     explicit AtomIterator(const gmx_mtop_t& mtop, int globalAtomNumber = 0);
128
129     //! Prefix increment.
130     AtomIterator& operator++();
131     //! Postfix increment.
132     AtomIterator operator++(int);
133
134     //! Equality comparison.
135     bool operator==(const AtomIterator& o) const;
136     //! Non-equal comparison.
137     bool operator!=(const AtomIterator& o) const;
138
139     //! Dereference operator. Returns proxy.
140     AtomProxy operator*() const { return { this }; }
141     //! Member of pointer operator.
142     ProxyPtr<AtomProxy> operator->() const { return { this }; }
143
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_;
159
160     friend class AtomProxy;
161 };
162
163 //! Range over all atoms of topology.
164 class AtomRange
165 {
166 public:
167     //! Default constructor.
168     explicit AtomRange(const gmx_mtop_t& mtop) : begin_(mtop), end_(mtop, mtop.natoms) {}
169     //! Iterator to begin of range.
170     AtomIterator& begin() { return begin_; }
171     //! Iterator to end of range.
172     AtomIterator& end() { return end_; }
173
174 private:
175     AtomIterator begin_, end_;
176 };
177
178 /* Abstract type for atom loop over atoms in all molecule blocks */
179 typedef struct gmx_mtop_atomloop_block* gmx_mtop_atomloop_block_t;
180
181 /* Initialize an atom loop over atoms in all molecule blocks the system.
182  */
183 gmx_mtop_atomloop_block_t gmx_mtop_atomloop_block_init(const gmx_mtop_t* mtop);
184
185 /* Loop to the next atom.
186  * When not at the end:
187  *   returns TRUE
188  *   sets the pointer atom to the t_atom struct of that atom
189  *   and return the number of molecules corresponding to this atom.
190  * When at the end, destroys aloop and returns FALSE.
191  * Use as:
192  * gmx_mtop_atomloop_block_t aloop;
193  * aloop = gmx_mtop_atomloop_block_init(mtop)
194  * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
195  *     ...
196  * }
197  */
198 gmx_bool gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop, const t_atom** atom, int* nmol);
199
200
201 /* Abstract type for ilist loop over all ilists */
202 typedef struct gmx_mtop_ilistloop* gmx_mtop_ilistloop_t;
203
204 /* Initialize an ilist loop over all molecule types in the system. */
205 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t* mtop);
206
207 /* Initialize an ilist loop over all molecule types in the system. */
208 gmx_mtop_ilistloop_t gmx_mtop_ilistloop_init(const gmx_mtop_t& mtop);
209
210 /* Loop to the next molecule,
211  * When not at the end:
212  *   returns a valid pointer to the next array ilist_mol[F_NRE],
213  *   writes the number of molecules for this ilist in *nmol.
214  * When at the end, destroys iloop and returns nullptr.
215  */
216 const InteractionLists* gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop, int* nmol);
217
218 /* Abstract type for ilist loop over all ilists of all molecules */
219 typedef struct gmx_mtop_ilistloop_all* gmx_mtop_ilistloop_all_t;
220
221 /* Initialize an ilist loop over all molecule types in the system.
222  * Only use this when you really need to loop over all molecules,
223  * i.e. when you use groups which might differ per molecule,
224  * otherwise use gmx_mtop_ilistloop.
225  */
226 gmx_mtop_ilistloop_all_t gmx_mtop_ilistloop_all_init(const gmx_mtop_t* mtop);
227
228 /* Loop to the next molecule,
229  * When not at the end:
230  *   returns a valid pointer to the next array ilist_mol[F_NRE],
231  *   writes the atom offset which should be added to iatoms in atnr_offset.
232  * When at the end, destroys iloop and returns nullptr.
233  */
234 const InteractionLists* gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop, int* atnr_offset);
235
236
237 /* Returns the total number of interactions in the system of type ftype */
238 int gmx_mtop_ftype_count(const gmx_mtop_t* mtop, int ftype);
239
240 /* Returns the total number of interactions in the system of type ftype */
241 int gmx_mtop_ftype_count(const gmx_mtop_t& mtop, int ftype);
242
243 /* Returns the total number of interactions in the system with all interaction flags that are set in \p if_flags set */
244 int gmx_mtop_interaction_count(const gmx_mtop_t& mtop, int unsigned if_flags);
245
246 /* Returns a single t_atoms struct for the whole system */
247 t_atoms gmx_mtop_global_atoms(const gmx_mtop_t* mtop);
248
249
250 /*! \brief
251  * Populate a 'local' topology for the whole system.
252  *
253  * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
254  * be sorted to the end.
255  *
256  * \param[in]     mtop                        The global topology used to populate the local one.
257  * \param[in,out] top                         New local topology populated from global \p mtop.
258  * \param[in]     freeEnergyInteractionsAtEnd If free energy interactions will be sorted.
259  */
260 void gmx_mtop_generate_local_top(const gmx_mtop_t& mtop, gmx_localtop_t* top, bool freeEnergyInteractionsAtEnd);
261
262
263 /*!\brief Creates and returns a struct with begin/end atom indices of all molecules
264  *
265  * \param[in] mtop  The global topology
266  * \returns A RangePartitioning object with numBlocks() equal to the number
267  * of molecules and atom indices such that molecule m contains atoms a with:
268  * index[m] <= a < index[m+1].
269  */
270 gmx::RangePartitioning gmx_mtop_molecules(const gmx_mtop_t& mtop);
271
272
273 /* Converts a gmx_mtop_t struct to t_topology.
274  *
275  * If the lifetime of the returned topology should be longer than that
276  * of mtop, your need to pass freeMtop==true.
277  * If freeMTop == true, memory related to mtop will be freed so that done_top()
278  * on the result value will free all memory.
279  * If freeMTop == false, mtop and the return value will share some of their
280  * memory, and there is currently no way to consistently free all the memory.
281  */
282 t_topology gmx_mtop_t_to_t_topology(gmx_mtop_t* mtop, bool freeMTop);
283
284 /*! \brief Get vector of atoms indices from topology
285  *
286  * This function returns the indices of all particles with type
287  * eptAtom, that is shells, vsites etc. are left out.
288  * \param[in]  mtop Molecular topology
289  * \returns Vector that will be filled with the atom indices
290  */
291 std::vector<int> get_atom_index(const gmx_mtop_t* mtop);
292
293 /*! \brief Converts a t_atoms struct to an mtop struct
294  *
295  * All pointers contained in \p atoms will be copied into \p mtop.
296  * Note that this will produce one moleculetype encompassing the whole system.
297  *
298  * \param[in]  symtab  The symbol table
299  * \param[in]  name    Pointer to the name for the topology
300  * \param[in]  atoms   The atoms to convert
301  * \param[out] mtop    The molecular topology output containing atoms.
302  */
303 void convertAtomsToMtop(t_symtab* symtab, char** name, t_atoms* atoms, gmx_mtop_t* mtop);
304
305 #endif