Gromacs  2025-dev-20241003-bd59e46
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Public Member Functions | Static Public Member Functions | Friends
gmx::ColvarsForceProvider Class Referencefinal

#include <gromacs/applied_forces/colvars/colvarsforceprovider.h>

+ Inheritance diagram for gmx::ColvarsForceProvider:
+ Collaboration diagram for gmx::ColvarsForceProvider:

Description

Implements IForceProvider for colvars. Override the ColvarProxyGromacs generic class for the communication.

Public Member Functions

 ColvarsForceProvider (const std::string &colvarsConfigString, t_atoms atoms, PbcType pbcType, const MDLogger *logger, const std::map< std::string, std::string > &inputStrings, real ensembleTemperature, int seed, LocalAtomSetManager *localAtomSetManager, const t_commrec *cr, double simulationTimeStep, const std::vector< RVec > &colvarsCoords, const std::string &outputPrefix, const ColvarsForceProviderState &state)
 Construct ColvarsForceProvider from its parameters. More...
 
void calculateForces (const ForceProviderInput &forceProviderInput, ForceProviderOutput *forceProviderOutput) override
 Calculate colvars forces. More...
 
void writeCheckpointData (MDModulesWriteCheckpointData checkpointWriting, const std::string &moduleName)
 Write internal colvars data to checkpoint file. More...
 
void processAtomsRedistributedSignal (const MDModulesAtomsRedistributedSignal &atomsRedistributedSignal)
 Process atomsRedistributedSignal notification during mdrun. More...
 
void add_energy (cvm::real energy) override
 From colvarproxy. More...
 
- Public Member Functions inherited from gmx::ColvarProxyGromacs
 ColvarProxyGromacs (const std::string &colvarsConfigString, t_atoms atoms, PbcType pbcType, const MDLogger *logger, bool doParsing, const std::map< std::string, std::string > &inputStrings, real ensembleTemperature, int seed)
 Construct ColvarProxyGromacs from its parameters. More...
 
void updateAtomProperties (int index)
 Update colvars topology of one atom mass and charge from the GROMACS topology.
 
cvm::real rand_gaussian () override
 From colvarproxy. More...
 
void log (std::string const &message) override
 Print a message to the main log.
 
void error (std::string const &message) override
 Print a message to the main log and let GROMACS handle the error.
 
int backup_file (char const *filename) override
 Rename the given file, before overwriting it.
 
int set_unit_system (std::string const &unitsIn, bool colvarsDefined) override
 Request to set the units used internally by Colvars.
 
int init_atom (int atomNumber) override
 Initialize colvars atom from GROMACS topology.
 
int check_atom_id (int atomNumber) override
 Check if atom belongs to the global index of atoms. More...
 
cvm::rvector position_distance (cvm::atom_pos const &pos1, cvm::atom_pos const &pos2) const override
 Compute the minimum distance with respect to the PBC between 2 atoms.
 

Static Public Member Functions

static void addVirialTerm (matrix vir, const rvec &f, const gmx::RVec &x)
 Compute virial tensor for position r and force f, and add to matrix vir.
 

Friends

class cvm::atom
 

Additional Inherited Members

- Protected Attributes inherited from gmx::ColvarProxyGromacs
t_atoms gmxAtoms_
 Atoms topology.
 
PbcType pbcType_
 Box infos.
 
t_pbc gmxPbc_
 
const MDLoggerlogger_ = nullptr
 
bool doParsing_
 Activate or not the parsing of the Colvars config file.
 
DefaultRandomEngine rng_
 
TabulatedNormalDistribution normalDistribution_
 

Constructor & Destructor Documentation

gmx::ColvarsForceProvider::ColvarsForceProvider ( const std::string &  colvarsConfigString,
t_atoms  atoms,
PbcType  pbcType,
const MDLogger logger,
const std::map< std::string, std::string > &  inputStrings,
real  ensembleTemperature,
int  seed,
LocalAtomSetManager localAtomSetManager,
const t_commrec *  cr,
double  simulationTimeStep,
const std::vector< RVec > &  colvarsCoords,
const std::string &  outputPrefix,
const ColvarsForceProviderState state 
)

Construct ColvarsForceProvider from its parameters.

Parameters
[in]colvarsConfigStringContent of the colvars input file.
[in]atomsAtoms topology
[in]pbcTypePeriodic boundary conditions
[in]loggerGROMACS logger instance
[in]inputStringsInput files stored as string in the TPR's KVT
[in]ensembleTemperaturethe constant ensemble temperature
[in]seedThe colvars seed for random number generator
[in]localAtomSetManagerAtom Manager to retrieve Colvars index atoms
[in]crCommunication Record
[in]simulationTimeStepThe simulation time step
[in]colvarsCoordsThe colvars atoms coordinates retrived from the TPR's KVT
[in]outputPrefixThe prefix for output colvars files
[in]stateThe state of colvars force provider to be written in the checkpoint

Member Function Documentation

void gmx::ColvarsForceProvider::add_energy ( cvm::real  energy)
override

From colvarproxy.

add energy to the total count of bias energy biasEnergy

Parameters
[in]energythe value of energy to add
void gmx::ColvarsForceProvider::calculateForces ( const ForceProviderInput forceProviderInput,
ForceProviderOutput forceProviderOutput 
)
overridevirtual

Calculate colvars forces.

Parameters
[in]forceProviderInputinput for force provider
[out]forceProviderOutputoutput for force provider

Implements gmx::IForceProvider.

void gmx::ColvarsForceProvider::processAtomsRedistributedSignal ( const MDModulesAtomsRedistributedSignal atomsRedistributedSignal)

Process atomsRedistributedSignal notification during mdrun.

Parameters
[in]atomsRedistributedSignalsignal recieved
void gmx::ColvarsForceProvider::writeCheckpointData ( MDModulesWriteCheckpointData  checkpointWriting,
const std::string &  moduleName 
)

Write internal colvars data to checkpoint file.

Parameters
[in]checkpointWritingenables writing to the Key-Value-Tree that is used for storing the checkpoint information
[in]moduleNamenames the module that is checkpointing this force-provider
Note
The provided state to checkpoint has to change if checkpointing is moved before the force provider call in the MD-loop.

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