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::BFS Class Reference

Basis Function Sampling Algorithm. More...

#include <BasisFunc.h>

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

Public Member Functions

 BFS (const MPI_Comm &world, const MPI_Comm &comm, Grid< unsigned int > *h, Grid< std::vector< double >> *f, Grid< double > *b, const std::vector< BasisFunction * > &functions, 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 double temperature, const double tol, const double weight, bool converge)
 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 SetIteration (const int iter)
 Set the current iteration. More...
 
void SetBasis (const std::vector< double > &coeff, std::vector< double > &unbias)
 Set the values for the basis. More...
 
 ~BFS ()
 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 BFSBuild (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 ProjectBias (const CVList &cvs, const double beta)
 The functions which calculates the updated bias and coefficients and then stores them. More...
 
void InBounds (const CVList &cvs)
 A function which checks to see if the CVs are still in bounds. More...
 
void PrintBias (const CVList &cvs, const double beta)
 Prints the current bias to a defined file from the JSON. More...
 

Private Attributes

Grid< unsigned int > * h_
 Grid of visited states. More...
 
Grid< double > * b_
 Stored bias potential. More...
 
Grid< std::vector< double > > * f_
 Stored gradients. More...
 
std::vector< double > unbias_
 The biased histogram of states. More...
 
std::vector< double > coeffArr_
 The coefficient array for restart runs.
 
BasisEvaluator evaluator_
 The basis evaluator class. More...
 
std::vector< double > restraint_
 Spring constants for restrained system. More...
 
std::vector< double > boundUp_
 Upper position of the spring restraint.
 
std::vector< double > boundLow_
 Lower position of the spring restraint.
 
unsigned int cyclefreq_
 Frequency of coefficient updates.
 
unsigned int cyclefreqTmp_
 Temporary cyclefrequency that ensures sweeps are of proper length for restrained simulations.
 
unsigned int step_
 Step counter for the cyclefrequency.
 
unsigned int mpiid_
 The node that the current system belongs to, primarily for printing and debugging.
 
double weight_
 Weighting for potentially faster sampling. More...
 
double temperature_
 Self-defined temperature. More...
 
double tol_
 The tolerance criteria for the system .
 
bool bounds_
 A variable to check to see if the simulation is in bounds or not.
 
bool convergeExit_
 A check to see if you want the system to end when it reaches the convergence criteria.
 
std::string bnme_
 
unsigned int iteration_
 Iteration counter.
 

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

Basis Function Sampling Algorithm.

Implementation of the Basis Function Sampling Method based on [5].

Definition at line 38 of file BasisFunc.h.

Constructor & Destructor Documentation

◆ BFS()

SSAGES::BFS::BFS ( const MPI_Comm &  world,
const MPI_Comm &  comm,
Grid< unsigned int > *  h,
Grid< std::vector< double >> *  f,
Grid< double > *  b,
const std::vector< BasisFunction * > &  functions,
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 double  temperature,
const double  tol,
const double  weight,
bool  converge 
)
inline

Constructor.

Parameters
worldMPI global communicator.
commMPI local communicator.
hGrid containing histogram.
fGrid containing histogram.
bGrid containing histogram.
functionsVector of basis functions.
restraintRestraint spring constants.
boundUpUpper bounds of restraint springs.
boundLowLower bounds of restraint springs.
cyclefreqCycle frequency.
frequencyFrequency with which this Method is applied.
bnmeBasis file name.
temperatureAutomatically set temperature.
tolThreshold for tolerance criterion.
weightWeight for improved sampling.
convergeIf True quit on convergence.

Constructs an instance of the Basis function sampling method. The coefficients describe the basis projection of the system. This is updated once every cyclefreq_. For now, only the Legendre polynomial is implemented. Others will be added later.

Definition at line 181 of file BasisFunc.h.

196  :
197  Method(frequency, world, comm),
198  h_(h), b_(b), f_(f), unbias_(), coeffArr_(), evaluator_(functions),
199  restraint_(restraint), boundUp_(boundUp), boundLow_(boundLow),
200  cyclefreq_(cyclefreq), mpiid_(0), weight_(weight),
201  temperature_(temperature), tol_(tol),
202  convergeExit_(converge), bnme_(bnme), iteration_(0)
203  {
204  }
double weight_
Weighting for potentially faster sampling.
Definition: BasisFunc.h:110
bool convergeExit_
A check to see if you want the system to end when it reaches the convergence criteria.
Definition: BasisFunc.h:127
std::vector< double > restraint_
Spring constants for restrained system.
Definition: BasisFunc.h:84
double temperature_
Self-defined temperature.
Definition: BasisFunc.h:118
Grid< unsigned int > * h_
Grid of visited states.
Definition: BasisFunc.h:46
Grid< std::vector< double > > * f_
Stored gradients.
Definition: BasisFunc.h:58
BasisEvaluator evaluator_
The basis evaluator class.
Definition: BasisFunc.h:76
Grid< double > * b_
Stored bias potential.
Definition: BasisFunc.h:52
std::string bnme_
Definition: BasisFunc.h:151
std::vector< double > coeffArr_
The coefficient array for restart runs.
Definition: BasisFunc.h:69
std::vector< double > boundUp_
Upper position of the spring restraint.
Definition: BasisFunc.h:87
double tol_
The tolerance criteria for the system .
Definition: BasisFunc.h:121
unsigned int cyclefreq_
Frequency of coefficient updates.
Definition: BasisFunc.h:93
unsigned int iteration_
Iteration counter.
Definition: BasisFunc.h:154
std::vector< double > boundLow_
Lower position of the spring restraint.
Definition: BasisFunc.h:90
unsigned int mpiid_
The node that the current system belongs to, primarily for printing and debugging.
Definition: BasisFunc.h:102
std::vector< double > unbias_
The biased histogram of states.
Definition: BasisFunc.h:66
Method(unsigned int frequency, const MPI_Comm &world, const MPI_Comm &comm)
Constructor.
Definition: Method.h:61

Member Function Documentation

◆ Build()

BFS * SSAGES::BFS::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 308 of file BasisFunc.cpp.

312  {
313  ObjectRequirement validator;
314  Value schema;
315  CharReaderBuilder rbuilder;
316  CharReader* reader = rbuilder.newCharReader();
317 
318  reader->parse(JsonSchema::BFSMethod.c_str(),
319  JsonSchema::BFSMethod.c_str() + JsonSchema::BFSMethod.size(),
320  &schema, nullptr);
321  validator.Parse(schema, path);
322 
323  //Validate Inputs
324  validator.Validate(json, path);
325  if(validator.HasErrors())
326  throw BuildException(validator.GetErrors());
327 
328  std::vector<double> restrCV(0);
329  for(auto& restr : json["CV_restraint_spring_constants"])
330  restrCV.push_back(restr.asDouble());
331 
332  std::vector<double> boundLow(0);
333  for(auto& bndl : json["CV_restraint_minimums"])
334  boundLow.push_back(bndl.asDouble());
335 
336  std::vector<double> boundUp(0);
337  for(auto& bndu : json["CV_restraint_maximums"])
338  boundUp.push_back(bndu.asDouble());
339 
340  auto cyclefreq = json.get("cycle_frequency", 100000).asInt();
341  auto freq = json.get("frequency", 1).asInt();
342  auto wght = json.get("weight", 1.0).asDouble();
343  auto bnme = json.get("basis_filename", "").asString();
344  auto temp = json.get("temperature", 0.0).asDouble();
345  auto tol = json.get("tolerance", 1e-6).asDouble();
346  auto conv = json.get("convergence_exit", false).asBool();
347 
348  Grid<unsigned int> *h = Grid<unsigned int>::BuildGrid(
349  json.get("grid", Json::Value()) );
350 
351  Grid<std::vector<double>> *f = Grid<std::vector<double>>::BuildGrid(
352  json.get("grid", Json::Value()) );
353 
354  Grid<double> *b = Grid<double>::BuildGrid(
355  json.get("grid", Json::Value()) );
356 
357  size_t ii = 0;
358  std::vector<BasisFunction*> functions;
359  for(auto& m : json["basis_functions"])
360  {
361  auto *bf = BasisFunction::Build(m, path, b->GetNumPoints(ii));
362  functions.push_back(bf);
363  ii++;
364  }
365 
366  // If restraints were not specified then set them equal to zero
367  if (restrCV.size() == 0) {
368  // Choose simplest array that should have equivalent size
369  for(size_t i = 0; i<functions.size(); i++)
370  restrCV.push_back(0); // Push back value of 0
371  }
372 
373  // This also needs to be the same for upper bounds and lower bounds
374  if (boundLow.size() == 0) {
375  // Choose simplest array that should have equivalent size
376  for(size_t i = 0; i<functions.size(); i++)
377  boundLow.push_back(0); // Push back value of 0
378  }
379 
380  // This also needs to be the same for upper bounds and lower bounds
381  if (boundUp.size() == 0) {
382  // Choose simplest array that should have equivalent size
383  for(size_t i = 0; i<functions.size(); i++)
384  boundUp.push_back(0); // Push back value of 0
385  }
386 
387  // Check to make sure that the arrays for boundaries and restraints are of the same length
388  if (boundUp.size() != boundLow.size() && boundLow.size() != restrCV.size()) {
389  std::cerr<<"ERROR: Number of restraint boundaries do not match number of spring constants."<<std::endl;
390  std::cerr<<"Boundary size is " << boundUp.size() << " boundary low size is " << boundLow.size() << " restraints size is" << restrCV.size() << std::endl;
391  MPI_Abort(world, EXIT_FAILURE);
392  }
393 
394  auto* m = new BFS(world, comm, h, f, b, functions, restrCV, boundUp, boundLow,
395  cyclefreq, freq, bnme, temp, tol, wght,
396  conv);
397 
398  if(json.isMember("iteration"))
399  m->SetIteration(json.get("iteration",0).asInt());
400 
401  // Check if previously saved grids exist. If so, check that data match and load grids.
402  std::ifstream restrFile;
403  restrFile.open("restart"+bnme+".out");
404  if(restrFile.is_open())
405  {
406  std::string line;
407  std::vector<double> coeff;
408  std::vector<double> unbias;
409  bool coeffFlag = false;
410  bool basisFlag = false;
411  std::cout << "Attempting to load data from a previous run of BFS." << std::endl;
412 
413  while(!std::getline(restrFile,line).eof())
414  {
415  if(line.find("COEFFICIENTS") != std::string::npos) {
416  std::getline(restrFile,line);
417  coeffFlag = true;
418  }
419  else if (line.find("HISTOGRAM") != std::string::npos) {
420  std::getline(restrFile,line);
421  basisFlag = true;
422  coeffFlag = false;
423  }
424  if(coeffFlag) {
425  coeff.push_back(std::stod(line));
426  }
427  else if (basisFlag) {
428  unbias.push_back(std::stod(line));
429  }
430  }
431  m->SetBasis(coeff, unbias);
432  restrFile.close();
433  }
434 
435  return m;
436  }
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
BFS(const MPI_Comm &world, const MPI_Comm &comm, Grid< unsigned int > *h, Grid< std::vector< double >> *f, Grid< double > *b, const std::vector< BasisFunction * > &functions, 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 double temperature, const double tol, const double weight, bool converge)
Constructor.
Definition: BasisFunc.h:181
static BasisFunction * Build(const Json::Value &json, const std::string &path, unsigned int nbins)
Build BasisFunction from JSON value.
Definition: Basis.cpp:9
static Grid< T > * BuildGrid(const Json::Value &json)
Set up the grid.
Definition: Grid.h:127

References Json::Requirement::GetErrors(), SSAGES::GridBase< T >::GetNumPoints(), Json::Requirement::HasErrors(), Json::ObjectRequirement::Parse(), and Json::ObjectRequirement::Validate().

Here is the call graph for this function:

◆ InBounds()

void SSAGES::BFS::InBounds ( const CVList cvs)
private

A function which checks to see if the CVs are still in bounds.

Parameters
cvsList of CVs to check

Definition at line 272 of file BasisFunc.cpp.

273  {
274  std::vector<double> x (cvs.size(),0);
275  //This is calculating the derivatives for the bias force
276  for (size_t j = 0; j < cvs.size(); ++j)
277  {
278  x[j] = cvs[j]->GetValue();
279  double min = h_->GetLower(j);
280  double max = h_->GetUpper(j);
281 
282  if(!h_->GetPeriodic(j))
283  {
284  // In order to prevent the index for the histogram from going out of bounds a check is in place
285  if(x[j] > max && bounds_)
286  {
287  //std::cout<<"WARNING: CV is above the maximum boundary."<<std::endl;
288  //std::cout<<"Statistics will not be gathered during this interval"<<std::endl;
289  bounds_ = false;
290  break;
291  }
292  else if(x[j] < min && bounds_)
293  {
294  //std::cout<<"WARNING: CV is below the minimum boundary."<<std::endl;
295  //std::cout<<"Statistics will not be gathered during this interval"<<std::endl;
296  bounds_ = false;
297  break;
298  }
299  else if(x[j] < max && x[j] > min && !bounds_)
300  {
301  //std::cout<<"CV has returned in between bounds. Run is resuming"<<std::endl;
302  bounds_ = true;
303  }
304  }
305  }
306  }
bool bounds_
A variable to check to see if the simulation is in bounds or not.
Definition: BasisFunc.h:124
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
const std::vector< bool > & GetPeriodic() const
Return the periodicity of the Grid.
Definition: GridBase.h:295

◆ PostIntegration()

void SSAGES::BFS::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 81 of file BasisFunc.cpp.

82  {
83  auto cvs = cvmanager.GetCVs(cvmask_);
84  std::vector<double> x(cvs.size(),0);
85 
86  /*The binned cv space is updated at every step
87  *After a certain number of steps has been passed, the system updates a
88  *bias projection based on the visited histogram states
89  */
90  for(size_t i = 0; i < cvs.size(); ++i)
91  {
92  x[i] = cvs[i]->GetValue();
93  }
94 
95  // The histogram is updated based on the index
96  InBounds(cvs);
97  if (bounds_) {h_->at(x) += 1;}
98  else {cyclefreq_++;}
99 
100  // Update the basis projection after the requisite number of steps
101  int steps = snapshot->GetIteration()-step_;
102  if(steps % cyclefreq_ == 0) {
103  step_ = snapshot->GetIteration();
105  double beta;
106  beta = 1.0 / (temperature_ * snapshot->GetKb());
107 
108  // For systems with poorly defined temperature (ie: 1 particle) the
109  // user needs to define their own temperature. This is a hack that
110  // will be removed in future versions.
111 
112  iteration_+= 1;
113  ProjectBias(cvs,beta);
114  std::cout<<"Node: ["<<mpiid_<<"]"<<std::setw(10)<<"\tSweep: "<<iteration_<<std::endl;
115  }
116 
117  // This gets the bias force from the grid
118  std::vector<double> bias_grad(cvs.size(),0);
119 
120  if (bounds_)
121  {
122  auto& bias = f_->at(x);
123  for (size_t ii = 0; ii<bias.size(); ii++)
124  bias_grad[ii] = bias[ii];
125  }
126  else
127  {
128  // This is where the wall potentials are going to be thrown into the method if the system is not a periodic CV
129  for(size_t j = 0; j < cvs.size(); ++j)
130  {
131  if(!h_->GetPeriodic(j))
132  {
133  if(x[j] > boundUp_[j])
134  bias_grad[j] = -restraint_[j] * (x[j] - boundUp_[j]);
135  else if(x[j] < boundLow_[j])
136  bias_grad[j] = restraint_[j] * (x[j] - boundLow_[j]);
137  }
138  }
139  }
140 
141  // Take each CV and add its biased forces to the atoms using the chain rule
142  auto& forces = snapshot->GetForces();
143  auto& virial = snapshot->GetVirial();
144  for(size_t i = 0; i < cvs.size(); ++i)
145  {
146  auto& grad = cvs[i]->GetGradient();
147  auto& boxgrad = cvs[i]->GetBoxGradient();
148 
149  /* Update the forces in snapshot by adding in the force bias from each
150  *CV to each atom based on the gradient of the CV.
151  */
152  for (size_t j = 0; j < forces.size(); ++j)
153  forces[j] += bias_grad[i]*grad[j];
154 
155  virial -= bias_grad[i]*boxgrad;
156  }
157  }
unsigned int cyclefreqTmp_
Temporary cyclefrequency that ensures sweeps are of proper length for restrained simulations.
Definition: BasisFunc.h:96
void ProjectBias(const CVList &cvs, const double beta)
The functions which calculates the updated bias and coefficients and then stores them.
Definition: BasisFunc.cpp:166
void InBounds(const CVList &cvs)
A function which checks to see if the CVs are still in bounds.
Definition: BasisFunc.cpp:272
unsigned int step_
Step counter for the cyclefrequency.
Definition: BasisFunc.h:99
const T & at(const std::vector< int > &indices) const
Access Grid element read-only.
Definition: GridBase.h:546
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(), SSAGES::Snapshot::GetKb(), and SSAGES::Snapshot::GetVirial().

Here is the call graph for this function:

◆ PostSimulation()

void SSAGES::BFS::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 160 of file BasisFunc.cpp.

161  {
162  std::cout<<"Run has finished"<<std::endl;
163  }

◆ PreSimulation()

void SSAGES::BFS::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 36 of file BasisFunc.cpp.

37  {
38  auto cvs = cvmanager.GetCVs(cvmask_);
39  // For print statements and file I/O, the walker IDs are used
40  mpiid_ = snapshot->GetWalkerID();
41 
42  // Make sure the iteration index is set correctly
43  iteration_ = 0;
44  step_ = 0;
45 
46  // Set the temporary sweep frequency to the same as the input
48 
49  // There are a few error messages / checks that are in place with
50  // defining CVs and grids
51  if(h_->GetDimension() != cvs.size())
52  {
53  std::cerr<<"ERROR: Grid dimensions doesn't match number of CVS."<<std::endl;
54  MPI_Abort(world_, EXIT_FAILURE);
55  }
56 
57  // This is to check for non-periodic bounds. It comes into play in the update bias function
58  bounds_ = true;
59  unbias_.resize(h_->size(),0);
60  std::fill(unbias_.begin(),unbias_.end(),1.0);
61 
62  // Initialize the mapping for the hist function
63  for (auto &val : *h_) {
64  val = 0;
65  }
66 
67  // Reset the values in the force grid
68  for (auto &val : *f_) {
69  for (size_t i = 0; i <val.size(); i++) {
70  val[i] = 0;
71  }
72  }
73 
74  // Reset the values in the bias potential grid
75  for (auto &val : *b_) {
76  val = 0;
77  }
78  }
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
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46

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

Here is the call graph for this function:

◆ PrintBias()

void SSAGES::BFS::PrintBias ( const CVList cvs,
const double  beta 
)
private

Prints the current bias to a defined file from the JSON.

Parameters
cvsList of collective variables.
betaScale parameter.

Definition at line 227 of file BasisFunc.cpp.

228  {
229  // The filenames will have a standard name, with a user-defined suffix
230  std::string filename1 = "basis"+bnme_+".out";
231  std::string filename2 = "restart"+bnme_+".out";
232  std::ofstream basisout;
233  std::ofstream coeffout;
234  basisout.precision(5);
235  coeffout.precision(5);
236  basisout.open(filename1.c_str());
237  coeffout.open(filename2.c_str());
238 
239  // The CV values, PMF projection, PMF, and biased histogram are output for the user
240  basisout << "The information stored in this file is the output of a Basis Function Sampling run" << std::endl;
241  basisout << "CV Values" << std::setw(15*cvs.size()) << "Basis Set Bias" << std::setw(15) << "PMF Estimate" << std::endl;
242  coeffout << "The information stored in this file is for the purpose of restarting simulations with BFS" << std::endl;
243  coeffout << "***COEFFICIENTS***" << std::endl;
244 
245  size_t j = 0;
246  for(Grid<double>::iterator it = b_->begin(); it != b_->end(); ++it, ++j)
247  {
248  for(size_t k = 0; k < cvs.size(); ++k)
249  {
250  // Evaluate the CV values for printing purposes
251  basisout << it.coordinate(k) << std::setw(15);
252  }
253  basisout << -(*it) << std::setw(15) <<
254  -1.0/beta*log(unbias_[j]) <<
255  std::endl;
256  }
257 
258  std::vector<double> coeff = evaluator_.GetCoeff();
259  for (auto& val : coeff) {
260  coeffout << val << std::endl;
261  }
262  coeffout << "***HISTOGRAM***" << std::endl;
263  for (auto& val : unbias_) {
264  coeffout << val << std::endl;
265  }
266 
267  basisout.close();
268  coeffout.close();
269  }
std::vector< double > GetCoeff(void)
Gets the coefficient array.
Definition: Basis.h:394
iterator begin()
Return iterator at first grid point.
Definition: Grid.h:527
GridIterator< double > iterator
Custom iterator over a grid.
Definition: Grid.h:515
iterator end()
Return iterator after last valid grid point.
Definition: Grid.h:540

References SSAGES::Grid< T >::begin().

Here is the call graph for this function:

◆ ProjectBias()

void SSAGES::BFS::ProjectBias ( const CVList cvs,
const double  beta 
)
private

The functions which calculates the updated bias and coefficients and then stores them.

Parameters
cvsList of CVs onto which the bias is projected
betaThermodynamic beta (1/kT)

Definition at line 166 of file BasisFunc.cpp.

167  {
168  double sum = 0.0;
169 
170  // For multiple walkers, the struct is unpacked
171  Grid<unsigned int> histlocal(*h_);
172 
173  // Summed between all walkers
174  MPI_Allreduce(histlocal.data(), h_->data(), h_->size(), MPI_INT, MPI_SUM, world_);
175 
176  h_->syncGrid();
177 
178  // Construct the biased histogram
179  size_t i = 0;
180  for (Grid<unsigned int>::iterator it2 = h_->begin(); it2 != h_->end(); ++it2, ++i)
181  {
182  /* The evaluation of the biased histogram which projects the histogram to the
183  * current bias of CV space.
184  */
185  unbias_[i] += (*it2) * exp(beta * b_->at(it2.coordinates()));
186  }
187 
188  // Reset histogram
189  for (auto &val : *h_) {
190  val = 0;
191  }
192 
193  // Create the log array and send that to the integrator
194  std::vector<double> z (unbias_.size(),0);
195  for (i = 0; i < unbias_.size(); i++)
196  // Make sure the log of the biased histogram is a number
197  z[i] = 1.0/beta*log(unbias_[i] * weight_);
198 
199  //Update the coefficients and determine the difference from the previous iteration
200  sum = evaluator_.UpdateCoeff(z,h_);
202  //Update both the gradient and the bias on the grids
204 
205  if(IsMasterRank(world_)) {
206  // Write coeff at this step, but only one walker
207  std::cout<<"Coefficient difference is: " <<sum<<std::endl;
208  PrintBias(cvs,beta);
209  }
210 
211  // The convergence tolerance and whether the user wants to exit are incorporated here
212  if(sum < tol_)
213  {
214  std::cout<<"System has converged"<<std::endl;
215  if(convergeExit_)
216  {
217  std::cout<<"User has elected to exit. System is now exiting"<<std::endl;
218  exit(EXIT_SUCCESS);
219  }
220  }
221  }
void PrintBias(const CVList &cvs, const double beta)
Prints the current bias to a defined file from the JSON.
Definition: BasisFunc.cpp:227
void UpdateBias(Grid< double > *bias, Grid< std::vector< double >> *grad)
Outputs the basis projection at a specific coordinate.
Definition: Basis.cpp:154
double UpdateCoeff(const std::vector< double > &array, Grid< unsigned int > *hist)
Definition: Basis.cpp:194
static bool IsMasterRank(const MPI_Comm &comm)
Check if current processor is master.
T * data()
Get pointer to the internal data storage vector.
Definition: GridBase.h:338
void syncGrid()
Sync the grid.
Definition: GridBase.h:146

References SSAGES::Grid< T >::begin(), and SSAGES::GridBase< T >::data().

Here is the call graph for this function:

◆ SetBasis()

void SSAGES::BFS::SetBasis ( const std::vector< double > &  coeff,
std::vector< double > &  unbias 
)
inline

Set the values for the basis.

Parameters
coeffList of coefficients.
unbiasList of unbiased values.

This function is used to set starting values at the beginning of a run. For example when continuing from a restart value.

Definition at line 235 of file BasisFunc.h.

236  {
237  evaluator_.SetCoeff(coeff);
238  unbias_ = unbias;
239  }
void SetCoeff(const std::vector< double > &coeff)
Set the coefficient array in the event of a restart run.
Definition: Basis.h:408

References evaluator_, SSAGES::BasisEvaluator::SetCoeff(), and unbias_.

Here is the call graph for this function:

◆ SetIteration()

void SSAGES::BFS::SetIteration ( const int  iter)
inline

Set the current iteration.

Parameters
iterNew value for the current iteration.

This function is used to set the current iteration, for example when continuing from a restart.

Definition at line 222 of file BasisFunc.h.

223  {
224  iteration_ = iter;
225  }

References iteration_.

Member Data Documentation

◆ b_

Grid<double>* SSAGES::BFS::b_
private

Stored bias potential.

The sum of basis functions that adds up to the bias potential of the surface

Definition at line 52 of file BasisFunc.h.

◆ bnme_

std::string SSAGES::BFS::bnme_
private

The option to name both the basis and coefficient files will be given Basis filename

Definition at line 151 of file BasisFunc.h.

◆ evaluator_

BasisEvaluator SSAGES::BFS::evaluator_
private

The basis evaluator class.

A class which holds the basis functions and updates the coefficients at each cycle

Definition at line 76 of file BasisFunc.h.

Referenced by SetBasis().

◆ f_

Grid<std::vector<double> >* SSAGES::BFS::f_
private

Stored gradients.

The gradients corresponding to each point on the grid

Definition at line 58 of file BasisFunc.h.

◆ h_

Grid<unsigned int>* SSAGES::BFS::h_
private

Grid of visited states.

Grid is stored locally.

Definition at line 46 of file BasisFunc.h.

◆ restraint_

std::vector<double> SSAGES::BFS::restraint_
private

Spring constants for restrained system.

The system uses this to determine if the system is to be restrained on the defined interval. The user inputs the spring constants if the system is not periodic.

Definition at line 84 of file BasisFunc.h.

◆ temperature_

double SSAGES::BFS::temperature_
private

Self-defined temperature.

In the case of the MD engine using a poorly defined temperature, the option to throw it into the method is available. If not provided it takes the value from the engine.

Definition at line 118 of file BasisFunc.h.

◆ unbias_

std::vector<double> SSAGES::BFS::unbias_
private

The biased histogram of states.

The biased histogram of states has the form hist_*exp(phi*beta), where phi is the bias potential and beta is the inverse of the temperature. It is defined globally.

Definition at line 66 of file BasisFunc.h.

Referenced by SetBasis().

◆ weight_

double SSAGES::BFS::weight_
private

Weighting for potentially faster sampling.

Weighting can be used to potentially sample faster, however it can cause the simulation to explode. By default this value will be set to 1.0

Definition at line 110 of file BasisFunc.h.


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