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::ParticleCoordinateCV Class Reference

Collective variable on a particle coordinate. More...

#include <ParticleCoordinateCV.h>

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

Public Member Functions

 ParticleCoordinateCV (const Label &atomids, Dimension dim)
 Constructor. More...
 
void Initialize (const Snapshot &snapshot) override
 Initialize necessary variables. More...
 
void Evaluate (const Snapshot &snapshot) override
 Evaluate the CV. More...
 
double GetPeriodicValue (double Location) const override
 Return value taking periodic boundary conditions into account. More...
 
double GetDifference (const double Location) const override
 Return difference considering periodic boundaries. 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...
 
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...
 

Static Public Member Functions

static ParticleCoordinateCVBuild (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_
 IDs of atoms of interest.
 
Dimension dim_
 Index of dimension.
 

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 on a particle coordinate.

This will return the value of either the x, y, or z coordinate, depending on the user specification for a defined particle, which is a collection of one or more atoms.

Definition at line 40 of file ParticleCoordinateCV.h.

Constructor & Destructor Documentation

◆ ParticleCoordinateCV()

SSAGES::ParticleCoordinateCV::ParticleCoordinateCV ( const Label atomids,
Dimension  dim 
)
inline

Constructor.

Parameters
atomidsAtom ID's of interest.
dimIndex of dimension.

Construct a particle coordinate CV. The atomids specify a vector of the atom ID's of interest, and index specifies the dimension to report (x,y,z).

Todo:
Bounds needs to be an input.

Definition at line 60 of file ParticleCoordinateCV.h.

60  :
61  atomids_(atomids), dim_(dim)
62  {}
Dimension dim_
Index of dimension.
Label atomids_
IDs of atoms of interest.

Member Function Documentation

◆ Build()

static ParticleCoordinateCV* SSAGES::ParticleCoordinateCV::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 168 of file ParticleCoordinateCV.h.

169  {
170  Json::ObjectRequirement validator;
171  Json::Value schema;
172  Json::CharReaderBuilder rbuilder;
173  Json::CharReader* reader = rbuilder.newCharReader();
174 
175  reader->parse(JsonSchema::ParticleCoordinateCV.c_str(),
176  JsonSchema::ParticleCoordinateCV.c_str() + JsonSchema::ParticleCoordinateCV.size(),
177  &schema, nullptr);
178  validator.Parse(schema, path);
179 
180  // Validate inputs.
181  validator.Validate(json, path);
182  if(validator.HasErrors())
183  throw BuildException(validator.GetErrors());
184 
185  Label atomids;
186  for(auto& id : json["atom_ids"])
187  atomids.push_back(id.asInt());
188 
189  auto indextype = json.get("dimension","x").asString();
190 
191  Dimension index;
192  if(indextype == "x")
193  index = Dimension::x;
194  else if(indextype == "y")
195  index = Dimension::y;
196  else if(indextype == "z")
197  index = Dimension::z;
198  else
199  throw BuildException({"Could not obtain ParticleCoordinate dimension specified."});
200 
201  return new ParticleCoordinateCV(atomids, index);
202  }
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
ParticleCoordinateCV(const Label &atomids, Dimension dim)
Constructor.
Dimension
Enum for dimension.
Definition: types.h:54
std::vector< int > Label
List of integers.
Definition: types.h:48

References 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::ParticleCoordinateCV::Evaluate ( const Snapshot snapshot)
inlineoverride

Evaluate the CV.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 97 of file ParticleCoordinateCV.h.

98  {
99  // Get local atom indices and compute COM.
100  std::vector<int> idx;
101  snapshot.GetLocalIndices(atomids_, &idx);
102 
103  // Get data from snapshot.
104  auto n = snapshot.GetNumAtoms();
105  const auto& masses = snapshot.GetMasses();
106 
107  // Initialize gradient.
108  std::fill(grad_.begin(), grad_.end(), Vector3{0,0,0});
109  grad_.resize(n, Vector3{0,0,0});
110 
111  // Compute total and center of mass.
112  auto masstot = snapshot.TotalMass(idx);
113  Vector3 com = snapshot.CenterOfMass(idx);
114 
115  // Assign CV value.
116  switch(dim_)
117  {
118  case Dimension::x:
119  val_ = com[0];
120  break;
121  case Dimension::y:
122  val_ = com[1];
123  break;
124  case Dimension::z:
125  val_ = com[2];
126  break;
127  }
128 
129  // Assign gradient to appropriate atoms.
130  for(auto& id : idx)
131  {
132  grad_[id][0] = (dim_ == Dimension::x) ? masses[id]/masstot : 0;
133  grad_[id][1] = (dim_ == Dimension::y) ? masses[id]/masstot : 0;
134  grad_[id][2] = (dim_ == Dimension::z) ? masses[id]/masstot : 0;
135  }
136  }
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 atomids_, SSAGES::Snapshot::CenterOfMass(), dim_, SSAGES::Snapshot::GetLocalIndices(), SSAGES::Snapshot::GetMasses(), SSAGES::Snapshot::GetNumAtoms(), SSAGES::CollectiveVariable::grad_, and SSAGES::Snapshot::TotalMass().

Here is the call graph for this function:

◆ GetDifference()

double SSAGES::ParticleCoordinateCV::GetDifference ( const double  Location) const
inlineoverridevirtual

Return difference considering periodic boundaries.

Parameters
LocationCalculate difference of CV value to this location.
Returns
Direct difference.

As the AtomCoordinate CV does not consider periodic boundary conditions, the difference between the current value of the CV and another value is always the direct difference.

Reimplemented from SSAGES::CollectiveVariable.

Definition at line 162 of file ParticleCoordinateCV.h.

163  {
164  return val_ - Location;
165  }

References SSAGES::CollectiveVariable::val_.

◆ GetPeriodicValue()

double SSAGES::ParticleCoordinateCV::GetPeriodicValue ( double  Location) const
inlineoverridevirtual

Return value taking periodic boundary conditions into account.

Parameters
LocationGet wrapped value of this location.
Returns
Input value

The AtomCoordinate CV does not consider periodic boundary conditions. Thus, this function always returns the input value.

Reimplemented from SSAGES::CollectiveVariable.

Definition at line 147 of file ParticleCoordinateCV.h.

148  {
149  return Location;
150  }

◆ Initialize()

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

Initialize necessary variables.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 68 of file ParticleCoordinateCV.h.

69  {
70  using std::to_string;
71 
72  auto n = atomids_.size();
73 
74  // Make sure atom ID's are on at least one processor.
75  std::vector<int> found(n, 0);
76  for(size_t i = 0; i < n; ++i)
77  {
78  if(snapshot.GetLocalIndex(atomids_[i]) != -1)
79  found[i] = 1;
80  }
81 
82  MPI_Allreduce(MPI_IN_PLACE, found.data(), n, MPI_INT, MPI_SUM, snapshot.GetCommunicator());
83  unsigned ntot = std::accumulate(found.begin(), found.end(), 0, std::plus<int>());
84  if(ntot != n)
85  throw BuildException({
86  "ParticleCoordinateCV: Expected to find " +
87  to_string(n) +
88  " atoms, but only found " +
89  to_string(ntot) + "."
90  });
91  }

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

Here is the call graph for this function:

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