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