SSAGES  0.9.3
Software Suite for Advanced General Ensemble Simulations
ABF.h
1 
22 #pragma once
23 
24 #include "Method.h"
25 #include <fstream>
26 #include <vector>
27 #include <string>
28 #include "Grids/Grid.h"
29 
30 
31 namespace SSAGES
32 {
34 
40  class ABF : public Method
41  {
42  private:
44 
49 
51 
57 
62  std::vector<Grid<double>*> F_;
63 
65 
70  std::vector<Grid<double>*> Fworld_;
71 
73 
81  std::vector<std::vector<double>> restraint_;
82 
84  std::vector<bool> isperiodic_;
85 
87 
94  std::vector<std::vector<double>> periodicboundaries_;
95 
98  int min_;
99 
101  Eigen::VectorXd wdotp1_;
102 
104  Eigen::VectorXd wdotp2_;
105 
107  Eigen::VectorXd Fold_;
108 
111 
113  std::vector<Vector3> biases_;
114 
116  unsigned int dim_;
117 
119  std::ofstream worldout_;
120 
123  std::string filename_;
124 
126  std::string Nworld_filename_;
127 
129  std::string Fworld_filename_;
130 
132 
140  std::vector<std::vector<double>> histdetails_;
141 
143 
147 
149 
155  double unitconv_;
156 
158 
163  void CalcBiasForce(const Snapshot* snapshot, const CVList& cvs, const std::vector<double> &cvVals);
164 
166  void WriteData();
167 
169  double timestep_;
170 
172  unsigned int iteration_;
173 
175  Eigen::VectorXd mass_;
176 
178 
183  bool boundsCheck(const std::vector<double> &CVs);
184 
185 
186  public:
188 
213  ABF(const MPI_Comm& world,
214  const MPI_Comm& comm,
215  Grid<int> *N,
216  Grid<int> *Nworld,
217  std::vector<Grid<double>*> F,
218  std::vector<Grid<double>*> Fworld,
219  std::vector<std::vector<double>> restraint,
220  std::vector<bool> isperiodic,
221  std::vector<std::vector<double>> periodicboundaries,
222  double min,
223  bool massweigh,
224  std::string filename,
225  std::string Nworld_filename,
226  std::string Fworld_filename,
227  const std::vector<std::vector<double>>& histdetails,
228  int FBackupInterv,
229  double unitconv,
230  double timestep,
231  unsigned int frequency) :
232  Method(frequency, world, comm), N_(N), Nworld_(Nworld), F_(F), Fworld_(Fworld), restraint_(restraint), isperiodic_(isperiodic), periodicboundaries_(periodicboundaries),
233  min_(min), wdotp1_(), wdotp2_(), Fold_(), massweigh_(massweigh),
234  biases_(), dim_(0), filename_(filename), Nworld_filename_(Nworld_filename), Fworld_filename_(Fworld_filename),
235  histdetails_(histdetails),
236  FBackupInterv_(FBackupInterv), unitconv_(unitconv), timestep_(timestep),
237  iteration_(0)
238  {
239  }
240 
242  void PreSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
243 
245  void PostIntegration(Snapshot* snapshot, const class CVManager& cvmanager) override;
246 
248  void PostSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
249 
251 
254  void SetIteration(const int iter)
255  {
256  iteration_ = iter;
257  }
258 
260  static ABF* Build(const Json::Value& json,
261  const MPI_Comm& world,
262  const MPI_Comm& comm,
263  const std::string& path);
264 
266  ~ABF() {}
267  };
268 }
269 
Adaptive Biasing Force Algorithm.
Definition: ABF.h:41
double unitconv_
Unit Conversion Constant from W dot P to Force.
Definition: ABF.h:155
Eigen::VectorXd Fold_
To hold last iterations F_ value for removing bias.
Definition: ABF.h:107
void PostSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call prior to simulation initiation.
Definition: ABF.cpp:171
std::vector< std::vector< double > > histdetails_
Histogram details.
Definition: ABF.h:140
Eigen::VectorXd wdotp2_
To hold second to last iteration wdotp value for numerical derivative.
Definition: ABF.h:104
void PreSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call prior to simulation initiation.
Definition: ABF.cpp:39
double timestep_
Timestep of integration.
Definition: ABF.h:169
Grid< int > * N_
To store number of local hits at a given CV bin.
Definition: ABF.h:48
std::vector< Grid< double > * > F_
To store running total of the local walker.
Definition: ABF.h:62
std::string Nworld_filename_
Nworld print out filename, for restarts.
Definition: ABF.h:126
std::ofstream worldout_
Output stream for F/N world data.
Definition: ABF.h:119
ABF(const MPI_Comm &world, const MPI_Comm &comm, Grid< int > *N, Grid< int > *Nworld, std::vector< Grid< double > * > F, std::vector< Grid< double > * > Fworld, std::vector< std::vector< double >> restraint, std::vector< bool > isperiodic, std::vector< std::vector< double >> periodicboundaries, double min, bool massweigh, std::string filename, std::string Nworld_filename, std::string Fworld_filename, const std::vector< std::vector< double >> &histdetails, int FBackupInterv, double unitconv, double timestep, unsigned int frequency)
Constructor.
Definition: ABF.h:213
Eigen::VectorXd wdotp1_
To hold last iteration wdotp value for numerical derivative.
Definition: ABF.h:101
std::vector< Grid< double > * > Fworld_
Will hold the global total, synced across walkers at every time step.
Definition: ABF.h:70
std::vector< std::vector< double > > restraint_
Information for a harmonic restraint to keep CV in the region of interest.
Definition: ABF.h:81
std::vector< std::vector< double > > periodicboundaries_
Holds periodic boundaries of CVs.
Definition: ABF.h:94
std::string filename_
Definition: ABF.h:123
int min_
Definition: ABF.h:98
unsigned int dim_
Number of CVs in system.
Definition: ABF.h:116
static ABF * Build(const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path)
Build a derived method from JSON node.
Definition: ABF.cpp:297
void CalcBiasForce(const Snapshot *snapshot, const CVList &cvs, const std::vector< double > &cvVals)
Computes the bias force.
Definition: ABF.cpp:176
unsigned int iteration_
Iteration counter.
Definition: ABF.h:172
std::vector< bool > isperiodic_
For each CV, holds whether that CV has periodic boundaries or not.
Definition: ABF.h:84
~ABF()
Destructor.
Definition: ABF.h:266
bool massweigh_
Mass weighing of bias enabled/disabled.
Definition: ABF.h:110
void WriteData()
Writes out data to file.
Definition: ABF.cpp:232
void SetIteration(const int iter)
Set iteration of the method.
Definition: ABF.h:254
Eigen::VectorXd mass_
Mass vector. Empty unless required.
Definition: ABF.h:175
int FBackupInterv_
Integer to hold F estimate backup information.
Definition: ABF.h:146
void PostIntegration(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call post integration.
Definition: ABF.cpp:72
bool boundsCheck(const std::vector< double > &CVs)
Checks whether the local walker is within CV bounds.
Definition: ABF.cpp:284
Grid< int > * Nworld_
To store number of global hits at a given CV bin.
Definition: ABF.h:55
std::vector< Vector3 > biases_
Biases applied to atoms each timestep.
Definition: ABF.h:113
std::string Fworld_filename_
Fworld print out filename, for restarts.
Definition: ABF.h:129
Collective variable manager.
Definition: CVManager.h:43
Interface for Method implementations.
Definition: Method.h:44
Class containing a snapshot of the current simulation in time.
Definition: Snapshot.h:48
std::vector< CollectiveVariable * > CVList
List of Collective Variables.
Definition: types.h:51