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