2 * This file is part of the GROMACS molecular simulation package.
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.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * Implements functions in indexutil.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gromacs/selection/indexutil.h"
51 #include "gromacs/topology/block.h"
52 #include "gromacs/topology/index.h"
53 #include "gromacs/topology/topology.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/smalloc.h"
57 /********************************************************************
58 * gmx_ana_indexgrps_t functions
59 ********************************************************************/
62 * Stores a set of index groups.
64 struct gmx_ana_indexgrps_t
66 //! Initializes an empty set of groups.
67 explicit gmx_ana_indexgrps_t(int nr) : nr(nr), g(NULL)
72 ~gmx_ana_indexgrps_t()
74 for (int i = 0; i < nr; ++i)
76 gmx_ana_index_deinit(&g[i]);
81 /** Number of index groups. */
83 /** Array of index groups. */
86 std::vector<std::string> names;
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.
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.
103 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
106 t_blocka *block = NULL;
111 block = init_index(fnm, &names);
115 block = new_blocka();
116 analyse(&top->atoms, block, &names, FALSE, FALSE);
120 *g = new gmx_ana_indexgrps_t(0);
126 *g = new gmx_ana_indexgrps_t(block->nr);
127 for (int i = 0; i < block->nr; ++i)
129 gmx_ana_index_t *grp = &(*g)->g[i];
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)
135 grp->index[j] = block->a[block->index[i]+j];
137 grp->nalloc_index = grp->isize;
138 (*g)->names.push_back(names[i]);
143 for (int i = 0; i < block->nr; ++i)
152 for (int i = 0; i < block->nr; ++i)
162 * \param[in] g Index groups structure.
164 * The pointer \p g is invalid after the call.
167 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
173 * \param[out] g Index group structure.
174 * \returns true if \p g is empty, i.e., has 0 index groups.
177 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g)
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.
187 * The returned pointer should not be freed.
190 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n)
192 if (n < 0 || n >= g->nr)
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.
207 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
208 gmx_ana_indexgrps_t *src, int n)
211 if (n < 0 || n >= src->nr)
217 if (destName != NULL)
219 *destName = src->names[n];
221 gmx_ana_index_copy(dest, &src->g[n], true);
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.
232 * Uses the Gromacs routine find_group() to find the actual group;
233 * the comparison is case-insensitive.
236 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
237 gmx_ana_indexgrps_t *src,
243 snew(names, src->nr);
244 for (int i = 0; i < src->nr; ++i)
246 names[i] = src->names[i].c_str();
248 int n = find_group(const_cast<char *>(name), src->nr,
249 const_cast<char **>(names));
257 return gmx_ana_indexgrps_extract(dest, destName, src, n);
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).
267 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn)
269 for (int i = 0; i < g->nr; ++i)
271 fprintf(fp, " Group %2d \"%s\" ", i, g->names[i].c_str());
272 gmx_ana_index_dump(fp, &g->g[i], maxn);
276 /********************************************************************
277 * gmx_ana_index_t functions
278 ********************************************************************/
281 * \param[in,out] g Index group structure.
282 * \param[in] isize Maximum number of atoms to reserve space for.
285 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
287 if (g->nalloc_index < isize)
289 srenew(g->index, isize);
290 g->nalloc_index = isize;
295 * \param[in,out] g Index group structure.
297 * Resizes the memory allocated for holding the indices such that the
298 * current contents fit.
301 gmx_ana_index_squeeze(gmx_ana_index_t *g)
303 srenew(g->index, g->isize);
304 g->nalloc_index = g->isize;
308 * \param[out] g Output structure.
310 * Any contents of \p g are discarded without freeing.
313 gmx_ana_index_clear(gmx_ana_index_t *g)
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())
327 * No copy if \p index is made.
330 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, int nalloc)
334 g->nalloc_index = nalloc;
338 * \param[out] g Output structure.
339 * \param[in] natoms Number of atoms.
342 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms)
347 snew(g->index, natoms);
348 for (i = 0; i < natoms; ++i)
352 g->nalloc_index = natoms;
356 * \param[in] g Index group structure.
358 * The pointer \p g is not freed.
361 gmx_ana_index_deinit(gmx_ana_index_t *g)
363 if (g->nalloc_index > 0)
367 gmx_ana_index_clear(g);
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.
376 * A deep copy of the name is only made if \p bAlloc is true.
379 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
381 dest->isize = src->isize;
386 snew(dest->index, dest->isize);
387 dest->nalloc_index = dest->isize;
389 std::memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
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).
399 gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int maxn)
403 fprintf(fp, "(%d atoms)", g->isize);
408 if (maxn >= 0 && n > maxn)
412 for (j = 0; j < n; ++j)
414 fprintf(fp, " %d", g->index[j]+1);
425 gmx_ana_index_get_max_index(gmx_ana_index_t *g)
433 return *std::max_element(g->index, g->index + g->isize);
438 * \param[in] g Index group to check.
439 * \returns true if the index group is sorted and has no duplicates,
443 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
447 for (i = 0; i < g->isize-1; ++i)
449 if (g->index[i+1] <= g->index[i])
458 gmx_ana_index_check_range(gmx_ana_index_t *g, int natoms)
460 for (int i = 0; i < g->isize; ++i)
462 if (g->index[i] < 0 || g->index[i] >= natoms)
470 /********************************************************************
472 ********************************************************************/
474 /** Helper function for gmx_ana_index_sort(). */
476 cmp_atomid(const void *a, const void *b)
478 if (*(atom_id *)a < *(atom_id *)b)
482 if (*(atom_id *)a > *(atom_id *)b)
490 * \param[in,out] g Index group to be sorted.
493 gmx_ana_index_sort(gmx_ana_index_t *g)
495 std::qsort(g->index, g->isize, sizeof(*g->index), cmp_atomid);
499 * \param[in] a Index group to check.
500 * \param[in] b Index group to check.
501 * \returns true if \p a and \p b are equal, false otherwise.
504 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
508 if (a->isize != b->isize)
512 for (i = 0; i < a->isize; ++i)
514 if (a->index[i] != b->index[i])
523 * \param[in] a Index group to check against.
524 * \param[in] b Index group to check.
525 * \returns true if \p b is contained in \p a,
528 * If the elements are not in the same order in both groups, the function
529 * fails. However, the groups do not need to be sorted.
532 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
536 for (i = j = 0; j < b->isize; ++i, ++j)
538 while (i < a->isize && a->index[i] != b->index[j])
551 * \param[out] dest Output index group (the intersection of \p a and \p b).
552 * \param[in] a First index group.
553 * \param[in] b Second index group.
555 * \p dest can be the same as \p a or \p b.
558 gmx_ana_index_intersection(gmx_ana_index_t *dest,
559 gmx_ana_index_t *a, gmx_ana_index_t *b)
563 for (i = j = k = 0; i < a->isize && j < b->isize; ++i)
565 while (j < b->isize && b->index[j] < a->index[i])
569 if (j < b->isize && b->index[j] == a->index[i])
571 dest->index[k++] = b->index[j++];
578 * \param[out] dest Output index group (the difference \p a - \p b).
579 * \param[in] a First index group.
580 * \param[in] b Second index group.
582 * \p dest can equal \p a, but not \p b.
585 gmx_ana_index_difference(gmx_ana_index_t *dest,
586 gmx_ana_index_t *a, gmx_ana_index_t *b)
590 for (i = j = k = 0; i < a->isize; ++i)
592 while (j < b->isize && b->index[j] < a->index[i])
596 if (j == b->isize || b->index[j] != a->index[i])
598 dest->index[k++] = a->index[i];
605 * \param[in] a First index group.
606 * \param[in] b Second index group.
607 * \returns Size of the difference \p a - \p b.
610 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
614 for (i = j = k = 0; i < a->isize; ++i)
616 while (j < b->isize && b->index[j] < a->index[i])
620 if (j == b->isize || b->index[j] != a->index[i])
629 * \param[out] dest1 Output group 1 (will equal \p g).
630 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
631 * \param[in] src Group to be partitioned.
632 * \param[in] g One partition.
634 * \pre \p g is a subset of \p src and both sets are sorted
635 * \pre \p dest1 has allocated storage to store \p src
636 * \post \p dest1 == \p g
637 * \post \p dest2 == \p src - \p g
639 * No storage should be allocated for \p dest2; after the call,
640 * \p dest2->index points to the memory allocated for \p dest1
641 * (to a part that is not used by \p dest1).
643 * The calculation can be performed in-place by setting \p dest1 equal to
647 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
648 gmx_ana_index_t *src, gmx_ana_index_t *g)
652 dest2->index = dest1->index + g->isize;
653 dest2->isize = src->isize - g->isize;
654 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
656 while (j >= 0 && src->index[j] != g->index[i])
658 dest2->index[k--] = src->index[j--];
663 dest2->index[k--] = src->index[j--];
665 gmx_ana_index_copy(dest1, g, false);
669 * \param[out] dest Output index group (the union of \p a and \p b).
670 * \param[in] a First index group.
671 * \param[in] b Second index group.
673 * \p a and \p b can have common items.
674 * \p dest can equal \p a or \p b.
676 * \see gmx_ana_index_merge()
679 gmx_ana_index_union(gmx_ana_index_t *dest,
680 gmx_ana_index_t *a, gmx_ana_index_t *b)
685 dsize = gmx_ana_index_difference_size(b, a);
688 dest->isize = a->isize + dsize;
689 for (k = dest->isize - 1; k >= 0; k--)
691 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
693 dest->index[k] = b->index[j--];
697 if (j >= 0 && a->index[i] == b->index[j])
701 dest->index[k] = a->index[i--];
707 * \param[out] dest Output index group (the union of \p a and \p b).
708 * \param[in] a First index group.
709 * \param[in] b Second index group.
711 * \p a and \p b should not have common items.
712 * \p dest can equal \p a or \p b.
714 * \see gmx_ana_index_union()
717 gmx_ana_index_merge(gmx_ana_index_t *dest,
718 gmx_ana_index_t *a, gmx_ana_index_t *b)
724 dest->isize = a->isize + b->isize;
725 for (k = dest->isize - 1; k >= 0; k--)
727 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
729 dest->index[k] = b->index[j--];
733 dest->index[k] = a->index[i--];
738 /********************************************************************
739 * gmx_ana_indexmap_t and related things
740 ********************************************************************/
743 * \param[in,out] t Output block.
744 * \param[in] top Topology structure
745 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
747 * \param[in] g Index group
748 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
749 * \param[in] type Type of partitioning to make.
750 * \param[in] bComplete
751 * If true, the index group is expanded to include any residue/molecule
752 * (depending on \p type) that is partially contained in the group.
753 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
755 * \p m should have been initialized somehow (calloc() is enough) unless
756 * \p type is INDEX_UNKNOWN.
757 * \p g should be sorted.
760 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
761 e_index_t type, bool bComplete)
766 if (type == INDEX_UNKNOWN)
779 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
781 if (type != INDEX_RES && type != INDEX_MOL)
785 /* Allocate memory for the atom array and fill it unless we are using
790 /* We may allocate some extra memory here because we don't know in
791 * advance how much will be needed. */
792 if (t->nalloc_a < top->atoms.nr)
794 srenew(t->a, top->atoms.nr);
795 t->nalloc_a = top->atoms.nr;
801 if (t->nalloc_a < g->isize)
803 srenew(t->a, g->isize);
804 t->nalloc_a = g->isize;
806 std::memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
809 /* Allocate memory for the block index. We don't know in advance
810 * how much will be needed, so we allocate some extra and free it in the
812 if (t->nalloc_index < g->isize + 1)
814 srenew(t->index, g->isize + 1);
815 t->nalloc_index = g->isize + 1;
819 j = 0; /* j is used by residue completion for the first atom not stored */
821 for (i = 0; i < g->isize; ++i)
824 /* Find the ID number of the atom/residue/molecule corresponding to
832 id = top->atoms.atom[ai].resind;
835 while (ai >= top->mols.index[id+1])
840 case INDEX_UNKNOWN: /* Should not occur */
845 /* If this is the first atom in a new block, initialize the block. */
850 /* For completion, we first set the start of the block. */
851 t->index[t->nr++] = t->nra;
852 /* And then we find all the atoms that should be included. */
856 while (top->atoms.atom[j].resind != id)
860 while (j < top->atoms.nr && top->atoms.atom[j].resind == id)
868 for (j = top->mols.index[id]; j < top->mols.index[id+1]; ++j)
874 default: /* Should not be reached */
875 GMX_RELEASE_ASSERT(false, "Unreachable code was reached");
881 /* If not using completion, simply store the start of the block. */
882 t->index[t->nr++] = i;
887 /* Set the end of the last block */
888 t->index[t->nr] = t->nra;
889 /* Free any unnecessary memory */
890 srenew(t->index, t->nr+1);
891 t->nalloc_index = t->nr+1;
894 srenew(t->a, t->nra);
895 t->nalloc_a = t->nra;
900 * \param[in] g Index group to check.
901 * \param[in] b Block data to check against.
902 * \returns true if \p g consists of one or more complete blocks from \p b,
905 * The atoms in \p g are assumed to be sorted.
908 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b)
913 /* Each round in the loop matches one block */
916 /* Find the block that begins with the first unmatched atom */
917 while (bi < b->nr && b->index[bi] != g->index[i])
921 /* If not found, or if too large, return */
922 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
926 /* Check that the block matches the index */
927 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
929 if (g->index[i] != j)
934 /* Move the search to the next block */
941 * \param[in] g Index group to check.
942 * \param[in] b Block data to check against.
943 * \returns true if \p g consists of one or more complete blocks from \p b,
946 * The atoms in \p g and \p b->a are assumed to be in the same order.
949 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
954 /* Each round in the loop matches one block */
957 /* Find the block that begins with the first unmatched atom */
958 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
962 /* If not found, or if too large, return */
963 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
967 /* Check that the block matches the index */
968 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
970 if (b->a[j] != g->index[i])
975 /* Move the search to the next block */
982 * \param[in] g Index group to check.
983 * \param[in] type Block data to check against.
984 * \param[in] top Topology data.
985 * \returns true if \p g consists of one or more complete elements of type
986 * \p type, false otherwise.
988 * \p g is assumed to be sorted, otherwise may return false negatives.
990 * If \p type is \ref INDEX_ATOM, the return value is always true.
991 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
995 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
998 // TODO: Consider whether unsorted groups need to be supported better.
1014 for (i = 0; i < g->isize; ++i)
1017 id = top->atoms.atom[ai].resind;
1020 if (ai > 0 && top->atoms.atom[ai-1].resind == id)
1024 if (i > 0 && g->index[i-1] < top->atoms.nr - 1
1025 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1032 if (g->index[i-1] < top->atoms.nr - 1
1033 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1041 return gmx_ana_index_has_full_blocks(g, &top->mols);
1047 * \param[out] m Output structure.
1049 * Any contents of \p m are discarded without freeing.
1052 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1054 m->type = INDEX_UNKNOWN;
1058 m->mapb.index = NULL;
1059 m->mapb.nalloc_index = 0;
1062 m->mapb.nalloc_a = 0;
1068 m->b.nalloc_index = 0;
1074 * \param[in,out] m Mapping structure.
1075 * \param[in] nr Maximum number of blocks to reserve space for.
1076 * \param[in] isize Maximum number of atoms to reserve space for.
1079 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1081 if (m->mapb.nalloc_index < nr + 1)
1083 srenew(m->refid, nr);
1084 srenew(m->mapid, nr);
1085 srenew(m->orgid, nr);
1086 srenew(m->mapb.index, nr + 1);
1087 srenew(m->b.index, nr + 1);
1088 m->mapb.nalloc_index = nr + 1;
1089 m->b.nalloc_index = nr + 1;
1091 if (m->b.nalloc_a < isize)
1093 srenew(m->b.a, isize);
1094 m->b.nalloc_a = isize;
1099 * \param[in,out] m Mapping structure to initialize.
1100 * \param[in] g Index group to map
1101 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1102 * \param[in] top Topology structure
1103 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1104 * \param[in] type Type of mapping to construct.
1106 * Initializes a new index group mapping.
1107 * The index group provided to gmx_ana_indexmap_update() should always be a
1108 * subset of the \p g given here.
1110 * \p m should have been initialized somehow (calloc() is enough).
1113 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1114 t_topology *top, e_index_t type)
1119 gmx_ana_index_make_block(&m->b, top, g, type, false);
1120 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1121 for (i = mi = 0; i < m->b.nr; ++i)
1123 ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1130 m->orgid[i] = top->atoms.atom[ii].resind;
1133 while (top->mols.index[mi+1] <= ii)
1145 for (i = 0; i < m->b.nr; ++i)
1148 m->mapid[i] = m->orgid[i];
1150 m->mapb.nr = m->b.nr;
1151 m->mapb.nra = m->b.nra;
1153 std::memcpy(m->mapb.index, m->b.index, (m->b.nr+1)*sizeof(*(m->mapb.index)));
1158 * \param[in,out] m Mapping structure to initialize.
1159 * \param[in] b Block information to use for data.
1161 * Frees some memory that is not necessary for static index group mappings.
1162 * Internal pointers are set to point to data in \p b; it is the responsibility
1163 * of the caller to ensure that the block information matches the contents of
1165 * After this function has been called, the index group provided to
1166 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1168 * This function breaks modularity of the index group mapping interface in an
1169 * ugly way, but allows reducing memory usage of static selections by a
1170 * significant amount.
1173 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1176 sfree(m->mapb.index);
1179 m->mapb.nalloc_index = 0;
1180 m->mapb.nalloc_a = 0;
1181 m->b.nalloc_index = 0;
1183 m->mapid = m->orgid;
1184 m->mapb.index = b->index;
1186 m->b.index = b->index;
1191 * \param[in,out] dest Destination data structure.
1192 * \param[in] src Source mapping.
1193 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1194 * copy is made; otherwise, only variable parts are copied, and no memory
1197 * \p dest should have been initialized somehow (calloc() is enough).
1200 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1204 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1205 dest->type = src->type;
1206 dest->b.nr = src->b.nr;
1207 dest->b.nra = src->b.nra;
1208 std::memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1209 std::memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1210 std::memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1212 dest->mapb.nr = src->mapb.nr;
1213 dest->mapb.nra = src->mapb.nra;
1214 if (src->mapb.nalloc_a > 0)
1218 snew(dest->mapb.a, src->mapb.nalloc_a);
1219 dest->mapb.nalloc_a = src->mapb.nalloc_a;
1221 std::memcpy(dest->mapb.a, src->mapb.a, dest->mapb.nra*sizeof(*dest->mapb.a));
1225 dest->mapb.a = src->mapb.a;
1227 std::memcpy(dest->refid, src->refid, dest->mapb.nr*sizeof(*dest->refid));
1228 std::memcpy(dest->mapid, src->mapid, dest->mapb.nr*sizeof(*dest->mapid));
1229 std::memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1230 dest->bStatic = src->bStatic;
1234 * Helper function to set the source atoms in an index map.
1236 * \param[in,out] m Mapping structure.
1237 * \param[in] isize Number of atoms in the \p index array.
1238 * \param[in] index List of atoms.
1241 set_atoms(gmx_ana_indexmap_t *m, int isize, int *index)
1243 m->mapb.nra = isize;
1244 if (m->mapb.nalloc_a == 0)
1250 for (int i = 0; i < isize; ++i)
1252 m->mapb.a[i] = index[i];
1258 * \param[in,out] m Mapping structure.
1259 * \param[in] g Current index group.
1260 * \param[in] bMaskOnly true if the unused blocks should be masked with
1261 * -1 instead of removing them.
1263 * Updates the index group mapping with the new index group \p g.
1265 * \see gmx_ana_indexmap_t
1268 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1273 /* Process the simple cases first */
1274 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1278 if (m->type == INDEX_ALL)
1280 set_atoms(m, g->isize, g->index);
1283 m->mapb.index[1] = g->isize;
1287 /* Reset the reference IDs and mapping if necessary */
1288 const bool bToFull = (g->isize == m->b.nra);
1289 const bool bWasFull = (m->mapb.nra == m->b.nra);
1290 if (bToFull || bMaskOnly)
1294 for (bj = 0; bj < m->b.nr; ++bj)
1301 for (bj = 0; bj < m->b.nr; ++bj)
1303 m->mapid[bj] = m->orgid[bj];
1305 for (bj = 0; bj <= m->b.nr; ++bj)
1307 m->mapb.index[bj] = m->b.index[bj];
1310 set_atoms(m, m->b.nra, m->b.a);
1311 m->mapb.nr = m->b.nr;
1313 /* Exit immediately if the group is static */
1322 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1324 /* Find the next atom in the block */
1325 while (m->b.a[j] != g->index[i])
1329 /* Mark blocks that did not contain any atoms */
1330 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1332 m->refid[bj++] = -1;
1334 /* Advance the block index if we have reached the next block */
1335 if (m->b.index[bj] <= j)
1340 /* Mark the last blocks as not accessible */
1341 while (bj < m->b.nr)
1343 m->refid[bj++] = -1;
1348 set_atoms(m, g->isize, g->index);
1349 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1351 /* Find the next atom in the block */
1352 while (m->b.a[j] != g->index[i])
1356 /* If we have reached a new block, add it */
1357 if (m->b.index[bj+1] <= j)
1359 /* Skip any blocks in between */
1360 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1365 m->mapid[bi] = m->orgid[bj];
1366 m->mapb.index[bi] = i;
1370 /* Update the number of blocks */
1371 m->mapb.index[bi] = g->isize;
1378 * \param[in,out] m Mapping structure to free.
1380 * All the memory allocated for the mapping structure is freed, and
1381 * the pointers set to NULL.
1382 * The pointer \p m is not freed.
1385 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1388 if (m->mapid != m->orgid)
1392 if (m->mapb.nalloc_index > 0)
1394 sfree(m->mapb.index);
1396 if (m->mapb.nalloc_a > 0)
1401 if (m->b.nalloc_index > 0)
1405 if (m->b.nalloc_a > 0)
1409 gmx_ana_indexmap_clear(m);