Gromacs
2020.4
|
#include "gmxpre.h"
#include "indexutil.h"
#include <cstdlib>
#include <cstring>
#include <algorithm>
#include <numeric>
#include <string>
#include <vector>
#include "gromacs/topology/block.h"
#include "gromacs/topology/index.h"
#include "gromacs/topology/mtop_lookup.h"
#include "gromacs/topology/mtop_util.h"
#include "gromacs/topology/topology.h"
#include "gromacs/utility/exceptions.h"
#include "gromacs/utility/gmxassert.h"
#include "gromacs/utility/smalloc.h"
#include "gromacs/utility/stringutil.h"
#include "gromacs/utility/textwriter.h"
Implements functions in indexutil.h.
Classes | |
struct | gmx_ana_indexgrps_t |
Stores a set of index groups. More... | |
Functions | |
void | gmx_ana_indexgrps_init (gmx_ana_indexgrps_t **g, gmx_mtop_t *top, const char *fnm) |
Reads index groups from a file or constructs them from topology. More... | |
void | gmx_ana_indexgrps_free (gmx_ana_indexgrps_t *g) |
Frees memory allocated for index groups. More... | |
bool | gmx_ana_indexgrps_extract (gmx_ana_index_t *dest, std::string *destName, gmx_ana_indexgrps_t *src, int n) |
Extracts a single index group. More... | |
bool | gmx_ana_indexgrps_find (gmx_ana_index_t *dest, std::string *destName, gmx_ana_indexgrps_t *src, const char *name) |
Finds and extracts a single index group by name. More... | |
void | gmx_ana_indexgrps_print (gmx::TextWriter *writer, gmx_ana_indexgrps_t *g, int maxn) |
Writes out a list of index groups. More... | |
void | gmx_ana_index_reserve (gmx_ana_index_t *g, int isize) |
Reserves memory to store an index group of size isize . More... | |
void | gmx_ana_index_squeeze (gmx_ana_index_t *g) |
Frees any memory not necessary to hold the current contents. More... | |
void | gmx_ana_index_clear (gmx_ana_index_t *g) |
Initializes an empty index group. More... | |
void | gmx_ana_index_set (gmx_ana_index_t *g, int isize, int *index, int nalloc) |
Constructs a gmx_ana_index_t from given values. More... | |
void | gmx_ana_index_init_simple (gmx_ana_index_t *g, int natoms) |
Creates a simple index group from the first to the natoms'th atom. More... | |
void | gmx_ana_index_deinit (gmx_ana_index_t *g) |
Frees memory allocated for an index group. More... | |
void | gmx_ana_index_copy (gmx_ana_index_t *dest, gmx_ana_index_t *src, bool bAlloc) |
Copies a gmx_ana_index_t . More... | |
void | gmx_ana_index_dump (gmx::TextWriter *writer, gmx_ana_index_t *g, int maxn) |
Writes out the contents of a index group. More... | |
int | gmx_ana_index_get_max_index (gmx_ana_index_t *g) |
Returns maximum atom index that appears in an index group. More... | |
bool | gmx_ana_index_check_sorted (gmx_ana_index_t *g) |
Checks whether an index group is sorted. More... | |
bool | gmx_ana_index_check_range (gmx_ana_index_t *g, int natoms) |
Checks whether an index group has atoms from a defined range. More... | |
void | gmx_ana_index_sort (gmx_ana_index_t *g) |
Sorts the indices within an index group. More... | |
void | gmx_ana_index_remove_duplicates (gmx_ana_index_t *g) |
Removes duplicates from a sorted index group. More... | |
bool | gmx_ana_index_equals (gmx_ana_index_t *a, gmx_ana_index_t *b) |
Checks whether two index groups are equal. More... | |
bool | gmx_ana_index_contains (gmx_ana_index_t *a, gmx_ana_index_t *b) |
Checks whether a sorted index group contains another sorted index group. More... | |
void | gmx_ana_index_intersection (gmx_ana_index_t *dest, gmx_ana_index_t *a, gmx_ana_index_t *b) |
Calculates the intersection between two sorted index groups. More... | |
void | gmx_ana_index_difference (gmx_ana_index_t *dest, gmx_ana_index_t *a, gmx_ana_index_t *b) |
Calculates the set difference between two sorted index groups. More... | |
int | gmx_ana_index_difference_size (gmx_ana_index_t *a, gmx_ana_index_t *b) |
Calculates the size of the difference between two sorted index groups. More... | |
void | gmx_ana_index_partition (gmx_ana_index_t *dest1, gmx_ana_index_t *dest2, gmx_ana_index_t *src, gmx_ana_index_t *g) |
Calculates the intersection and the difference in one call. More... | |
void | gmx_ana_index_union (gmx_ana_index_t *dest, gmx_ana_index_t *a, gmx_ana_index_t *b) |
Calculates the union of two sorted index groups. More... | |
void | gmx_ana_index_union_unsorted (gmx_ana_index_t *dest, gmx_ana_index_t *a, gmx_ana_index_t *b) |
Calculates the union of two index groups, where the second group may not be sorted. More... | |
void | gmx_ana_index_merge (gmx_ana_index_t *dest, gmx_ana_index_t *a, gmx_ana_index_t *b) |
Merges two distinct sorted index groups. More... | |
static bool | next_group_index (int atomIndex, const gmx_mtop_t *top, e_index_t type, int *id) |
Helper for splitting a sequence of atom indices into groups. More... | |
void | gmx_ana_index_make_block (t_blocka *t, const gmx_mtop_t *top, gmx_ana_index_t *g, e_index_t type, bool bComplete) |
Partition a group based on topology information. More... | |
bool | gmx_ana_index_has_full_blocks (const gmx_ana_index_t *g, const gmx::RangePartitioning *b) |
Checks whether a group consists of full blocks. More... | |
bool | gmx_ana_index_has_full_ablocks (gmx_ana_index_t *g, t_blocka *b) |
Checks whether a group consists of full blocks. More... | |
static bool | is_at_residue_boundary (const gmx_mtop_t *top, int a, int *molb) |
Returns if an atom is at a residue boundary. More... | |
bool | gmx_ana_index_has_complete_elems (gmx_ana_index_t *g, e_index_t type, const gmx_mtop_t *top) |
Checks whether a group consists of full residues/molecules. More... | |
void | gmx_ana_indexmap_clear (gmx_ana_indexmap_t *m) |
Initializes an empty index group mapping. More... | |
void | gmx_ana_indexmap_reserve (gmx_ana_indexmap_t *m, int nr, int isize) |
Reserves memory for an index group mapping. More... | |
void | gmx_ana_indexmap_init (gmx_ana_indexmap_t *m, gmx_ana_index_t *g, const gmx_mtop_t *top, e_index_t type) |
Initializes an index group mapping. More... | |
int | gmx_ana_indexmap_init_orgid_group (gmx_ana_indexmap_t *m, const gmx_mtop_t *top, e_index_t type) |
Initializes orgid entries based on topology grouping. More... | |
void | gmx_ana_indexmap_set_static (gmx_ana_indexmap_t *m, t_blocka *b) |
Sets an index group mapping to be static. More... | |
void | gmx_ana_indexmap_copy (gmx_ana_indexmap_t *dest, gmx_ana_indexmap_t *src, bool bFirst) |
Makes a deep copy of an index group mapping. More... | |
static void | set_atoms (gmx_ana_indexmap_t *m, int isize, int *index) |
Helper function to set the source atoms in an index map. More... | |
void | gmx_ana_indexmap_update (gmx_ana_indexmap_t *m, gmx_ana_index_t *g, bool bMaskOnly) |
Updates an index group mapping. More... | |
void | gmx_ana_indexmap_deinit (gmx_ana_indexmap_t *m) |
Frees memory allocated for index group mapping. More... | |
bool gmx_ana_index_check_range | ( | gmx_ana_index_t * | g, |
int | natoms | ||
) |
Checks whether an index group has atoms from a defined range.
[in] | g | Index group to check. |
[in] | natoms | Largest atom number allowed. |
natoms
(i.e., no atoms over natoms
are referenced). bool gmx_ana_index_check_sorted | ( | gmx_ana_index_t * | g | ) |
Checks whether an index group is sorted.
[in] | g | Index group to check. |
void gmx_ana_index_clear | ( | gmx_ana_index_t * | g | ) |
Initializes an empty index group.
[out] | g | Output structure. |
Any contents of g
are discarded without freeing.
bool gmx_ana_index_contains | ( | gmx_ana_index_t * | a, |
gmx_ana_index_t * | b | ||
) |
Checks whether a sorted index group contains another sorted index group.
[in] | a | Index group to check against. |
[in] | b | Index group to check. |
b
is contained in a
, false otherwise.If the elements are not in the same order in both groups, the function fails. However, the groups do not need to be sorted.
void gmx_ana_index_copy | ( | gmx_ana_index_t * | dest, |
gmx_ana_index_t * | src, | ||
bool | bAlloc | ||
) |
Copies a gmx_ana_index_t
.
[out] | dest | Destination index group. |
[in] | src | Source index group. |
[in] | bAlloc | If true, memory is allocated at dest ; otherwise, it is assumed that enough memory has been allocated for index. |
void gmx_ana_index_deinit | ( | gmx_ana_index_t * | g | ) |
Frees memory allocated for an index group.
[in] | g | Index group structure. |
The pointer g
is not freed.
void gmx_ana_index_difference | ( | gmx_ana_index_t * | dest, |
gmx_ana_index_t * | a, | ||
gmx_ana_index_t * | b | ||
) |
Calculates the set difference between two sorted index groups.
[out] | dest | Output index group (the difference a - b ). |
[in] | a | First index group. |
[in] | b | Second index group. |
dest
can equal a
, but not b
.
int gmx_ana_index_difference_size | ( | gmx_ana_index_t * | a, |
gmx_ana_index_t * | b | ||
) |
Calculates the size of the difference between two sorted index groups.
[in] | a | First index group. |
[in] | b | Second index group. |
a
- b
. void gmx_ana_index_dump | ( | gmx::TextWriter * | writer, |
gmx_ana_index_t * | g, | ||
int | maxn | ||
) |
Writes out the contents of a index group.
[in] | writer | Writer to use for output. |
[in] | g | Index group to print. |
[in] | maxn | Maximum number of indices to print (-1 = print all). |
bool gmx_ana_index_equals | ( | gmx_ana_index_t * | a, |
gmx_ana_index_t * | b | ||
) |
Checks whether two index groups are equal.
[in] | a | Index group to check. |
[in] | b | Index group to check. |
a
and b
are equal, false otherwise. int gmx_ana_index_get_max_index | ( | gmx_ana_index_t * | g | ) |
Returns maximum atom index that appears in an index group.
[in] | g | Index group to query. |
g
, or zero if g
is empty. bool gmx_ana_index_has_complete_elems | ( | gmx_ana_index_t * | g, |
e_index_t | type, | ||
const gmx_mtop_t * | top | ||
) |
Checks whether a group consists of full residues/molecules.
[in] | g | Index group to check. |
[in] | type | Block data to check against. |
[in] | top | Topology data. |
g
consists of one or more complete elements of type type
, false otherwise.g
is assumed to be sorted, otherwise may return false negatives.
If type
is INDEX_ATOM, the return value is always true. If type
is INDEX_UNKNOWN or INDEX_ALL, the return value is always false.
bool gmx_ana_index_has_full_ablocks | ( | gmx_ana_index_t * | g, |
t_blocka * | b | ||
) |
Checks whether a group consists of full blocks.
[in] | g | Index group to check. |
[in] | b | Block data to check against. |
g
consists of one or more complete blocks from b
, false otherwise.The atoms in g
and b->a
are assumed to be in the same order.
bool gmx_ana_index_has_full_blocks | ( | const gmx_ana_index_t * | g, |
const gmx::RangePartitioning * | b | ||
) |
Checks whether a group consists of full blocks.
[in] | g | Index group to check. |
[in] | b | Block data to check against. |
g
consists of one or more complete blocks from b
, false otherwise.The atoms in g
are assumed to be sorted.
void gmx_ana_index_init_simple | ( | gmx_ana_index_t * | g, |
int | natoms | ||
) |
Creates a simple index group from the first to the natoms'th
atom.
[out] | g | Output structure. |
[in] | natoms | Number of atoms. |
void gmx_ana_index_intersection | ( | gmx_ana_index_t * | dest, |
gmx_ana_index_t * | a, | ||
gmx_ana_index_t * | b | ||
) |
Calculates the intersection between two sorted index groups.
[out] | dest | Output index group (the intersection of a and b ). |
[in] | a | First index group. |
[in] | b | Second index group. |
dest
can be the same as a
or b
.
void gmx_ana_index_make_block | ( | t_blocka * | t, |
const gmx_mtop_t * | top, | ||
gmx_ana_index_t * | g, | ||
e_index_t | type, | ||
bool | bComplete | ||
) |
Partition a group based on topology information.
[in,out] | t | Output block. |
[in] | top | Topology structure (only used if type is INDEX_RES or INDEX_MOL, can be NULL otherwise). |
[in] | g | Index group (can be NULL if type is INDEX_UNKNOWN). |
[in] | type | Type of partitioning to make. |
[in] | bComplete | If true, the index group is expanded to include any residue/molecule (depending on type ) that is partially contained in the group. If type is not INDEX_RES or INDEX_MOL, this has no effect. |
m
should have been initialized somehow (calloc() is enough). g
should be sorted.
void gmx_ana_index_merge | ( | gmx_ana_index_t * | dest, |
gmx_ana_index_t * | a, | ||
gmx_ana_index_t * | b | ||
) |
Merges two distinct sorted index groups.
[out] | dest | Output index group (the union of a and b ). |
[in] | a | First index group. |
[in] | b | Second index group. |
a
and b
should not have common items. dest
can equal a
or b
.
void gmx_ana_index_partition | ( | gmx_ana_index_t * | dest1, |
gmx_ana_index_t * | dest2, | ||
gmx_ana_index_t * | src, | ||
gmx_ana_index_t * | g | ||
) |
Calculates the intersection and the difference in one call.
[out] | dest1 | Output group 1 (will equal g ). |
[out] | dest2 | Output group 2 (will equal src - g ). |
[in] | src | Group to be partitioned. |
[in] | g | One partition. |
g
is a subset of src
and both sets are sorted dest1
has allocated storage to store src
dest1
== g
dest2
== src
- g
No storage should be allocated for dest2
; after the call, dest2->index
points to the memory allocated for dest1
(to a part that is not used by dest1
).
The calculation can be performed in-place by setting dest1
equal to src
.
void gmx_ana_index_remove_duplicates | ( | gmx_ana_index_t * | g | ) |
Removes duplicates from a sorted index group.
[in,out] | g | Index group to be processed. |
void gmx_ana_index_reserve | ( | gmx_ana_index_t * | g, |
int | isize | ||
) |
Reserves memory to store an index group of size isize
.
[in,out] | g | Index group structure. |
[in] | isize | Maximum number of atoms to reserve space for. |
void gmx_ana_index_set | ( | gmx_ana_index_t * | g, |
int | isize, | ||
int * | index, | ||
int | nalloc | ||
) |
Constructs a gmx_ana_index_t
from given values.
[out] | g | Output structure. |
[in] | isize | Number of atoms in the new group. |
[in] | index | Array of isize atoms (can be NULL if isize is 0). |
[in] | nalloc | Number of elements allocated for index (if 0, index is not freed in gmx_ana_index_deinit()) |
No copy if index
is made.
void gmx_ana_index_sort | ( | gmx_ana_index_t * | g | ) |
Sorts the indices within an index group.
[in,out] | g | Index group to be sorted. |
void gmx_ana_index_squeeze | ( | gmx_ana_index_t * | g | ) |
Frees any memory not necessary to hold the current contents.
[in,out] | g | Index group structure. |
Resizes the memory allocated for holding the indices such that the current contents fit.
void gmx_ana_index_union | ( | gmx_ana_index_t * | dest, |
gmx_ana_index_t * | a, | ||
gmx_ana_index_t * | b | ||
) |
Calculates the union of two sorted index groups.
[out] | dest | Output index group (the union of a and b ). |
[in] | a | First index group. |
[in] | b | Second index group. |
a
and b
can have common items. dest
can equal a
or b
.
void gmx_ana_index_union_unsorted | ( | gmx_ana_index_t * | dest, |
gmx_ana_index_t * | a, | ||
gmx_ana_index_t * | b | ||
) |
Calculates the union of two index groups, where the second group may not be sorted.
[out] | dest | Output index group (the union of a and b ). |
[in] | a | First index group (must be sorted). |
[in] | b | Second index group. |
a
and b
can have common items. dest
can equal a
or b
.
bool gmx_ana_indexgrps_extract | ( | gmx_ana_index_t * | dest, |
std::string * | destName, | ||
gmx_ana_indexgrps_t * | src, | ||
int | n | ||
) |
Extracts a single index group.
[out] | dest | Output structure. |
[out] | destName | Receives the name of the group if found. |
[in] | src | Input index groups. |
[in] | n | Number of the group to extract. |
n
is a valid group in src
, false otherwise. bool gmx_ana_indexgrps_find | ( | gmx_ana_index_t * | dest, |
std::string * | destName, | ||
gmx_ana_indexgrps_t * | src, | ||
const char * | name | ||
) |
Finds and extracts a single index group by name.
[out] | dest | Output structure. |
[out] | destName | Receives the name of the group if found. |
[in] | src | Input index groups. |
[in] | name | Name (or part of the name) of the group to extract. |
name
is a valid group in src
, false otherwise.Uses the Gromacs routine find_group() to find the actual group; the comparison is case-insensitive.
void gmx_ana_indexgrps_free | ( | gmx_ana_indexgrps_t * | g | ) |
Frees memory allocated for index groups.
[in] | g | Index groups structure. |
The pointer g
is invalid after the call.
void gmx_ana_indexgrps_init | ( | gmx_ana_indexgrps_t ** | g, |
gmx_mtop_t * | top, | ||
const char * | fnm | ||
) |
Reads index groups from a file or constructs them from topology.
[out] | g | Index group structure. |
[in] | top | Topology structure. |
[in] | fnm | File name for the index file. Memory is automatically allocated. |
One or both of top
or fnm
can be NULL. If top
is NULL, an index file is required and the groups are read from the file (uses Gromacs routine init_index()). If fnm
is NULL, default groups are constructed based on the topology (uses Gromacs routine analyse()). If both are null, the index group structure is initialized empty.
void gmx_ana_indexgrps_print | ( | gmx::TextWriter * | writer, |
gmx_ana_indexgrps_t * | g, | ||
int | maxn | ||
) |
Writes out a list of index groups.
[in] | writer | Writer to use for output. |
[in] | g | Index groups to print. |
[in] | maxn | Maximum number of indices to print (-1 = print all, 0 = print only names). |
void gmx_ana_indexmap_clear | ( | gmx_ana_indexmap_t * | m | ) |
Initializes an empty index group mapping.
[out] | m | Output structure. |
Any contents of m
are discarded without freeing.
void gmx_ana_indexmap_copy | ( | gmx_ana_indexmap_t * | dest, |
gmx_ana_indexmap_t * | src, | ||
bool | bFirst | ||
) |
Makes a deep copy of an index group mapping.
[in,out] | dest | Destination data structure. |
[in] | src | Source mapping. |
[in] | bFirst | If true, memory is allocated for dest and a full copy is made; otherwise, only variable parts are copied, and no memory is allocated. |
dest
should have been initialized somehow (calloc() is enough).
void gmx_ana_indexmap_deinit | ( | gmx_ana_indexmap_t * | m | ) |
Frees memory allocated for index group mapping.
[in,out] | m | Mapping structure to free. |
All the memory allocated for the mapping structure is freed, and the pointers set to NULL. The pointer m
is not freed.
void gmx_ana_indexmap_init | ( | gmx_ana_indexmap_t * | m, |
gmx_ana_index_t * | g, | ||
const gmx_mtop_t * | top, | ||
e_index_t | type | ||
) |
Initializes an index group mapping.
[in,out] | m | Mapping structure to initialize. |
[in] | g | Index group to map (can be NULL if type is INDEX_UNKNOWN). |
[in] | top | Topology structure (can be NULL if type is not INDEX_RES or INDEX_MOL). |
[in] | type | Type of mapping to construct. |
Initializes a new index group mapping. The index group provided to gmx_ana_indexmap_update() should always be a subset of the g
given here.
m
should have been initialized somehow (calloc() is enough).
int gmx_ana_indexmap_init_orgid_group | ( | gmx_ana_indexmap_t * | m, |
const gmx_mtop_t * | top, | ||
e_index_t | type | ||
) |
Initializes orgid
entries based on topology grouping.
[in,out] | m | Mapping structure to use/initialize. |
[in] | top | Topology structure (can be NULL if not required for type ). |
[in] | type | Type of groups to use. |
type
that were present in m
. InconsistentInputError | if the blocks in m do not have a unique group (e.g., contain atoms from multiple residues with type == INDEX_RES ). |
By default, the gmx_ana_indexmap_t::orgid fields are initialized to atom/residue/molecule indices from the topology (see documentation for the struct for more details). This function can be used to set the field to a zero-based group index instead. The first block will always get orgid[0] = 0
, and all following blocks that belong to the same residue/molecule (type
) will get the same index. Each time a block does not belong to the same group, it will get the next available number. If type == INDEX_ATOM
, the orgid
field is initialized as 0, 1, ..., independent of whether the blocks are single atoms or not.
Strong exception safety guarantee.
void gmx_ana_indexmap_reserve | ( | gmx_ana_indexmap_t * | m, |
int | nr, | ||
int | isize | ||
) |
Reserves memory for an index group mapping.
[in,out] | m | Mapping structure. |
[in] | nr | Maximum number of blocks to reserve space for. |
[in] | isize | Maximum number of atoms to reserve space for. |
void gmx_ana_indexmap_set_static | ( | gmx_ana_indexmap_t * | m, |
t_blocka * | b | ||
) |
Sets an index group mapping to be static.
[in,out] | m | Mapping structure to initialize. |
[in] | b | Block information to use for data. |
Frees some memory that is not necessary for static index group mappings. Internal pointers are set to point to data in b
; it is the responsibility of the caller to ensure that the block information matches the contents of the mapping. After this function has been called, the index group provided to gmx_ana_indexmap_update() should always be the same as g
given here.
This function breaks modularity of the index group mapping interface in an ugly way, but allows reducing memory usage of static selections by a significant amount.
void gmx_ana_indexmap_update | ( | gmx_ana_indexmap_t * | m, |
gmx_ana_index_t * | g, | ||
bool | bMaskOnly | ||
) |
Updates an index group mapping.
[in,out] | m | Mapping structure. |
[in] | g | Current index group. |
[in] | bMaskOnly | true if the unused blocks should be masked with -1 instead of removing them. |
Updates the index group mapping with the new index group g
.
|
static |
Returns if an atom is at a residue boundary.
[in] | top | Topology data. |
[in] | a | Atom index to check, should be -1 <= a < top->natoms. |
[in,out] | molb | The molecule block of atom a |
a
and a
+ 1 are in different residues, false otherwise.
|
static |
Helper function to set the source atoms in an index map.
[in,out] | m | Mapping structure. |
[in] | isize | Number of atoms in the index array. |
[in] | index | List of atoms. |