d2e15d775b02a6178383c88c17b8025a34f3c14c
[alexxy/gromacs.git] / include / mtop_util.h
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
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.
14
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.
19  * 
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.
26  * 
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.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Gromacs Runs On Most of All Computer Systems
34  */
35 #include "typedefs.h"
36
37 #ifdef __cplusplus
38 extern "C" {
39 #endif
40
41 /* Should be called after generating or reading mtop,
42  * to set some compute intesive variables to avoid
43  * N^2 operations later on.
44  */
45 void
46 gmx_mtop_finalize(gmx_mtop_t *mtop);
47
48
49 /* Returns the total number of charge groups in mtop */
50 int
51 ncg_mtop(const gmx_mtop_t *mtop);
52
53
54 /* Returns a pointer to the t_atom struct belonging to atnr_global.
55  * This can be an expensive operation, so if possible use
56  * one of the atom loop constructs below.
57  */
58 void
59 gmx_mtop_atomnr_to_atom(const gmx_mtop_t *mtop,int atnr_global,
60                         t_atom **atom);
61
62
63 /* Returns a pointer to the molecule interaction array ilist_mol[F_NRE]
64  * and the local atom number in the molecule belonging to atnr_global.
65  */
66 void
67 gmx_mtop_atomnr_to_ilist(const gmx_mtop_t *mtop,int atnr_global,
68                          t_ilist **ilist_mol,int *atnr_offset);
69
70
71 /* Returns the molecule block index
72  * and the molecule number in the block
73  * and the atom number offset for the atom indices in moltype
74  * belonging to atnr_global.
75  */
76 void
77 gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_t *mtop,int atnr_global,
78                                 int *molb,int *molnr,int *atnr_mol);
79
80
81 /* Returns atom name, global resnr and residue name  of atom atnr_global */
82 void
83 gmx_mtop_atominfo_global(const gmx_mtop_t *mtop,int atnr_global,
84                          char **atomname,int *resnr,char **resname);
85
86
87 /* Abstract type for atom loop over all atoms */
88 typedef struct gmx_mtop_atomloop_all *gmx_mtop_atomloop_all_t;
89
90 /* Initialize an atom loop over all atoms in the system.
91  * The order of the atoms will be as in the state struct.
92  * Only use this when you really need to loop over all atoms,
93  * i.e. when you use groups which might differ per molecule,
94  * otherwise use gmx_mtop_atomloop_block.
95  */
96 gmx_mtop_atomloop_all_t
97 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop);
98
99 /* Loop to the next atom.
100  * When not at the end:
101  *   returns TRUE and at_global,
102  *   writes the global atom number in *at_global
103  *   and sets the pointer atom to the t_atom struct of that atom.
104  * When at the end, destroys aloop and returns FALSE.
105  * Use as:
106  * gmx_mtop_atomloop_all_t aloop;
107  * aloop = gmx_mtop_atomloop_all_init(mtop)
108  * while (gmx_mtop_atomloop_all_next(aloop,&at_global,&atom)) {
109  *     ...
110  * }
111  */
112 bool
113 gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
114                            int *at_global,t_atom **atom);
115
116 /* Return the atomname, the residue number and residue name
117  * of the current atom in the loop.
118  */
119 void
120 gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
121                             char **atomname,int *resnr,char **resname);
122
123 /* Return the a pointer to the moltype struct of the current atom
124  * in the loop and the atom number in the molecule.
125  */
126 void
127 gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
128                               gmx_moltype_t **moltype,int *at_mol);
129
130
131 /* Abstract type for atom loop over atoms in all molecule blocks */
132 typedef struct gmx_mtop_atomloop_block *gmx_mtop_atomloop_block_t;
133
134 /* Initialize an atom loop over atoms in all molecule blocks the system.
135  */
136 gmx_mtop_atomloop_block_t
137 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop);
138
139 /* Loop to the next atom.
140  * When not at the end:
141  *   returns TRUE
142  *   sets the pointer atom to the t_atom struct of that atom
143  *   and return the number of molecules corresponding to this atom.
144  * When at the end, destroys aloop and returns FALSE.
145  * Use as:
146  * gmx_mtop_atomloop_block_t aloop;
147  * aloop = gmx_mtop_atomloop_block_init(mtop)
148  * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
149  *     ...
150  * }
151  */
152 bool
153 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
154                              t_atom **atom,int *nmol);
155
156
157 /* Abstract type for ilist loop over all ilists */
158 typedef struct gmx_mtop_ilistloop *gmx_mtop_ilistloop_t;
159
160 /* Initialize an ilist loop over all molecule types in the system. */
161 gmx_mtop_ilistloop_t
162 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop);
163
164
165 /* Loop to the next molecule,
166  * When not at the end:
167  *   returns TRUE and a pointer to the next array ilist_mol[F_NRE],
168  *   writes the number of molecules for this ilist in *nmol.
169  * When at the end, destroys iloop and returns FALSE.
170  */
171 bool
172 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
173                         t_ilist **ilist_mol,int *nmol);
174
175
176 /* Abstract type for ilist loop over all ilists of all molecules */
177 typedef struct gmx_mtop_ilistloop_all *gmx_mtop_ilistloop_all_t;
178
179 /* Initialize an ilist loop over all molecule types in the system.
180  * Only use this when you really need to loop over all molecules,
181  * i.e. when you use groups which might differ per molecule,
182  * otherwise use gmx_mtop_ilistloop.
183  */
184 gmx_mtop_ilistloop_all_t
185 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop);
186
187 /* Loop to the next molecule,
188  * When not at the end:
189  *   returns TRUE and a pointer to the next array ilist_mol[F_NRE],
190  *   writes the atom offset which should be added to iatoms in atnr_offset.
191  * When at the end, destroys iloop and returns FALSE.
192  */
193 bool
194 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
195                             t_ilist **ilist_mol,int *atnr_offset);
196
197
198 /* Returns the total number of interactions in the system of type ftype */
199 int
200 gmx_mtop_ftype_count(const gmx_mtop_t *mtop,int ftype);
201
202
203 /* Returns a charge group index for the whole system */
204 t_block
205 gmx_mtop_global_cgs(const gmx_mtop_t *mtop);
206
207
208 /* Returns a single t_atoms struct for the whole system */ 
209 t_atoms
210 gmx_mtop_global_atoms(const gmx_mtop_t *mtop);
211
212
213 /* Make all charge groups the size of one atom.
214  * When bKeepSingleMolCG==TRUE keep charge groups for molecules
215  * that consist of a single charge group.
216  */
217 void
218 gmx_mtop_make_atomic_charge_groups(gmx_mtop_t *mtop,bool bKeepSingleMolCG);
219
220
221 /* Generate a 'local' topology for the whole system.
222  * When ir!=NULL the free energy interactions will be sorted to the end.
223  */
224 gmx_localtop_t *
225 gmx_mtop_generate_local_top(const gmx_mtop_t *mtop,const t_inputrec *ir);
226
227
228 /* Converts a gmx_mtop_t struct to t_topology.
229  * All memory relating only to mtop will be freed.
230  */
231 t_topology
232 gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop);
233
234 #ifdef __cplusplus
235 }
236 #endif
237