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"
50 #include "gromacs/topology/block.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/smalloc.h"
56 /********************************************************************
57 * gmx_ana_indexgrps_t functions
58 ********************************************************************/
61 * Stores a set of index groups.
63 struct gmx_ana_indexgrps_t
65 //! Initializes an empty set of groups.
66 explicit gmx_ana_indexgrps_t(int nr) : nr(nr), g(NULL)
71 ~gmx_ana_indexgrps_t()
73 for (int i = 0; i < nr; ++i)
75 gmx_ana_index_deinit(&g[i]);
80 /** Number of index groups. */
82 /** Array of index groups. */
85 std::vector<std::string> names;
89 * \param[out] g Index group structure.
90 * \param[in] top Topology structure.
91 * \param[in] fnm File name for the index file.
92 * Memory is automatically allocated.
94 * One or both of \p top or \p fnm can be NULL.
95 * If \p top is NULL, an index file is required and the groups are read
96 * from the file (uses Gromacs routine init_index()).
97 * If \p fnm is NULL, default groups are constructed based on the
98 * topology (uses Gromacs routine analyse()).
99 * If both are null, the index group structure is initialized empty.
102 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
105 t_blocka *block = NULL;
110 block = init_index(fnm, &names);
114 block = new_blocka();
115 analyse(&top->atoms, block, &names, FALSE, FALSE);
119 *g = new gmx_ana_indexgrps_t(0);
125 *g = new gmx_ana_indexgrps_t(block->nr);
126 for (int i = 0; i < block->nr; ++i)
128 gmx_ana_index_t *grp = &(*g)->g[i];
130 grp->isize = block->index[i+1] - block->index[i];
131 snew(grp->index, grp->isize);
132 for (int j = 0; j < grp->isize; ++j)
134 grp->index[j] = block->a[block->index[i]+j];
136 grp->nalloc_index = grp->isize;
137 (*g)->names.push_back(names[i]);
142 for (int i = 0; i < block->nr; ++i)
151 for (int i = 0; i < block->nr; ++i)
161 * \param[in] g Index groups structure.
163 * The pointer \p g is invalid after the call.
166 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
172 * \param[out] g Index group structure.
173 * \returns true if \p g is empty, i.e., has 0 index groups.
176 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g)
182 * \param[in] g Index groups structure.
183 * \param[in] n Index group number to get.
184 * \returns Pointer to the \p n'th index group in \p g.
186 * The returned pointer should not be freed.
189 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n)
191 if (n < 0 || n >= g->nr)
199 * \param[out] dest Output structure.
200 * \param[out] destName Receives the name of the group if found.
201 * \param[in] src Input index groups.
202 * \param[in] n Number of the group to extract.
203 * \returns true if \p n is a valid group in \p src, false otherwise.
206 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
207 gmx_ana_indexgrps_t *src, int n)
210 if (n < 0 || n >= src->nr)
216 if (destName != NULL)
218 *destName = src->names[n];
220 gmx_ana_index_copy(dest, &src->g[n], true);
225 * \param[out] dest Output structure.
226 * \param[out] destName Receives the name of the group if found.
227 * \param[in] src Input index groups.
228 * \param[in] name Name (or part of the name) of the group to extract.
229 * \returns true if \p name is a valid group in \p src, false otherwise.
231 * Uses the Gromacs routine find_group() to find the actual group;
232 * the comparison is case-insensitive.
235 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
236 gmx_ana_indexgrps_t *src,
242 snew(names, src->nr);
243 for (int i = 0; i < src->nr; ++i)
245 names[i] = src->names[i].c_str();
247 int n = find_group(const_cast<char *>(name), src->nr,
248 const_cast<char **>(names));
256 return gmx_ana_indexgrps_extract(dest, destName, src, n);
260 * \param[in] fp Where to print the output.
261 * \param[in] g Index groups to print.
262 * \param[in] maxn Maximum number of indices to print
263 * (-1 = print all, 0 = print only names).
266 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn)
268 for (int i = 0; i < g->nr; ++i)
270 fprintf(fp, " Group %2d \"%s\" ", i, g->names[i].c_str());
271 gmx_ana_index_dump(fp, &g->g[i], maxn);
275 /********************************************************************
276 * gmx_ana_index_t functions
277 ********************************************************************/
280 * \param[in,out] g Index group structure.
281 * \param[in] isize Maximum number of atoms to reserve space for.
284 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
286 if (g->nalloc_index < isize)
288 srenew(g->index, isize);
289 g->nalloc_index = isize;
294 * \param[in,out] g Index group structure.
296 * Resizes the memory allocated for holding the indices such that the
297 * current contents fit.
300 gmx_ana_index_squeeze(gmx_ana_index_t *g)
302 srenew(g->index, g->isize);
303 g->nalloc_index = g->isize;
307 * \param[out] g Output structure.
309 * Any contents of \p g are discarded without freeing.
312 gmx_ana_index_clear(gmx_ana_index_t *g)
320 * \param[out] g Output structure.
321 * \param[in] isize Number of atoms in the new group.
322 * \param[in] index Array of \p isize atoms (can be NULL if \p isize is 0).
323 * \param[in] nalloc Number of elements allocated for \p index
324 * (if 0, \p index is not freed in gmx_ana_index_deinit())
326 * No copy if \p index is made.
329 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, int nalloc)
333 g->nalloc_index = nalloc;
337 * \param[out] g Output structure.
338 * \param[in] natoms Number of atoms.
341 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms)
346 snew(g->index, natoms);
347 for (i = 0; i < natoms; ++i)
351 g->nalloc_index = natoms;
355 * \param[in] g Index group structure.
357 * The pointer \p g is not freed.
360 gmx_ana_index_deinit(gmx_ana_index_t *g)
362 if (g->nalloc_index > 0)
366 gmx_ana_index_clear(g);
370 * \param[out] dest Destination index group.
371 * \param[in] src Source index group.
372 * \param[in] bAlloc If true, memory is allocated at \p dest; otherwise,
373 * it is assumed that enough memory has been allocated for index.
375 * A deep copy of the name is only made if \p bAlloc is true.
378 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
380 dest->isize = src->isize;
385 snew(dest->index, dest->isize);
386 dest->nalloc_index = dest->isize;
388 std::memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
393 * \param[in] fp Where to print the output.
394 * \param[in] g Index group to print.
395 * \param[in] maxn Maximum number of indices to print (-1 = print all).
398 gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int maxn)
402 fprintf(fp, "(%d atoms)", g->isize);
407 if (maxn >= 0 && n > maxn)
411 for (j = 0; j < n; ++j)
413 fprintf(fp, " %d", g->index[j]+1);
424 * \param[in] g Index group to check.
425 * \returns true if the index group is sorted and has no duplicates,
429 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
433 for (i = 0; i < g->isize-1; ++i)
435 if (g->index[i+1] <= g->index[i])
443 /********************************************************************
445 ********************************************************************/
447 /** Helper function for gmx_ana_index_sort(). */
449 cmp_atomid(const void *a, const void *b)
451 if (*(atom_id *)a < *(atom_id *)b)
455 if (*(atom_id *)a > *(atom_id *)b)
463 * \param[in,out] g Index group to be sorted.
466 gmx_ana_index_sort(gmx_ana_index_t *g)
468 std::qsort(g->index, g->isize, sizeof(*g->index), cmp_atomid);
472 * \param[in] a Index group to check.
473 * \param[in] b Index group to check.
474 * \returns true if \p a and \p b are equal, false otherwise.
477 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
481 if (a->isize != b->isize)
485 for (i = 0; i < a->isize; ++i)
487 if (a->index[i] != b->index[i])
496 * \param[in] a Index group to check against.
497 * \param[in] b Index group to check.
498 * \returns true if \p b is contained in \p a,
501 * If the elements are not in the same order in both groups, the function
502 * fails. However, the groups do not need to be sorted.
505 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
509 for (i = j = 0; j < b->isize; ++i, ++j)
511 while (i < a->isize && a->index[i] != b->index[j])
524 * \param[out] dest Output index group (the intersection of \p a and \p b).
525 * \param[in] a First index group.
526 * \param[in] b Second index group.
528 * \p dest can be the same as \p a or \p b.
531 gmx_ana_index_intersection(gmx_ana_index_t *dest,
532 gmx_ana_index_t *a, gmx_ana_index_t *b)
536 for (i = j = k = 0; i < a->isize && j < b->isize; ++i)
538 while (j < b->isize && b->index[j] < a->index[i])
542 if (j < b->isize && b->index[j] == a->index[i])
544 dest->index[k++] = b->index[j++];
551 * \param[out] dest Output index group (the difference \p a - \p b).
552 * \param[in] a First index group.
553 * \param[in] b Second index group.
555 * \p dest can equal \p a, but not \p b.
558 gmx_ana_index_difference(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; ++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++] = a->index[i];
578 * \param[in] a First index group.
579 * \param[in] b Second index group.
580 * \returns Size of the difference \p a - \p b.
583 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
587 for (i = j = k = 0; i < a->isize; ++i)
589 while (j < b->isize && b->index[j] < a->index[i])
593 if (j == b->isize || b->index[j] != a->index[i])
602 * \param[out] dest1 Output group 1 (will equal \p g).
603 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
604 * \param[in] src Group to be partitioned.
605 * \param[in] g One partition.
607 * \pre \p g is a subset of \p src and both sets are sorted
608 * \pre \p dest1 has allocated storage to store \p src
609 * \post \p dest1 == \p g
610 * \post \p dest2 == \p src - \p g
612 * No storage should be allocated for \p dest2; after the call,
613 * \p dest2->index points to the memory allocated for \p dest1
614 * (to a part that is not used by \p dest1).
616 * The calculation can be performed in-place by setting \p dest1 equal to
620 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
621 gmx_ana_index_t *src, gmx_ana_index_t *g)
625 dest2->index = dest1->index + g->isize;
626 dest2->isize = src->isize - g->isize;
627 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
629 while (j >= 0 && src->index[j] != g->index[i])
631 dest2->index[k--] = src->index[j--];
636 dest2->index[k--] = src->index[j--];
638 gmx_ana_index_copy(dest1, g, false);
642 * \param[out] dest Output index group (the union of \p a and \p b).
643 * \param[in] a First index group.
644 * \param[in] b Second index group.
646 * \p a and \p b can have common items.
647 * \p dest can equal \p a or \p b.
649 * \see gmx_ana_index_merge()
652 gmx_ana_index_union(gmx_ana_index_t *dest,
653 gmx_ana_index_t *a, gmx_ana_index_t *b)
658 dsize = gmx_ana_index_difference_size(b, a);
661 dest->isize = a->isize + dsize;
662 for (k = dest->isize - 1; k >= 0; k--)
664 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
666 dest->index[k] = b->index[j--];
670 if (j >= 0 && a->index[i] == b->index[j])
674 dest->index[k] = a->index[i--];
680 * \param[out] dest Output index group (the union of \p a and \p b).
681 * \param[in] a First index group.
682 * \param[in] b Second index group.
684 * \p a and \p b should not have common items.
685 * \p dest can equal \p a or \p b.
687 * \see gmx_ana_index_union()
690 gmx_ana_index_merge(gmx_ana_index_t *dest,
691 gmx_ana_index_t *a, gmx_ana_index_t *b)
697 dest->isize = a->isize + b->isize;
698 for (k = dest->isize - 1; k >= 0; k--)
700 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
702 dest->index[k] = b->index[j--];
706 dest->index[k] = a->index[i--];
711 /********************************************************************
712 * gmx_ana_indexmap_t and related things
713 ********************************************************************/
716 * \param[in,out] t Output block.
717 * \param[in] top Topology structure
718 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
720 * \param[in] g Index group
721 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
722 * \param[in] type Type of partitioning to make.
723 * \param[in] bComplete
724 * If true, the index group is expanded to include any residue/molecule
725 * (depending on \p type) that is partially contained in the group.
726 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
728 * \p m should have been initialized somehow (calloc() is enough) unless
729 * \p type is INDEX_UNKNOWN.
730 * \p g should be sorted.
733 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
734 e_index_t type, bool bComplete)
739 if (type == INDEX_UNKNOWN)
752 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
754 if (type != INDEX_RES && type != INDEX_MOL)
758 /* Allocate memory for the atom array and fill it unless we are using
763 /* We may allocate some extra memory here because we don't know in
764 * advance how much will be needed. */
765 if (t->nalloc_a < top->atoms.nr)
767 srenew(t->a, top->atoms.nr);
768 t->nalloc_a = top->atoms.nr;
774 if (t->nalloc_a < g->isize)
776 srenew(t->a, g->isize);
777 t->nalloc_a = g->isize;
779 std::memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
782 /* Allocate memory for the block index. We don't know in advance
783 * how much will be needed, so we allocate some extra and free it in the
785 if (t->nalloc_index < g->isize + 1)
787 srenew(t->index, g->isize + 1);
788 t->nalloc_index = g->isize + 1;
792 j = 0; /* j is used by residue completion for the first atom not stored */
794 for (i = 0; i < g->isize; ++i)
797 /* Find the ID number of the atom/residue/molecule corresponding to
805 id = top->atoms.atom[ai].resind;
808 while (ai >= top->mols.index[id+1])
813 case INDEX_UNKNOWN: /* Should not occur */
818 /* If this is the first atom in a new block, initialize the block. */
823 /* For completion, we first set the start of the block. */
824 t->index[t->nr++] = t->nra;
825 /* And then we find all the atoms that should be included. */
829 while (top->atoms.atom[j].resind != id)
833 while (j < top->atoms.nr && top->atoms.atom[j].resind == id)
841 for (j = top->mols.index[id]; j < top->mols.index[id+1]; ++j)
847 default: /* Should not be reached */
848 GMX_RELEASE_ASSERT(false, "Unreachable code was reached");
854 /* If not using completion, simply store the start of the block. */
855 t->index[t->nr++] = i;
860 /* Set the end of the last block */
861 t->index[t->nr] = t->nra;
862 /* Free any unnecessary memory */
863 srenew(t->index, t->nr+1);
864 t->nalloc_index = t->nr+1;
867 srenew(t->a, t->nra);
868 t->nalloc_a = t->nra;
873 * \param[in] g Index group to check.
874 * \param[in] b Block data to check against.
875 * \returns true if \p g consists of one or more complete blocks from \p b,
878 * The atoms in \p g are assumed to be sorted.
881 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b)
886 /* Each round in the loop matches one block */
889 /* Find the block that begins with the first unmatched atom */
890 while (bi < b->nr && b->index[bi] != g->index[i])
894 /* If not found, or if too large, return */
895 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
899 /* Check that the block matches the index */
900 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
902 if (g->index[i] != j)
907 /* Move the search to the next block */
914 * \param[in] g Index group to check.
915 * \param[in] b Block data to check against.
916 * \returns true if \p g consists of one or more complete blocks from \p b,
919 * The atoms in \p g and \p b->a are assumed to be in the same order.
922 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
927 /* Each round in the loop matches one block */
930 /* Find the block that begins with the first unmatched atom */
931 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
935 /* If not found, or if too large, return */
936 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
940 /* Check that the block matches the index */
941 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
943 if (b->a[j] != g->index[i])
948 /* Move the search to the next block */
955 * \param[in] g Index group to check.
956 * \param[in] type Block data to check against.
957 * \param[in] top Topology data.
958 * \returns true if \p g consists of one or more complete elements of type
959 * \p type, false otherwise.
961 * \p g is assumed to be sorted, otherwise may return false negatives.
963 * If \p type is \ref INDEX_ATOM, the return value is always true.
964 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
968 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
971 // TODO: Consider whether unsorted groups need to be supported better.
987 for (i = 0; i < g->isize; ++i)
990 id = top->atoms.atom[ai].resind;
993 if (ai > 0 && top->atoms.atom[ai-1].resind == id)
997 if (i > 0 && g->index[i-1] < top->atoms.nr - 1
998 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1005 if (g->index[i-1] < top->atoms.nr - 1
1006 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1014 return gmx_ana_index_has_full_blocks(g, &top->mols);
1020 * \param[out] m Output structure.
1022 * Any contents of \p m are discarded without freeing.
1025 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1027 m->type = INDEX_UNKNOWN;
1031 m->mapb.index = NULL;
1032 m->mapb.nalloc_index = 0;
1035 m->mapb.nalloc_a = 0;
1041 m->b.nalloc_index = 0;
1047 * \param[in,out] m Mapping structure.
1048 * \param[in] nr Maximum number of blocks to reserve space for.
1049 * \param[in] isize Maximum number of atoms to reserve space for.
1052 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1054 if (m->mapb.nalloc_index < nr + 1)
1056 srenew(m->refid, nr);
1057 srenew(m->mapid, nr);
1058 srenew(m->orgid, nr);
1059 srenew(m->mapb.index, nr + 1);
1060 srenew(m->b.index, nr + 1);
1061 m->mapb.nalloc_index = nr + 1;
1062 m->b.nalloc_index = nr + 1;
1064 if (m->b.nalloc_a < isize)
1066 srenew(m->b.a, isize);
1067 m->b.nalloc_a = isize;
1072 * \param[in,out] m Mapping structure to initialize.
1073 * \param[in] g Index group to map
1074 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1075 * \param[in] top Topology structure
1076 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1077 * \param[in] type Type of mapping to construct.
1079 * Initializes a new index group mapping.
1080 * The index group provided to gmx_ana_indexmap_update() should always be a
1081 * subset of the \p g given here.
1083 * \p m should have been initialized somehow (calloc() is enough).
1086 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1087 t_topology *top, e_index_t type)
1092 gmx_ana_index_make_block(&m->b, top, g, type, false);
1093 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1094 for (i = mi = 0; i < m->b.nr; ++i)
1096 ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1103 m->orgid[i] = top->atoms.atom[ii].resind;
1106 while (top->mols.index[mi+1] <= ii)
1118 for (i = 0; i < m->b.nr; ++i)
1121 m->mapid[i] = m->orgid[i];
1123 m->mapb.nr = m->b.nr;
1124 m->mapb.nra = m->b.nra;
1126 std::memcpy(m->mapb.index, m->b.index, (m->b.nr+1)*sizeof(*(m->mapb.index)));
1131 * \param[in,out] m Mapping structure to initialize.
1132 * \param[in] b Block information to use for data.
1134 * Frees some memory that is not necessary for static index group mappings.
1135 * Internal pointers are set to point to data in \p b; it is the responsibility
1136 * of the caller to ensure that the block information matches the contents of
1138 * After this function has been called, the index group provided to
1139 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1141 * This function breaks modularity of the index group mapping interface in an
1142 * ugly way, but allows reducing memory usage of static selections by a
1143 * significant amount.
1146 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1149 sfree(m->mapb.index);
1152 m->mapb.nalloc_index = 0;
1153 m->mapb.nalloc_a = 0;
1154 m->b.nalloc_index = 0;
1156 m->mapid = m->orgid;
1157 m->mapb.index = b->index;
1159 m->b.index = b->index;
1164 * \param[in,out] dest Destination data structure.
1165 * \param[in] src Source mapping.
1166 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1167 * copy is made; otherwise, only variable parts are copied, and no memory
1170 * \p dest should have been initialized somehow (calloc() is enough).
1173 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1177 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1178 dest->type = src->type;
1179 dest->b.nr = src->b.nr;
1180 dest->b.nra = src->b.nra;
1181 std::memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1182 std::memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1183 std::memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1185 dest->mapb.nr = src->mapb.nr;
1186 dest->mapb.nra = src->mapb.nra;
1187 if (src->mapb.nalloc_a > 0)
1191 snew(dest->mapb.a, src->mapb.nalloc_a);
1192 dest->mapb.nalloc_a = src->mapb.nalloc_a;
1194 std::memcpy(dest->mapb.a, src->mapb.a, dest->mapb.nra*sizeof(*dest->mapb.a));
1198 dest->mapb.a = src->mapb.a;
1200 std::memcpy(dest->refid, src->refid, dest->mapb.nr*sizeof(*dest->refid));
1201 std::memcpy(dest->mapid, src->mapid, dest->mapb.nr*sizeof(*dest->mapid));
1202 std::memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1203 dest->bStatic = src->bStatic;
1207 * Helper function to set the source atoms in an index map.
1209 * \param[in,out] m Mapping structure.
1210 * \param[in] isize Number of atoms in the \p index array.
1211 * \param[in] index List of atoms.
1214 set_atoms(gmx_ana_indexmap_t *m, int isize, int *index)
1216 m->mapb.nra = isize;
1217 if (m->mapb.nalloc_a == 0)
1223 for (int i = 0; i < isize; ++i)
1225 m->mapb.a[i] = index[i];
1231 * \param[in,out] m Mapping structure.
1232 * \param[in] g Current index group.
1233 * \param[in] bMaskOnly true if the unused blocks should be masked with
1234 * -1 instead of removing them.
1236 * Updates the index group mapping with the new index group \p g.
1238 * \see gmx_ana_indexmap_t
1241 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1246 /* Process the simple cases first */
1247 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1251 if (m->type == INDEX_ALL)
1253 set_atoms(m, g->isize, g->index);
1256 m->mapb.index[1] = g->isize;
1260 /* Reset the reference IDs and mapping if necessary */
1261 const bool bToFull = (g->isize == m->b.nra);
1262 const bool bWasFull = (m->mapb.nra == m->b.nra);
1263 if (bToFull || bMaskOnly)
1267 for (bj = 0; bj < m->b.nr; ++bj)
1274 for (bj = 0; bj < m->b.nr; ++bj)
1276 m->mapid[bj] = m->orgid[bj];
1278 for (bj = 0; bj <= m->b.nr; ++bj)
1280 m->mapb.index[bj] = m->b.index[bj];
1283 set_atoms(m, m->b.nra, m->b.a);
1284 m->mapb.nr = m->b.nr;
1286 /* Exit immediately if the group is static */
1295 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1297 /* Find the next atom in the block */
1298 while (m->b.a[j] != g->index[i])
1302 /* Mark blocks that did not contain any atoms */
1303 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1305 m->refid[bj++] = -1;
1307 /* Advance the block index if we have reached the next block */
1308 if (m->b.index[bj] <= j)
1313 /* Mark the last blocks as not accessible */
1314 while (bj < m->b.nr)
1316 m->refid[bj++] = -1;
1321 set_atoms(m, g->isize, g->index);
1322 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1324 /* Find the next atom in the block */
1325 while (m->b.a[j] != g->index[i])
1329 /* If we have reached a new block, add it */
1330 if (m->b.index[bj+1] <= j)
1332 /* Skip any blocks in between */
1333 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1338 m->mapid[bi] = m->orgid[bj];
1339 m->mapb.index[bi] = i;
1343 /* Update the number of blocks */
1344 m->mapb.index[bi] = g->isize;
1351 * \param[in,out] m Mapping structure to free.
1353 * All the memory allocated for the mapping structure is freed, and
1354 * the pointers set to NULL.
1355 * The pointer \p m is not freed.
1358 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1361 if (m->mapid != m->orgid)
1365 if (m->mapb.nalloc_index > 0)
1367 sfree(m->mapb.index);
1369 if (m->mapb.nalloc_a > 0)
1374 if (m->b.nalloc_index > 0)
1378 if (m->b.nalloc_a > 0)
1382 gmx_ana_indexmap_clear(m);