Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / selection / indexutil.cpp
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
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.
13
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.
18  *
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.
25  *
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.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \internal \file
32  * \brief
33  * Implements functions in indexutil.h.
34  *
35  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36  * \ingroup module_selection
37  */
38 #include "gromacs/legacyheaders/index.h"
39 #include "gromacs/legacyheaders/gmx_fatal.h"
40 #include "gromacs/legacyheaders/smalloc.h"
41 #include "gromacs/legacyheaders/string2.h"
42 #include "gromacs/legacyheaders/typedefs.h"
43
44 #include "gromacs/selection/indexutil.h"
45
46 /********************************************************************
47  * gmx_ana_indexgrps_t functions
48  ********************************************************************/
49
50 /*! \internal \brief
51  * Stores a set of index groups.
52  */
53 struct gmx_ana_indexgrps_t
54 {
55     /** Number of index groups. */
56     int                 nr;
57     /** Array of index groups. */
58     gmx_ana_index_t    *g;
59 };
60
61 /*!
62  * \param[out] g     Index group structure.
63  * \param[in]  ngrps Number of groups for which memory is allocated.
64  */
65 void
66 gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps)
67 {
68     snew(*g, 1);
69     (*g)->nr = ngrps;
70     snew((*g)->g,    ngrps);
71 }
72
73 /*!
74  * \param[out] g     Index group structure.
75  * \param[in]  top   Topology structure.
76  * \param[in]  fnm   File name for the index file.
77  *   Memory is automatically allocated.
78  *
79  * One or both of \p top or \p fnm can be NULL.
80  * If \p top is NULL, an index file is required and the groups are read
81  * from the file (uses Gromacs routine init_index()).
82  * If \p fnm is NULL, default groups are constructed based on the
83  * topology (uses Gromacs routine analyse()).
84  * If both are null, the index group structure is initialized empty.
85  */
86 void
87 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
88                        const char *fnm)
89 {
90     t_blocka *block = NULL;
91     char    **names = NULL;
92     int       i, j;
93
94     if (fnm)
95     {
96         block = init_index(fnm, &names);
97     }
98     else if (top)
99     {
100         block = new_blocka();
101         analyse(&top->atoms, block, &names, FALSE, FALSE);
102     }
103     else
104     {
105         snew(*g, 1);
106         (*g)->nr    = 0;
107         (*g)->g     = NULL;
108         return;
109     }
110
111     gmx_ana_indexgrps_alloc(g, block->nr);
112     for (i = 0; i < block->nr; ++i)
113     {
114         gmx_ana_index_t *grp = &(*g)->g[i];
115
116         grp->isize = block->index[i+1] - block->index[i];
117         snew(grp->index, grp->isize);
118         for (j = 0; j < grp->isize; ++j)
119         {
120             grp->index[j] = block->a[block->index[i]+j];
121         }
122         grp->name = names[i];
123         grp->nalloc_index = grp->isize;
124     }
125
126     done_blocka(block);
127     sfree(block);
128     sfree(names);
129 }
130
131 /*!
132  * \param[in] g  Index groups structure.
133  *
134  * The pointer \p g is invalid after the call.
135  */
136 void
137 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
138 {
139     int  i;
140
141     if (g->nr == 0)
142     {
143         sfree(g);
144         return;
145     }
146     for (i = 0; i < g->nr; ++i)
147     {
148         gmx_ana_index_deinit(&g->g[i]);
149     }
150     sfree(g->g);
151     g->nr    = 0;
152     g->g     = NULL;
153     sfree(g);
154 }
155
156 /*!
157  * \param[out] dest Destination index groups.
158  * \param[in]  src  Source index groups.
159  *
160  * A deep copy is made for all fields, including the group names.
161  */
162 void
163 gmx_ana_indexgrps_clone(gmx_ana_indexgrps_t **dest, gmx_ana_indexgrps_t *src)
164 {
165     int g;
166
167     gmx_ana_indexgrps_alloc(dest, src->nr);
168     for (g = 0; g < src->nr; ++g)
169     {
170         gmx_ana_index_copy(&(*dest)->g[g], &src->g[g], true);
171     }
172 }
173
174 /*!
175  * \param[out] g     Index group structure.
176  * \returns    true if \p g is empty, i.e., has 0 index groups.
177  */
178 bool
179 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g)
180 {
181     return g->nr == 0;
182 }
183
184 /*!
185  * \param[in]  g     Index groups structure.
186  * \param[in]  n     Index group number to get.
187  * \returns    Pointer to the \p n'th index group in \p g.
188  *
189  * The returned pointer should not be freed.
190  */
191 gmx_ana_index_t *
192 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n)
193 {
194     if (n < 0 || n >= g->nr)
195     {
196         return NULL;
197     }
198     return &g->g[n];
199 }
200
201 /*!
202  * \param[out] dest Output structure.
203  * \param[in]  src  Input index groups.
204  * \param[in]  n    Number of the group to extract.
205  * \returns true if \p n is a valid group in \p src, false otherwise.
206  */
207 bool
208 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, int n)
209 {
210     if (n < 0 || n >= src->nr)
211     {
212         dest->isize = 0;
213         return false;
214     }
215
216     gmx_ana_index_copy(dest, &src->g[n], true);
217     return true;
218 }
219
220 /*!
221  * \param[out] dest Output structure.
222  * \param[in]  src  Input index groups.
223  * \param[in]  name Name (or part of the name) of the group to extract.
224  * \returns true if \p name is a valid group in \p src, false otherwise.
225  *
226  * Uses the Gromacs routine find_group() to find the actual group;
227  * the comparison is case-insensitive.
228  */
229 bool
230 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, char *name)
231 {
232     int    i;
233     char **names;
234
235     snew(names, src->nr);
236     for (i = 0; i < src->nr; ++i)
237     {
238         names[i] = src->g[i].name;
239     }
240     i = find_group(name, src->nr, names);
241     sfree(names);
242     if (i == NOTSET)
243     {
244         dest->isize = 0;
245         return false;
246     }
247
248     return gmx_ana_indexgrps_extract(dest, src, i);
249 }
250
251 /*!
252  * \param[in]  fp     Where to print the output.
253  * \param[in]  g      Index groups to print.
254  * \param[in]  maxn   Maximum number of indices to print
255  *      (-1 = print all, 0 = print only names).
256  */
257 void
258 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn)
259 {
260     int  i;
261
262     for (i = 0; i < g->nr; ++i)
263     {
264         fprintf(fp, " %2d: ", i);
265         gmx_ana_index_dump(fp, &g->g[i], i, maxn);
266     }
267 }
268
269 /********************************************************************
270  * gmx_ana_index_t functions
271  ********************************************************************/
272
273 /*!
274  * \param[in,out] g      Index group structure.
275  * \param[in]     isize  Maximum number of atoms to reserve space for.
276  */
277 void
278 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
279 {
280     if (g->nalloc_index < isize)
281     {
282         srenew(g->index, isize);
283         g->nalloc_index = isize;
284     }
285 }
286
287 /*!
288  * \param[in,out] g      Index group structure.
289  *
290  * Resizes the memory allocated for holding the indices such that the
291  * current contents fit.
292  */
293 void
294 gmx_ana_index_squeeze(gmx_ana_index_t *g)
295 {
296     srenew(g->index, g->isize);
297     g->nalloc_index = g->isize;
298 }
299
300 /*!
301  * \param[out] g      Output structure.
302  *
303  * Any contents of \p g are discarded without freeing.
304  */
305 void
306 gmx_ana_index_clear(gmx_ana_index_t *g)
307 {
308     g->isize        = 0;
309     g->index        = NULL;
310     g->name         = NULL;
311     g->nalloc_index = 0;
312 }
313
314 /*!
315  * \param[out] g      Output structure.
316  * \param[in]  isize  Number of atoms in the new group.
317  * \param[in]  index  Array of \p isize atoms (can be NULL if \p isize is 0).
318  * \param[in]  name   Name for the new group (can be NULL).
319  * \param[in]  nalloc Number of elements allocated for \p index
320  *   (if 0, \p index is not freed in gmx_ana_index_deinit())
321  *
322  * No copy if \p index is made.
323  */
324 void
325 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, char *name,
326                   int nalloc)
327 {
328     g->isize        = isize;
329     g->index        = index;
330     g->name         = name;
331     g->nalloc_index = nalloc;
332 }
333
334 /*!
335  * \param[out] g      Output structure.
336  * \param[in]  natoms Number of atoms.
337  * \param[in]  name   Name for the new group (can be NULL).
338  */
339 void
340 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms, char *name)
341 {
342     int  i;
343
344     g->isize = natoms;
345     snew(g->index, natoms);
346     for (i = 0; i < natoms; ++i)
347     {
348         g->index[i] = i;
349     }
350     g->name = name;
351     g->nalloc_index = natoms;
352 }
353
354 /*!
355  * \param[in] g  Index group structure.
356  *
357  * The pointer \p g is not freed.
358  */
359 void
360 gmx_ana_index_deinit(gmx_ana_index_t *g)
361 {
362     if (g->nalloc_index > 0)
363     {
364         sfree(g->index);
365     }
366     sfree(g->name);
367     gmx_ana_index_clear(g);
368 }
369
370 /*!
371  * \param[out] dest   Destination index group.
372  * \param[in]  src    Source index group.
373  * \param[in]  bAlloc If true, memory is allocated at \p dest; otherwise,
374  *   it is assumed that enough memory has been allocated for index.
375  *
376  * A deep copy of the name is only made if \p bAlloc is true.
377  */
378 void
379 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
380 {
381     dest->isize = src->isize;
382     if (dest->isize > 0)
383     {
384         if (bAlloc)
385         {
386             snew(dest->index, dest->isize);
387             dest->nalloc_index = dest->isize;
388         }
389         memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
390     }
391     if (bAlloc && src->name)
392     {
393         dest->name = strdup(src->name);
394     }
395     else if (bAlloc || src->name)
396     {
397         dest->name = src->name;
398     }
399 }
400
401 /*!
402  * \param[in]  fp     Where to print the output.
403  * \param[in]  g      Index group to print.
404  * \param[in]  i      Group number to use if the name is NULL.
405  * \param[in]  maxn   Maximum number of indices to print (-1 = print all).
406  */
407 void
408 gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int i, int maxn)
409 {
410     int  j, n;
411
412     if (g->name)
413     {
414         fprintf(fp, "\"%s\"", g->name);
415     }
416     else
417     {
418         fprintf(fp, "Group %d", i+1);
419     }
420     fprintf(fp, " (%d atoms)", g->isize);
421     if (maxn != 0)
422     {
423         fprintf(fp, ":");
424         n = g->isize;
425         if (maxn >= 0 && n > maxn)
426         {
427             n = maxn;
428         }
429         for (j = 0; j < n; ++j)
430         {
431             fprintf(fp, " %d", g->index[j]+1);
432         }
433         if (n < g->isize)
434         {
435             fprintf(fp, " ...");
436         }
437     }
438     fprintf(fp, "\n");
439 }
440
441 /*!
442  * \param[in]  g      Input index group.
443  * \param[in]  natoms Number of atoms to check against.
444  *
445  * If any atom index in the index group is less than zero or >= \p natoms,
446  * gmx_fatal() is called.
447  */
448 void
449 gmx_ana_index_check(gmx_ana_index_t *g, int natoms)
450 {
451     int  j;
452
453     for (j = 0; j < g->isize; ++j)
454     {
455         if (g->index[j] >= natoms)
456         {
457             gmx_fatal(FARGS,"Atom index (%d) in index group %s (%d atoms) "
458               "larger than number of atoms in trajectory (%d atoms)",
459               g->index[j], g->name, g->isize, natoms);
460         }
461         else if (g->index[j] < 0)
462         {
463             gmx_fatal(FARGS,"Atom index (%d) in index group %s (%d atoms) "
464               "is less than zero",
465               g->index[j], g->name, g->isize);
466         }
467     }
468 }
469
470 /*!
471  * \param[in]  g      Index group to check.
472  * \returns    true if the index group is sorted and has no duplicates,
473  *   false otherwise.
474  */
475 bool
476 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
477 {
478     int  i;
479
480     for (i = 0; i < g->isize-1; ++i)
481     {
482         if (g->index[i+1] <= g->index[i])
483         {
484             return false;
485         }
486     }
487     return true;
488 }
489
490 /********************************************************************
491  * Set operations
492  ********************************************************************/
493
494 /** Helper function for gmx_ana_index_sort(). */
495 static int
496 cmp_atomid(const void *a, const void *b)
497 {
498     if (*(atom_id *)a < *(atom_id *)b) return -1;
499     if (*(atom_id *)a > *(atom_id *)b) return 1;
500     return 0;
501 }
502
503 /*!
504  * \param[in,out] g  Index group to be sorted.
505  */
506 void
507 gmx_ana_index_sort(gmx_ana_index_t *g)
508 {
509     qsort(g->index, g->isize, sizeof(*g->index), cmp_atomid);
510 }
511
512 /*!
513  * \param[in]  a      Index group to check.
514  * \param[in]  b      Index group to check.
515  * \returns    true if \p a and \p b are equal, false otherwise.
516  */
517 bool
518 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
519 {
520     int  i;
521
522     if (a->isize != b->isize)
523     {
524         return false;
525     }
526     for (i = 0; i < a->isize; ++i)
527     {
528         if (a->index[i] != b->index[i])
529         {
530             return false;
531         }
532     }
533     return true;
534 }
535
536 /*!
537  * \param[in]  a      Index group to check against.
538  * \param[in]  b      Index group to check.
539  * \returns    true if \p b is contained in \p a,
540  *   false otherwise.
541  *
542  * If the elements are not in the same order in both groups, the function
543  * fails. However, the groups do not need to be sorted.
544  */
545 bool
546 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
547 {
548     int  i, j;
549
550     for (i = j = 0; j < b->isize; ++i, ++j) {
551         while (i < a->isize && a->index[i] != b->index[j])
552         {
553             ++i;
554         }
555         if (i == a->isize)
556         {
557             return false;
558         }
559     }
560     return true;
561 }
562
563 /*!
564  * \param[out] dest Output index group (the intersection of \p a and \p b).
565  * \param[in]  a    First index group.
566  * \param[in]  b    Second index group.
567  *
568  * \p dest can be the same as \p a or \p b.
569  */
570 void
571 gmx_ana_index_intersection(gmx_ana_index_t *dest,
572                            gmx_ana_index_t *a, gmx_ana_index_t *b)
573 {
574     int i, j, k;
575
576     for (i = j = k = 0; i < a->isize && j < b->isize; ++i) {
577         while (j < b->isize && b->index[j] < a->index[i])
578         {
579             ++j;
580         }
581         if (j < b->isize && b->index[j] == a->index[i])
582         {
583             dest->index[k++] = b->index[j++];
584         }
585     }
586     dest->isize = k;
587 }
588
589 /*!
590  * \param[out] dest Output index group (the difference \p a - \p b).
591  * \param[in]  a    First index group.
592  * \param[in]  b    Second index group.
593  *
594  * \p dest can equal \p a, but not \p b.
595  */
596 void
597 gmx_ana_index_difference(gmx_ana_index_t *dest,
598                          gmx_ana_index_t *a, gmx_ana_index_t *b)
599 {
600     int i, j, k;
601
602     for (i = j = k = 0; i < a->isize; ++i)
603     {
604         while (j < b->isize && b->index[j] < a->index[i])
605         {
606             ++j;
607         }
608         if (j == b->isize || b->index[j] != a->index[i])
609         {
610             dest->index[k++] = a->index[i];
611         }
612     }
613     dest->isize = k;
614 }
615
616 /*!
617  * \param[in]  a    First index group.
618  * \param[in]  b    Second index group.
619  * \returns    Size of the difference \p a - \p b.
620  */
621 int
622 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
623 {
624     int i, j, k;
625
626     for (i = j = k = 0; i < a->isize; ++i)
627     {
628         while (j < b->isize && b->index[j] < a->index[i])
629         {
630             ++j;
631         }
632         if (j == b->isize || b->index[j] != a->index[i])
633         {
634             ++k;
635         }
636     }
637     return k;
638 }
639
640 /*!
641  * \param[out] dest1 Output group 1 (will equal \p g).
642  * \param[out] dest2 Output group 2 (will equal \p src - \p g).
643  * \param[in]  src   Group to be partitioned.
644  * \param[in]  g     One partition.
645  *
646  * \pre \p g is a subset of \p src and both sets are sorted
647  * \pre \p dest1 has allocated storage to store \p src
648  * \post \p dest1 == \p g
649  * \post \p dest2 == \p src - \p g
650  *
651  * No storage should be allocated for \p dest2; after the call,
652  * \p dest2->index points to the memory allocated for \p dest1
653  * (to a part that is not used by \p dest1).
654  *
655  * The calculation can be performed in-place by setting \p dest1 equal to
656  * \p src.
657  */
658 void
659 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
660                         gmx_ana_index_t *src, gmx_ana_index_t *g)
661 {
662     int i, j, k;
663
664     dest2->index = dest1->index + g->isize;
665     dest2->isize = src->isize - g->isize;
666     for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
667     {
668         while (j >= 0 && src->index[j] != g->index[i])
669         {
670             dest2->index[k--] = src->index[j--];
671         }
672     }
673     while (j >= 0)
674     {
675         dest2->index[k--] = src->index[j--];
676     }
677     gmx_ana_index_copy(dest1, g, false);
678 }
679
680 /*!
681  * \param[out] dest Output index group (the union of \p a and \p b).
682  * \param[in]  a    First index group.
683  * \param[in]  b    Second index group.
684  *
685  * \p a and \p b can have common items.
686  * \p dest can equal \p a or \p b.
687  *
688  * \see gmx_ana_index_merge()
689  */
690 void
691 gmx_ana_index_union(gmx_ana_index_t *dest,
692                     gmx_ana_index_t *a, gmx_ana_index_t *b)
693 {
694     int dsize;
695     int i, j, k;
696
697     dsize = gmx_ana_index_difference_size(b, a);
698     i = a->isize - 1;
699     j = b->isize - 1;
700     dest->isize = a->isize + dsize;
701     for (k = dest->isize - 1; k >= 0; k--)
702     {
703         if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
704         {
705             dest->index[k] = b->index[j--];
706         }
707         else
708         {
709             if (j >= 0 && a->index[i] == b->index[j])
710             {
711                 --j;
712             }
713             dest->index[k] = a->index[i--];
714         }
715     }
716 }
717
718 /*!
719  * \param[out] dest Output index group (the union of \p a and \p b).
720  * \param[in]  a    First index group.
721  * \param[in]  b    Second index group.
722  *
723  * \p a and \p b should not have common items.
724  * \p dest can equal \p a or \p b.
725  *
726  * \see gmx_ana_index_union()
727  */
728 void
729 gmx_ana_index_merge(gmx_ana_index_t *dest,
730                     gmx_ana_index_t *a, gmx_ana_index_t *b)
731 {
732     int i, j, k;
733
734     i = a->isize - 1;
735     j = b->isize - 1;
736     dest->isize = a->isize + b->isize;
737     for (k = dest->isize - 1; k >= 0; k--)
738     {
739         if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
740         {
741             dest->index[k] = b->index[j--];
742         }
743         else
744         {
745             dest->index[k] = a->index[i--];
746         }
747     }
748 }
749
750 /********************************************************************
751  * gmx_ana_indexmap_t and related things
752  ********************************************************************/
753
754 /*!
755  * \param[in,out] t    Output block.
756  * \param[in]  top  Topology structure
757  *   (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
758  *   otherwise).
759  * \param[in]  g    Index group
760  *   (can be NULL if \p type is \ref INDEX_UNKNOWN).
761  * \param[in]  type Type of partitioning to make.
762  * \param[in]  bComplete
763  *   If true, the index group is expanded to include any residue/molecule
764  *   (depending on \p type) that is partially contained in the group.
765  *   If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
766  *
767  * \p m should have been initialized somehow (calloc() is enough) unless
768  * \p type is INDEX_UNKNOWN.
769  * \p g should be sorted.
770  */
771 void
772 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
773                          e_index_t type, bool bComplete)
774 {
775     int      i, j, ai;
776     int      id, cur;
777
778     if (type == INDEX_UNKNOWN)
779     {
780         t->nr           = 1;
781         snew(t->index, 2);
782         t->nalloc_index = 2;
783         t->index[0]     = 0;
784         t->index[1]     = 0;
785         t->nra          = 0;
786         t->a            = NULL;
787         t->nalloc_a     = 0;
788         return;
789     }
790
791     /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
792      * off otherwise. */
793     if (type != INDEX_RES && type != INDEX_MOL)
794     {
795         bComplete = false;
796     }
797     /* Allocate memory for the atom array and fill it unless we are using
798      * completion. */
799     if (bComplete)
800     {
801         t->nra = 0;
802         /* We may allocate some extra memory here because we don't know in
803          * advance how much will be needed. */
804         if (t->nalloc_a < top->atoms.nr)
805         {
806             srenew(t->a, top->atoms.nr);
807             t->nalloc_a = top->atoms.nr;
808         }
809     }
810     else
811     {
812         t->nra      = g->isize;
813         if (t->nalloc_a < g->isize)
814         {
815             srenew(t->a, g->isize);
816             t->nalloc_a = g->isize;
817         }
818         memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
819     }
820
821     /* Allocate memory for the block index. We don't know in advance
822      * how much will be needed, so we allocate some extra and free it in the
823      * end. */
824     if (t->nalloc_index < g->isize + 1)
825     {
826         srenew(t->index, g->isize + 1);
827         t->nalloc_index = g->isize + 1;
828     }
829     /* Clear counters */
830     t->nr = 0;
831     j = 0; /* j is used by residue completion for the first atom not stored */
832     id = cur = -1;
833     for (i = 0; i < g->isize; ++i)
834     {
835         ai = g->index[i];
836         /* Find the ID number of the atom/residue/molecule corresponding to
837          * atom ai. */
838         switch (type)
839         {
840             case INDEX_ATOM:
841                 id = ai;
842                 break;
843             case INDEX_RES:
844                 id = top->atoms.atom[ai].resind;
845                 break;
846             case INDEX_MOL:
847                 while (ai >= top->mols.index[id+1])
848                 {
849                     id++;
850                 }
851                 break;
852             case INDEX_UNKNOWN: /* Should not occur */
853             case INDEX_ALL:
854                 id = 0;
855                 break;
856         }
857         /* If this is the first atom in a new block, initialize the block. */
858         if (id != cur)
859         {
860             if (bComplete)
861             {
862                 /* For completion, we first set the start of the block. */
863                 t->index[t->nr++] = t->nra;
864                 /* And then we find all the atoms that should be included. */
865                 switch (type)
866                 {
867                     case INDEX_RES:
868                         while (top->atoms.atom[j].resind != id)
869                         {
870                             ++j;
871                         }
872                         while (j < top->atoms.nr && top->atoms.atom[j].resind == id)
873                         {
874                             t->a[t->nra++] = j;
875                             ++j;
876                         }
877                         break;
878
879                     case INDEX_MOL:
880                         for (j = top->mols.index[id]; j < top->mols.index[id+1]; ++j)
881                         {
882                             t->a[t->nra++] = j;
883                         }
884                         break;
885
886                     default: /* Should not be reached */
887                         gmx_bug("internal error");
888                         break;
889                 }
890             }
891             else
892             {
893                 /* If not using completion, simply store the start of the block. */
894                 t->index[t->nr++] = i;
895             }
896             cur = id;
897         }
898     }
899     /* Set the end of the last block */
900     t->index[t->nr] = t->nra;
901     /* Free any unnecessary memory */
902     srenew(t->index, t->nr+1);
903     t->nalloc_index = t->nr+1;
904     if (bComplete)
905     {
906         srenew(t->a, t->nra);
907         t->nalloc_a = t->nra;
908     }
909 }
910
911 /*!
912  * \param[in] g   Index group to check.
913  * \param[in] b   Block data to check against.
914  * \returns   true if \p g consists of one or more complete blocks from \p b,
915  *   false otherwise.
916  *
917  * The atoms in \p g are assumed to be sorted.
918  */
919 bool
920 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b)
921 {
922     int  i, j, bi;
923
924     i = bi = 0;
925     /* Each round in the loop matches one block */
926     while (i < g->isize)
927     {
928         /* Find the block that begins with the first unmatched atom */
929         while (bi < b->nr && b->index[bi] != g->index[i])
930         {
931             ++bi;
932         }
933         /* If not found, or if too large, return */
934         if (bi == b->nr || i + b->index[bi+1] -  b->index[bi] > g->isize)
935         {
936             return false;
937         }
938         /* Check that the block matches the index */
939         for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
940         {
941             if (g->index[i] != j)
942             {
943                 return false;
944             }
945         }
946         /* Move the search to the next block */
947         ++bi;
948     }
949     return true;
950 }
951
952 /*!
953  * \param[in] g   Index group to check.
954  * \param[in] b   Block data to check against.
955  * \returns   true if \p g consists of one or more complete blocks from \p b,
956  *   false otherwise.
957  *
958  * The atoms in \p g and \p b->a are assumed to be in the same order.
959  */
960 bool
961 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
962 {
963     int  i, j, bi;
964
965     i = bi = 0;
966     /* Each round in the loop matches one block */
967     while (i < g->isize)
968     {
969         /* Find the block that begins with the first unmatched atom */
970         while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
971         {
972             ++bi;
973         }
974         /* If not found, or if too large, return */
975         if (bi == b->nr || i + b->index[bi+1] -  b->index[bi] > g->isize)
976         {
977             return false;
978         }
979         /* Check that the block matches the index */
980         for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
981         {
982             if (b->a[j] != g->index[i])
983             {
984                 return false;
985             }
986         }
987         /* Move the search to the next block */
988         ++bi;
989     }
990     return true;
991 }
992
993 /*!
994  * \param[in] g     Index group to check.
995  * \param[in] type  Block data to check against.
996  * \param[in] top   Topology data.
997  * \returns   true if \p g consists of one or more complete elements of type
998  *   \p type, false otherwise.
999  *
1000  * If \p type is \ref INDEX_ATOM, the return value is always true.
1001  * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
1002  * always false.
1003  */
1004 bool
1005 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
1006                                  t_topology *top)
1007 {
1008     switch (type)
1009     {
1010         case INDEX_UNKNOWN:
1011         case INDEX_ALL:
1012             return false;
1013
1014         case INDEX_ATOM:
1015             return true;
1016
1017         case INDEX_RES:
1018         {
1019             int      i, ai;
1020             int      id, prev;
1021
1022             prev = -1;
1023             for (i = 0; i < g->isize; ++i)
1024             {
1025                 ai = g->index[i];
1026                 id = top->atoms.atom[ai].resind;
1027                 if (id != prev)
1028                 {
1029                     if (ai > 0 && top->atoms.atom[ai-1].resind == id)
1030                     {
1031                         return false;
1032                     }
1033                     if (i > 0 && g->index[i-1] < top->atoms.nr - 1
1034                         && top->atoms.atom[g->index[i-1]+1].resind == prev)
1035                     {
1036                         return false;
1037                     }
1038                 }
1039                 prev = id;
1040             }
1041             if (g->index[i-1] < top->atoms.nr - 1
1042                 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1043             {
1044                 return false;
1045             }
1046             break;
1047         }
1048
1049         case INDEX_MOL:
1050             return gmx_ana_index_has_full_blocks(g, &top->mols);
1051     }
1052     return true;
1053 }
1054
1055 /*!
1056  * \param[out] m      Output structure.
1057  *
1058  * Any contents of \p m are discarded without freeing.
1059  */
1060 void
1061 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1062 {
1063     m->type              = INDEX_UNKNOWN;
1064     m->nr                = 0;
1065     m->refid             = NULL;
1066     m->mapid             = NULL;
1067     m->mapb.nr           = 0;
1068     m->mapb.index        = NULL;
1069     m->mapb.nalloc_index = 0;
1070     m->orgid             = NULL;
1071     m->b.nr              = 0;
1072     m->b.index           = NULL;
1073     m->b.nra             = 0;
1074     m->b.a               = NULL;
1075     m->b.nalloc_index    = 0;
1076     m->b.nalloc_a        = 0;
1077     m->bStatic           = true;
1078     m->bMapStatic        = true;
1079 }
1080
1081 /*!
1082  * \param[in,out] m      Mapping structure.
1083  * \param[in]     nr     Maximum number of blocks to reserve space for.
1084  * \param[in]     isize  Maximum number of atoms to reserve space for.
1085  */
1086 void
1087 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1088 {
1089     if (m->mapb.nalloc_index < nr + 1)
1090     {
1091         srenew(m->refid,      nr);
1092         srenew(m->mapid,      nr);
1093         srenew(m->orgid,      nr);
1094         srenew(m->mapb.index, nr + 1);
1095         srenew(m->b.index,    nr + 1);
1096         m->mapb.nalloc_index = nr + 1;
1097         m->b.nalloc_index    = nr + 1;
1098     }
1099     if (m->b.nalloc_a < isize)
1100     {
1101         srenew(m->b.a,        isize);
1102         m->b.nalloc_a = isize;
1103     }
1104 }
1105
1106 /*!
1107  * \param[in,out] m    Mapping structure to initialize.
1108  * \param[in]     g    Index group to map
1109  *   (can be NULL if \p type is \ref INDEX_UNKNOWN).
1110  * \param[in]     top  Topology structure
1111  *   (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1112  * \param[in]     type Type of mapping to construct.
1113  *
1114  * Initializes a new index group mapping.
1115  * The index group provided to gmx_ana_indexmap_update() should always be a
1116  * subset of the \p g given here.
1117  *
1118  * \p m should have been initialized somehow (calloc() is enough).
1119  */
1120 void
1121 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1122                       t_topology *top, e_index_t type)
1123 {
1124     int      i, ii, mi;
1125
1126     m->type   = type;
1127     gmx_ana_index_make_block(&m->b, top, g, type, false);
1128     gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1129     m->nr = m->b.nr;
1130     for (i = mi = 0; i < m->nr; ++i)
1131     {
1132         ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1133         switch (type)
1134         {
1135             case INDEX_ATOM:
1136                 m->orgid[i] = ii;
1137                 break;
1138             case INDEX_RES:
1139                 m->orgid[i] = top->atoms.atom[ii].resind;
1140                 break;
1141             case INDEX_MOL:
1142                 while (top->mols.index[mi+1] <= ii)
1143                 {
1144                     ++mi;
1145                 }
1146                 m->orgid[i] = mi;
1147                 break;
1148             case INDEX_ALL:
1149             case INDEX_UNKNOWN:
1150                 m->orgid[i] = 0;
1151                 break;
1152         }
1153     }
1154     for (i = 0; i < m->nr; ++i)
1155     {
1156         m->refid[i] = i;
1157         m->mapid[i] = m->orgid[i];
1158     }
1159     m->mapb.nr = m->nr;
1160     memcpy(m->mapb.index, m->b.index, (m->nr+1)*sizeof(*(m->mapb.index)));
1161     m->bStatic    = true;
1162     m->bMapStatic = true;
1163 }
1164
1165 /*!
1166  * \param[in,out] m    Mapping structure to initialize.
1167  * \param[in]     b    Block information to use for data.
1168  *
1169  * Frees some memory that is not necessary for static index group mappings.
1170  * Internal pointers are set to point to data in \p b; it is the responsibility
1171  * of the caller to ensure that the block information matches the contents of
1172  * the mapping.
1173  * After this function has been called, the index group provided to
1174  * gmx_ana_indexmap_update() should always be the same as \p g given here.
1175  *
1176  * This function breaks modularity of the index group mapping interface in an
1177  * ugly way, but allows reducing memory usage of static selections by a
1178  * significant amount.
1179  */
1180 void
1181 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1182 {
1183     sfree(m->mapid);
1184     m->mapid = m->orgid;
1185     sfree(m->b.index);
1186     m->b.nalloc_index = 0;
1187     m->b.index = b->index;
1188     sfree(m->mapb.index);
1189     m->mapb.nalloc_index = 0;
1190     m->mapb.index = m->b.index;
1191     sfree(m->b.a);
1192     m->b.nalloc_a = 0;
1193     m->b.a = b->a;
1194 }
1195
1196 /*!
1197  * \param[in,out] dest Destination data structure.
1198  * \param[in]     src  Source mapping.
1199  * \param[in]     bFirst If true, memory is allocated for \p dest and a full
1200  *   copy is made; otherwise, only variable parts are copied, and no memory
1201  *   is allocated.
1202  *
1203  * \p dest should have been initialized somehow (calloc() is enough).
1204  */
1205 void
1206 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1207 {
1208     if (bFirst)
1209     {
1210         gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1211         dest->type       = src->type;
1212         dest->b.nr       = src->b.nr;
1213         dest->b.nra      = src->b.nra;
1214         memcpy(dest->orgid,      src->orgid,      dest->b.nr*sizeof(*dest->orgid));
1215         memcpy(dest->b.index,    src->b.index,   (dest->b.nr+1)*sizeof(*dest->b.index));
1216         memcpy(dest->b.a,        src->b.a,        dest->b.nra*sizeof(*dest->b.a));
1217     }
1218     dest->nr         = src->nr;
1219     dest->mapb.nr    = src->mapb.nr;
1220     memcpy(dest->refid,      src->refid,      dest->nr*sizeof(*dest->refid));
1221     memcpy(dest->mapid,      src->mapid,      dest->nr*sizeof(*dest->mapid));
1222     memcpy(dest->mapb.index, src->mapb.index,(dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1223     dest->bStatic    = src->bStatic;
1224     dest->bMapStatic = src->bMapStatic;
1225 }
1226
1227 /*!
1228  * \param[in,out] m         Mapping structure.
1229  * \param[in]     g         Current index group.
1230  * \param[in]     bMaskOnly true if the unused blocks should be masked with
1231  *   -1 instead of removing them.
1232  *
1233  * Updates the index group mapping with the new index group \p g.
1234  *
1235  * \see gmx_ana_indexmap_t
1236  */
1237 void
1238 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1239                         bool bMaskOnly)
1240 {
1241     int i, j, bi, bj;
1242     bool bStatic;
1243
1244     /* Process the simple cases first */
1245     if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1246     {
1247         return;
1248     }
1249     if (m->type == INDEX_ALL)
1250     {
1251         if (m->b.nr > 0)
1252         {
1253             m->mapb.index[1] = g->isize;
1254         }
1255         return;
1256     }
1257     /* Reset the reference IDs and mapping if necessary */
1258     bStatic = (g->isize == m->b.nra && m->nr == m->b.nr);
1259     if (bStatic || bMaskOnly)
1260     {
1261         if (!m->bStatic)
1262         {
1263             for (bj = 0; bj < m->b.nr; ++bj)
1264             {
1265                 m->refid[bj] = bj;
1266             }
1267         }
1268         if (!m->bMapStatic)
1269         {
1270             for (bj = 0; bj < m->b.nr; ++bj)
1271             {
1272                 m->mapid[bj] = m->orgid[bj];
1273             }
1274             for (bj = 0; bj <= m->b.nr; ++bj)
1275             {
1276                 m->mapb.index[bj] = m->b.index[bj];
1277             }
1278             m->bMapStatic = true;
1279         }
1280     }
1281     /* Exit immediately if the group is static */
1282     if (bStatic)
1283     {
1284         m->bStatic = true;
1285         return;
1286     }
1287
1288     if (bMaskOnly)
1289     {
1290         m->nr = m->b.nr;
1291         for (i = j = bj = 0; i < g->isize; ++i, ++j)
1292         {
1293             /* Find the next atom in the block */
1294             while (m->b.a[j] != g->index[i])
1295             {
1296                 ++j;
1297             }
1298             /* Mark blocks that did not contain any atoms */
1299             while (bj < m->b.nr && m->b.index[bj+1] <= j)
1300             {
1301                 m->refid[bj++] = -1;
1302             }
1303             /* Advance the block index if we have reached the next block */
1304             if (m->b.index[bj] <= j)
1305             {
1306                 ++bj;
1307             }
1308         }
1309         /* Mark the last blocks as not accessible */
1310         while (bj < m->b.nr)
1311         {
1312             m->refid[bj++] = -1;
1313         }
1314     }
1315     else
1316     {
1317         for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1318         {
1319             /* Find the next atom in the block */
1320             while (m->b.a[j] != g->index[i])
1321             {
1322                 ++j;
1323             }
1324             /* If we have reached a new block, add it */
1325             if (m->b.index[bj+1] <= j)
1326             {
1327                 /* Skip any blocks in between */
1328                 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1329                 {
1330                     ++bj;
1331                 }
1332                 m->refid[bi] = bj;
1333                 m->mapid[bi] = m->orgid[bj];
1334                 m->mapb.index[bi] = i;
1335                 bi++;
1336             }
1337         }
1338         /* Update the number of blocks */
1339         m->mapb.index[bi] = g->isize;
1340         m->nr             = bi;
1341         m->bMapStatic     = false;
1342     }
1343     m->mapb.nr = m->nr;
1344     m->bStatic = false;
1345 }
1346
1347 /*!
1348  * \param[in,out] m         Mapping structure to free.
1349  *
1350  * All the memory allocated for the mapping structure is freed, and
1351  * the pointers set to NULL.
1352  * The pointer \p m is not freed.
1353  */
1354 void
1355 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1356 {
1357     sfree(m->refid);
1358     if (m->mapid != m->orgid)
1359     {
1360         sfree(m->mapid);
1361     }
1362     if (m->mapb.nalloc_index > 0)
1363     {
1364         sfree(m->mapb.index);
1365     }
1366     sfree(m->orgid);
1367     if (m->b.nalloc_index > 0)
1368     {
1369         sfree(m->b.index);
1370     }
1371     if (m->b.nalloc_a > 0)
1372     {
1373         sfree(m->b.a);
1374     }
1375     gmx_ana_indexmap_clear(m);
1376 }