2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013, 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"
49 #include "gromacs/legacyheaders/index.h"
50 #include "gromacs/legacyheaders/gmx_fatal.h"
51 #include "gromacs/legacyheaders/smalloc.h"
52 #include "gromacs/legacyheaders/typedefs.h"
55 /********************************************************************
56 * gmx_ana_indexgrps_t functions
57 ********************************************************************/
60 * Stores a set of index groups.
62 struct gmx_ana_indexgrps_t
64 //! Initializes an empty set of groups.
65 explicit gmx_ana_indexgrps_t(int nr) : nr(nr), g(NULL)
70 ~gmx_ana_indexgrps_t()
72 for (int i = 0; i < nr; ++i)
74 gmx_ana_index_deinit(&g[i]);
79 /** Number of index groups. */
81 /** Array of index groups. */
84 std::vector<std::string> names;
88 * \param[out] g Index group structure.
89 * \param[in] top Topology structure.
90 * \param[in] fnm File name for the index file.
91 * Memory is automatically allocated.
93 * One or both of \p top or \p fnm can be NULL.
94 * If \p top is NULL, an index file is required and the groups are read
95 * from the file (uses Gromacs routine init_index()).
96 * If \p fnm is NULL, default groups are constructed based on the
97 * topology (uses Gromacs routine analyse()).
98 * If both are null, the index group structure is initialized empty.
101 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
104 t_blocka *block = NULL;
109 block = init_index(fnm, &names);
113 block = new_blocka();
114 analyse(&top->atoms, block, &names, FALSE, FALSE);
118 *g = new gmx_ana_indexgrps_t(0);
124 *g = new gmx_ana_indexgrps_t(block->nr);
125 for (int i = 0; i < block->nr; ++i)
127 gmx_ana_index_t *grp = &(*g)->g[i];
129 grp->isize = block->index[i+1] - block->index[i];
130 snew(grp->index, grp->isize);
131 for (int j = 0; j < grp->isize; ++j)
133 grp->index[j] = block->a[block->index[i]+j];
135 grp->nalloc_index = grp->isize;
136 (*g)->names.push_back(names[i]);
141 for (int i = 0; i < block->nr; ++i)
150 for (int i = 0; i < block->nr; ++i)
160 * \param[in] g Index groups structure.
162 * The pointer \p g is invalid after the call.
165 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
171 * \param[out] g Index group structure.
172 * \returns true if \p g is empty, i.e., has 0 index groups.
175 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g)
181 * \param[in] g Index groups structure.
182 * \param[in] n Index group number to get.
183 * \returns Pointer to the \p n'th index group in \p g.
185 * The returned pointer should not be freed.
188 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n)
190 if (n < 0 || n >= g->nr)
198 * \param[out] dest Output structure.
199 * \param[out] destName Receives the name of the group if found.
200 * \param[in] src Input index groups.
201 * \param[in] n Number of the group to extract.
202 * \returns true if \p n is a valid group in \p src, false otherwise.
205 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
206 gmx_ana_indexgrps_t *src, int n)
209 if (n < 0 || n >= src->nr)
215 if (destName != NULL)
217 *destName = src->names[n];
219 gmx_ana_index_copy(dest, &src->g[n], true);
224 * \param[out] dest Output structure.
225 * \param[out] destName Receives the name of the group if found.
226 * \param[in] src Input index groups.
227 * \param[in] name Name (or part of the name) of the group to extract.
228 * \returns true if \p name is a valid group in \p src, false otherwise.
230 * Uses the Gromacs routine find_group() to find the actual group;
231 * the comparison is case-insensitive.
234 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
235 gmx_ana_indexgrps_t *src,
241 snew(names, src->nr);
242 for (int i = 0; i < src->nr; ++i)
244 names[i] = src->names[i].c_str();
246 int n = find_group(const_cast<char *>(name), src->nr,
247 const_cast<char **>(names));
255 return gmx_ana_indexgrps_extract(dest, destName, src, n);
259 * \param[in] fp Where to print the output.
260 * \param[in] g Index groups to print.
261 * \param[in] maxn Maximum number of indices to print
262 * (-1 = print all, 0 = print only names).
265 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn)
267 for (int i = 0; i < g->nr; ++i)
269 fprintf(fp, " Group %2d \"%s\" ", i + 1, g->names[i].c_str());
270 gmx_ana_index_dump(fp, &g->g[i], maxn);
274 /********************************************************************
275 * gmx_ana_index_t functions
276 ********************************************************************/
279 * \param[in,out] g Index group structure.
280 * \param[in] isize Maximum number of atoms to reserve space for.
283 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
285 if (g->nalloc_index < isize)
287 srenew(g->index, isize);
288 g->nalloc_index = isize;
293 * \param[in,out] g Index group structure.
295 * Resizes the memory allocated for holding the indices such that the
296 * current contents fit.
299 gmx_ana_index_squeeze(gmx_ana_index_t *g)
301 srenew(g->index, g->isize);
302 g->nalloc_index = g->isize;
306 * \param[out] g Output structure.
308 * Any contents of \p g are discarded without freeing.
311 gmx_ana_index_clear(gmx_ana_index_t *g)
319 * \param[out] g Output structure.
320 * \param[in] isize Number of atoms in the new group.
321 * \param[in] index Array of \p isize atoms (can be NULL if \p isize is 0).
322 * \param[in] nalloc Number of elements allocated for \p index
323 * (if 0, \p index is not freed in gmx_ana_index_deinit())
325 * No copy if \p index is made.
328 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, int nalloc)
332 g->nalloc_index = nalloc;
336 * \param[out] g Output structure.
337 * \param[in] natoms Number of atoms.
340 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms)
345 snew(g->index, natoms);
346 for (i = 0; i < natoms; ++i)
350 g->nalloc_index = natoms;
354 * \param[in] g Index group structure.
356 * The pointer \p g is not freed.
359 gmx_ana_index_deinit(gmx_ana_index_t *g)
361 if (g->nalloc_index > 0)
365 gmx_ana_index_clear(g);
369 * \param[out] dest Destination index group.
370 * \param[in] src Source index group.
371 * \param[in] bAlloc If true, memory is allocated at \p dest; otherwise,
372 * it is assumed that enough memory has been allocated for index.
374 * A deep copy of the name is only made if \p bAlloc is true.
377 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
379 dest->isize = src->isize;
384 snew(dest->index, dest->isize);
385 dest->nalloc_index = dest->isize;
387 std::memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
392 * \param[in] fp Where to print the output.
393 * \param[in] g Index group to print.
394 * \param[in] maxn Maximum number of indices to print (-1 = print all).
397 gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int maxn)
401 fprintf(fp, "(%d atoms)", g->isize);
406 if (maxn >= 0 && n > maxn)
410 for (j = 0; j < n; ++j)
412 fprintf(fp, " %d", g->index[j]+1);
423 * \param[in] g Index group to check.
424 * \returns true if the index group is sorted and has no duplicates,
428 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
432 for (i = 0; i < g->isize-1; ++i)
434 if (g->index[i+1] <= g->index[i])
442 /********************************************************************
444 ********************************************************************/
446 /** Helper function for gmx_ana_index_sort(). */
448 cmp_atomid(const void *a, const void *b)
450 if (*(atom_id *)a < *(atom_id *)b)
454 if (*(atom_id *)a > *(atom_id *)b)
462 * \param[in,out] g Index group to be sorted.
465 gmx_ana_index_sort(gmx_ana_index_t *g)
467 qsort(g->index, g->isize, sizeof(*g->index), cmp_atomid);
471 * \param[in] a Index group to check.
472 * \param[in] b Index group to check.
473 * \returns true if \p a and \p b are equal, false otherwise.
476 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
480 if (a->isize != b->isize)
484 for (i = 0; i < a->isize; ++i)
486 if (a->index[i] != b->index[i])
495 * \param[in] a Index group to check against.
496 * \param[in] b Index group to check.
497 * \returns true if \p b is contained in \p a,
500 * If the elements are not in the same order in both groups, the function
501 * fails. However, the groups do not need to be sorted.
504 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
508 for (i = j = 0; j < b->isize; ++i, ++j)
510 while (i < a->isize && a->index[i] != b->index[j])
523 * \param[out] dest Output index group (the intersection of \p a and \p b).
524 * \param[in] a First index group.
525 * \param[in] b Second index group.
527 * \p dest can be the same as \p a or \p b.
530 gmx_ana_index_intersection(gmx_ana_index_t *dest,
531 gmx_ana_index_t *a, gmx_ana_index_t *b)
535 for (i = j = k = 0; i < a->isize && j < b->isize; ++i)
537 while (j < b->isize && b->index[j] < a->index[i])
541 if (j < b->isize && b->index[j] == a->index[i])
543 dest->index[k++] = b->index[j++];
550 * \param[out] dest Output index group (the difference \p a - \p b).
551 * \param[in] a First index group.
552 * \param[in] b Second index group.
554 * \p dest can equal \p a, but not \p b.
557 gmx_ana_index_difference(gmx_ana_index_t *dest,
558 gmx_ana_index_t *a, gmx_ana_index_t *b)
562 for (i = j = k = 0; i < a->isize; ++i)
564 while (j < b->isize && b->index[j] < a->index[i])
568 if (j == b->isize || b->index[j] != a->index[i])
570 dest->index[k++] = a->index[i];
577 * \param[in] a First index group.
578 * \param[in] b Second index group.
579 * \returns Size of the difference \p a - \p b.
582 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
586 for (i = j = k = 0; i < a->isize; ++i)
588 while (j < b->isize && b->index[j] < a->index[i])
592 if (j == b->isize || b->index[j] != a->index[i])
601 * \param[out] dest1 Output group 1 (will equal \p g).
602 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
603 * \param[in] src Group to be partitioned.
604 * \param[in] g One partition.
606 * \pre \p g is a subset of \p src and both sets are sorted
607 * \pre \p dest1 has allocated storage to store \p src
608 * \post \p dest1 == \p g
609 * \post \p dest2 == \p src - \p g
611 * No storage should be allocated for \p dest2; after the call,
612 * \p dest2->index points to the memory allocated for \p dest1
613 * (to a part that is not used by \p dest1).
615 * The calculation can be performed in-place by setting \p dest1 equal to
619 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
620 gmx_ana_index_t *src, gmx_ana_index_t *g)
624 dest2->index = dest1->index + g->isize;
625 dest2->isize = src->isize - g->isize;
626 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
628 while (j >= 0 && src->index[j] != g->index[i])
630 dest2->index[k--] = src->index[j--];
635 dest2->index[k--] = src->index[j--];
637 gmx_ana_index_copy(dest1, g, false);
641 * \param[out] dest Output index group (the union of \p a and \p b).
642 * \param[in] a First index group.
643 * \param[in] b Second index group.
645 * \p a and \p b can have common items.
646 * \p dest can equal \p a or \p b.
648 * \see gmx_ana_index_merge()
651 gmx_ana_index_union(gmx_ana_index_t *dest,
652 gmx_ana_index_t *a, gmx_ana_index_t *b)
657 dsize = gmx_ana_index_difference_size(b, a);
660 dest->isize = a->isize + dsize;
661 for (k = dest->isize - 1; k >= 0; k--)
663 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
665 dest->index[k] = b->index[j--];
669 if (j >= 0 && a->index[i] == b->index[j])
673 dest->index[k] = a->index[i--];
679 * \param[out] dest Output index group (the union of \p a and \p b).
680 * \param[in] a First index group.
681 * \param[in] b Second index group.
683 * \p a and \p b should not have common items.
684 * \p dest can equal \p a or \p b.
686 * \see gmx_ana_index_union()
689 gmx_ana_index_merge(gmx_ana_index_t *dest,
690 gmx_ana_index_t *a, gmx_ana_index_t *b)
696 dest->isize = a->isize + b->isize;
697 for (k = dest->isize - 1; k >= 0; k--)
699 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
701 dest->index[k] = b->index[j--];
705 dest->index[k] = a->index[i--];
710 /********************************************************************
711 * gmx_ana_indexmap_t and related things
712 ********************************************************************/
715 * \param[in,out] t Output block.
716 * \param[in] top Topology structure
717 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
719 * \param[in] g Index group
720 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
721 * \param[in] type Type of partitioning to make.
722 * \param[in] bComplete
723 * If true, the index group is expanded to include any residue/molecule
724 * (depending on \p type) that is partially contained in the group.
725 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
727 * \p m should have been initialized somehow (calloc() is enough) unless
728 * \p type is INDEX_UNKNOWN.
729 * \p g should be sorted.
732 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
733 e_index_t type, bool bComplete)
738 if (type == INDEX_UNKNOWN)
751 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
753 if (type != INDEX_RES && type != INDEX_MOL)
757 /* Allocate memory for the atom array and fill it unless we are using
762 /* We may allocate some extra memory here because we don't know in
763 * advance how much will be needed. */
764 if (t->nalloc_a < top->atoms.nr)
766 srenew(t->a, top->atoms.nr);
767 t->nalloc_a = top->atoms.nr;
773 if (t->nalloc_a < g->isize)
775 srenew(t->a, g->isize);
776 t->nalloc_a = g->isize;
778 std::memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
781 /* Allocate memory for the block index. We don't know in advance
782 * how much will be needed, so we allocate some extra and free it in the
784 if (t->nalloc_index < g->isize + 1)
786 srenew(t->index, g->isize + 1);
787 t->nalloc_index = g->isize + 1;
791 j = 0; /* j is used by residue completion for the first atom not stored */
793 for (i = 0; i < g->isize; ++i)
796 /* Find the ID number of the atom/residue/molecule corresponding to
804 id = top->atoms.atom[ai].resind;
807 while (ai >= top->mols.index[id+1])
812 case INDEX_UNKNOWN: /* Should not occur */
817 /* If this is the first atom in a new block, initialize the block. */
822 /* For completion, we first set the start of the block. */
823 t->index[t->nr++] = t->nra;
824 /* And then we find all the atoms that should be included. */
828 while (top->atoms.atom[j].resind != id)
832 while (j < top->atoms.nr && top->atoms.atom[j].resind == id)
840 for (j = top->mols.index[id]; j < top->mols.index[id+1]; ++j)
846 default: /* Should not be reached */
847 gmx_bug("internal error");
853 /* If not using completion, simply store the start of the block. */
854 t->index[t->nr++] = i;
859 /* Set the end of the last block */
860 t->index[t->nr] = t->nra;
861 /* Free any unnecessary memory */
862 srenew(t->index, t->nr+1);
863 t->nalloc_index = t->nr+1;
866 srenew(t->a, t->nra);
867 t->nalloc_a = t->nra;
872 * \param[in] g Index group to check.
873 * \param[in] b Block data to check against.
874 * \returns true if \p g consists of one or more complete blocks from \p b,
877 * The atoms in \p g are assumed to be sorted.
880 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b)
885 /* Each round in the loop matches one block */
888 /* Find the block that begins with the first unmatched atom */
889 while (bi < b->nr && b->index[bi] != g->index[i])
893 /* If not found, or if too large, return */
894 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
898 /* Check that the block matches the index */
899 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
901 if (g->index[i] != j)
906 /* Move the search to the next block */
913 * \param[in] g Index group to check.
914 * \param[in] b Block data to check against.
915 * \returns true if \p g consists of one or more complete blocks from \p b,
918 * The atoms in \p g and \p b->a are assumed to be in the same order.
921 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
926 /* Each round in the loop matches one block */
929 /* Find the block that begins with the first unmatched atom */
930 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
934 /* If not found, or if too large, return */
935 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
939 /* Check that the block matches the index */
940 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
942 if (b->a[j] != g->index[i])
947 /* Move the search to the next block */
954 * \param[in] g Index group to check.
955 * \param[in] type Block data to check against.
956 * \param[in] top Topology data.
957 * \returns true if \p g consists of one or more complete elements of type
958 * \p type, false otherwise.
960 * \p g is assumed to be sorted, otherwise may return false negatives.
962 * If \p type is \ref INDEX_ATOM, the return value is always true.
963 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
967 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
970 // TODO: Consider whether unsorted groups need to be supported better.
986 for (i = 0; i < g->isize; ++i)
989 id = top->atoms.atom[ai].resind;
992 if (ai > 0 && top->atoms.atom[ai-1].resind == id)
996 if (i > 0 && g->index[i-1] < top->atoms.nr - 1
997 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1004 if (g->index[i-1] < top->atoms.nr - 1
1005 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1013 return gmx_ana_index_has_full_blocks(g, &top->mols);
1019 * \param[out] m Output structure.
1021 * Any contents of \p m are discarded without freeing.
1024 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1026 m->type = INDEX_UNKNOWN;
1030 m->mapb.index = NULL;
1031 m->mapb.nalloc_index = 0;
1034 m->mapb.nalloc_a = 0;
1040 m->b.nalloc_index = 0;
1046 * \param[in,out] m Mapping structure.
1047 * \param[in] nr Maximum number of blocks to reserve space for.
1048 * \param[in] isize Maximum number of atoms to reserve space for.
1051 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1053 if (m->mapb.nalloc_index < nr + 1)
1055 srenew(m->refid, nr);
1056 srenew(m->mapid, nr);
1057 srenew(m->orgid, nr);
1058 srenew(m->mapb.index, nr + 1);
1059 srenew(m->b.index, nr + 1);
1060 m->mapb.nalloc_index = nr + 1;
1061 m->b.nalloc_index = nr + 1;
1063 if (m->b.nalloc_a < isize)
1065 srenew(m->b.a, isize);
1066 m->b.nalloc_a = isize;
1071 * \param[in,out] m Mapping structure to initialize.
1072 * \param[in] g Index group to map
1073 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1074 * \param[in] top Topology structure
1075 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1076 * \param[in] type Type of mapping to construct.
1078 * Initializes a new index group mapping.
1079 * The index group provided to gmx_ana_indexmap_update() should always be a
1080 * subset of the \p g given here.
1082 * \p m should have been initialized somehow (calloc() is enough).
1085 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1086 t_topology *top, e_index_t type)
1091 gmx_ana_index_make_block(&m->b, top, g, type, false);
1092 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1093 for (i = mi = 0; i < m->b.nr; ++i)
1095 ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1102 m->orgid[i] = top->atoms.atom[ii].resind;
1105 while (top->mols.index[mi+1] <= ii)
1117 for (i = 0; i < m->b.nr; ++i)
1120 m->mapid[i] = m->orgid[i];
1122 m->mapb.nr = m->b.nr;
1123 m->mapb.nra = m->b.nra;
1125 std::memcpy(m->mapb.index, m->b.index, (m->b.nr+1)*sizeof(*(m->mapb.index)));
1130 * \param[in,out] m Mapping structure to initialize.
1131 * \param[in] b Block information to use for data.
1133 * Frees some memory that is not necessary for static index group mappings.
1134 * Internal pointers are set to point to data in \p b; it is the responsibility
1135 * of the caller to ensure that the block information matches the contents of
1137 * After this function has been called, the index group provided to
1138 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1140 * This function breaks modularity of the index group mapping interface in an
1141 * ugly way, but allows reducing memory usage of static selections by a
1142 * significant amount.
1145 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1148 sfree(m->mapb.index);
1151 m->mapb.nalloc_index = 0;
1152 m->mapb.nalloc_a = 0;
1153 m->b.nalloc_index = 0;
1155 m->mapid = m->orgid;
1156 m->mapb.index = b->index;
1158 m->b.index = b->index;
1163 * \param[in,out] dest Destination data structure.
1164 * \param[in] src Source mapping.
1165 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1166 * copy is made; otherwise, only variable parts are copied, and no memory
1169 * \p dest should have been initialized somehow (calloc() is enough).
1172 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1176 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1177 dest->type = src->type;
1178 dest->b.nr = src->b.nr;
1179 dest->b.nra = src->b.nra;
1180 std::memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1181 std::memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1182 std::memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1184 dest->mapb.nr = src->mapb.nr;
1185 dest->mapb.nra = src->mapb.nra;
1186 if (src->mapb.nalloc_a > 0)
1190 snew(dest->mapb.a, src->mapb.nalloc_a);
1191 dest->mapb.nalloc_a = src->mapb.nalloc_a;
1193 std::memcpy(dest->mapb.a, src->mapb.a, dest->mapb.nra*sizeof(*dest->mapb.a));
1197 dest->mapb.a = src->mapb.a;
1199 std::memcpy(dest->refid, src->refid, dest->mapb.nr*sizeof(*dest->refid));
1200 std::memcpy(dest->mapid, src->mapid, dest->mapb.nr*sizeof(*dest->mapid));
1201 std::memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1202 dest->bStatic = src->bStatic;
1206 * Helper function to set the source atoms in an index map.
1208 * \param[in,out] m Mapping structure.
1209 * \param[in] isize Number of atoms in the \p index array.
1210 * \param[in] index List of atoms.
1213 set_atoms(gmx_ana_indexmap_t *m, int isize, int *index)
1215 m->mapb.nra = isize;
1216 if (m->mapb.nalloc_a == 0)
1222 for (int i = 0; i < isize; ++i)
1224 m->mapb.a[i] = index[i];
1230 * \param[in,out] m Mapping structure.
1231 * \param[in] g Current index group.
1232 * \param[in] bMaskOnly true if the unused blocks should be masked with
1233 * -1 instead of removing them.
1235 * Updates the index group mapping with the new index group \p g.
1237 * \see gmx_ana_indexmap_t
1240 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1245 /* Process the simple cases first */
1246 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1250 if (m->type == INDEX_ALL)
1252 set_atoms(m, g->isize, g->index);
1255 m->mapb.index[1] = g->isize;
1259 /* Reset the reference IDs and mapping if necessary */
1260 const bool bToFull = (g->isize == m->b.nra);
1261 const bool bWasFull = (m->mapb.nra == m->b.nra);
1262 if (bToFull || bMaskOnly)
1266 for (bj = 0; bj < m->b.nr; ++bj)
1273 for (bj = 0; bj < m->b.nr; ++bj)
1275 m->mapid[bj] = m->orgid[bj];
1277 for (bj = 0; bj <= m->b.nr; ++bj)
1279 m->mapb.index[bj] = m->b.index[bj];
1282 set_atoms(m, m->b.nra, m->b.a);
1283 m->mapb.nr = m->b.nr;
1285 /* Exit immediately if the group is static */
1294 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1296 /* Find the next atom in the block */
1297 while (m->b.a[j] != g->index[i])
1301 /* Mark blocks that did not contain any atoms */
1302 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1304 m->refid[bj++] = -1;
1306 /* Advance the block index if we have reached the next block */
1307 if (m->b.index[bj] <= j)
1312 /* Mark the last blocks as not accessible */
1313 while (bj < m->b.nr)
1315 m->refid[bj++] = -1;
1320 set_atoms(m, g->isize, g->index);
1321 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1323 /* Find the next atom in the block */
1324 while (m->b.a[j] != g->index[i])
1328 /* If we have reached a new block, add it */
1329 if (m->b.index[bj+1] <= j)
1331 /* Skip any blocks in between */
1332 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1337 m->mapid[bi] = m->orgid[bj];
1338 m->mapb.index[bi] = i;
1342 /* Update the number of blocks */
1343 m->mapb.index[bi] = g->isize;
1350 * \param[in,out] m Mapping structure to free.
1352 * All the memory allocated for the mapping structure is freed, and
1353 * the pointers set to NULL.
1354 * The pointer \p m is not freed.
1357 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1360 if (m->mapid != m->orgid)
1364 if (m->mapb.nalloc_index > 0)
1366 sfree(m->mapb.index);
1368 if (m->mapb.nalloc_a > 0)
1373 if (m->b.nalloc_index > 0)
1377 if (m->b.nalloc_a > 0)
1381 gmx_ana_indexmap_clear(m);