ee2cf93de0c426b30ffda17abe9d82d95ef8eba7
[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 "visibility.h"
36 #include "typedefs.h"
37
38 #ifdef __cplusplus
39 extern "C" {
40 #endif
41
42 /* Should be called after generating or reading mtop,
43  * to set some compute intesive variables to avoid
44  * N^2 operations later on.
45  */
46 GMX_LIBGMX_EXPORT
47 void
48 gmx_mtop_finalize(gmx_mtop_t *mtop);
49
50
51 /* Returns the total number of charge groups in mtop */
52 GMX_LIBGMX_EXPORT
53 int
54 ncg_mtop(const gmx_mtop_t *mtop);
55
56 /* Removes the charge groups, i.e. makes single atom charge groups, in mtop */
57 GMX_LIBGMX_EXPORT
58 void gmx_mtop_remove_chargegroups(gmx_mtop_t *mtop);
59
60
61 /* Abstract data type for looking up atoms by global atom number */
62 typedef struct gmx_mtop_atomlookup *gmx_mtop_atomlookup_t;
63
64 /* Initialize atom lookup by global atom number */
65 GMX_LIBGMX_EXPORT
66 gmx_mtop_atomlookup_t
67 gmx_mtop_atomlookup_init(const gmx_mtop_t *mtop);
68
69 /* As gmx_mtop_atomlookup_init, but optimized for atoms involved in settle */
70 GMX_LIBGMX_EXPORT
71 gmx_mtop_atomlookup_t
72 gmx_mtop_atomlookup_settle_init(const gmx_mtop_t *mtop);
73
74 /* Destroy a gmx_mtop_atomlookup_t data structure */
75 GMX_LIBGMX_EXPORT
76 void
77 gmx_mtop_atomlookup_destroy(gmx_mtop_atomlookup_t alook);
78
79
80 /* Returns a pointer to the t_atom struct belonging to atnr_global.
81  * This can be an expensive operation, so if possible use
82  * one of the atom loop constructs below.
83  */
84 GMX_LIBGMX_EXPORT
85 void
86 gmx_mtop_atomnr_to_atom(const gmx_mtop_atomlookup_t alook,
87                         int atnr_global,
88                         t_atom **atom);
89
90
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.
93  */
94 GMX_LIBGMX_EXPORT
95 void
96 gmx_mtop_atomnr_to_ilist(const gmx_mtop_atomlookup_t alook,
97                          int atnr_global,
98                          t_ilist **ilist_mol,int *atnr_offset);
99
100
101 /* Returns the molecule block index
102  * and the molecule number in the block
103  * and the atom number offset for the atom indices in moltype
104  * belonging to atnr_global.
105  */
106 GMX_LIBGMX_EXPORT
107 void
108 gmx_mtop_atomnr_to_molblock_ind(const gmx_mtop_atomlookup_t alook,
109                                 int atnr_global,
110                                 int *molb,int *molnr,int *atnr_mol);
111
112
113 /* Returns atom name, global resnr and residue name  of atom atnr_global */
114 GMX_LIBGMX_EXPORT
115 void
116 gmx_mtop_atominfo_global(const gmx_mtop_t *mtop,int atnr_global,
117                          char **atomname,int *resnr,char **resname);
118
119
120 /* Abstract type for atom loop over all atoms */
121 typedef struct gmx_mtop_atomloop_all *gmx_mtop_atomloop_all_t;
122
123 /* Initialize an atom loop over all atoms in the system.
124  * The order of the atoms will be as in the state struct.
125  * Only use this when you really need to loop over all atoms,
126  * i.e. when you use groups which might differ per molecule,
127  * otherwise use gmx_mtop_atomloop_block.
128  */
129 GMX_LIBGMX_EXPORT
130 gmx_mtop_atomloop_all_t
131 gmx_mtop_atomloop_all_init(const gmx_mtop_t *mtop);
132
133 /* Loop to the next atom.
134  * When not at the end:
135  *   returns TRUE and at_global,
136  *   writes the global atom number in *at_global
137  *   and sets the pointer atom to the t_atom struct of that atom.
138  * When at the end, destroys aloop and returns FALSE.
139  * Use as:
140  * gmx_mtop_atomloop_all_t aloop;
141  * aloop = gmx_mtop_atomloop_all_init(mtop)
142  * while (gmx_mtop_atomloop_all_next(aloop,&at_global,&atom)) {
143  *     ...
144  * }
145  */
146 GMX_LIBGMX_EXPORT
147 gmx_bool
148 gmx_mtop_atomloop_all_next(gmx_mtop_atomloop_all_t aloop,
149                            int *at_global,t_atom **atom);
150
151 /* Return the atomname, the residue number and residue name
152  * of the current atom in the loop.
153  */
154 void
155 gmx_mtop_atomloop_all_names(gmx_mtop_atomloop_all_t aloop,
156                             char **atomname,int *resnr,char **resname);
157
158 /* Return the a pointer to the moltype struct of the current atom
159  * in the loop and the atom number in the molecule.
160  */
161 void
162 gmx_mtop_atomloop_all_moltype(gmx_mtop_atomloop_all_t aloop,
163                               gmx_moltype_t **moltype,int *at_mol);
164
165
166 /* Abstract type for atom loop over atoms in all molecule blocks */
167 typedef struct gmx_mtop_atomloop_block *gmx_mtop_atomloop_block_t;
168
169 /* Initialize an atom loop over atoms in all molecule blocks the system.
170  */
171 GMX_LIBGMX_EXPORT
172 gmx_mtop_atomloop_block_t
173 gmx_mtop_atomloop_block_init(const gmx_mtop_t *mtop);
174
175 /* Loop to the next atom.
176  * When not at the end:
177  *   returns TRUE
178  *   sets the pointer atom to the t_atom struct of that atom
179  *   and return the number of molecules corresponding to this atom.
180  * When at the end, destroys aloop and returns FALSE.
181  * Use as:
182  * gmx_mtop_atomloop_block_t aloop;
183  * aloop = gmx_mtop_atomloop_block_init(mtop)
184  * while (gmx_mtop_atomloop_block_next(aloop,&atom,&nmol)) {
185  *     ...
186  * }
187  */
188 GMX_LIBGMX_EXPORT
189 gmx_bool
190 gmx_mtop_atomloop_block_next(gmx_mtop_atomloop_block_t aloop,
191                              t_atom **atom,int *nmol);
192
193
194 /* Abstract type for ilist loop over all ilists */
195 typedef struct gmx_mtop_ilistloop *gmx_mtop_ilistloop_t;
196
197 /* Initialize an ilist loop over all molecule types in the system. */
198 GMX_LIBGMX_EXPORT
199 gmx_mtop_ilistloop_t
200 gmx_mtop_ilistloop_init(const gmx_mtop_t *mtop);
201
202
203 /* Loop to the next molecule,
204  * When not at the end:
205  *   returns TRUE and a pointer to the next array ilist_mol[F_NRE],
206  *   writes the number of molecules for this ilist in *nmol.
207  * When at the end, destroys iloop and returns FALSE.
208  */
209 GMX_LIBGMX_EXPORT
210 gmx_bool
211 gmx_mtop_ilistloop_next(gmx_mtop_ilistloop_t iloop,
212                         t_ilist **ilist_mol,int *nmol);
213
214
215 /* Abstract type for ilist loop over all ilists of all molecules */
216 typedef struct gmx_mtop_ilistloop_all *gmx_mtop_ilistloop_all_t;
217
218 /* Initialize an ilist loop over all molecule types in the system.
219  * Only use this when you really need to loop over all molecules,
220  * i.e. when you use groups which might differ per molecule,
221  * otherwise use gmx_mtop_ilistloop.
222  */
223 GMX_LIBGMX_EXPORT
224 gmx_mtop_ilistloop_all_t
225 gmx_mtop_ilistloop_all_init(const gmx_mtop_t *mtop);
226
227 /* Loop to the next molecule,
228  * When not at the end:
229  *   returns TRUE and a pointer to the next array ilist_mol[F_NRE],
230  *   writes the atom offset which should be added to iatoms in atnr_offset.
231  * When at the end, destroys iloop and returns FALSE.
232  */
233 GMX_LIBGMX_EXPORT
234 gmx_bool
235 gmx_mtop_ilistloop_all_next(gmx_mtop_ilistloop_all_t iloop,
236                             t_ilist **ilist_mol,int *atnr_offset);
237
238
239 /* Returns the total number of interactions in the system of type ftype */
240 GMX_LIBGMX_EXPORT
241 int
242 gmx_mtop_ftype_count(const gmx_mtop_t *mtop,int ftype);
243
244
245 /* Returns a charge group index for the whole system */
246 GMX_LIBGMX_EXPORT
247 t_block
248 gmx_mtop_global_cgs(const gmx_mtop_t *mtop);
249
250
251 /* Returns a single t_atoms struct for the whole system */ 
252 GMX_LIBGMX_EXPORT
253 t_atoms
254 gmx_mtop_global_atoms(const gmx_mtop_t *mtop);
255
256
257 /* Make all charge groups the size of one atom.
258  * When bKeepSingleMolCG==TRUE keep charge groups for molecules
259  * that consist of a single charge group.
260  */
261 void
262 gmx_mtop_make_atomic_charge_groups(gmx_mtop_t *mtop,gmx_bool bKeepSingleMolCG);
263
264
265 /* Generate a 'local' topology for the whole system.
266  * When ir!=NULL the free energy interactions will be sorted to the end.
267  */
268 GMX_LIBGMX_EXPORT
269 gmx_localtop_t *
270 gmx_mtop_generate_local_top(const gmx_mtop_t *mtop,const t_inputrec *ir);
271
272
273 /* Converts a gmx_mtop_t struct to t_topology.
274  * All memory relating only to mtop will be freed.
275  */
276 GMX_LIBGMX_EXPORT
277 t_topology
278 gmx_mtop_t_to_t_topology(gmx_mtop_t *mtop);
279
280 #ifdef __cplusplus
281 }
282 #endif
283