SSAGES  0.9.3
Software Suite for Advanced General Ensemble Simulations
CFF.h
1 
21 #pragma once
22 
23 #include "Method.h"
24 #include "Grids/Grid.h"
25 #include "nnet/nnet.h"
26 
27 namespace SSAGES
28 {
30 
36  class CFF : public Method
37  {
38  private:
40  Eigen::VectorXi topol_;
41 
43  unsigned int sweep_, nsweep_;
44 
46  unsigned int citers_;
47 
49  nnet::neural_net net_;
50 
52  nnet::neural_net net2_;
53 
55  double timestep_;
56 
58  double pweight_, weight_;
59 
61  double temp_, kbt_;
62 
64  std::vector<Grid<double>*> F_;
65 
67  std::vector<Grid<double>*> Fworld_;
68 
70  Eigen::VectorXd Fold_;
71 
73  Eigen::VectorXd wdotp1_;
74 
76  Eigen::VectorXd wdotp2_;
77 
80 
83 
85  Eigen::ArrayXi Nworld_;
86 
89 
91  double ratio_;
92 
94  Eigen::MatrixXd hist_, bias_;
95 
97  std::vector<double> lowerb_, upperb_;
98 
100  std::vector<double> lowerk_, upperk_;
101 
103  std::string outfile_;
104 
107 
109  void TrainNetwork();
110 
112  void WriteBias();
113 
115  int dim_;
116 
118  double unitconv_;
119 
122  int min_;
123 
125  Eigen::MatrixXd force_to_val_ratio_;
126 
127  public:
129 
152  const MPI_Comm& world,
153  const MPI_Comm& comm,
154  const Eigen::VectorXi& topol,
155  Grid<Eigen::VectorXd>* fgrid,
156  Grid<unsigned int>* hgrid,
157  Grid<double>* ugrid,
158  std::vector<Grid<double>*> F,
159  std::vector<Grid<double>*> Fworld,
160  const std::vector<double>& lowerb,
161  const std::vector<double>& upperb,
162  const std::vector<double>& lowerk,
163  const std::vector<double>& upperk,
164  double temperature,
165  double unitconv,
166  double timestep,
167  double weight,
168  unsigned int nsweep,
169  int min
170  ) :
171  Method(1, world, comm), topol_(topol), sweep_(0), nsweep_(nsweep), citers_(0),
172  net_(topol), net2_(topol), timestep_(timestep), pweight_(1.), weight_(weight),
173  temp_(temperature), kbt_(), F_(F), Fworld_(Fworld), fgrid_(fgrid),
174  hgrid_(hgrid), ugrid_(ugrid), hist_(), bias_(), lowerb_(lowerb),
175  upperb_(upperb), lowerk_(lowerk), upperk_(upperk), outfile_("CFF.out"),
176  overwrite_(true), unitconv_(unitconv), min_(min)
177  {
178  // Create histogram grid matrix.
179  hist_.resize(hgrid_->size(), hgrid_->GetDimension());
180 
181  // Fill it up.
182  int i = 0;
183  for(auto it = hgrid_->begin(); it != hgrid_->end(); ++it)
184  {
185  auto coord = it.coordinates();
186  auto n = coord.size();
187  for(decltype(n) j = 0; j < n; ++j)
188  hist_(i, j) = coord[j];
189  ++i;
190  }
191 
192  // Initialize free energy surface vector.
193  bias_.resize(hgrid_->size(), 1);
194  net_.forward_pass(hist_);
195  net2_.forward_pass(hist_);
196  bias_.array() = net_.get_activation().col(0).array();
197  }
198 
200  void PreSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
201 
203  void PostIntegration(Snapshot* snapshot, const class CVManager& cvmanager) override;
204 
206  void PostSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
207 
209  void SetPrevWeight(double h)
210  {
211  pweight_ = h;
212  }
213 
215  void SetOutput(const std::string& outfile)
216  {
217  outfile_ = outfile;
218  }
219 
221  void SetOutputOverwrite(bool overwrite)
222  {
223  overwrite_ = overwrite;
224  }
225 
227  void SetConvergeIters(unsigned int citers)
228  {
229  citers_ = citers;
230  }
231 
233  void SetMaxIters(unsigned int iters)
234  {
235  auto params = net_.get_train_params();
236  params.max_iter = iters;
237  net_.set_train_params(params);
238 
239  auto params2 = net2_.get_train_params();
240  params2.max_iter = iters;
241  net2_.set_train_params(params2);
242  }
243 
245  void SetMinLoss(double loss)
246  {
247  auto params = net_.get_train_params();
248  params.min_loss = loss;
249  net_.set_train_params(params);
250 
251  auto params2 = net2_.get_train_params();
252  params2.min_loss = loss;
253  net2_.set_train_params(params2);
254  }
255 
257  void SetRatio(double trainratio)
258  {
259  auto params = net_.get_train_params();
260  params.ratio = trainratio;
261  net_.set_train_params(params);
262 
263  auto params2 = net2_.get_train_params();
264  params2.ratio = trainratio;
265  net2_.set_train_params(params2);
266  }
267 
269  static CFF* Build(
270  const Json::Value& json,
271  const MPI_Comm& world,
272  const MPI_Comm& comm,
273  const std::string& path
274  );
275 
278  {
279  delete fgrid_;
280  delete hgrid_;
281  }
282  };
283 }
Combined Force Frequency (CFF) Algorithm.
Definition: CFF.h:37
unsigned int citers_
Number of iterations after which we turn on full weight.
Definition: CFF.h:46
void SetRatio(double trainratio)
Set training ratio for gradient vs value.
Definition: CFF.h:257
int min_
Definition: CFF.h:122
bool overwrite_
Overwrite outputs.
Definition: CFF.h:106
nnet::neural_net net_
Neural network trained using visit frequency.
Definition: CFF.h:49
unsigned int sweep_
Number of iterations per sweep.
Definition: CFF.h:43
Grid< unsigned int > * hgrid_
Histogram grid.
Definition: CFF.h:82
void SetConvergeIters(unsigned int citers)
Set number of iterations after which we turn on full weight.
Definition: CFF.h:227
double temp_
System temperature and energy units.
Definition: CFF.h:61
~CFF()
Destructor.
Definition: CFF.h:277
double timestep_
Timestep.
Definition: CFF.h:55
Eigen::VectorXi topol_
Neural network topology.
Definition: CFF.h:40
double ratio_
Hold parameters to adjust ratio of 1st and 2nd neural networks (freq vs force-based).
Definition: CFF.h:91
void SetMinLoss(double loss)
Set minimum loss function value (should be zero for production).
Definition: CFF.h:245
nnet::neural_net net2_
Neural network trained on both visit frequency and force.
Definition: CFF.h:52
std::string outfile_
Output filename.
Definition: CFF.h:103
Eigen::ArrayXi Nworld_
Histogram grid for that sotres the number of global hits.
Definition: CFF.h:85
void SetPrevWeight(double h)
Set previous history weight.
Definition: CFF.h:209
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.
Definition: CFF.cpp:407
Eigen::VectorXd Fold_
To hold the last iterations F_ value for removing bias.
Definition: CFF.h:70
std::vector< Grid< double > * > Fworld_
Generalized force grid that stors the global total.
Definition: CFF.h:67
Eigen::VectorXd wdotp2_
To hold second to last iteration wdotp value for numerical derivative.
Definition: CFF.h:76
void SetMaxIters(unsigned int iters)
Set maximum number of training iterations per sweep.
Definition: CFF.h:233
std::vector< double > lowerb_
Bounds.
Definition: CFF.h:97
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.
Definition: CFF.h:151
double pweight_
Previous and current histogram weight.
Definition: CFF.h:58
Grid< double > * ugrid_
Unbiased histogram grid.
Definition: CFF.h:88
void PostIntegration(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call post integration.
Definition: CFF.cpp:79
void WriteBias()
Writes out the bias and CFF output files.
Definition: CFF.cpp:369
Grid< Eigen::VectorXd > * fgrid_
Force grid.
Definition: CFF.h:79
std::vector< double > lowerk_
Bound restraints.
Definition: CFF.h:100
Eigen::MatrixXd force_to_val_ratio_
To hold booleans for training neural network only in specific region for net2_.
Definition: CFF.h:125
void PreSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call prior to simulation initiation.
Definition: CFF.cpp:41
Eigen::MatrixXd hist_
Eigen matrices of grids.
Definition: CFF.h:94
std::vector< Grid< double > * > F_
Generalized force grid that stores total of the local walker.
Definition: CFF.h:64
void PostSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call post simulation.
Definition: CFF.cpp:227
void TrainNetwork()
Trains the neural network.
Definition: CFF.cpp:232
int dim_
Number of CVs in system.
Definition: CFF.h:115
double unitconv_
Unit conversion from mass*velocity/time to force.
Definition: CFF.h:118
void SetOutputOverwrite(bool overwrite)
Set overwrite flag on output file.
Definition: CFF.h:221
Eigen::VectorXd wdotp1_
To hold last iteration wdotp value for numerical derivative.
Definition: CFF.h:73
void SetOutput(const std::string &outfile)
Set name of output file.
Definition: CFF.h:215
Collective variable manager.
Definition: CVManager.h:43
size_t size() const
Get the size of the internal storage vector.
Definition: GridBase.h:326
size_t GetDimension() const
Get the dimension.
Definition: GridBase.h:195
iterator begin()
Return iterator at first grid point.
Definition: Grid.h:527
iterator end()
Return iterator after last valid grid point.
Definition: Grid.h:540
Interface for Method implementations.
Definition: Method.h:44
Class containing a snapshot of the current simulation in time.
Definition: Snapshot.h:48