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"
54 #include "gromacs/topology/block.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/topology/mtop_lookup.h"
57 #include "gromacs/topology/mtop_util.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/exceptions.h"
60 #include "gromacs/utility/gmxassert.h"
61 #include "gromacs/utility/smalloc.h"
62 #include "gromacs/utility/stringutil.h"
63 #include "gromacs/utility/textwriter.h"
68 IndexGroupsAndNames::IndexGroupsAndNames(
69 const t_blocka
&indexGroup
, ArrayRef
<char const * const> groupNames
)
70 : indexGroup_
{indexGroup
}
72 std::copy(groupNames
.begin(), groupNames
.end(), std::back_inserter(groupNames_
));
73 GMX_ASSERT(indexGroup_
.nr
== ssize(groupNames
),
74 "Number of groups must match number of group names.");
77 bool IndexGroupsAndNames::containsGroupName(const std::string
&groupName
) const
79 return std::any_of(std::begin(groupNames_
), std::end(groupNames_
),
80 [&groupName
](const std::string
&name
){return equalCaseInsensitive(groupName
, name
); });
83 std::vector
<index
> IndexGroupsAndNames::indices(const std::string
&groupName
) const
85 if (!containsGroupName(groupName
))
87 GMX_THROW(InconsistentInputError(std::string("Group ") + groupName
+
88 " referenced in the .mdp file was not found in the index file.\n"
89 "Group names must match either [moleculetype] names or custom index group\n"
90 "names, in which case you must supply an index file to the '-n' option\n"
93 const auto groupNamePosition
= std::find_if(std::begin(groupNames_
), std::end(groupNames_
),
94 [&groupName
](const std::string
&name
){return equalCaseInsensitive(groupName
, name
); });
95 const auto groupIndex
= std::distance(std::begin(groupNames_
), groupNamePosition
);
96 const auto groupSize
= indexGroup_
.index
[groupIndex
+1] - indexGroup_
.index
[groupIndex
];
97 std::vector
<index
> groupIndices(groupSize
);
98 const auto startingIndex
= indexGroup_
.index
[groupIndex
];
99 std::iota(std::begin(groupIndices
), std::end(groupIndices
), startingIndex
);
100 std::transform(std::begin(groupIndices
), std::end(groupIndices
), std::begin(groupIndices
),
101 [blockLookup
= indexGroup_
.a
](auto i
){return blockLookup
[i
]; });
107 /********************************************************************
108 * gmx_ana_indexgrps_t functions
109 ********************************************************************/
112 * Stores a set of index groups.
114 struct gmx_ana_indexgrps_t
116 //! Initializes an empty set of groups.
117 explicit gmx_ana_indexgrps_t(int nr
)
122 ~gmx_ana_indexgrps_t()
124 for (auto &indexGrp
: g
)
126 gmx_ana_index_deinit(&indexGrp
);
131 /** Array of index groups. */
132 std::vector
<gmx_ana_index_t
> g
;
134 std::vector
<std::string
> names
;
138 * \param[out] g Index group structure.
139 * \param[in] top Topology structure.
140 * \param[in] fnm File name for the index file.
141 * Memory is automatically allocated.
143 * One or both of \p top or \p fnm can be NULL.
144 * If \p top is NULL, an index file is required and the groups are read
145 * from the file (uses Gromacs routine init_index()).
146 * If \p fnm is NULL, default groups are constructed based on the
147 * topology (uses Gromacs routine analyse()).
148 * If both are null, the index group structure is initialized empty.
151 gmx_ana_indexgrps_init(gmx_ana_indexgrps_t
**g
, gmx_mtop_t
*top
,
154 t_blocka
*block
= nullptr;
155 char **names
= nullptr;
159 block
= init_index(fnm
, &names
);
163 block
= new_blocka();
164 // TODO: Propagate mtop further.
165 t_atoms atoms
= gmx_mtop_global_atoms(top
);
166 analyse(&atoms
, block
, &names
, FALSE
, FALSE
);
171 *g
= new gmx_ana_indexgrps_t(0);
177 *g
= new gmx_ana_indexgrps_t(block
->nr
);
178 for (int i
= 0; i
< block
->nr
; ++i
)
180 gmx_ana_index_t
*grp
= &(*g
)->g
[i
];
182 grp
->isize
= block
->index
[i
+1] - block
->index
[i
];
183 snew(grp
->index
, grp
->isize
);
184 for (int j
= 0; j
< grp
->isize
; ++j
)
186 grp
->index
[j
] = block
->a
[block
->index
[i
]+j
];
188 grp
->nalloc_index
= grp
->isize
;
189 (*g
)->names
.emplace_back(names
[i
]);
194 for (int i
= 0; i
< block
->nr
; ++i
)
203 for (int i
= 0; i
< block
->nr
; ++i
)
213 * \param[in] g Index groups structure.
215 * The pointer \p g is invalid after the call.
218 gmx_ana_indexgrps_free(gmx_ana_indexgrps_t
*g
)
225 * \param[out] dest Output structure.
226 * \param[out] destName Receives the name of the group if found.
227 * \param[in] src Input index groups.
228 * \param[in] n Number of the group to extract.
229 * \returns true if \p n is a valid group in \p src, false otherwise.
232 gmx_ana_indexgrps_extract(gmx_ana_index_t
*dest
, std::string
*destName
,
233 gmx_ana_indexgrps_t
*src
, int n
)
236 if (n
< 0 || n
>= gmx::index(src
->g
.size()))
242 if (destName
!= nullptr)
244 *destName
= src
->names
[n
];
246 gmx_ana_index_copy(dest
, &src
->g
[n
], true);
251 * \param[out] dest Output structure.
252 * \param[out] destName Receives the name of the group if found.
253 * \param[in] src Input index groups.
254 * \param[in] name Name (or part of the name) of the group to extract.
255 * \returns true if \p name is a valid group in \p src, false otherwise.
257 * Uses the Gromacs routine find_group() to find the actual group;
258 * the comparison is case-insensitive.
261 gmx_ana_indexgrps_find(gmx_ana_index_t
*dest
, std::string
*destName
,
262 gmx_ana_indexgrps_t
*src
,
268 snew(names
, src
->g
.size());
269 for (size_t i
= 0; i
< src
->g
.size(); ++i
)
271 names
[i
] = src
->names
[i
].c_str();
273 int n
= find_group(const_cast<char *>(name
), src
->g
.size(),
274 const_cast<char **>(names
));
282 return gmx_ana_indexgrps_extract(dest
, destName
, src
, n
);
286 * \param[in] writer Writer to use for output.
287 * \param[in] g Index groups to print.
288 * \param[in] maxn Maximum number of indices to print
289 * (-1 = print all, 0 = print only names).
292 gmx_ana_indexgrps_print(gmx::TextWriter
*writer
, gmx_ana_indexgrps_t
*g
, int maxn
)
294 for (gmx::index i
= 0; i
< gmx::ssize(g
->g
); ++i
)
296 writer
->writeString(gmx::formatString(" Group %2zd \"%s\" ",
297 i
, g
->names
[i
].c_str()));
298 gmx_ana_index_dump(writer
, &g
->g
[i
], maxn
);
302 /********************************************************************
303 * gmx_ana_index_t functions
304 ********************************************************************/
307 * \param[in,out] g Index group structure.
308 * \param[in] isize Maximum number of atoms to reserve space for.
311 gmx_ana_index_reserve(gmx_ana_index_t
*g
, int isize
)
313 if (g
->nalloc_index
< isize
)
315 srenew(g
->index
, isize
);
316 g
->nalloc_index
= isize
;
321 * \param[in,out] g Index group structure.
323 * Resizes the memory allocated for holding the indices such that the
324 * current contents fit.
327 gmx_ana_index_squeeze(gmx_ana_index_t
*g
)
329 srenew(g
->index
, g
->isize
);
330 g
->nalloc_index
= g
->isize
;
334 * \param[out] g Output structure.
336 * Any contents of \p g are discarded without freeing.
339 gmx_ana_index_clear(gmx_ana_index_t
*g
)
347 * \param[out] g Output structure.
348 * \param[in] isize Number of atoms in the new group.
349 * \param[in] index Array of \p isize atoms (can be NULL if \p isize is 0).
350 * \param[in] nalloc Number of elements allocated for \p index
351 * (if 0, \p index is not freed in gmx_ana_index_deinit())
353 * No copy if \p index is made.
356 gmx_ana_index_set(gmx_ana_index_t
*g
, int isize
, int *index
, int nalloc
)
360 g
->nalloc_index
= nalloc
;
364 * \param[out] g Output structure.
365 * \param[in] natoms Number of atoms.
368 gmx_ana_index_init_simple(gmx_ana_index_t
*g
, int natoms
)
373 snew(g
->index
, natoms
);
374 for (i
= 0; i
< natoms
; ++i
)
378 g
->nalloc_index
= natoms
;
382 * \param[in] g Index group structure.
384 * The pointer \p g is not freed.
387 gmx_ana_index_deinit(gmx_ana_index_t
*g
)
389 if (g
->nalloc_index
> 0)
393 gmx_ana_index_clear(g
);
397 * \param[out] dest Destination index group.
398 * \param[in] src Source index group.
399 * \param[in] bAlloc If true, memory is allocated at \p dest; otherwise,
400 * it is assumed that enough memory has been allocated for index.
403 gmx_ana_index_copy(gmx_ana_index_t
*dest
, gmx_ana_index_t
*src
, bool bAlloc
)
405 dest
->isize
= src
->isize
;
408 snew(dest
->index
, dest
->isize
);
409 dest
->nalloc_index
= dest
->isize
;
413 std::memcpy(dest
->index
, src
->index
, dest
->isize
*sizeof(*dest
->index
));
418 * \param[in] writer Writer to use for output.
419 * \param[in] g Index group to print.
420 * \param[in] maxn Maximum number of indices to print (-1 = print all).
423 gmx_ana_index_dump(gmx::TextWriter
*writer
, gmx_ana_index_t
*g
, int maxn
)
425 writer
->writeString(gmx::formatString("(%d atoms)", g
->isize
));
428 writer
->writeString(":");
430 if (maxn
>= 0 && n
> maxn
)
434 for (int j
= 0; j
< n
; ++j
)
436 writer
->writeString(gmx::formatString(" %d", g
->index
[j
]+1));
440 writer
->writeString(" ...");
443 writer
->ensureLineBreak();
447 gmx_ana_index_get_max_index(gmx_ana_index_t
*g
)
455 return *std::max_element(g
->index
, g
->index
+ g
->isize
);
460 * \param[in] g Index group to check.
461 * \returns true if the index group is sorted and has no duplicates,
465 gmx_ana_index_check_sorted(gmx_ana_index_t
*g
)
469 for (i
= 0; i
< g
->isize
-1; ++i
)
471 if (g
->index
[i
+1] <= g
->index
[i
])
480 gmx_ana_index_check_range(gmx_ana_index_t
*g
, int natoms
)
482 for (int i
= 0; i
< g
->isize
; ++i
)
484 if (g
->index
[i
] < 0 || g
->index
[i
] >= natoms
)
492 /********************************************************************
494 ********************************************************************/
497 * \param[in,out] g Index group to be sorted.
500 gmx_ana_index_sort(gmx_ana_index_t
*g
)
502 std::sort(g
->index
, g
->index
+g
->isize
);
506 gmx_ana_index_remove_duplicates(gmx_ana_index_t
*g
)
509 for (int i
= 0; i
< g
->isize
; ++i
)
511 if (i
== 0 || g
->index
[i
-1] != g
->index
[i
])
513 g
->index
[j
] = g
->index
[i
];
521 * \param[in] a Index group to check.
522 * \param[in] b Index group to check.
523 * \returns true if \p a and \p b are equal, false otherwise.
526 gmx_ana_index_equals(gmx_ana_index_t
*a
, gmx_ana_index_t
*b
)
530 if (a
->isize
!= b
->isize
)
534 for (i
= 0; i
< a
->isize
; ++i
)
536 if (a
->index
[i
] != b
->index
[i
])
545 * \param[in] a Index group to check against.
546 * \param[in] b Index group to check.
547 * \returns true if \p b is contained in \p a,
550 * If the elements are not in the same order in both groups, the function
551 * fails. However, the groups do not need to be sorted.
554 gmx_ana_index_contains(gmx_ana_index_t
*a
, gmx_ana_index_t
*b
)
558 for (i
= j
= 0; j
< b
->isize
; ++i
, ++j
)
560 while (i
< a
->isize
&& a
->index
[i
] != b
->index
[j
])
573 * \param[out] dest Output index group (the intersection of \p a and \p b).
574 * \param[in] a First index group.
575 * \param[in] b Second index group.
577 * \p dest can be the same as \p a or \p b.
580 gmx_ana_index_intersection(gmx_ana_index_t
*dest
,
581 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
)
585 for (i
= j
= k
= 0; i
< a
->isize
&& j
< b
->isize
; ++i
)
587 while (j
< b
->isize
&& b
->index
[j
] < a
->index
[i
])
591 if (j
< b
->isize
&& b
->index
[j
] == a
->index
[i
])
593 dest
->index
[k
++] = b
->index
[j
++];
600 * \param[out] dest Output index group (the difference \p a - \p b).
601 * \param[in] a First index group.
602 * \param[in] b Second index group.
604 * \p dest can equal \p a, but not \p b.
607 gmx_ana_index_difference(gmx_ana_index_t
*dest
,
608 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
)
612 for (i
= j
= k
= 0; i
< a
->isize
; ++i
)
614 while (j
< b
->isize
&& b
->index
[j
] < a
->index
[i
])
618 if (j
== b
->isize
|| b
->index
[j
] != a
->index
[i
])
620 dest
->index
[k
++] = a
->index
[i
];
627 * \param[in] a First index group.
628 * \param[in] b Second index group.
629 * \returns Size of the difference \p a - \p b.
632 gmx_ana_index_difference_size(gmx_ana_index_t
*a
, gmx_ana_index_t
*b
)
636 for (i
= j
= k
= 0; i
< a
->isize
; ++i
)
638 while (j
< b
->isize
&& b
->index
[j
] < a
->index
[i
])
642 if (j
== b
->isize
|| b
->index
[j
] != a
->index
[i
])
651 * \param[out] dest1 Output group 1 (will equal \p g).
652 * \param[out] dest2 Output group 2 (will equal \p src - \p g).
653 * \param[in] src Group to be partitioned.
654 * \param[in] g One partition.
656 * \pre \p g is a subset of \p src and both sets are sorted
657 * \pre \p dest1 has allocated storage to store \p src
658 * \post \p dest1 == \p g
659 * \post \p dest2 == \p src - \p g
661 * No storage should be allocated for \p dest2; after the call,
662 * \p dest2->index points to the memory allocated for \p dest1
663 * (to a part that is not used by \p dest1).
665 * The calculation can be performed in-place by setting \p dest1 equal to
669 gmx_ana_index_partition(gmx_ana_index_t
*dest1
, gmx_ana_index_t
*dest2
,
670 gmx_ana_index_t
*src
, gmx_ana_index_t
*g
)
674 dest2
->index
= dest1
->index
+ g
->isize
;
675 dest2
->isize
= src
->isize
- g
->isize
;
676 for (i
= g
->isize
-1, j
= src
->isize
-1, k
= dest2
->isize
-1; i
>= 0; --i
, --j
)
678 while (j
>= 0 && src
->index
[j
] != g
->index
[i
])
680 dest2
->index
[k
--] = src
->index
[j
--];
685 dest2
->index
[k
--] = src
->index
[j
--];
687 gmx_ana_index_copy(dest1
, g
, false);
691 * \param[out] dest Output index group (the union of \p a and \p b).
692 * \param[in] a First index group.
693 * \param[in] b Second index group.
695 * \p a and \p b can have common items.
696 * \p dest can equal \p a or \p b.
698 * \see gmx_ana_index_merge()
701 gmx_ana_index_union(gmx_ana_index_t
*dest
,
702 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
)
707 dsize
= gmx_ana_index_difference_size(b
, a
);
710 dest
->isize
= a
->isize
+ dsize
;
711 for (k
= dest
->isize
- 1; k
>= 0; k
--)
713 if (i
< 0 || (j
>= 0 && a
->index
[i
] < b
->index
[j
]))
715 dest
->index
[k
] = b
->index
[j
--];
719 if (j
>= 0 && a
->index
[i
] == b
->index
[j
])
723 dest
->index
[k
] = a
->index
[i
--];
729 gmx_ana_index_union_unsorted(gmx_ana_index_t
*dest
,
730 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
)
732 if (gmx_ana_index_check_sorted(b
))
734 gmx_ana_index_union(dest
, a
, b
);
739 gmx_ana_index_copy(&tmp
, b
, true);
740 gmx_ana_index_sort(&tmp
);
741 gmx_ana_index_remove_duplicates(&tmp
);
742 gmx_ana_index_union(dest
, a
, &tmp
);
743 gmx_ana_index_deinit(&tmp
);
748 * \param[out] dest Output index group (the union of \p a and \p b).
749 * \param[in] a First index group.
750 * \param[in] b Second index group.
752 * \p a and \p b should not have common items.
753 * \p dest can equal \p a or \p b.
755 * \see gmx_ana_index_union()
758 gmx_ana_index_merge(gmx_ana_index_t
*dest
,
759 gmx_ana_index_t
*a
, gmx_ana_index_t
*b
)
765 dest
->isize
= a
->isize
+ b
->isize
;
766 for (k
= dest
->isize
- 1; k
>= 0; k
--)
768 if (i
< 0 || (j
>= 0 && a
->index
[i
] < b
->index
[j
]))
770 dest
->index
[k
] = b
->index
[j
--];
774 dest
->index
[k
] = a
->index
[i
--];
779 /********************************************************************
780 * gmx_ana_indexmap_t and related things
781 ********************************************************************/
784 * Helper for splitting a sequence of atom indices into groups.
786 * \param[in] atomIndex Index of the next atom in the sequence.
787 * \param[in] top Topology structure.
788 * \param[in] type Type of group to split into.
789 * \param[in,out] id Variable to receive the next group id.
790 * \returns `true` if \p atomIndex starts a new group in the sequence, i.e.,
791 * if \p *id was changed.
793 * \p *id should be initialized to `-1` before first call of this function, and
794 * then each atom index in the sequence passed to the function in turn.
796 * \ingroup module_selection
799 next_group_index(int atomIndex
, const gmx_mtop_t
*top
,
800 e_index_t type
, int *id
)
810 int resind
, molb
= 0;
811 mtopGetAtomAndResidueName(top
, atomIndex
, &molb
, nullptr, nullptr, nullptr, &resind
);
818 *id
= mtopGetMoleculeIndex(top
, atomIndex
, &molb
);
830 * \param[in,out] t Output block.
831 * \param[in] top Topology structure
832 * (only used if \p type is \ref INDEX_RES or \ref INDEX_MOL, can be NULL
834 * \param[in] g Index group
835 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
836 * \param[in] type Type of partitioning to make.
837 * \param[in] bComplete
838 * If true, the index group is expanded to include any residue/molecule
839 * (depending on \p type) that is partially contained in the group.
840 * If \p type is not INDEX_RES or INDEX_MOL, this has no effect.
842 * \p m should have been initialized somehow (calloc() is enough).
843 * \p g should be sorted.
846 gmx_ana_index_make_block(t_blocka
*t
, const gmx_mtop_t
*top
, gmx_ana_index_t
*g
,
847 e_index_t type
, bool bComplete
)
849 if (type
== INDEX_UNKNOWN
)
863 // TODO: Check callers and either check these there as well, or turn these
865 GMX_RELEASE_ASSERT(top
!= nullptr || (type
!= INDEX_RES
&& type
!= INDEX_MOL
),
866 "Topology must be provided for residue or molecule blocks");
867 GMX_RELEASE_ASSERT(type
!= INDEX_MOL
|| top
->haveMoleculeIndices
,
868 "Molecule information must be present for molecule blocks");
870 /* bComplete only does something for INDEX_RES or INDEX_MOL, so turn it
872 if (type
!= INDEX_RES
&& type
!= INDEX_MOL
)
876 /* Allocate memory for the atom array and fill it unless we are using
881 /* We may allocate some extra memory here because we don't know in
882 * advance how much will be needed. */
883 if (t
->nalloc_a
< top
->natoms
)
885 srenew(t
->a
, top
->natoms
);
886 t
->nalloc_a
= top
->natoms
;
892 if (t
->nalloc_a
< g
->isize
)
894 srenew(t
->a
, g
->isize
);
895 t
->nalloc_a
= g
->isize
;
897 std::memcpy(t
->a
, g
->index
, g
->isize
*sizeof(*(t
->a
)));
900 /* Allocate memory for the block index. We don't know in advance
901 * how much will be needed, so we allocate some extra and free it in the
903 if (t
->nalloc_index
< g
->isize
+ 1)
905 srenew(t
->index
, g
->isize
+ 1);
906 t
->nalloc_index
= g
->isize
+ 1;
912 for (int i
= 0; i
< g
->isize
; ++i
)
914 const int ai
= g
->index
[i
];
915 /* Find the ID number of the atom/residue/molecule corresponding to
917 if (next_group_index(ai
, top
, type
, &id
))
919 /* If this is the first atom in a new block, initialize the block. */
922 /* For completion, we first set the start of the block. */
923 t
->index
[t
->nr
++] = t
->nra
;
924 /* And then we find all the atoms that should be included. */
930 mtopGetMolblockIndex(top
, ai
, &molb
, &molnr
, &atnr_mol
);
931 const t_atoms
&mol_atoms
= top
->moltype
[top
->molblock
[molb
].type
].atoms
;
932 int last_atom
= atnr_mol
+ 1;
933 const int currentResid
= mol_atoms
.atom
[atnr_mol
].resind
;
934 while (last_atom
< mol_atoms
.nr
935 && mol_atoms
.atom
[last_atom
].resind
== currentResid
)
939 int first_atom
= atnr_mol
- 1;
940 while (first_atom
>= 0
941 && mol_atoms
.atom
[first_atom
].resind
== currentResid
)
945 const MoleculeBlockIndices
&molBlock
= top
->moleculeBlockIndices
[molb
];
946 int first_mol_atom
= molBlock
.globalAtomStart
;
947 first_mol_atom
+= molnr
*molBlock
.numAtomsPerMolecule
;
948 first_atom
= first_mol_atom
+ first_atom
+ 1;
949 last_atom
= first_mol_atom
+ last_atom
- 1;
950 for (int j
= first_atom
; j
<= last_atom
; ++j
)
959 mtopGetMolblockIndex(top
, ai
, &molb
, &molnr
, &atnr_mol
);
960 const MoleculeBlockIndices
&blockIndices
= top
->moleculeBlockIndices
[molb
];
961 const int atomStart
= blockIndices
.globalAtomStart
+ (id
- blockIndices
.moleculeIndexStart
)*blockIndices
.numAtomsPerMolecule
;
962 for (int j
= 0; j
< blockIndices
.numAtomsPerMolecule
; ++j
)
964 t
->a
[t
->nra
++] = atomStart
+ j
;
968 default: /* Should not be reached */
969 GMX_RELEASE_ASSERT(false, "Unreachable code was reached");
975 /* If not using completion, simply store the start of the block. */
976 t
->index
[t
->nr
++] = i
;
980 /* Set the end of the last block */
981 t
->index
[t
->nr
] = t
->nra
;
982 /* Free any unnecessary memory */
983 srenew(t
->index
, t
->nr
+1);
984 t
->nalloc_index
= t
->nr
+1;
987 srenew(t
->a
, t
->nra
);
988 t
->nalloc_a
= t
->nra
;
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 are assumed to be sorted.
1001 gmx_ana_index_has_full_blocks(const gmx_ana_index_t
*g
,
1002 const gmx::RangePartitioning
*b
)
1007 /* Each round in the loop matches one block */
1008 while (i
< g
->isize
)
1010 /* Find the block that begins with the first unmatched atom */
1011 while (bi
< b
->numBlocks() && *b
->block(bi
).begin() != g
->index
[i
])
1015 /* If not found, or if too large, return */
1016 if (bi
== b
->numBlocks() || i
+ b
->block(bi
).size() > g
->isize
)
1020 /* Check that the block matches the index */
1021 for (j
= *b
->block(bi
).begin(); j
< *b
->block(bi
).end(); ++j
, ++i
)
1023 if (g
->index
[i
] != j
)
1028 /* Move the search to the next block */
1035 * \param[in] g Index group to check.
1036 * \param[in] b Block data to check against.
1037 * \returns true if \p g consists of one or more complete blocks from \p b,
1040 * The atoms in \p g and \p b->a are assumed to be in the same order.
1043 gmx_ana_index_has_full_ablocks(gmx_ana_index_t
*g
, t_blocka
*b
)
1048 /* Each round in the loop matches one block */
1049 while (i
< g
->isize
)
1051 /* Find the block that begins with the first unmatched atom */
1052 while (bi
< b
->nr
&& b
->a
[b
->index
[bi
]] != g
->index
[i
])
1056 /* If not found, or if too large, return */
1057 if (bi
== b
->nr
|| i
+ b
->index
[bi
+1] - b
->index
[bi
] > g
->isize
)
1061 /* Check that the block matches the index */
1062 for (j
= b
->index
[bi
]; j
< b
->index
[bi
+1]; ++j
, ++i
)
1064 if (b
->a
[j
] != g
->index
[i
])
1069 /* Move the search to the next block */
1076 * \brief Returns if an atom is at a residue boundary.
1078 * \param[in] top Topology data.
1079 * \param[in] a Atom index to check, should be -1 <= \p a < top->natoms.
1080 * \param[in,out] molb The molecule block of atom a
1081 * \returns true if atoms \p a and \p a + 1 are in different residues, false otherwise.
1083 static bool is_at_residue_boundary(const gmx_mtop_t
*top
, int a
, int *molb
)
1085 if (a
== -1 || a
+ 1 == top
->natoms
)
1090 mtopGetAtomAndResidueName(top
, a
, molb
,
1091 nullptr, nullptr, nullptr, &resindA
);
1093 mtopGetAtomAndResidueName(top
, a
+ 1, molb
,
1094 nullptr, nullptr, nullptr, &resindAPlusOne
);
1095 return resindAPlusOne
!= resindA
;
1099 * \param[in] g Index group to check.
1100 * \param[in] type Block data to check against.
1101 * \param[in] top Topology data.
1102 * \returns true if \p g consists of one or more complete elements of type
1103 * \p type, false otherwise.
1105 * \p g is assumed to be sorted, otherwise may return false negatives.
1107 * If \p type is \ref INDEX_ATOM, the return value is always true.
1108 * If \p type is \ref INDEX_UNKNOWN or \ref INDEX_ALL, the return value is
1112 gmx_ana_index_has_complete_elems(gmx_ana_index_t
*g
, e_index_t type
,
1113 const gmx_mtop_t
*top
)
1120 // TODO: Consider whether unsorted groups need to be supported better.
1134 for (int i
= 0; i
< g
->isize
; ++i
)
1136 const int a
= g
->index
[i
];
1137 // Check if a is consecutive or on a residue boundary
1140 if (!is_at_residue_boundary(top
, aPrev
, &molb
))
1144 if (!is_at_residue_boundary(top
, a
- 1, &molb
))
1151 GMX_ASSERT(g
->isize
> 0, "We return above when isize=0");
1152 const int a
= g
->index
[g
->isize
- 1];
1153 if (!is_at_residue_boundary(top
, a
, &molb
))
1162 auto molecules
= gmx_mtop_molecules(*top
);
1163 return gmx_ana_index_has_full_blocks(g
, &molecules
);
1170 * \param[out] m Output structure.
1172 * Any contents of \p m are discarded without freeing.
1175 gmx_ana_indexmap_clear(gmx_ana_indexmap_t
*m
)
1177 m
->type
= INDEX_UNKNOWN
;
1181 m
->mapb
.index
= nullptr;
1182 m
->mapb
.nalloc_index
= 0;
1184 m
->mapb
.a
= nullptr;
1185 m
->mapb
.nalloc_a
= 0;
1188 m
->b
.index
= nullptr;
1191 m
->b
.nalloc_index
= 0;
1197 * \param[in,out] m Mapping structure.
1198 * \param[in] nr Maximum number of blocks to reserve space for.
1199 * \param[in] isize Maximum number of atoms to reserve space for.
1202 gmx_ana_indexmap_reserve(gmx_ana_indexmap_t
*m
, int nr
, int isize
)
1204 if (m
->mapb
.nalloc_index
< nr
+ 1)
1206 srenew(m
->refid
, nr
);
1207 srenew(m
->mapid
, nr
);
1208 srenew(m
->orgid
, nr
);
1209 srenew(m
->mapb
.index
, nr
+ 1);
1210 srenew(m
->b
.index
, nr
+ 1);
1211 m
->mapb
.nalloc_index
= nr
+ 1;
1212 m
->b
.nalloc_index
= nr
+ 1;
1214 if (m
->b
.nalloc_a
< isize
)
1216 srenew(m
->b
.a
, isize
);
1217 m
->b
.nalloc_a
= isize
;
1222 * \param[in,out] m Mapping structure to initialize.
1223 * \param[in] g Index group to map
1224 * (can be NULL if \p type is \ref INDEX_UNKNOWN).
1225 * \param[in] top Topology structure
1226 * (can be NULL if \p type is not \ref INDEX_RES or \ref INDEX_MOL).
1227 * \param[in] type Type of mapping to construct.
1229 * Initializes a new index group mapping.
1230 * The index group provided to gmx_ana_indexmap_update() should always be a
1231 * subset of the \p g given here.
1233 * \p m should have been initialized somehow (calloc() is enough).
1236 gmx_ana_indexmap_init(gmx_ana_indexmap_t
*m
, gmx_ana_index_t
*g
,
1237 const gmx_mtop_t
*top
, e_index_t type
)
1240 gmx_ana_index_make_block(&m
->b
, top
, g
, type
, false);
1241 gmx_ana_indexmap_reserve(m
, m
->b
.nr
, m
->b
.nra
);
1243 for (int i
= 0; i
< m
->b
.nr
; ++i
)
1245 const int ii
= (type
== INDEX_UNKNOWN
? 0 : m
->b
.a
[m
->b
.index
[i
]]);
1246 next_group_index(ii
, top
, type
, &id
);
1251 m
->mapb
.nr
= m
->b
.nr
;
1252 m
->mapb
.nra
= m
->b
.nra
;
1254 std::memcpy(m
->mapb
.index
, m
->b
.index
, (m
->b
.nr
+1)*sizeof(*(m
->mapb
.index
)));
1259 gmx_ana_indexmap_init_orgid_group(gmx_ana_indexmap_t
*m
, const gmx_mtop_t
*top
,
1262 GMX_RELEASE_ASSERT(m
->bStatic
,
1263 "Changing original IDs is not supported after starting "
1264 "to use the mapping");
1265 GMX_RELEASE_ASSERT(top
!= nullptr || (type
!= INDEX_RES
&& type
!= INDEX_MOL
),
1266 "Topology must be provided for residue or molecule blocks");
1267 // Check that all atoms in each block belong to the same group.
1268 // This is a separate loop for better error handling (no state is modified
1269 // if there is an error.
1270 if (type
== INDEX_RES
|| type
== INDEX_MOL
)
1273 for (int i
= 0; i
< m
->b
.nr
; ++i
)
1275 const int ii
= m
->b
.a
[m
->b
.index
[i
]];
1276 if (next_group_index(ii
, top
, type
, &id
))
1278 for (int j
= m
->b
.index
[i
] + 1; j
< m
->b
.index
[i
+1]; ++j
)
1280 if (next_group_index(m
->b
.a
[j
], top
, type
, &id
))
1282 std::string
message("Grouping into residues/molecules is ambiguous");
1283 GMX_THROW(gmx::InconsistentInputError(message
));
1289 // Do a second loop, where things are actually set.
1292 for (int i
= 0; i
< m
->b
.nr
; ++i
)
1294 const int ii
= (type
== INDEX_UNKNOWN
? 0 : m
->b
.a
[m
->b
.index
[i
]]);
1295 if (next_group_index(ii
, top
, type
, &id
))
1299 m
->mapid
[i
] = group
;
1300 m
->orgid
[i
] = group
;
1302 // Count also the last group.
1308 * \param[in,out] m Mapping structure to initialize.
1309 * \param[in] b Block information to use for data.
1311 * Frees some memory that is not necessary for static index group mappings.
1312 * Internal pointers are set to point to data in \p b; it is the responsibility
1313 * of the caller to ensure that the block information matches the contents of
1315 * After this function has been called, the index group provided to
1316 * gmx_ana_indexmap_update() should always be the same as \p g given here.
1318 * This function breaks modularity of the index group mapping interface in an
1319 * ugly way, but allows reducing memory usage of static selections by a
1320 * significant amount.
1323 gmx_ana_indexmap_set_static(gmx_ana_indexmap_t
*m
, t_blocka
*b
)
1326 sfree(m
->mapb
.index
);
1329 m
->mapb
.nalloc_index
= 0;
1330 m
->mapb
.nalloc_a
= 0;
1331 m
->b
.nalloc_index
= 0;
1333 m
->mapid
= m
->orgid
;
1334 m
->mapb
.index
= b
->index
;
1336 m
->b
.index
= b
->index
;
1341 * \param[in,out] dest Destination data structure.
1342 * \param[in] src Source mapping.
1343 * \param[in] bFirst If true, memory is allocated for \p dest and a full
1344 * copy is made; otherwise, only variable parts are copied, and no memory
1347 * \p dest should have been initialized somehow (calloc() is enough).
1350 gmx_ana_indexmap_copy(gmx_ana_indexmap_t
*dest
, gmx_ana_indexmap_t
*src
, bool bFirst
)
1354 gmx_ana_indexmap_reserve(dest
, src
->b
.nr
, src
->b
.nra
);
1355 dest
->type
= src
->type
;
1356 dest
->b
.nr
= src
->b
.nr
;
1357 dest
->b
.nra
= src
->b
.nra
;
1358 std::memcpy(dest
->orgid
, src
->orgid
, dest
->b
.nr
*sizeof(*dest
->orgid
));
1359 std::memcpy(dest
->b
.index
, src
->b
.index
, (dest
->b
.nr
+1)*sizeof(*dest
->b
.index
));
1360 std::memcpy(dest
->b
.a
, src
->b
.a
, dest
->b
.nra
*sizeof(*dest
->b
.a
));
1362 dest
->mapb
.nr
= src
->mapb
.nr
;
1363 dest
->mapb
.nra
= src
->mapb
.nra
;
1364 if (src
->mapb
.nalloc_a
> 0)
1368 snew(dest
->mapb
.a
, src
->mapb
.nalloc_a
);
1369 dest
->mapb
.nalloc_a
= src
->mapb
.nalloc_a
;
1371 std::memcpy(dest
->mapb
.a
, src
->mapb
.a
, dest
->mapb
.nra
*sizeof(*dest
->mapb
.a
));
1375 dest
->mapb
.a
= src
->mapb
.a
;
1377 std::memcpy(dest
->refid
, src
->refid
, dest
->mapb
.nr
*sizeof(*dest
->refid
));
1378 std::memcpy(dest
->mapid
, src
->mapid
, dest
->mapb
.nr
*sizeof(*dest
->mapid
));
1379 std::memcpy(dest
->mapb
.index
, src
->mapb
.index
, (dest
->mapb
.nr
+1)*sizeof(*dest
->mapb
.index
));
1380 dest
->bStatic
= src
->bStatic
;
1384 * Helper function to set the source atoms in an index map.
1386 * \param[in,out] m Mapping structure.
1387 * \param[in] isize Number of atoms in the \p index array.
1388 * \param[in] index List of atoms.
1391 set_atoms(gmx_ana_indexmap_t
*m
, int isize
, int *index
)
1393 m
->mapb
.nra
= isize
;
1394 if (m
->mapb
.nalloc_a
== 0)
1400 for (int i
= 0; i
< isize
; ++i
)
1402 m
->mapb
.a
[i
] = index
[i
];
1408 * \param[in,out] m Mapping structure.
1409 * \param[in] g Current index group.
1410 * \param[in] bMaskOnly true if the unused blocks should be masked with
1411 * -1 instead of removing them.
1413 * Updates the index group mapping with the new index group \p g.
1415 * \see gmx_ana_indexmap_t
1418 gmx_ana_indexmap_update(gmx_ana_indexmap_t
*m
, gmx_ana_index_t
*g
,
1423 /* Process the simple cases first */
1424 if (m
->type
== INDEX_UNKNOWN
&& m
->b
.nra
== 0)
1428 if (m
->type
== INDEX_ALL
)
1430 set_atoms(m
, g
->isize
, g
->index
);
1433 m
->mapb
.index
[1] = g
->isize
;
1437 /* Reset the reference IDs and mapping if necessary */
1438 const bool bToFull
= (g
->isize
== m
->b
.nra
);
1439 const bool bWasFull
= (m
->mapb
.nra
== m
->b
.nra
);
1440 if (bToFull
|| bMaskOnly
)
1444 for (bj
= 0; bj
< m
->b
.nr
; ++bj
)
1451 for (bj
= 0; bj
< m
->b
.nr
; ++bj
)
1453 m
->mapid
[bj
] = m
->orgid
[bj
];
1455 for (bj
= 0; bj
<= m
->b
.nr
; ++bj
)
1457 m
->mapb
.index
[bj
] = m
->b
.index
[bj
];
1460 set_atoms(m
, m
->b
.nra
, m
->b
.a
);
1461 m
->mapb
.nr
= m
->b
.nr
;
1463 /* Exit immediately if the group is static */
1472 for (i
= j
= bj
= 0; i
< g
->isize
; ++i
, ++j
)
1474 /* Find the next atom in the block */
1475 while (m
->b
.a
[j
] != g
->index
[i
])
1479 /* Mark blocks that did not contain any atoms */
1480 while (bj
< m
->b
.nr
&& m
->b
.index
[bj
+1] <= j
)
1482 m
->refid
[bj
++] = -1;
1484 /* Advance the block index if we have reached the next block */
1485 if (m
->b
.index
[bj
] <= j
)
1490 /* Mark the last blocks as not accessible */
1491 while (bj
< m
->b
.nr
)
1493 m
->refid
[bj
++] = -1;
1498 set_atoms(m
, g
->isize
, g
->index
);
1499 for (i
= j
= bi
= 0, bj
= -1; i
< g
->isize
; ++i
)
1501 /* Find the next atom in the block */
1502 while (m
->b
.a
[j
] != g
->index
[i
])
1506 /* If we have reached a new block, add it */
1507 if (m
->b
.index
[bj
+1] <= j
)
1509 /* Skip any blocks in between */
1510 while (bj
< m
->b
.nr
&& m
->b
.index
[bj
+1] <= j
)
1515 m
->mapid
[bi
] = m
->orgid
[bj
];
1516 m
->mapb
.index
[bi
] = i
;
1520 /* Update the number of blocks */
1521 m
->mapb
.index
[bi
] = g
->isize
;
1528 * \param[in,out] m Mapping structure to free.
1530 * All the memory allocated for the mapping structure is freed, and
1531 * the pointers set to NULL.
1532 * The pointer \p m is not freed.
1535 gmx_ana_indexmap_deinit(gmx_ana_indexmap_t
*m
)
1538 if (m
->mapid
!= m
->orgid
)
1542 if (m
->mapb
.nalloc_index
> 0)
1544 sfree(m
->mapb
.index
);
1546 if (m
->mapb
.nalloc_a
> 0)
1551 if (m
->b
.nalloc_index
> 0)
1555 if (m
->b
.nalloc_a
> 0)
1559 gmx_ana_indexmap_clear(m
);