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