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
44 #include "gromacs/selection/indexutil.h"
53 #include "gromacs/topology/block.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/gmxassert.h"
57 #include "gromacs/utility/smalloc.h"
59 /********************************************************************
60 * gmx_ana_indexgrps_t functions
61 ********************************************************************/
64 * Stores a set of index groups.
66 struct gmx_ana_indexgrps_t
68 //! Initializes an empty set of groups.
69 explicit gmx_ana_indexgrps_t(int nr) : nr(nr), g(NULL)
74 ~gmx_ana_indexgrps_t()
76 for (int i = 0; i < nr; ++i)
78 gmx_ana_index_deinit(&g[i]);
83 /** Number of index groups. */
85 /** Array of index groups. */
88 std::vector<std::string> names;
92 * \param[out] g Index group structure.
93 * \param[in] top Topology structure.
94 * \param[in] fnm File name for the index file.
95 * Memory is automatically allocated.
97 * One or both of \p top or \p fnm can be NULL.
98 * If \p top is NULL, an index file is required and the groups are read
99 * from the file (uses Gromacs routine init_index()).
100 * If \p fnm is NULL, default groups are constructed based on the
101 * topology (uses Gromacs routine analyse()).
102 * If both are null, the index group structure is initialized empty.
105 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
108 t_blocka *block = NULL;
113 block = init_index(fnm, &names);
117 block = new_blocka();
118 analyse(&top->atoms, block, &names, FALSE, FALSE);
122 *g = new gmx_ana_indexgrps_t(0);
128 *g = new gmx_ana_indexgrps_t(block->nr);
129 for (int i = 0; i < block->nr; ++i)
131 gmx_ana_index_t *grp = &(*g)->g[i];
133 grp->isize = block->index[i+1] - block->index[i];
134 snew(grp->index, grp->isize);
135 for (int j = 0; j < grp->isize; ++j)
137 grp->index[j] = block->a[block->index[i]+j];
139 grp->nalloc_index = grp->isize;
140 (*g)->names.push_back(names[i]);
145 for (int i = 0; i < block->nr; ++i)
154 for (int i = 0; i < block->nr; ++i)
164 * \param[in] g Index groups structure.
166 * The pointer \p g is invalid after the call.
169 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
175 * \param[out] g Index group structure.
176 * \returns true if \p g is empty, i.e., has 0 index groups.
179 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g)
185 * \param[in] g Index groups structure.
186 * \param[in] n Index group number to get.
187 * \returns Pointer to the \p n'th index group in \p g.
189 * The returned pointer should not be freed.
192 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n)
194 if (n < 0 || n >= g->nr)
202 * \param[out] dest Output structure.
203 * \param[out] destName Receives the name of the group if found.
204 * \param[in] src Input index groups.
205 * \param[in] n Number of the group to extract.
206 * \returns true if \p n is a valid group in \p src, false otherwise.
209 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
210 gmx_ana_indexgrps_t *src, int n)
213 if (n < 0 || n >= src->nr)
219 if (destName != NULL)
221 *destName = src->names[n];
223 gmx_ana_index_copy(dest, &src->g[n], true);
228 * \param[out] dest Output structure.
229 * \param[out] destName Receives the name of the group if found.
230 * \param[in] src Input index groups.
231 * \param[in] name Name (or part of the name) of the group to extract.
232 * \returns true if \p name is a valid group in \p src, false otherwise.
234 * Uses the Gromacs routine find_group() to find the actual group;
235 * the comparison is case-insensitive.
238 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
239 gmx_ana_indexgrps_t *src,
245 snew(names, src->nr);
246 for (int i = 0; i < src->nr; ++i)
248 names[i] = src->names[i].c_str();
250 int n = find_group(const_cast<char *>(name), src->nr,
251 const_cast<char **>(names));
259 return gmx_ana_indexgrps_extract(dest, destName, src, n);
263 * \param[in] fp Where to print the output.
264 * \param[in] g Index groups to print.
265 * \param[in] maxn Maximum number of indices to print
266 * (-1 = print all, 0 = print only names).
269 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn)
271 for (int i = 0; i < g->nr; ++i)
273 fprintf(fp, " Group %2d \"%s\" ", i, g->names[i].c_str());
274 gmx_ana_index_dump(fp, &g->g[i], maxn);
278 /********************************************************************
279 * gmx_ana_index_t functions
280 ********************************************************************/
283 * \param[in,out] g Index group structure.
284 * \param[in] isize Maximum number of atoms to reserve space for.
287 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
289 if (g->nalloc_index < isize)
291 srenew(g->index, isize);
292 g->nalloc_index = isize;
297 * \param[in,out] g Index group structure.
299 * Resizes the memory allocated for holding the indices such that the
300 * current contents fit.
303 gmx_ana_index_squeeze(gmx_ana_index_t *g)
305 srenew(g->index, g->isize);
306 g->nalloc_index = g->isize;
310 * \param[out] g Output structure.
312 * Any contents of \p g are discarded without freeing.
315 gmx_ana_index_clear(gmx_ana_index_t *g)
323 * \param[out] g Output structure.
324 * \param[in] isize Number of atoms in the new group.
325 * \param[in] index Array of \p isize atoms (can be NULL if \p isize is 0).
326 * \param[in] nalloc Number of elements allocated for \p index
327 * (if 0, \p index is not freed in gmx_ana_index_deinit())
329 * No copy if \p index is made.
332 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, int nalloc)
336 g->nalloc_index = nalloc;
340 * \param[out] g Output structure.
341 * \param[in] natoms Number of atoms.
344 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms)
349 snew(g->index, natoms);
350 for (i = 0; i < natoms; ++i)
354 g->nalloc_index = natoms;
358 * \param[in] g Index group structure.
360 * The pointer \p g is not freed.
363 gmx_ana_index_deinit(gmx_ana_index_t *g)
365 if (g->nalloc_index > 0)
369 gmx_ana_index_clear(g);
373 * \param[out] dest Destination index group.
374 * \param[in] src Source index group.
375 * \param[in] bAlloc If true, memory is allocated at \p dest; otherwise,
376 * it is assumed that enough memory has been allocated for index.
378 * A deep copy of the name is only made if \p bAlloc is true.
381 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
383 dest->isize = src->isize;
388 snew(dest->index, dest->isize);
389 dest->nalloc_index = dest->isize;
391 std::memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
396 * \param[in] fp Where to print the output.
397 * \param[in] g Index group to print.
398 * \param[in] maxn Maximum number of indices to print (-1 = print all).
401 gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int maxn)
405 fprintf(fp, "(%d atoms)", g->isize);
410 if (maxn >= 0 && n > maxn)
414 for (j = 0; j < n; ++j)
416 fprintf(fp, " %d", g->index[j]+1);
427 gmx_ana_index_get_max_index(gmx_ana_index_t *g)
435 return *std::max_element(g->index, g->index + g->isize);
440 * \param[in] g Index group to check.
441 * \returns true if the index group is sorted and has no duplicates,
445 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
449 for (i = 0; i < g->isize-1; ++i)
451 if (g->index[i+1] <= g->index[i])
460 gmx_ana_index_check_range(gmx_ana_index_t *g, int natoms)
462 for (int i = 0; i < g->isize; ++i)
464 if (g->index[i] < 0 || g->index[i] >= natoms)
472 /********************************************************************
474 ********************************************************************/
476 /** Helper function for gmx_ana_index_sort(). */
478 cmp_atomid(const void *a, const void *b)
480 if (*(atom_id *)a < *(atom_id *)b)
484 if (*(atom_id *)a > *(atom_id *)b)
492 * \param[in,out] g Index group to be sorted.
495 gmx_ana_index_sort(gmx_ana_index_t *g)
497 std::qsort(g->index, g->isize, sizeof(*g->index), cmp_atomid);
501 * \param[in] a Index group to check.
502 * \param[in] b Index group to check.
503 * \returns true if \p a and \p b are equal, false otherwise.
506 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
510 if (a->isize != b->isize)
514 for (i = 0; i < a->isize; ++i)
516 if (a->index[i] != b->index[i])
525 * \param[in] a Index group to check against.
526 * \param[in] b Index group to check.
527 * \returns true if \p b is contained in \p a,
530 * If the elements are not in the same order in both groups, the function
531 * fails. However, the groups do not need to be sorted.
534 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
538 for (i = j = 0; j < b->isize; ++i, ++j)
540 while (i < a->isize && a->index[i] != b->index[j])
553 * \param[out] dest Output index group (the intersection of \p a and \p b).
554 * \param[in] a First index group.
555 * \param[in] b Second index group.
557 * \p dest can be the same as \p a or \p b.
560 gmx_ana_index_intersection(gmx_ana_index_t *dest,
561 gmx_ana_index_t *a, gmx_ana_index_t *b)
565 for (i = j = k = 0; i < a->isize && j < b->isize; ++i)
567 while (j < b->isize && b->index[j] < a->index[i])
571 if (j < b->isize && b->index[j] == a->index[i])
573 dest->index[k++] = b->index[j++];
580 * \param[out] dest Output index group (the difference \p a - \p b).
581 * \param[in] a First index group.
582 * \param[in] b Second index group.
584 * \p dest can equal \p a, but not \p b.
587 gmx_ana_index_difference(gmx_ana_index_t *dest,
588 gmx_ana_index_t *a, gmx_ana_index_t *b)
592 for (i = j = k = 0; i < a->isize; ++i)
594 while (j < b->isize && b->index[j] < a->index[i])
598 if (j == b->isize || b->index[j] != a->index[i])
600 dest->index[k++] = a->index[i];
607 * \param[in] a First index group.
608 * \param[in] b Second index group.
609 * \returns Size of the difference \p a - \p b.
612 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
616 for (i = j = k = 0; i < a->isize; ++i)
618 while (j < b->isize && b->index[j] < a->index[i])
622 if (j == b->isize || b->index[j] != a->index[i])
631 * \param[out] dest1 Output group 1 (will equal \p g).
632 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
633 * \param[in] src Group to be partitioned.
634 * \param[in] g One partition.
636 * \pre \p g is a subset of \p src and both sets are sorted
637 * \pre \p dest1 has allocated storage to store \p src
638 * \post \p dest1 == \p g
639 * \post \p dest2 == \p src - \p g
641 * No storage should be allocated for \p dest2; after the call,
642 * \p dest2->index points to the memory allocated for \p dest1
643 * (to a part that is not used by \p dest1).
645 * The calculation can be performed in-place by setting \p dest1 equal to
649 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
650 gmx_ana_index_t *src, gmx_ana_index_t *g)
654 dest2->index = dest1->index + g->isize;
655 dest2->isize = src->isize - g->isize;
656 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
658 while (j >= 0 && src->index[j] != g->index[i])
660 dest2->index[k--] = src->index[j--];
665 dest2->index[k--] = src->index[j--];
667 gmx_ana_index_copy(dest1, g, false);
671 * \param[out] dest Output index group (the union of \p a and \p b).
672 * \param[in] a First index group.
673 * \param[in] b Second index group.
675 * \p a and \p b can have common items.
676 * \p dest can equal \p a or \p b.
678 * \see gmx_ana_index_merge()
681 gmx_ana_index_union(gmx_ana_index_t *dest,
682 gmx_ana_index_t *a, gmx_ana_index_t *b)
687 dsize = gmx_ana_index_difference_size(b, a);
690 dest->isize = a->isize + dsize;
691 for (k = dest->isize - 1; k >= 0; k--)
693 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
695 dest->index[k] = b->index[j--];
699 if (j >= 0 && a->index[i] == b->index[j])
703 dest->index[k] = a->index[i--];
709 * \param[out] dest Output index group (the union of \p a and \p b).
710 * \param[in] a First index group.
711 * \param[in] b Second index group.
713 * \p a and \p b should not have common items.
714 * \p dest can equal \p a or \p b.
716 * \see gmx_ana_index_union()
719 gmx_ana_index_merge(gmx_ana_index_t *dest,
720 gmx_ana_index_t *a, gmx_ana_index_t *b)
726 dest->isize = a->isize + b->isize;
727 for (k = dest->isize - 1; k >= 0; k--)
729 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
731 dest->index[k] = b->index[j--];
735 dest->index[k] = a->index[i--];
740 /********************************************************************
741 * gmx_ana_indexmap_t and related things
742 ********************************************************************/
745 * \param[in,out] t Output block.
746 * \param[in] top Topology structure
747 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
749 * \param[in] g Index group
750 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
751 * \param[in] type Type of partitioning to make.
752 * \param[in] bComplete
753 * If true, the index group is expanded to include any residue/molecule
754 * (depending on \p type) that is partially contained in the group.
755 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
757 * \p m should have been initialized somehow (calloc() is enough) unless
758 * \p type is INDEX_UNKNOWN.
759 * \p g should be sorted.
762 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
763 e_index_t type, bool bComplete)
768 if (type == INDEX_UNKNOWN)
781 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
783 if (type != INDEX_RES && type != INDEX_MOL)
787 /* Allocate memory for the atom array and fill it unless we are using
792 /* We may allocate some extra memory here because we don't know in
793 * advance how much will be needed. */
794 if (t->nalloc_a < top->atoms.nr)
796 srenew(t->a, top->atoms.nr);
797 t->nalloc_a = top->atoms.nr;
803 if (t->nalloc_a < g->isize)
805 srenew(t->a, g->isize);
806 t->nalloc_a = g->isize;
808 std::memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
811 /* Allocate memory for the block index. We don't know in advance
812 * how much will be needed, so we allocate some extra and free it in the
814 if (t->nalloc_index < g->isize + 1)
816 srenew(t->index, g->isize + 1);
817 t->nalloc_index = g->isize + 1;
821 j = 0; /* j is used by residue completion for the first atom not stored */
823 for (i = 0; i < g->isize; ++i)
826 /* Find the ID number of the atom/residue/molecule corresponding to
834 id = top->atoms.atom[ai].resind;
837 while (ai >= top->mols.index[id+1])
842 case INDEX_UNKNOWN: /* Should not occur */
847 /* If this is the first atom in a new block, initialize the block. */
852 /* For completion, we first set the start of the block. */
853 t->index[t->nr++] = t->nra;
854 /* And then we find all the atoms that should be included. */
858 while (top->atoms.atom[j].resind != id)
862 while (j < top->atoms.nr && top->atoms.atom[j].resind == id)
870 for (j = top->mols.index[id]; j < top->mols.index[id+1]; ++j)
876 default: /* Should not be reached */
877 GMX_RELEASE_ASSERT(false, "Unreachable code was reached");
883 /* If not using completion, simply store the start of the block. */
884 t->index[t->nr++] = i;
889 /* Set the end of the last block */
890 t->index[t->nr] = t->nra;
891 /* Free any unnecessary memory */
892 srenew(t->index, t->nr+1);
893 t->nalloc_index = t->nr+1;
896 srenew(t->a, t->nra);
897 t->nalloc_a = t->nra;
902 * \param[in] g Index group to check.
903 * \param[in] b Block data to check against.
904 * \returns true if \p g consists of one or more complete blocks from \p b,
907 * The atoms in \p g are assumed to be sorted.
910 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b)
915 /* Each round in the loop matches one block */
918 /* Find the block that begins with the first unmatched atom */
919 while (bi < b->nr && b->index[bi] != g->index[i])
923 /* If not found, or if too large, return */
924 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
928 /* Check that the block matches the index */
929 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
931 if (g->index[i] != j)
936 /* Move the search to the next block */
943 * \param[in] g Index group to check.
944 * \param[in] b Block data to check against.
945 * \returns true if \p g consists of one or more complete blocks from \p b,
948 * The atoms in \p g and \p b->a are assumed to be in the same order.
951 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
956 /* Each round in the loop matches one block */
959 /* Find the block that begins with the first unmatched atom */
960 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
964 /* If not found, or if too large, return */
965 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
969 /* Check that the block matches the index */
970 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
972 if (b->a[j] != g->index[i])
977 /* Move the search to the next block */
984 * \param[in] g Index group to check.
985 * \param[in] type Block data to check against.
986 * \param[in] top Topology data.
987 * \returns true if \p g consists of one or more complete elements of type
988 * \p type, false otherwise.
990 * \p g is assumed to be sorted, otherwise may return false negatives.
992 * If \p type is \ref INDEX_ATOM, the return value is always true.
993 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
997 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
1000 // TODO: Consider whether unsorted groups need to be supported better.
1016 for (i = 0; i < g->isize; ++i)
1019 id = top->atoms.atom[ai].resind;
1022 if (ai > 0 && top->atoms.atom[ai-1].resind == id)
1026 if (i > 0 && g->index[i-1] < top->atoms.nr - 1
1027 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1034 if (g->index[i-1] < top->atoms.nr - 1
1035 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1043 return gmx_ana_index_has_full_blocks(g, &top->mols);
1049 * \param[out] m Output structure.
1051 * Any contents of \p m are discarded without freeing.
1054 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1056 m->type = INDEX_UNKNOWN;
1060 m->mapb.index = NULL;
1061 m->mapb.nalloc_index = 0;
1064 m->mapb.nalloc_a = 0;
1070 m->b.nalloc_index = 0;
1076 * \param[in,out] m Mapping structure.
1077 * \param[in] nr Maximum number of blocks to reserve space for.
1078 * \param[in] isize Maximum number of atoms to reserve space for.
1081 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1083 if (m->mapb.nalloc_index < nr + 1)
1085 srenew(m->refid, nr);
1086 srenew(m->mapid, nr);
1087 srenew(m->orgid, nr);
1088 srenew(m->mapb.index, nr + 1);
1089 srenew(m->b.index, nr + 1);
1090 m->mapb.nalloc_index = nr + 1;
1091 m->b.nalloc_index = nr + 1;
1093 if (m->b.nalloc_a < isize)
1095 srenew(m->b.a, isize);
1096 m->b.nalloc_a = isize;
1101 * \param[in,out] m Mapping structure to initialize.
1102 * \param[in] g Index group to map
1103 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1104 * \param[in] top Topology structure
1105 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1106 * \param[in] type Type of mapping to construct.
1108 * Initializes a new index group mapping.
1109 * The index group provided to gmx_ana_indexmap_update() should always be a
1110 * subset of the \p g given here.
1112 * \p m should have been initialized somehow (calloc() is enough).
1115 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1116 t_topology *top, e_index_t type)
1121 gmx_ana_index_make_block(&m->b, top, g, type, false);
1122 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1123 for (i = mi = 0; i < m->b.nr; ++i)
1125 ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1132 m->orgid[i] = top->atoms.atom[ii].resind;
1135 while (top->mols.index[mi+1] <= ii)
1147 for (i = 0; i < m->b.nr; ++i)
1150 m->mapid[i] = m->orgid[i];
1152 m->mapb.nr = m->b.nr;
1153 m->mapb.nra = m->b.nra;
1155 std::memcpy(m->mapb.index, m->b.index, (m->b.nr+1)*sizeof(*(m->mapb.index)));
1160 * \param[in,out] m Mapping structure to initialize.
1161 * \param[in] b Block information to use for data.
1163 * Frees some memory that is not necessary for static index group mappings.
1164 * Internal pointers are set to point to data in \p b; it is the responsibility
1165 * of the caller to ensure that the block information matches the contents of
1167 * After this function has been called, the index group provided to
1168 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1170 * This function breaks modularity of the index group mapping interface in an
1171 * ugly way, but allows reducing memory usage of static selections by a
1172 * significant amount.
1175 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1178 sfree(m->mapb.index);
1181 m->mapb.nalloc_index = 0;
1182 m->mapb.nalloc_a = 0;
1183 m->b.nalloc_index = 0;
1185 m->mapid = m->orgid;
1186 m->mapb.index = b->index;
1188 m->b.index = b->index;
1193 * \param[in,out] dest Destination data structure.
1194 * \param[in] src Source mapping.
1195 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1196 * copy is made; otherwise, only variable parts are copied, and no memory
1199 * \p dest should have been initialized somehow (calloc() is enough).
1202 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1206 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1207 dest->type = src->type;
1208 dest->b.nr = src->b.nr;
1209 dest->b.nra = src->b.nra;
1210 std::memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1211 std::memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1212 std::memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1214 dest->mapb.nr = src->mapb.nr;
1215 dest->mapb.nra = src->mapb.nra;
1216 if (src->mapb.nalloc_a > 0)
1220 snew(dest->mapb.a, src->mapb.nalloc_a);
1221 dest->mapb.nalloc_a = src->mapb.nalloc_a;
1223 std::memcpy(dest->mapb.a, src->mapb.a, dest->mapb.nra*sizeof(*dest->mapb.a));
1227 dest->mapb.a = src->mapb.a;
1229 std::memcpy(dest->refid, src->refid, dest->mapb.nr*sizeof(*dest->refid));
1230 std::memcpy(dest->mapid, src->mapid, dest->mapb.nr*sizeof(*dest->mapid));
1231 std::memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1232 dest->bStatic = src->bStatic;
1236 * Helper function to set the source atoms in an index map.
1238 * \param[in,out] m Mapping structure.
1239 * \param[in] isize Number of atoms in the \p index array.
1240 * \param[in] index List of atoms.
1243 set_atoms(gmx_ana_indexmap_t *m, int isize, int *index)
1245 m->mapb.nra = isize;
1246 if (m->mapb.nalloc_a == 0)
1252 for (int i = 0; i < isize; ++i)
1254 m->mapb.a[i] = index[i];
1260 * \param[in,out] m Mapping structure.
1261 * \param[in] g Current index group.
1262 * \param[in] bMaskOnly true if the unused blocks should be masked with
1263 * -1 instead of removing them.
1265 * Updates the index group mapping with the new index group \p g.
1267 * \see gmx_ana_indexmap_t
1270 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1275 /* Process the simple cases first */
1276 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1280 if (m->type == INDEX_ALL)
1282 set_atoms(m, g->isize, g->index);
1285 m->mapb.index[1] = g->isize;
1289 /* Reset the reference IDs and mapping if necessary */
1290 const bool bToFull = (g->isize == m->b.nra);
1291 const bool bWasFull = (m->mapb.nra == m->b.nra);
1292 if (bToFull || bMaskOnly)
1296 for (bj = 0; bj < m->b.nr; ++bj)
1303 for (bj = 0; bj < m->b.nr; ++bj)
1305 m->mapid[bj] = m->orgid[bj];
1307 for (bj = 0; bj <= m->b.nr; ++bj)
1309 m->mapb.index[bj] = m->b.index[bj];
1312 set_atoms(m, m->b.nra, m->b.a);
1313 m->mapb.nr = m->b.nr;
1315 /* Exit immediately if the group is static */
1324 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1326 /* Find the next atom in the block */
1327 while (m->b.a[j] != g->index[i])
1331 /* Mark blocks that did not contain any atoms */
1332 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1334 m->refid[bj++] = -1;
1336 /* Advance the block index if we have reached the next block */
1337 if (m->b.index[bj] <= j)
1342 /* Mark the last blocks as not accessible */
1343 while (bj < m->b.nr)
1345 m->refid[bj++] = -1;
1350 set_atoms(m, g->isize, g->index);
1351 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1353 /* Find the next atom in the block */
1354 while (m->b.a[j] != g->index[i])
1358 /* If we have reached a new block, add it */
1359 if (m->b.index[bj+1] <= j)
1361 /* Skip any blocks in between */
1362 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1367 m->mapid[bi] = m->orgid[bj];
1368 m->mapb.index[bi] = i;
1372 /* Update the number of blocks */
1373 m->mapb.index[bi] = g->isize;
1380 * \param[in,out] m Mapping structure to free.
1382 * All the memory allocated for the mapping structure is freed, and
1383 * the pointers set to NULL.
1384 * The pointer \p m is not freed.
1387 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1390 if (m->mapid != m->orgid)
1394 if (m->mapb.nalloc_index > 0)
1396 sfree(m->mapb.index);
1398 if (m->mapb.nalloc_a > 0)
1403 if (m->b.nalloc_index > 0)
1407 if (m->b.nalloc_a > 0)
1411 gmx_ana_indexmap_clear(m);