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)
502 if (*(atom_id *)a > *(atom_id *)b)
510 * \param[in,out] g Index group to be sorted.
513 gmx_ana_index_sort(gmx_ana_index_t *g)
515 qsort(g->index, g->isize, sizeof(*g->index), cmp_atomid);
519 * \param[in] a Index group to check.
520 * \param[in] b Index group to check.
521 * \returns true if \p a and \p b are equal, false otherwise.
524 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
528 if (a->isize != b->isize)
532 for (i = 0; i < a->isize; ++i)
534 if (a->index[i] != b->index[i])
543 * \param[in] a Index group to check against.
544 * \param[in] b Index group to check.
545 * \returns true if \p b is contained in \p a,
548 * If the elements are not in the same order in both groups, the function
549 * fails. However, the groups do not need to be sorted.
552 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
556 for (i = j = 0; j < b->isize; ++i, ++j)
558 while (i < a->isize && a->index[i] != b->index[j])
571 * \param[out] dest Output index group (the intersection of \p a and \p b).
572 * \param[in] a First index group.
573 * \param[in] b Second index group.
575 * \p dest can be the same as \p a or \p b.
578 gmx_ana_index_intersection(gmx_ana_index_t *dest,
579 gmx_ana_index_t *a, gmx_ana_index_t *b)
583 for (i = j = k = 0; i < a->isize && j < b->isize; ++i)
585 while (j < b->isize && b->index[j] < a->index[i])
589 if (j < b->isize && b->index[j] == a->index[i])
591 dest->index[k++] = b->index[j++];
598 * \param[out] dest Output index group (the difference \p a - \p b).
599 * \param[in] a First index group.
600 * \param[in] b Second index group.
602 * \p dest can equal \p a, but not \p b.
605 gmx_ana_index_difference(gmx_ana_index_t *dest,
606 gmx_ana_index_t *a, gmx_ana_index_t *b)
610 for (i = j = k = 0; i < a->isize; ++i)
612 while (j < b->isize && b->index[j] < a->index[i])
616 if (j == b->isize || b->index[j] != a->index[i])
618 dest->index[k++] = a->index[i];
625 * \param[in] a First index group.
626 * \param[in] b Second index group.
627 * \returns Size of the difference \p a - \p b.
630 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
634 for (i = j = k = 0; i < a->isize; ++i)
636 while (j < b->isize && b->index[j] < a->index[i])
640 if (j == b->isize || b->index[j] != a->index[i])
649 * \param[out] dest1 Output group 1 (will equal \p g).
650 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
651 * \param[in] src Group to be partitioned.
652 * \param[in] g One partition.
654 * \pre \p g is a subset of \p src and both sets are sorted
655 * \pre \p dest1 has allocated storage to store \p src
656 * \post \p dest1 == \p g
657 * \post \p dest2 == \p src - \p g
659 * No storage should be allocated for \p dest2; after the call,
660 * \p dest2->index points to the memory allocated for \p dest1
661 * (to a part that is not used by \p dest1).
663 * The calculation can be performed in-place by setting \p dest1 equal to
667 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
668 gmx_ana_index_t *src, gmx_ana_index_t *g)
672 dest2->index = dest1->index + g->isize;
673 dest2->isize = src->isize - g->isize;
674 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
676 while (j >= 0 && src->index[j] != g->index[i])
678 dest2->index[k--] = src->index[j--];
683 dest2->index[k--] = src->index[j--];
685 gmx_ana_index_copy(dest1, g, false);
689 * \param[out] dest Output index group (the union of \p a and \p b).
690 * \param[in] a First index group.
691 * \param[in] b Second index group.
693 * \p a and \p b can have common items.
694 * \p dest can equal \p a or \p b.
696 * \see gmx_ana_index_merge()
699 gmx_ana_index_union(gmx_ana_index_t *dest,
700 gmx_ana_index_t *a, gmx_ana_index_t *b)
705 dsize = gmx_ana_index_difference_size(b, a);
708 dest->isize = a->isize + dsize;
709 for (k = dest->isize - 1; k >= 0; k--)
711 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
713 dest->index[k] = b->index[j--];
717 if (j >= 0 && a->index[i] == b->index[j])
721 dest->index[k] = a->index[i--];
727 * \param[out] dest Output index group (the union of \p a and \p b).
728 * \param[in] a First index group.
729 * \param[in] b Second index group.
731 * \p a and \p b should not have common items.
732 * \p dest can equal \p a or \p b.
734 * \see gmx_ana_index_union()
737 gmx_ana_index_merge(gmx_ana_index_t *dest,
738 gmx_ana_index_t *a, gmx_ana_index_t *b)
744 dest->isize = a->isize + b->isize;
745 for (k = dest->isize - 1; k >= 0; k--)
747 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
749 dest->index[k] = b->index[j--];
753 dest->index[k] = a->index[i--];
758 /********************************************************************
759 * gmx_ana_indexmap_t and related things
760 ********************************************************************/
763 * \param[in,out] t Output block.
764 * \param[in] top Topology structure
765 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
767 * \param[in] g Index group
768 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
769 * \param[in] type Type of partitioning to make.
770 * \param[in] bComplete
771 * If true, the index group is expanded to include any residue/molecule
772 * (depending on \p type) that is partially contained in the group.
773 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
775 * \p m should have been initialized somehow (calloc() is enough) unless
776 * \p type is INDEX_UNKNOWN.
777 * \p g should be sorted.
780 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
781 e_index_t type, bool bComplete)
786 if (type == INDEX_UNKNOWN)
799 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
801 if (type != INDEX_RES && type != INDEX_MOL)
805 /* Allocate memory for the atom array and fill it unless we are using
810 /* We may allocate some extra memory here because we don't know in
811 * advance how much will be needed. */
812 if (t->nalloc_a < top->atoms.nr)
814 srenew(t->a, top->atoms.nr);
815 t->nalloc_a = top->atoms.nr;
821 if (t->nalloc_a < g->isize)
823 srenew(t->a, g->isize);
824 t->nalloc_a = g->isize;
826 memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
829 /* Allocate memory for the block index. We don't know in advance
830 * how much will be needed, so we allocate some extra and free it in the
832 if (t->nalloc_index < g->isize + 1)
834 srenew(t->index, g->isize + 1);
835 t->nalloc_index = g->isize + 1;
839 j = 0; /* j is used by residue completion for the first atom not stored */
841 for (i = 0; i < g->isize; ++i)
844 /* Find the ID number of the atom/residue/molecule corresponding to
852 id = top->atoms.atom[ai].resind;
855 while (ai >= top->mols.index[id+1])
860 case INDEX_UNKNOWN: /* Should not occur */
865 /* If this is the first atom in a new block, initialize the block. */
870 /* For completion, we first set the start of the block. */
871 t->index[t->nr++] = t->nra;
872 /* And then we find all the atoms that should be included. */
876 while (top->atoms.atom[j].resind != id)
880 while (j < top->atoms.nr && top->atoms.atom[j].resind == id)
888 for (j = top->mols.index[id]; j < top->mols.index[id+1]; ++j)
894 default: /* Should not be reached */
895 gmx_bug("internal error");
901 /* If not using completion, simply store the start of the block. */
902 t->index[t->nr++] = i;
907 /* Set the end of the last block */
908 t->index[t->nr] = t->nra;
909 /* Free any unnecessary memory */
910 srenew(t->index, t->nr+1);
911 t->nalloc_index = t->nr+1;
914 srenew(t->a, t->nra);
915 t->nalloc_a = t->nra;
920 * \param[in] g Index group to check.
921 * \param[in] b Block data to check against.
922 * \returns true if \p g consists of one or more complete blocks from \p b,
925 * The atoms in \p g are assumed to be sorted.
928 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b)
933 /* Each round in the loop matches one block */
936 /* Find the block that begins with the first unmatched atom */
937 while (bi < b->nr && b->index[bi] != g->index[i])
941 /* If not found, or if too large, return */
942 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
946 /* Check that the block matches the index */
947 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
949 if (g->index[i] != j)
954 /* Move the search to the next block */
961 * \param[in] g Index group to check.
962 * \param[in] b Block data to check against.
963 * \returns true if \p g consists of one or more complete blocks from \p b,
966 * The atoms in \p g and \p b->a are assumed to be in the same order.
969 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
974 /* Each round in the loop matches one block */
977 /* Find the block that begins with the first unmatched atom */
978 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
982 /* If not found, or if too large, return */
983 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
987 /* Check that the block matches the index */
988 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
990 if (b->a[j] != g->index[i])
995 /* Move the search to the next block */
1002 * \param[in] g Index group to check.
1003 * \param[in] type Block data to check against.
1004 * \param[in] top Topology data.
1005 * \returns true if \p g consists of one or more complete elements of type
1006 * \p type, false otherwise.
1008 * If \p type is \ref INDEX_ATOM, the return value is always true.
1009 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
1013 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
1031 for (i = 0; i < g->isize; ++i)
1034 id = top->atoms.atom[ai].resind;
1037 if (ai > 0 && top->atoms.atom[ai-1].resind == id)
1041 if (i > 0 && g->index[i-1] < top->atoms.nr - 1
1042 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1049 if (g->index[i-1] < top->atoms.nr - 1
1050 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1058 return gmx_ana_index_has_full_blocks(g, &top->mols);
1064 * \param[out] m Output structure.
1066 * Any contents of \p m are discarded without freeing.
1069 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1071 m->type = INDEX_UNKNOWN;
1076 m->mapb.index = NULL;
1077 m->mapb.nalloc_index = 0;
1083 m->b.nalloc_index = 0;
1086 m->bMapStatic = true;
1090 * \param[in,out] m Mapping structure.
1091 * \param[in] nr Maximum number of blocks to reserve space for.
1092 * \param[in] isize Maximum number of atoms to reserve space for.
1095 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1097 if (m->mapb.nalloc_index < nr + 1)
1099 srenew(m->refid, nr);
1100 srenew(m->mapid, nr);
1101 srenew(m->orgid, nr);
1102 srenew(m->mapb.index, nr + 1);
1103 srenew(m->b.index, nr + 1);
1104 m->mapb.nalloc_index = nr + 1;
1105 m->b.nalloc_index = nr + 1;
1107 if (m->b.nalloc_a < isize)
1109 srenew(m->b.a, isize);
1110 m->b.nalloc_a = isize;
1115 * \param[in,out] m Mapping structure to initialize.
1116 * \param[in] g Index group to map
1117 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1118 * \param[in] top Topology structure
1119 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1120 * \param[in] type Type of mapping to construct.
1122 * Initializes a new index group mapping.
1123 * The index group provided to gmx_ana_indexmap_update() should always be a
1124 * subset of the \p g given here.
1126 * \p m should have been initialized somehow (calloc() is enough).
1129 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1130 t_topology *top, e_index_t type)
1135 gmx_ana_index_make_block(&m->b, top, g, type, false);
1136 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1138 for (i = mi = 0; i < m->nr; ++i)
1140 ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1147 m->orgid[i] = top->atoms.atom[ii].resind;
1150 while (top->mols.index[mi+1] <= ii)
1162 for (i = 0; i < m->nr; ++i)
1165 m->mapid[i] = m->orgid[i];
1168 memcpy(m->mapb.index, m->b.index, (m->nr+1)*sizeof(*(m->mapb.index)));
1170 m->bMapStatic = true;
1174 * \param[in,out] m Mapping structure to initialize.
1175 * \param[in] b Block information to use for data.
1177 * Frees some memory that is not necessary for static index group mappings.
1178 * Internal pointers are set to point to data in \p b; it is the responsibility
1179 * of the caller to ensure that the block information matches the contents of
1181 * After this function has been called, the index group provided to
1182 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1184 * This function breaks modularity of the index group mapping interface in an
1185 * ugly way, but allows reducing memory usage of static selections by a
1186 * significant amount.
1189 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1192 m->mapid = m->orgid;
1194 m->b.nalloc_index = 0;
1195 m->b.index = b->index;
1196 sfree(m->mapb.index);
1197 m->mapb.nalloc_index = 0;
1198 m->mapb.index = m->b.index;
1205 * \param[in,out] dest Destination data structure.
1206 * \param[in] src Source mapping.
1207 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1208 * copy is made; otherwise, only variable parts are copied, and no memory
1211 * \p dest should have been initialized somehow (calloc() is enough).
1214 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1218 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1219 dest->type = src->type;
1220 dest->b.nr = src->b.nr;
1221 dest->b.nra = src->b.nra;
1222 memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1223 memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1224 memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1227 dest->mapb.nr = src->mapb.nr;
1228 memcpy(dest->refid, src->refid, dest->nr*sizeof(*dest->refid));
1229 memcpy(dest->mapid, src->mapid, dest->nr*sizeof(*dest->mapid));
1230 memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1231 dest->bStatic = src->bStatic;
1232 dest->bMapStatic = src->bMapStatic;
1236 * \param[in,out] m Mapping structure.
1237 * \param[in] g Current index group.
1238 * \param[in] bMaskOnly true if the unused blocks should be masked with
1239 * -1 instead of removing them.
1241 * Updates the index group mapping with the new index group \p g.
1243 * \see gmx_ana_indexmap_t
1246 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1252 /* Process the simple cases first */
1253 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1257 if (m->type == INDEX_ALL)
1261 m->mapb.index[1] = g->isize;
1265 /* Reset the reference IDs and mapping if necessary */
1266 bStatic = (g->isize == m->b.nra && m->nr == m->b.nr);
1267 if (bStatic || bMaskOnly)
1271 for (bj = 0; bj < m->b.nr; ++bj)
1278 for (bj = 0; bj < m->b.nr; ++bj)
1280 m->mapid[bj] = m->orgid[bj];
1282 for (bj = 0; bj <= m->b.nr; ++bj)
1284 m->mapb.index[bj] = m->b.index[bj];
1286 m->bMapStatic = true;
1289 /* Exit immediately if the group is static */
1299 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1301 /* Find the next atom in the block */
1302 while (m->b.a[j] != g->index[i])
1306 /* Mark blocks that did not contain any atoms */
1307 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1309 m->refid[bj++] = -1;
1311 /* Advance the block index if we have reached the next block */
1312 if (m->b.index[bj] <= j)
1317 /* Mark the last blocks as not accessible */
1318 while (bj < m->b.nr)
1320 m->refid[bj++] = -1;
1325 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1327 /* Find the next atom in the block */
1328 while (m->b.a[j] != g->index[i])
1332 /* If we have reached a new block, add it */
1333 if (m->b.index[bj+1] <= j)
1335 /* Skip any blocks in between */
1336 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1341 m->mapid[bi] = m->orgid[bj];
1342 m->mapb.index[bi] = i;
1346 /* Update the number of blocks */
1347 m->mapb.index[bi] = g->isize;
1349 m->bMapStatic = false;
1356 * \param[in,out] m Mapping structure to free.
1358 * All the memory allocated for the mapping structure is freed, and
1359 * the pointers set to NULL.
1360 * The pointer \p m is not freed.
1363 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1366 if (m->mapid != m->orgid)
1370 if (m->mapb.nalloc_index > 0)
1372 sfree(m->mapb.index);
1375 if (m->b.nalloc_index > 0)
1379 if (m->b.nalloc_a > 0)
1383 gmx_ana_indexmap_clear(m);