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