2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
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.
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.
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.
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.
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.
36 * \brief API for handling index files and index groups.
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.
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.
57 * \author Teemu Murtola <teemu.murtola@gmail.com>
59 * \ingroup module_selection
61 #ifndef GMX_SELECTION_INDEXUTIL_H
62 #define GMX_SELECTION_INDEXUTIL_H
68 #include "gromacs/topology/block.h"
77 /** Stores a set of index groups. */
78 struct gmx_ana_indexgrps_t
;
81 * Specifies the type of index partition or index mapping in several contexts.
83 * \see gmx_ana_index_make_block(), gmx_ana_indexmap_init()
87 INDEX_UNKNOWN
, /**< Unknown index type.*/
88 INDEX_ATOM
, /**< Each atom in a separate block.*/
89 INDEX_RES
, /**< Each residue in a separate block.*/
90 INDEX_MOL
, /**< Each molecule in a separate block.*/
91 INDEX_ALL
/**< All atoms in a single block.*/
95 * Stores a single index group.
97 struct gmx_ana_index_t
99 /** Number of atoms. */
101 /** List of atoms. */
103 /** Number of items allocated for \p index. */
108 * Data structure for calculating index group mappings.
110 struct gmx_ana_indexmap_t
112 /** Type of the mapping. */
115 * Current reference IDs.
117 * This array provides a mapping from the current index group (last given
118 * to gmx_ana_indexmap_update()) to the blocks in \p b, i.e., the
119 * original index group used in gmx_ana_indexmap_init().
120 * The mapping is zero-based.
121 * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), the indices
122 * for blocks not present in the current group are set to -1, otherwise
123 * they are removed completely and the \p nr field updated.
127 * Current mapped IDs.
129 * This array provides IDs for the current index group. Instead of a
130 * zero-based mapping that \p refid provides, the values from the \p orgid
131 * array are used, thus allowing the mapping to be customized.
132 * In other words, `mapid[i] = orgid[refid[i]]`.
133 * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), this array
138 * Mapped block structure.
140 * A block structure that corresponds to the current index group.
141 * \c mapb.nra and \c mapb.a correspond to the last mapped index group.
146 * Customizable ID numbers for the original blocks.
148 * This array has \p b.nr elements, each defining an original ID number for
149 * a block in \p b (i.e., in the original group passed to
150 * gmx_ana_indexmap_init()).
151 * These are initialized in gmx_ana_indexmap_init() based on the type:
152 * - \ref INDEX_ATOM : the atom indices
153 * - \ref INDEX_RES : the residue indices
154 * - \ref INDEX_MOL : the molecule indices
156 * All the above numbers are zero-based.
157 * After gmx_ana_indexmap_init(), the caller is free to change these values
158 * if the above are not appropriate.
159 * The mapped values can be read through \p mapid.
164 * Block data that defines the mapping (internal use only).
166 * The data is initialized by gmx_ana_indexmap_init() and is not changed
168 * Hence, it cannot be directly applied to the index group passed to
169 * gmx_ana_indexmap_update() unless \p bMaskOnly was specified or the
170 * index group is identical to the one provided to gmx_ana_indexmap_init().
174 * true if the current reference IDs are for the whole group (internal use only).
176 * This is used internally to optimize the evaluation such that
177 * gmx_ana_indexmap_update() does not take any time if the group is
184 /*! \name Functions for handling gmx_ana_indexgrps_t
187 /** Reads index groups from a file or constructs them from topology. */
189 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t
**g
, gmx_mtop_t
*top
,
191 /** Frees memory allocated for index groups. */
193 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t
*g
);
194 /** Returns true if the index group structure is emtpy. */
196 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t
*g
);
198 /** Returns a pointer to an index group. */
200 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t
*g
, int n
);
201 /** Extracts a single index group. */
203 gmx_ana_indexgrps_extract(gmx_ana_index_t
*dest
, std::string
*destName
,
204 gmx_ana_indexgrps_t
*src
, int n
);
205 /** Finds and extracts a single index group by name. */
207 gmx_ana_indexgrps_find(gmx_ana_index_t
*dest
, std::string
*destName
,
208 gmx_ana_indexgrps_t
*src
, const char *name
);
210 /** Writes out a list of index groups. */
212 gmx_ana_indexgrps_print(gmx::TextWriter
*writer
, gmx_ana_indexgrps_t
*g
, int maxn
);
215 /*! \name Functions for handling gmx_ana_index_t
218 /** Reserves memory to store an index group of size \p isize. */
220 gmx_ana_index_reserve(gmx_ana_index_t
*g
, int isize
);
221 /** Frees any memory not necessary to hold the current contents. */
223 gmx_ana_index_squeeze(gmx_ana_index_t
*g
);
224 /** Initializes an empty index group. */
226 gmx_ana_index_clear(gmx_ana_index_t
*g
);
227 /** Constructs a \c gmx_ana_index_t from given values. */
229 gmx_ana_index_set(gmx_ana_index_t
*g
, int isize
, int *index
, int nalloc
);
230 /** Creates a simple index group from the first to the \p natoms'th atom. */
232 gmx_ana_index_init_simple(gmx_ana_index_t
*g
, int natoms
);
233 /** Frees memory allocated for an index group. */
235 gmx_ana_index_deinit(gmx_ana_index_t
*g
);
236 /** Copies a \c gmx_ana_index_t. */
238 gmx_ana_index_copy(gmx_ana_index_t
*dest
, gmx_ana_index_t
*src
, bool bAlloc
);
240 /** Writes out the contents of a index group. */
242 gmx_ana_index_dump(gmx::TextWriter
*writer
, gmx_ana_index_t
*g
, int maxn
);
245 * Returns maximum atom index that appears in an index group.
247 * \param[in] g Index group to query.
248 * \returns Largest atom index that appears in \p g, or zero if \p g is empty.
251 gmx_ana_index_get_max_index(gmx_ana_index_t
*g
);
252 /** Checks whether an index group is sorted. */
254 gmx_ana_index_check_sorted(gmx_ana_index_t
*g
);
256 * Checks whether an index group has atoms from a defined range.
258 * \param[in] g Index group to check.
259 * \param[in] natoms Largest atom number allowed.
260 * \returns true if all atoms in the index group are in the
261 * range 0 to \p natoms (i.e., no atoms over \p natoms are referenced).
264 gmx_ana_index_check_range(gmx_ana_index_t
*g
, int natoms
);
267 /*! \name Functions for set operations on gmx_ana_index_t
270 /** Sorts the indices within an index group. */
272 gmx_ana_index_sort(gmx_ana_index_t
*g
);
274 * Removes duplicates from a sorted index group.
276 * \param[in,out] g Index group to be processed.
279 gmx_ana_index_remove_duplicates(gmx_ana_index_t
*g
);
280 /** Checks whether two index groups are equal. */
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. */
285 gmx_ana_index_contains(gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
287 /** Calculates the intersection between two sorted index groups. */
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. */
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. */
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. */
300 gmx_ana_index_union(gmx_ana_index_t
*dest
,
301 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
303 * Calculates the union of two index groups, where the second group may not be sorted.
305 * \param[out] dest Output index group (the union of \p a and \p b).
306 * \param[in] a First index group (must be sorted).
307 * \param[in] b Second index group.
309 * \p a and \p b can have common items.
310 * \p dest can equal \p a or \p b.
313 gmx_ana_index_union_unsorted(gmx_ana_index_t
*dest
,
314 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
315 /** Merges two distinct sorted index groups. */
317 gmx_ana_index_merge(gmx_ana_index_t
*dest
,
318 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
);
319 /** Calculates the intersection and the difference in one call. */
321 gmx_ana_index_partition(gmx_ana_index_t
*dest1
, gmx_ana_index_t
*dest2
,
322 gmx_ana_index_t
*src
, gmx_ana_index_t
*g
);
325 /*! \name Functions for handling gmx_ana_indexmap_t and related things
328 /** Partition a group based on topology information. */
330 gmx_ana_index_make_block(t_blocka
*t
, const gmx_mtop_t
*top
, gmx_ana_index_t
*g
,
331 e_index_t type
, bool bComplete
);
332 /** Checks whether a group consists of full blocks. */
334 gmx_ana_index_has_full_blocks(const gmx_ana_index_t
*g
, const t_block
*b
);
335 /** Checks whether a group consists of full blocks. */
337 gmx_ana_index_has_full_ablocks(gmx_ana_index_t
*g
, t_blocka
*b
);
338 /** Checks whether a group consists of full residues/molecules. */
340 gmx_ana_index_has_complete_elems(gmx_ana_index_t
*g
, e_index_t type
,
341 const gmx_mtop_t
*top
);
343 /** Initializes an empty index group mapping. */
345 gmx_ana_indexmap_clear(gmx_ana_indexmap_t
*m
);
346 /** Reserves memory for an index group mapping. */
348 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t
*m
, int nr
, int isize
);
349 /** Initializes an index group mapping. */
351 gmx_ana_indexmap_init(gmx_ana_indexmap_t
*m
, gmx_ana_index_t
*g
,
352 const gmx_mtop_t
*top
, e_index_t type
);
354 * Initializes `orgid` entries based on topology grouping.
356 * \param[in,out] m Mapping structure to use/initialize.
357 * \param[in] top Topology structure
358 * (can be NULL if not required for \p type).
359 * \param[in] type Type of groups to use.
360 * \returns The number of groups of type \p type that were present in \p m.
361 * \throws InconsistentInputError if the blocks in \p m do not have a unique
362 * group (e.g., contain atoms from multiple residues with `type == INDEX_RES`).
364 * By default, the gmx_ana_indexmap_t::orgid fields are initialized to
365 * atom/residue/molecule indices from the topology (see documentation for the
366 * struct for more details).
367 * This function can be used to set the field to a zero-based group index
368 * instead. The first block will always get `orgid[0] = 0`, and all following
369 * blocks that belong to the same residue/molecule (\p type) will get the same
370 * index. Each time a block does not belong to the same group, it will get the
371 * next available number.
372 * If `type == INDEX_ATOM`, the `orgid` field is initialized as 0, 1, ...,
373 * independent of whether the blocks are single atoms or not.
375 * Strong exception safety guarantee.
378 gmx_ana_indexmap_init_orgid_group(gmx_ana_indexmap_t
*m
, const gmx_mtop_t
*top
,
380 /** Sets an index group mapping to be static. */
382 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t
*m
, t_blocka
*b
);
383 /** Frees memory allocated for index group mapping. */
385 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t
*m
);
386 /** Makes a deep copy of an index group mapping. */
388 gmx_ana_indexmap_copy(gmx_ana_indexmap_t
*dest
, gmx_ana_indexmap_t
*src
, bool bFirst
);
389 /** Updates an index group mapping. */
391 gmx_ana_indexmap_update(gmx_ana_indexmap_t
*m
, gmx_ana_index_t
*g
, bool bMaskOnly
);