Gromacs  2024.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Public Member Functions | Static Public Member Functions
gmx::Awh Class Reference

#include <gromacs/applied_forces/awh/awh.h>

Description

Coupling of the accelerated weight histogram method (AWH) with the system.

AWH calculates the free energy along order parameters of the system. Free energy barriers are overcome by adaptively tuning a bias potential along the order parameter such that the biased distribution along the parameter converges toward a chosen target distribution.

The Awh class takes care of the coupling between the system and the AWH bias(es). The Awh class contains one or more BiasCoupledToSystem objects. The BiasCoupledToSystem class takes care of the reaction coordinate input and force output for the single Bias object it containts.

Todo:
Update parameter reading and checkpointing, when general C++ framework is ready.

Public Member Functions

 Awh (FILE *fplog, const t_inputrec &inputRecord, const t_commrec *commRecord, const gmx_multisim_t *multiSimRecord, const AwhParams &awhParams, const std::string &biasInitFilename, pull_t *pull_work, int numLambdaStates, int lambdaState)
 Construct an AWH at the start of a simulation. More...
 
real applyBiasForcesAndUpdateBias (PbcType pbcType, ArrayRef< const double > neighborLambdaEnergies, ArrayRef< const double > neighborLambdaDhdl, const matrix box, double t, int64_t step, gmx_wallcycle *wallcycle, FILE *fplog)
 Peform an AWH update, to be called every MD step. More...
 
void updateHistory (AwhHistory *awhHistory) const
 Update the AWH history in preparation for writing to checkpoint file. More...
 
std::shared_ptr< AwhHistory > initHistoryFromState () const
 Allocate and initialize an AWH history with the given AWH state. More...
 
void restoreStateFromHistory (const AwhHistory *awhHistory)
 Restore the AWH state from the given history. More...
 
void writeToEnergyFrame (int64_t step, t_enxframe *fr)
 Fills the AWH data block of an energy frame with data at certain steps. More...
 
int fepLambdaState () const
 Get the current free energy lambda state. More...
 
bool needForeignEnergyDifferences (int64_t step) const
 Returns if foreign energy differences are required during this step. More...
 
bool hasFepLambdaDimension () const
 Returns true if AWH has a bias with a free energy lambda state dimension at all.
 

Static Public Member Functions

static const char * externalPotentialString ()
 Returns string "AWH" for registering AWH as an external potential provider with the pull module.
 
static void registerAwhWithPull (const AwhParams &awhParams, pull_t *pull_work)
 Register the AWH biased coordinates with pull. More...
 

Constructor & Destructor Documentation

gmx::Awh::Awh ( FILE *  fplog,
const t_inputrec &  inputRecord,
const t_commrec *  commRecord,
const gmx_multisim_t multiSimRecord,
const AwhParams awhParams,
const std::string &  biasInitFilename,
pull_t *  pull_work,
int  numLambdaStates,
int  lambdaState 
)

Construct an AWH at the start of a simulation.

AWH will here also register itself with the pull struct as the potential provider for the pull coordinates given as AWH coordinates in the user input. This allows AWH to later apply the bias force to these coordinate in Awh::applyBiasForcesAndUpdateBias.

Parameters
[in,out]fplogGeneral output file, normally md.log, can be nullptr.
[in]inputRecordGeneral input parameters (as set up by grompp).
[in]commRecordStruct for communication, can be nullptr.
[in]multiSimRecordMulti-sim handler
[in]awhParamsAWH input parameters, consistent with the relevant parts of inputRecord (as set up by grompp).
[in]biasInitFilenameName of file to read PMF and target from.
[in,out]pull_workPointer to a pull struct which AWH will couple to, has to be initialized, is assumed not to change during the lifetime of the Awh object.
[in]numLambdaStatesThe number of free energy lambda states.
[in]lambdaStateThe initial free energy lambda state of the system.
Exceptions
InvalidInputErrorIf there is user input (or combinations thereof) that is not allowed.

Member Function Documentation

real gmx::Awh::applyBiasForcesAndUpdateBias ( PbcType  pbcType,
ArrayRef< const double >  neighborLambdaEnergies,
ArrayRef< const double >  neighborLambdaDhdl,
const matrix  box,
double  t,
int64_t  step,
gmx_wallcycle *  wallcycle,
FILE *  fplog 
)

Peform an AWH update, to be called every MD step.

An update has two tasks: apply the bias force and improve the bias and the free energy estimate that AWH keeps internally.

For the first task, AWH retrieves the pull coordinate values from the pull struct. With these, the bias potential and forces are calculated. The bias force together with the atom forces and virial are passed on to pull which applies the bias force to the atoms. This is done at every step.

Secondly, coordinate values are regularly sampled and kept by AWH. Convergence of the bias and free energy estimate is achieved by updating the AWH bias state after a certain number of samples has been collected.

Note
Requires that pull_potential() from pull.h has been called first since AWH needs the current coordinate values (the pull code checks for this). Requires that pull_apply_forces() is called after this to take into account the updated pull forces AWH might be applied to.
Parameters
[in]pbcTypeType of periodic boundary conditions.
[in]neighborLambdaEnergiesAn array containing the energy of the system in neighboring lambdas. The array is of length numLambdas+1, where numLambdas is the number of free energy lambda states. Element 0 in the array is the energy of the current state and elements 1..numLambdas contain the energy of the system in the neighboring lambda states (also including the current state). When there are no free energy lambda state dimensions this can be empty.
[in]neighborLambdaDhdlAn array containing the dHdL at the neighboring lambda points. The array is of length numLambdas+1, where numLambdas is the number of free energy lambda states. Element 0 in the array is the dHdL of the current state and elements 1..numLambdas contain the dHdL of the system in the neighboring lambda states (also including the current state). When there are no free energy lambda state dimensions this can be empty.
[in]boxBox vectors.
[in]tTime.
[in]stepThe current MD step.
[in,out]wallcycleWallcycle counter, can be nullptr.
[in,out]fplogGeneral output file, normally md.log, can be nullptr.
Returns
the potential energy for the bias.
int gmx::Awh::fepLambdaState ( ) const
inline

Get the current free energy lambda state.

Returns
The value of lambdaState_.
std::shared_ptr< AwhHistory > gmx::Awh::initHistoryFromState ( ) const

Allocate and initialize an AWH history with the given AWH state.

This function should be called at the start of a new simulation at least on the main rank. Note that only constant data will be initialized here. History data is set by Awh::updateHistory.

Returns
a shared pointer to the AWH history on the rank that does I/O, nullptr otherwise.
bool gmx::Awh::needForeignEnergyDifferences ( int64_t  step) const

Returns if foreign energy differences are required during this step.

Parameters
[in]stepThe current MD step.
void gmx::Awh::registerAwhWithPull ( const AwhParams awhParams,
pull_t *  pull_work 
)
static

Register the AWH biased coordinates with pull.

This function is public because it needs to be called by grompp (and is otherwise only called by Awh()). Pull requires all external potentials to register themselves before the end of pre-processing and before the first MD step. If this has not happened, pull with throw an error.

Parameters
[in]awhParamsThe AWH parameters.
[in,out]pull_workPull struct which AWH will register the bias into.
void gmx::Awh::restoreStateFromHistory ( const AwhHistory *  awhHistory)

Restore the AWH state from the given history.

Should be called on all ranks (for internal MPI broadcast). Should pass a point to an AwhHistory on the main rank that is compatible with the AWH setup in this simulation. Will throw an exception if it is not compatible.

Parameters
[in]awhHistoryAWH history to restore from.
void gmx::Awh::updateHistory ( AwhHistory *  awhHistory) const

Update the AWH history in preparation for writing to checkpoint file.

Should be called at least on the main rank at checkpoint steps.

Should be called with a valid awhHistory (is checked).

Parameters
[in,out]awhHistoryAWH history to set.
void gmx::Awh::writeToEnergyFrame ( int64_t  step,
t_enxframe *  fr 
)

Fills the AWH data block of an energy frame with data at certain steps.

Parameters
[in]stepThe current MD step.
[in,out]frEnergy data frame.

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