SSAGES
0.9.3
Software Suite for Advanced General Ensemble Simulations
|
Combined Force Frequency (CFF) Algorithm. More...
#include <CFF.h>
Public Member Functions | |
CFF (const MPI_Comm &world, const MPI_Comm &comm, const Eigen::VectorXi &topol, Grid< Eigen::VectorXd > *fgrid, Grid< unsigned int > *hgrid, Grid< double > *ugrid, std::vector< Grid< double > * > F, std::vector< Grid< double > * > Fworld, const std::vector< double > &lowerb, const std::vector< double > &upperb, const std::vector< double > &lowerk, const std::vector< double > &upperk, double temperature, double unitconv, double timestep, double weight, unsigned int nsweep, int min) | |
Constructor. More... | |
void | PreSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override |
Method call prior to simulation initiation. More... | |
void | PostIntegration (Snapshot *snapshot, const class CVManager &cvmanager) override |
Method call post integration. More... | |
void | PostSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override |
Method call post simulation. More... | |
void | SetPrevWeight (double h) |
Set previous history weight. | |
void | SetOutput (const std::string &outfile) |
Set name of output file. | |
void | SetOutputOverwrite (bool overwrite) |
Set overwrite flag on output file. | |
void | SetConvergeIters (unsigned int citers) |
Set number of iterations after which we turn on full weight. | |
void | SetMaxIters (unsigned int iters) |
Set maximum number of training iterations per sweep. | |
void | SetMinLoss (double loss) |
Set minimum loss function value (should be zero for production). | |
void | SetRatio (double trainratio) |
Set training ratio for gradient vs value. | |
~CFF () | |
Destructor. | |
Public Member Functions inherited from SSAGES::Method | |
Method (unsigned int frequency, const MPI_Comm &world, const MPI_Comm &comm) | |
Constructor. More... | |
void | SetCVMask (const std::vector< unsigned int > &mask) |
Sets the collective variable mask. More... | |
virtual | ~Method () |
Destructor. | |
Public Member Functions inherited from SSAGES::EventListener | |
EventListener (unsigned int frequency) | |
Constructor. More... | |
unsigned int | GetFrequency () const |
Get frequency of event listener. More... | |
virtual | ~EventListener () |
Destructor. | |
Static Public Member Functions | |
static CFF * | Build (const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path) |
Build a derived method from JSON node. More... | |
Static Public Member Functions inherited from SSAGES::Method | |
static Method * | BuildMethod (const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path) |
Build a derived method from JSON node. More... | |
Static Public Member Functions inherited from SSAGES::EventListener | |
static unsigned int | GetWalkerID (const MPI_Comm &world, const MPI_Comm &comm) |
Get walker ID number of specified communicator. More... | |
static unsigned int | GetNumWalkers (const MPI_Comm &world, const MPI_Comm &comm) |
Get total number of walkers in the simulation. More... | |
static bool | IsMasterRank (const MPI_Comm &comm) |
Check if current processor is master. More... | |
Private Member Functions | |
void | TrainNetwork () |
Trains the neural network. More... | |
void | WriteBias () |
Writes out the bias and CFF output files. | |
Private Attributes | |
Eigen::VectorXi | topol_ |
Neural network topology. | |
unsigned int | sweep_ |
Number of iterations per sweep. | |
unsigned int | nsweep_ |
unsigned int | citers_ |
Number of iterations after which we turn on full weight. | |
nnet::neural_net | net_ |
Neural network trained using visit frequency. | |
nnet::neural_net | net2_ |
Neural network trained on both visit frequency and force. | |
double | timestep_ |
Timestep. | |
double | pweight_ |
Previous and current histogram weight. | |
double | weight_ |
double | temp_ |
System temperature and energy units. | |
double | kbt_ |
std::vector< Grid< double > * > | F_ |
Generalized force grid that stores total of the local walker. | |
std::vector< Grid< double > * > | Fworld_ |
Generalized force grid that stors the global total. | |
Eigen::VectorXd | Fold_ |
To hold the last iterations F_ value for removing bias. | |
Eigen::VectorXd | wdotp1_ |
To hold last iteration wdotp value for numerical derivative. | |
Eigen::VectorXd | wdotp2_ |
To hold second to last iteration wdotp value for numerical derivative. | |
Grid< Eigen::VectorXd > * | fgrid_ |
Force grid. | |
Grid< unsigned int > * | hgrid_ |
Histogram grid. | |
Eigen::ArrayXi | Nworld_ |
Histogram grid for that sotres the number of global hits. | |
Grid< double > * | ugrid_ |
Unbiased histogram grid. | |
double | ratio_ |
Hold parameters to adjust ratio of 1st and 2nd neural networks (freq vs force-based). | |
Eigen::MatrixXd | hist_ |
Eigen matrices of grids. | |
Eigen::MatrixXd | bias_ |
std::vector< double > | lowerb_ |
Bounds. | |
std::vector< double > | upperb_ |
std::vector< double > | lowerk_ |
Bound restraints. | |
std::vector< double > | upperk_ |
std::string | outfile_ |
Output filename. | |
bool | overwrite_ |
Overwrite outputs. | |
int | dim_ |
Number of CVs in system. | |
double | unitconv_ |
Unit conversion from mass*velocity/time to force. | |
int | min_ |
Eigen::MatrixXd | force_to_val_ratio_ |
To hold booleans for training neural network only in specific region for net2_. | |
Additional Inherited Members | |
Protected Attributes inherited from SSAGES::Method | |
mxx::comm | world_ |
Global MPI communicator. | |
mxx::comm | comm_ |
Local MPI communicator. | |
std::vector< unsigned int > | cvmask_ |
Mask which identifies which CVs to act on. | |
Combined Force Frequency (CFF) Algorithm.
Implementation of the CFF algorithm based on Sevgen et al. J. Chem. Theory Comput. 2020, 16, 3, 1448-1455.
|
inline |
Constructor.
world | MPI global communicator. |
comm | MPI local communicator. |
topol | Topology of network. |
fgrid | Grid containing biasing forces. |
hgrid | Grid containing histogram. |
ugrid | Grid containing unbiased histogram. |
F | Vector of grids holding local raw generalized force totals per bin per CV. |
Fworld | Vector of grids holding global raw generalized force totals per bin per CV. |
lowerb | Lower bounds for CVs. |
upperb | Upper bounds for CVs. |
lowerk | Lower bound restraints for CVs. |
upperk | Upper bound restraints for CVs. |
temperature | Temperature of the simulation. |
unitconv | Unit conversion from d(momentum)/d(time) to force. |
timestep | Simulation time step. |
weight | Relative weight of the statistics in sweep. |
nsweep | Number of iterations in the sweep. |
min | Number of counts for scaling back force biasing |
Constructs an instance of Combined Force Frequency method.
Definition at line 151 of file CFF.h.
References SSAGES::Grid< T >::begin(), SSAGES::Grid< T >::end(), SSAGES::GridBase< T >::GetDimension(), hgrid_, hist_, net2_, net_, and SSAGES::GridBase< T >::size().
|
static |
Build a derived method from JSON node.
json | JSON Value containing all input information. |
world | MPI global communicator. |
comm | MPI local communicator. |
path | Path for JSON path specification. |
This function builds a registered method from a JSON node. The difference between this function and "Build" is that this automatically determines the appropriate derived type based on the JSON node information.
Definition at line 407 of file CFF.cpp.
References Json::Requirement::GetErrors(), Json::Requirement::HasErrors(), Json::ObjectRequirement::Parse(), and Json::ObjectRequirement::Validate().
|
overridevirtual |
Method call post integration.
Post-integration hook.
snapshot | Pointer to the simulation snapshot. |
cvmanager | Collective variable manager. |
This function will be called after each integration step.
First, information from current snapshot is retrieved and stored in v-integration hook. Then, coordinates in CV space are determined. Then, for each CV, biasing force is calculated. Then, neural networks are trained if called. Then, information is printed out if called. Finally, bias is applied using either genralized force (for the first sweep) or neural networks.
Implements SSAGES::Method.
Definition at line 79 of file CFF.cpp.
References SSAGES::CVManager::GetCVs(), SSAGES::Snapshot::GetForces(), SSAGES::Snapshot::GetIteration(), SSAGES::Snapshot::GetMasses(), SSAGES::Snapshot::GetNumAtoms(), SSAGES::Snapshot::GetVelocities(), and SSAGES::Snapshot::GetVirial().
|
overridevirtual |
Method call post simulation.
Post-simulation hook.
snapshot | Pointer to the simulation snapshot. |
cvmanager | Collective variable manager. |
This function will be called after the end of the simulation run.
Implements SSAGES::Method.
|
overridevirtual |
Method call prior to simulation initiation.
Pre-simulation hook.
snapshot | Pointer to the simulation snapshot. |
cvmanager | Collective variable manager. |
This function will be called before the simulation is started.
Initialize biasing forces and histogram.
Implements SSAGES::Method.
Definition at line 41 of file CFF.cpp.
References SSAGES::CVManager::GetCVs(), and SSAGES::Snapshot::GetKb().
|
private |
Trains the neural network.
Train neural networks.
Definition at line 232 of file CFF.cpp.
|
private |