SSAGES  0.9.3
Software Suite for Advanced General Ensemble Simulations
Meta.h
1 
22 #pragma once
23 
24 #include "Method.h"
25 #include "Grids/Grid.h"
26 #include <fstream>
27 #include <vector>
28 
29 namespace SSAGES
30 {
32 
37  struct Hill
38  {
40  std::vector<double> center;
41 
43  std::vector<double> width;
44 
46  double height;
47 
49 
54  Hill(const std::vector<double>& center,
55  const std::vector<double>& sigma,
56  double height) :
57  center(center), width(sigma), height(height)
58  {}
59  };
60 
62 
68  class Meta : public Method
69  {
70  private:
72  std::vector<Hill> hills_;
73 
75  double height_;
76 
78  std::vector<double> widths_;
79 
82  std::vector<double> derivatives_, tder_, dx_;
84 
86  unsigned int hillfreq_;
87 
90 
93  std::vector<double> upperb_, lowerb_;
95 
98  std::vector<double> upperk_, lowerk_;
100 
102 
106  void AddHill(const CVList& cvs, int iteration);
107 
109 
112  void CalcBiasForce(const CVList& cvs);
113 
115 
119  void PrintHill(const Hill& hill, int iteration);
120 
122  std::ofstream hillsout_;
123 
124  public:
126 
143  Meta(const MPI_Comm& world,
144  const MPI_Comm& comm,
145  double height,
146  const std::vector<double>& widths,
147  const std::vector<double>& lowerb,
148  const std::vector<double>& upperb,
149  const std::vector<double>& lowerk,
150  const std::vector<double>& upperk,
151  Grid<Vector>* grid,
152  unsigned int hillfreq,
153  unsigned int frequency) :
154  Method(frequency, world, comm), hills_(), height_(height), widths_(widths),
155  derivatives_(0), tder_(0), dx_(0), hillfreq_(hillfreq), grid_(grid),
156  upperb_(upperb), lowerb_(lowerb), upperk_(upperk), lowerk_(lowerk)
157  {
158  }
159 
161  void PreSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
162 
164  void PostIntegration(Snapshot* snapshot, const class CVManager& cvmanager) override;
165 
167  void PostSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
168 
170 
177  void LoadHills(const std::string& filename);
178 
180  static Meta* Build(const Json::Value& json,
181  const MPI_Comm& world,
182  const MPI_Comm& comm,
183  const std::string& path);
184 
186  ~Meta() {}
187  };
188 }
189 
Collective variable manager.
Definition: CVManager.h:43
"Vanilla" multi-dimensional Metadynamics.
Definition: Meta.h:69
void PostSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call post simulation.
Definition: Meta.cpp:131
~Meta()
Destructor.
Definition: Meta.h:186
std::vector< double > upperk_
Definition: Meta.h:98
void PostIntegration(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call post integration.
Definition: Meta.cpp:100
Meta(const MPI_Comm &world, const MPI_Comm &comm, double height, const std::vector< double > &widths, const std::vector< double > &lowerb, const std::vector< double > &upperb, const std::vector< double > &lowerk, const std::vector< double > &upperk, Grid< Vector > *grid, unsigned int hillfreq, unsigned int frequency)
Constructor.
Definition: Meta.h:143
std::vector< double > upperb_
Definition: Meta.h:93
static Meta * 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: Meta.cpp:332
unsigned int hillfreq_
Frequency of new hills.
Definition: Meta.h:86
std::vector< Hill > hills_
Hills.
Definition: Meta.h:72
std::vector< double > derivatives_
Definition: Meta.h:82
void PreSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Method call prior to simulation initiation.
Definition: Meta.cpp:62
void LoadHills(const std::string &filename)
Load hills from file.
Definition: Meta.cpp:136
void CalcBiasForce(const CVList &cvs)
Computes the bias force.
Definition: Meta.cpp:253
Grid< Vector > * grid_
CV Grid.
Definition: Meta.h:89
std::vector< double > widths_
Hill widths.
Definition: Meta.h:78
std::ofstream hillsout_
Output stream for hill data.
Definition: Meta.h:122
void PrintHill(const Hill &hill, int iteration)
Prints the new hill to file.
Definition: Meta.cpp:236
double height_
Hill height.
Definition: Meta.h:75
void AddHill(const CVList &cvs, int iteration)
Adds a new hill.
Definition: Meta.cpp:168
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
Multidimensional hill.
Definition: Meta.h:38
std::vector< double > center
Hill center.
Definition: Meta.h:40
Hill(const std::vector< double > &center, const std::vector< double > &sigma, double height)
Constructs a multidimensional Hill (Gaussian)
Definition: Meta.h:54
std::vector< double > width
Hill width.
Definition: Meta.h:43
double height
Hill height.
Definition: Meta.h:46