SSAGES  0.9.3
Software Suite for Advanced General Ensemble Simulations
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
SSAGES::StringMethod Class Referenceabstract

String base class for FTS, Swarm, and elastic band. More...

#include <StringMethod.h>

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

Public Member Functions

 StringMethod (const MPI_Comm &world, const MPI_Comm &comm, const std::vector< double > &centers, unsigned int maxiterations, const std::vector< double > cvspring, unsigned int frequency)
 Constructor. More...
 
void PreSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override
 Method call prior to simulation initiation. More...
 
virtual void PostIntegration (Snapshot *snapshot, const class CVManager &cvmanager) override=0
 Method call post integration. More...
 
void PostSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override
 Method call post simulation. More...
 
void SetTolerance (std::vector< double > tol)
 Set the tolerance for quitting method. More...
 
void SetSendRecvNeighbors ()
 Communicate neighbor lists over MPI.
 
virtual ~StringMethod ()
 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 StringMethodBuild (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...
 

Protected Member Functions

virtual void StringUpdate ()=0
 Updates the position of the string.
 
double distance (const std::vector< double > &x, const std::vector< double > &y) const
 Helper function for calculating distances. More...
 
void PrintString (const CVList &CV)
 Prints the CV positions to file. More...
 
void GatherNeighbors (std::vector< double > *lcv0, std::vector< double > *ucv0)
 Gather neighbors over MPI. More...
 
void StringReparam (double alpha_star)
 Reparameterize the string. More...
 
void UpdateWorldString (const CVList &cvs)
 Update the world string over MPI. More...
 
bool TolCheck () const
 Check whether tolerance criteria have been met. More...
 
bool CheckEnd (const CVList &CV)
 Check if method reached one of the exit criteria. More...
 

Protected Attributes

std::vector< double > centers_
 CV starting location values.
 
std::vector< double > newcenters_
 CV starting location values.
 
std::vector< std::vector< double > > worldstring_
 The world's strings centers for each CV. More...
 
int mpiid_
 The node this belongs to.
 
std::vector< double > tol_
 Tolerance criteria for determining when to stop string (default 0 if no tolerance criteria)
 
int numnodes_
 Number of nodes on a string.
 
unsigned int maxiterator_
 Maximum cap on number of string method iterations performed.
 
std::vector< double > cvspring_
 Vector of spring constants.
 
unsigned int iterator_
 The local method iterator.
 
unsigned int iteration_
 The global method iteration.
 
std::ofstream stringout_
 Output stream for string data.
 
int sendneigh_
 Neighbor to send info to.
 
int recneigh_
 Neighbor to gain info from.
 
std::vector< std::vector< double > > prev_positions_
 Store positions for starting trajectories.
 
std::vector< std::vector< double > > prev_velocities_
 Store velocities for starting trajectories.
 
std::vector< std::vector< int > > prev_IDs_
 Store atom IDs for starting trajectories.
 
- 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

String base class for FTS, Swarm, and elastic band.

Implementation of a multi-walker finite string method with hard wall voronoi cells and running block averages.

Definition at line 38 of file StringMethod.h.

Constructor & Destructor Documentation

◆ StringMethod()

SSAGES::StringMethod::StringMethod ( const MPI_Comm &  world,
const MPI_Comm &  comm,
const std::vector< double > &  centers,
unsigned int  maxiterations,
const std::vector< double >  cvspring,
unsigned int  frequency 
)
inline

Constructor.

Parameters
worldMPI global communicator.
commMPI local communicator.
centersList of centers.
maxiterationsMaximum number of iterations.
cvspringSpring constants for cvs.
frequencyFrequency with which this method is invoked.

Definition at line 176 of file StringMethod.h.

181  :
182  Method(frequency, world, comm), centers_(centers),
183  maxiterator_(maxiterations),
184  cvspring_(cvspring), iterator_(1), iteration_(0)
185  {
186  newcenters_.resize(centers_.size(), 0);
187  }
Method(unsigned int frequency, const MPI_Comm &world, const MPI_Comm &comm)
Constructor.
Definition: Method.h:61
std::vector< double > cvspring_
Vector of spring constants.
Definition: StringMethod.h:67
std::vector< double > centers_
CV starting location values.
Definition: StringMethod.h:43
unsigned int iteration_
The global method iteration.
Definition: StringMethod.h:73
unsigned int iterator_
The local method iterator.
Definition: StringMethod.h:70
unsigned int maxiterator_
Maximum cap on number of string method iterations performed.
Definition: StringMethod.h:64
std::vector< double > newcenters_
CV starting location values.
Definition: StringMethod.h:46

References centers_, and newcenters_.

Member Function Documentation

◆ Build()

StringMethod * SSAGES::StringMethod::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 220 of file StringMethod.cpp.

224  {
225  ObjectRequirement validator;
226  Value schema;
227  CharReaderBuilder rbuilder;
228  CharReader* reader = rbuilder.newCharReader();
229 
230  StringMethod* m = nullptr;
231 
232  reader->parse(JsonSchema::StringMethod.c_str(),
233  JsonSchema::StringMethod.c_str() + JsonSchema::StringMethod.size(),
234  &schema, nullptr);
235  validator.Parse(schema, path);
236 
237  // Validate inputs.
238  validator.Validate(json, path);
239  if(validator.HasErrors())
240  throw BuildException(validator.GetErrors());
241 
242  unsigned int wid = GetWalkerID(world, comm);
243  std::vector<double> centers;
244  for(auto& s : json["centers"][wid])
245  centers.push_back(s.asDouble());
246 
247  auto maxiterator = json.get("max_iterations", 0).asInt();
248 
249  std::vector<double> ksprings;
250  for(auto& s : json["ksprings"])
251  ksprings.push_back(s.asDouble());
252 
253  auto freq = json.get("frequency", 1).asInt();
254 
255  // Get stringmethod flavor.
256  std::string flavor = json.get("flavor", "none").asString();
257  if(flavor == "ElasticBand")
258  {
259  reader->parse(JsonSchema::ElasticBandMethod.c_str(),
260  JsonSchema::ElasticBandMethod.c_str() + JsonSchema::ElasticBandMethod.size(),
261  &schema, nullptr);
262  validator.Parse(schema, path);
263 
264  // Validate inputs.
265  validator.Validate(json, path);
266  if(validator.HasErrors())
267  throw BuildException(validator.GetErrors());
268 
269  auto eqsteps = json.get("equilibration_steps", 20).asInt();
270  auto evsteps = json.get("evolution_steps", 5).asInt();
271  auto stringspring = json.get("kstring", 10.0).asDouble();
272  auto isteps = json.get("block_iterations", 100).asInt();
273  auto tau = json.get("time_step", 0.1).asDouble();
274 
275  m = new ElasticBand(world, comm, centers,
276  maxiterator, isteps,
277  tau, ksprings, eqsteps,
278  evsteps, stringspring, freq);
279 
280  if(json.isMember("tolerance"))
281  {
282  std::vector<double> tol;
283  for(auto& s : json["tolerance"])
284  tol.push_back(s.asDouble());
285 
286  m->SetTolerance(tol);
287  }
288  }
289  else if(flavor == "FTS")
290  {
291  reader->parse(JsonSchema::FTSMethod.c_str(),
292  JsonSchema::FTSMethod.c_str() + JsonSchema::FTSMethod.size(),
293  &schema, nullptr);
294  validator.Parse(schema, path);
295 
296  // Validate inputs.
297  validator.Validate(json, path);
298  if(validator.HasErrors())
299  throw BuildException(validator.GetErrors());
300 
301  auto isteps = json.get("block_iterations", 2000).asInt();
302  auto tau = json.get("time_step", 0.1).asDouble();
303  auto kappa = json.get("kappa", 0.1).asDouble();
304  auto springiter = json.get("umbrella_iterations",2000).asDouble();
305  m = new FiniteTempString(world, comm, centers,
306  maxiterator, isteps,
307  tau, ksprings, kappa,
308  springiter, freq);
309 
310  if(json.isMember("tolerance"))
311  {
312  std::vector<double> tol;
313  for(auto& s : json["tolerance"])
314  tol.push_back(s.asDouble());
315 
316  m->SetTolerance(tol);
317  }
318  }
319  else if(flavor == "SWARM")
320  {
321  reader->parse(JsonSchema::SwarmMethod.c_str(),
322  JsonSchema::SwarmMethod.c_str() + JsonSchema::SwarmMethod.size(),
323  &schema, nullptr);
324  validator.Parse(schema, path);
325 
326  //Validate input
327  validator.Validate(json, path);
328  if(validator.HasErrors())
329  throw BuildException(validator.GetErrors());
330 
331  auto InitialSteps = json.get("initial_steps", 2500).asInt();
332  auto HarvestLength = json.get("harvest_length", 10).asInt();
333  auto NumberTrajectories = json.get("number_of_trajectories", 250).asInt();
334  auto SwarmLength = json.get("swarm_length", 20).asInt();
335 
336  m = new Swarm(world, comm, centers, maxiterator, ksprings, freq, InitialSteps, HarvestLength, NumberTrajectories, SwarmLength);
337 
338  if(json.isMember("tolerance"))
339  {
340  std::vector<double> tol;
341  for(auto& s : json["tolerance"])
342  tol.push_back(s.asDouble());
343 
344  m->SetTolerance(tol);
345  }
346  }
347 
348  return m;
349  }
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 unsigned int GetWalkerID(const MPI_Comm &world, const MPI_Comm &comm)
Get walker ID number of specified communicator.
Definition: EventListener.h:89
StringMethod(const MPI_Comm &world, const MPI_Comm &comm, const std::vector< double > &centers, unsigned int maxiterations, const std::vector< double > cvspring, unsigned int frequency)
Constructor.
Definition: StringMethod.h:176

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

Here is the call graph for this function:

◆ CheckEnd()

bool SSAGES::StringMethod::CheckEnd ( const CVList CV)
protected

Check if method reached one of the exit criteria.

Parameters
CVlist of CVs.
Returns
True if one of the exit criteria is met.

The string method exits if either the maximum number of iteration has been reached or if all CVs are within the given tolerance thresholds.

Definition at line 128 of file StringMethod.cpp.

129  {
131  {
132  std::cout << "System has reached max string method iterations (" << maxiterator_ << ") as specified in the input file(s)." << std::endl;
133  std::cout << "Exiting now" << std::endl;
134  PrintString(CV); //Ensure that the system prints out if it's about to exit
135  MPI_Abort(world_, EXIT_FAILURE);
136  }
137 
138  int local_tolvalue = TolCheck();
139 
140  MPI_Allreduce(MPI_IN_PLACE, &local_tolvalue, 1, MPI_INT, MPI_LAND, world_);
141 
142  if(local_tolvalue)
143  {
144  std::cout << "System has converged within tolerance criteria. Exiting now" << std::endl;
145  PrintString(CV); //Ensure that the system prints out if it's about to exit
146  MPI_Abort(world_, EXIT_FAILURE);
147  }
148 
149  return true;
150  }
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
bool TolCheck() const
Check whether tolerance criteria have been met.
Definition: StringMethod.h:140
void PrintString(const CVList &CV)
Prints the CV positions to file.

Referenced by SSAGES::ElasticBand::PostIntegration(), SSAGES::FiniteTempString::PostIntegration(), and SSAGES::Swarm::PostIntegration().

Here is the caller graph for this function:

◆ distance()

double SSAGES::StringMethod::distance ( const std::vector< double > &  x,
const std::vector< double > &  y 
) const
inlineprotected

Helper function for calculating distances.

Parameters
xList of coordinates.
yList of coordinates.
Returns
Sum of distances between the x-values and y-values.

Definition at line 102 of file StringMethod.h.

103  {
104  double distance = 0;
105  for (size_t i = 0; i < x.size(); i++)
106  distance += pow((x[i] - y[i]),2);
107 
108  return sqrt(distance);
109  }
double distance(const std::vector< double > &x, const std::vector< double > &y) const
Helper function for calculating distances.
Definition: StringMethod.h:102

Referenced by SSAGES::FiniteTempString::StringUpdate(), and SSAGES::Swarm::StringUpdate().

Here is the caller graph for this function:

◆ GatherNeighbors()

void SSAGES::StringMethod::GatherNeighbors ( std::vector< double > *  lcv0,
std::vector< double > *  ucv0 
)
protected

Gather neighbors over MPI.

Parameters
lcv0Pointer to store array of lower CV values.
ucv0Pointer to store array of upper CV values.

Definition at line 53 of file StringMethod.cpp.

54  {
55  MPI_Status status;
56 
57  if(IsMasterRank(comm_))
58  {
59  MPI_Sendrecv(&centers_[0], centers_.size(), MPI_DOUBLE, sendneigh_, 1234,
60  &(*lcv0)[0], centers_.size(), MPI_DOUBLE, recneigh_, 1234,
61  world_, &status);
62 
63  MPI_Sendrecv(&centers_[0], centers_.size(), MPI_DOUBLE, recneigh_, 4321,
64  &(*ucv0)[0], centers_.size(), MPI_DOUBLE, sendneigh_, 4321,
65  world_, &status);
66  }
67 
68  MPI_Bcast(&(*lcv0)[0],centers_.size(),MPI_DOUBLE,0,comm_);
69  MPI_Bcast(&(*ucv0)[0],centers_.size(),MPI_DOUBLE,0,comm_);
70  }
static bool IsMasterRank(const MPI_Comm &comm)
Check if current processor is master.
mxx::comm comm_
Local MPI communicator.
Definition: Method.h:47
int recneigh_
Neighbor to gain info from.
Definition: StringMethod.h:82
int sendneigh_
Neighbor to send info to.
Definition: StringMethod.h:79

Referenced by SSAGES::ElasticBand::StringUpdate(), SSAGES::FiniteTempString::StringUpdate(), and SSAGES::Swarm::StringUpdate().

Here is the caller graph for this function:

◆ PostIntegration()

virtual void SSAGES::StringMethod::PostIntegration ( Snapshot snapshot,
const class CVManager cvmanager 
)
overridepure virtual

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.

Implemented in SSAGES::Swarm, SSAGES::FiniteTempString, and SSAGES::ElasticBand.

◆ PostSimulation()

void SSAGES::StringMethod::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 169 of file StringMethod.cpp.

170  {
171  stringout_.close();
172  }
std::ofstream stringout_
Output stream for string data.
Definition: StringMethod.h:76

◆ PreSimulation()

void SSAGES::StringMethod::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 152 of file StringMethod.cpp.

153  {
154  char file[1024];
155  mpiid_ = snapshot->GetWalkerID();
156  sprintf(file, "node-%04d.log",mpiid_);
157  stringout_.open(file);
158 
159  auto cvs = cvmanager.GetCVs(cvmask_);
161  worldstring_.resize(numnodes_);
162  for(auto& w : worldstring_)
163  w.resize(centers_.size());
164 
165  UpdateWorldString(cvs);
166  PrintString(cvs);
167  }
std::vector< unsigned int > cvmask_
Mask which identifies which CVs to act on.
Definition: Method.h:50
std::vector< std::vector< double > > worldstring_
The world's strings centers for each CV.
Definition: StringMethod.h:52
int mpiid_
The node this belongs to.
Definition: StringMethod.h:55
void UpdateWorldString(const CVList &cvs)
Update the world string over MPI.
void SetSendRecvNeighbors()
Communicate neighbor lists over MPI.
int numnodes_
Number of nodes on a string.
Definition: StringMethod.h:61

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

Here is the call graph for this function:

◆ PrintString()

void SSAGES::StringMethod::PrintString ( const CVList CV)
protected

Prints the CV positions to file.

Parameters
CVCollective Variables to be printed

Definition at line 38 of file StringMethod.cpp.

39  {
40  if(IsMasterRank(comm_))
41  {
42  //Write node, iteration, centers of the string and current CV value to output file
43  stringout_.precision(8);
44  stringout_ << mpiid_ << " " << iteration_ << " ";
45 
46  for(size_t i = 0; i < centers_.size(); i++)
47  stringout_ << worldstring_[mpiid_][i] << " " << CV[i]->GetValue() << " ";
48 
49  stringout_ << std::endl;
50  }
51  }

Referenced by SSAGES::ElasticBand::PostIntegration(), SSAGES::FiniteTempString::PostIntegration(), and SSAGES::Swarm::PostIntegration().

Here is the caller graph for this function:

◆ SetTolerance()

void SSAGES::StringMethod::SetTolerance ( std::vector< double >  tol)
inline

Set the tolerance for quitting method.

Parameters
tolList of tolerances for each CV.

Set the tolerance levels until which the method should run. The method quit, when the tolerance level is reached for all CVs.

Definition at line 205 of file StringMethod.h.

206  {
207  for(auto& t : tol)
208  tol_.push_back(t);
209  }
std::vector< double > tol_
Tolerance criteria for determining when to stop string (default 0 if no tolerance criteria)
Definition: StringMethod.h:58

References tol_.

Referenced by Build().

Here is the caller graph for this function:

◆ StringReparam()

void SSAGES::StringMethod::StringReparam ( double  alpha_star)
protected

Reparameterize the string.

Parameters
alpha_starFactor for reparametrization.

Definition at line 72 of file StringMethod.cpp.

73  {
74  std::vector<double> alpha_star_vector(numnodes_,0.0);
75 
76  //Reparameterization
77  //Alpha star is the uneven mesh, approximated as linear distance between points
78  if(IsMasterRank(comm_))
79  alpha_star_vector[mpiid_] = mpiid_ == 0 ? 0 : alpha_star;
80 
81  //Gather each alpha_star into a vector
82  MPI_Allreduce(MPI_IN_PLACE, &alpha_star_vector[0], numnodes_, MPI_DOUBLE, MPI_SUM, world_);
83 
84  for(size_t i = 1; i < alpha_star_vector.size(); i++)
85  alpha_star_vector[i] += alpha_star_vector[i-1];
86 
87  for(size_t i = 1; i < alpha_star_vector.size(); i++)
88  alpha_star_vector[i] /= alpha_star_vector[numnodes_ - 1];
89 
90  tk::spline spl; //Cubic spline interpolation
91 
92  for(size_t i = 0; i < centers_.size(); i++)
93  {
94  std::vector<double> cvs_new(numnodes_, 0.0);
95 
96  if(IsMasterRank(comm_))
97  cvs_new[mpiid_] = centers_[i];
98 
99  MPI_Allreduce(MPI_IN_PLACE, &cvs_new[0], numnodes_, MPI_DOUBLE, MPI_SUM, world_);
100 
101  spl.set_points(alpha_star_vector, cvs_new);
102  centers_[i] = spl(mpiid_/(numnodes_ - 1.0));
103  }
104  }

Referenced by SSAGES::FiniteTempString::StringUpdate(), and SSAGES::Swarm::StringUpdate().

Here is the caller graph for this function:

◆ TolCheck()

bool SSAGES::StringMethod::TolCheck ( ) const
inlineprotected

Check whether tolerance criteria have been met.

Returns
Boolean for whether tolerance criteria have been met.

Definition at line 140 of file StringMethod.h.

141  {
142  if(tol_.size() == 0)
143  return false;
144 
145  for(size_t i = 0; i < centers_.size(); i++)
146  {
147  if(fabs(centers_[i] - worldstring_[mpiid_][i]) > tol_[i])
148  {
149  return false;
150  }
151  }
152 
153  return true;
154  }

References centers_, mpiid_, tol_, and worldstring_.

◆ UpdateWorldString()

void SSAGES::StringMethod::UpdateWorldString ( const CVList cvs)
protected

Update the world string over MPI.

Parameters
cvsList of CVs.

Definition at line 106 of file StringMethod.cpp.

107  {
108  for(size_t i = 0; i < centers_.size(); i++)
109  {
110  std::vector<double> cvs_new(numnodes_, 0.0);
111 
112  if(IsMasterRank(comm_))
113  {
114  cvs_new[mpiid_] = centers_[i];
115  }
116 
117  MPI_Allreduce(MPI_IN_PLACE, &cvs_new[0], numnodes_, MPI_DOUBLE, MPI_SUM, world_);
118 
119  for(int j = 0; j < numnodes_; j++)
120  {
121  worldstring_[j][i] = cvs_new[j];
122  //Represent worldstring in periodic space
123  worldstring_[j][i] = cvs[i]->GetPeriodicValue(worldstring_[j][i]);
124  }
125  }
126  }

Referenced by SSAGES::ElasticBand::PostIntegration(), SSAGES::FiniteTempString::PostIntegration(), and SSAGES::Swarm::PostIntegration().

Here is the caller graph for this function:

Member Data Documentation

◆ worldstring_

std::vector<std::vector<double> > SSAGES::StringMethod::worldstring_
protected

The world's strings centers for each CV.

worldstring_[node#][cv#]

Definition at line 52 of file StringMethod.h.

Referenced by SSAGES::FiniteTempString::InCell(), and TolCheck().


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