SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
BasisFunc.h
1 
21 #pragma once
22 
23 #include "Method.h"
24 #include "Grids/Histogram.h"
25 #include <fstream>
26 #include <vector>
27 
28 namespace SSAGES
29 {
30 
32 
37  struct Map
38  {
40  double value;
41 
43  std::vector<int> map;
44 
46 
50  Map(const std::vector<int>& map,
51  double value) :
52  value(value), map(map)
53  {}
54  };
55 
57 
63  struct BasisLUT
64  {
66  std::vector<double> values;
67 
69  std::vector<double> derivs;
70 
72 
76  BasisLUT(const std::vector<double>& values,
77  const std::vector<double>& derivs) :
78  values(values), derivs(derivs)
79  {}
80  };
81 
83 
89  class Basis : public Method
90  {
91  private:
92 
94 
98 
100 
104  std::vector<Map> coeff_;
105 
107 
112  std::vector<double> unbias_;
113 
115  std::vector<double> coeff_arr_;
116 
118  std::vector<BasisLUT> LUT_;
119 
121 
125  std::vector<double> derivatives_;
126 
128  std::vector<unsigned int> polyords_;
129 
131 
136  std::vector<double> restraint_;
137 
139  std::vector<double> boundUp_;
140 
142  std::vector<double> boundLow_;
143 
145  unsigned int cyclefreq_;
146 
148  unsigned int mpiid_;
149 
151 
156  double weight_;
157 
159 
164  double temperature_;
165 
167  double tol_;
168 
170  bool bounds_;
171 
174 
176 
180  void UpdateBias(const CVList& cvs, const double);
181 
183 
186  void CalcBiasForce(const CVList& cvs);
187 
189 
193  void PrintBias(const CVList& cvs, const double beta);
194 
196 
199  void BasisInit(const CVList& cvs);
200 
202  std::ofstream basisout_;
203 
205  std::ofstream coeffout_;
206 
209  std::string bnme_;
210 
212  std::string cnme_;
213 
216 
217  public:
219 
240  Basis(const MPI_Comm& world,
241  const MPI_Comm& comm,
242  Histogram<int> *hist,
243  const std::vector<unsigned int>& polyord,
244  const std::vector<double>& restraint,
245  const std::vector<double>& boundUp,
246  const std::vector<double>& boundLow,
247  unsigned int cyclefreq,
248  unsigned int frequency,
249  const std::string bnme,
250  const std::string cnme,
251  const double temperature,
252  const double tol,
253  const double weight,
254  bool converge) :
255  Method(frequency, world, comm), hist_(hist),
256  coeff_(), unbias_(), coeff_arr_(), LUT_(), derivatives_(), polyords_(polyord),
257  restraint_(restraint), boundUp_(boundUp), boundLow_(boundLow),
258  cyclefreq_(cyclefreq), mpiid_(0), weight_(weight),
259  temperature_(temperature), tol_(tol),
260  converge_exit_(converge), bnme_(bnme), cnme_(cnme), iteration_(0)
261  {
262  }
263 
265 
269  void PreSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
270 
272 
276  void PostIntegration(Snapshot* snapshot, const class CVManager& cvmanager) override;
277 
279 
283  void PostSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
284 
286 
292  void SetIteration(const int iter)
293  {
294  iteration_ = iter;
295  }
296 
298 
305  void SetBasis(const std::vector<double>&coeff, std::vector<double>&unbias)
306  {
307  coeff_arr_ = coeff;
308  unbias_ = unbias;
309  }
310 
312  static Basis* Build(const Json::Value& json,
313  const MPI_Comm& world,
314  const MPI_Comm& comm,
315  const std::string& path);
316 
319  {
320  delete hist_;
321  }
322  };
323 }
324 
double value
The coefficient value.
Definition: BasisFunc.h:40
Histogram< int > * hist_
Histogram of visited states.
Definition: BasisFunc.h:97
std::vector< double > derivatives_
Derivatives of the bias potential.
Definition: BasisFunc.h:125
void PostSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Post-simulation hook.
Definition: BasisFunc.cpp:190
Collective variable manager.
Definition: CVManager.h:40
BasisLUT(const std::vector< double > &values, const std::vector< double > &derivs)
Constructor.
Definition: BasisFunc.h:76
unsigned int cyclefreq_
Frequency of coefficient updates.
Definition: BasisFunc.h:145
std::vector< CollectiveVariable * > CVList
List of Collective Variables.
Definition: types.h:51
uint iteration_
Iteration counter.
Definition: BasisFunc.h:215
void PrintBias(const CVList &cvs, const double beta)
Prints the current bias to a defined file from the JSON.
Definition: BasisFunc.cpp:365
std::vector< double > boundUp_
Upper position of the spring restraint.
Definition: BasisFunc.h:139
void SetIteration(const int iter)
Set the current iteration.
Definition: BasisFunc.h:292
Class containing a snapshot of the current simulation in time.
Definition: Snapshot.h:43
Basis Function Sampling Algorithm.
Definition: BasisFunc.h:89
Map(const std::vector< int > &map, double value)
Constructor.
Definition: BasisFunc.h:50
Interface for Method implementations.
Definition: Method.h:43
void UpdateBias(const CVList &cvs, const double)
Updates the bias projection of the PMF.
Definition: BasisFunc.cpp:244
std::vector< BasisLUT > LUT_
The Basis set lookup table, also defined globally.
Definition: BasisFunc.h:118
void PreSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Pre-simulation hook.
Definition: BasisFunc.cpp:36
std::vector< double > coeff_arr_
The coefficient array for restart runs.
Definition: BasisFunc.h:115
Map for histogram and coefficients.
Definition: BasisFunc.h:37
void SetBasis(const std::vector< double > &coeff, std::vector< double > &unbias)
Set the values for the basis.
Definition: BasisFunc.h:305
std::string bnme_
Definition: BasisFunc.h:209
void BasisInit(const CVList &cvs)
Initializes the look up tables for polynomials.
Definition: BasisFunc.cpp:198
unsigned int mpiid_
The node that the current system belongs to, primarily for printing and debugging.
Definition: BasisFunc.h:148
double tol_
The tolerance criteria for the system .
Definition: BasisFunc.h:167
Look-up table for basis functions.
Definition: BasisFunc.h:63
std::ofstream basisout_
Output stream for basis projection data.
Definition: BasisFunc.h:202
void PostIntegration(Snapshot *snapshot, const class CVManager &cvmanager) override
Post-integration hook.
Definition: BasisFunc.cpp:124
std::vector< double > restraint_
Spring constants for restrained system.
Definition: BasisFunc.h:136
static Basis * Build(const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path)
Definition: BasisFunc.cpp:521
std::vector< double > values
The values of the basis sets.
Definition: BasisFunc.h:66
std::ofstream coeffout_
Output stream for coefficients (for reading purposes)
Definition: BasisFunc.h:205
~Basis()
Destructor.
Definition: BasisFunc.h:318
std::vector< double > derivs
The values of the derivatives of the basis sets.
Definition: BasisFunc.h:69
void CalcBiasForce(const CVList &cvs)
Computes the bias force.
Definition: BasisFunc.cpp:440
std::vector< unsigned int > polyords_
The order of the basis polynomials.
Definition: BasisFunc.h:128
std::vector< Map > coeff_
Globally located coefficient values.
Definition: BasisFunc.h:104
double temperature_
Self-defined temperature.
Definition: BasisFunc.h:164
Basis(const MPI_Comm &world, const MPI_Comm &comm, Histogram< int > *hist, const std::vector< unsigned int > &polyord, 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 std::string cnme, const double temperature, const double tol, const double weight, bool converge)
Constructor.
Definition: BasisFunc.h:240
double weight_
Weighting for potentially faster sampling.
Definition: BasisFunc.h:156
bool converge_exit_
A check to see if you want the system to end when it reaches the convergence criteria.
Definition: BasisFunc.h:173
std::vector< int > map
The mapping in an array of integers.
Definition: BasisFunc.h:43
bool bounds_
A variable to check to see if the simulation is in bounds or not.
Definition: BasisFunc.h:170
std::vector< double > boundLow_
Lower position of the spring restraint.
Definition: BasisFunc.h:142
std::vector< double > unbias_
The biased histogram of states.
Definition: BasisFunc.h:112
std::string cnme_
Coefficient filename.
Definition: BasisFunc.h:212