SSAGES  0.1
A MetaDynamics Package
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Groups Pages
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, bool fixx, bool fixy, bool fixz)
 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.
 
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.
 
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)
 
- 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 fix_
 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

SSAGES::ParticlePositionCV::ParticlePositionCV ( const Label atomids,
const Vector3 position,
bool  fixx,
bool  fixy,
bool  fixz 
)
inline

Constructor.

Parameters
atomidsVector of atom ID's.
positionPoint in (1,2,3)D space for the distance calculation.
fixxIf False, constrain x dimension.
fixyIf False, constrain y dimension.
fixzIf False, constrain z dimension.

Construct a particle position CV. If a dimension is constrained, this dimension does not contribute to the distance calculation. For example, if the y- and z-dimension are constrained, the CV calculates only the distance in x-direction.

Todo:
Bounds needs to be an input.

Definition at line 67 of file ParticlePositionCV.h.

67  :
68  atomids_(atomids), point_(position), fix_{fixx, fixy, fixz}
69  {
70  }
Vector3 point_
Target point in space.
Bool3 fix_
Each dimension determines if a constraint is applied by the CV.
Label atomids_
Vector of atom ids of interest.

Member Function Documentation

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

Evaluate the CV.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 105 of file ParticlePositionCV.h.

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

106  {
107  // Get local atom indices and compute COM.
108  std::vector<int> idx;
109  snapshot.GetLocalIndices(atomids_, &idx);
110 
111  // Get data from snapshot.
112  auto n = snapshot.GetNumAtoms();
113  const auto& masses = snapshot.GetMasses();
114 
115  // Initialize gradient.
116  std::fill(grad_.begin(), grad_.end(), Vector3{0,0,0});
117  grad_.resize(n, Vector3{0,0,0});
118 
119  // Compute total and center of mass.
120  auto masstot = snapshot.TotalMass(idx);
121  Vector3 com = snapshot.CenterOfMass(idx, masstot);
122 
123  // Compute difference between point and account for requested
124  // dimensions only.
125  Vector3 dx = snapshot.ApplyMinimumImage(com - point_).cwiseProduct(fix_.cast<double>());
126  val_ = dx.norm();
127 
128  // If distance is zero, we have nothing to do.
129  if(val_ == 0)
130  return;
131 
132  for(auto& id : idx)
133  grad_[id] = dx/val_*masses[id]/masstot;
134  }
Vector3 point_
Target point in space.
Bool3 fix_
Each dimension determines if a constraint is applied by the CV.
std::vector< Vector3 > grad_
Gradient vector dCv/dxi.
Label atomids_
Vector of atom ids of interest.
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33
double val_
Current value of CV.

Here is the call graph for this function:

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

Initialize necessary variables.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 76 of file ParticlePositionCV.h.

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

77  {
78  using std::to_string;
79 
80  auto n = atomids_.size();
81 
82  // Make sure atom ID's are on at least one processor.
83  std::vector<int> found(n, 0);
84  for(size_t i = 0; i < n; ++i)
85  {
86  if(snapshot.GetLocalIndex(atomids_[i]) != -1)
87  found[i] = 1;
88  }
89 
90  MPI_Allreduce(MPI_IN_PLACE, found.data(), n, MPI_INT, MPI_SUM, snapshot.GetCommunicator());
91  unsigned ntot = std::accumulate(found.begin(), found.end(), 0, std::plus<int>());
92  if(ntot != n)
93  throw BuildException({
94  "ParticlePositionCV: Expected to find " +
95  to_string(n) +
96  " atoms, but only found " +
97  to_string(ntot) + "."
98  });
99  }
Label atomids_
Vector of atom ids of interest.

Here is the call graph for this function:


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