23 #include "CollectiveVariable.h"
24 #include "Validator/ObjectRequirement.h"
25 #include "Drivers/DriverException.h"
89 std::vector<int> found1(n1, 0), found2(n2, 0);
90 for(
size_t i = 0; i < n1; ++i)
96 for(
size_t i = 0; i < n2; ++i)
102 MPI_Allreduce(MPI_IN_PLACE, found1.data(), n1, MPI_INT, MPI_SUM, snapshot.
GetCommunicator());
103 MPI_Allreduce(MPI_IN_PLACE, found2.data(), n2, MPI_INT, MPI_SUM, snapshot.
GetCommunicator());
105 unsigned ntot1 = std::accumulate(found1.begin(), found1.end(), 0, std::plus<int>());
106 unsigned ntot2 = std::accumulate(found2.begin(), found2.end(), 0, std::plus<int>());
109 "ParticleSeparationCV: Expected to find " +
111 " atoms in particle 1, but only found " +
112 to_string(ntot1) +
"."
117 "ParticleSeparationCV: Expected to find " +
119 " atoms in particle 2, but only found " +
120 to_string(ntot2) +
"."
131 std::vector<int> idx1, idx2;
137 const auto& masses = snapshot.
GetMasses();
171 Json::CharReaderBuilder rbuilder;
172 Json::CharReader* reader = rbuilder.newCharReader();
174 reader->parse(JsonSchema::ParticleSeparationCV.c_str(),
175 JsonSchema::ParticleSeparationCV.c_str() + JsonSchema::ParticleSeparationCV.size(),
177 validator.
Parse(schema, path);
184 std::vector<int> group1, group2;
186 for(
auto& s : json[
"group1"])
187 group1.push_back(s.asInt());
189 for(
auto& s : json[
"group2"])
190 group2.push_back(s.asInt());
193 if(json.isMember(
"dimension"))
195 auto fixx = json[
"dimension"][0].asBool();
196 auto fixy = json[
"dimension"][1].asBool();
197 auto fixz = json[
"dimension"][2].asBool();
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.
bool HasErrors()
Check if errors have occured.
Exception to be thrown when building the Driver fails.
Abstract class for a collective variable.
Matrix3 boxgrad_
Gradient w.r.t box vectors dCv/dHij.
std::vector< Vector3 > grad_
Gradient vector dCv/dxi.
double val_
Current value of CV.
Collective variable on the distance between two particles' centers of mass.
Bool3 dim_
Each dimension determines if it is applied by the CV.
static ParticleSeparationCV * Build(const Json::Value &json, const std::string &path)
Set up collective variable.
Label group1_
Atom ID's of group 1.
void Initialize(const Snapshot &snapshot) override
Initialize necessary variables.
Label group2_
Atom ID's of group 2.
ParticleSeparationCV(const Label &group1, const Label &group2, bool dimx, bool dimy, bool dimz)
Constructor.
ParticleSeparationCV(const Label &group1, const Label &group2)
Constructor.
void Evaluate(const Snapshot &snapshot) override
Evaluate the CV.
Class containing a snapshot of the current simulation in time.
unsigned GetNumAtoms() const
Get number of atoms in this snapshot.
void GetLocalIndices(const Label &ids, Label *indices) const
const mxx::comm & GetCommunicator() const
Get communicator for walker.
int GetLocalIndex(int id) const
Gets the local atom index corresponding to an atom ID.
double TotalMass(const Label &indices) const
Compute the total mass of a group of particles based on index.
const std::vector< double > & GetMasses() const
Const access to the particle masses.
Vector3 CenterOfMass(const Label &indices, bool mass_weight=true) const
Compute center of mass of a group of atoms based on index.
void ApplyMinimumImage(Vector3 *v) const
Apply minimum image to a vector.
Eigen::Vector3d Vector3
Three-dimensional vector.
Eigen::Matrix< bool, 3, 1 > Bool3
Three-dimensional boolean.
std::vector< int > Label
List of integers.