Gromacs  2026.0-dev-20241212-74b8831
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Public Member Functions
gmx::UpdateGroupsCog Class Reference

#include <gromacs/mdlib/updategroupscog.h>

Description

Class for managing and computing centers of geometry of update groups.

Public Member Functions

 UpdateGroupsCog (const gmx_mtop_t &mtop, gmx::ArrayRef< const gmx::RangePartitioning > updateGroupingsPerMoleculeType, real temperature)
 Constructor. More...
 
void addCogs (ArrayRef< const int > globalAtomIndices, ArrayRef< const RVec > coordinates, ArrayRef< const int > numAtomsPerNbnxmGridColumn)
 Compute centers of geometry for supplied coordinates. More...
 
int numCogs () const
 Returns the number of centers of geometry currently stored.
 
RVeccog (int cogIndex)
 Returns a reference to a center of geometry. More...
 
bool cogIsFillerParticle (int cogIndex) const
 Returns whether a COG entry consists of a (single) filler particle. More...
 
int cogIndex (int atomIndex) const
 Returns the COG index given an atom index. More...
 
const RVeccogForAtom (int atomIndex) const
 Return a reference to a center of geometry given an atom index. More...
 
void clear ()
 Clear the lists of stored center of geometry coordinates.
 
real maxUpdateGroupRadius () const
 Returns the maximum radius over all update groups.
 

Constructor & Destructor Documentation

gmx::UpdateGroupsCog::UpdateGroupsCog ( const gmx_mtop_t &  mtop,
gmx::ArrayRef< const gmx::RangePartitioning updateGroupingsPerMoleculeType,
real  temperature 
)

Constructor.

The temperature is used for computing the maximum update group radius.

Parameters
[in]mtopThe global topology
[in]updateGroupingsPerMoleculeTypeList of update groups for each molecule type in mtop
[in]temperatureThe maximum reference temperature, pass -1 when unknown or not applicable

Member Function Documentation

void gmx::UpdateGroupsCog::addCogs ( ArrayRef< const int >  globalAtomIndices,
ArrayRef< const RVec coordinates,
ArrayRef< const int >  numAtomsPerNbnxmGridColumn 
)

Compute centers of geometry for supplied coordinates.

Coordinates are processed starting after the last index processed in the previous call to addCogs(), unless clear() was called last, in which case processing starts at 0. Coordinates are processed until globalAtomIndices.size().

Parameters
[in]globalAtomIndicesList of global atom indices for the atoms belonging to coordinates
[in]coordinatesList of coordinates to be processed, processing might not start at 0 (see above)
[in]numAtomsPerNbnxmGridColumnA list of the number of atoms per NBNxM grid column, used for OpenMP parallelization, as update groups are never split over columns; can be empty
RVec& gmx::UpdateGroupsCog::cog ( int  cogIndex)
inline

Returns a reference to a center of geometry.

Parameters
[in]cogIndexThe COG requested, should be; 0 <= cogIndex < cogs_.size()
const RVec& gmx::UpdateGroupsCog::cogForAtom ( int  atomIndex) const
inline

Return a reference to a center of geometry given an atom index.

Parameters
[in]atomIndexThe atom index that maps to the COG requested
int gmx::UpdateGroupsCog::cogIndex ( int  atomIndex) const
inline

Returns the COG index given an atom index.

Parameters
[in]atomIndexThe local atom index that maps to the COG
bool gmx::UpdateGroupsCog::cogIsFillerParticle ( int  cogIndex) const
inline

Returns whether a COG entry consists of a (single) filler particle.

Parameters
[in]cogIndexThe COG requested, should be; 0 <= cogIndex < cogs_.size()

The documentation for this class was generated from the following files: