2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017,2018, 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
44 #include "indexutil.h"
53 #include "gromacs/topology/block.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/topology/mtop_lookup.h"
56 #include "gromacs/topology/mtop_util.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/utility/exceptions.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/smalloc.h"
61 #include "gromacs/utility/stringutil.h"
62 #include "gromacs/utility/textwriter.h"
64 /********************************************************************
65 * gmx_ana_indexgrps_t functions
66 ********************************************************************/
69 * Stores a set of index groups.
71 struct gmx_ana_indexgrps_t
73 //! Initializes an empty set of groups.
74 explicit gmx_ana_indexgrps_t(int nr) : nr(nr), g(nullptr)
79 ~gmx_ana_indexgrps_t()
81 for (int i = 0; i < nr; ++i)
83 gmx_ana_index_deinit(&g[i]);
88 /** Number of index groups. */
90 /** Array of index groups. */
93 std::vector<std::string> names;
97 * \param[out] g Index group structure.
98 * \param[in] top Topology structure.
99 * \param[in] fnm File name for the index file.
100 * Memory is automatically allocated.
102 * One or both of \p top or \p fnm can be NULL.
103 * If \p top is NULL, an index file is required and the groups are read
104 * from the file (uses Gromacs routine init_index()).
105 * If \p fnm is NULL, default groups are constructed based on the
106 * topology (uses Gromacs routine analyse()).
107 * If both are null, the index group structure is initialized empty.
110 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, gmx_mtop_t *top,
113 t_blocka *block = nullptr;
114 char **names = nullptr;
118 block = init_index(fnm, &names);
122 block = new_blocka();
123 // TODO: Propagate mtop further.
124 t_atoms atoms = gmx_mtop_global_atoms(top);
125 analyse(&atoms, block, &names, FALSE, FALSE);
130 *g = new gmx_ana_indexgrps_t(0);
136 *g = new gmx_ana_indexgrps_t(block->nr);
137 for (int i = 0; i < block->nr; ++i)
139 gmx_ana_index_t *grp = &(*g)->g[i];
141 grp->isize = block->index[i+1] - block->index[i];
142 snew(grp->index, grp->isize);
143 for (int j = 0; j < grp->isize; ++j)
145 grp->index[j] = block->a[block->index[i]+j];
147 grp->nalloc_index = grp->isize;
148 (*g)->names.emplace_back(names[i]);
153 for (int i = 0; i < block->nr; ++i)
162 for (int i = 0; i < block->nr; ++i)
172 * \param[in] g Index groups structure.
174 * The pointer \p g is invalid after the call.
177 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
183 * \param[out] g Index group structure.
184 * \returns true if \p g is empty, i.e., has 0 index groups.
187 gmx_ana_indexgrps_is_empty(gmx_ana_indexgrps_t *g)
193 * \param[in] g Index groups structure.
194 * \param[in] n Index group number to get.
195 * \returns Pointer to the \p n'th index group in \p g.
197 * The returned pointer should not be freed.
200 gmx_ana_indexgrps_get_grp(gmx_ana_indexgrps_t *g, int n)
202 if (n < 0 || n >= g->nr)
210 * \param[out] dest Output structure.
211 * \param[out] destName Receives the name of the group if found.
212 * \param[in] src Input index groups.
213 * \param[in] n Number of the group to extract.
214 * \returns true if \p n is a valid group in \p src, false otherwise.
217 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
218 gmx_ana_indexgrps_t *src, int n)
221 if (n < 0 || n >= src->nr)
227 if (destName != nullptr)
229 *destName = src->names[n];
231 gmx_ana_index_copy(dest, &src->g[n], true);
236 * \param[out] dest Output structure.
237 * \param[out] destName Receives the name of the group if found.
238 * \param[in] src Input index groups.
239 * \param[in] name Name (or part of the name) of the group to extract.
240 * \returns true if \p name is a valid group in \p src, false otherwise.
242 * Uses the Gromacs routine find_group() to find the actual group;
243 * the comparison is case-insensitive.
246 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
247 gmx_ana_indexgrps_t *src,
253 snew(names, src->nr);
254 for (int i = 0; i < src->nr; ++i)
256 names[i] = src->names[i].c_str();
258 int n = find_group(const_cast<char *>(name), src->nr,
259 const_cast<char **>(names));
267 return gmx_ana_indexgrps_extract(dest, destName, src, n);
271 * \param[in] writer Writer to use for output.
272 * \param[in] g Index groups to print.
273 * \param[in] maxn Maximum number of indices to print
274 * (-1 = print all, 0 = print only names).
277 gmx_ana_indexgrps_print(gmx::TextWriter *writer, gmx_ana_indexgrps_t *g, int maxn)
279 for (int i = 0; i < g->nr; ++i)
281 writer->writeString(gmx::formatString(" Group %2d \"%s\" ",
282 i, g->names[i].c_str()));
283 gmx_ana_index_dump(writer, &g->g[i], maxn);
287 /********************************************************************
288 * gmx_ana_index_t functions
289 ********************************************************************/
292 * \param[in,out] g Index group structure.
293 * \param[in] isize Maximum number of atoms to reserve space for.
296 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
298 if (g->nalloc_index < isize)
300 srenew(g->index, isize);
301 g->nalloc_index = isize;
306 * \param[in,out] g Index group structure.
308 * Resizes the memory allocated for holding the indices such that the
309 * current contents fit.
312 gmx_ana_index_squeeze(gmx_ana_index_t *g)
314 srenew(g->index, g->isize);
315 g->nalloc_index = g->isize;
319 * \param[out] g Output structure.
321 * Any contents of \p g are discarded without freeing.
324 gmx_ana_index_clear(gmx_ana_index_t *g)
332 * \param[out] g Output structure.
333 * \param[in] isize Number of atoms in the new group.
334 * \param[in] index Array of \p isize atoms (can be NULL if \p isize is 0).
335 * \param[in] nalloc Number of elements allocated for \p index
336 * (if 0, \p index is not freed in gmx_ana_index_deinit())
338 * No copy if \p index is made.
341 gmx_ana_index_set(gmx_ana_index_t *g, int isize, int *index, int nalloc)
345 g->nalloc_index = nalloc;
349 * \param[out] g Output structure.
350 * \param[in] natoms Number of atoms.
353 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms)
358 snew(g->index, natoms);
359 for (i = 0; i < natoms; ++i)
363 g->nalloc_index = natoms;
367 * \param[in] g Index group structure.
369 * The pointer \p g is not freed.
372 gmx_ana_index_deinit(gmx_ana_index_t *g)
374 if (g->nalloc_index > 0)
378 gmx_ana_index_clear(g);
382 * \param[out] dest Destination index group.
383 * \param[in] src Source index group.
384 * \param[in] bAlloc If true, memory is allocated at \p dest; otherwise,
385 * it is assumed that enough memory has been allocated for index.
388 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
390 dest->isize = src->isize;
393 snew(dest->index, dest->isize);
394 dest->nalloc_index = dest->isize;
398 std::memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
403 * \param[in] writer Writer to use for output.
404 * \param[in] g Index group to print.
405 * \param[in] maxn Maximum number of indices to print (-1 = print all).
408 gmx_ana_index_dump(gmx::TextWriter *writer, gmx_ana_index_t *g, int maxn)
410 writer->writeString(gmx::formatString("(%d atoms)", g->isize));
413 writer->writeString(":");
415 if (maxn >= 0 && n > maxn)
419 for (int j = 0; j < n; ++j)
421 writer->writeString(gmx::formatString(" %d", g->index[j]+1));
425 writer->writeString(" ...");
428 writer->ensureLineBreak();
432 gmx_ana_index_get_max_index(gmx_ana_index_t *g)
440 return *std::max_element(g->index, g->index + g->isize);
445 * \param[in] g Index group to check.
446 * \returns true if the index group is sorted and has no duplicates,
450 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
454 for (i = 0; i < g->isize-1; ++i)
456 if (g->index[i+1] <= g->index[i])
465 gmx_ana_index_check_range(gmx_ana_index_t *g, int natoms)
467 for (int i = 0; i < g->isize; ++i)
469 if (g->index[i] < 0 || g->index[i] >= natoms)
477 /********************************************************************
479 ********************************************************************/
482 * \param[in,out] g Index group to be sorted.
485 gmx_ana_index_sort(gmx_ana_index_t *g)
487 std::sort(g->index, g->index+g->isize);
491 gmx_ana_index_remove_duplicates(gmx_ana_index_t *g)
494 for (int i = 0; i < g->isize; ++i)
496 if (i == 0 || g->index[i-1] != g->index[i])
498 g->index[j] = g->index[i];
506 * \param[in] a Index group to check.
507 * \param[in] b Index group to check.
508 * \returns true if \p a and \p b are equal, false otherwise.
511 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
515 if (a->isize != b->isize)
519 for (i = 0; i < a->isize; ++i)
521 if (a->index[i] != b->index[i])
530 * \param[in] a Index group to check against.
531 * \param[in] b Index group to check.
532 * \returns true if \p b is contained in \p a,
535 * If the elements are not in the same order in both groups, the function
536 * fails. However, the groups do not need to be sorted.
539 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
543 for (i = j = 0; j < b->isize; ++i, ++j)
545 while (i < a->isize && a->index[i] != b->index[j])
558 * \param[out] dest Output index group (the intersection of \p a and \p b).
559 * \param[in] a First index group.
560 * \param[in] b Second index group.
562 * \p dest can be the same as \p a or \p b.
565 gmx_ana_index_intersection(gmx_ana_index_t *dest,
566 gmx_ana_index_t *a, gmx_ana_index_t *b)
570 for (i = j = k = 0; i < a->isize && j < b->isize; ++i)
572 while (j < b->isize && b->index[j] < a->index[i])
576 if (j < b->isize && b->index[j] == a->index[i])
578 dest->index[k++] = b->index[j++];
585 * \param[out] dest Output index group (the difference \p a - \p b).
586 * \param[in] a First index group.
587 * \param[in] b Second index group.
589 * \p dest can equal \p a, but not \p b.
592 gmx_ana_index_difference(gmx_ana_index_t *dest,
593 gmx_ana_index_t *a, gmx_ana_index_t *b)
597 for (i = j = k = 0; i < a->isize; ++i)
599 while (j < b->isize && b->index[j] < a->index[i])
603 if (j == b->isize || b->index[j] != a->index[i])
605 dest->index[k++] = a->index[i];
612 * \param[in] a First index group.
613 * \param[in] b Second index group.
614 * \returns Size of the difference \p a - \p b.
617 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
621 for (i = j = k = 0; i < a->isize; ++i)
623 while (j < b->isize && b->index[j] < a->index[i])
627 if (j == b->isize || b->index[j] != a->index[i])
636 * \param[out] dest1 Output group 1 (will equal \p g).
637 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
638 * \param[in] src Group to be partitioned.
639 * \param[in] g One partition.
641 * \pre \p g is a subset of \p src and both sets are sorted
642 * \pre \p dest1 has allocated storage to store \p src
643 * \post \p dest1 == \p g
644 * \post \p dest2 == \p src - \p g
646 * No storage should be allocated for \p dest2; after the call,
647 * \p dest2->index points to the memory allocated for \p dest1
648 * (to a part that is not used by \p dest1).
650 * The calculation can be performed in-place by setting \p dest1 equal to
654 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
655 gmx_ana_index_t *src, gmx_ana_index_t *g)
659 dest2->index = dest1->index + g->isize;
660 dest2->isize = src->isize - g->isize;
661 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
663 while (j >= 0 && src->index[j] != g->index[i])
665 dest2->index[k--] = src->index[j--];
670 dest2->index[k--] = src->index[j--];
672 gmx_ana_index_copy(dest1, g, false);
676 * \param[out] dest Output index group (the union of \p a and \p b).
677 * \param[in] a First index group.
678 * \param[in] b Second index group.
680 * \p a and \p b can have common items.
681 * \p dest can equal \p a or \p b.
683 * \see gmx_ana_index_merge()
686 gmx_ana_index_union(gmx_ana_index_t *dest,
687 gmx_ana_index_t *a, gmx_ana_index_t *b)
692 dsize = gmx_ana_index_difference_size(b, a);
695 dest->isize = a->isize + dsize;
696 for (k = dest->isize - 1; k >= 0; k--)
698 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
700 dest->index[k] = b->index[j--];
704 if (j >= 0 && a->index[i] == b->index[j])
708 dest->index[k] = a->index[i--];
714 gmx_ana_index_union_unsorted(gmx_ana_index_t *dest,
715 gmx_ana_index_t *a, gmx_ana_index_t *b)
717 if (gmx_ana_index_check_sorted(b))
719 gmx_ana_index_union(dest, a, b);
724 gmx_ana_index_copy(&tmp, b, true);
725 gmx_ana_index_sort(&tmp);
726 gmx_ana_index_remove_duplicates(&tmp);
727 gmx_ana_index_union(dest, a, &tmp);
728 gmx_ana_index_deinit(&tmp);
733 * \param[out] dest Output index group (the union of \p a and \p b).
734 * \param[in] a First index group.
735 * \param[in] b Second index group.
737 * \p a and \p b should not have common items.
738 * \p dest can equal \p a or \p b.
740 * \see gmx_ana_index_union()
743 gmx_ana_index_merge(gmx_ana_index_t *dest,
744 gmx_ana_index_t *a, gmx_ana_index_t *b)
750 dest->isize = a->isize + b->isize;
751 for (k = dest->isize - 1; k >= 0; k--)
753 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
755 dest->index[k] = b->index[j--];
759 dest->index[k] = a->index[i--];
764 /********************************************************************
765 * gmx_ana_indexmap_t and related things
766 ********************************************************************/
769 * Helper for splitting a sequence of atom indices into groups.
771 * \param[in] atomIndex Index of the next atom in the sequence.
772 * \param[in] top Topology structure.
773 * \param[in] type Type of group to split into.
774 * \param[in,out] id Variable to receive the next group id.
775 * \returns `true` if \p atomIndex starts a new group in the sequence, i.e.,
776 * if \p *id was changed.
778 * \p *id should be initialized to `-1` before first call of this function, and
779 * then each atom index in the sequence passed to the function in turn.
781 * \ingroup module_selection
784 next_group_index(int atomIndex, const gmx_mtop_t *top,
785 e_index_t type, int *id)
795 int resind, molb = 0;
796 mtopGetAtomAndResidueName(top, atomIndex, &molb, nullptr, nullptr, nullptr, &resind);
803 *id = mtopGetMoleculeIndex(top, atomIndex, &molb);
815 * \param[in,out] t Output block.
816 * \param[in] top Topology structure
817 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
819 * \param[in] g Index group
820 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
821 * \param[in] type Type of partitioning to make.
822 * \param[in] bComplete
823 * If true, the index group is expanded to include any residue/molecule
824 * (depending on \p type) that is partially contained in the group.
825 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
827 * \p m should have been initialized somehow (calloc() is enough).
828 * \p g should be sorted.
831 gmx_ana_index_make_block(t_blocka *t, const gmx_mtop_t *top, gmx_ana_index_t *g,
832 e_index_t type, bool bComplete)
834 if (type == INDEX_UNKNOWN)
848 // TODO: Check callers and either check these there as well, or turn these
850 GMX_RELEASE_ASSERT(top != nullptr || (type != INDEX_RES && type != INDEX_MOL),
851 "Topology must be provided for residue or molecule blocks");
852 GMX_RELEASE_ASSERT(type != INDEX_MOL || top->haveMoleculeIndices,
853 "Molecule information must be present for molecule blocks");
855 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
857 if (type != INDEX_RES && type != INDEX_MOL)
861 /* Allocate memory for the atom array and fill it unless we are using
866 /* We may allocate some extra memory here because we don't know in
867 * advance how much will be needed. */
868 if (t->nalloc_a < top->natoms)
870 srenew(t->a, top->natoms);
871 t->nalloc_a = top->natoms;
877 if (t->nalloc_a < g->isize)
879 srenew(t->a, g->isize);
880 t->nalloc_a = g->isize;
882 std::memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
885 /* Allocate memory for the block index. We don't know in advance
886 * how much will be needed, so we allocate some extra and free it in the
888 if (t->nalloc_index < g->isize + 1)
890 srenew(t->index, g->isize + 1);
891 t->nalloc_index = g->isize + 1;
897 for (int i = 0; i < g->isize; ++i)
899 const int ai = g->index[i];
900 /* Find the ID number of the atom/residue/molecule corresponding to
902 if (next_group_index(ai, top, type, &id))
904 /* If this is the first atom in a new block, initialize the block. */
907 /* For completion, we first set the start of the block. */
908 t->index[t->nr++] = t->nra;
909 /* And then we find all the atoms that should be included. */
915 mtopGetMolblockIndex(top, ai, &molb, &molnr, &atnr_mol);
916 const t_atoms &mol_atoms = top->moltype[top->molblock[molb].type].atoms;
917 int last_atom = atnr_mol + 1;
918 while (last_atom < mol_atoms.nr
919 && mol_atoms.atom[last_atom].resind == id)
923 int first_atom = atnr_mol - 1;
924 while (first_atom >= 0
925 && mol_atoms.atom[first_atom].resind == id)
929 int first_mol_atom = top->moleculeBlockIndices[molb].globalAtomStart;
930 first_mol_atom += molnr*top->moleculeBlockIndices[molb].numAtomsPerMolecule;
931 first_atom = first_mol_atom + first_atom + 1;
932 last_atom = first_mol_atom + last_atom - 1;
933 for (int j = first_atom; j <= last_atom; ++j)
942 while (molb + 1 < top->molblock.size() && id >= top->moleculeBlockIndices[molb].moleculeIndexStart)
946 const MoleculeBlockIndices &blockIndices = top->moleculeBlockIndices[molb];
947 const int atomStart = blockIndices.globalAtomStart + (id - blockIndices.moleculeIndexStart)*blockIndices.numAtomsPerMolecule;
948 for (int j = 0; j < blockIndices.numAtomsPerMolecule; ++j)
950 t->a[t->nra++] = atomStart + j;
954 default: /* Should not be reached */
955 GMX_RELEASE_ASSERT(false, "Unreachable code was reached");
961 /* If not using completion, simply store the start of the block. */
962 t->index[t->nr++] = i;
966 /* Set the end of the last block */
967 t->index[t->nr] = t->nra;
968 /* Free any unnecessary memory */
969 srenew(t->index, t->nr+1);
970 t->nalloc_index = t->nr+1;
973 srenew(t->a, t->nra);
974 t->nalloc_a = t->nra;
979 * \param[in] g Index group to check.
980 * \param[in] b Block data to check against.
981 * \returns true if \p g consists of one or more complete blocks from \p b,
984 * The atoms in \p g are assumed to be sorted.
987 gmx_ana_index_has_full_blocks(const gmx_ana_index_t *g,
988 const gmx::RangePartitioning *b)
993 /* Each round in the loop matches one block */
996 /* Find the block that begins with the first unmatched atom */
997 while (bi < b->numBlocks() && b->block(bi).begin() != g->index[i])
1001 /* If not found, or if too large, return */
1002 if (bi == b->numBlocks() || i + b->block(bi).size() > g->isize)
1006 /* Check that the block matches the index */
1007 for (j = b->block(bi).begin(); j < b->block(bi).end(); ++j, ++i)
1009 if (g->index[i] != j)
1014 /* Move the search to the next block */
1021 * \param[in] g Index group to check.
1022 * \param[in] b Block data to check against.
1023 * \returns true if \p g consists of one or more complete blocks from \p b,
1026 * The atoms in \p g and \p b->a are assumed to be in the same order.
1029 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
1034 /* Each round in the loop matches one block */
1035 while (i < g->isize)
1037 /* Find the block that begins with the first unmatched atom */
1038 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
1042 /* If not found, or if too large, return */
1043 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
1047 /* Check that the block matches the index */
1048 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
1050 if (b->a[j] != g->index[i])
1055 /* Move the search to the next block */
1062 * \brief Returns if an atom is at a residue boundary.
1064 * \param[in] top Topology data.
1065 * \param[in] a Atom index to check, should be -1 <= \p a < top->natoms.
1066 * \param[in,out] molb The molecule block of atom a
1067 * \returns true if atoms \p a and \p a + 1 are in different residues, false otherwise.
1069 static bool is_at_residue_boundary(const gmx_mtop_t *top, int a, int *molb)
1071 if (a == -1 || a + 1 == top->natoms)
1076 mtopGetAtomAndResidueName(top, a, molb,
1077 nullptr, nullptr, nullptr, &resindA);
1079 mtopGetAtomAndResidueName(top, a + 1, molb,
1080 nullptr, nullptr, nullptr, &resindAPlusOne);
1081 return resindAPlusOne != resindA;
1085 * \param[in] g Index group to check.
1086 * \param[in] type Block data to check against.
1087 * \param[in] top Topology data.
1088 * \returns true if \p g consists of one or more complete elements of type
1089 * \p type, false otherwise.
1091 * \p g is assumed to be sorted, otherwise may return false negatives.
1093 * If \p type is \ref INDEX_ATOM, the return value is always true.
1094 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
1098 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
1099 const gmx_mtop_t *top)
1106 // TODO: Consider whether unsorted groups need to be supported better.
1120 for (int i = 0; i < g->isize; ++i)
1122 const int a = g->index[i];
1123 // Check if a is consecutive or on a residue boundary
1126 if (!is_at_residue_boundary(top, aPrev, &molb))
1130 if (!is_at_residue_boundary(top, a - 1, &molb))
1137 GMX_ASSERT(g->isize > 0, "We return above when isize=0");
1138 const int a = g->index[g->isize - 1];
1139 if (!is_at_residue_boundary(top, a, &molb))
1148 auto molecules = gmx_mtop_molecules(*top);
1149 return gmx_ana_index_has_full_blocks(g, &molecules);
1156 * \param[out] m Output structure.
1158 * Any contents of \p m are discarded without freeing.
1161 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1163 m->type = INDEX_UNKNOWN;
1167 m->mapb.index = nullptr;
1168 m->mapb.nalloc_index = 0;
1170 m->mapb.a = nullptr;
1171 m->mapb.nalloc_a = 0;
1174 m->b.index = nullptr;
1177 m->b.nalloc_index = 0;
1183 * \param[in,out] m Mapping structure.
1184 * \param[in] nr Maximum number of blocks to reserve space for.
1185 * \param[in] isize Maximum number of atoms to reserve space for.
1188 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1190 if (m->mapb.nalloc_index < nr + 1)
1192 srenew(m->refid, nr);
1193 srenew(m->mapid, nr);
1194 srenew(m->orgid, nr);
1195 srenew(m->mapb.index, nr + 1);
1196 srenew(m->b.index, nr + 1);
1197 m->mapb.nalloc_index = nr + 1;
1198 m->b.nalloc_index = nr + 1;
1200 if (m->b.nalloc_a < isize)
1202 srenew(m->b.a, isize);
1203 m->b.nalloc_a = isize;
1208 * \param[in,out] m Mapping structure to initialize.
1209 * \param[in] g Index group to map
1210 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1211 * \param[in] top Topology structure
1212 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1213 * \param[in] type Type of mapping to construct.
1215 * Initializes a new index group mapping.
1216 * The index group provided to gmx_ana_indexmap_update() should always be a
1217 * subset of the \p g given here.
1219 * \p m should have been initialized somehow (calloc() is enough).
1222 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1223 const gmx_mtop_t *top, e_index_t type)
1226 gmx_ana_index_make_block(&m->b, top, g, type, false);
1227 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1229 for (int i = 0; i < m->b.nr; ++i)
1231 const int ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1232 next_group_index(ii, top, type, &id);
1237 m->mapb.nr = m->b.nr;
1238 m->mapb.nra = m->b.nra;
1240 std::memcpy(m->mapb.index, m->b.index, (m->b.nr+1)*sizeof(*(m->mapb.index)));
1245 gmx_ana_indexmap_init_orgid_group(gmx_ana_indexmap_t *m, const gmx_mtop_t *top,
1248 GMX_RELEASE_ASSERT(m->bStatic,
1249 "Changing original IDs is not supported after starting "
1250 "to use the mapping");
1251 GMX_RELEASE_ASSERT(top != nullptr || (type != INDEX_RES && type != INDEX_MOL),
1252 "Topology must be provided for residue or molecule blocks");
1253 // Check that all atoms in each block belong to the same group.
1254 // This is a separate loop for better error handling (no state is modified
1255 // if there is an error.
1256 if (type == INDEX_RES || type == INDEX_MOL)
1259 for (int i = 0; i < m->b.nr; ++i)
1261 const int ii = m->b.a[m->b.index[i]];
1262 if (next_group_index(ii, top, type, &id))
1264 for (int j = m->b.index[i] + 1; j < m->b.index[i+1]; ++j)
1266 if (next_group_index(m->b.a[j], top, type, &id))
1268 std::string message("Grouping into residues/molecules is ambiguous");
1269 GMX_THROW(gmx::InconsistentInputError(message));
1275 // Do a second loop, where things are actually set.
1278 for (int i = 0; i < m->b.nr; ++i)
1280 const int ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1281 if (next_group_index(ii, top, type, &id))
1285 m->mapid[i] = group;
1286 m->orgid[i] = group;
1288 // Count also the last group.
1294 * \param[in,out] m Mapping structure to initialize.
1295 * \param[in] b Block information to use for data.
1297 * Frees some memory that is not necessary for static index group mappings.
1298 * Internal pointers are set to point to data in \p b; it is the responsibility
1299 * of the caller to ensure that the block information matches the contents of
1301 * After this function has been called, the index group provided to
1302 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1304 * This function breaks modularity of the index group mapping interface in an
1305 * ugly way, but allows reducing memory usage of static selections by a
1306 * significant amount.
1309 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1312 sfree(m->mapb.index);
1315 m->mapb.nalloc_index = 0;
1316 m->mapb.nalloc_a = 0;
1317 m->b.nalloc_index = 0;
1319 m->mapid = m->orgid;
1320 m->mapb.index = b->index;
1322 m->b.index = b->index;
1327 * \param[in,out] dest Destination data structure.
1328 * \param[in] src Source mapping.
1329 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1330 * copy is made; otherwise, only variable parts are copied, and no memory
1333 * \p dest should have been initialized somehow (calloc() is enough).
1336 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1340 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1341 dest->type = src->type;
1342 dest->b.nr = src->b.nr;
1343 dest->b.nra = src->b.nra;
1344 std::memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1345 std::memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1346 std::memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1348 dest->mapb.nr = src->mapb.nr;
1349 dest->mapb.nra = src->mapb.nra;
1350 if (src->mapb.nalloc_a > 0)
1354 snew(dest->mapb.a, src->mapb.nalloc_a);
1355 dest->mapb.nalloc_a = src->mapb.nalloc_a;
1357 std::memcpy(dest->mapb.a, src->mapb.a, dest->mapb.nra*sizeof(*dest->mapb.a));
1361 dest->mapb.a = src->mapb.a;
1363 std::memcpy(dest->refid, src->refid, dest->mapb.nr*sizeof(*dest->refid));
1364 std::memcpy(dest->mapid, src->mapid, dest->mapb.nr*sizeof(*dest->mapid));
1365 std::memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1366 dest->bStatic = src->bStatic;
1370 * Helper function to set the source atoms in an index map.
1372 * \param[in,out] m Mapping structure.
1373 * \param[in] isize Number of atoms in the \p index array.
1374 * \param[in] index List of atoms.
1377 set_atoms(gmx_ana_indexmap_t *m, int isize, int *index)
1379 m->mapb.nra = isize;
1380 if (m->mapb.nalloc_a == 0)
1386 for (int i = 0; i < isize; ++i)
1388 m->mapb.a[i] = index[i];
1394 * \param[in,out] m Mapping structure.
1395 * \param[in] g Current index group.
1396 * \param[in] bMaskOnly true if the unused blocks should be masked with
1397 * -1 instead of removing them.
1399 * Updates the index group mapping with the new index group \p g.
1401 * \see gmx_ana_indexmap_t
1404 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1409 /* Process the simple cases first */
1410 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1414 if (m->type == INDEX_ALL)
1416 set_atoms(m, g->isize, g->index);
1419 m->mapb.index[1] = g->isize;
1423 /* Reset the reference IDs and mapping if necessary */
1424 const bool bToFull = (g->isize == m->b.nra);
1425 const bool bWasFull = (m->mapb.nra == m->b.nra);
1426 if (bToFull || bMaskOnly)
1430 for (bj = 0; bj < m->b.nr; ++bj)
1437 for (bj = 0; bj < m->b.nr; ++bj)
1439 m->mapid[bj] = m->orgid[bj];
1441 for (bj = 0; bj <= m->b.nr; ++bj)
1443 m->mapb.index[bj] = m->b.index[bj];
1446 set_atoms(m, m->b.nra, m->b.a);
1447 m->mapb.nr = m->b.nr;
1449 /* Exit immediately if the group is static */
1458 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1460 /* Find the next atom in the block */
1461 while (m->b.a[j] != g->index[i])
1465 /* Mark blocks that did not contain any atoms */
1466 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1468 m->refid[bj++] = -1;
1470 /* Advance the block index if we have reached the next block */
1471 if (m->b.index[bj] <= j)
1476 /* Mark the last blocks as not accessible */
1477 while (bj < m->b.nr)
1479 m->refid[bj++] = -1;
1484 set_atoms(m, g->isize, g->index);
1485 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1487 /* Find the next atom in the block */
1488 while (m->b.a[j] != g->index[i])
1492 /* If we have reached a new block, add it */
1493 if (m->b.index[bj+1] <= j)
1495 /* Skip any blocks in between */
1496 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1501 m->mapid[bi] = m->orgid[bj];
1502 m->mapb.index[bi] = i;
1506 /* Update the number of blocks */
1507 m->mapb.index[bi] = g->isize;
1514 * \param[in,out] m Mapping structure to free.
1516 * All the memory allocated for the mapping structure is freed, and
1517 * the pointers set to NULL.
1518 * The pointer \p m is not freed.
1521 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1524 if (m->mapid != m->orgid)
1528 if (m->mapb.nalloc_index > 0)
1530 sfree(m->mapb.index);
1532 if (m->mapb.nalloc_a > 0)
1537 if (m->b.nalloc_index > 0)
1541 if (m->b.nalloc_a > 0)
1545 gmx_ana_indexmap_clear(m);