918bb3dcfdbd96077bd0d2a89a6ee7e4cb577572
[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, 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 #include "gromacs/utility/gmxassert.h"
47
48 struct gmx_localtop_t;
49 struct t_atom;
50 struct t_atoms;
51 struct t_block;
52 struct t_ilist;
53 struct t_symtab;
54
55 /* Should be called after generating or reading mtop,
56  * to set some compute intesive variables to avoid
57  * N^2 operations later on.
58  */
59 void
60 gmx_mtop_finalize(gmx_mtop_t *mtop);
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
67 gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[]);
68
69 /* Returns the total number of charge groups in mtop */
70 int
71 ncg_mtop(const gmx_mtop_t *mtop);
72
73 /*!\brief Returns the total number of molecules in mtop
74  *
75  * \param[in] mtop  The global topology
76  */
77 int gmx_mtop_num_molecules(const gmx_mtop_t &mtop);
78
79 /* Returns the total number of residues in mtop. */
80 int gmx_mtop_nres(const gmx_mtop_t *mtop);
81
82 /* Removes the charge groups, i.e. makes single atom charge groups, in mtop */
83 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop);
84
85 /* Abstract type for atom loop over all atoms */
86 typedef struct gmx_mtop_atomloop_all *gmx_mtop_atomloop_all_t;
87
88 /* Initialize an atom loop over all atoms in the system.
89  * The order of the atoms will be as in the state struct.
90  * Only use this when you really need to loop over all atoms,
91  * i.e. when you use groups which might differ per molecule,
92  * otherwise use gmx_mtop_atomloop_block.
93  */
94 gmx_mtop_atomloop_all_t
95 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop);
96
97 /* Loop to the next atom.
98  * When not at the end:
99  *   returns TRUE and at_global,
100  *   writes the global atom number in *at_global
101  *   and sets the pointer atom to the t_atom struct of that atom.
102  * When at the end, destroys aloop and returns FALSE.
103  * Use as:
104  * gmx_mtop_atomloop_all_t aloop;
105  * aloop = gmx_mtop_atomloop_all_init(mtop)
106  * while (gmx_mtop_atomloop_all_next(aloop,&at_global,&atom)) {
107  *     ...
108  * }
109  */
110 gmx_bool
111 gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
112                            int *at_global, const t_atom **atom);
113
114 /* Return the atomname, the residue number and residue name
115  * of the current atom in the loop.
116  */
117 void
118 gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
119                             char **atomname, int *resnr, char **resname);
120
121 /* Return the a pointer to the moltype struct of the current atom
122  * in the loop and the atom number in the molecule.
123  */
124 void
125 gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
126                               gmx_moltype_t **moltype, int *at_mol);
127
128
129 /* Abstract type for atom loop over atoms in all molecule blocks */
130 typedef struct gmx_mtop_atomloop_block *gmx_mtop_atomloop_block_t;
131
132 /* Initialize an atom loop over atoms in all molecule blocks the system.
133  */
134 gmx_mtop_atomloop_block_t
135 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop);
136
137 /* Loop to the next atom.
138  * When not at the end:
139  *   returns TRUE
140  *   sets the pointer atom to the t_atom struct of that atom
141  *   and return the number of molecules corresponding to this atom.
142  * When at the end, destroys aloop and returns FALSE.
143  * Use as:
144  * gmx_mtop_atomloop_block_t aloop;
145  * aloop = gmx_mtop_atomloop_block_init(mtop)
146  * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
147  *     ...
148  * }
149  */
150 gmx_bool
151 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
152                              const t_atom **atom, int *nmol);
153
154
155 /* Abstract type for ilist loop over all ilists */
156 typedef struct gmx_mtop_ilistloop *gmx_mtop_ilistloop_t;
157
158 /* Initialize an ilist loop over all molecule types in the system. */
159 gmx_mtop_ilistloop_t
160 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop);
161
162
163 /* Loop to the next molecule,
164  * When not at the end:
165  *   returns TRUE and a pointer to the next array ilist_mol[F_NRE],
166  *   writes the number of molecules for this ilist in *nmol.
167  * When at the end, destroys iloop and returns FALSE.
168  */
169 gmx_bool
170 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
171                         t_ilist **ilist_mol, int *nmol);
172
173
174 /* Abstract type for ilist loop over all ilists of all molecules */
175 typedef struct gmx_mtop_ilistloop_all *gmx_mtop_ilistloop_all_t;
176
177 /* Initialize an ilist loop over all molecule types in the system.
178  * Only use this when you really need to loop over all molecules,
179  * i.e. when you use groups which might differ per molecule,
180  * otherwise use gmx_mtop_ilistloop.
181  */
182 gmx_mtop_ilistloop_all_t
183 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop);
184
185 /* Loop to the next molecule,
186  * When not at the end:
187  *   returns TRUE and a pointer to the next array ilist_mol[F_NRE],
188  *   writes the atom offset which should be added to iatoms in atnr_offset.
189  * When at the end, destroys iloop and returns FALSE.
190  */
191 gmx_bool
192 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
193                             t_ilist **ilist_mol, int *atnr_offset);
194
195
196 /* Returns the total number of interactions in the system of type ftype */
197 int
198 gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype);
199
200
201 /* Returns a charge group index for the whole system */
202 t_block
203 gmx_mtop_global_cgs(const gmx_mtop_t *mtop);
204
205
206 /* Returns a single t_atoms struct for the whole system */
207 t_atoms
208 gmx_mtop_global_atoms(const gmx_mtop_t *mtop);
209
210
211 /* Generate a 'local' topology for the whole system.
212  * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
213  * be sorted to the end.
214  */
215 gmx_localtop_t *
216 gmx_mtop_generate_local_top(const gmx_mtop_t *mtop, bool freeEnergyInteractionsAtEnd);
217
218
219 /*!\brief Creates and returns a struct with begin/end atom indices of all molecules
220  *
221  * \param[in] mtop  The global topology
222  * \returns A struct of type BlockRanges with numBlocks() equal to the number
223  * of molecules and atom indices such that molecule m contains atoms a with:
224  * index[m] <= a < index[m+1].
225  */
226 gmx::BlockRanges gmx_mtop_molecules(const gmx_mtop_t &mtop);
227
228
229 /* Converts a gmx_mtop_t struct to t_topology.
230  *
231  * If freeMTop == true, memory related to mtop will be freed so that done_top()
232  * on the result value will free all memory.
233  * If freeMTop == false, mtop and the return value will share some of their
234  * memory, and there is currently no way to consistently free all the memory.
235  */
236 t_topology
237 gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop, bool freeMTop);
238
239 /*! \brief Get vector of atoms indices from topology
240  *
241  * This function returns the indices of all particles with type
242  * eptAtom, that is shells, vsites etc. are left out.
243  * \param[in]  mtop Molecular topology
244  * \returns Vector that will be filled with the atom indices
245  */
246 std::vector<size_t> get_atom_index(const gmx_mtop_t *mtop);
247
248 /*! \brief Converts a t_atoms struct to an mtop struct
249  *
250  * All pointers contained in \p atoms will be copied into \p mtop.
251  * Note that this will produce one moleculetype encompassing the whole system.
252  *
253  * \param[in]  symtab  The symbol table
254  * \param[in]  name    Pointer to the name for the topology
255  * \param[in]  atoms   The atoms to convert
256  * \param[out] mtop    The molecular topology output containing atoms.
257  */
258 void
259 convertAtomsToMtop(t_symtab    *symtab,
260                    char       **name,
261                    t_atoms     *atoms,
262                    gmx_mtop_t  *mtop);
263
264 #endif