2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2009, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief API for handling index files and index groups.
41 * The API contains functions and data structures for handling index
42 * files more conveniently than as several separate variables.
43 * In addition to basic functions for initializing the data structures and
44 * making copies, functions are provided for performing (most) set operations
45 * on sorted index groups.
46 * There is also a function for partitioning a index group based on
47 * topology information such as residues or molecules.
48 * Finally, there is a set of functions for constructing mappings between
49 * an index group and its subgroups such.
50 * These can be used with dynamic index group in calculations if one
51 * needs to have a unique ID for each possible atom/residue/molecule in the
52 * selection, e.g., for analysis of dynamics or for look-up tables.
54 * Mostly, these functions are used internally by the library and the
56 * However, some of the checking functions can be useful in user code to
57 * check the validity of input groups.
58 * Also, the mapping functions are useful when dealing with dynamic index
63 #include "visibility.h"
71 /** Stores a set of index groups. */
72 typedef struct gmx_ana_indexgrps_t gmx_ana_indexgrps_t
;
75 * Specifies the type of index partition or index mapping in several contexts.
77 * \see gmx_ana_index_make_block(), gmx_ana_indexmap_init()
81 INDEX_UNKNOWN
, /**< Unknown index type.*/
82 INDEX_ATOM
, /**< Each atom in a separate block.*/
83 INDEX_RES
, /**< Each residue in a separate block.*/
84 INDEX_MOL
, /**< Each molecule in a separate block.*/
85 INDEX_ALL
/**< All atoms in a single block.*/
89 * Stores a single index group.
91 typedef struct gmx_ana_index_t
93 /** Number of atoms. */
99 /** Number of items allocated for \p index. */
104 * Data structure for calculating index group mappings.
106 typedef struct gmx_ana_indexmap_t
108 /** Type of the mapping. */
111 * Current number of mapped values.
113 * This is the current number of values in the \p refid and \p mapid
115 * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), this
116 * is always equal to \p b.nr, i.e., the number of blocks in the
117 * original index group.
121 * Current reference IDs.
123 * This array provides a mapping from the current index group (last given
124 * to gmx_ana_indexmap_update()) to the blocks in \p b, i.e., the
125 * original index group used in gmx_ana_indexmap_init().
126 * The mapping is zero-based.
127 * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), the indices
128 * for blocks not present in the current group are set to -1, otherwise
129 * they are removed completely and the \p nr field updated.
133 * Current mapped IDs.
135 * This array provides an arbitrary mapping from the current index group
136 * to the original index group. Instead of a zero-based mapping, the
137 * values from the \p orgid array are used. That is,
138 * \c mapid[i]=orgid[refid[i]].
139 * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), this array
144 * Mapped block structure.
146 * A block structure that corresponds to the current index group.
151 * Arbitrary ID numbers for the blocks.
153 * This array has \p b.nr elements, each defining an ID number for a
155 * These are initialized in gmx_ana_indexmap_init() based on the type:
156 * - \ref INDEX_ATOM : the atom indices
157 * - \ref INDEX_RES : the residue numbers
158 * - \ref INDEX_MOL : the molecule numbers
160 * All the above numbers are zero-based.
161 * After gmx_ana_indexmap_init(), the user is free to change these values
162 * if the above are not appropriate.
163 * The mapped values can be read through \p mapid.
168 * Block data that defines the mapping (internal use only).
170 * The data is initialized by gmx_ana_indexmap_init() and is not changed
172 * Hence, it cannot be directly applied to the index group passed to
173 * gmx_ana_indexmap_update() unless \p bMaskOnly was specified or the
174 * index group is identical to the one provided to gmx_ana_indexmap_init().
178 * TRUE if the current reference IDs are for the whole group (internal use only).
180 * This is used internally to optimize the evaluation such that
181 * gmx_ana_indexmap_update() does not take any time if the group is
186 * TRUE if the current mapping is for the whole group (internal use only).
188 * This is used internally to optimize the evaluation such that
189 * gmx_ana_indexmap_update() does not take any time if the group is
193 } gmx_ana_indexmap_t
;
196 /*! \name Functions for handling gmx_ana_indexgrps_t
199 /** Allocate memory for index groups. */
201 gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t
**g
, int ngrps
);
202 /** Initializes index groups from arrays. */
204 gmx_ana_indexgrps_set(gmx_ana_indexgrps_t
**g
, int ngrps
, int *isize
,
205 atom_id
**index
, char **name
, gmx_bool bFree
);
206 /** Reads index groups from a file or constructs them from topology. */
208 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t
**g
, t_topology
*top
,
210 /** Ask user to select index groups, possibly constructing groups from
213 gmx_ana_indexgrps_get(gmx_ana_indexgrps_t
**g
, t_topology
*top
,
214 const char *fnm
, int ngrps
);
215 /** Ask user to select index groups from those specified in a file. */
217 gmx_ana_indexgrps_rd(gmx_ana_indexgrps_t
**g
, const char *fnm
, int ngrps
);
218 /** Frees memory allocated for index groups. */
220 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t
*g
);
221 /** Create a deep copy of \c gmx_ana_indexgrps_t. */
223 gmx_ana_indexgrps_clone(gmx_ana_indexgrps_t
**dest
, gmx_ana_indexgrps_t
*src
);
224 /** Returns TRUE if the index group structure is emtpy. */
226 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t
*g
);
228 /** Returns a pointer to an index group. */
230 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t
*g
, int n
);
231 /** Extracts a single index group. */
233 gmx_ana_indexgrps_extract(gmx_ana_index_t
*dest
, gmx_ana_indexgrps_t
*src
, int n
);
234 /** Finds and extracts a single index group by name. */
236 gmx_ana_indexgrps_find(gmx_ana_index_t
*dest
, gmx_ana_indexgrps_t
*src
, char *name
);
238 /** Writes out a list of index groups. */
240 gmx_ana_indexgrps_print(gmx_ana_indexgrps_t
*g
, int maxn
);
243 /*! \name Functions for handling gmx_ana_index_t
246 /** Reserves memory to store an index group of size \p isize. */
248 gmx_ana_index_reserve(gmx_ana_index_t
*g
, int isize
);
249 /** Frees any memory not necessary to hold the current contents. */
251 gmx_ana_index_squeeze(gmx_ana_index_t
*g
);
252 /** Initializes an empty index group. */
254 gmx_ana_index_clear(gmx_ana_index_t
*g
);
255 /** Constructs a \c gmx_ana_index_t from given values. */
257 gmx_ana_index_set(gmx_ana_index_t
*g
, int isize
, atom_id
*index
, char *name
,
259 /** Creates a simple index group from the first to the \p natoms'th atom. */
261 gmx_ana_index_init_simple(gmx_ana_index_t
*g
, int natoms
, char *name
);
262 /** Frees memory allocated for an index group. */
264 gmx_ana_index_deinit(gmx_ana_index_t
*g
);
265 /** Copies a \c gmx_ana_index_t. */
267 gmx_ana_index_copy(gmx_ana_index_t
*dest
, gmx_ana_index_t
*src
, gmx_bool bAlloc
);
269 /** Writes out the contents of a index group. */
271 gmx_ana_index_dump(gmx_ana_index_t
*g
, int i
, int maxn
);
273 /** Checks whether all indices are between 0 and \p natoms. */
275 gmx_ana_index_check(gmx_ana_index_t
*g
, int natoms
);
276 /** Checks whether an index group is sorted. */
278 gmx_ana_index_check_sorted(gmx_ana_index_t
*g
);
281 /*! \name Functions for set operations on gmx_ana_index_t
284 /** Sorts the indices within an index group. */
286 gmx_ana_index_sort(gmx_ana_index_t
*g
);
287 /** Checks whether two index groups are equal. */
289 gmx_ana_index_equals(gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
290 /** Checks whether a sorted index group contains another sorted index group. */
292 gmx_ana_index_contains(gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
294 /** Calculates the intersection between two sorted index groups. */
296 gmx_ana_index_intersection(gmx_ana_index_t
*dest
,
297 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
298 /** Calculates the set difference between two sorted index groups. */
300 gmx_ana_index_difference(gmx_ana_index_t
*dest
,
301 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
302 /** Calculates the size of the difference between two sorted index groups. */
304 gmx_ana_index_difference_size(gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
305 /** Calculates the union of two sorted index groups. */
307 gmx_ana_index_union(gmx_ana_index_t
*dest
,
308 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
309 /** Merges two distinct sorted index groups. */
311 gmx_ana_index_merge(gmx_ana_index_t
*dest
,
312 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
313 /** Calculates the intersection and the difference in one call. */
315 gmx_ana_index_partition(gmx_ana_index_t
*dest1
, gmx_ana_index_t
*dest2
,
316 gmx_ana_index_t
*src
, gmx_ana_index_t
*g
);
319 /*! \name Functions for handling gmx_ana_indexmap_t and related things
322 /** Partition a group based on topology information. */
324 gmx_ana_index_make_block(t_blocka
*t
, t_topology
*top
, gmx_ana_index_t
*g
,
325 e_index_t type
, gmx_bool bComplete
);
326 /** Checks whether a group consists of full blocks. */
328 gmx_ana_index_has_full_blocks(gmx_ana_index_t
*g
, t_block
*b
);
329 /** Checks whether a group consists of full blocks. */
331 gmx_ana_index_has_full_ablocks(gmx_ana_index_t
*g
, t_blocka
*b
);
332 /** Checks whether a group consists of full residues/molecules. */
334 gmx_ana_index_has_complete_elems(gmx_ana_index_t
*g
, e_index_t type
, t_topology
*top
);
336 /** Initializes an empty index group mapping. */
338 gmx_ana_indexmap_clear(gmx_ana_indexmap_t
*m
);
339 /** Reserves memory for an index group mapping. */
341 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t
*m
, int nr
, int isize
);
342 /** Initializes an index group mapping. */
345 gmx_ana_indexmap_init(gmx_ana_indexmap_t
*m
, gmx_ana_index_t
*g
,
346 t_topology
*top
, e_index_t type
);
347 /** Sets an index group mapping to be static. */
349 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t
*m
, t_blocka
*b
);
350 /** Frees memory allocated for index group mapping. */
352 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t
*m
);
353 /** Makes a deep copy of an index group mapping. */
355 gmx_ana_indexmap_copy(gmx_ana_indexmap_t
*dest
, gmx_ana_indexmap_t
*src
, gmx_bool bFirst
);
356 /** Updates an index group mapping. */
359 gmx_ana_indexmap_update(gmx_ana_indexmap_t
*m
, gmx_ana_index_t
*g
, gmx_bool bMaskOnly
);