SSAGES  0.9.3
Software Suite for Advanced General Ensemble Simulations
Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
SSAGES::Meta Class Reference

"Vanilla" multi-dimensional Metadynamics. More...

#include <Meta.h>

Inheritance diagram for SSAGES::Meta:
Inheritance graph
[legend]

Public Member Functions

 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. More...
 
void PreSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override
 Method call prior to simulation initiation. More...
 
void PostIntegration (Snapshot *snapshot, const class CVManager &cvmanager) override
 Method call post integration. More...
 
void PostSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override
 Method call post simulation. More...
 
void LoadHills (const std::string &filename)
 Load hills from file. More...
 
 ~Meta ()
 Destructor.
 
- Public Member Functions inherited from SSAGES::Method
 Method (unsigned int frequency, const MPI_Comm &world, const MPI_Comm &comm)
 Constructor. More...
 
void SetCVMask (const std::vector< unsigned int > &mask)
 Sets the collective variable mask. More...
 
virtual ~Method ()
 Destructor.
 
- Public Member Functions inherited from SSAGES::EventListener
 EventListener (unsigned int frequency)
 Constructor. More...
 
unsigned int GetFrequency () const
 Get frequency of event listener. More...
 
virtual ~EventListener ()
 Destructor.
 

Static Public Member Functions

static MetaBuild (const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path)
 Build a derived method from JSON node. More...
 
- Static Public Member Functions inherited from SSAGES::Method
static MethodBuildMethod (const Json::Value &json, const MPI_Comm &world, const MPI_Comm &comm, const std::string &path)
 Build a derived method from JSON node. More...
 
- Static Public Member Functions inherited from SSAGES::EventListener
static unsigned int GetWalkerID (const MPI_Comm &world, const MPI_Comm &comm)
 Get walker ID number of specified communicator. More...
 
static unsigned int GetNumWalkers (const MPI_Comm &world, const MPI_Comm &comm)
 Get total number of walkers in the simulation. More...
 
static bool IsMasterRank (const MPI_Comm &comm)
 Check if current processor is master. More...
 

Private Member Functions

void AddHill (const CVList &cvs, int iteration)
 Adds a new hill. More...
 
void CalcBiasForce (const CVList &cvs)
 Computes the bias force. More...
 
void PrintHill (const Hill &hill, int iteration)
 Prints the new hill to file. More...
 

Private Attributes

std::vector< Hillhills_
 Hills.
 
double height_
 Hill height.
 
std::vector< double > widths_
 Hill widths.
 
unsigned int hillfreq_
 Frequency of new hills.
 
Grid< Vector > * grid_
 CV Grid.
 
std::ofstream hillsout_
 Output stream for hill data.
 
std::vector< double > derivatives_
 
std::vector< double > tder_
 
std::vector< double > dx_
 
std::vector< double > upperb_
 
std::vector< double > lowerb_
 
std::vector< double > upperk_
 
std::vector< double > lowerk_
 

Additional Inherited Members

- Protected Attributes inherited from SSAGES::Method
mxx::comm world_
 Global MPI communicator.
 
mxx::comm comm_
 Local MPI communicator.
 
std::vector< unsigned int > cvmask_
 Mask which identifies which CVs to act on.
 

Detailed Description

"Vanilla" multi-dimensional Metadynamics.

Implementation of a "vanilla" multi-dimensional Metadynamics method with no bells and whistles.

Definition at line 68 of file Meta.h.

Constructor & Destructor Documentation

◆ Meta()

SSAGES::Meta::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 
)
inline

Constructor.

Parameters
worldMPI global communicator.
commMPI local communicator.
heightHeight of the hills to be deposited.
widthsWidth of the hills to be deposited.
lowerbLower bounds for CVs.
upperbUpper bounds for CVs.
lowerkLower bound restraints for CVs.
upperkUpper bound restraints for CVs.
gridGrid containing bias.
hillfreqFrequency of depositing hills.
frequencyFrequency of invoking this method.

Constructs an instance of Metadynamics method.

Note
The size of "widths" should be commensurate with the number of CV's expected.

Definition at line 143 of file Meta.h.

153  :
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  }
std::vector< double > upperk_
Definition: Meta.h:98
std::vector< double > upperb_
Definition: Meta.h:93
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
Grid< Vector > * grid_
CV Grid.
Definition: Meta.h:89
std::vector< double > widths_
Hill widths.
Definition: Meta.h:78
double height_
Hill height.
Definition: Meta.h:75
Method(unsigned int frequency, const MPI_Comm &world, const MPI_Comm &comm)
Constructor.
Definition: Method.h:61

Member Function Documentation

◆ AddHill()

void SSAGES::Meta::AddHill ( const CVList cvs,
int  iteration 
)
private

Adds a new hill.

Parameters
cvsList of CVs.
iterationCurrent iteration.

Definition at line 168 of file Meta.cpp.

169  {
170  size_t n = cvs.size();
171 
172  // Assume we have the same number of procs per walker.
173  int nwalkers = world_.size()/comm_.size();
174 
175  // We need to exchange CV values across the walkers
176  // and to each proc on a walker.
177  std::vector<double> cvals(n*nwalkers, 0);
178 
179  if(IsMasterRank(comm_))
180  {
181  for(size_t i = 0, j = GetWalkerID(world_,comm_)*n; i < n; ++i,++j)
182  cvals[j] = cvs[i]->GetValue();
183  }
184 
185  // Reduce across all processors and add hills.
186  MPI_Allreduce(MPI_IN_PLACE, cvals.data(), n*nwalkers, MPI_DOUBLE, MPI_SUM, world_);
187 
188  for(size_t i = 0; i < n*nwalkers; i += n)
189  {
190  std::vector<double> cval(cvals.begin() + i, cvals.begin() + i + n);
191  hills_.emplace_back(cval, widths_, height_);
192 
193  // Write hill to file.
194  if(IsMasterRank(world_))
195  PrintHill(hills_.back(), iteration);
196  }
197 
198  // If grid is defined, add bias onto grid.
199  if(grid_ != nullptr)
200  {
201  std::vector<double> dx(n, 0.0), df(n, 1.0);
202  auto& hill = hills_.back();
203  for(auto it = grid_->begin(); it != grid_->end(); ++it)
204  {
205  auto& val = *it;
206  auto coord = it.coordinates();
207 
208  // Compute difference between grid point and current val.
209  for(size_t i = 0; i < n; ++i)
210  {
211  dx[i] = -cvs[i]->GetDifference(coord[i]);
212  df[i] = 1.;
213  }
214 
215  // Compute derivative.
216  for(size_t i = 0; i < n; ++i)
217  {
218  for(size_t j = 0; j < n; ++j)
219  {
220  if(j != i)
221  df[i] *= gaussian(dx[j], hill.width[j]);
222  else
223  df[i] *= gaussianDerv(dx[j], hill.width[j]);
224  }
225  }
226 
227  // Add to grid.
228  for(size_t i = 0; i < n; ++i)
229  val[i] += height_*df[i];
230  }
231  }
232  }
static bool IsMasterRank(const MPI_Comm &comm)
Check if current processor is master.
static unsigned int GetWalkerID(const MPI_Comm &world, const MPI_Comm &comm)
Get walker ID number of specified communicator.
Definition: EventListener.h:89
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
void PrintHill(const Hill &hill, int iteration)
Prints the new hill to file.
Definition: Meta.cpp:236
mxx::comm comm_
Local MPI communicator.
Definition: Method.h:47
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
double gaussianDerv(double dx, double sigma)
Evaluate Gaussian derivative. Helper function.
Definition: Meta.cpp:54
double gaussian(double dx, double sigma)
Evlauate Gaussian. Helper function.
Definition: Meta.cpp:42

References SSAGES::gaussian(), and SSAGES::gaussianDerv().

Here is the call graph for this function:

◆ Build()

Meta * SSAGES::Meta::Build ( const Json::Value &  json,
const MPI_Comm &  world,
const MPI_Comm &  comm,
const std::string &  path 
)
static

Build a derived method from JSON node.

Parameters
jsonJSON Value containing all input information.
worldMPI global communicator.
commMPI local communicator.
pathPath for JSON path specification.
Returns
Pointer to the Method built. nullptr if an unknown error occurred.

This function builds a registered method from a JSON node. The difference between this function and "Build" is that this automatically determines the appropriate derived type based on the JSON node information.

Note
Object lifetime is the caller's responsibility.

Definition at line 332 of file Meta.cpp.

336  {
337  ObjectRequirement validator;
338  Value schema;
339  CharReaderBuilder rbuilder;
340  CharReader* reader = rbuilder.newCharReader();
341 
342  reader->parse(JsonSchema::MetadynamicsMethod.c_str(),
343  JsonSchema::MetadynamicsMethod.c_str() + JsonSchema::MetadynamicsMethod.size(),
344  &schema, nullptr);
345  validator.Parse(schema, path);
346 
347  // Validate inputs.
348  validator.Validate(json, path);
349  if(validator.HasErrors())
350  throw BuildException(validator.GetErrors());
351 
352  std::vector<double> widths;
353  for(auto& s : json["widths"])
354  widths.push_back(s.asDouble());
355 
356  std::vector<double> lowerb, upperb, lowerk, upperk;
357  Grid<Vector>* grid = nullptr;
358  if(json.isMember("grid"))
359  grid = Grid<Vector>::BuildGrid(json.get("grid", Json::Value()));
360  else if(!json.isMember("lower_bounds") || !json.isMember("upper_bounds"))
361  throw BuildException({
362  "#/Method/Metadynamics: Both upper_bounds and lower_bounds "
363  "must be defined if grid is not being used."});
364 
365  // Assume all vectors are the same size.
366  for(int i = 0; i < static_cast<int>(json["lower_bound_restraints"].size()); ++i)
367  {
368  lowerk.push_back(json["lower_bound_restraints"][i].asDouble());
369  upperk.push_back(json["upper_bound_restraints"][i].asDouble());
370  lowerb.push_back(json["lower_bounds"][i].asDouble());
371  upperb.push_back(json["upper_bounds"][i].asDouble());
372  }
373 
374  auto height = json.get("height", 1.0).asDouble();
375  auto hillfreq = json.get("hill_frequency", 1).asInt();
376  auto freq = json.get("frequency", 1).asInt();
377 
378  auto* m = new Meta(
379  world, comm, height, widths,
380  lowerb, upperb, lowerk, upperk,
381  grid, hillfreq, freq
382  );
383 
384  if(json.isMember("load_hills"))
385  m->LoadHills(json["load_hills"].asString());
386 
387  return m;
388  }
Requirements on an object.
virtual void Parse(Value json, const std::string &path) override
Parse JSON value to generate Requirement(s).
virtual void Validate(const Value &json, const std::string &path) override
Validate JSON value.
std::vector< std::string > GetErrors()
Get list of error messages.
Definition: Requirement.h:92
bool HasErrors()
Check if errors have occured.
Definition: Requirement.h:86
static Grid< T > * BuildGrid(const Json::Value &json)
Set up the grid.
Definition: Grid.h:127
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

References Json::Requirement::GetErrors(), Json::Requirement::HasErrors(), Json::ObjectRequirement::Parse(), and Json::ObjectRequirement::Validate().

Here is the call graph for this function:

◆ CalcBiasForce()

void SSAGES::Meta::CalcBiasForce ( const CVList cvs)
private

Computes the bias force.

Parameters
cvsList of CVs.

Definition at line 253 of file Meta.cpp.

254  {
255  // Reset bias and derivatives.
256  double bias = 0.;
257  auto n = cvs.size();
258 
259  // Reset vectors.
260  std::fill(derivatives_.begin(), derivatives_.end(), 0);
261 
262  // Look up and apply grid bias.
263  if(grid_ != nullptr)
264  {
265  bool inbounds = true;
266  std::vector<double> val(n, 0.);
267  for(size_t i = 0; i < n; ++i)
268  {
269  val[i] = cvs[i]->GetValue();
270  if(val[i] < grid_->GetLower(i) || val[i] > grid_->GetUpper(i))
271  inbounds = false;
272  }
273 
274  if(inbounds)
275  {
276  auto frc = (*grid_)[val];
277  for(size_t i = 0; i < n; ++i)
278  derivatives_[i] = frc[i];
279  }
280  else
281  {
282  if(IsMasterRank(comm_))
283  {
284  std::cerr << "Metadynamics: out of bounds ( ";
285  for(auto& v : val)
286  std::cerr << v << " ";
287  std::cerr << ")" << std::endl;
288  }
289  }
290  }
291  else
292  {
293  // Loop through hills and calculate the bias force.
294  for(auto& hill : hills_)
295  {
296  auto tbias = 1.;
297  std::fill(tder_.begin(), tder_.end(), 1.0);
298  std::fill(dx_.begin(), dx_.end(), 1.0);
299 
300  for(size_t i = 0; i < n; ++i)
301  {
302  dx_[i] = cvs[i]->GetDifference(hill.center[i]);
303  tbias *= gaussian(dx_[i], hill.width[i]);
304  }
305 
306  for(size_t i = 0; i < n; ++i)
307  for(size_t j = 0; j < n; ++j)
308  {
309  if(j != i)
310  tder_[i] *= gaussian(dx_[j], hill.width[j]);
311  else
312  tder_[i] *= gaussianDerv(dx_[j], hill.width[j]);
313  }
314 
315  bias += height_ * tbias;
316  for(size_t i = 0; i < n; ++i)
317  derivatives_[i] += height_*tder_[i];
318  }
319  }
320 
321  // Restraints.
322  for(size_t i = 0; i < n; ++i)
323  {
324  auto cval = cvs[i]->GetValue();
325  if(cval < lowerb_[i])
326  derivatives_[i] += lowerk_[i]*(cval - lowerb_[i]);
327  else if(cval > upperb_[i])
328  derivatives_[i] += upperk_[i]*(cval - upperb_[i]);
329  }
330  }
const std::vector< double > GetLower() const
Return the lower edges of the Grid.
Definition: GridBase.h:231
const std::vector< double > GetUpper() const
Return the upper edges of the Grid.
Definition: GridBase.h:262

References SSAGES::gaussian(), and SSAGES::gaussianDerv().

Here is the call graph for this function:

◆ LoadHills()

void SSAGES::Meta::LoadHills ( const std::string &  filename)

Load hills from file.

Parameters
filenameFile name containing hills.
Note
File format must match the output written by this method and the dimensionality of the hills must match the initialized Meta class.

Definition at line 136 of file Meta.cpp.

137  {
138  std::ifstream file(filename);
139  std::string line;
140 
141  auto dim = widths_.size();
142  double iteration = 0, height = 0;
143  std::vector<double> width(dim, 0.), center(dim, 0);
144 
145  // First line is a comment.
146  std::getline(file, line);
147  while(std::getline(file, line))
148  {
149  std::istringstream iss(line);
150  iss >> iteration; // Not really using this.
151 
152  // Load centers.
153  for(size_t i = 0; i < dim; ++i)
154  iss >> center[i];
155 
156  // Load widths.
157  for(size_t i = 0; i < dim; ++i)
158  iss >> width[i];
159 
160  // Load height.
161  iss >> height;
162 
163  hills_.emplace_back(center, width, height);
164  }
165  }

◆ PostIntegration()

void SSAGES::Meta::PostIntegration ( Snapshot snapshot,
const class CVManager cvmanager 
)
overridevirtual

Method call post integration.

Parameters
snapshotPointer to the simulation snapshot.
cvmanagerCollective variable manager.

This function will be called after each integration step.

Implements SSAGES::Method.

Definition at line 100 of file Meta.cpp.

101  {
102  auto cvs = cvmanager.GetCVs(cvmask_);
103  // Add hills when needed.
104  if(snapshot->GetIteration() % hillfreq_ == 0)
105  AddHill(cvs, snapshot->GetIteration());
106 
107  // Always calculate the current bias.
108  CalcBiasForce(cvs);
109 
110  // Take each CV and add its biased forces to the atoms
111  // using the chain rule.
112  auto& forces = snapshot->GetForces();
113  auto& virial = snapshot->GetVirial();
114 
115  for(size_t i = 0; i < cvs.size(); ++i)
116  {
117  auto& grad = cvs[i]->GetGradient();
118  auto& boxgrad = cvs[i]->GetBoxGradient();
119 
120  // Update the forces in snapshot by adding in the force bias from each
121  // CV to each atom based on the gradient of the CV.
122  for (size_t j = 0; j < forces.size(); ++j)
123  for(size_t k = 0; k < 3; ++k)
124  forces[j][k] -= derivatives_[i]*grad[j][k];
125 
126  virial += derivatives_[i]*boxgrad;
127  }
128  }
void CalcBiasForce(const CVList &cvs)
Computes the bias force.
Definition: Meta.cpp:253
void AddHill(const CVList &cvs, int iteration)
Adds a new hill.
Definition: Meta.cpp:168
std::vector< unsigned int > cvmask_
Mask which identifies which CVs to act on.
Definition: Method.h:50

References SSAGES::CVManager::GetCVs(), SSAGES::Snapshot::GetForces(), SSAGES::Snapshot::GetIteration(), and SSAGES::Snapshot::GetVirial().

Here is the call graph for this function:

◆ PostSimulation()

void SSAGES::Meta::PostSimulation ( Snapshot snapshot,
const class CVManager cvmanager 
)
overridevirtual

Method call post simulation.

Parameters
snapshotPointer to the simulation snapshot.
cvmanagerCollective variable manager.

This function will be called after the end of the simulation run.

Implements SSAGES::Method.

Definition at line 131 of file Meta.cpp.

132  {
133  }

◆ PreSimulation()

void SSAGES::Meta::PreSimulation ( Snapshot snapshot,
const class CVManager cvmanager 
)
overridevirtual

Method call prior to simulation initiation.

Parameters
snapshotPointer to the simulation snapshot.
cvmanagerCollective variable manager.

This function will be called before the simulation is started.

Implements SSAGES::Method.

Definition at line 62 of file Meta.cpp.

63  {
64  auto cvs = cvmanager.GetCVs(cvmask_);
65  // Write ouput file header.
66  if(IsMasterRank(world_))
67  {
68  hillsout_.open("hills.out");
69  hillsout_ << "#Iteration ";
70 
71  for(size_t i = 0; i < cvs.size(); ++i)
72  hillsout_ << "center." << i << " ";
73 
74  for(size_t i = 0; i < cvs.size(); ++i)
75  hillsout_ << "sigma." << i << " ";
76 
77  hillsout_ << "height" << std::endl;
78 
79  hillsout_.close();
80  }
81 
82  // Initialize grid to zero.
83  if(grid_ != nullptr)
84  {
85  Vector vec = Vector::Zero(cvs.size());
86  std::fill(grid_->begin(), grid_->end(), vec);
87  }
88 
89  auto n = snapshot->GetTargetIterations();
90  n = n ? n : 1e5; // Pre-allocate at least something.
91 
92  hills_.reserve(n+1);
93  widths_.reserve(n+1);
94  derivatives_.resize(cvs.size());
95  tder_.resize(cvs.size());
96  dx_.resize(cvs.size());
97  }
std::ofstream hillsout_
Output stream for hill data.
Definition: Meta.h:122
Eigen::VectorXd Vector
Arbitrary length vector.
Definition: types.h:30

References SSAGES::CVManager::GetCVs(), and SSAGES::Snapshot::GetTargetIterations().

Here is the call graph for this function:

◆ PrintHill()

void SSAGES::Meta::PrintHill ( const Hill hill,
int  iteration 
)
private

Prints the new hill to file.

Parameters
hillHill to be printed.
iterationCurrent iteration.

Definition at line 236 of file Meta.cpp.

237  {
238  hillsout_.open("hills.out", std::fstream::app);
239 
240  hillsout_ << iteration << " ";
241  hillsout_.precision(8);
242 
243  for(auto& cv : hill.center)
244  hillsout_ << cv << " ";
245 
246  for(auto& w : hill.width)
247  hillsout_ << w << " ";
248 
249  hillsout_ << height_ << std::endl;
250  hillsout_.close();
251  }

References SSAGES::Hill::center, and SSAGES::Hill::width.

Member Data Documentation

◆ derivatives_

std::vector<double> SSAGES::Meta::derivatives_
private

Derivatives and temporary storage vectors.

Definition at line 82 of file Meta.h.

◆ upperb_

std::vector<double> SSAGES::Meta::upperb_
private

Bounds

Definition at line 93 of file Meta.h.

◆ upperk_

std::vector<double> SSAGES::Meta::upperk_
private

Bound restraints.

Definition at line 98 of file Meta.h.


The documentation for this class was generated from the following files: