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

Collective variable to calculate angle. More...

#include <AngleCV.h>

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

Public Member Functions

 AngleCV (int atomid1, int atomid2, int atomid3)
 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.
 
virtual ~CollectiveVariable ()
 Destructor.
 
virtual void Initialize (const class Snapshot &)
 Initialize CV. More...
 
virtual void Evaluate (const class Snapshot &)=0
 Evaluate CV. More...
 
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. More...
 
const std::array< double, 2 > & GetBoundaries ()
 Get CV boundaries. More...
 
virtual double GetDifference (double location) const
 

Static Public Member Functions

static AngleCVBuild (const Json::Value &json, const std::string &path)
 Set up collective variable. More...
 
- 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 atomids_
 Vector of 3 atom ID's of interest.
 

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

Collective variable to calculate angle.

Definition at line 35 of file AngleCV.h.

Constructor & Destructor Documentation

◆ AngleCV()

SSAGES::AngleCV::AngleCV ( int  atomid1,
int  atomid2,
int  atomid3 
)
inline

Constructor.

Parameters
atomid1ID of the first atom defining the angle.
atomid2ID of the second atom defining the angle.
atomid3ID of the third atom defining the angle.

Construct an dihedral CV.

Todo:
Bounds needs to be an input and periodic boundary conditions

Definition at line 53 of file AngleCV.h.

53  :
54  atomids_({atomid1, atomid2, atomid3})
55  {
56  bounds_ = {{0,M_PI}};
57  }
Label atomids_
Vector of 3 atom ID's of interest.
Definition: AngleCV.h:39
std::array< double, 2 > bounds_
Bounds on CV.

Referenced by Build().

Here is the caller graph for this function:

Member Function Documentation

◆ Build()

static AngleCV* SSAGES::AngleCV::Build ( const Json::Value &  json,
const std::string &  path 
)
inlinestatic

Set up collective variable.

Parameters
jsonJSON input value.
pathPath for JSON path specification.
Returns
Pointer to the CV built by this function. nullptr in case of unknown error.

Builds a CV from a JSON node. Returns a pointer to the built cv. If an unknown error is encountered, this function will return a nullptr, but generally it will throw a BuildException on failure.

Warning
Object lifetime is the caller's responsibility.

Definition at line 137 of file AngleCV.h.

138  {
139  Json::ObjectRequirement validator;
140  Json::Value schema;
141  Json::CharReaderBuilder rbuilder;
142  Json::CharReader* reader = rbuilder.newCharReader();
143 
144  reader->parse(JsonSchema::AngleCV.c_str(),
145  JsonSchema::AngleCV.c_str() + JsonSchema::AngleCV.size(),
146  &schema, nullptr);
147  validator.Parse(schema, path);
148 
149  // Validate inputs.
150  validator.Validate(json, path);
151  if(validator.HasErrors())
152  throw BuildException(validator.GetErrors());
153 
154  std::vector<int> atomids;
155  for(auto& s : json["atom_ids"])
156  atomids.push_back(s.asInt());
157 
158  return new AngleCV(atomids[0], atomids[1], atomids[2]);
159  }
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
AngleCV(int atomid1, int atomid2, int atomid3)
Constructor.
Definition: AngleCV.h:53

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

Referenced by SSAGES::CollectiveVariable::BuildCV().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Evaluate()

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

Evaluate the CV.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 85 of file AngleCV.h.

86  {
87  // Get data from snapshot.
88  auto n = snapshot.GetNumAtoms();
89  const auto& pos = snapshot.GetPositions();
90  auto& comm = snapshot.GetCommunicator();
91 
92  // Initialize gradient.
93  std::fill(grad_.begin(), grad_.end(), Vector3{0,0,0});
94  grad_.resize(n, Vector3{0,0,0});
95 
96  Vector3 xi{0, 0, 0}, xj{0, 0, 0}, xk{0, 0, 0};
97 
98  auto iindex = snapshot.GetLocalIndex(atomids_[0]);
99  auto jindex = snapshot.GetLocalIndex(atomids_[1]);
100  auto kindex = snapshot.GetLocalIndex(atomids_[2]);
101 
102  if(iindex != -1) xi = pos[iindex];
103  if(jindex != -1) xj = pos[jindex];
104  if(kindex != -1) xk = pos[kindex];
105 
106  // By performing a reduce, we actually collect all. This can
107  // be converted to a more intelligent allgater on rank then bcast.
108  MPI_Allreduce(MPI_IN_PLACE, xi.data(), 3, MPI_DOUBLE, MPI_SUM, comm);
109  MPI_Allreduce(MPI_IN_PLACE, xj.data(), 3, MPI_DOUBLE, MPI_SUM, comm);
110  MPI_Allreduce(MPI_IN_PLACE, xk.data(), 3, MPI_DOUBLE, MPI_SUM, comm);
111 
112  // Two vectors
113  auto rij = snapshot.ApplyMinimumImage(xi - xj);
114  auto rkj = snapshot.ApplyMinimumImage(xk - xj);
115 
116  auto dotP = rij.dot(rkj);
117  auto nrij = rij.norm();
118  auto nrkj = rkj.norm();
119  auto dotPnn = dotP/(nrij*nrkj);
120 
121  val_ = acos(dotPnn);
122 
123  // Calculate gradients
124  auto prefactor = -1.0/(sqrt(1. - dotPnn*dotPnn)*nrij*nrkj);
125 
126  Vector3 gradi{0, 0, 0}, gradk{0, 0, 0};
127  if(iindex != -1) gradi = prefactor*(rkj - dotP*rij/(nrij*nrij));
128  if(kindex != -1) gradk = prefactor*(rij - dotP*rkj/(nrkj*nrkj));
129  MPI_Allreduce(MPI_IN_PLACE, gradi.data(), 3, MPI_DOUBLE, MPI_SUM, comm);
130  MPI_Allreduce(MPI_IN_PLACE, gradk.data(), 3, MPI_DOUBLE, MPI_SUM, comm);
131  if(iindex != -1) grad_[iindex] = gradi;
132  if(kindex != -1) grad_[kindex] = gradk;
133  if(jindex != -1) grad_[jindex] = -gradi - gradk;
134  }
std::vector< Vector3 > grad_
Gradient vector dCv/dxi.
double val_
Current value of CV.
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33

References SSAGES::Snapshot::ApplyMinimumImage(), atomids_, SSAGES::Snapshot::GetCommunicator(), SSAGES::Snapshot::GetLocalIndex(), SSAGES::Snapshot::GetNumAtoms(), SSAGES::Snapshot::GetPositions(), SSAGES::CollectiveVariable::grad_, and SSAGES::CollectiveVariable::val_.

Here is the call graph for this function:

◆ Initialize()

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

Initialize necessary variables.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 63 of file AngleCV.h.

64  {
65  using std::to_string;
66 
67  std::vector<int> found;
68  snapshot.GetLocalIndices(atomids_, &found);
69  int nfound = found.size();
70  MPI_Allreduce(MPI_IN_PLACE, &nfound, 1, MPI_INT, MPI_SUM, snapshot.GetCommunicator());
71 
72  if(nfound != 3)
73  throw BuildException({
74  "TorsionalCV: Expected to find " +
75  to_string(3) +
76  " atoms, but only found " +
77  to_string(nfound) + "."
78  });
79  }

References atomids_, SSAGES::Snapshot::GetCommunicator(), and SSAGES::Snapshot::GetLocalIndices().

Here is the call graph for this function:

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