Add conserved quantity for Berendsen P-couple
[gromacs.git] / src / gromacs / topology / mtop_util.h
blob366cd13e9eaacb191ac982e41a97010eeb14d525
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/topology/topology.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/gmxassert.h"
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;
55 /* Should be called after generating or reading mtop,
56 * to set some compute intesive variables to avoid
57 * N^2 operations later on.
59 void
60 gmx_mtop_finalize(gmx_mtop_t *mtop);
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.
66 void
67 gmx_mtop_count_atomtypes(const gmx_mtop_t *mtop, int state, int typecount[]);
69 /* Returns the total number of charge groups in mtop */
70 int
71 ncg_mtop(const gmx_mtop_t *mtop);
73 /* Returns the total number of residues in mtop. */
74 int gmx_mtop_nres(const gmx_mtop_t *mtop);
76 /* Removes the charge groups, i.e. makes single atom charge groups, in mtop */
77 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop);
79 /* Abstract type for atom loop over all atoms */
80 typedef struct gmx_mtop_atomloop_all *gmx_mtop_atomloop_all_t;
82 /* Initialize an atom loop over all atoms in the system.
83 * The order of the atoms will be as in the state struct.
84 * Only use this when you really need to loop over all atoms,
85 * i.e. when you use groups which might differ per molecule,
86 * otherwise use gmx_mtop_atomloop_block.
88 gmx_mtop_atomloop_all_t
89 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop);
91 /* Loop to the next atom.
92 * When not at the end:
93 * returns TRUE and at_global,
94 * writes the global atom number in *at_global
95 * and sets the pointer atom to the t_atom struct of that atom.
96 * When at the end, destroys aloop and returns FALSE.
97 * Use as:
98 * gmx_mtop_atomloop_all_t aloop;
99 * aloop = gmx_mtop_atomloop_all_init(mtop)
100 * while (gmx_mtop_atomloop_all_next(aloop,&at_global,&atom)) {
101 * ...
104 gmx_bool
105 gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
106 int *at_global, const t_atom **atom);
108 /* Return the atomname, the residue number and residue name
109 * of the current atom in the loop.
111 void
112 gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
113 char **atomname, int *resnr, char **resname);
115 /* Return the a pointer to the moltype struct of the current atom
116 * in the loop and the atom number in the molecule.
118 void
119 gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
120 gmx_moltype_t **moltype, int *at_mol);
123 /* Abstract type for atom loop over atoms in all molecule blocks */
124 typedef struct gmx_mtop_atomloop_block *gmx_mtop_atomloop_block_t;
126 /* Initialize an atom loop over atoms in all molecule blocks the system.
128 gmx_mtop_atomloop_block_t
129 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop);
131 /* Loop to the next atom.
132 * When not at the end:
133 * returns TRUE
134 * sets the pointer atom to the t_atom struct of that atom
135 * and return the number of molecules corresponding to this atom.
136 * When at the end, destroys aloop and returns FALSE.
137 * Use as:
138 * gmx_mtop_atomloop_block_t aloop;
139 * aloop = gmx_mtop_atomloop_block_init(mtop)
140 * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
141 * ...
144 gmx_bool
145 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
146 const t_atom **atom, int *nmol);
149 /* Abstract type for ilist loop over all ilists */
150 typedef struct gmx_mtop_ilistloop *gmx_mtop_ilistloop_t;
152 /* Initialize an ilist loop over all molecule types in the system. */
153 gmx_mtop_ilistloop_t
154 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop);
157 /* Loop to the next molecule,
158 * When not at the end:
159 * returns TRUE and a pointer to the next array ilist_mol[F_NRE],
160 * writes the number of molecules for this ilist in *nmol.
161 * When at the end, destroys iloop and returns FALSE.
163 gmx_bool
164 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
165 t_ilist **ilist_mol, int *nmol);
168 /* Abstract type for ilist loop over all ilists of all molecules */
169 typedef struct gmx_mtop_ilistloop_all *gmx_mtop_ilistloop_all_t;
171 /* Initialize an ilist loop over all molecule types in the system.
172 * Only use this when you really need to loop over all molecules,
173 * i.e. when you use groups which might differ per molecule,
174 * otherwise use gmx_mtop_ilistloop.
176 gmx_mtop_ilistloop_all_t
177 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop);
179 /* Loop to the next molecule,
180 * When not at the end:
181 * returns TRUE and a pointer to the next array ilist_mol[F_NRE],
182 * writes the atom offset which should be added to iatoms in atnr_offset.
183 * When at the end, destroys iloop and returns FALSE.
185 gmx_bool
186 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
187 t_ilist **ilist_mol, int *atnr_offset);
190 /* Returns the total number of interactions in the system of type ftype */
192 gmx_mtop_ftype_count(const gmx_mtop_t *mtop, int ftype);
195 /* Returns a charge group index for the whole system */
196 t_block
197 gmx_mtop_global_cgs(const gmx_mtop_t *mtop);
200 /* Returns a single t_atoms struct for the whole system */
201 t_atoms
202 gmx_mtop_global_atoms(const gmx_mtop_t *mtop);
205 /* Generate a 'local' topology for the whole system.
206 * When freeEnergyInteractionsAtEnd == true, the free energy interactions will
207 * be sorted to the end.
209 gmx_localtop_t *
210 gmx_mtop_generate_local_top(const gmx_mtop_t *mtop, bool freeEnergyInteractionsAtEnd);
213 /* Converts a gmx_mtop_t struct to t_topology.
215 * If freeMTop == true, memory related to mtop will be freed so that done_top()
216 * on the result value will free all memory.
217 * If freeMTop == false, mtop and the return value will share some of their
218 * memory, and there is currently no way to consistently free all the memory.
220 t_topology
221 gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop, bool freeMTop);
223 /*! \brief Get vector of atoms indices from topology
225 * This function returns the indices of all particles with type
226 * eptAtom, that is shells, vsites etc. are left out.
227 * \param[in] mtop Molecular topology
228 * \returns Vector that will be filled with the atom indices
230 std::vector<size_t> get_atom_index(const gmx_mtop_t *mtop);
232 /*! \brief Converts a t_atoms struct to an mtop struct
234 * All pointers contained in \p atoms will be copied into \p mtop.
235 * Note that this will produce one moleculetype encompassing the whole system.
237 * \param[in] symtab The symbol table
238 * \param[in] name Pointer to the name for the topology
239 * \param[in] atoms The atoms to convert
240 * \param[out] mtop The molecular topology output containing atoms.
242 void
243 convertAtomsToMtop(t_symtab *symtab,
244 char **name,
245 t_atoms *atoms,
246 gmx_mtop_t *mtop);
248 #endif