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

Finite Temperature Spring Method. More...

#include <FiniteTempString.h>

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

Public Member Functions

 FiniteTempString (const MPI_Comm &world, const MPI_Comm &comm, const std::vector< double > &centers, unsigned int maxiterations, unsigned int blockiterations, double tau, const std::vector< double > cvspring, double kappa, unsigned int springiter, unsigned int frequency)
 Constructor. More...
 
void PostIntegration (Snapshot *snapshot, const class CVManager &cvmanager) override
 Method call post integration. More...
 
 ~FiniteTempString ()
 Destructor.
 
- Public Member Functions inherited from 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)
 Constructor. More...
 
void PreSimulation (Snapshot *snapshot, const class CVManager &cvmanager) override
 Method call prior to simulation initiation. 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.
 

Private Member Functions

bool InCell (const CVList &cvs) const
 Checks if CV is in Voronoi cell. More...
 
void StringUpdate () override
 Updates the string according to the FTS method.
 

Private Attributes

double kappa_
 String modification parameter.
 
unsigned int blockiterations_
 Number of steps to block average the CV's postions over.
 
double tau_
 Time step of string change.
 
unsigned int min_num_umbrella_steps_
 Minimum number of steps to apply umbrella sampling.
 
int run_umbrella_
 Flag to run umbrella or not during post-integration.
 
unsigned int umbrella_iter_
 Iterator that keeps track of umbrella iterations.
 
std::vector< double > prev_CVs_
 Stores the last positions of the CVs.
 
bool reset_for_umbrella
 Flag for whether a system was to run umbrella sampling before checking against other systems.
 

Additional Inherited Members

- Static Public Member Functions inherited from SSAGES::StringMethod
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 inherited from SSAGES::StringMethod
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 inherited from SSAGES::StringMethod
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

Finite Temperature Spring Method.

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

Definition at line 35 of file FiniteTempString.h.

Constructor & Destructor Documentation

◆ FiniteTempString()

SSAGES::FiniteTempString::FiniteTempString ( const MPI_Comm &  world,
const MPI_Comm &  comm,
const std::vector< double > &  centers,
unsigned int  maxiterations,
unsigned int  blockiterations,
double  tau,
const std::vector< double >  cvspring,
double  kappa,
unsigned int  springiter,
unsigned int  frequency 
)
inline

Constructor.

Parameters
worldMPI global communicator.
commMPI local communicator.
centersList of centers.
maxiterationsMaximum number of iterations.
blockiterationsNumber of iterations per block averaging.
tauValue of tau (default: 0.1).
cvspringSpring constants for cvs.
kappaValue of kappa (default: 0.1).
springiterMinimum number of umbrella steps.
frequencyFrequency with which this method is invoked.

Constructs an instance of Finite String method.

Store positions for starting trajectories

Store velocities for starting trajectories

Definition at line 90 of file FiniteTempString.h.

99  :
100  StringMethod(world, comm, centers, maxiterations, cvspring, frequency),
101  kappa_(kappa), blockiterations_(blockiterations), tau_(tau),
102  min_num_umbrella_steps_(springiter), run_umbrella_(true),
104  {
106  prev_positions_.resize(1);
107 
109  prev_velocities_.resize(1);
110 
111  prev_IDs_.resize(1);
112 
113  }
unsigned int min_num_umbrella_steps_
Minimum number of steps to apply umbrella sampling.
unsigned int blockiterations_
Number of steps to block average the CV's postions over.
double tau_
Time step of string change.
bool reset_for_umbrella
Flag for whether a system was to run umbrella sampling before checking against other systems.
int run_umbrella_
Flag to run umbrella or not during post-integration.
double kappa_
String modification parameter.
unsigned int umbrella_iter_
Iterator that keeps track of umbrella iterations.
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
std::vector< std::vector< double > > prev_velocities_
Store velocities for starting trajectories.
Definition: StringMethod.h:88
std::vector< std::vector< double > > prev_positions_
Store positions for starting trajectories.
Definition: StringMethod.h:85
std::vector< std::vector< int > > prev_IDs_
Store atom IDs for starting trajectories.
Definition: StringMethod.h:91

References SSAGES::StringMethod::prev_IDs_, SSAGES::StringMethod::prev_positions_, and SSAGES::StringMethod::prev_velocities_.

Member Function Documentation

◆ InCell()

bool SSAGES::FiniteTempString::InCell ( const CVList cvs) const
private

Checks if CV is in Voronoi cell.

Parameters
cvsList of CVs to check
Returns
Boolean if all CVs are in the Voronoi cell

Definition at line 32 of file FiniteTempString.cpp.

33  {
34  std::vector<double> dists (numnodes_, 0);
35 
36  // Record the difference between all cvs and all nodes
37  for (int i = 0; i < numnodes_; i++)
38  for(size_t j = 0; j < cvs.size(); j++)
39  dists[i]+= pow(cvs[j]->GetDifference(worldstring_[i][j]),2);
40 
41  if(std::min_element(dists.begin(), dists.end()) - dists.begin() == mpiid_)
42  return true;
43 
44  return false;
45  }
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
int numnodes_
Number of nodes on a string.
Definition: StringMethod.h:61

References SSAGES::StringMethod::mpiid_, SSAGES::StringMethod::numnodes_, and SSAGES::StringMethod::worldstring_.

Referenced by PostIntegration().

Here is the caller graph for this function:

◆ PostIntegration()

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

Method call post integration.

Post-integration hook.

Parameters
snapshotPointer to the simulation snapshot.
cvmanagerCollective variable manager.

This function will be called after each integration step.

Implements SSAGES::StringMethod.

Definition at line 48 of file FiniteTempString.cpp.

49  {
50  auto cvs = cvmanager.GetCVs(cvmask_);
51  auto& forces = snapshot->GetForces();
52  bool insidecell;
53 
55  {
56  //If the system was not going to run the umbrella anymore, reset its position
57  for(auto& force : forces)
58  {
59  force.setZero();
60  }
61 
62  auto& Pos = snapshot->GetPositions();
63 
64  for(size_t i = 0; i < prev_IDs_[0].size(); i++)
65  {
66  auto localindex = snapshot->GetLocalIndex(prev_IDs_[0][i]);
67  if(localindex!= -1)
68  {
69  Pos[localindex][0] = prev_positions_[0][i*3];
70  Pos[localindex][1] = prev_positions_[0][i*3 + 1];
71  Pos[localindex][2] = prev_positions_[0][i*3 + 2];
72  }
73  }
74  }
75  for(auto& cv : cvs)
76  {
77  //Trigger a rebuild of the CVs since we reset the positions
78  cv->Evaluate(*snapshot);
79  }
80  insidecell = InCell(cvs);
81 
82  MPI_Allreduce(MPI_IN_PLACE, &run_umbrella_, 1, MPI_INT, MPI_LOR, world_); //Compares all run_umbrellas on all processors - if any is true, all are set to true
83 
84  if(run_umbrella_)
85  {
87  {
88  //After reaching the minimum umbrella sampling checkpoint, evaluate whether to keep going (per individual image)
89  if(insidecell)
90  {
91  //Stop umbrella sampling if inside cell
92  run_umbrella_ = false;
93  reset_for_umbrella = true;
94 
95  //This node is done initializing; so store this snapshot
96  prev_positions_[0] = snapshot->SerializePositions();
97  prev_IDs_[0] = snapshot->SerializeIDs();
98  }
99  umbrella_iter_ = 1;
100  }
101  else
102  {
103  if(!reset_for_umbrella)
104  {
105  //If not at the minimum step checkpoint AND not a system which was already done, add restraining force
106  for(size_t i = 0; i < cvs.size(); i++)
107  {
108  // Get current cv and gradient
109  auto& cv = cvs[i];
110  auto& grad = cv->GetGradient();
111 
112  // Compute dV/dCV
113  auto D = cvspring_[i]*(cv->GetDifference(centers_[i]));
114 
115  // Update forces
116  for(size_t j = 0; j < forces.size(); j++)
117  forces[j] -= D*grad[j];
118  }
119  umbrella_iter_++; //Progress toward checkpoint
120  }
121  else
122  {
123  //End umbrella sampling of this node if it was already initialized
124  run_umbrella_ = false;
125  }
126  }
127  return; //Don't do any sampling after restraining
128  }
129  else
130  {
131  //If run umbrella was false (only possible if every node was done) then set reset to false for next string method iteration (only proceed past umbrella sampling if no images require further umbrella sampling)
132  reset_for_umbrella = false;
133  }
134 
135  if(!insidecell)
136  {
137  //If last step took the system outside the Voronoi cell, reset the position and zero the force to reset the trajectory
138  for(auto& force : forces)
139  force.setZero();
140 
141  auto& Pos = snapshot->GetPositions();
142 
143  for(size_t i = 0; i < prev_IDs_[0].size(); i++)
144  {
145  auto localindex = snapshot->GetLocalIndex(prev_IDs_[0][i]);
146  if(localindex!= -1)
147  {
148  Pos[localindex][0] = prev_positions_[0][i*3];
149  Pos[localindex][1] = prev_positions_[0][i*3 + 1];
150  Pos[localindex][2] = prev_positions_[0][i*3 + 2];
151  }
152  }
153 
154  // Calculate running averages for each CV at each node based on previous CV
155  for(size_t i = 0; i < newcenters_.size(); i++)
156  {
159  }
160  }
161  else
162  {
163  // Calculate running averages for each CV at each node
164  for(size_t i = 0; i < newcenters_.size(); i++)
165  {
166  newcenters_[i] = newcenters_[i] * (iteration_ * blockiterations_ + iterator_ - 1) + cvs[i]->GetMinimumImage(centers_[i]);
168  }
169 
170  //Store info about this step in prev_vectors
171  prev_CVs_.clear();
172  for(size_t i = 0; i < centers_.size(); i++)
173  {
174  prev_CVs_.push_back(cvs[i]->GetMinimumImage(centers_[i]));
175  }
176  prev_positions_[0] = snapshot->SerializePositions();
177  prev_IDs_[0] = snapshot->SerializeIDs();
178  }
179 
180  // Update the string, every blockiterations_ string method iterations
181  if(iterator_ % blockiterations_ == 0)
182  {
183  MPI_Barrier(world_);
184  StringUpdate();
185  CheckEnd(cvs);
186  MPI_Barrier(world_);
187  UpdateWorldString(cvs);
188  PrintString(cvs);
189 
190  iterator_ = 0; //Number of steps within current method iteration
191  iteration_++; //Increment string method iteration
192 
193  if(!InCell(cvs))
194  {
195  run_umbrella_ = true;
196  reset_for_umbrella = false;
197  }
198  MPI_Allreduce(MPI_IN_PLACE, &run_umbrella_, 1, MPI_INT, MPI_LOR, world_); //Compares all run_umbrellas on all processors - if any is true, all are set to true
199  }
200 
201  iterator_++; //Iterate number of steps within current method iteration
202  }
std::vector< double > prev_CVs_
Stores the last positions of the CVs.
bool InCell(const CVList &cvs) const
Checks if CV is in Voronoi cell.
void StringUpdate() override
Updates the string according to the FTS method.
std::vector< unsigned int > cvmask_
Mask which identifies which CVs to act on.
Definition: Method.h:50
mxx::comm world_
Global MPI communicator.
Definition: Method.h:46
std::vector< double > cvspring_
Vector of spring constants.
Definition: StringMethod.h:67
void UpdateWorldString(const CVList &cvs)
Update the world string over MPI.
std::vector< double > centers_
CV starting location values.
Definition: StringMethod.h:43
bool CheckEnd(const CVList &CV)
Check if method reached one of the exit criteria.
unsigned int iteration_
The global method iteration.
Definition: StringMethod.h:73
void PrintString(const CVList &CV)
Prints the CV positions to file.
unsigned int iterator_
The local method iterator.
Definition: StringMethod.h:70
std::vector< double > newcenters_
CV starting location values.
Definition: StringMethod.h:46

References blockiterations_, SSAGES::StringMethod::centers_, SSAGES::StringMethod::CheckEnd(), SSAGES::Method::cvmask_, SSAGES::StringMethod::cvspring_, SSAGES::CVManager::GetCVs(), SSAGES::Snapshot::GetForces(), SSAGES::Snapshot::GetLocalIndex(), SSAGES::Snapshot::GetPositions(), InCell(), SSAGES::StringMethod::iteration_, SSAGES::StringMethod::iterator_, min_num_umbrella_steps_, SSAGES::StringMethod::newcenters_, prev_CVs_, SSAGES::StringMethod::prev_IDs_, SSAGES::StringMethod::prev_positions_, SSAGES::StringMethod::PrintString(), reset_for_umbrella, run_umbrella_, SSAGES::Snapshot::SerializeIDs(), SSAGES::Snapshot::SerializePositions(), StringUpdate(), umbrella_iter_, SSAGES::StringMethod::UpdateWorldString(), and SSAGES::Method::world_.

Here is the call graph for this function:

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