SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
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 {
33 
35 
42  class ABF : public Method
43  {
44  private:
46 
51  std::vector<Grid<double>*> F_;
52 
54 
59  std::vector<Grid<double>*> Fworld_;
60 
62 
67 
69 
74 
76 
84  std::vector<std::vector<double>> restraint_;
85 
87  std::vector<bool> isperiodic_;
88 
90 
97  std::vector<std::vector<double>> periodicboundaries_;
98 
101  int min_;
102 
104  Eigen::VectorXd wdotp1_;
105 
107  Eigen::VectorXd wdotp2_;
108 
110  Eigen::VectorXd Fold_;
111 
114 
116  std::vector<Vector3> biases_;
117 
119  unsigned int dim_;
120 
122  std::ofstream worldout_;
123 
126  std::string filename_;
127 
129  std::string Nworld_filename_;
130 
132  std::string Fworld_filename_;
133 
135 
143  std::vector<std::vector<double>> histdetails_;
144 
146 
150 
152 
158  double unitconv_;
159 
161 
166  void CalcBiasForce(const Snapshot* snapshot, const CVList& cvs, const std::vector<double> &cvVals);
167 
169  void WriteData();
170 
172  double timestep_;
173 
176 
178  Eigen::VectorXd mass_;
179 
181 
186  bool boundsCheck(const std::vector<double> &CVs);
187 
188 
189  public:
191 
216  ABF(const MPI_Comm& world,
217  const MPI_Comm& comm,
218  Grid<int> *N,
219  Grid<int> *Nworld,
220  std::vector<Grid<double>*> F,
221  std::vector<Grid<double>*> Fworld,
222  std::vector<std::vector<double>> restraint,
223  std::vector<bool> isperiodic,
224  std::vector<std::vector<double>> periodicboundaries,
225  double min,
226  bool massweigh,
227  std::string filename,
228  std::string Nworld_filename,
229  std::string Fworld_filename,
230  const std::vector<std::vector<double>>& histdetails,
231  int FBackupInterv,
232  double unitconv,
233  double timestep,
234  unsigned int frequency) :
235  Method(frequency, world, comm), N_(N), Nworld_(Nworld), F_(F), Fworld_(Fworld), restraint_(restraint), isperiodic_(isperiodic), periodicboundaries_(periodicboundaries),
236  min_(min), wdotp1_(), wdotp2_(), Fold_(), massweigh_(massweigh),
237  biases_(), dim_(0), filename_(filename), Nworld_filename_(Nworld_filename), Fworld_filename_(Fworld_filename),
238  histdetails_(histdetails),
239  FBackupInterv_(FBackupInterv), unitconv_(unitconv), timestep_(timestep),
240  iteration_(0)
241  {
242  }
243 
245 
249  void PreSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
250 
252 
256  void PostIntegration(Snapshot* snapshot, const class CVManager& cvmanager) override;
257 
259 
263  void PostSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
264 
266 
269  void SetIteration(const int iter)
270  {
271  iteration_ = iter;
272  }
273 
275  static ABF* Build(const Json::Value& json,
276  const MPI_Comm& world,
277  const MPI_Comm& comm,
278  const std::string& path);
279 
281  ~ABF() {}
282  };
283 }
284 
static ABF * Build(const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path)
Definition: ABF.cpp:311
void PostSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Post-simulation hook.
Definition: ABF.cpp:174
std::string filename_
Definition: ABF.h:126
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:216
Grid< int > * N_
To store number of local hits at a given CV bin.
Definition: ABF.h:66
Collective variable manager.
Definition: CVManager.h:40
Eigen::VectorXd mass_
Mass vector. Empty unless required.
Definition: ABF.h:178
std::vector< CollectiveVariable * > CVList
List of Collective Variables.
Definition: types.h:51
std::string Fworld_filename_
Fworld print out filename, for restarts.
Definition: ABF.h:132
~ABF()
Destructor.
Definition: ABF.h:281
Class containing a snapshot of the current simulation in time.
Definition: Snapshot.h:43
Interface for Method implementations.
Definition: Method.h:43
std::ofstream worldout_
Output stream for F/N world data.
Definition: ABF.h:122
void PreSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Pre-simulation hook.
Definition: ABF.cpp:39
void WriteData()
Writes out data to file.
Definition: ABF.cpp:235
double unitconv_
Unit Conversion Constant from W dot P to Force.
Definition: ABF.h:158
std::vector< std::vector< double > > histdetails_
Histogram details.
Definition: ABF.h:143
std::vector< Grid< double > * > F_
To store running total of the local walker.
Definition: ABF.h:51
bool massweigh_
Mass weighing of bias enabled/disabled.
Definition: ABF.h:113
Eigen::VectorXd Fold_
To hold last iterations F_ value for removing bias.
Definition: ABF.h:110
std::vector< bool > isperiodic_
For each CV, holds whether that CV has periodic boundaries or not.
Definition: ABF.h:87
Adaptive Biasing Force Algorithm.
Definition: ABF.h:42
void SetIteration(const int iter)
Set iteration of the method.
Definition: ABF.h:269
std::vector< Vector3 > biases_
Biases applied to atoms each timestep.
Definition: ABF.h:116
std::vector< std::vector< double > > restraint_
Information for a harmonic restraint to keep CV in the region of interest.
Definition: ABF.h:84
void CalcBiasForce(const Snapshot *snapshot, const CVList &cvs, const std::vector< double > &cvVals)
Computes the bias force.
Definition: ABF.cpp:179
std::string Nworld_filename_
Nworld print out filename, for restarts.
Definition: ABF.h:129
unsigned int dim_
Number of CVs in system.
Definition: ABF.h:119
std::vector< Grid< double > * > Fworld_
Will hold the global total, synced across walkers at every time step.
Definition: ABF.h:59
bool boundsCheck(const std::vector< double > &CVs)
Checks whether the local walker is within CV bounds.
Definition: ABF.cpp:298
void PostIntegration(Snapshot *snapshot, const class CVManager &cvmanager) override
Post-integration hook.
Definition: ABF.cpp:76
int min_
Definition: ABF.h:101
int FBackupInterv_
Integer to hold F estimate backup information.
Definition: ABF.h:149
Grid< int > * Nworld_
To store number of global hits at a given CV bin.
Definition: ABF.h:73
Eigen::VectorXd wdotp2_
To hold second to last iteration wdotp value for numerical derivative.
Definition: ABF.h:107
double timestep_
Timestep of integration.
Definition: ABF.h:172
uint iteration_
Iteration counter.
Definition: ABF.h:175
std::vector< std::vector< double > > periodicboundaries_
Holds periodic boundaries of CVs.
Definition: ABF.h:97
Eigen::VectorXd wdotp1_
To hold last iteration wdotp value for numerical derivative.
Definition: ABF.h:104