SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
SSAGES::PairwiseCV Class Reference

Generalized collective variable based on pairwise properties of atoms. More...

#include <PairwiseCV.h>

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

Public Member Functions

 PairwiseCV (const Label &group1, const Label &group2, PairwiseKernel *pk)
 Constructor. More...
 
void Initialize (const Snapshot &snapshot) override
 Initialize necessary variables. More...
 
void Evaluate (const Snapshot &snapshot) override
 Evaluate the CV. More...
 
- Public Member Functions inherited from SSAGES::CollectiveVariable
 CollectiveVariable ()
 Constructor.
 
double GetValue () const
 Get current value of the CV. More...
 
virtual double GetMinimumImage (double) const
 Returns the minimum image of a CV based on the input location. More...
 
virtual double GetPeriodicValue (double location) const
 Apply periodic boundaries to a given value. More...
 
const std::vector< Vector3 > & GetGradient () const
 Get current gradient of the CV. More...
 
const Matrix3GetBoxGradient () const
 Get gradient contribution to box.
 
const std::array< double, 2 > & GetBoundaries ()
 Get CV boundaries. More...
 
virtual double GetDifference (double location) const
 

Static Public Member Functions

static PairwiseCVBuild (const Json::Value &json, const std::string &path)
 
- Static Public Member Functions inherited from SSAGES::CollectiveVariable
static CollectiveVariableBuildCV (const Json::Value &json, const std::string &path)
 Set up collective variable. More...
 

Private Attributes

Label group1_
 IDs of the first group of atoms.
 
Label group2_
 IDs of the second group of atoms.
 
PairwiseKernelpk_
 Pairwise kernel function used for CV.
 

Additional Inherited Members

- Protected Attributes inherited from SSAGES::CollectiveVariable
std::vector< Vector3grad_
 Gradient vector dCv/dxi.
 
Matrix3 boxgrad_
 Gradient w.r.t box vectors dCv/dHij.
 
double val_
 Current value of CV.
 
std::array< double, 2 > bounds_
 Bounds on CV.
 

Detailed Description

Generalized collective variable based on pairwise properties of atoms.

Collective variable on pairwise properties between two groups of atoms. To ensure generality of usage, there are various pairwise kernel functions from which to choose.

Definition at line 40 of file PairwiseCV.h.

Constructor & Destructor Documentation

SSAGES::PairwiseCV::PairwiseCV ( const Label group1,
const Label group2,
PairwiseKernel pk 
)
inline

Constructor.

Parameters
group1IDs of the first group of atoms.
group2IDs of the second group of atoms.

Construct a PairwiseCV.

Definition at line 56 of file PairwiseCV.h.

56  :
57  group1_(group1), group2_(group2), pk_(pk)
58  {
59  }
Label group1_
IDs of the first group of atoms.
Definition: PairwiseCV.h:43
PairwiseKernel * pk_
Pairwise kernel function used for CV.
Definition: PairwiseCV.h:45
Label group2_
IDs of the second group of atoms.
Definition: PairwiseCV.h:44

Member Function Documentation

void SSAGES::PairwiseCV::Evaluate ( const Snapshot snapshot)
inlineoverride

Evaluate the CV.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 111 of file PairwiseCV.h.

References SSAGES::Snapshot::ApplyMinimumImage(), SSAGES::CollectiveVariable::boxgrad_, SSAGES::PairwiseKernel::Evaluate(), SSAGES::Snapshot::GetAtomIDs(), SSAGES::Snapshot::GetCommunicator(), SSAGES::Snapshot::GetLocalIndex(), SSAGES::Snapshot::GetLocalIndices(), SSAGES::Snapshot::GetNumAtoms(), SSAGES::Snapshot::GetPositions(), SSAGES::CollectiveVariable::grad_, group1_, group2_, pk_, and SSAGES::CollectiveVariable::val_.

112  {
113  // Get local atom indices.
114  std::vector<int> idx1, idx2;
115  snapshot.GetLocalIndices(group1_, &idx1);
116  snapshot.GetLocalIndices(group2_, &idx2);
117 
118  // Get data from snapshot.
119  auto n = snapshot.GetNumAtoms();
120  auto& atomids = snapshot.GetAtomIDs();
121  auto& positions = snapshot.GetPositions();
122  auto& comm = snapshot.GetCommunicator();
123 
124  // Initialize gradient.
125  std::fill(grad_.begin(), grad_.end(), Vector3{0,0,0});
126  grad_.resize(n, Vector3{0,0,0});
127  boxgrad_ = Matrix3::Zero();
128 
129  // Fill eigen matrix with coordinate data.
130  std::vector<double> pos1(3*idx1.size()), pos2(3*idx2.size());
131  std::vector<int> id1(idx1.size()), id2(idx2.size());
132  for(size_t i = 0; i < idx1.size(); ++i)
133  {
134  pos1[3*i+0] = positions[idx1[i]][0];
135  pos1[3*i+1] = positions[idx1[i]][1];
136  pos1[3*i+2] = positions[idx1[i]][2];
137 
138  id1[i] = atomids[idx1[i]];
139  }
140 
141  for(size_t i = 0; i < idx2.size(); ++i)
142  {
143  pos2[3*i+0] = positions[idx2[i]][0];
144  pos2[3*i+1] = positions[idx2[i]][1];
145  pos2[3*i+2] = positions[idx2[i]][2];
146 
147  id2[i] = atomids[idx2[i]];
148  }
149 
150  // Gather across all procs.
151  pos1 = mxx::allgatherv(pos1.data(), pos1.size(), comm);
152  pos2 = mxx::allgatherv(pos2.data(), pos2.size(), comm);
153  id1 = mxx::allgatherv(id1.data(), id1.size(), comm);
154  id2 = mxx::allgatherv(id2.data(), id2.size(), comm);
155 
156  // Create Eigen map for ease of use.
157  using Map = Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, 3, Eigen::RowMajor>>;
158  Map gpos1(pos1.data(), group1_.size(), 3), gpos2(pos2.data(), group2_.size(), 3);
159 
160  val_ = 0;
161  double df = 0.;
162  // Compute gradient and value.
163  for(size_t i = 0; i < group1_.size(); ++i)
164  {
165  auto t1 = id1[i];
166  const auto& pi = gpos1.row(i);
167  for(size_t j = 0; j < group2_.size(); ++j)
168  {
169  auto t2 = id2[j];
170  const auto& pj = gpos2.row(j);
171 
172  // Skip identical pairs.
173  if(t1 == t2)
174  continue;
175 
176  // Compute distance and pairwise function.
177  Vector3 rij = pi - pj;
178  snapshot.ApplyMinimumImage(&rij);
179  auto r = rij.norm();
180  val_ += pk_->Evaluate(r, df);
181 
182  // Get local indices and sum gradients.
183  auto lid1 = snapshot.GetLocalIndex(t1);
184  if(lid1 != -1)
185  grad_[lid1] += df*rij/r;
186 
187  auto lid2 = snapshot.GetLocalIndex(t2);
188  if(lid2 != -1)
189  grad_[lid2] -= df*rij/r;
190 
191  // Only sum box gradient on a single proc.
192  if(comm.rank() == 0)
193  boxgrad_ += df*rij/r*rij.transpose();
194  }
195  }
196  }
Label group1_
IDs of the first group of atoms.
Definition: PairwiseCV.h:43
std::vector< Vector3 > grad_
Gradient vector dCv/dxi.
virtual double Evaluate(double rij, double &df) const =0
Evaluate the pairwise kernel function.
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33
double val_
Current value of CV.
Matrix3 boxgrad_
Gradient w.r.t box vectors dCv/dHij.
PairwiseKernel * pk_
Pairwise kernel function used for CV.
Definition: PairwiseCV.h:45
Label group2_
IDs of the second group of atoms.
Definition: PairwiseCV.h:44

Here is the call graph for this function:

void SSAGES::PairwiseCV::Initialize ( const Snapshot snapshot)
inlineoverride

Initialize necessary variables.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 65 of file PairwiseCV.h.

References SSAGES::Snapshot::GetCommunicator(), SSAGES::Snapshot::GetLocalIndex(), group1_, and group2_.

66  {
67  using std::to_string;
68 
69  auto n1 = group1_.size(), n2 = group2_.size();
70 
71  // Make sure atom ID's are on at least one processor.
72  std::vector<int> found1(n1, 0), found2(n2, 0);
73  for(size_t i = 0; i < n1; ++i)
74  {
75  if(snapshot.GetLocalIndex(group1_[i]) != -1)
76  found1[i] = 1;
77  }
78 
79  for(size_t i = 0; i < n2; ++i)
80  {
81  if(snapshot.GetLocalIndex(group2_[i]) != -1)
82  found2[i] = 1;
83  }
84 
85  MPI_Allreduce(MPI_IN_PLACE, found1.data(), n1, MPI_INT, MPI_SUM, snapshot.GetCommunicator());
86  MPI_Allreduce(MPI_IN_PLACE, found2.data(), n2, MPI_INT, MPI_SUM, snapshot.GetCommunicator());
87 
88  unsigned ntot1 = std::accumulate(found1.begin(), found1.end(), 0, std::plus<int>());
89  unsigned ntot2 = std::accumulate(found2.begin(), found2.end(), 0, std::plus<int>());
90  if(ntot1 != n1)
91  throw BuildException({
92  "PairwiseCV: Expected to find " +
93  to_string(n1) +
94  " atoms in group 1, but only found " +
95  to_string(ntot1) + "."
96  });
97 
98  if(ntot2 != n2)
99  throw BuildException({
100  "PairwiseCV: Expected to find " +
101  to_string(n2) +
102  " atoms in group 2, but only found " +
103  to_string(ntot2) + "."
104  });
105  }
Label group1_
IDs of the first group of atoms.
Definition: PairwiseCV.h:43
Label group2_
IDs of the second group of atoms.
Definition: PairwiseCV.h:44

Here is the call graph for this function:


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