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