SSAGES
0.9.3
Software Suite for Advanced General Ensemble Simulations
|
Collective variable on the distance between two particles' centers of mass. More...
#include <ParticleSeparationCV.h>
Public Member Functions | |
ParticleSeparationCV (const Label &group1, const Label &group2) | |
Constructor. More... | |
ParticleSeparationCV (const Label &group1, const Label &group2, 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 Matrix3 & | GetBoxGradient () 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 ParticleSeparationCV * | Build (const Json::Value &json, const std::string &path) |
Set up collective variable. More... | |
Static Public Member Functions inherited from SSAGES::CollectiveVariable | |
static CollectiveVariable * | BuildCV (const Json::Value &json, const std::string &path) |
Set up collective variable. More... | |
Private Attributes | |
Label | group1_ |
Atom ID's of group 1. | |
Label | group2_ |
Atom ID's of group 2. | |
Bool3 | dim_ |
Each dimension determines if it is applied by the CV. | |
Additional Inherited Members | |
Protected Attributes inherited from SSAGES::CollectiveVariable | |
std::vector< Vector3 > | grad_ |
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. | |
Collective variable on the distance between two particles' centers of mass.
Collective variable on two particle positions. This CV will return the distance between two specific atom groups of the simulation.
Definition at line 39 of file ParticleSeparationCV.h.
|
inline |
Constructor.
group1 | Vector of atom ID's of the first particle. |
group2 | Vector of atom ID's of the second particle. |
Construct a paarticle separation CV.
Definition at line 58 of file ParticleSeparationCV.h.
Referenced by Build().
|
inline |
Constructor.
group1 | Vector of atom ID's of the first particle. |
group2 | Vector of atom ID's of the second particle. |
dimx | If True , include x dimension. |
dimy | If True , include y dimension. |
dimz | If True , include z dimension. |
Construct a particle separation CV.
Definition at line 74 of file ParticleSeparationCV.h.
|
inlinestatic |
Set up collective variable.
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.
Definition at line 167 of file ParticleSeparationCV.h.
References Json::Requirement::GetErrors(), Json::Requirement::HasErrors(), Json::ObjectRequirement::Parse(), ParticleSeparationCV(), and Json::ObjectRequirement::Validate().
Referenced by SSAGES::CollectiveVariable::BuildCV().
|
inlineoverride |
Evaluate the CV.
snapshot | Current simulation snapshot. |
Definition at line 128 of file ParticleSeparationCV.h.
References SSAGES::Snapshot::ApplyMinimumImage(), SSAGES::CollectiveVariable::boxgrad_, SSAGES::Snapshot::CenterOfMass(), dim_, SSAGES::Snapshot::GetLocalIndices(), SSAGES::Snapshot::GetMasses(), SSAGES::Snapshot::GetNumAtoms(), SSAGES::CollectiveVariable::grad_, group1_, group2_, SSAGES::Snapshot::TotalMass(), and SSAGES::CollectiveVariable::val_.
|
inlineoverride |
Initialize necessary variables.
snapshot | Current simulation snapshot. |
Definition at line 82 of file ParticleSeparationCV.h.
References SSAGES::Snapshot::GetCommunicator(), SSAGES::Snapshot::GetLocalIndex(), group1_, and group2_.