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

Collective variable on an particle position. More...

#include <ParticlePositionCV.h>

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

Public Member Functions

 ParticlePositionCV (const Label &atomids, const Vector3 &position)
 Constructor. More...
 
 ParticlePositionCV (const Label &atomids, const Vector3 &position, bool dimx, bool dimy, bool dimz)
 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 ParticlePositionCVBuild (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 atom ids of interest.
 
Vector3 point_
 Target point in space.
 
Bool3 dim_
 Each dimension determines if a constraint is applied by the 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

Collective variable on an particle position.

This CV will return the distance of a particle (collection of atoms) from a particular point in (1,2,3)-dimensional space.

Definition at line 39 of file ParticlePositionCV.h.

Constructor & Destructor Documentation

◆ ParticlePositionCV() [1/2]

SSAGES::ParticlePositionCV::ParticlePositionCV ( const Label atomids,
const Vector3 position 
)
inline

Constructor.

Parameters
atomidsVector of atom ID's.
positionPoint in (1,2,3)D space for the distance calculation.

Construct a particle position CV.

Todo:
Bounds needs to be an input.

Definition at line 61 of file ParticlePositionCV.h.

61  :
62  atomids_(atomids), point_(position), dim_{true, true, true}
63  {
64  }
Label atomids_
Vector of atom ids of interest.
Bool3 dim_
Each dimension determines if a constraint is applied by the CV.
Vector3 point_
Target point in space.

Referenced by Build().

Here is the caller graph for this function:

◆ ParticlePositionCV() [2/2]

SSAGES::ParticlePositionCV::ParticlePositionCV ( const Label atomids,
const Vector3 position,
bool  dimx,
bool  dimy,
bool  dimz 
)
inline

Constructor.

Parameters
atomidsVector of atom ID's.
positionPoint in (1,2,3)D space for the distance calculation.
dimxIf True, include x dimension.
dimyIf True, include y dimension.
dimzIf True, include z dimension.

Construct a particle position CV.

Todo:
Bounds needs to be an input.

Definition at line 78 of file ParticlePositionCV.h.

78  :
79  atomids_(atomids), point_(position), dim_{dimx, dimy, dimz}
80  {
81  }

Member Function Documentation

◆ Build()

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

149  {
150  Json::ObjectRequirement validator;
151  Json::Value schema;
152  Json::CharReaderBuilder rbuilder;
153  Json::CharReader* reader = rbuilder.newCharReader();
154 
155  reader->parse(JsonSchema::ParticlePositionCV.c_str(),
156  JsonSchema::ParticlePositionCV.c_str() + JsonSchema::ParticlePositionCV.size(),
157  &schema, nullptr);
158  validator.Parse(schema, path);
159 
160  // Validate inputs.
161  validator.Validate(json, path);
162  if(validator.HasErrors())
163  throw BuildException(validator.GetErrors());
164 
165  Label atomids;
166  for(auto& id : json["atom_ids"])
167  atomids.push_back(id.asInt());
168 
169  Vector3 position;
170  position[0] = json["position"][0].asDouble();
171  position[1] = json["position"][1].asDouble();
172  position[2] = json["position"][2].asDouble();
173 
175 
176  if(json.isMember("dimension"))
177  {
178  auto dimx = json["dimension"][0].asBool();
179  auto dimy = json["dimension"][1].asBool();
180  auto dimz = json["dimension"][2].asBool();
181  c = new ParticlePositionCV(atomids, position, dimx, dimy, dimz);
182  }
183  else
184  c = new ParticlePositionCV(atomids, position);
185 
186  return c;
187  }
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
ParticlePositionCV(const Label &atomids, const Vector3 &position)
Constructor.
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33
std::vector< int > Label
List of integers.
Definition: types.h:48

References Json::Requirement::GetErrors(), Json::Requirement::HasErrors(), Json::ObjectRequirement::Parse(), ParticlePositionCV(), 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::ParticlePositionCV::Evaluate ( const Snapshot snapshot)
inlineoverride

Evaluate the CV.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 116 of file ParticlePositionCV.h.

117  {
118  // Get local atom indices and compute COM.
119  std::vector<int> idx;
120  snapshot.GetLocalIndices(atomids_, &idx);
121 
122  // Get data from snapshot.
123  auto n = snapshot.GetNumAtoms();
124  const auto& masses = snapshot.GetMasses();
125 
126  // Initialize gradient.
127  std::fill(grad_.begin(), grad_.end(), Vector3{0,0,0});
128  grad_.resize(n, Vector3{0,0,0});
129 
130  // Compute total and center of mass.
131  auto masstot = snapshot.TotalMass(idx);
132  Vector3 com = snapshot.CenterOfMass(idx);
133 
134  // Compute difference between point and account for requested
135  // dimensions only.
136  Vector3 dx = snapshot.ApplyMinimumImage(com - point_).cwiseProduct(dim_.cast<double>());
137  val_ = dx.norm();
138 
139  // If distance is zero, we have nothing to do.
140  if(val_ == 0)
141  return;
142 
143  for(auto& id : idx)
144  grad_[id] = dx/val_*masses[id]/masstot;
145  }
std::vector< Vector3 > grad_
Gradient vector dCv/dxi.
double val_
Current value of CV.

References SSAGES::Snapshot::ApplyMinimumImage(), atomids_, SSAGES::Snapshot::CenterOfMass(), dim_, SSAGES::Snapshot::GetLocalIndices(), SSAGES::Snapshot::GetMasses(), SSAGES::Snapshot::GetNumAtoms(), SSAGES::CollectiveVariable::grad_, point_, SSAGES::Snapshot::TotalMass(), and SSAGES::CollectiveVariable::val_.

Here is the call graph for this function:

◆ Initialize()

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

Initialize necessary variables.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 87 of file ParticlePositionCV.h.

88  {
89  using std::to_string;
90 
91  auto n = atomids_.size();
92 
93  // Make sure atom ID's are on at least one processor.
94  std::vector<int> found(n, 0);
95  for(size_t i = 0; i < n; ++i)
96  {
97  if(snapshot.GetLocalIndex(atomids_[i]) != -1)
98  found[i] = 1;
99  }
100 
101  MPI_Allreduce(MPI_IN_PLACE, found.data(), n, MPI_INT, MPI_SUM, snapshot.GetCommunicator());
102  unsigned ntot = std::accumulate(found.begin(), found.end(), 0, std::plus<int>());
103  if(ntot != n)
104  throw BuildException({
105  "ParticlePositionCV: Expected to find " +
106  to_string(n) +
107  " atoms, but only found " +
108  to_string(ntot) + "."
109  });
110  }

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: