SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
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
 Pre-simulation hook.
 
virtual void PostIntegration (Snapshot *snapshot, const class CVManager &cvmanager) override=0
 Post-integration hook.
 
void PostSimulation (Snapshot *, const class 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 (uint frequency, const MPI_Comm &world, const MPI_Comm &comm)
 Constructor. More...
 
void SetCVMask (const std::vector< uint > &mask)
 Sets the collective variable mask.
 
virtual ~Method ()
 Destructor.
 
- Public Member Functions inherited from SSAGES::EventListener
 EventListener (uint frequency)
 Constructor. More...
 
uint 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...
 

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.
 
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 has been met.
 
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.
 
uint 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< uint > 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

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 170 of file StringMethod.h.

References centers_, and newcenters_.

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

Member Function Documentation

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 215 of file StringMethod.cpp.

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

219  {
220  ObjectRequirement validator;
221  Value schema;
222  Reader reader;
223 
224  StringMethod* m = nullptr;
225 
226  reader.parse(JsonSchema::StringMethod, schema);
227  validator.Parse(schema, path);
228 
229  // Validate inputs.
230  validator.Validate(json, path);
231  if(validator.HasErrors())
232  throw BuildException(validator.GetErrors());
233 
234  unsigned int wid = mxx::comm(world).rank()/mxx::comm(comm).size();
235  std::vector<double> centers;
236  for(auto& s : json["centers"][wid])
237  centers.push_back(s.asDouble());
238 
239  auto maxiterator = json.get("max_iterations", 0).asInt();
240 
241  std::vector<double> ksprings;
242  for(auto& s : json["ksprings"])
243  ksprings.push_back(s.asDouble());
244 
245  auto freq = json.get("frequency", 1).asInt();
246 
247  // Get stringmethod flavor.
248  std::string flavor = json.get("flavor", "none").asString();
249  if(flavor == "ElasticBand")
250  {
251  reader.parse(JsonSchema::ElasticBandMethod, schema);
252  validator.Parse(schema, path);
253 
254  // Validate inputs.
255  validator.Validate(json, path);
256  if(validator.HasErrors())
257  throw BuildException(validator.GetErrors());
258 
259  auto eqsteps = json.get("equilibration_steps", 20).asInt();
260  auto evsteps = json.get("evolution_steps", 5).asInt();
261  auto stringspring = json.get("kstring", 10.0).asDouble();
262  auto isteps = json.get("block_iterations", 100).asInt();
263  auto tau = json.get("time_step", 0.1).asDouble();
264 
265  m = new ElasticBand(world, comm, centers,
266  maxiterator, isteps,
267  tau, ksprings, eqsteps,
268  evsteps, stringspring, freq);
269 
270  if(json.isMember("tolerance"))
271  {
272  std::vector<double> tol;
273  for(auto& s : json["tolerance"])
274  tol.push_back(s.asDouble());
275 
276  m->SetTolerance(tol);
277  }
278  }
279  else if(flavor == "FTS")
280  {
281  reader.parse(JsonSchema::FTSMethod, schema);
282  validator.Parse(schema, path);
283 
284  // Validate inputs.
285  validator.Validate(json, path);
286  if(validator.HasErrors())
287  throw BuildException(validator.GetErrors());
288 
289  auto isteps = json.get("block_iterations", 2000).asInt();
290  auto tau = json.get("time_step", 0.1).asDouble();
291  auto kappa = json.get("kappa", 0.1).asDouble();
292  auto springiter = json.get("umbrella_iterations",2000).asDouble();
293  m = new FiniteTempString(world, comm, centers,
294  maxiterator, isteps,
295  tau, ksprings, kappa,
296  springiter, freq);
297 
298  if(json.isMember("tolerance"))
299  {
300  std::vector<double> tol;
301  for(auto& s : json["tolerance"])
302  tol.push_back(s.asDouble());
303 
304  m->SetTolerance(tol);
305  }
306  }
307  else if(flavor == "SWARM")
308  {
309  reader.parse(JsonSchema::SwarmMethod, schema);
310  validator.Parse(schema, path);
311 
312  //Validate input
313  validator.Validate(json, path);
314  if(validator.HasErrors())
315  throw BuildException(validator.GetErrors());
316 
317  auto InitialSteps = json.get("initial_steps", 2500).asInt();
318  auto HarvestLength = json.get("harvest_length", 10).asInt();
319  auto NumberTrajectories = json.get("number_of_trajectories", 250).asInt();
320  auto SwarmLength = json.get("swarm_length", 20).asInt();
321 
322  m = new Swarm(world, comm, centers, maxiterator, ksprings, freq, InitialSteps, HarvestLength, NumberTrajectories, SwarmLength);
323 
324  if(json.isMember("tolerance"))
325  {
326  std::vector<double> tol;
327  for(auto& s : json["tolerance"])
328  tol.push_back(s.asDouble());
329 
330  m->SetTolerance(tol);
331  }
332  }
333 
334  return m;
335  }
bool HasErrors()
Check if errors have occured.
Definition: Requirement.h:86
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:170
virtual void Parse(Value json, const std::string &path) override
Parse JSON value to generate Requirement(s).
std::vector< std::string > GetErrors()
Get list of error messages.
Definition: Requirement.h:92
Requirements on an object.
virtual void Validate(const Value &json, const std::string &path) override
Validate JSON value.

Here is the call graph for this function:

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.

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

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  }
bool TolCheck() const
Check whether tolerance criteria has been met.
Definition: StringMethod.h:134
uint iteration_
The global method iteration.
Definition: StringMethod.h:73
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
unsigned int maxiterator_
Maximum cap on number of string method iterations performed.
Definition: StringMethod.h:64
void PrintString(const CVList &CV)
Prints the CV positions to file.

Here is the caller graph for this function:

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.

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

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

Here is the caller graph for this function:

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.

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

54  {
55  MPI_Status status;
56 
57  if(comm_.rank() == 0)
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  }
int sendneigh_
Neighbor to send info to.
Definition: StringMethod.h:79
std::vector< double > centers_
CV starting location values.
Definition: StringMethod.h:43
mxx::comm comm_
Local MPI communicator.
Definition: Method.h:47
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
int recneigh_
Neighbor to gain info from.
Definition: StringMethod.h:82

Here is the caller graph for this function:

void SSAGES::StringMethod::PostSimulation ( Snapshot snapshot,
const class CVManager cvmanager 
)
inlineoverridevirtual

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 190 of file StringMethod.h.

References stringout_.

191  {
192  stringout_.close();
193  }
std::ofstream stringout_
Output stream for string data.
Definition: StringMethod.h:76
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 202 of file StringMethod.h.

References tol_.

Referenced by Build().

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

Here is the caller graph for this function:

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

Reparameterize the string.

Parameters
alpha_starFactor for reparametrization.

Definition at line 72 of file StringMethod.cpp.

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

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(comm_.rank()==0)
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(comm_.rank() == 0)
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  }
std::vector< double > centers_
CV starting location values.
Definition: StringMethod.h:43
mxx::comm comm_
Local MPI communicator.
Definition: Method.h:47
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
int numnodes_
Number of nodes on a string.
Definition: StringMethod.h:61
int mpiid_
The node this belongs to.
Definition: StringMethod.h:55

Here is the caller graph for this function:

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.

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

107  {
108  for(size_t i = 0; i < centers_.size(); i++)
109  {
110  std::vector<double> cvs_new(numnodes_, 0.0);
111 
112  if(comm_.rank() == 0)
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  }
std::vector< double > centers_
CV starting location values.
Definition: StringMethod.h:43
mxx::comm comm_
Local MPI communicator.
Definition: Method.h:47
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
int numnodes_
Number of nodes on a string.
Definition: StringMethod.h:61
int mpiid_
The node this belongs to.
Definition: StringMethod.h:55
std::vector< std::vector< double > > worldstring_
The world's strings centers for each CV.
Definition: StringMethod.h:52

Here is the caller graph for this function:

Member Data Documentation

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: