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) 2013, 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.
43 /* Should be called after generating or reading mtop,
44 * to set some compute intesive variables to avoid
45 * N^2 operations later on.
48 gmx_mtop_finalize(gmx_mtop_t
*mtop
);
50 /* Counts the number of atoms of each type. State should be 0 for
51 * state A and 1 for state B types. typecount should have at
52 * least mtop->ffparams.atnr elements.
55 gmx_mtop_count_atomtypes(const gmx_mtop_t
*mtop
, int state
, int typecount
[]);
57 /* Returns the total number of charge groups in mtop */
59 ncg_mtop(const gmx_mtop_t
*mtop
);
61 /* Removes the charge groups, i.e. makes single atom charge groups, in mtop */
62 void gmx_mtop_remove_chargegroups(gmx_mtop_t
*mtop
);
65 /* Abstract data type for looking up atoms by global atom number */
66 typedef struct gmx_mtop_atomlookup
*gmx_mtop_atomlookup_t
;
68 /* Initialize atom lookup by global atom number */
70 gmx_mtop_atomlookup_init(const gmx_mtop_t
*mtop
);
72 /* As gmx_mtop_atomlookup_init, but optimized for atoms involved in settle */
74 gmx_mtop_atomlookup_settle_init(const gmx_mtop_t
*mtop
);
76 /* Destroy a gmx_mtop_atomlookup_t data structure */
78 gmx_mtop_atomlookup_destroy(gmx_mtop_atomlookup_t alook
);
81 /* Returns a pointer to the t_atom struct belonging to atnr_global.
82 * This can be an expensive operation, so if possible use
83 * one of the atom loop constructs below.
86 gmx_mtop_atomnr_to_atom(const gmx_mtop_atomlookup_t alook
,
91 /* Returns a pointer to the molecule interaction array ilist_mol[F_NRE]
92 * and the local atom number in the molecule belonging to atnr_global.
95 gmx_mtop_atomnr_to_ilist(const gmx_mtop_atomlookup_t alook
,
97 t_ilist
**ilist_mol
, int *atnr_offset
);
100 /* Returns the molecule block index
101 * and the molecule number in the block
102 * and the atom number offset for the atom indices in moltype
103 * belonging to atnr_global.
106 gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_atomlookup_t alook
,
108 int *molb
, int *molnr
, int *atnr_mol
);
111 /* Returns atom name, global resnr and residue name of atom atnr_global */
113 gmx_mtop_atominfo_global(const gmx_mtop_t
*mtop
, int atnr_global
,
114 char **atomname
, int *resnr
, char **resname
);
117 /* Abstract type for atom loop over all atoms */
118 typedef struct gmx_mtop_atomloop_all
*gmx_mtop_atomloop_all_t
;
120 /* Initialize an atom loop over all atoms in the system.
121 * The order of the atoms will be as in the state struct.
122 * Only use this when you really need to loop over all atoms,
123 * i.e. when you use groups which might differ per molecule,
124 * otherwise use gmx_mtop_atomloop_block.
126 gmx_mtop_atomloop_all_t
127 gmx_mtop_atomloop_all_init(const gmx_mtop_t
*mtop
);
129 /* Loop to the next atom.
130 * When not at the end:
131 * returns TRUE and at_global,
132 * writes the global atom number in *at_global
133 * and sets the pointer atom to the t_atom struct of that atom.
134 * When at the end, destroys aloop and returns FALSE.
136 * gmx_mtop_atomloop_all_t aloop;
137 * aloop = gmx_mtop_atomloop_all_init(mtop)
138 * while (gmx_mtop_atomloop_all_next(aloop,&at_global,&atom)) {
143 gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop
,
144 int *at_global
, t_atom
**atom
);
146 /* Return the atomname, the residue number and residue name
147 * of the current atom in the loop.
150 gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop
,
151 char **atomname
, int *resnr
, char **resname
);
153 /* Return the a pointer to the moltype struct of the current atom
154 * in the loop and the atom number in the molecule.
157 gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop
,
158 gmx_moltype_t
**moltype
, int *at_mol
);
161 /* Abstract type for atom loop over atoms in all molecule blocks */
162 typedef struct gmx_mtop_atomloop_block
*gmx_mtop_atomloop_block_t
;
164 /* Initialize an atom loop over atoms in all molecule blocks the system.
166 gmx_mtop_atomloop_block_t
167 gmx_mtop_atomloop_block_init(const gmx_mtop_t
*mtop
);
169 /* Loop to the next atom.
170 * When not at the end:
172 * sets the pointer atom to the t_atom struct of that atom
173 * and return the number of molecules corresponding to this atom.
174 * When at the end, destroys aloop and returns FALSE.
176 * gmx_mtop_atomloop_block_t aloop;
177 * aloop = gmx_mtop_atomloop_block_init(mtop)
178 * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
183 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop
,
184 t_atom
**atom
, int *nmol
);
187 /* Abstract type for ilist loop over all ilists */
188 typedef struct gmx_mtop_ilistloop
*gmx_mtop_ilistloop_t
;
190 /* Initialize an ilist loop over all molecule types in the system. */
192 gmx_mtop_ilistloop_init(const gmx_mtop_t
*mtop
);
195 /* Loop to the next molecule,
196 * When not at the end:
197 * returns TRUE and a pointer to the next array ilist_mol[F_NRE],
198 * writes the number of molecules for this ilist in *nmol.
199 * When at the end, destroys iloop and returns FALSE.
202 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop
,
203 t_ilist
**ilist_mol
, int *nmol
);
206 /* Abstract type for ilist loop over all ilists of all molecules */
207 typedef struct gmx_mtop_ilistloop_all
*gmx_mtop_ilistloop_all_t
;
209 /* Initialize an ilist loop over all molecule types in the system.
210 * Only use this when you really need to loop over all molecules,
211 * i.e. when you use groups which might differ per molecule,
212 * otherwise use gmx_mtop_ilistloop.
214 gmx_mtop_ilistloop_all_t
215 gmx_mtop_ilistloop_all_init(const gmx_mtop_t
*mtop
);
217 /* Loop to the next molecule,
218 * When not at the end:
219 * returns TRUE and a pointer to the next array ilist_mol[F_NRE],
220 * writes the atom offset which should be added to iatoms in atnr_offset.
221 * When at the end, destroys iloop and returns FALSE.
224 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop
,
225 t_ilist
**ilist_mol
, int *atnr_offset
);
228 /* Returns the total number of interactions in the system of type ftype */
230 gmx_mtop_ftype_count(const gmx_mtop_t
*mtop
, int ftype
);
233 /* Returns a charge group index for the whole system */
235 gmx_mtop_global_cgs(const gmx_mtop_t
*mtop
);
238 /* Returns a single t_atoms struct for the whole system */
240 gmx_mtop_global_atoms(const gmx_mtop_t
*mtop
);
243 /* Make all charge groups the size of one atom.
244 * When bKeepSingleMolCG==TRUE keep charge groups for molecules
245 * that consist of a single charge group.
248 gmx_mtop_make_atomic_charge_groups(gmx_mtop_t
*mtop
, gmx_bool bKeepSingleMolCG
);
251 /* Generate a 'local' topology for the whole system.
252 * When ir!=NULL the free energy interactions will be sorted to the end.
255 gmx_mtop_generate_local_top(const gmx_mtop_t
*mtop
, const t_inputrec
*ir
);
258 /* Converts a gmx_mtop_t struct to t_topology.
259 * All memory relating only to mtop will be freed.
262 gmx_mtop_t_to_t_topology(gmx_mtop_t
*mtop
);