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,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.
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
44 #include "gromacs/topology/topology.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/gmxassert.h"
48 struct gmx_localtop_t
;
53 enum struct GmxQmmmMode
;
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.
59 /* Should be called after generating or reading mtop,
60 * to set some compute intesive variables to avoid
61 * N^2 operations later on.
64 gmx_mtop_finalize(gmx_mtop_t
*mtop
);
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.
71 gmx_mtop_count_atomtypes(const gmx_mtop_t
*mtop
, int state
, int typecount
[]);
73 /* Returns the total number of charge groups in mtop */
75 ncg_mtop(const gmx_mtop_t
*mtop
);
77 /*!\brief Returns the total number of molecules in mtop
79 * \param[in] mtop The global topology
81 int gmx_mtop_num_molecules(const gmx_mtop_t
&mtop
);
83 /* Returns the total number of residues in mtop. */
84 int gmx_mtop_nres(const gmx_mtop_t
*mtop
);
86 /* Removes the charge groups, i.e. makes single atom charge groups, in mtop */
87 void gmx_mtop_remove_chargegroups(gmx_mtop_t
*mtop
);
89 /* Abstract type for atom loop over all atoms */
90 typedef struct gmx_mtop_atomloop_all
*gmx_mtop_atomloop_all_t
;
92 /* Initialize an atom loop over all atoms in the system.
93 * The order of the atoms will be as in the state struct.
94 * Only use this when you really need to loop over all atoms,
95 * i.e. when you use groups which might differ per molecule,
96 * otherwise use gmx_mtop_atomloop_block.
98 gmx_mtop_atomloop_all_t
99 gmx_mtop_atomloop_all_init(const gmx_mtop_t
*mtop
);
101 /* Loop to the next atom.
102 * When not at the end:
103 * returns TRUE and at_global,
104 * writes the global atom number in *at_global
105 * and sets the pointer atom to the t_atom struct of that atom.
106 * When at the end, destroys aloop and returns FALSE.
108 * gmx_mtop_atomloop_all_t aloop;
109 * aloop = gmx_mtop_atomloop_all_init(mtop)
110 * while (gmx_mtop_atomloop_all_next(aloop,&at_global,&atom)) {
115 gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop
,
116 int *at_global
, const t_atom
**atom
);
118 /* Return the atomname, the residue number and residue name
119 * of the current atom in the loop.
122 gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop
,
123 char **atomname
, int *resnr
, char **resname
);
125 /* Return the a pointer to the moltype struct of the current atom
126 * in the loop and the atom number in the molecule.
129 gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop
,
130 const gmx_moltype_t
**moltype
, int *at_mol
);
133 /* Abstract type for atom loop over atoms in all molecule blocks */
134 typedef struct gmx_mtop_atomloop_block
*gmx_mtop_atomloop_block_t
;
136 /* Initialize an atom loop over atoms in all molecule blocks the system.
138 gmx_mtop_atomloop_block_t
139 gmx_mtop_atomloop_block_init(const gmx_mtop_t
*mtop
);
141 /* Loop to the next atom.
142 * When not at the end:
144 * sets the pointer atom to the t_atom struct of that atom
145 * and return the number of molecules corresponding to this atom.
146 * When at the end, destroys aloop and returns FALSE.
148 * gmx_mtop_atomloop_block_t aloop;
149 * aloop = gmx_mtop_atomloop_block_init(mtop)
150 * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
155 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop
,
156 const t_atom
**atom
, int *nmol
);
159 /* Abstract type for ilist loop over all ilists */
160 typedef struct gmx_mtop_ilistloop
*gmx_mtop_ilistloop_t
;
162 /* Initialize an ilist loop over all molecule types in the system. */
164 gmx_mtop_ilistloop_init(const gmx_mtop_t
*mtop
);
166 /* Initialize an ilist loop over all molecule types in the system. */
168 gmx_mtop_ilistloop_init(const gmx_mtop_t
&mtop
);
170 /* Loop to the next molecule,
171 * When not at the end:
172 * returns a valid pointer to the next array ilist_mol[F_NRE],
173 * writes the number of molecules for this ilist in *nmol.
174 * When at the end, destroys iloop and returns nullptr.
176 const InteractionLists
*
177 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop
,
180 /* Abstract type for ilist loop over all ilists of all molecules */
181 typedef struct gmx_mtop_ilistloop_all
*gmx_mtop_ilistloop_all_t
;
183 /* Initialize an ilist loop over all molecule types in the system.
184 * Only use this when you really need to loop over all molecules,
185 * i.e. when you use groups which might differ per molecule,
186 * otherwise use gmx_mtop_ilistloop.
188 gmx_mtop_ilistloop_all_t
189 gmx_mtop_ilistloop_all_init(const gmx_mtop_t
*mtop
);
191 /* Loop to the next molecule,
192 * When not at the end:
193 * returns a valid pointer to the next array ilist_mol[F_NRE],
194 * writes the atom offset which should be added to iatoms in atnr_offset.
195 * When at the end, destroys iloop and returns nullptr.
197 const InteractionLists
*
198 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop
,
202 /* Returns the total number of interactions in the system of type ftype */
204 gmx_mtop_ftype_count(const gmx_mtop_t
*mtop
, int ftype
);
206 /* Returns the total number of interactions in the system of type ftype */
208 gmx_mtop_ftype_count(const gmx_mtop_t
&mtop
, int ftype
);
210 /* Returns a charge group index for the whole system */
212 gmx_mtop_global_cgs(const gmx_mtop_t
*mtop
);
215 /* Returns a single t_atoms struct for the whole system */
217 gmx_mtop_global_atoms(const gmx_mtop_t
*mtop
);
220 /* Generate a 'local' topology for the whole system.
221 * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
222 * be sorted to the end.
225 gmx_mtop_generate_local_top(const gmx_mtop_t
*mtop
, bool freeEnergyInteractionsAtEnd
);
228 /*!\brief Creates and returns a struct with begin/end atom indices of all molecules
230 * \param[in] mtop The global topology
231 * \returns A RangePartitioning object with numBlocks() equal to the number
232 * of molecules and atom indices such that molecule m contains atoms a with:
233 * index[m] <= a < index[m+1].
235 gmx::RangePartitioning
gmx_mtop_molecules(const gmx_mtop_t
&mtop
);
238 /* Converts a gmx_mtop_t struct to t_topology.
240 * If the lifetime of the returned topology should be longer than that
241 * of mtop, your need to pass freeMtop==true.
242 * If freeMTop == true, memory related to mtop will be freed so that done_top()
243 * on the result value will free all memory.
244 * If freeMTop == false, mtop and the return value will share some of their
245 * memory, and there is currently no way to consistently free all the memory.
248 gmx_mtop_t_to_t_topology(gmx_mtop_t
*mtop
, bool freeMTop
);
250 /*! \brief Get vector of atoms indices from topology
252 * This function returns the indices of all particles with type
253 * eptAtom, that is shells, vsites etc. are left out.
254 * \param[in] mtop Molecular topology
255 * \returns Vector that will be filled with the atom indices
257 std::vector
<int> get_atom_index(const gmx_mtop_t
*mtop
);
259 /*! \brief Converts a t_atoms struct to an mtop struct
261 * All pointers contained in \p atoms will be copied into \p mtop.
262 * Note that this will produce one moleculetype encompassing the whole system.
264 * \param[in] symtab The symbol table
265 * \param[in] name Pointer to the name for the topology
266 * \param[in] atoms The atoms to convert
267 * \param[out] mtop The molecular topology output containing atoms.
270 convertAtomsToMtop(t_symtab
*symtab
,