Gromacs  2021-beta3-UNCHECKED
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Functions
groupcoord.h File Reference
#include <stdio.h>
#include "gromacs/math/vectypes.h"
#include "gromacs/utility/basedefinitions.h"
+ Include dependency graph for groupcoord.h:

Description

Assemble atomic positions of a (small) subset of atoms and distribute to all nodes.

This file contains functions to assemble the positions of a subset of the atoms and to do operations on it like determining the center of mass, or doing translations and rotations. These functions are useful when a subset of the positions needs to be compared to some set of reference positions, as e.g. done for essential dynamics.

Functions

void dd_make_local_group_indices (const gmx_ga2la_t *ga2la, int nr, int anrs[], int *nr_loc, int *anrs_loc[], int *nalloc_loc, int coll_ind[])
 Select local atoms of a group. More...
 
void communicate_group_positions (const t_commrec *cr, rvec *xcoll, ivec *shifts, ivec *extra_shifts, gmx_bool bNS, const rvec *x_loc, int nr, int nr_loc, const int *anrs_loc, const int *coll_ind, rvec *xcoll_old, const matrix box)
 Assemble local positions into a collective array present on all nodes. More...
 
void get_center (rvec x[], real weight[], int nr, rvec center)
 Calculates the center of the positions x locally. More...
 
double get_sum_of_positions (rvec x[], real weight[], int nr, dvec dsumvec)
 Calculates the sum of the positions x locally. More...
 
void get_center_comm (const t_commrec *cr, rvec x_loc[], real weight_loc[], int nr_loc, int nr_group, rvec center)
 Calculates the global center of all local arrays x. More...
 
void translate_x (rvec x[], int nr, const rvec transvec)
 Translate positions. More...
 
void rotate_x (rvec x[], int nr, matrix rmat)
 Rotate positions. More...
 

Function Documentation

void communicate_group_positions ( const t_commrec *  cr,
rvec *  xcoll,
ivec *  shifts,
ivec *  extra_shifts,
gmx_bool  bNS,
const rvec *  x_loc,
int  nr,
int  nr_loc,
const int *  anrs_loc,
const int *  coll_ind,
rvec *  xcoll_old,
const matrix  box 
)

Assemble local positions into a collective array present on all nodes.

Communicate the positions of the group's atoms such that every node has all of them. Unless running on huge number of cores, this is not a big performance impact as long as the collective subset [0..nr] is kept small. The atom indices are retrieved from anrs_loc[0..nr_loc]. If you call the routine for the serial case, provide an array coll_ind[i] = i for i in 1..nr.

If shifts != NULL, the PBC representation of each atom is chosen such that a continuous trajectory results. Therefore, if the group is whole at the start of the simulation, it will always stay whole. If shifts = NULL, the group positions are not made whole again, but assembled and distributed to all nodes. The variables marked "optional" are not used in that case.

Parameters
[in]crPointer to MPI communication data.
[out]xcollCollective array of positions, identical on all nodes after this routine has been called.
[in,out]shiftsCollective array of shifts for xcoll, needed to make the group whole. This array remembers the shifts since the start of the simulation (where the group is whole) and must therefore not be changed outside of this routine! If NULL, the group will not be made whole and the optional variables are ignored.
[out]extra_shiftsExtra shifts since last time step, only needed as buffer variable [0..nr] (optional).
[in]bNSNeighbor searching / domain re-decomposition has been performed at the begin of this time step such that the shifts have changed and need to be updated (optional).
[in]x_locPointer to the local atom positions this node has.
[in]nrTotal number of atoms in the group.
[in]nr_locNumber of group atoms on the local node.
[in]anrs_locArray of the local atom indices.
[in]coll_indThis array of size nr stores for each local atom where it belongs in the collective array so that the local contributions can be gmx_summed. It is provided by dd_make_local_group_indices.
[in,out]xcoll_oldPositions from the last time step, used to make the group whole (optional).
[in]boxSimulation box matrix, needed to shift xcoll such that the group becomes whole (optional).
void dd_make_local_group_indices ( const gmx_ga2la_t ga2la,
int  nr,
int  anrs[],
int *  nr_loc,
int *  anrs_loc[],
int *  nalloc_loc,
int  coll_ind[] 
)

Select local atoms of a group.

Selects the indices of local atoms of a group and stores them in anrs_loc[0..nr_loc]. If you need the positions of the group's atoms on all nodes, provide a coll_ind[0..nr] array and pass it on to communicate_group_positions. Thus the collective array will always have the same atom order (ascending indices).

Parameters
[in]ga2laGlobal to local atom index conversion data.
[in]nrThe total number of atoms that the group contains.
[in]anrsThe global atom number of the group's atoms.
[out]nr_locThe number of group atoms present on the local node.
[out]anrs_locThe local atom numbers of the group.
[in,out]nalloc_locLocal allocation size of anrs_loc array.
[out]coll_indIf not NULL this array must be of size nr. It stores for each local atom where it belongs in the global (collective) array such that it can be gmx_summed in the communicate_group_positions routine.
void get_center ( rvec  x[],
real  weight[],
int  nr,
rvec  center 
)

Calculates the center of the positions x locally.

Calculates the center of mass (if masses are given in the weight array) or the geometrical center (if NULL is passed as weight).

Parameters
[in]xPositions.
[in]weightCan be NULL or an array of weights. If masses are given as weights, the COM is calculated.
[in]nrNumber of positions and weights if present.
[out]centerThe (weighted) center of the positions.
void get_center_comm ( const t_commrec *  cr,
rvec  x_loc[],
real  weight_loc[],
int  nr_loc,
int  nr_group,
rvec  center 
)

Calculates the global center of all local arrays x.

Get the center from local positions [0..nr_loc], this involves communication. Not that the positions must already have the correct PBC representation. Use this routine if no collective coordinates are assembled from which the center could be calculated without communication.

Parameters
[in]crPointer to MPI communication data.
[in]x_locArray of local positions [0..nr_loc].
[in]weight_locArray of local weights, these are the masses if the center of mass is to be calculated.
[in]nr_locThe number of positions on the local node.
[in]nr_groupThe number of positions in the whole group. Since this is known anyway, we do not need to communicate and sum nr_loc if we pass it over.
[out]centerThe (weighted) center of all x_loc from all the nodes.
double get_sum_of_positions ( rvec  x[],
real  weight[],
int  nr,
dvec  dsumvec 
)

Calculates the sum of the positions x locally.

Calculates the (weighted) sum of position vectors and returns the sum of weights, which is needed when local contributions shall be summed to a global weighted center.

Parameters
[in]xArray of positions.
[in]weightCan be NULL or an array of weights.
[in]nrNumber of positions and weights if present.
[out]dsumvecThe (weighted) sum of the positions.
Returns
Sum of weights.
void rotate_x ( rvec  x[],
int  nr,
matrix  rmat 
)

Rotate positions.

Rotate the positions with the rotation matrix.

Parameters
[in,out]xArray of positions.
[in]nrNumber of entries in the position array.
[in]rmatRotation matrix to operate on all positions.
void translate_x ( rvec  x[],
int  nr,
const rvec  transvec 
)

Translate positions.

Add a translation vector to the positions x.

Parameters
[in,out]xArray of positions.
[in]nrNumber of entries in the position array.
[in]transvecTranslation vector to be added to all positions.