2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017,2018,2019, 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)
79 ~gmx_ana_indexgrps_t()
81 for (auto &indexGrp : g)
83 gmx_ana_index_deinit(&indexGrp);
88 /** Array of index groups. */
89 std::vector<gmx_ana_index_t> g;
91 std::vector<std::string> names;
95 * \param[out] g Index group structure.
96 * \param[in] top Topology structure.
97 * \param[in] fnm File name for the index file.
98 * Memory is automatically allocated.
100 * One or both of \p top or \p fnm can be NULL.
101 * If \p top is NULL, an index file is required and the groups are read
102 * from the file (uses Gromacs routine init_index()).
103 * If \p fnm is NULL, default groups are constructed based on the
104 * topology (uses Gromacs routine analyse()).
105 * If both are null, the index group structure is initialized empty.
108 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t **g, gmx_mtop_t *top,
111 t_blocka *block = nullptr;
112 char **names = nullptr;
116 block = init_index(fnm, &names);
120 block = new_blocka();
121 // TODO: Propagate mtop further.
122 t_atoms atoms = gmx_mtop_global_atoms(top);
123 analyse(&atoms, block, &names, FALSE, FALSE);
128 *g = new gmx_ana_indexgrps_t(0);
134 *g = new gmx_ana_indexgrps_t(block->nr);
135 for (int i = 0; i < block->nr; ++i)
137 gmx_ana_index_t *grp = &(*g)->g[i];
139 grp->isize = block->index[i+1] - block->index[i];
140 snew(grp->index, grp->isize);
141 for (int j = 0; j < grp->isize; ++j)
143 grp->index[j] = block->a[block->index[i]+j];
145 grp->nalloc_index = grp->isize;
146 (*g)->names.emplace_back(names[i]);
151 for (int i = 0; i < block->nr; ++i)
160 for (int i = 0; i < block->nr; ++i)
170 * \param[in] g Index groups structure.
172 * The pointer \p g is invalid after the call.
175 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t *g)
182 * \param[out] dest Output structure.
183 * \param[out] destName Receives the name of the group if found.
184 * \param[in] src Input index groups.
185 * \param[in] n Number of the group to extract.
186 * \returns true if \p n is a valid group in \p src, false otherwise.
189 gmx_ana_indexgrps_extract(gmx_ana_index_t *dest, std::string *destName,
190 gmx_ana_indexgrps_t *src, int n)
193 if (n < 0 || n >= gmx::index(src->g.size()))
199 if (destName != nullptr)
201 *destName = src->names[n];
203 gmx_ana_index_copy(dest, &src->g[n], true);
208 * \param[out] dest Output structure.
209 * \param[out] destName Receives the name of the group if found.
210 * \param[in] src Input index groups.
211 * \param[in] name Name (or part of the name) of the group to extract.
212 * \returns true if \p name is a valid group in \p src, false otherwise.
214 * Uses the Gromacs routine find_group() to find the actual group;
215 * the comparison is case-insensitive.
218 gmx_ana_indexgrps_find(gmx_ana_index_t *dest, std::string *destName,
219 gmx_ana_indexgrps_t *src,
225 snew(names, src->g.size());
226 for (size_t i = 0; i < src->g.size(); ++i)
228 names[i] = src->names[i].c_str();
230 int n = find_group(const_cast<char *>(name), src->g.size(),
231 const_cast<char **>(names));
239 return gmx_ana_indexgrps_extract(dest, destName, src, n);
243 * \param[in] writer Writer to use for output.
244 * \param[in] g Index groups to print.
245 * \param[in] maxn Maximum number of indices to print
246 * (-1 = print all, 0 = print only names).
249 gmx_ana_indexgrps_print(gmx::TextWriter *writer, gmx_ana_indexgrps_t *g, int maxn)
251 for (int i = 0; i < gmx::ssize(g->g); ++i)
253 writer->writeString(gmx::formatString(" Group %2d \"%s\" ",
254 i, g->names[i].c_str()));
255 gmx_ana_index_dump(writer, &g->g[i], maxn);
259 /********************************************************************
260 * gmx_ana_index_t functions
261 ********************************************************************/
264 * \param[in,out] g Index group structure.
265 * \param[in] isize Maximum number of atoms to reserve space for.
268 gmx_ana_index_reserve(gmx_ana_index_t *g, int isize)
270 if (g->nalloc_index < isize)
272 srenew(g->index, isize);
273 g->nalloc_index = isize;
278 * \param[in,out] g Index group structure.
280 * Resizes the memory allocated for holding the indices such that the
281 * current contents fit.
284 gmx_ana_index_squeeze(gmx_ana_index_t *g)
286 srenew(g->index, g->isize);
287 g->nalloc_index = g->isize;
291 * \param[out] g Output structure.
293 * Any contents of \p g are discarded without freeing.
296 gmx_ana_index_clear(gmx_ana_index_t *g)
304 * \param[out] g Output structure.
305 * \param[in] isize Number of atoms in the new group.
306 * \param[in] index Array of \p isize atoms (can be NULL if \p isize is 0).
307 * \param[in] nalloc Number of elements allocated for \p index
308 * (if 0, \p index is not freed in gmx_ana_index_deinit())
310 * No copy if \p index is made.
313 gmx_ana_index_set(gmx_ana_index_t *g, int isize, int *index, int nalloc)
317 g->nalloc_index = nalloc;
321 * \param[out] g Output structure.
322 * \param[in] natoms Number of atoms.
325 gmx_ana_index_init_simple(gmx_ana_index_t *g, int natoms)
330 snew(g->index, natoms);
331 for (i = 0; i < natoms; ++i)
335 g->nalloc_index = natoms;
339 * \param[in] g Index group structure.
341 * The pointer \p g is not freed.
344 gmx_ana_index_deinit(gmx_ana_index_t *g)
346 if (g->nalloc_index > 0)
350 gmx_ana_index_clear(g);
354 * \param[out] dest Destination index group.
355 * \param[in] src Source index group.
356 * \param[in] bAlloc If true, memory is allocated at \p dest; otherwise,
357 * it is assumed that enough memory has been allocated for index.
360 gmx_ana_index_copy(gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc)
362 dest->isize = src->isize;
365 snew(dest->index, dest->isize);
366 dest->nalloc_index = dest->isize;
370 std::memcpy(dest->index, src->index, dest->isize*sizeof(*dest->index));
375 * \param[in] writer Writer to use for output.
376 * \param[in] g Index group to print.
377 * \param[in] maxn Maximum number of indices to print (-1 = print all).
380 gmx_ana_index_dump(gmx::TextWriter *writer, gmx_ana_index_t *g, int maxn)
382 writer->writeString(gmx::formatString("(%d atoms)", g->isize));
385 writer->writeString(":");
387 if (maxn >= 0 && n > maxn)
391 for (int j = 0; j < n; ++j)
393 writer->writeString(gmx::formatString(" %d", g->index[j]+1));
397 writer->writeString(" ...");
400 writer->ensureLineBreak();
404 gmx_ana_index_get_max_index(gmx_ana_index_t *g)
412 return *std::max_element(g->index, g->index + g->isize);
417 * \param[in] g Index group to check.
418 * \returns true if the index group is sorted and has no duplicates,
422 gmx_ana_index_check_sorted(gmx_ana_index_t *g)
426 for (i = 0; i < g->isize-1; ++i)
428 if (g->index[i+1] <= g->index[i])
437 gmx_ana_index_check_range(gmx_ana_index_t *g, int natoms)
439 for (int i = 0; i < g->isize; ++i)
441 if (g->index[i] < 0 || g->index[i] >= natoms)
449 /********************************************************************
451 ********************************************************************/
454 * \param[in,out] g Index group to be sorted.
457 gmx_ana_index_sort(gmx_ana_index_t *g)
459 std::sort(g->index, g->index+g->isize);
463 gmx_ana_index_remove_duplicates(gmx_ana_index_t *g)
466 for (int i = 0; i < g->isize; ++i)
468 if (i == 0 || g->index[i-1] != g->index[i])
470 g->index[j] = g->index[i];
478 * \param[in] a Index group to check.
479 * \param[in] b Index group to check.
480 * \returns true if \p a and \p b are equal, false otherwise.
483 gmx_ana_index_equals(gmx_ana_index_t *a, gmx_ana_index_t *b)
487 if (a->isize != b->isize)
491 for (i = 0; i < a->isize; ++i)
493 if (a->index[i] != b->index[i])
502 * \param[in] a Index group to check against.
503 * \param[in] b Index group to check.
504 * \returns true if \p b is contained in \p a,
507 * If the elements are not in the same order in both groups, the function
508 * fails. However, the groups do not need to be sorted.
511 gmx_ana_index_contains(gmx_ana_index_t *a, gmx_ana_index_t *b)
515 for (i = j = 0; j < b->isize; ++i, ++j)
517 while (i < a->isize && a->index[i] != b->index[j])
530 * \param[out] dest Output index group (the intersection of \p a and \p b).
531 * \param[in] a First index group.
532 * \param[in] b Second index group.
534 * \p dest can be the same as \p a or \p b.
537 gmx_ana_index_intersection(gmx_ana_index_t *dest,
538 gmx_ana_index_t *a, gmx_ana_index_t *b)
542 for (i = j = k = 0; i < a->isize && j < b->isize; ++i)
544 while (j < b->isize && b->index[j] < a->index[i])
548 if (j < b->isize && b->index[j] == a->index[i])
550 dest->index[k++] = b->index[j++];
557 * \param[out] dest Output index group (the difference \p a - \p b).
558 * \param[in] a First index group.
559 * \param[in] b Second index group.
561 * \p dest can equal \p a, but not \p b.
564 gmx_ana_index_difference(gmx_ana_index_t *dest,
565 gmx_ana_index_t *a, gmx_ana_index_t *b)
569 for (i = j = k = 0; i < a->isize; ++i)
571 while (j < b->isize && b->index[j] < a->index[i])
575 if (j == b->isize || b->index[j] != a->index[i])
577 dest->index[k++] = a->index[i];
584 * \param[in] a First index group.
585 * \param[in] b Second index group.
586 * \returns Size of the difference \p a - \p b.
589 gmx_ana_index_difference_size(gmx_ana_index_t *a, gmx_ana_index_t *b)
593 for (i = j = k = 0; i < a->isize; ++i)
595 while (j < b->isize && b->index[j] < a->index[i])
599 if (j == b->isize || b->index[j] != a->index[i])
608 * \param[out] dest1 Output group 1 (will equal \p g).
609 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
610 * \param[in] src Group to be partitioned.
611 * \param[in] g One partition.
613 * \pre \p g is a subset of \p src and both sets are sorted
614 * \pre \p dest1 has allocated storage to store \p src
615 * \post \p dest1 == \p g
616 * \post \p dest2 == \p src - \p g
618 * No storage should be allocated for \p dest2; after the call,
619 * \p dest2->index points to the memory allocated for \p dest1
620 * (to a part that is not used by \p dest1).
622 * The calculation can be performed in-place by setting \p dest1 equal to
626 gmx_ana_index_partition(gmx_ana_index_t *dest1, gmx_ana_index_t *dest2,
627 gmx_ana_index_t *src, gmx_ana_index_t *g)
631 dest2->index = dest1->index + g->isize;
632 dest2->isize = src->isize - g->isize;
633 for (i = g->isize-1, j = src->isize-1, k = dest2->isize-1; i >= 0; --i, --j)
635 while (j >= 0 && src->index[j] != g->index[i])
637 dest2->index[k--] = src->index[j--];
642 dest2->index[k--] = src->index[j--];
644 gmx_ana_index_copy(dest1, g, false);
648 * \param[out] dest Output index group (the union of \p a and \p b).
649 * \param[in] a First index group.
650 * \param[in] b Second index group.
652 * \p a and \p b can have common items.
653 * \p dest can equal \p a or \p b.
655 * \see gmx_ana_index_merge()
658 gmx_ana_index_union(gmx_ana_index_t *dest,
659 gmx_ana_index_t *a, gmx_ana_index_t *b)
664 dsize = gmx_ana_index_difference_size(b, a);
667 dest->isize = a->isize + dsize;
668 for (k = dest->isize - 1; k >= 0; k--)
670 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
672 dest->index[k] = b->index[j--];
676 if (j >= 0 && a->index[i] == b->index[j])
680 dest->index[k] = a->index[i--];
686 gmx_ana_index_union_unsorted(gmx_ana_index_t *dest,
687 gmx_ana_index_t *a, gmx_ana_index_t *b)
689 if (gmx_ana_index_check_sorted(b))
691 gmx_ana_index_union(dest, a, b);
696 gmx_ana_index_copy(&tmp, b, true);
697 gmx_ana_index_sort(&tmp);
698 gmx_ana_index_remove_duplicates(&tmp);
699 gmx_ana_index_union(dest, a, &tmp);
700 gmx_ana_index_deinit(&tmp);
705 * \param[out] dest Output index group (the union of \p a and \p b).
706 * \param[in] a First index group.
707 * \param[in] b Second index group.
709 * \p a and \p b should not have common items.
710 * \p dest can equal \p a or \p b.
712 * \see gmx_ana_index_union()
715 gmx_ana_index_merge(gmx_ana_index_t *dest,
716 gmx_ana_index_t *a, gmx_ana_index_t *b)
722 dest->isize = a->isize + b->isize;
723 for (k = dest->isize - 1; k >= 0; k--)
725 if (i < 0 || (j >= 0 && a->index[i] < b->index[j]))
727 dest->index[k] = b->index[j--];
731 dest->index[k] = a->index[i--];
736 /********************************************************************
737 * gmx_ana_indexmap_t and related things
738 ********************************************************************/
741 * Helper for splitting a sequence of atom indices into groups.
743 * \param[in] atomIndex Index of the next atom in the sequence.
744 * \param[in] top Topology structure.
745 * \param[in] type Type of group to split into.
746 * \param[in,out] id Variable to receive the next group id.
747 * \returns `true` if \p atomIndex starts a new group in the sequence, i.e.,
748 * if \p *id was changed.
750 * \p *id should be initialized to `-1` before first call of this function, and
751 * then each atom index in the sequence passed to the function in turn.
753 * \ingroup module_selection
756 next_group_index(int atomIndex, const gmx_mtop_t *top,
757 e_index_t type, int *id)
767 int resind, molb = 0;
768 mtopGetAtomAndResidueName(top, atomIndex, &molb, nullptr, nullptr, nullptr, &resind);
775 *id = mtopGetMoleculeIndex(top, atomIndex, &molb);
787 * \param[in,out] t Output block.
788 * \param[in] top Topology structure
789 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
791 * \param[in] g Index group
792 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
793 * \param[in] type Type of partitioning to make.
794 * \param[in] bComplete
795 * If true, the index group is expanded to include any residue/molecule
796 * (depending on \p type) that is partially contained in the group.
797 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
799 * \p m should have been initialized somehow (calloc() is enough).
800 * \p g should be sorted.
803 gmx_ana_index_make_block(t_blocka *t, const gmx_mtop_t *top, gmx_ana_index_t *g,
804 e_index_t type, bool bComplete)
806 if (type == INDEX_UNKNOWN)
820 // TODO: Check callers and either check these there as well, or turn these
822 GMX_RELEASE_ASSERT(top != nullptr || (type != INDEX_RES && type != INDEX_MOL),
823 "Topology must be provided for residue or molecule blocks");
824 GMX_RELEASE_ASSERT(type != INDEX_MOL || top->haveMoleculeIndices,
825 "Molecule information must be present for molecule blocks");
827 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
829 if (type != INDEX_RES && type != INDEX_MOL)
833 /* Allocate memory for the atom array and fill it unless we are using
838 /* We may allocate some extra memory here because we don't know in
839 * advance how much will be needed. */
840 if (t->nalloc_a < top->natoms)
842 srenew(t->a, top->natoms);
843 t->nalloc_a = top->natoms;
849 if (t->nalloc_a < g->isize)
851 srenew(t->a, g->isize);
852 t->nalloc_a = g->isize;
854 std::memcpy(t->a, g->index, g->isize*sizeof(*(t->a)));
857 /* Allocate memory for the block index. We don't know in advance
858 * how much will be needed, so we allocate some extra and free it in the
860 if (t->nalloc_index < g->isize + 1)
862 srenew(t->index, g->isize + 1);
863 t->nalloc_index = g->isize + 1;
869 for (int i = 0; i < g->isize; ++i)
871 const int ai = g->index[i];
872 /* Find the ID number of the atom/residue/molecule corresponding to
874 if (next_group_index(ai, top, type, &id))
876 /* If this is the first atom in a new block, initialize the block. */
879 /* For completion, we first set the start of the block. */
880 t->index[t->nr++] = t->nra;
881 /* And then we find all the atoms that should be included. */
887 mtopGetMolblockIndex(top, ai, &molb, &molnr, &atnr_mol);
888 const t_atoms &mol_atoms = top->moltype[top->molblock[molb].type].atoms;
889 int last_atom = atnr_mol + 1;
890 while (last_atom < mol_atoms.nr
891 && mol_atoms.atom[last_atom].resind == id)
895 int first_atom = atnr_mol - 1;
896 while (first_atom >= 0
897 && mol_atoms.atom[first_atom].resind == id)
901 int first_mol_atom = top->moleculeBlockIndices[molb].globalAtomStart;
902 first_mol_atom += molnr*top->moleculeBlockIndices[molb].numAtomsPerMolecule;
903 first_atom = first_mol_atom + first_atom + 1;
904 last_atom = first_mol_atom + last_atom - 1;
905 for (int j = first_atom; j <= last_atom; ++j)
914 while (molb + 1 < top->molblock.size() && id >= top->moleculeBlockIndices[molb].moleculeIndexStart)
918 const MoleculeBlockIndices &blockIndices = top->moleculeBlockIndices[molb];
919 const int atomStart = blockIndices.globalAtomStart + (id - blockIndices.moleculeIndexStart)*blockIndices.numAtomsPerMolecule;
920 for (int j = 0; j < blockIndices.numAtomsPerMolecule; ++j)
922 t->a[t->nra++] = atomStart + j;
926 default: /* Should not be reached */
927 GMX_RELEASE_ASSERT(false, "Unreachable code was reached");
933 /* If not using completion, simply store the start of the block. */
934 t->index[t->nr++] = i;
938 /* Set the end of the last block */
939 t->index[t->nr] = t->nra;
940 /* Free any unnecessary memory */
941 srenew(t->index, t->nr+1);
942 t->nalloc_index = t->nr+1;
945 srenew(t->a, t->nra);
946 t->nalloc_a = t->nra;
951 * \param[in] g Index group to check.
952 * \param[in] b Block data to check against.
953 * \returns true if \p g consists of one or more complete blocks from \p b,
956 * The atoms in \p g are assumed to be sorted.
959 gmx_ana_index_has_full_blocks(const gmx_ana_index_t *g,
960 const gmx::RangePartitioning *b)
965 /* Each round in the loop matches one block */
968 /* Find the block that begins with the first unmatched atom */
969 while (bi < b->numBlocks() && b->block(bi).begin() != g->index[i])
973 /* If not found, or if too large, return */
974 if (bi == b->numBlocks() || i + b->block(bi).size() > g->isize)
978 /* Check that the block matches the index */
979 for (j = b->block(bi).begin(); j < b->block(bi).end(); ++j, ++i)
981 if (g->index[i] != j)
986 /* Move the search to the next block */
993 * \param[in] g Index group to check.
994 * \param[in] b Block data to check against.
995 * \returns true if \p g consists of one or more complete blocks from \p b,
998 * The atoms in \p g and \p b->a are assumed to be in the same order.
1001 gmx_ana_index_has_full_ablocks(gmx_ana_index_t *g, t_blocka *b)
1006 /* Each round in the loop matches one block */
1007 while (i < g->isize)
1009 /* Find the block that begins with the first unmatched atom */
1010 while (bi < b->nr && b->a[b->index[bi]] != g->index[i])
1014 /* If not found, or if too large, return */
1015 if (bi == b->nr || i + b->index[bi+1] - b->index[bi] > g->isize)
1019 /* Check that the block matches the index */
1020 for (j = b->index[bi]; j < b->index[bi+1]; ++j, ++i)
1022 if (b->a[j] != g->index[i])
1027 /* Move the search to the next block */
1034 * \brief Returns if an atom is at a residue boundary.
1036 * \param[in] top Topology data.
1037 * \param[in] a Atom index to check, should be -1 <= \p a < top->natoms.
1038 * \param[in,out] molb The molecule block of atom a
1039 * \returns true if atoms \p a and \p a + 1 are in different residues, false otherwise.
1041 static bool is_at_residue_boundary(const gmx_mtop_t *top, int a, int *molb)
1043 if (a == -1 || a + 1 == top->natoms)
1048 mtopGetAtomAndResidueName(top, a, molb,
1049 nullptr, nullptr, nullptr, &resindA);
1051 mtopGetAtomAndResidueName(top, a + 1, molb,
1052 nullptr, nullptr, nullptr, &resindAPlusOne);
1053 return resindAPlusOne != resindA;
1057 * \param[in] g Index group to check.
1058 * \param[in] type Block data to check against.
1059 * \param[in] top Topology data.
1060 * \returns true if \p g consists of one or more complete elements of type
1061 * \p type, false otherwise.
1063 * \p g is assumed to be sorted, otherwise may return false negatives.
1065 * If \p type is \ref INDEX_ATOM, the return value is always true.
1066 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
1070 gmx_ana_index_has_complete_elems(gmx_ana_index_t *g, e_index_t type,
1071 const gmx_mtop_t *top)
1078 // TODO: Consider whether unsorted groups need to be supported better.
1092 for (int i = 0; i < g->isize; ++i)
1094 const int a = g->index[i];
1095 // Check if a is consecutive or on a residue boundary
1098 if (!is_at_residue_boundary(top, aPrev, &molb))
1102 if (!is_at_residue_boundary(top, a - 1, &molb))
1109 GMX_ASSERT(g->isize > 0, "We return above when isize=0");
1110 const int a = g->index[g->isize - 1];
1111 if (!is_at_residue_boundary(top, a, &molb))
1120 auto molecules = gmx_mtop_molecules(*top);
1121 return gmx_ana_index_has_full_blocks(g, &molecules);
1128 * \param[out] m Output structure.
1130 * Any contents of \p m are discarded without freeing.
1133 gmx_ana_indexmap_clear(gmx_ana_indexmap_t *m)
1135 m->type = INDEX_UNKNOWN;
1139 m->mapb.index = nullptr;
1140 m->mapb.nalloc_index = 0;
1142 m->mapb.a = nullptr;
1143 m->mapb.nalloc_a = 0;
1146 m->b.index = nullptr;
1149 m->b.nalloc_index = 0;
1155 * \param[in,out] m Mapping structure.
1156 * \param[in] nr Maximum number of blocks to reserve space for.
1157 * \param[in] isize Maximum number of atoms to reserve space for.
1160 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t *m, int nr, int isize)
1162 if (m->mapb.nalloc_index < nr + 1)
1164 srenew(m->refid, nr);
1165 srenew(m->mapid, nr);
1166 srenew(m->orgid, nr);
1167 srenew(m->mapb.index, nr + 1);
1168 srenew(m->b.index, nr + 1);
1169 m->mapb.nalloc_index = nr + 1;
1170 m->b.nalloc_index = nr + 1;
1172 if (m->b.nalloc_a < isize)
1174 srenew(m->b.a, isize);
1175 m->b.nalloc_a = isize;
1180 * \param[in,out] m Mapping structure to initialize.
1181 * \param[in] g Index group to map
1182 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1183 * \param[in] top Topology structure
1184 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1185 * \param[in] type Type of mapping to construct.
1187 * Initializes a new index group mapping.
1188 * The index group provided to gmx_ana_indexmap_update() should always be a
1189 * subset of the \p g given here.
1191 * \p m should have been initialized somehow (calloc() is enough).
1194 gmx_ana_indexmap_init(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1195 const gmx_mtop_t *top, e_index_t type)
1198 gmx_ana_index_make_block(&m->b, top, g, type, false);
1199 gmx_ana_indexmap_reserve(m, m->b.nr, m->b.nra);
1201 for (int i = 0; i < m->b.nr; ++i)
1203 const int ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1204 next_group_index(ii, top, type, &id);
1209 m->mapb.nr = m->b.nr;
1210 m->mapb.nra = m->b.nra;
1212 std::memcpy(m->mapb.index, m->b.index, (m->b.nr+1)*sizeof(*(m->mapb.index)));
1217 gmx_ana_indexmap_init_orgid_group(gmx_ana_indexmap_t *m, const gmx_mtop_t *top,
1220 GMX_RELEASE_ASSERT(m->bStatic,
1221 "Changing original IDs is not supported after starting "
1222 "to use the mapping");
1223 GMX_RELEASE_ASSERT(top != nullptr || (type != INDEX_RES && type != INDEX_MOL),
1224 "Topology must be provided for residue or molecule blocks");
1225 // Check that all atoms in each block belong to the same group.
1226 // This is a separate loop for better error handling (no state is modified
1227 // if there is an error.
1228 if (type == INDEX_RES || type == INDEX_MOL)
1231 for (int i = 0; i < m->b.nr; ++i)
1233 const int ii = m->b.a[m->b.index[i]];
1234 if (next_group_index(ii, top, type, &id))
1236 for (int j = m->b.index[i] + 1; j < m->b.index[i+1]; ++j)
1238 if (next_group_index(m->b.a[j], top, type, &id))
1240 std::string message("Grouping into residues/molecules is ambiguous");
1241 GMX_THROW(gmx::InconsistentInputError(message));
1247 // Do a second loop, where things are actually set.
1250 for (int i = 0; i < m->b.nr; ++i)
1252 const int ii = (type == INDEX_UNKNOWN ? 0 : m->b.a[m->b.index[i]]);
1253 if (next_group_index(ii, top, type, &id))
1257 m->mapid[i] = group;
1258 m->orgid[i] = group;
1260 // Count also the last group.
1266 * \param[in,out] m Mapping structure to initialize.
1267 * \param[in] b Block information to use for data.
1269 * Frees some memory that is not necessary for static index group mappings.
1270 * Internal pointers are set to point to data in \p b; it is the responsibility
1271 * of the caller to ensure that the block information matches the contents of
1273 * After this function has been called, the index group provided to
1274 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1276 * This function breaks modularity of the index group mapping interface in an
1277 * ugly way, but allows reducing memory usage of static selections by a
1278 * significant amount.
1281 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t *m, t_blocka *b)
1284 sfree(m->mapb.index);
1287 m->mapb.nalloc_index = 0;
1288 m->mapb.nalloc_a = 0;
1289 m->b.nalloc_index = 0;
1291 m->mapid = m->orgid;
1292 m->mapb.index = b->index;
1294 m->b.index = b->index;
1299 * \param[in,out] dest Destination data structure.
1300 * \param[in] src Source mapping.
1301 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1302 * copy is made; otherwise, only variable parts are copied, and no memory
1305 * \p dest should have been initialized somehow (calloc() is enough).
1308 gmx_ana_indexmap_copy(gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst)
1312 gmx_ana_indexmap_reserve(dest, src->b.nr, src->b.nra);
1313 dest->type = src->type;
1314 dest->b.nr = src->b.nr;
1315 dest->b.nra = src->b.nra;
1316 std::memcpy(dest->orgid, src->orgid, dest->b.nr*sizeof(*dest->orgid));
1317 std::memcpy(dest->b.index, src->b.index, (dest->b.nr+1)*sizeof(*dest->b.index));
1318 std::memcpy(dest->b.a, src->b.a, dest->b.nra*sizeof(*dest->b.a));
1320 dest->mapb.nr = src->mapb.nr;
1321 dest->mapb.nra = src->mapb.nra;
1322 if (src->mapb.nalloc_a > 0)
1326 snew(dest->mapb.a, src->mapb.nalloc_a);
1327 dest->mapb.nalloc_a = src->mapb.nalloc_a;
1329 std::memcpy(dest->mapb.a, src->mapb.a, dest->mapb.nra*sizeof(*dest->mapb.a));
1333 dest->mapb.a = src->mapb.a;
1335 std::memcpy(dest->refid, src->refid, dest->mapb.nr*sizeof(*dest->refid));
1336 std::memcpy(dest->mapid, src->mapid, dest->mapb.nr*sizeof(*dest->mapid));
1337 std::memcpy(dest->mapb.index, src->mapb.index, (dest->mapb.nr+1)*sizeof(*dest->mapb.index));
1338 dest->bStatic = src->bStatic;
1342 * Helper function to set the source atoms in an index map.
1344 * \param[in,out] m Mapping structure.
1345 * \param[in] isize Number of atoms in the \p index array.
1346 * \param[in] index List of atoms.
1349 set_atoms(gmx_ana_indexmap_t *m, int isize, int *index)
1351 m->mapb.nra = isize;
1352 if (m->mapb.nalloc_a == 0)
1358 for (int i = 0; i < isize; ++i)
1360 m->mapb.a[i] = index[i];
1366 * \param[in,out] m Mapping structure.
1367 * \param[in] g Current index group.
1368 * \param[in] bMaskOnly true if the unused blocks should be masked with
1369 * -1 instead of removing them.
1371 * Updates the index group mapping with the new index group \p g.
1373 * \see gmx_ana_indexmap_t
1376 gmx_ana_indexmap_update(gmx_ana_indexmap_t *m, gmx_ana_index_t *g,
1381 /* Process the simple cases first */
1382 if (m->type == INDEX_UNKNOWN && m->b.nra == 0)
1386 if (m->type == INDEX_ALL)
1388 set_atoms(m, g->isize, g->index);
1391 m->mapb.index[1] = g->isize;
1395 /* Reset the reference IDs and mapping if necessary */
1396 const bool bToFull = (g->isize == m->b.nra);
1397 const bool bWasFull = (m->mapb.nra == m->b.nra);
1398 if (bToFull || bMaskOnly)
1402 for (bj = 0; bj < m->b.nr; ++bj)
1409 for (bj = 0; bj < m->b.nr; ++bj)
1411 m->mapid[bj] = m->orgid[bj];
1413 for (bj = 0; bj <= m->b.nr; ++bj)
1415 m->mapb.index[bj] = m->b.index[bj];
1418 set_atoms(m, m->b.nra, m->b.a);
1419 m->mapb.nr = m->b.nr;
1421 /* Exit immediately if the group is static */
1430 for (i = j = bj = 0; i < g->isize; ++i, ++j)
1432 /* Find the next atom in the block */
1433 while (m->b.a[j] != g->index[i])
1437 /* Mark blocks that did not contain any atoms */
1438 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1440 m->refid[bj++] = -1;
1442 /* Advance the block index if we have reached the next block */
1443 if (m->b.index[bj] <= j)
1448 /* Mark the last blocks as not accessible */
1449 while (bj < m->b.nr)
1451 m->refid[bj++] = -1;
1456 set_atoms(m, g->isize, g->index);
1457 for (i = j = bi = 0, bj = -1; i < g->isize; ++i)
1459 /* Find the next atom in the block */
1460 while (m->b.a[j] != g->index[i])
1464 /* If we have reached a new block, add it */
1465 if (m->b.index[bj+1] <= j)
1467 /* Skip any blocks in between */
1468 while (bj < m->b.nr && m->b.index[bj+1] <= j)
1473 m->mapid[bi] = m->orgid[bj];
1474 m->mapb.index[bi] = i;
1478 /* Update the number of blocks */
1479 m->mapb.index[bi] = g->isize;
1486 * \param[in,out] m Mapping structure to free.
1488 * All the memory allocated for the mapping structure is freed, and
1489 * the pointers set to NULL.
1490 * The pointer \p m is not freed.
1493 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t *m)
1496 if (m->mapid != m->orgid)
1500 if (m->mapb.nalloc_index > 0)
1502 sfree(m->mapb.index);
1504 if (m->mapb.nalloc_a > 0)
1509 if (m->b.nalloc_index > 0)
1513 if (m->b.nalloc_a > 0)
1517 gmx_ana_indexmap_clear(m);