Gromacs  2020.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Enumerations
#include <cstdio>
#include <string>
#include <vector>
#include "gromacs/topology/block.h"
#include "gromacs/utility/arrayref.h"
+ Include dependency graph for indexutil.h:
+ This graph shows which files directly or indirectly include this file:

Description

API for handling index files and index groups.

The API contains functions and data structures for handling index files more conveniently than as several separate variables. In addition to basic functions for initializing the data structures and making copies, functions are provided for performing (most) set operations on sorted index groups. There is also a function for partitioning a index group based on topology information such as residues or molecules. Finally, there is a set of functions for constructing mappings between an index group and its subgroups such. These can be used with dynamic index group in calculations if one needs to have a unique ID for each possible atom/residue/molecule in the selection, e.g., for analysis of dynamics or for look-up tables.

Mostly, these functions are used internally by the selection engine, but it is necessary to use some of these functions in order to provide external index groups to a gmx::SelectionCollection. Some of the checking functions can be useful outside the selection engine to check the validity of input groups.

Author
Teemu Murtola teemu.nosp@m..mur.nosp@m.tola@.nosp@m.gmai.nosp@m.l.com

Classes

class  gmx::IndexGroupsAndNames
 Bundle index groups with their names. More...
 
struct  gmx_ana_index_t
 Stores a single index group. More...
 
struct  gmx_ana_indexmap_t
 Data structure for calculating index group mappings. More...
 

Enumerations

enum  e_index_t {
  INDEX_UNKNOWN, INDEX_ATOM, INDEX_RES, INDEX_MOL,
  INDEX_ALL
}
 Specifies the type of index partition or index mapping in several contexts. More...
 

Functions

Functions for handling gmx_ana_indexgrps_t
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_is_empty (gmx_ana_indexgrps_t *g)
 Returns true if the index group structure is emtpy. More...
 
gmx_ana_index_tgmx_ana_indexgrps_get_grp (gmx_ana_indexgrps_t *g, int n)
 Returns a pointer to an index group. 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...
 
Functions for handling gmx_ana_index_t
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...
 
Functions for set operations on gmx_ana_index_t
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_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...
 
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...
 
Functions for handling gmx_ana_indexmap_t and related things
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...
 
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_deinit (gmx_ana_indexmap_t *m)
 Frees memory allocated for index group mapping. 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...
 
void gmx_ana_indexmap_update (gmx_ana_indexmap_t *m, gmx_ana_index_t *g, bool bMaskOnly)
 Updates an index group mapping. More...
 

Enumeration Type Documentation

enum e_index_t

Specifies the type of index partition or index mapping in several contexts.

See Also
gmx_ana_index_make_block(), gmx_ana_indexmap_init()
Enumerator
INDEX_UNKNOWN 

Unknown index type.

INDEX_ATOM 

Each atom in a separate block.

INDEX_RES 

Each residue in a separate block.

INDEX_MOL 

Each molecule in a separate block.

INDEX_ALL 

All atoms in a single block.

Function Documentation

bool gmx_ana_index_check_range ( gmx_ana_index_t g,
int  natoms 
)

Checks whether an index group has atoms from a defined range.

Parameters
[in]gIndex group to check.
[in]natomsLargest atom number allowed.
Returns
true if all atoms in the index group are in the range 0 to 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.

Parameters
[in]gIndex group to check.
Returns
true if the index group is sorted and has no duplicates, false otherwise.
void gmx_ana_index_clear ( gmx_ana_index_t g)

Initializes an empty index group.

Parameters
[out]gOutput 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.

Parameters
[in]aIndex group to check against.
[in]bIndex group to check.
Returns
true if 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.

Parameters
[out]destDestination index group.
[in]srcSource index group.
[in]bAllocIf 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.

Parameters
[in]gIndex 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.

Parameters
[out]destOutput index group (the difference a - b).
[in]aFirst index group.
[in]bSecond 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.

Parameters
[in]aFirst index group.
[in]bSecond index group.
Returns
Size of the difference 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.

Parameters
[in]writerWriter to use for output.
[in]gIndex group to print.
[in]maxnMaximum 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.

Parameters
[in]aIndex group to check.
[in]bIndex group to check.
Returns
true if 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.

Parameters
[in]gIndex group to query.
Returns
Largest atom index that appears in 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.

Parameters
[in]gIndex group to check.
[in]typeBlock data to check against.
[in]topTopology data.
Returns
true if 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.

Parameters
[in]gIndex group to check.
[in]bBlock data to check against.
Returns
true if 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.

Parameters
[in]gIndex group to check.
[in]bBlock data to check against.
Returns
true if 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.

Parameters
[out]gOutput structure.
[in]natomsNumber 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.

Parameters
[out]destOutput index group (the intersection of a and b).
[in]aFirst index group.
[in]bSecond 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.

Parameters
[in,out]tOutput block.
[in]topTopology structure (only used if type is INDEX_RES or INDEX_MOL, can be NULL otherwise).
[in]gIndex group (can be NULL if type is INDEX_UNKNOWN).
[in]typeType of partitioning to make.
[in]bCompleteIf 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.

Parameters
[out]destOutput index group (the union of a and b).
[in]aFirst index group.
[in]bSecond index group.

a and b should not have common items. dest can equal a or b.

See Also
gmx_ana_index_union()
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.

Parameters
[out]dest1Output group 1 (will equal g).
[out]dest2Output group 2 (will equal src - g).
[in]srcGroup to be partitioned.
[in]gOne partition.
Precondition
g is a subset of src and both sets are sorted
dest1 has allocated storage to store src
Postcondition
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.

Parameters
[in,out]gIndex 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.

Parameters
[in,out]gIndex group structure.
[in]isizeMaximum 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.

Parameters
[out]gOutput structure.
[in]isizeNumber of atoms in the new group.
[in]indexArray of isize atoms (can be NULL if isize is 0).
[in]nallocNumber 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.

Parameters
[in,out]gIndex group to be sorted.
void gmx_ana_index_squeeze ( gmx_ana_index_t g)

Frees any memory not necessary to hold the current contents.

Parameters
[in,out]gIndex 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.

Parameters
[out]destOutput index group (the union of a and b).
[in]aFirst index group.
[in]bSecond index group.

a and b can have common items. dest can equal a or b.

See Also
gmx_ana_index_merge()
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.

Parameters
[out]destOutput index group (the union of a and b).
[in]aFirst index group (must be sorted).
[in]bSecond 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.

Parameters
[out]destOutput structure.
[out]destNameReceives the name of the group if found.
[in]srcInput index groups.
[in]nNumber of the group to extract.
Returns
true if 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.

Parameters
[out]destOutput structure.
[out]destNameReceives the name of the group if found.
[in]srcInput index groups.
[in]nameName (or part of the name) of the group to extract.
Returns
true if 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.

Parameters
[in]gIndex groups structure.

The pointer g is invalid after the call.

gmx_ana_index_t* gmx_ana_indexgrps_get_grp ( gmx_ana_indexgrps_t g,
int  n 
)

Returns a pointer to an index group.

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.

Parameters
[out]gIndex group structure.
[in]topTopology structure.
[in]fnmFile 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.

bool gmx_ana_indexgrps_is_empty ( gmx_ana_indexgrps_t g)

Returns true if the index group structure is emtpy.

void gmx_ana_indexgrps_print ( gmx::TextWriter writer,
gmx_ana_indexgrps_t g,
int  maxn 
)

Writes out a list of index groups.

Parameters
[in]writerWriter to use for output.
[in]gIndex groups to print.
[in]maxnMaximum 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.

Parameters
[out]mOutput 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.

Parameters
[in,out]destDestination data structure.
[in]srcSource mapping.
[in]bFirstIf 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.

Parameters
[in,out]mMapping 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.

Parameters
[in,out]mMapping structure to initialize.
[in]gIndex group to map (can be NULL if type is INDEX_UNKNOWN).
[in]topTopology structure (can be NULL if type is not INDEX_RES or INDEX_MOL).
[in]typeType 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.

Parameters
[in,out]mMapping structure to use/initialize.
[in]topTopology structure (can be NULL if not required for type).
[in]typeType of groups to use.
Returns
The number of groups of type type that were present in m.
Exceptions
InconsistentInputErrorif 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.

Parameters
[in,out]mMapping structure.
[in]nrMaximum number of blocks to reserve space for.
[in]isizeMaximum 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.

Parameters
[in,out]mMapping structure to initialize.
[in]bBlock 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.

Parameters
[in,out]mMapping structure.
[in]gCurrent index group.
[in]bMaskOnlytrue if the unused blocks should be masked with -1 instead of removing them.

Updates the index group mapping with the new index group g.

See Also
gmx_ana_indexmap_t