SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
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 
81  std::vector<double> derivatives_, tder_, dx_;
82 
84  unsigned int hillfreq_;
85 
88 
90  std::vector<double> upperb_, lowerb_;
91 
93  std::vector<double> upperk_, lowerk_;
94 
96 
100  void AddHill(const CVList& cvs, int iteration);
101 
103 
106  void CalcBiasForce(const CVList& cvs);
107 
109 
113  void PrintHill(const Hill& hill, int interation);
114 
116  std::ofstream hillsout_;
117 
118  public:
120 
132  Meta(const MPI_Comm& world,
133  const MPI_Comm& comm,
134  double height,
135  const std::vector<double>& widths,
136  const std::vector<double>& lowerb,
137  const std::vector<double>& upperb,
138  const std::vector<double>& lowerk,
139  const std::vector<double>& upperk,
140  Grid<Vector>* grid,
141  unsigned int hillfreq,
142  unsigned int frequency) :
143  Method(frequency, world, comm), hills_(), height_(height), widths_(widths),
144  derivatives_(0), tder_(0), dx_(0), hillfreq_(hillfreq), grid_(grid),
145  upperb_(upperb), lowerb_(lowerb), upperk_(upperk), lowerk_(lowerk)
146  {
147  }
148 
150 
154  void PreSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
155 
157 
161  void PostIntegration(Snapshot* snapshot, const class CVManager& cvmanager) override;
162 
164 
168  void PostSimulation(Snapshot* snapshot, const class CVManager& cvmanager) override;
169 
171 
178  void LoadHills(const std::string& filename);
179 
181  static Meta* Build(const Json::Value& json,
182  const MPI_Comm& world,
183  const MPI_Comm& comm,
184  const std::string& path);
185 
187  ~Meta() {}
188  };
189 }
190 
Collective variable manager.
Definition: CVManager.h:40
std::vector< CollectiveVariable * > CVList
List of Collective Variables.
Definition: types.h:51
~Meta()
Destructor.
Definition: Meta.h:187
static Meta * Build(const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path)
Definition: Meta.cpp:332
void PreSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Pre-simulation hook.
Definition: Meta.cpp:62
void LoadHills(const std::string &filename)
Load hills from file.
Definition: Meta.cpp:136
Class containing a snapshot of the current simulation in time.
Definition: Snapshot.h:43
unsigned int hillfreq_
Frequency of new hills.
Definition: Meta.h:84
Interface for Method implementations.
Definition: Method.h:43
void PrintHill(const Hill &hill, int interation)
Prints the new hill to file.
Definition: Meta.cpp:236
void PostSimulation(Snapshot *snapshot, const class CVManager &cvmanager) override
Post-simulation hook.
Definition: Meta.cpp:131
std::vector< double > upperb_
Bounds.
Definition: Meta.h:90
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:132
void CalcBiasForce(const CVList &cvs)
Computes the bias force.
Definition: Meta.cpp:253
std::vector< Hill > hills_
Hills.
Definition: Meta.h:72
"Vanilla" multi-dimensional Metadynamics.
Definition: Meta.h:68
std::vector< double > upperk_
Bound restraints.
Definition: Meta.h:93
Multidimensional hill.
Definition: Meta.h:37
std::vector< double > center
Hill center.
Definition: Meta.h:40
std::vector< double > derivatives_
Derivatives and temporary storage vectors.
Definition: Meta.h:81
Grid< Vector > * grid_
CV Grid.
Definition: Meta.h:87
double height
Hill height.
Definition: Meta.h:46
Hill(const std::vector< double > &center, const std::vector< double > &sigma, double height)
Constructs a multidimensional Hill (Gaussian)
Definition: Meta.h:54
void AddHill(const CVList &cvs, int iteration)
Adds a new hill.
Definition: Meta.cpp:168
std::vector< double > width
Hill width.
Definition: Meta.h:43
std::vector< double > widths_
Hill widths.
Definition: Meta.h:78
void PostIntegration(Snapshot *snapshot, const class CVManager &cvmanager) override
Post-integration hook.
Definition: Meta.cpp:100
std::ofstream hillsout_
Output stream for hill data.
Definition: Meta.h:116
double height_
Hill height.
Definition: Meta.h:75