A bit more modular selection position handling.
[gromacs/qmmm-gamess-us.git] / src / gmxlib / trajana / indexutil.c
blob95351f4fb6da40c1f1843fcfd4a060bf669029ea
1 /*
3 * This source code is part of
5 * G R O M A C S
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
31 /*! \internal \file
32 * \brief Implementation of functions in indexutil.h.
34 #ifdef HAVE_CONFIG_H
35 #include <config.h>
36 #endif
38 #include <index.h>
39 #include <smalloc.h>
40 #include <string2.h>
41 #include <typedefs.h>
43 #include <indexutil.h>
45 /********************************************************************
46 * gmx_ana_indexgrps_t functions
47 ********************************************************************/
49 /*! \internal \brief
50 * Stores a set of index groups.
52 struct gmx_ana_indexgrps_t
54 /** Number of index groups. */
55 int nr;
56 /** Array of index groups. */
57 gmx_ana_index_t *g;
60 /*!
61 * \param[out] g Index group structure.
62 * \param[in] ngrps Number of groups for which memory is allocated.
64 void
65 gmx_ana_indexgrps_alloc(gmx_ana_indexgrps_t **g, int ngrps)
67 snew(*g, 1);
68 (*g)->nr = ngrps;
69 snew((*g)->g, ngrps);
72 /*!
73 * \param[out] g Index group structure.
74 * \param[in] ngrps Number of index groups.
75 * \param[in] isize Array of index group sizes.
76 * \param[in] index Array of pointers to indices of each group.
77 * \param[in] name Array of names of the groups.
78 * \param[in] bFree If TRUE, the \p isize, \p index and \p name arrays
79 * are freed after they have been copied.
81 void
82 gmx_ana_indexgrps_set(gmx_ana_indexgrps_t **g, int ngrps, int *isize,
83 atom_id **index, char **name, bool bFree)
85 int i;
87 gmx_ana_indexgrps_alloc(g, ngrps);
88 for (i = 0; i < ngrps; ++i)
90 gmx_ana_index_set(&(*g)->g[i], isize[i], index[i], name[i], isize[i]);
92 if (bFree)
94 sfree(isize);
95 sfree(index);
96 sfree(name);
101 * \param[out] g Index group structure.
102 * \param[in] top Topology structure.
103 * \param[in] fnm File name for the index file.
104 * Memory is automatically allocated.
106 * One or both of \p top or \p fnm can be NULL.
107 * If \p top is NULL, an index file is required and the groups are read
108 * from the file (uses Gromacs routine init_index()).
109 * If \p fnm is NULL, default groups are constructed based on the
110 * topology (uses Gromacs routine analyse()).
111 * If both are null, the index group structure is initialized empty.
113 void
114 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, t_topology *top,
115 const char *fnm)
117 t_blocka *block = NULL;
118 char **names = NULL;
119 int i, j;
121 if (fnm)
123 block = init_index(fnm, &names);
125 else if (top)
127 block = new_blocka();
128 analyse(&top->atoms, block, &names, FALSE, FALSE);
130 else
132 snew(*g, 1);
133 (*g)->nr = 0;
134 (*g)->g = NULL;
135 return;
138 gmx_ana_indexgrps_alloc(g, block->nr);
139 for (i = 0; i < block->nr; ++i)
141 gmx_ana_index_t *grp = &(*g)->g[i];
143 grp->isize = block->index[i+1] - block->index[i];
144 snew(grp->index, grp->isize);
145 for (j = 0; j < grp->isize; ++j)
147 grp->index[j] = block->a[block->index[i]+j];
149 grp->name = names[i];
150 grp->nalloc_index = grp->isize;
153 done_blocka(block);
154 sfree(block);
155 sfree(names);
159 * \param[out] g Index group structure.
160 * \param[in] top Topology structure.
161 * \param[in] fnm File name for the index file.
162 * \param[in] ngrps Number of required groups.
163 * Memory is automatically allocated.
165 * One of \p top or \p fnm can be NULL, but not both.
166 * If \p top is NULL, an index file is required and the groups are read
167 * from the file (uses Gromacs routine rd_index()).
168 * If \p fnm is NULL, default groups are constructed based on the
169 * topology (uses Gromacs routine get_index()).
171 void
172 gmx_ana_indexgrps_get(gmx_ana_indexgrps_t **g, t_topology *top,
173 const char *fnm, int ngrps)
175 int *isize;
176 atom_id **index;
177 char **name;
179 snew(isize, ngrps);
180 snew(index, ngrps);
181 snew(name, ngrps);
182 if (!top)
184 rd_index(fnm, ngrps, isize, index, name);
186 else
188 get_index(&(top->atoms), fnm, ngrps, isize, index, name);
190 gmx_ana_indexgrps_set(g, ngrps, isize, index, name, TRUE);
194 * \param[out] g Index group structure.
195 * \param[in] fnm File name for the index file.
196 * \param[in] ngrps Number of required groups.
197 * Memory is automatically allocated.
199 * This is a convenience function for calling the Gromacs routine
200 * rd_index().
202 void
203 gmx_ana_indexgrps_rd(gmx_ana_indexgrps_t **g, const char *fnm, int ngrps)
205 gmx_ana_indexgrps_get(g, NULL, fnm, ngrps);
209 * \param[in] g Index groups structure.
211 * The pointer \p g is invalid after the call.
213 void
214 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
216 int i;
218 if (g->nr == 0)
220 sfree(g);
221 return;
223 for (i = 0; i < g->nr; ++i)
225 gmx_ana_index_deinit(&g->g[i]);
227 sfree(g->g);
228 g->nr = 0;
229 g->g = NULL;
230 sfree(g);
234 * \param[out] dest Destination index groups.
235 * \param[in] src Source index groups.
237 * A deep copy is made for all fields, including the group names.
239 void
240 gmx_ana_indexgrps_clone(gmx_ana_indexgrps_t **dest, gmx_ana_indexgrps_t *src)
242 int g;
244 gmx_ana_indexgrps_alloc(dest, src->nr);
245 for (g = 0; g < src->nr; ++g)
247 gmx_ana_index_copy(&(*dest)->g[g], &src->g[g], TRUE);
252 * \param[out] g Index group structure.
253 * \returns TRUE if \p g is empty, i.e., has 0 index groups.
255 bool
256 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g)
258 return g->nr == 0;
262 * \param[in] g Index groups structure.
263 * \param[in] n Index group number to get.
264 * \returns Pointer to the \p n'th index group in \p g.
266 * The returned pointer should not be freed.
268 gmx_ana_index_t *
269 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n)
271 if (n < 0 || n >= g->nr)
273 return NULL;
275 return &g->g[n];
279 * \param[out] dest Output structure.
280 * \param[in] src Input index groups.
281 * \param[in] n Number of the group to extract.
282 * \returns TRUE if \p n is a valid group in \p src, FALSE otherwise.
284 bool
285 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, int n)
287 if (n < 0 || n >= src->nr)
289 dest->isize = 0;
290 return FALSE;
293 gmx_ana_index_copy(dest, &src->g[n], TRUE);
294 return TRUE;
298 * \param[out] dest Output structure.
299 * \param[in] src Input index groups.
300 * \param[in] name Name (or part of the name) of the group to extract.
301 * \returns TRUE if \p name is a valid group in \p src, FALSE otherwise.
303 * Uses the Gromacs routine find_group() to find the actual group;
304 * the comparison is case-insensitive.
306 bool
307 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, gmx_ana_indexgrps_t *src, char *name)
309 int i;
310 char **names;
312 snew(names, src->nr);
313 for (i = 0; i < src->nr; ++i)
315 names[i] = src->g[i].name;
317 i = find_group(name, src->nr, names);
318 sfree(names);
319 if (i == NOTSET)
321 dest->isize = 0;
322 return FALSE;
325 return gmx_ana_indexgrps_extract(dest, src, i);
329 * \param[in] g Index groups to print.
330 * \param[in] maxn Maximum number of indices to print
331 * (-1 = print all, 0 = print only names).
333 void
334 gmx_ana_indexgrps_print(gmx_ana_indexgrps_t *g, int maxn)
336 int i;
338 for (i = 0; i < g->nr; ++i)
340 fprintf(stderr, " %2d: ", i);
341 gmx_ana_index_dump(&g->g[i], i, maxn);
345 /********************************************************************
346 * gmx_ana_index_t functions
347 ********************************************************************/
350 * \param[in,out] g Index group structure.
351 * \param[in] isize Maximum number of atoms to reserve space for.
353 void
354 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
356 if (g->nalloc_index < isize)
358 srenew(g->index, isize);
359 g->nalloc_index = isize;
364 * \param[in,out] g Index group structure.
366 * Resizes the memory allocated for holding the indices such that the
367 * current contents fit.
369 void
370 gmx_ana_index_squeeze(gmx_ana_index_t *g)
372 srenew(g->index, g->isize);
373 g->nalloc_index = g->isize;
377 * \param[out] g Output structure.
379 * Any contents of \p g are discarded without freeing.
381 void
382 gmx_ana_index_clear(gmx_ana_index_t *g)
384 g->isize = 0;
385 g->index = NULL;
386 g->name = NULL;
387 g->nalloc_index = 0;
391 * \param[out] g Output structure.
392 * \param[in] isize Number of atoms in the new group.
393 * \param[in] index Array of \p isize atoms (can be NULL if \p isize is 0).
394 * \param[in] name Name for the new group (can be NULL).
395 * \param[in] nalloc Number of elements allocated for \p index
396 * (if 0, \p index is not freed in gmx_ana_index_deinit())
398 * No copy if \p index is made.
400 void
401 gmx_ana_index_set(gmx_ana_index_t *g, int isize, atom_id *index, char *name,
402 int nalloc)
404 g->isize = isize;
405 g->index = index;
406 g->name = name;
407 g->nalloc_index = nalloc;
411 * \param[out] g Output structure.
412 * \param[in] natoms Number of atoms.
413 * \param[in] name Name for the new group (can be NULL).
415 void
416 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms, char *name)
418 int i;
420 g->isize = natoms;
421 snew(g->index, natoms);
422 for (i = 0; i < natoms; ++i)
424 g->index[i] = i;
426 g->name = name;
427 g->nalloc_index = natoms;
431 * \param[in] g Index group structure.
433 * The pointer \p g is not freed.
435 void
436 gmx_ana_index_deinit(gmx_ana_index_t *g)
438 if (g->nalloc_index > 0)
440 sfree(g->index);
442 sfree(g->name);
443 gmx_ana_index_clear(g);
447 * \param[out] dest Destination index group.
448 * \param[in] src Source index group.
449 * \param[in] bAlloc If TRUE, memory is allocated at \p dest; otherwise,
450 * it is assumed that enough memory has been allocated for index.
452 * A deep copy of the name is only made if \p bAlloc is TRUE.
454 void
455 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
457 dest->isize = src->isize;
458 if (dest->isize > 0)
460 if (bAlloc)
462 snew(dest->index, dest->isize);
463 dest->nalloc_index = dest->isize;
465 memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
467 if (bAlloc && src->name)
469 dest->name = strdup(src->name);
471 else if (bAlloc || src->name)
473 dest->name = src->name;
478 * \param[in] g Index group to print.
479 * \param[in] i Group number to use if the name is NULL.
480 * \param[in] maxn Maximum number of indices to print (-1 = print all).
482 void
483 gmx_ana_index_dump(gmx_ana_index_t *g, int i, int maxn)
485 int j, n;
487 if (g->name)
489 fprintf(stderr, "\"%s\"", g->name);
491 else
493 fprintf(stderr, "Group %d", i+1);
495 fprintf(stderr, " (%d atoms)", g->isize);
496 if (maxn != 0)
498 fprintf(stderr, ":");
499 n = g->isize;
500 if (maxn >= 0 && n > maxn)
502 n = maxn;
504 for (j = 0; j < n; ++j)
506 fprintf(stderr, " %d", g->index[j]+1);
508 if (n < g->isize)
510 fprintf(stderr, " ...");
513 fprintf(stderr, "\n");
517 * \param[in] g Input index group.
518 * \param[in] natoms Number of atoms to check against.
520 * If any atom index in the index group is less than zero or >= \p natoms,
521 * gmx_fatal() is called.
523 void
524 gmx_ana_index_check(gmx_ana_index_t *g, int natoms)
526 int j;
528 for (j = 0; j < g->isize; ++j)
530 if (g->index[j] >= natoms)
532 gmx_fatal(FARGS,"Atom index (%d) in index group %s (%d atoms) "
533 "larger than number of atoms in trajectory (%d atoms)",
534 g->index[j], g->name, g->isize, natoms);
536 else if (g->index[j] < 0)
538 gmx_fatal(FARGS,"Atom index (%d) in index group %s (%d atoms) "
539 "is less than zero",
540 g->index[j], g->name, g->isize);
546 * \param[in] g Index group to check.
547 * \returns TRUE if the index group is sorted and has no duplicates,
548 * FALSE otherwise.
550 bool
551 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
553 int i;
555 for (i = 0; i < g->isize-1; ++i)
557 if (g->index[i+1] <= g->index[i])
559 return FALSE;
562 return TRUE;
565 /********************************************************************
566 * Set operations
567 ********************************************************************/
569 /** Helper function for gmx_ana_index_sort(). */
570 static int
571 cmp_atomid(const void *a, const void *b)
573 if (*(atom_id *)a < *(atom_id *)b) return -1;
574 if (*(atom_id *)a > *(atom_id *)b) return 1;
575 return 0;
579 * \param[in,out] g Index group to be sorted.
581 void
582 gmx_ana_index_sort(gmx_ana_index_t *g)
584 qsort(g->index, g->isize, sizeof(*g->index), cmp_atomid);
588 * \param[in] a Index group to check.
589 * \param[in] b Index group to check.
590 * \returns TRUE if \p a and \p b are equal, FALSE otherwise.
592 bool
593 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
595 int i;
597 if (a->isize != b->isize)
599 return FALSE;
601 for (i = 0; i < a->isize; ++i)
603 if (a->index[i] != b->index[i])
605 return FALSE;
608 return TRUE;
612 * \param[in] a Index group to check against.
613 * \param[in] b Index group to check.
614 * \returns TRUE if \p b is contained in \p a,
615 * FALSE otherwise.
617 * If the elements are not in the same order in both groups, the function
618 * fails. However, the groups do not need to be sorted.
620 bool
621 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
623 int i, j;
625 for (i = j = 0; j < b->isize; ++i, ++j) {
626 while (i < a->isize && a->index[i] != b->index[j])
628 ++i;
630 if (i == a->isize)
632 return FALSE;
635 return TRUE;
639 * \param[out] dest Output index group (the intersection of \p a and \p b).
640 * \param[in] a First index group.
641 * \param[in] b Second index group.
643 * \p dest can be the same as \p a or \p b.
645 void
646 gmx_ana_index_intersection(gmx_ana_index_t *dest,
647 gmx_ana_index_t *a, gmx_ana_index_t *b)
649 int i, j, k;
651 for (i = j = k = 0; i < a->isize && j < b->isize; ++i) {
652 while (j < b->isize && b->index[j] < a->index[i])
654 ++j;
656 if (j < b->isize && b->index[j] == a->index[i])
658 dest->index[k++] = b->index[j++];
661 dest->isize = k;
665 * \param[out] dest Output index group (the difference \p a - \p b).
666 * \param[in] a First index group.
667 * \param[in] b Second index group.
669 * \p dest can equal \p a, but not \p b.
671 void
672 gmx_ana_index_difference(gmx_ana_index_t *dest,
673 gmx_ana_index_t *a, gmx_ana_index_t *b)
675 int i, j, k;
677 for (i = j = k = 0; i < a->isize; ++i)
679 while (j < b->isize && b->index[j] < a->index[i])
681 ++j;
683 if (j == b->isize || b->index[j] != a->index[i])
685 dest->index[k++] = a->index[i];
688 dest->isize = k;
692 * \param[in] a First index group.
693 * \param[in] b Second index group.
694 * \returns Size of the difference \p a - \p b.
697 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
699 int i, j, k;
701 for (i = j = k = 0; i < a->isize; ++i)
703 while (j < b->isize && b->index[j] < a->index[i])
705 ++j;
707 if (j == b->isize || b->index[j] != a->index[i])
709 ++k;
712 return k;
716 * \param[out] dest1 Output group 1 (will equal \p g).
717 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
718 * \param[in] src Group to be partitioned.
719 * \param[in] g One partition.
721 * \pre \p g is a subset of \p src and both sets are sorted
722 * \pre \p dest1 has allocated storage to store \p src
723 * \post \p dest1 == \p g
724 * \post \p dest2 == \p src - \p g
726 * No storage should be allocated for \p dest2; after the call,
727 * \p dest2->index points to the memory allocated for \p dest1
728 * (to a part that is not used by \p dest1).
730 * The calculation can be performed in-place by setting \p dest1 equal to
731 * \p src.
733 void
734 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
735 gmx_ana_index_t *src, gmx_ana_index_t *g)
738 int i, j, k;
740 dest2->index = dest1->index + g->isize;
741 dest2->isize = src->isize - g->isize;
742 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
744 while (j >= 0 && src->index[j] != g->index[i])
746 dest2->index[k--] = src->index[j--];
749 while (j >= 0)
751 dest2->index[k--] = src->index[j--];
753 gmx_ana_index_copy(dest1, g, FALSE);
757 * \param[out] dest Output index group (the union of \p a and \p b).
758 * \param[in] a First index group.
759 * \param[in] b Second index group.
761 * \p a and \p b can have common items.
762 * \p dest can equal \p a or \p b.
764 * \see gmx_ana_index_merge()
766 void
767 gmx_ana_index_union(gmx_ana_index_t *dest,
768 gmx_ana_index_t *a, gmx_ana_index_t *b)
770 int dsize;
771 int i, j, k;
773 dsize = gmx_ana_index_difference_size(b, a);
774 i = a->isize - 1;
775 j = b->isize - 1;
776 dest->isize = a->isize + dsize;
777 for (k = dest->isize - 1; k >= 0; k--)
779 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
781 dest->index[k] = b->index[j--];
783 else
785 if (j >= 0 && a->index[i] == b->index[j])
787 --j;
789 dest->index[k] = a->index[i--];
795 * \param[out] dest Output index group (the union of \p a and \p b).
796 * \param[in] a First index group.
797 * \param[in] b Second index group.
799 * \p a and \p b should not have common items.
800 * \p dest can equal \p a or \p b.
802 * \see gmx_ana_index_union()
804 void
805 gmx_ana_index_merge(gmx_ana_index_t *dest,
806 gmx_ana_index_t *a, gmx_ana_index_t *b)
808 int i, j, k;
810 i = a->isize - 1;
811 j = b->isize - 1;
812 dest->isize = a->isize + b->isize;
813 for (k = dest->isize - 1; k >= 0; k--)
815 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
817 dest->index[k] = b->index[j--];
819 else
821 dest->index[k] = a->index[i--];
826 /********************************************************************
827 * gmx_ana_indexmap_t and related things
828 ********************************************************************/
831 * \param[in,out] t Output block.
832 * \param[in] top Topology structure
833 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
834 * otherwise).
835 * \param[in] g Index group
836 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
837 * \param[in] type Type of partitioning to make.
838 * \param[in] bComplete
839 * If TRUE, the index group is expanded to include any residue/molecule
840 * (depending on \p type) that is partially contained in the group.
841 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
843 * \p m should have been initialized somehow (calloc() is enough) unless
844 * \p type is INDEX_UNKNOWN.
845 * \p g should be sorted.
847 void
848 gmx_ana_index_make_block(t_blocka *t, t_topology *top, gmx_ana_index_t *g,
849 e_index_t type, bool bComplete)
851 int i, j, ai;
852 int id, cur;
854 if (type == INDEX_UNKNOWN)
856 t->nr = 1;
857 snew(t->index, 2);
858 t->nalloc_index = 2;
859 t->index[0] = 0;
860 t->index[1] = 0;
861 t->nra = 0;
862 t->a = NULL;
863 t->nalloc_a = 0;
864 return;
867 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
868 * off otherwise. */
869 if (type != INDEX_RES && type != INDEX_MOL)
871 bComplete = FALSE;
873 /* Allocate memory for the atom array and fill it unless we are using
874 * completion. */
875 if (bComplete)
877 t->nra = 0;
878 /* We may allocate some extra memory here because we don't know in
879 * advance how much will be needed. */
880 if (t->nalloc_a < top->atoms.nr)
882 srenew(t->a, top->atoms.nr);
883 t->nalloc_a = top->atoms.nr;
886 else
888 t->nra = g->isize;
889 if (t->nalloc_a < g->isize)
891 srenew(t->a, g->isize);
892 t->nalloc_a = g->isize;
894 memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
897 /* Allocate memory for the block index. We don't know in advance
898 * how much will be needed, so we allocate some extra and free it in the
899 * end. */
900 if (t->nalloc_index < g->isize + 1)
902 srenew(t->index, g->isize + 1);
903 t->nalloc_index = g->isize + 1;
905 /* Clear counters */
906 t->nr = 0;
907 j = 0; /* j is used by residue completion for the first atom not stored */
908 id = cur = -1;
909 for (i = 0; i < g->isize; ++i)
911 ai = g->index[i];
912 /* Find the ID number of the atom/residue/molecule corresponding to
913 * atom ai. */
914 switch (type)
916 case INDEX_ATOM:
917 id = ai;
918 break;
919 case INDEX_RES:
920 id = top->atoms.atom[ai].resind;
921 break;
922 case INDEX_MOL:
923 while (ai >= top->mols.index[id+1])
925 id++;
927 break;
928 case INDEX_UNKNOWN: /* Should not occur */
929 case INDEX_ALL:
930 id = 0;
931 break;
933 /* If this is the first atom in a new block, initialize the block. */
934 if (id != cur)
936 if (bComplete)
938 /* For completion, we first set the start of the block. */
939 t->index[t->nr++] = t->nra;
940 /* And then we find all the atoms that should be included. */
941 switch (type)
943 case INDEX_RES:
944 while (top->atoms.atom[j].resind != id)
946 ++j;
948 while (j < top->atoms.nr && top->atoms.atom[j].resind == id)
950 t->a[t->nra++] = j;
951 ++j;
953 break;
955 case INDEX_MOL:
956 for (j = top->mols.index[id]; j < top->mols.index[id+1]; ++j)
958 t->a[t->nra++] = j;
960 break;
962 default: /* Should not be reached */
963 gmx_bug("internal error");
964 break;
967 else
969 /* If not using completion, simply store the start of the block. */
970 t->index[t->nr++] = i;
972 cur = id;
975 /* Set the end of the last block */
976 t->index[t->nr] = t->nra;
977 /* Free any unnecessary memory */
978 srenew(t->index, t->nr+1);
979 t->nalloc_index = t->nr+1;
980 if (bComplete)
982 srenew(t->a, t->nra);
983 t->nalloc_a = t->nra;
988 * \param[in] g Index group to check.
989 * \param[in] b Block data to check against.
990 * \returns TRUE if \p g consists of one or more complete blocks from \p b,
991 * FALSE otherwise.
993 * The atoms in \p g are assumed to be sorted.
995 bool
996 gmx_ana_index_has_full_blocks(gmx_ana_index_t *g, t_block *b)
998 int i, j, bi;
1000 i = bi = 0;
1001 /* Each round in the loop matches one block */
1002 while (i < g->isize)
1004 /* Find the block that begins with the first unmatched atom */
1005 while (bi < b->nr && b->index[bi] != g->index[i])
1007 ++bi;
1009 /* If not found, or if too large, return */
1010 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
1012 return FALSE;
1014 /* Check that the block matches the index */
1015 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
1017 if (g->index[i] != j)
1019 return FALSE;
1022 /* Move the search to the next block */
1023 ++bi;
1025 return TRUE;
1029 * \param[in] g Index group to check.
1030 * \param[in] b Block data to check against.
1031 * \returns TRUE if \p g consists of one or more complete blocks from \p b,
1032 * FALSE otherwise.
1034 * The atoms in \p g and \p b->a are assumed to be in the same order.
1036 bool
1037 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
1039 int i, j, bi;
1041 i = bi = 0;
1042 /* Each round in the loop matches one block */
1043 while (i < g->isize)
1045 /* Find the block that begins with the first unmatched atom */
1046 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
1048 ++bi;
1050 /* If not found, or if too large, return */
1051 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
1053 return FALSE;
1055 /* Check that the block matches the index */
1056 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
1058 if (b->a[j] != g->index[i])
1060 return FALSE;
1063 /* Move the search to the next block */
1064 ++bi;
1066 return TRUE;
1070 * \param[in] g Index group to check.
1071 * \param[in] type Block data to check against.
1072 * \param[in] top Topology data.
1073 * \returns TRUE if \p g consists of one or more complete elements of type
1074 * \p type, FALSE otherwise.
1076 * If \p type is \ref INDEX_ATOM, the return value is always TRUE.
1077 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
1078 * always FALSE.
1080 bool
1081 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
1082 t_topology *top)
1084 switch (type)
1086 case INDEX_UNKNOWN:
1087 case INDEX_ALL:
1088 return FALSE;
1090 case INDEX_ATOM:
1091 return TRUE;
1093 case INDEX_RES:
1095 int i, ai;
1096 int id, prev;
1098 prev = -1;
1099 for (i = 0; i < g->isize; ++i)
1101 ai = g->index[i];
1102 id = top->atoms.atom[ai].resind;
1103 if (id != prev)
1105 if (ai > 0 && top->atoms.atom[ai-1].resind == id)
1107 return FALSE;
1109 if (i > 0 && g->index[i-1] < top->atoms.nr - 1
1110 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1112 return FALSE;
1115 prev = id;
1117 if (g->index[i-1] < top->atoms.nr - 1
1118 && top->atoms.atom[g->index[i-1]+1].resind == prev)
1120 return FALSE;
1122 break;
1125 case INDEX_MOL:
1126 return gmx_ana_index_has_full_blocks(g, &top->mols);
1128 return TRUE;
1132 * \param[out] m Output structure.
1134 * Any contents of \p m are discarded without freeing.
1136 void
1137 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1139 m->type = INDEX_UNKNOWN;
1140 m->nr = 0;
1141 m->refid = NULL;
1142 m->mapid = NULL;
1143 m->mapb.nr = 0;
1144 m->mapb.index = NULL;
1145 m->mapb.nalloc_index = 0;
1146 m->orgid = NULL;
1147 m->b.nr = 0;
1148 m->b.index = NULL;
1149 m->b.nra = 0;
1150 m->b.a = NULL;
1151 m->b.nalloc_index = 0;
1152 m->b.nalloc_a = 0;
1153 m->bStatic = TRUE;
1154 m->bMapStatic = TRUE;
1158 * \param[in,out] m Mapping structure.
1159 * \param[in] nr Maximum number of blocks to reserve space for.
1160 * \param[in] isize Maximum number of atoms to reserve space for.
1162 void
1163 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1165 if (m->mapb.nalloc_index < nr + 1)
1167 srenew(m->refid, nr);
1168 srenew(m->mapid, nr);
1169 srenew(m->orgid, nr);
1170 srenew(m->mapb.index, nr + 1);
1171 srenew(m->b.index, nr + 1);
1172 m->mapb.nalloc_index = nr + 1;
1173 m->b.nalloc_index = nr + 1;
1175 if (m->b.nalloc_a < isize)
1177 srenew(m->b.a, isize);
1178 m->b.nalloc_a = isize;
1183 * \param[in,out] m Mapping structure to initialize.
1184 * \param[in] g Index group to map
1185 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1186 * \param[in] top Topology structure
1187 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1188 * \param[in] type Type of mapping to construct.
1190 * Initializes a new index group mapping.
1191 * The index group provided to gmx_ana_indexmap_update() should always be a
1192 * subset of the \p g given here.
1194 * \p m should have been initialized somehow (calloc() is enough).
1196 void
1197 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1198 t_topology *top, e_index_t type)
1200 int i, ii, mi;
1202 m->type = type;
1203 gmx_ana_index_make_block(&m->b, top, g, type, FALSE);
1204 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1205 m->nr = m->b.nr;
1206 for (i = mi = 0; i < m->nr; ++i)
1208 ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1209 switch (type)
1211 case INDEX_ATOM:
1212 m->orgid[i] = ii;
1213 break;
1214 case INDEX_RES:
1215 m->orgid[i] = top->atoms.atom[ii].resind;
1216 break;
1217 case INDEX_MOL:
1218 while (top->mols.index[mi+1] <= ii)
1220 ++mi;
1222 m->orgid[i] = mi;
1223 break;
1224 case INDEX_ALL:
1225 case INDEX_UNKNOWN:
1226 m->orgid[i] = 0;
1227 break;
1230 for (i = 0; i < m->nr; ++i)
1232 m->refid[i] = i;
1233 m->mapid[i] = m->orgid[i];
1235 m->mapb.nr = m->nr;
1236 memcpy(m->mapb.index, m->b.index, (m->nr+1)*sizeof(*(m->mapb.index)));
1237 m->bStatic = TRUE;
1238 m->bMapStatic = TRUE;
1242 * \param[in,out] dest Destination data structure.
1243 * \param[in] src Source mapping.
1244 * \param[in] bFirst If TRUE, memory is allocated for \p dest and a full
1245 * copy is made; otherwise, only variable parts are copied, and no memory
1246 * is allocated.
1248 * \p dest should have been initialized somehow (calloc() is enough).
1250 void
1251 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1253 if (bFirst)
1255 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1256 dest->type = src->type;
1257 dest->b.nr = src->b.nr;
1258 dest->b.nra = src->b.nra;
1259 memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1260 memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1261 memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1263 dest->nr = src->nr;
1264 dest->mapb.nr = src->mapb.nr;
1265 memcpy(dest->refid, src->refid, dest->nr*sizeof(*dest->refid));
1266 memcpy(dest->mapid, src->mapid, dest->nr*sizeof(*dest->mapid));
1267 memcpy(dest->mapb.index, src->mapb.index,(dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1268 dest->bStatic = src->bStatic;
1269 dest->bMapStatic = src->bMapStatic;
1273 * \param[in,out] m Mapping structure.
1274 * \param[in] g Current index group.
1275 * \param[in] bMaskOnly TRUE if the unused blocks should be masked with
1276 * -1 instead of removing them.
1278 * Updates the index group mapping with the new index group \p g.
1280 * \see gmx_ana_indexmap_t
1282 void
1283 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1284 bool bMaskOnly)
1286 int i, j, bi, bj;
1287 bool bStatic;
1289 /* Process the simple cases first */
1290 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1292 return;
1294 if (m->type == INDEX_ALL)
1296 m->mapb.index[1] = g->isize;
1297 return;
1299 /* Reset the reference IDs and mapping if necessary */
1300 bStatic = (g->isize == m->b.nra && m->nr == m->b.nr);
1301 if (bStatic || bMaskOnly)
1303 if (!m->bStatic)
1305 for (bj = 0; bj < m->b.nr; ++bj)
1307 m->refid[bj] = bj;
1310 if (!m->bMapStatic)
1312 for (bj = 0; bj < m->b.nr; ++bj)
1314 m->mapid[bj] = m->orgid[bj];
1316 for (bj = 0; bj <= m->b.nr; ++bj)
1318 m->mapb.index[bj] = m->b.index[bj];
1320 m->bMapStatic = TRUE;
1323 /* Exit immediately if the group is static */
1324 if (bStatic)
1326 m->bStatic = TRUE;
1327 return;
1330 if (bMaskOnly)
1332 m->nr = m->b.nr;
1333 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1335 /* Find the next atom in the block */
1336 while (m->b.a[j] != g->index[i])
1338 ++j;
1340 /* Mark blocks that did not contain any atoms */
1341 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1343 m->refid[bj++] = -1;
1345 /* Advance the block index if we have reached the next block */
1346 if (m->b.index[bj] <= j)
1348 ++bj;
1351 /* Mark the last blocks as not accessible */
1352 while (bj < m->b.nr)
1354 m->refid[bj++] = -1;
1357 else
1359 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1361 /* Find the next atom in the block */
1362 while (m->b.a[j] != g->index[i])
1364 ++j;
1366 /* If we have reached a new block, add it */
1367 if (m->b.index[bj+1] <= j)
1369 /* Skip any blocks in between */
1370 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1372 ++bj;
1374 m->refid[bi] = bj;
1375 m->mapid[bi] = m->orgid[bj];
1376 m->mapb.index[bi] = i;
1377 bi++;
1380 /* Update the number of blocks */
1381 m->mapb.index[bi] = g->isize;
1382 m->nr = bi;
1383 m->bMapStatic = FALSE;
1385 m->mapb.nr = m->nr;
1386 m->bStatic = FALSE;
1390 * \param[in,out] m Mapping structure to free.
1392 * All the memory allocated for the mapping structure is freed, and
1393 * the pointers set to NULL.
1394 * The pointer \p m is not freed.
1396 void
1397 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1399 sfree(m->refid);
1400 sfree(m->mapid);
1401 sfree(m->mapb.index);
1402 sfree(m->orgid);
1403 sfree(m->b.index);
1404 sfree(m->b.a);
1405 gmx_ana_indexmap_clear(m);