SSAGES
0.9.3
Software Suite for Advanced General Ensemble Simulations
|
Abstract class for a collective variable. More...
#include <CollectiveVariable.h>
Public Member Functions | |
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 CollectiveVariable * | BuildCV (const Json::Value &json, const std::string &path) |
Set up collective variable. More... | |
Protected Attributes | |
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. | |
Abstract class for a collective variable.
Definition at line 40 of file CollectiveVariable.h.
|
static |
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 43 of file CollectiveVariable.cpp.
References SSAGES::AlphaRMSDCV::Build(), SSAGES::AngleCV::Build(), SSAGES::ANNCV::Build(), SSAGES::AntiBetaRMSDCV::Build(), SSAGES::BoxVolumeCV::Build(), SSAGES::GyrationTensorCV::Build(), SSAGES::PairwiseCV::Build(), SSAGES::ParallelBetaRMSDCV::Build(), SSAGES::ParticleCoordinateCV::Build(), SSAGES::ParticlePositionCV::Build(), SSAGES::ParticleSeparationCV::Build(), SSAGES::RMSDCV::Build(), SSAGES::RouseModeCV::Build(), and SSAGES::TorsionalCV::Build().
Referenced by SSAGES::ResourceHandler::Build().
|
pure virtual |
Evaluate CV.
Evaluation of the CV. This function is called by the Hook in the post-integration phase every iteration. The CV should compute its value and gradient, storing them in a local private variable.
|
inline |
Get CV boundaries.
Returns the boundaries of the CV. These represent the bounds within which the CV is expected to be constrained. There is no requirement on the method to respect the values returned here.
Definition at line 151 of file CollectiveVariable.h.
|
inline |
Get gradient contribution to box.
Returns the gradient of the CV with respect to the box.
Definition at line 138 of file CollectiveVariable.h.
|
inlinevirtual |
Get difference between current CV value and a given value, taking periodic boundaries into account.
location | Value whose distance from the current CV value should be calculated. |
Returns the difference betwen the current cv value and Location: (value_ - Location) respecting periodic boundary conditions of the CV, if the CV has periodic boundary conditions. For example Torsional angle has boundaries at pi and -pi, in which the difference beteen the angles is 0 not 2pi
Reimplemented in SSAGES::TorsionalCV, and SSAGES::ParticleCoordinateCV.
Definition at line 168 of file CollectiveVariable.h.
|
inline |
Get current gradient of the CV.
Returns the current value of the CV gradient. This should be an n length vector, where n is the number of atoms in the snapshot. Each element in the vector is the derivative of the CV with respect to the atom's coordinate (dCV/dxi).
Definition at line 127 of file CollectiveVariable.h.
|
inlinevirtual |
Returns the minimum image of a CV based on the input location.
Takes the input location and applies the periodic boundary conditions to return a minimum image of the CV.
Reimplemented in SSAGES::TorsionalCV.
Definition at line 98 of file CollectiveVariable.h.
|
inlinevirtual |
Apply periodic boundaries to a given value.
location | Value to which the periodic boundaries should be applied. |
Takes location and applies periodic boundaries of the CV on it and returns a correct value. Example would be torsional angle which has bounds at pi and -pi. If location = 2pi, GetPeriodicValue(location) would return 0.
Reimplemented in SSAGES::TorsionalCV, and SSAGES::ParticleCoordinateCV.
Definition at line 113 of file CollectiveVariable.h.
|
inline |
Get current value of the CV.
Returns the current value of the CV which has been computed before via the call to CollectiveVariable::Evaluate().
Definition at line 86 of file CollectiveVariable.h.
|
inlinevirtual |
Initialize CV.
Initialization of the CV. This is an optional method and is called during the pre-simulation phase of the hook. It is typically used to allocate/reserve memory.
Definition at line 69 of file CollectiveVariable.h.