3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gromacs Runs On Most of All Computer Systems
42 /* Returns the total number of charge groups in mtop */
44 ncg_mtop(const gmx_mtop_t
*mtop
);
47 /* Returns a pointer to the t_atom struct belonging to atnr_global.
48 * This can be an expensive operation, so if possible use
49 * one of the atom loop constructs below.
52 gmx_mtop_atomnr_to_atom(const gmx_mtop_t
*mtop
,int atnr_global
,
56 /* Returns a pointer to the molecule interaction array ilist_mol[F_NRE]
57 * and the local atom number in the molecule belonging to atnr_global.
60 gmx_mtop_atomnr_to_ilist(const gmx_mtop_t
*mtop
,int atnr_global
,
61 t_ilist
**ilist_mol
,int *atnr_offset
);
64 /* Returns the molecule block index
65 * and the molecule number in the block
66 * and the atom number offset for the atom indices in moltype
67 * belonging to atnr_global.
70 gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_t
*mtop
,int atnr_global
,
71 int *molb
,int *molnr
,int *atnr_mol
);
74 /* Returns atom name, global resnr and residue name of atom atnr_global */
76 gmx_mtop_atominfo_global(const gmx_mtop_t
*mtop
,int atnr_global
,
77 char **atomname
,int *resnr
,char **resname
);
80 /* Abstract type for atom loop over all atoms */
81 typedef struct gmx_mtop_atomloop_all
*gmx_mtop_atomloop_all_t
;
83 /* Initialize an atom loop over all atoms in the system.
84 * The order of the atoms will be as in the state struct.
85 * Only use this when you really need to loop over all atoms,
86 * i.e. when you use groups which might differ per molecule,
87 * otherwise use gmx_mtop_atomloop_block.
89 extern gmx_mtop_atomloop_all_t
90 gmx_mtop_atomloop_all_init(const gmx_mtop_t
*mtop
);
92 /* Loop to the next atom.
93 * When not at the end:
94 * returns TRUE and at_global,
95 * writes the global atom number in *at_global
96 * and sets the pointer atom to the t_atom struct of that atom.
97 * When at the end, destroys aloop and returns FALSE.
99 * gmx_mtop_atomloop_all_t aloop;
100 * aloop = gmx_mtop_atomloop_all_init(mtop)
101 * while (gmx_mtop_atomloop_all_next(aloop,&at_global,&atom)) {
106 gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop
,
107 int *at_global
,t_atom
**atom
);
109 /* Return the atomname, the residue number and residue name
110 * of the current atom in the loop.
113 gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop
,
114 char **atomname
,int *resnr
,char **resname
);
116 /* Return the a pointer to the moltype struct of the current atom
117 * in the loop and the atom number in the molecule.
120 gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop
,
121 gmx_moltype_t
**moltype
,int *at_mol
);
124 /* Abstract type for atom loop over atoms in all molecule blocks */
125 typedef struct gmx_mtop_atomloop_block
*gmx_mtop_atomloop_block_t
;
127 /* Initialize an atom loop over atoms in all molecule blocks the system.
129 extern gmx_mtop_atomloop_block_t
130 gmx_mtop_atomloop_block_init(const gmx_mtop_t
*mtop
);
132 /* Loop to the next atom.
133 * When not at the end:
135 * sets the pointer atom to the t_atom struct of that atom
136 * and return the number of molecules corresponding to this atom.
137 * When at the end, destroys aloop and returns FALSE.
139 * gmx_mtop_atomloop_block_t aloop;
140 * aloop = gmx_mtop_atomloop_block_init(mtop)
141 * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
146 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop
,
147 t_atom
**atom
,int *nmol
);
150 /* Abstract type for ilist loop over all ilists */
151 typedef struct gmx_mtop_ilistloop
*gmx_mtop_ilistloop_t
;
153 /* Initialize an ilist loop over all molecule types in the system. */
154 extern gmx_mtop_ilistloop_t
155 gmx_mtop_ilistloop_init(const gmx_mtop_t
*mtop
);
158 /* Loop to the next molecule,
159 * When not at the end:
160 * returns TRUE and a pointer to the next array ilist_mol[F_NRE],
161 * writes the number of molecules for this ilist in *nmol.
162 * When at the end, destroys iloop and returns FALSE.
165 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop
,
166 t_ilist
**ilist_mol
,int *nmol
);
169 /* Abstract type for ilist loop over all ilists of all molecules */
170 typedef struct gmx_mtop_ilistloop_all
*gmx_mtop_ilistloop_all_t
;
172 /* Initialize an ilist loop over all molecule types in the system.
173 * Only use this when you really need to loop over all molecules,
174 * i.e. when you use groups which might differ per molecule,
175 * otherwise use gmx_mtop_ilistloop.
177 extern gmx_mtop_ilistloop_all_t
178 gmx_mtop_ilistloop_all_init(const gmx_mtop_t
*mtop
);
180 /* Loop to the next molecule,
181 * When not at the end:
182 * returns TRUE and a pointer to the next array ilist_mol[F_NRE],
183 * writes the atom offset which should be added to iatoms in atnr_offset.
184 * When at the end, destroys iloop and returns FALSE.
187 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop
,
188 t_ilist
**ilist_mol
,int *atnr_offset
);
191 /* Returns the total number of interactions in the system of type ftype */
193 gmx_mtop_ftype_count(const gmx_mtop_t
*mtop
,int ftype
);
196 /* Returns a charge group index for the whole system */
198 gmx_mtop_global_cgs(const gmx_mtop_t
*mtop
);
201 /* Returns a single t_atoms struct for the whole system */
203 gmx_mtop_global_atoms(const gmx_mtop_t
*mtop
);
206 /* Make all charge groups the size of one atom.
207 * When bKeepSingleMolCG==TRUE keep charge groups for molecules
208 * that consist of a single charge group.
211 gmx_mtop_make_atomic_charge_groups(gmx_mtop_t
*mtop
,bool bKeepSingleMolCG
);
214 /* Generate a 'local' topology for the whole system.
215 * When ir!=NULL the free energy interactions will be sorted to the end.
217 extern gmx_localtop_t
*
218 gmx_mtop_generate_local_top(const gmx_mtop_t
*mtop
,const t_inputrec
*ir
);
221 /* Converts a gmx_mtop_t struct to t_topology.
222 * All memory relating only to mtop will be freed.
225 gmx_mtop_t_to_t_topology(gmx_mtop_t
*mtop
);