Add readConfAndTopology
[gromacs.git] / src / gromacs / topology / mtop_util.h
blob77b1569d81a37263c47230d02643c799250474b6
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
37 #ifndef GMX_TOPOLOGY_MTOP_UTIL_H
38 #define GMX_TOPOLOGY_MTOP_UTIL_H
40 #include <cstddef>
42 #include <vector>
44 #include "gromacs/utility/basedefinitions.h"
46 struct gmx_localtop_t;
47 struct gmx_moltype_t;
48 struct gmx_mtop_t;
49 struct t_atom;
50 struct t_atoms;
51 struct t_block;
52 struct t_ilist;
53 struct t_symtab;
54 struct t_topology;
56 /* Should be called after generating or reading mtop,
57 * to set some compute intesive variables to avoid
58 * N^2 operations later on.
60 void
61 gmx_mtop_finalize(gmx_mtop_t *mtop);
63 /* Counts the number of atoms of each type. State should be 0 for
64 * state A and 1 for state B types. typecount should have at
65 * least mtop->ffparams.atnr elements.
67 void
68 gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[]);
70 /* Returns the total number of charge groups in mtop */
71 int
72 ncg_mtop(const gmx_mtop_t *mtop);
74 /* Removes the charge groups, i.e. makes single atom charge groups, in mtop */
75 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop);
78 /* Abstract data type for looking up atoms by global atom number */
79 typedef struct gmx_mtop_atomlookup *gmx_mtop_atomlookup_t;
81 /* Initialize atom lookup by global atom number */
82 gmx_mtop_atomlookup_t
83 gmx_mtop_atomlookup_init(const gmx_mtop_t *mtop);
85 /* As gmx_mtop_atomlookup_init, but optimized for atoms involved in settle */
86 gmx_mtop_atomlookup_t
87 gmx_mtop_atomlookup_settle_init(const gmx_mtop_t *mtop);
89 /* Destroy a gmx_mtop_atomlookup_t data structure */
90 void
91 gmx_mtop_atomlookup_destroy(gmx_mtop_atomlookup_t alook);
94 /* Returns a pointer to the t_atom struct belonging to atnr_global.
95 * This can be an expensive operation, so if possible use
96 * one of the atom loop constructs below.
98 void
99 gmx_mtop_atomnr_to_atom(const gmx_mtop_atomlookup_t alook,
100 int atnr_global,
101 t_atom **atom);
104 /* Returns a pointer to the molecule interaction array ilist_mol[F_NRE]
105 * and the local atom number in the molecule belonging to atnr_global.
107 void
108 gmx_mtop_atomnr_to_ilist(const gmx_mtop_atomlookup_t alook,
109 int atnr_global,
110 t_ilist **ilist_mol, int *atnr_offset);
113 /* Returns the molecule block index
114 * and the molecule number in the block
115 * and the atom number offset for the atom indices in moltype
116 * belonging to atnr_global.
118 void
119 gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_atomlookup_t alook,
120 int atnr_global,
121 int *molb, int *molnr, int *atnr_mol);
124 /* Returns atom name, global resnr and residue name of atom atnr_global */
125 void
126 gmx_mtop_atominfo_global(const gmx_mtop_t *mtop, int atnr_global,
127 char **atomname, int *resnr, char **resname);
130 /* Abstract type for atom loop over all atoms */
131 typedef struct gmx_mtop_atomloop_all *gmx_mtop_atomloop_all_t;
133 /* Initialize an atom loop over all atoms in the system.
134 * The order of the atoms will be as in the state struct.
135 * Only use this when you really need to loop over all atoms,
136 * i.e. when you use groups which might differ per molecule,
137 * otherwise use gmx_mtop_atomloop_block.
139 gmx_mtop_atomloop_all_t
140 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop);
142 /* Loop to the next atom.
143 * When not at the end:
144 * returns TRUE and at_global,
145 * writes the global atom number in *at_global
146 * and sets the pointer atom to the t_atom struct of that atom.
147 * When at the end, destroys aloop and returns FALSE.
148 * Use as:
149 * gmx_mtop_atomloop_all_t aloop;
150 * aloop = gmx_mtop_atomloop_all_init(mtop)
151 * while (gmx_mtop_atomloop_all_next(aloop,&at_global,&atom)) {
152 * ...
155 gmx_bool
156 gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
157 int *at_global, t_atom **atom);
159 /* Return the atomname, the residue number and residue name
160 * of the current atom in the loop.
162 void
163 gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
164 char **atomname, int *resnr, char **resname);
166 /* Return the a pointer to the moltype struct of the current atom
167 * in the loop and the atom number in the molecule.
169 void
170 gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
171 gmx_moltype_t **moltype, int *at_mol);
174 /* Abstract type for atom loop over atoms in all molecule blocks */
175 typedef struct gmx_mtop_atomloop_block *gmx_mtop_atomloop_block_t;
177 /* Initialize an atom loop over atoms in all molecule blocks the system.
179 gmx_mtop_atomloop_block_t
180 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop);
182 /* Loop to the next atom.
183 * When not at the end:
184 * returns TRUE
185 * sets the pointer atom to the t_atom struct of that atom
186 * and return the number of molecules corresponding to this atom.
187 * When at the end, destroys aloop and returns FALSE.
188 * Use as:
189 * gmx_mtop_atomloop_block_t aloop;
190 * aloop = gmx_mtop_atomloop_block_init(mtop)
191 * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
192 * ...
195 gmx_bool
196 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
197 t_atom **atom, int *nmol);
200 /* Abstract type for ilist loop over all ilists */
201 typedef struct gmx_mtop_ilistloop *gmx_mtop_ilistloop_t;
203 /* Initialize an ilist loop over all molecule types in the system. */
204 gmx_mtop_ilistloop_t
205 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop);
208 /* Loop to the next molecule,
209 * When not at the end:
210 * returns TRUE and a pointer to the next array ilist_mol[F_NRE],
211 * writes the number of molecules for this ilist in *nmol.
212 * When at the end, destroys iloop and returns FALSE.
214 gmx_bool
215 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
216 t_ilist **ilist_mol, int *nmol);
219 /* Abstract type for ilist loop over all ilists of all molecules */
220 typedef struct gmx_mtop_ilistloop_all *gmx_mtop_ilistloop_all_t;
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.
227 gmx_mtop_ilistloop_all_t
228 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop);
230 /* Loop to the next molecule,
231 * When not at the end:
232 * returns TRUE and a 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 FALSE.
236 gmx_bool
237 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
238 t_ilist **ilist_mol, int *atnr_offset);
241 /* Returns the total number of interactions in the system of type ftype */
243 gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype);
246 /* Returns a charge group index for the whole system */
247 t_block
248 gmx_mtop_global_cgs(const gmx_mtop_t *mtop);
251 /* Returns a single t_atoms struct for the whole system */
252 t_atoms
253 gmx_mtop_global_atoms(const gmx_mtop_t *mtop);
256 /* Generate a 'local' topology for the whole system.
257 * When feeEnergyInteractionsAtEnd == true, the free energy interactions will
258 * be sorted to the end.
260 gmx_localtop_t *
261 gmx_mtop_generate_local_top(const gmx_mtop_t *mtop, bool freeEnergyInteractionsAtEnd);
264 /* Converts a gmx_mtop_t struct to t_topology.
265 * All memory relating only to mtop will be freed.
267 t_topology
268 gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop);
270 /*! \brief Get vector of atoms indices from topology
272 * This function returns the indices of all particles with type
273 * eptAtom, that is shells, vsites etc. are left out.
274 * \param[in] mtop Molecular topology
275 * \returns Vector that will be filled with the atom indices
277 std::vector<size_t> get_atom_index(const gmx_mtop_t *mtop);
279 /*! \brief Converts a t_atoms struct to an mtop struct
281 * All pointers contained in \p atoms will be copied into \p mtop.
282 * Note that this will produce one moleculetype encompassing the whole system.
284 * \param[in] symtab The symbol table
285 * \param[in] name Pointer to the name for the topology
286 * \param[in] atoms The atoms to convert
287 * \param[out] mtop The molecular topology output containing atoms.
289 void
290 convertAtomsToMtop(t_symtab *symtab,
291 char **name,
292 t_atoms *atoms,
293 gmx_mtop_t *mtop);
295 #endif