SSAGES
0.9.3
Software Suite for Advanced General Ensemble Simulations
|
Basis Function Sampling Algorithm. More...
#include <BasisFunc.h>
Public Member Functions | |
BFS (const MPI_Comm &world, const MPI_Comm &comm, Grid< unsigned int > *h, Grid< std::vector< double >> *f, Grid< double > *b, const std::vector< BasisFunction * > &functions, const std::vector< double > &restraint, const std::vector< double > &boundUp, const std::vector< double > &boundLow, unsigned int cyclefreq, unsigned int frequency, const std::string bnme, const double temperature, const double tol, const double weight, bool converge) | |
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 | SetIteration (const int iter) |
Set the current iteration. More... | |
void | SetBasis (const std::vector< double > &coeff, std::vector< double > &unbias) |
Set the values for the basis. More... | |
~BFS () | |
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 BFS * | 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 | ProjectBias (const CVList &cvs, const double beta) |
The functions which calculates the updated bias and coefficients and then stores them. More... | |
void | InBounds (const CVList &cvs) |
A function which checks to see if the CVs are still in bounds. More... | |
void | PrintBias (const CVList &cvs, const double beta) |
Prints the current bias to a defined file from the JSON. More... | |
Private Attributes | |
Grid< unsigned int > * | h_ |
Grid of visited states. More... | |
Grid< double > * | b_ |
Stored bias potential. More... | |
Grid< std::vector< double > > * | f_ |
Stored gradients. More... | |
std::vector< double > | unbias_ |
The biased histogram of states. More... | |
std::vector< double > | coeffArr_ |
The coefficient array for restart runs. | |
BasisEvaluator | evaluator_ |
The basis evaluator class. More... | |
std::vector< double > | restraint_ |
Spring constants for restrained system. More... | |
std::vector< double > | boundUp_ |
Upper position of the spring restraint. | |
std::vector< double > | boundLow_ |
Lower position of the spring restraint. | |
unsigned int | cyclefreq_ |
Frequency of coefficient updates. | |
unsigned int | cyclefreqTmp_ |
Temporary cyclefrequency that ensures sweeps are of proper length for restrained simulations. | |
unsigned int | step_ |
Step counter for the cyclefrequency. | |
unsigned int | mpiid_ |
The node that the current system belongs to, primarily for printing and debugging. | |
double | weight_ |
Weighting for potentially faster sampling. More... | |
double | temperature_ |
Self-defined temperature. More... | |
double | tol_ |
The tolerance criteria for the system . | |
bool | bounds_ |
A variable to check to see if the simulation is in bounds or not. | |
bool | convergeExit_ |
A check to see if you want the system to end when it reaches the convergence criteria. | |
std::string | bnme_ |
unsigned int | iteration_ |
Iteration counter. | |
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. | |
Basis Function Sampling Algorithm.
Implementation of the Basis Function Sampling Method based on [5].
Definition at line 38 of file BasisFunc.h.
|
inline |
Constructor.
world | MPI global communicator. |
comm | MPI local communicator. |
h | Grid containing histogram. |
f | Grid containing histogram. |
b | Grid containing histogram. |
functions | Vector of basis functions. |
restraint | Restraint spring constants. |
boundUp | Upper bounds of restraint springs. |
boundLow | Lower bounds of restraint springs. |
cyclefreq | Cycle frequency. |
frequency | Frequency with which this Method is applied. |
bnme | Basis file name. |
temperature | Automatically set temperature. |
tol | Threshold for tolerance criterion. |
weight | Weight for improved sampling. |
converge | If True quit on convergence. |
Constructs an instance of the Basis function sampling method. The coefficients describe the basis projection of the system. This is updated once every cyclefreq_. For now, only the Legendre polynomial is implemented. Others will be added later.
Definition at line 181 of file BasisFunc.h.
|
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 308 of file BasisFunc.cpp.
References Json::Requirement::GetErrors(), SSAGES::GridBase< T >::GetNumPoints(), Json::Requirement::HasErrors(), Json::ObjectRequirement::Parse(), and Json::ObjectRequirement::Validate().
|
private |
A function which checks to see if the CVs are still in bounds.
cvs | List of CVs to check |
Definition at line 272 of file BasisFunc.cpp.
|
overridevirtual |
Method call post integration.
snapshot | Pointer to the simulation snapshot. |
cvmanager | Collective variable manager. |
This function will be called after each integration step.
Implements SSAGES::Method.
Definition at line 81 of file BasisFunc.cpp.
References SSAGES::CVManager::GetCVs(), SSAGES::Snapshot::GetForces(), SSAGES::Snapshot::GetIteration(), SSAGES::Snapshot::GetKb(), and SSAGES::Snapshot::GetVirial().
|
overridevirtual |
Method call post simulation.
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.
Definition at line 160 of file BasisFunc.cpp.
|
overridevirtual |
Method call prior to simulation initiation.
snapshot | Pointer to the simulation snapshot. |
cvmanager | Collective variable manager. |
This function will be called before the simulation is started.
Implements SSAGES::Method.
Definition at line 36 of file BasisFunc.cpp.
References SSAGES::CVManager::GetCVs(), and SSAGES::Snapshot::GetWalkerID().
|
private |
Prints the current bias to a defined file from the JSON.
cvs | List of collective variables. |
beta | Scale parameter. |
Definition at line 227 of file BasisFunc.cpp.
References SSAGES::Grid< T >::begin().
|
private |
The functions which calculates the updated bias and coefficients and then stores them.
cvs | List of CVs onto which the bias is projected |
beta | Thermodynamic beta (1/kT) |
Definition at line 166 of file BasisFunc.cpp.
References SSAGES::Grid< T >::begin(), and SSAGES::GridBase< T >::data().
|
inline |
Set the values for the basis.
coeff | List of coefficients. |
unbias | List of unbiased values. |
This function is used to set starting values at the beginning of a run. For example when continuing from a restart value.
Definition at line 235 of file BasisFunc.h.
References evaluator_, SSAGES::BasisEvaluator::SetCoeff(), and unbias_.
|
inline |
Set the current iteration.
iter | New value for the current iteration. |
This function is used to set the current iteration, for example when continuing from a restart.
Definition at line 222 of file BasisFunc.h.
References iteration_.
|
private |
Stored bias potential.
The sum of basis functions that adds up to the bias potential of the surface
Definition at line 52 of file BasisFunc.h.
|
private |
The option to name both the basis and coefficient files will be given Basis filename
Definition at line 151 of file BasisFunc.h.
|
private |
The basis evaluator class.
A class which holds the basis functions and updates the coefficients at each cycle
Definition at line 76 of file BasisFunc.h.
Referenced by SetBasis().
|
private |
Stored gradients.
The gradients corresponding to each point on the grid
Definition at line 58 of file BasisFunc.h.
|
private |
|
private |
Spring constants for restrained system.
The system uses this to determine if the system is to be restrained on the defined interval. The user inputs the spring constants if the system is not periodic.
Definition at line 84 of file BasisFunc.h.
|
private |
Self-defined temperature.
In the case of the MD engine using a poorly defined temperature, the option to throw it into the method is available. If not provided it takes the value from the engine.
Definition at line 118 of file BasisFunc.h.
|
private |
The biased histogram of states.
The biased histogram of states has the form hist_*exp(phi*beta), where phi is the bias potential and beta is the inverse of the temperature. It is defined globally.
Definition at line 66 of file BasisFunc.h.
Referenced by SetBasis().
|
private |
Weighting for potentially faster sampling.
Weighting can be used to potentially sample faster, however it can cause the simulation to explode. By default this value will be set to 1.0
Definition at line 110 of file BasisFunc.h.