Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / selection / indexutil.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  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11  * Copyright (c) 2001-2009, The GROMACS development team,
12  * check out http://www.gromacs.org for more information.
13
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  *
19  * If you want to redistribute modifications, please consider that
20  * scientific software is very special. Version control is crucial -
21  * bugs must be traceable. We will be happy to consider code for
22  * inclusion in the official distribution, but derived work must not
23  * be called official GROMACS. Details are found in the README & COPYING
24  * files - if they are missing, get the official version at www.gromacs.org.
25  *
26  * To help us fund GROMACS development, we humbly ask that you cite
27  * the papers on the package - you can find them in the top README file.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \file
32  * \brief API for handling index files and index groups.
33  *
34  * The API contains functions and data structures for handling index
35  * files more conveniently than as several separate variables.
36  * In addition to basic functions for initializing the data structures and
37  * making copies, functions are provided for performing (most) set operations
38  * on sorted index groups.
39  * There is also a function for partitioning a index group based on
40  * topology information such as residues or molecules.
41  * Finally, there is a set of functions for constructing mappings between
42  * an index group and its subgroups such.
43  * These can be used with dynamic index group in calculations if one
44  * needs to have a unique ID for each possible atom/residue/molecule in the
45  * selection, e.g., for analysis of dynamics or for look-up tables.
46  *
47  * Mostly, these functions are used internally by the selection engine, but
48  * it is necessary to use some of these functions in order to provide external
49  * index groups to a gmx::SelectionCollection.
50  * Some of the checking functions can be useful outside the selection engine to
51  * check the validity of input groups.
52  *
53  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
54  * \inlibraryapi
55  * \ingroup module_selection
56  */
57 #ifndef GMX_SELECTION_INDEXUTIL_H
58 #define GMX_SELECTION_INDEXUTIL_H
59
60 #include "../legacyheaders/typedefs.h"
61
62 /** Stores a set of index groups. */
63 typedef struct gmx_ana_indexgrps_t gmx_ana_indexgrps_t;
64
65 /*! \brief
66  * Specifies the type of index partition or index mapping in several contexts.
67  *
68  * \see gmx_ana_index_make_block(), gmx_ana_indexmap_init()
69  */
70 typedef enum
71 {
72     INDEX_UNKNOWN, /**< Unknown index type.*/
73     INDEX_ATOM,    /**< Each atom in a separate block.*/
74     INDEX_RES,     /**< Each residue in a separate block.*/
75     INDEX_MOL,     /**< Each molecule in a separate block.*/
76     INDEX_ALL      /**< All atoms in a single block.*/
77 } e_index_t;
78
79 /*! \brief
80  * Stores a single index group.
81  */
82 typedef struct gmx_ana_index_t
83 {
84     /** Number of atoms. */
85     int                 isize;
86     /** List of atoms. */
87     atom_id            *index;
88     /** Group name. */
89     char               *name;
90     /** Number of items allocated for \p index. */
91     int                 nalloc_index;
92 } gmx_ana_index_t;
93
94 /*! \brief
95  * Data structure for calculating index group mappings.
96  */
97 typedef struct gmx_ana_indexmap_t
98 {
99     /** Type of the mapping. */
100     e_index_t           type;
101     /*! \brief
102      * Current number of mapped values.
103      *
104      * This is the current number of values in the \p refid and \p mapid
105      * arrays.
106      * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), this
107      * is always equal to \p b.nr, i.e., the number of blocks in the
108      * original index group.
109      */
110     int                 nr;
111     /*! \brief
112      * Current reference IDs.
113      *
114      * This array provides a mapping from the current index group (last given
115      * to gmx_ana_indexmap_update()) to the blocks in \p b, i.e., the
116      * original index group used in gmx_ana_indexmap_init().
117      * The mapping is zero-based.
118      * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), the indices
119      * for blocks not present in the current group are set to -1, otherwise
120      * they are removed completely and the \p nr field updated.
121      */
122     int                *refid;
123     /*! \brief
124      * Current mapped IDs.
125      *
126      * This array provides an arbitrary mapping from the current index group
127      * to the original index group. Instead of a zero-based mapping, the
128      * values from the \p orgid array are used. That is,
129      * \c mapid[i]=orgid[refid[i]].
130      * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), this array
131      * equals \p orgid.
132      */
133     int                *mapid;
134     /*! \brief
135      * Mapped block structure.
136      *
137      * A block structure that corresponds to the current index group.
138      */
139     t_block             mapb;
140
141     /*! \brief
142      * Arbitrary ID numbers for the blocks.
143      *
144      * This array has \p b.nr elements, each defining an ID number for a
145      * block in \p b.
146      * These are initialized in gmx_ana_indexmap_init() based on the type:
147      *  - \ref INDEX_ATOM : the atom indices
148      *  - \ref INDEX_RES :  the residue numbers
149      *  - \ref INDEX_MOL :  the molecule numbers
150      *
151      * All the above numbers are zero-based.
152      * After gmx_ana_indexmap_init(), the user is free to change these values
153      * if the above are not appropriate.
154      * The mapped values can be read through \p mapid.
155      */
156     int                *orgid;
157
158     /*! \brief
159      * Block data that defines the mapping (internal use only).
160      *
161      * The data is initialized by gmx_ana_indexmap_init() and is not changed
162      * after that.
163      * Hence, it cannot be directly applied to the index group passed to
164      * gmx_ana_indexmap_update() unless \p bMaskOnly was specified or the
165      * index group is identical to the one provided to gmx_ana_indexmap_init().
166      */
167     t_blocka            b;
168     /*! \brief
169      * true if the current reference IDs are for the whole group (internal use only).
170      *
171      * This is used internally to optimize the evaluation such that
172      * gmx_ana_indexmap_update() does not take any time if the group is
173      * actually static.
174      */
175     bool                bStatic;
176     /*! \brief
177      * true if the current mapping is for the whole group (internal use only).
178      *
179      * This is used internally to optimize the evaluation such that
180      * gmx_ana_indexmap_update() does not take any time if the group is
181      * actually static.
182      */
183     bool                bMapStatic;
184 } gmx_ana_indexmap_t;
185
186
187 /*! \name Functions for handling gmx_ana_indexgrps_t
188  */
189 /*@{*/
190 /** Allocate memory for index groups. */
191 void
192 gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps);
193 /** Reads index groups from a file or constructs them from topology. */
194 void
195 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
196                        const char *fnm);
197 /** Frees memory allocated for index groups. */
198 void
199 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g);
200 /** Create a deep copy of \c gmx_ana_indexgrps_t. */
201 void
202 gmx_ana_indexgrps_clone(gmx_ana_indexgrps_t **dest, gmx_ana_indexgrps_t *src);
203 /** Returns true if the index group structure is emtpy. */
204 bool
205 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g);
206
207 /** Returns a pointer to an index group. */
208 gmx_ana_index_t *
209 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n);
210 /** Extracts a single index group. */
211 bool
212 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, int n);
213 /** Finds and extracts a single index group by name. */
214 bool
215 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, char *name);
216
217 /** Writes out a list of index groups. */
218 void
219 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn);
220 /*@}*/
221
222 /*! \name Functions for handling gmx_ana_index_t
223  */
224 /*@{*/
225 /** Reserves memory to store an index group of size \p isize. */
226 void
227 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize);
228 /** Frees any memory not necessary to hold the current contents. */
229 void
230 gmx_ana_index_squeeze(gmx_ana_index_t *g);
231 /** Initializes an empty index group. */
232 void
233 gmx_ana_index_clear(gmx_ana_index_t *g);
234 /** Constructs a \c gmx_ana_index_t from given values. */
235 void
236 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, char *name,
237                   int nalloc);
238 /** Creates a simple index group from the first to the \p natoms'th atom. */
239 void
240 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms, char *name);
241 /** Frees memory allocated for an index group. */
242 void
243 gmx_ana_index_deinit(gmx_ana_index_t *g);
244 /** Copies a \c gmx_ana_index_t. */
245 void
246 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc);
247
248 /** Writes out the contents of a index group. */
249 void
250 gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int i, int maxn);
251
252 /** Checks whether all indices are between 0 and \p natoms. */
253 void
254 gmx_ana_index_check(gmx_ana_index_t *g, int natoms);
255 /** Checks whether an index group is sorted. */
256 bool
257 gmx_ana_index_check_sorted(gmx_ana_index_t *g);
258 /*@}*/
259
260 /*! \name Functions for set operations on gmx_ana_index_t
261  */
262 /*@{*/
263 /** Sorts the indices within an index group. */
264 void
265 gmx_ana_index_sort(gmx_ana_index_t *g);
266 /** Checks whether two index groups are equal. */
267 bool
268 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b);
269 /** Checks whether a sorted index group contains another sorted index group. */
270 bool
271 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b);
272
273 /** Calculates the intersection between two sorted index groups. */
274 void
275 gmx_ana_index_intersection(gmx_ana_index_t *dest,
276                            gmx_ana_index_t *a, gmx_ana_index_t *b);
277 /** Calculates the set difference between two sorted index groups. */
278 void
279 gmx_ana_index_difference(gmx_ana_index_t *dest,
280                          gmx_ana_index_t *a, gmx_ana_index_t *b);
281 /** Calculates the size of the difference between two sorted index groups. */
282 int
283 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b);
284 /** Calculates the union of two sorted index groups. */
285 void
286 gmx_ana_index_union(gmx_ana_index_t *dest,
287                     gmx_ana_index_t *a, gmx_ana_index_t *b);
288 /** Merges two distinct sorted index groups. */
289 void
290 gmx_ana_index_merge(gmx_ana_index_t *dest,
291                     gmx_ana_index_t *a, gmx_ana_index_t *b);
292 /** Calculates the intersection and the difference in one call. */
293 void
294 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
295                         gmx_ana_index_t *src, gmx_ana_index_t *g);
296 /*@}*/
297
298 /*! \name Functions for handling gmx_ana_indexmap_t and related things
299  */
300 /*@{*/
301 /** Partition a group based on topology information. */
302 void
303 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
304                          e_index_t type, bool bComplete);
305 /** Checks whether a group consists of full blocks. */
306 bool
307 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b);
308 /** Checks whether a group consists of full blocks. */
309 bool
310 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b);
311 /** Checks whether a group consists of full residues/molecules. */
312 bool
313 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type, t_topology *top);
314
315 /** Initializes an empty index group mapping. */
316 void
317 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m);
318 /** Reserves memory for an index group mapping. */
319 void
320 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize);
321 /** Initializes an index group mapping. */
322 void
323 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
324                       t_topology *top, e_index_t type);
325 /** Sets an index group mapping to be static. */
326 void
327 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b);
328 /** Frees memory allocated for index group mapping. */
329 void
330 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m);
331 /** Makes a deep copy of an index group mapping. */
332 void
333 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst);
334 /** Updates an index group mapping. */
335 void
336 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g, bool bMaskOnly);
337 /*@}*/
338
339 #endif