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