Enforced rotation: made fit type controllable by mdp option, reduced number of output...
[gromacs/adressmacs.git] / include / indexutil.h
blob7aec7d14fa08921a2e806ae563472b228258c99f
1 /*
3 * This source code is part of
5 * G R O M A C S
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
29 * For more info, check our website at http://www.gromacs.org
31 /*! \file
32 * \brief API for handling index files and index groups.
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.
47 * Mostly, these functions are used internally by the library and the
48 * selection engine.
49 * However, some of the checking functions can be useful in user code to
50 * check the validity of input groups.
51 * Also, the mapping functions are useful when dealing with dynamic index
52 * groups.
54 #ifndef INDEXUTIL_H
55 #define INDEXUTIL_H
57 #include <typedefs.h>
59 #ifdef __cplusplus
60 extern "C"
62 #endif
64 /*! Stores a set of index groups.*/
65 typedef struct gmx_ana_indexgrps_t gmx_ana_indexgrps_t;
67 /*! \brief
68 * Specifies the type of index partition or index mapping in several contexts.
70 * \see gmx_ana_index_make_block(), gmx_ana_indexmap_init()
72 typedef enum
74 INDEX_UNKNOWN, /*!< Unknown index type.*/
75 INDEX_ATOM, /*!< Each atom in a separate block.*/
76 INDEX_RES, /*!< Each residue in a separate block.*/
77 INDEX_MOL, /*!< Each molecule in a separate block.*/
78 INDEX_ALL, /*!< All atoms in a single block.*/
79 } e_index_t;
81 /*! \brief
82 * Stores a single index group.
84 typedef struct gmx_ana_index_t
86 /*! Number of atoms.*/
87 int isize;
88 /*! List of atoms.*/
89 atom_id *index;
90 /*! Group name.*/
91 char *name;
92 /*! Number of items allocated for \p index.*/
93 int nalloc_index;
94 } gmx_ana_index_t;
96 /*! \brief
97 * Data structure for calculating index group mappings.
99 typedef struct gmx_ana_indexmap_t
101 /*! Type of the mapping.*/
102 e_index_t type;
103 /*! \brief
104 * Current number of mapped values.
106 * This is the current number of values in the \p refid and \p mapid
107 * arrays.
108 * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), this
109 * is always equal to \p b.nr, i.e., the number of blocks in the
110 * original index group.
112 int nr;
113 /*! \brief
114 * Current reference IDs.
116 * This array provides a mapping from the current index group (last given
117 * to gmx_ana_indexmap_update()) to the blocks in \p b, i.e., the
118 * original index group used in gmx_ana_indexmap_init().
119 * The mapping is zero-based.
120 * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), the indices
121 * for blocks not present in the current group are set to -1, otherwise
122 * they are removed completely and the \p nr field updated.
124 int *refid;
125 /*! \brief
126 * Current mapped IDs.
128 * This array provides an arbitrary mapping from the current index group
129 * to the original index group. Instead of a zero-based mapping, the
130 * values from the \p orgid array are used. That is,
131 * \c mapid[i]=orgid[refid[i]].
132 * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), this array
133 * equals \p orgid.
135 int *mapid;
136 /*! \brief
137 * Mapped block structure.
139 * A block structure that corresponds to the current index group.
141 t_block mapb;
143 /*! \brief
144 * Arbitrary ID numbers for the blocks.
146 * This array has \p b.nr elements, each defining an ID number for a
147 * block in \p b.
148 * These are initialized in gmx_ana_indexmap_init() based on the type:
149 * - \ref INDEX_ATOM : the atom indices
150 * - \ref INDEX_RES : the residue numbers
151 * - \ref INDEX_MOL : the molecule numbers
153 * All the above numbers are zero-based.
154 * After gmx_ana_indexmap_init(), the user is free to change these values
155 * if the above are not appropriate.
156 * The mapped values can be read through \p mapid.
158 int *orgid;
160 /*! \brief
161 * Block data that defines the mapping (internal use only).
163 * The data is initialized by gmx_ana_indexmap_init() and is not changed
164 * after that.
165 * Hence, it cannot be directly applied to the index group passed to
166 * gmx_ana_indexmap_update() unless \p bMaskOnly was specified or the
167 * index group is identical to the one provided to gmx_ana_indexmap_init().
169 t_blocka b;
170 /*! \brief
171 * TRUE if the current reference IDs are for the whole group (internal use only).
173 * This is used internally to optimize the evaluation such that
174 * gmx_ana_indexmap_update() does not take any time if the group is
175 * actually static.
177 bool bStatic;
178 /*! \brief
179 * TRUE if the current mapping is for the whole group (internal use only).
181 * This is used internally to optimize the evaluation such that
182 * gmx_ana_indexmap_update() does not take any time if the group is
183 * actually static.
185 bool bMapStatic;
186 } gmx_ana_indexmap_t;
189 /*! \name Functions for handling gmx_ana_indexgrps_t
191 /*@{*/
192 /*! Allocate memory for index groups.*/
193 extern void
194 gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps);
195 /*! Initializes index groups from arrays.*/
196 extern void
197 gmx_ana_indexgrps_set(gmx_ana_indexgrps_t **g, int ngrps, int *isize,
198 atom_id **index, char **name, bool bFree);
199 /*! Reads index groups from a file or constructs them from topology.*/
200 extern void
201 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
202 const char *fnm);
203 /*! Ask user to select index groups, possibly constructing groups from
204 topology.*/
205 extern void
206 gmx_ana_indexgrps_get(gmx_ana_indexgrps_t **g, t_topology *top,
207 const char *fnm, int ngrps);
208 /*! Ask user to select index groups from those specified in a file.*/
209 extern void
210 gmx_ana_indexgrps_rd(gmx_ana_indexgrps_t **g, const char *fnm, int ngrps);
211 /*! Frees memory allocated for index groups.*/
212 extern void
213 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g);
214 /*! Create a deep copy of \c gmx_ana_indexgrps_t.*/
215 extern void
216 gmx_ana_indexgrps_clone(gmx_ana_indexgrps_t **dest, gmx_ana_indexgrps_t *src);
217 /*! Returns TRUE if the index group structure is emtpy.*/
218 extern bool
219 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g);
221 /*! Returns a pointer to an index group.*/
222 extern gmx_ana_index_t *
223 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n);
224 /*! Extracts a single index group.*/
225 extern bool
226 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, int n);
227 /*! Finds and extracts a single index group by name.*/
228 extern bool
229 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, char *name);
231 /*! Writes out the contents of index groups.*/
232 extern void
233 gmx_ana_indexgrps_dump(gmx_ana_indexgrps_t *g);
234 /*@}*/
236 /*! \name Functions for handling gmx_ana_index_t
238 /*@{*/
239 /*! Reserves memory to store an index group of size \p isize.*/
240 extern void
241 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize);
242 /*! Frees any memory not necessary to hold the current contents.*/
243 extern void
244 gmx_ana_index_squeeze(gmx_ana_index_t *g);
245 /*! Initializes an empty index group.*/
246 extern void
247 gmx_ana_index_clear(gmx_ana_index_t *g);
248 /*! Constructs a \c gmx_ana_index_t from given values.*/
249 extern void
250 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, char *name,
251 int nalloc);
252 /*! Creates a simple index group from the first to the \p natoms'th atom.*/
253 extern void
254 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms, char *name);
255 /*! Frees memory allocated for an index group.*/
256 extern void
257 gmx_ana_index_deinit(gmx_ana_index_t *g);
258 /*! Copies a \c gmx_ana_index_t.*/
259 extern void
260 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc);
262 /*! Writes out the contents of a index group.*/
263 extern void
264 gmx_ana_index_dump(gmx_ana_index_t *g, int i);
266 /*! Checks whether all indices are between 0 and \p natoms.*/
267 extern void
268 gmx_ana_index_check(gmx_ana_index_t *g, int natoms);
269 /*! Checks whether an index group is sorted.*/
270 extern bool
271 gmx_ana_index_check_sorted(gmx_ana_index_t *g);
272 /*@}*/
274 /*! \name Functions for set operations on gmx_ana_index_t
276 /*@{*/
277 /*! Sorts the indices within an index group.*/
278 extern void
279 gmx_ana_index_sort(gmx_ana_index_t *g);
280 /*! Checks whether two index groups are equal.*/
281 extern bool
282 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b);
283 /*! Checks whether a sorted index group contains another sorted index group.*/
284 extern bool
285 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b);
287 /*! Calculates the intersection between two sorted index groups.*/
288 extern void
289 gmx_ana_index_intersection(gmx_ana_index_t *dest,
290 gmx_ana_index_t *a, gmx_ana_index_t *b);
291 /*! Calculates the set difference between two sorted index groups.*/
292 extern void
293 gmx_ana_index_difference(gmx_ana_index_t *dest,
294 gmx_ana_index_t *a, gmx_ana_index_t *b);
295 /*! Calculates the size of the difference between two sorted index groups.*/
296 extern int
297 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b);
298 /*! Calculates the union of two sorted index groups.*/
299 extern void
300 gmx_ana_index_union(gmx_ana_index_t *dest,
301 gmx_ana_index_t *a, gmx_ana_index_t *b);
302 /*! Merges two distinct sorted index groups.*/
303 extern void
304 gmx_ana_index_merge(gmx_ana_index_t *dest,
305 gmx_ana_index_t *a, gmx_ana_index_t *b);
306 /*! Calculates the intersection and the difference in one call.*/
307 extern void
308 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
309 gmx_ana_index_t *src, gmx_ana_index_t *g);
310 /*@}*/
312 /*! \name Functions for handling gmx_ana_indexmap_t and related things
314 /*@{*/
315 /*! Partition a group based on topology information.*/
316 extern void
317 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
318 e_index_t type, bool bComplete);
319 /*! Checks whether a group consists of full blocks.*/
320 extern bool
321 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b);
322 /*! Checks whether a group consists of full blocks.*/
323 extern bool
324 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b);
325 /*! Checks whether a group consists of full residues/molecules.*/
326 extern bool
327 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type, t_topology *top);
329 /*! Initializes an empty index group mapping.*/
330 extern void
331 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m);
332 /*! Reserves memory for an index group mapping.*/
333 extern void
334 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize);
335 /*! Initialize index group mapping.*/
336 extern void
337 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
338 t_topology *top, e_index_t type);
339 /*! Frees memory allocated for index group mapping.*/
340 extern void
341 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m);
342 /*! Makes a deep copy of an index group mapping.*/
343 extern void
344 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst);
345 /*! Updates an index group mapping.*/
346 extern void
347 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g, bool bMaskOnly);
348 /*@}*/
350 #ifdef __cplusplus
352 #endif
354 #endif