3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
33 * Implements functions in indexutil.h.
35 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36 * \ingroup module_selection
38 #include "gromacs/legacyheaders/index.h"
39 #include "gromacs/legacyheaders/gmx_fatal.h"
40 #include "gromacs/legacyheaders/smalloc.h"
41 #include "gromacs/legacyheaders/string2.h"
42 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/selection/indexutil.h"
46 /********************************************************************
47 * gmx_ana_indexgrps_t functions
48 ********************************************************************/
51 * Stores a set of index groups.
53 struct gmx_ana_indexgrps_t
55 /** Number of index groups. */
57 /** Array of index groups. */
62 * \param[out] g Index group structure.
63 * \param[in] ngrps Number of groups for which memory is allocated.
66 gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps)
74 * \param[out] g Index group structure.
75 * \param[in] top Topology structure.
76 * \param[in] fnm File name for the index file.
77 * Memory is automatically allocated.
79 * One or both of \p top or \p fnm can be NULL.
80 * If \p top is NULL, an index file is required and the groups are read
81 * from the file (uses Gromacs routine init_index()).
82 * If \p fnm is NULL, default groups are constructed based on the
83 * topology (uses Gromacs routine analyse()).
84 * If both are null, the index group structure is initialized empty.
87 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
90 t_blocka *block = NULL;
96 block = init_index(fnm, &names);
100 block = new_blocka();
101 analyse(&top->atoms, block, &names, FALSE, FALSE);
111 gmx_ana_indexgrps_alloc(g, block->nr);
112 for (i = 0; i < block->nr; ++i)
114 gmx_ana_index_t *grp = &(*g)->g[i];
116 grp->isize = block->index[i+1] - block->index[i];
117 snew(grp->index, grp->isize);
118 for (j = 0; j < grp->isize; ++j)
120 grp->index[j] = block->a[block->index[i]+j];
122 grp->name = names[i];
123 grp->nalloc_index = grp->isize;
132 * \param[in] g Index groups structure.
134 * The pointer \p g is invalid after the call.
137 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
146 for (i = 0; i < g->nr; ++i)
148 gmx_ana_index_deinit(&g->g[i]);
157 * \param[out] dest Destination index groups.
158 * \param[in] src Source index groups.
160 * A deep copy is made for all fields, including the group names.
163 gmx_ana_indexgrps_clone(gmx_ana_indexgrps_t **dest, gmx_ana_indexgrps_t *src)
167 gmx_ana_indexgrps_alloc(dest, src->nr);
168 for (g = 0; g < src->nr; ++g)
170 gmx_ana_index_copy(&(*dest)->g[g], &src->g[g], true);
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[in] src Input index groups.
204 * \param[in] n Number of the group to extract.
205 * \returns true if \p n is a valid group in \p src, false otherwise.
208 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, int n)
210 if (n < 0 || n >= src->nr)
216 gmx_ana_index_copy(dest, &src->g[n], true);
221 * \param[out] dest Output structure.
222 * \param[in] src Input index groups.
223 * \param[in] name Name (or part of the name) of the group to extract.
224 * \returns true if \p name is a valid group in \p src, false otherwise.
226 * Uses the Gromacs routine find_group() to find the actual group;
227 * the comparison is case-insensitive.
230 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, char *name)
235 snew(names, src->nr);
236 for (i = 0; i < src->nr; ++i)
238 names[i] = src->g[i].name;
240 i = find_group(name, src->nr, names);
248 return gmx_ana_indexgrps_extract(dest, src, i);
252 * \param[in] fp Where to print the output.
253 * \param[in] g Index groups to print.
254 * \param[in] maxn Maximum number of indices to print
255 * (-1 = print all, 0 = print only names).
258 gmx_ana_indexgrps_print(FILE *fp, gmx_ana_indexgrps_t *g, int maxn)
262 for (i = 0; i < g->nr; ++i)
264 fprintf(fp, " %2d: ", i);
265 gmx_ana_index_dump(fp, &g->g[i], i, maxn);
269 /********************************************************************
270 * gmx_ana_index_t functions
271 ********************************************************************/
274 * \param[in,out] g Index group structure.
275 * \param[in] isize Maximum number of atoms to reserve space for.
278 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
280 if (g->nalloc_index < isize)
282 srenew(g->index, isize);
283 g->nalloc_index = isize;
288 * \param[in,out] g Index group structure.
290 * Resizes the memory allocated for holding the indices such that the
291 * current contents fit.
294 gmx_ana_index_squeeze(gmx_ana_index_t *g)
296 srenew(g->index, g->isize);
297 g->nalloc_index = g->isize;
301 * \param[out] g Output structure.
303 * Any contents of \p g are discarded without freeing.
306 gmx_ana_index_clear(gmx_ana_index_t *g)
315 * \param[out] g Output structure.
316 * \param[in] isize Number of atoms in the new group.
317 * \param[in] index Array of \p isize atoms (can be NULL if \p isize is 0).
318 * \param[in] name Name for the new group (can be NULL).
319 * \param[in] nalloc Number of elements allocated for \p index
320 * (if 0, \p index is not freed in gmx_ana_index_deinit())
322 * No copy if \p index is made.
325 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, char *name,
331 g->nalloc_index = nalloc;
335 * \param[out] g Output structure.
336 * \param[in] natoms Number of atoms.
337 * \param[in] name Name for the new group (can be NULL).
340 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms, char *name)
345 snew(g->index, natoms);
346 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)
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 memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
391 if (bAlloc && src->name)
393 dest->name = strdup(src->name);
395 else if (bAlloc || src->name)
397 dest->name = src->name;
402 * \param[in] fp Where to print the output.
403 * \param[in] g Index group to print.
404 * \param[in] i Group number to use if the name is NULL.
405 * \param[in] maxn Maximum number of indices to print (-1 = print all).
408 gmx_ana_index_dump(FILE *fp, gmx_ana_index_t *g, int i, int maxn)
414 fprintf(fp, "\"%s\"", g->name);
418 fprintf(fp, "Group %d", i+1);
420 fprintf(fp, " (%d atoms)", g->isize);
425 if (maxn >= 0 && n > maxn)
429 for (j = 0; j < n; ++j)
431 fprintf(fp, " %d", g->index[j]+1);
442 * \param[in] g Input index group.
443 * \param[in] natoms Number of atoms to check against.
445 * If any atom index in the index group is less than zero or >= \p natoms,
446 * gmx_fatal() is called.
449 gmx_ana_index_check(gmx_ana_index_t *g, int natoms)
453 for (j = 0; j < g->isize; ++j)
455 if (g->index[j] >= natoms)
457 gmx_fatal(FARGS,"Atom index (%d) in index group %s (%d atoms) "
458 "larger than number of atoms in trajectory (%d atoms)",
459 g->index[j], g->name, g->isize, natoms);
461 else if (g->index[j] < 0)
463 gmx_fatal(FARGS,"Atom index (%d) in index group %s (%d atoms) "
465 g->index[j], g->name, g->isize);
471 * \param[in] g Index group to check.
472 * \returns true if the index group is sorted and has no duplicates,
476 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
480 for (i = 0; i < g->isize-1; ++i)
482 if (g->index[i+1] <= g->index[i])
490 /********************************************************************
492 ********************************************************************/
494 /** Helper function for gmx_ana_index_sort(). */
496 cmp_atomid(const void *a, const void *b)
498 if (*(atom_id *)a < *(atom_id *)b) return -1;
499 if (*(atom_id *)a > *(atom_id *)b) return 1;
504 * \param[in,out] g Index group to be sorted.
507 gmx_ana_index_sort(gmx_ana_index_t *g)
509 qsort(g->index, g->isize, sizeof(*g->index), cmp_atomid);
513 * \param[in] a Index group to check.
514 * \param[in] b Index group to check.
515 * \returns true if \p a and \p b are equal, false otherwise.
518 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
522 if (a->isize != b->isize)
526 for (i = 0; i < a->isize; ++i)
528 if (a->index[i] != b->index[i])
537 * \param[in] a Index group to check against.
538 * \param[in] b Index group to check.
539 * \returns true if \p b is contained in \p a,
542 * If the elements are not in the same order in both groups, the function
543 * fails. However, the groups do not need to be sorted.
546 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
550 for (i = j = 0; j < b->isize; ++i, ++j) {
551 while (i < a->isize && a->index[i] != b->index[j])
564 * \param[out] dest Output index group (the intersection of \p a and \p b).
565 * \param[in] a First index group.
566 * \param[in] b Second index group.
568 * \p dest can be the same as \p a or \p b.
571 gmx_ana_index_intersection(gmx_ana_index_t *dest,
572 gmx_ana_index_t *a, gmx_ana_index_t *b)
576 for (i = j = k = 0; i < a->isize && j < b->isize; ++i) {
577 while (j < b->isize && b->index[j] < a->index[i])
581 if (j < b->isize && b->index[j] == a->index[i])
583 dest->index[k++] = b->index[j++];
590 * \param[out] dest Output index group (the difference \p a - \p b).
591 * \param[in] a First index group.
592 * \param[in] b Second index group.
594 * \p dest can equal \p a, but not \p b.
597 gmx_ana_index_difference(gmx_ana_index_t *dest,
598 gmx_ana_index_t *a, gmx_ana_index_t *b)
602 for (i = j = k = 0; i < a->isize; ++i)
604 while (j < b->isize && b->index[j] < a->index[i])
608 if (j == b->isize || b->index[j] != a->index[i])
610 dest->index[k++] = a->index[i];
617 * \param[in] a First index group.
618 * \param[in] b Second index group.
619 * \returns Size of the difference \p a - \p b.
622 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
626 for (i = j = k = 0; i < a->isize; ++i)
628 while (j < b->isize && b->index[j] < a->index[i])
632 if (j == b->isize || b->index[j] != a->index[i])
641 * \param[out] dest1 Output group 1 (will equal \p g).
642 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
643 * \param[in] src Group to be partitioned.
644 * \param[in] g One partition.
646 * \pre \p g is a subset of \p src and both sets are sorted
647 * \pre \p dest1 has allocated storage to store \p src
648 * \post \p dest1 == \p g
649 * \post \p dest2 == \p src - \p g
651 * No storage should be allocated for \p dest2; after the call,
652 * \p dest2->index points to the memory allocated for \p dest1
653 * (to a part that is not used by \p dest1).
655 * The calculation can be performed in-place by setting \p dest1 equal to
659 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
660 gmx_ana_index_t *src, gmx_ana_index_t *g)
664 dest2->index = dest1->index + g->isize;
665 dest2->isize = src->isize - g->isize;
666 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
668 while (j >= 0 && src->index[j] != g->index[i])
670 dest2->index[k--] = src->index[j--];
675 dest2->index[k--] = src->index[j--];
677 gmx_ana_index_copy(dest1, g, false);
681 * \param[out] dest Output index group (the union of \p a and \p b).
682 * \param[in] a First index group.
683 * \param[in] b Second index group.
685 * \p a and \p b can have common items.
686 * \p dest can equal \p a or \p b.
688 * \see gmx_ana_index_merge()
691 gmx_ana_index_union(gmx_ana_index_t *dest,
692 gmx_ana_index_t *a, gmx_ana_index_t *b)
697 dsize = gmx_ana_index_difference_size(b, a);
700 dest->isize = a->isize + dsize;
701 for (k = dest->isize - 1; k >= 0; k--)
703 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
705 dest->index[k] = b->index[j--];
709 if (j >= 0 && a->index[i] == b->index[j])
713 dest->index[k] = a->index[i--];
719 * \param[out] dest Output index group (the union of \p a and \p b).
720 * \param[in] a First index group.
721 * \param[in] b Second index group.
723 * \p a and \p b should not have common items.
724 * \p dest can equal \p a or \p b.
726 * \see gmx_ana_index_union()
729 gmx_ana_index_merge(gmx_ana_index_t *dest,
730 gmx_ana_index_t *a, gmx_ana_index_t *b)
736 dest->isize = a->isize + b->isize;
737 for (k = dest->isize - 1; k >= 0; k--)
739 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
741 dest->index[k] = b->index[j--];
745 dest->index[k] = a->index[i--];
750 /********************************************************************
751 * gmx_ana_indexmap_t and related things
752 ********************************************************************/
755 * \param[in,out] t Output block.
756 * \param[in] top Topology structure
757 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
759 * \param[in] g Index group
760 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
761 * \param[in] type Type of partitioning to make.
762 * \param[in] bComplete
763 * If true, the index group is expanded to include any residue/molecule
764 * (depending on \p type) that is partially contained in the group.
765 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
767 * \p m should have been initialized somehow (calloc() is enough) unless
768 * \p type is INDEX_UNKNOWN.
769 * \p g should be sorted.
772 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
773 e_index_t type, bool bComplete)
778 if (type == INDEX_UNKNOWN)
791 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
793 if (type != INDEX_RES && type != INDEX_MOL)
797 /* Allocate memory for the atom array and fill it unless we are using
802 /* We may allocate some extra memory here because we don't know in
803 * advance how much will be needed. */
804 if (t->nalloc_a < top->atoms.nr)
806 srenew(t->a, top->atoms.nr);
807 t->nalloc_a = top->atoms.nr;
813 if (t->nalloc_a < g->isize)
815 srenew(t->a, g->isize);
816 t->nalloc_a = g->isize;
818 memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
821 /* Allocate memory for the block index. We don't know in advance
822 * how much will be needed, so we allocate some extra and free it in the
824 if (t->nalloc_index < g->isize + 1)
826 srenew(t->index, g->isize + 1);
827 t->nalloc_index = g->isize + 1;
831 j = 0; /* j is used by residue completion for the first atom not stored */
833 for (i = 0; i < g->isize; ++i)
836 /* Find the ID number of the atom/residue/molecule corresponding to
844 id = top->atoms.atom[ai].resind;
847 while (ai >= top->mols.index[id+1])
852 case INDEX_UNKNOWN: /* Should not occur */
857 /* If this is the first atom in a new block, initialize the block. */
862 /* For completion, we first set the start of the block. */
863 t->index[t->nr++] = t->nra;
864 /* And then we find all the atoms that should be included. */
868 while (top->atoms.atom[j].resind != id)
872 while (j < top->atoms.nr && top->atoms.atom[j].resind == id)
880 for (j = top->mols.index[id]; j < top->mols.index[id+1]; ++j)
886 default: /* Should not be reached */
887 gmx_bug("internal error");
893 /* If not using completion, simply store the start of the block. */
894 t->index[t->nr++] = i;
899 /* Set the end of the last block */
900 t->index[t->nr] = t->nra;
901 /* Free any unnecessary memory */
902 srenew(t->index, t->nr+1);
903 t->nalloc_index = t->nr+1;
906 srenew(t->a, t->nra);
907 t->nalloc_a = t->nra;
912 * \param[in] g Index group to check.
913 * \param[in] b Block data to check against.
914 * \returns true if \p g consists of one or more complete blocks from \p b,
917 * The atoms in \p g are assumed to be sorted.
920 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b)
925 /* Each round in the loop matches one block */
928 /* Find the block that begins with the first unmatched atom */
929 while (bi < b->nr && b->index[bi] != g->index[i])
933 /* If not found, or if too large, return */
934 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
938 /* Check that the block matches the index */
939 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
941 if (g->index[i] != j)
946 /* Move the search to the next block */
953 * \param[in] g Index group to check.
954 * \param[in] b Block data to check against.
955 * \returns true if \p g consists of one or more complete blocks from \p b,
958 * The atoms in \p g and \p b->a are assumed to be in the same order.
961 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
966 /* Each round in the loop matches one block */
969 /* Find the block that begins with the first unmatched atom */
970 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
974 /* If not found, or if too large, return */
975 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
979 /* Check that the block matches the index */
980 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
982 if (b->a[j] != g->index[i])
987 /* Move the search to the next block */
994 * \param[in] g Index group to check.
995 * \param[in] type Block data to check against.
996 * \param[in] top Topology data.
997 * \returns true if \p g consists of one or more complete elements of type
998 * \p type, false otherwise.
1000 * If \p type is \ref INDEX_ATOM, the return value is always true.
1001 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
1005 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
1023 for (i = 0; i < g->isize; ++i)
1026 id = top->atoms.atom[ai].resind;
1029 if (ai > 0 && top->atoms.atom[ai-1].resind == id)
1033 if (i > 0 && g->index[i-1] < top->atoms.nr - 1
1034 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1041 if (g->index[i-1] < top->atoms.nr - 1
1042 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1050 return gmx_ana_index_has_full_blocks(g, &top->mols);
1056 * \param[out] m Output structure.
1058 * Any contents of \p m are discarded without freeing.
1061 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1063 m->type = INDEX_UNKNOWN;
1068 m->mapb.index = NULL;
1069 m->mapb.nalloc_index = 0;
1075 m->b.nalloc_index = 0;
1078 m->bMapStatic = true;
1082 * \param[in,out] m Mapping structure.
1083 * \param[in] nr Maximum number of blocks to reserve space for.
1084 * \param[in] isize Maximum number of atoms to reserve space for.
1087 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1089 if (m->mapb.nalloc_index < nr + 1)
1091 srenew(m->refid, nr);
1092 srenew(m->mapid, nr);
1093 srenew(m->orgid, nr);
1094 srenew(m->mapb.index, nr + 1);
1095 srenew(m->b.index, nr + 1);
1096 m->mapb.nalloc_index = nr + 1;
1097 m->b.nalloc_index = nr + 1;
1099 if (m->b.nalloc_a < isize)
1101 srenew(m->b.a, isize);
1102 m->b.nalloc_a = isize;
1107 * \param[in,out] m Mapping structure to initialize.
1108 * \param[in] g Index group to map
1109 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1110 * \param[in] top Topology structure
1111 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1112 * \param[in] type Type of mapping to construct.
1114 * Initializes a new index group mapping.
1115 * The index group provided to gmx_ana_indexmap_update() should always be a
1116 * subset of the \p g given here.
1118 * \p m should have been initialized somehow (calloc() is enough).
1121 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1122 t_topology *top, e_index_t type)
1127 gmx_ana_index_make_block(&m->b, top, g, type, false);
1128 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1130 for (i = mi = 0; i < m->nr; ++i)
1132 ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1139 m->orgid[i] = top->atoms.atom[ii].resind;
1142 while (top->mols.index[mi+1] <= ii)
1154 for (i = 0; i < m->nr; ++i)
1157 m->mapid[i] = m->orgid[i];
1160 memcpy(m->mapb.index, m->b.index, (m->nr+1)*sizeof(*(m->mapb.index)));
1162 m->bMapStatic = true;
1166 * \param[in,out] m Mapping structure to initialize.
1167 * \param[in] b Block information to use for data.
1169 * Frees some memory that is not necessary for static index group mappings.
1170 * Internal pointers are set to point to data in \p b; it is the responsibility
1171 * of the caller to ensure that the block information matches the contents of
1173 * After this function has been called, the index group provided to
1174 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1176 * This function breaks modularity of the index group mapping interface in an
1177 * ugly way, but allows reducing memory usage of static selections by a
1178 * significant amount.
1181 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1184 m->mapid = m->orgid;
1186 m->b.nalloc_index = 0;
1187 m->b.index = b->index;
1188 sfree(m->mapb.index);
1189 m->mapb.nalloc_index = 0;
1190 m->mapb.index = m->b.index;
1197 * \param[in,out] dest Destination data structure.
1198 * \param[in] src Source mapping.
1199 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1200 * copy is made; otherwise, only variable parts are copied, and no memory
1203 * \p dest should have been initialized somehow (calloc() is enough).
1206 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1210 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1211 dest->type = src->type;
1212 dest->b.nr = src->b.nr;
1213 dest->b.nra = src->b.nra;
1214 memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1215 memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1216 memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1219 dest->mapb.nr = src->mapb.nr;
1220 memcpy(dest->refid, src->refid, dest->nr*sizeof(*dest->refid));
1221 memcpy(dest->mapid, src->mapid, dest->nr*sizeof(*dest->mapid));
1222 memcpy(dest->mapb.index, src->mapb.index,(dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1223 dest->bStatic = src->bStatic;
1224 dest->bMapStatic = src->bMapStatic;
1228 * \param[in,out] m Mapping structure.
1229 * \param[in] g Current index group.
1230 * \param[in] bMaskOnly true if the unused blocks should be masked with
1231 * -1 instead of removing them.
1233 * Updates the index group mapping with the new index group \p g.
1235 * \see gmx_ana_indexmap_t
1238 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1244 /* Process the simple cases first */
1245 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1249 if (m->type == INDEX_ALL)
1253 m->mapb.index[1] = g->isize;
1257 /* Reset the reference IDs and mapping if necessary */
1258 bStatic = (g->isize == m->b.nra && m->nr == m->b.nr);
1259 if (bStatic || bMaskOnly)
1263 for (bj = 0; bj < m->b.nr; ++bj)
1270 for (bj = 0; bj < m->b.nr; ++bj)
1272 m->mapid[bj] = m->orgid[bj];
1274 for (bj = 0; bj <= m->b.nr; ++bj)
1276 m->mapb.index[bj] = m->b.index[bj];
1278 m->bMapStatic = true;
1281 /* Exit immediately if the group is static */
1291 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1293 /* Find the next atom in the block */
1294 while (m->b.a[j] != g->index[i])
1298 /* Mark blocks that did not contain any atoms */
1299 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1301 m->refid[bj++] = -1;
1303 /* Advance the block index if we have reached the next block */
1304 if (m->b.index[bj] <= j)
1309 /* Mark the last blocks as not accessible */
1310 while (bj < m->b.nr)
1312 m->refid[bj++] = -1;
1317 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1319 /* Find the next atom in the block */
1320 while (m->b.a[j] != g->index[i])
1324 /* If we have reached a new block, add it */
1325 if (m->b.index[bj+1] <= j)
1327 /* Skip any blocks in between */
1328 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1333 m->mapid[bi] = m->orgid[bj];
1334 m->mapb.index[bi] = i;
1338 /* Update the number of blocks */
1339 m->mapb.index[bi] = g->isize;
1341 m->bMapStatic = false;
1348 * \param[in,out] m Mapping structure to free.
1350 * All the memory allocated for the mapping structure is freed, and
1351 * the pointers set to NULL.
1352 * The pointer \p m is not freed.
1355 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1358 if (m->mapid != m->orgid)
1362 if (m->mapb.nalloc_index > 0)
1364 sfree(m->mapb.index);
1367 if (m->b.nalloc_index > 0)
1371 if (m->b.nalloc_a > 0)
1375 gmx_ana_indexmap_clear(m);