Sort all includes in src/gromacs
[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,2013,2014, 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.
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 <cstdio>
65
66 #include <string>
67
68 #include "gromacs/legacyheaders/types/simple.h"
69 #include "gromacs/topology/block.h"
70
71 struct t_topology;
72
73 /** Stores a set of index groups. */
74 struct gmx_ana_indexgrps_t;
75
76 /*! \brief
77  * Specifies the type of index partition or index mapping in several contexts.
78  *
79  * \see gmx_ana_index_make_block(), gmx_ana_indexmap_init()
80  */
81 typedef enum
82 {
83     INDEX_UNKNOWN, /**< Unknown index type.*/
84     INDEX_ATOM,    /**< Each atom in a separate block.*/
85     INDEX_RES,     /**< Each residue in a separate block.*/
86     INDEX_MOL,     /**< Each molecule in a separate block.*/
87     INDEX_ALL      /**< All atoms in a single block.*/
88 } e_index_t;
89
90 /*! \brief
91  * Stores a single index group.
92  */
93 struct gmx_ana_index_t
94 {
95     /** Number of atoms. */
96     int                 isize;
97     /** List of atoms. */
98     atom_id            *index;
99     /** Number of items allocated for \p index. */
100     int                 nalloc_index;
101 };
102
103 /*! \brief
104  * Data structure for calculating index group mappings.
105  */
106 struct gmx_ana_indexmap_t
107 {
108     /** Type of the mapping. */
109     e_index_t           type;
110     /*! \brief
111      * Current reference IDs.
112      *
113      * This array provides a mapping from the current index group (last given
114      * to gmx_ana_indexmap_update()) to the blocks in \p b, i.e., the
115      * original index group used in gmx_ana_indexmap_init().
116      * The mapping is zero-based.
117      * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), the indices
118      * for blocks not present in the current group are set to -1, otherwise
119      * they are removed completely and the \p nr field updated.
120      */
121     int                *refid;
122     /*! \brief
123      * Current mapped IDs.
124      *
125      * This array provides an arbitrary mapping from the current index group
126      * to the original index group. Instead of a zero-based mapping, the
127      * values from the \p orgid array are used. That is,
128      * \c mapid[i]=orgid[refid[i]].
129      * If \p bMaskOnly is provided to gmx_ana_indexmap_update(), this array
130      * equals \p orgid.
131      */
132     int                *mapid;
133     /*! \brief
134      * Mapped block structure.
135      *
136      * A block structure that corresponds to the current index group.
137      * \c mapb.nra and \c mapb.a correspond to the last mapped index group.
138      */
139     t_blocka            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 };
177
178
179 /*! \name Functions for handling gmx_ana_indexgrps_t
180  */
181 /*@{*/
182 /** Reads index groups from a file or constructs them from topology. */
183 void
184 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
185                        const char *fnm);
186 /** Frees memory allocated for index groups. */
187 void
188 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g);
189 /** Returns true if the index group structure is emtpy. */
190 bool
191 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g);
192
193 /** Returns a pointer to an index group. */
194 gmx_ana_index_t *
195 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n);
196 /** Extracts a single index group. */
197 bool
198 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
199                           gmx_ana_indexgrps_t *src, int n);
200 /** Finds and extracts a single index group by name. */
201 bool
202 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
203                        gmx_ana_indexgrps_t *src, const char *name);
204
205 /** Writes out a list of index groups. */
206 void
207 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn);
208 /*@}*/
209
210 /*! \name Functions for handling gmx_ana_index_t
211  */
212 /*@{*/
213 /** Reserves memory to store an index group of size \p isize. */
214 void
215 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize);
216 /** Frees any memory not necessary to hold the current contents. */
217 void
218 gmx_ana_index_squeeze(gmx_ana_index_t *g);
219 /** Initializes an empty index group. */
220 void
221 gmx_ana_index_clear(gmx_ana_index_t *g);
222 /** Constructs a \c gmx_ana_index_t from given values. */
223 void
224 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, int nalloc);
225 /** Creates a simple index group from the first to the \p natoms'th atom. */
226 void
227 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms);
228 /** Frees memory allocated for an index group. */
229 void
230 gmx_ana_index_deinit(gmx_ana_index_t *g);
231 /** Copies a \c gmx_ana_index_t. */
232 void
233 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc);
234
235 /** Writes out the contents of a index group. */
236 void
237 gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int maxn);
238
239 /*! \brief
240  * Returns maximum atom index that appears in an index group.
241  *
242  * \param[in]  g      Index group to query.
243  * \returns    Largest atom index that appears in \p g, or zero if \p g is empty.
244  */
245 int
246 gmx_ana_index_get_max_index(gmx_ana_index_t *g);
247 /** Checks whether an index group is sorted. */
248 bool
249 gmx_ana_index_check_sorted(gmx_ana_index_t *g);
250 /*! \brief
251  * Checks whether an index group has atoms from a defined range.
252  *
253  * \param[in]  g      Index group to check.
254  * \param[in]  natoms Largest atom number allowed.
255  * \returns    true if all atoms in the index group are in the
256  *     range 0 to \p natoms (i.e., no atoms over \p natoms are referenced).
257  */
258 bool
259 gmx_ana_index_check_range(gmx_ana_index_t *g, int natoms);
260 /*@}*/
261
262 /*! \name Functions for set operations on gmx_ana_index_t
263  */
264 /*@{*/
265 /** Sorts the indices within an index group. */
266 void
267 gmx_ana_index_sort(gmx_ana_index_t *g);
268 /** Checks whether two index groups are equal. */
269 bool
270 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b);
271 /** Checks whether a sorted index group contains another sorted index group. */
272 bool
273 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b);
274
275 /** Calculates the intersection between two sorted index groups. */
276 void
277 gmx_ana_index_intersection(gmx_ana_index_t *dest,
278                            gmx_ana_index_t *a, gmx_ana_index_t *b);
279 /** Calculates the set difference between two sorted index groups. */
280 void
281 gmx_ana_index_difference(gmx_ana_index_t *dest,
282                          gmx_ana_index_t *a, gmx_ana_index_t *b);
283 /** Calculates the size of the difference between two sorted index groups. */
284 int
285 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b);
286 /** Calculates the union of two sorted index groups. */
287 void
288 gmx_ana_index_union(gmx_ana_index_t *dest,
289                     gmx_ana_index_t *a, gmx_ana_index_t *b);
290 /** Merges two distinct sorted index groups. */
291 void
292 gmx_ana_index_merge(gmx_ana_index_t *dest,
293                     gmx_ana_index_t *a, gmx_ana_index_t *b);
294 /** Calculates the intersection and the difference in one call. */
295 void
296 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
297                         gmx_ana_index_t *src, gmx_ana_index_t *g);
298 /*@}*/
299
300 /*! \name Functions for handling gmx_ana_indexmap_t and related things
301  */
302 /*@{*/
303 /** Partition a group based on topology information. */
304 void
305 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
306                          e_index_t type, bool bComplete);
307 /** Checks whether a group consists of full blocks. */
308 bool
309 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b);
310 /** Checks whether a group consists of full blocks. */
311 bool
312 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b);
313 /** Checks whether a group consists of full residues/molecules. */
314 bool
315 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type, t_topology *top);
316
317 /** Initializes an empty index group mapping. */
318 void
319 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m);
320 /** Reserves memory for an index group mapping. */
321 void
322 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize);
323 /** Initializes an index group mapping. */
324 void
325 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
326                       t_topology *top, e_index_t type);
327 /** Sets an index group mapping to be static. */
328 void
329 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b);
330 /** Frees memory allocated for index group mapping. */
331 void
332 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m);
333 /** Makes a deep copy of an index group mapping. */
334 void
335 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst);
336 /** Updates an index group mapping. */
337 void
338 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g, bool bMaskOnly);
339 /*@}*/
340
341 #endif