45#include "rnemd/VSS.hpp"
60#include "brains/Thermo.hpp"
61#include "io/Globals.hpp"
62#include "math/ConvexHull.hpp"
69#include "rnemd/RNEMD.hpp"
70#include "rnemd/RNEMDParameters.hpp"
71#include "types/FixedChargeAdapter.hpp"
72#include "types/FluctuatingChargeAdapter.hpp"
73#include "utils/Constants.hpp"
75namespace OpenMD::RNEMD {
78 RNEMD {info, forceMan} {
79 rnemdMethodLabel_ =
"VSS";
81 RNEMDParameters* rnemdParams = info->getSimParams()->getRNEMDParameters();
83 bool hasKineticFlux = rnemdParams->haveKineticFlux();
84 bool hasMomentumFlux = rnemdParams->haveMomentumFlux();
85 bool hasMomentumFluxVector = rnemdParams->haveMomentumFluxVector();
86 bool hasAngularMomentumFlux = rnemdParams->haveAngularMomentumFlux();
87 bool hasAngularMomentumFluxVector =
88 rnemdParams->haveAngularMomentumFluxVector();
90 bool methodFluxMismatch =
false;
91 bool hasCorrectFlux =
false;
93 switch (rnemdFluxType_) {
97 hasCorrectFlux = hasKineticFlux;
102 hasCorrectFlux = hasMomentumFlux;
107 hasCorrectFlux = hasAngularMomentumFlux;
110 hasCorrectFlux = hasMomentumFluxVector;
113 hasCorrectFlux = hasAngularMomentumFluxVector;
117 hasCorrectFlux = hasMomentumFlux && hasKineticFlux;
122 hasCorrectFlux = hasAngularMomentumFlux && hasKineticFlux;
125 hasCorrectFlux = hasMomentumFluxVector && hasKineticFlux;
128 hasCorrectFlux = hasAngularMomentumFluxVector && hasKineticFlux;
131 methodFluxMismatch =
true;
135 if (methodFluxMismatch) {
136 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
137 "RNEMD: The current method,\n"
139 "\tcannot be used with the current flux type, %s\n",
140 rnemdFluxTypeLabel_.c_str());
141 painCave.isFatal = 1;
142 painCave.severity = OPENMD_ERROR;
146 if (!hasCorrectFlux) {
147 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
148 "RNEMD: The current method, VSS, and flux type, %s,\n"
149 "\tdid not have the correct flux value specified. Options\n"
150 "\tinclude: kineticFlux, momentumFlux, angularMomentumFlux,\n"
151 "\tmomentumFluxVector, and angularMomentumFluxVector.\n",
152 rnemdFluxTypeLabel_.c_str());
153 painCave.isFatal = 1;
154 painCave.severity = OPENMD_ERROR;
158 if (hasKineticFlux) {
159 setKineticFlux(rnemdParams->getKineticFlux());
164 if (hasMomentumFluxVector) {
165 setMomentumFluxVector(rnemdParams->getMomentumFluxVector());
167 std::vector<RealType> momentumFluxVector(3);
169 if (hasMomentumFlux) {
170 RealType momentumFlux = rnemdParams->getMomentumFlux();
172 switch (rnemdFluxType_) {
175 momentumFluxVector[0] = momentumFlux;
179 momentumFluxVector[1] = momentumFlux;
182 momentumFluxVector[2] = momentumFlux;
189 setMomentumFluxVector(momentumFluxVector);
192 if (hasAngularMomentumFluxVector) {
193 setAngularMomentumFluxVector(rnemdParams->getAngularMomentumFluxVector());
195 std::vector<RealType> angularMomentumFluxVector(3);
197 if (hasAngularMomentumFlux) {
198 RealType angularMomentumFlux = rnemdParams->getAngularMomentumFlux();
200 switch (rnemdFluxType_) {
203 angularMomentumFluxVector[0] = angularMomentumFlux;
207 angularMomentumFluxVector[1] = angularMomentumFlux;
211 angularMomentumFluxVector[2] = angularMomentumFlux;
217 setAngularMomentumFluxVector(angularMomentumFluxVector);
223 if (!doRNEMD_)
return;
229 vector<StuntDouble*> hotBin, coldBin;
246 bool doLinearPart =
false;
247 bool doAngularPart =
false;
249 switch (rnemdFluxType_) {
267 doAngularPart =
true;
273 if (usePeriodicBoundaryConditions_)
276 doAngularPart =
true;
285 if (usePeriodicBoundaryConditions_) currentSnap_->wrapVector(pos);
292 hotBin.push_back(sd);
296 Lh += mass *
cross(rPos, vel);
297 Ih -= outProduct(rPos, rPos) * mass;
299 Ih(0, 0) += mass * r2;
300 Ih(1, 1) += mass * r2;
301 Ih(2, 2) += mass * r2;
303 if (rnemdFluxType_ == rnemdFullKE) {
311 Kh += angMom[j] * angMom[j] / I(j, j) +
312 angMom[k] * angMom[k] / I(k, k);
314 Kh += angMom[0] * angMom[0] / I(0, 0) +
315 angMom[1] * angMom[1] / I(1, 1) +
316 angMom[2] * angMom[2] / I(2, 2);
327 if (usePeriodicBoundaryConditions_) currentSnap_->wrapVector(pos);
334 coldBin.push_back(sd);
338 Lc += mass *
cross(rPos, vel);
339 Ic -= outProduct(rPos, rPos) * mass;
341 Ic(0, 0) += mass * r2;
342 Ic(1, 1) += mass * r2;
343 Ic(2, 2) += mass * r2;
345 if (rnemdFluxType_ == rnemdFullKE) {
353 Kc += angMom[j] * angMom[j] / I(j, j) +
354 angMom[k] * angMom[k] / I(k, k);
356 Kc += angMom[0] * angMom[0] / I(0, 0) +
357 angMom[1] * angMom[1] / I(1, 1) +
358 angMom[2] * angMom[2] / I(2, 2);
368 MPI_Allreduce(MPI_IN_PLACE, &Ph[0], 3, MPI_REALTYPE, MPI_SUM,
370 MPI_Allreduce(MPI_IN_PLACE, &Pc[0], 3, MPI_REALTYPE, MPI_SUM,
372 MPI_Allreduce(MPI_IN_PLACE, &Lh[0], 3, MPI_REALTYPE, MPI_SUM,
374 MPI_Allreduce(MPI_IN_PLACE, &Lc[0], 3, MPI_REALTYPE, MPI_SUM,
376 MPI_Allreduce(MPI_IN_PLACE, &Mh, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
377 MPI_Allreduce(MPI_IN_PLACE, &Kh, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
378 MPI_Allreduce(MPI_IN_PLACE, &Mc, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
379 MPI_Allreduce(MPI_IN_PLACE, &Kc, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
380 MPI_Allreduce(MPI_IN_PLACE, Ih.getArrayPointer(), 9, MPI_REALTYPE, MPI_SUM,
382 MPI_Allreduce(MPI_IN_PLACE, Ic.getArrayPointer(), 9, MPI_REALTYPE, MPI_SUM,
389 bool successfulExchange =
false;
390 if ((Mh > 0.0) && (Mc > 0.0)) {
393 ac = -momentumTarget_ / Mc + vc;
394 acrec = -momentumTarget_ / Mc;
400 bc = -(Ici * angularMomentumTarget_) + omegac;
403 RealType cNumerator = Kc - kineticTarget_;
404 if (doLinearPart) cNumerator -= 0.5 * Mc * ac.
lengthSquare();
406 if (doAngularPart) cNumerator -= 0.5 * (
dot(bc, Ic * bc));
408 RealType cDenominator = Kc;
410 if (doLinearPart) cDenominator -= 0.5 * Mc * vc.
lengthSquare();
412 if (doAngularPart) cDenominator -= 0.5 * (
dot(omegac, Ic * omegac));
414 if (cNumerator / cDenominator > 0.0) {
415 RealType c = sqrt(cNumerator / cDenominator);
417 if ((c > 0.9) && (c < 1.1)) {
420 ah = momentumTarget_ / Mh + vh;
421 ahrec = momentumTarget_ / Mh;
427 bh = (Ihi * angularMomentumTarget_) + omegah;
430 RealType hNumerator = Kh + kineticTarget_;
431 if (doLinearPart) hNumerator -= 0.5 * Mh * ah.
lengthSquare();
433 if (doAngularPart) hNumerator -= 0.5 * (
dot(bh, Ih * bh));
435 RealType hDenominator = Kh;
436 if (doLinearPart) hDenominator -= 0.5 * Mh * vh.
lengthSquare();
437 if (doAngularPart) hDenominator -= 0.5 * (
dot(omegah, Ih * omegah));
439 if (hNumerator / hDenominator > 0.0) {
440 RealType h = sqrt(hNumerator / hDenominator);
442 if ((h > 0.9) && (h < 1.1)) {
443 vector<StuntDouble*>::iterator sdi;
447 for (sdi = coldBin.begin(); sdi != coldBin.end(); ++sdi) {
448 if (doLinearPart) vel = ((*sdi)->getVel() - vc) * c + ac;
450 rPos = (*sdi)->getPos() - coordinateOrigin_;
451 vel = ((*sdi)->getVel() -
cross(omegac, rPos)) * c +
457 if (rnemdFluxType_ == rnemdFullKE) {
458 if ((*sdi)->isDirectional()) {
459 Vector3d angMom = (*sdi)->getJ() * c;
460 (*sdi)->setJ(angMom);
465 for (sdi = hotBin.begin(); sdi != hotBin.end(); ++sdi) {
466 if (doLinearPart) vel = ((*sdi)->getVel() - vh) * h + ah;
468 rPos = (*sdi)->getPos() - coordinateOrigin_;
469 vel = ((*sdi)->getVel() -
cross(omegah, rPos)) * h +
475 if (rnemdFluxType_ == rnemdFullKE) {
476 if ((*sdi)->isDirectional()) {
477 Vector3d angMom = (*sdi)->getJ() * h;
478 (*sdi)->setJ(angMom);
483 successfulExchange =
true;
484 kineticExchange_ += kineticTarget_;
485 momentumExchange_ += momentumTarget_;
486 angularMomentumExchange_ += angularMomentumTarget_;
493 if (successfulExchange !=
true) {
494 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
495 "VSS exchange NOT performed - roots that solve\n"
496 "\tthe constraint equations may not exist or there may be\n"
497 "\tno selected objects in one or both slabs.\n");
498 painCave.isFatal = 0;
499 painCave.severity = OPENMD_INFO;
ForceManager is responsible for calculating both the short range (bonded) interactions and long range...
StuntDouble * nextSelected(int &i)
Finds the next selected StuntDouble in the selection.
StuntDouble * beginSelected(int &i)
Finds the first selected StuntDouble in the selection.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
SquareMatrix3< Real > inverse() const
Sets the value of this matrix to the inverse of itself.
"Don't move, or you're dead! Stand up! Captain, we've got them!"
Vector3d getVel()
Returns the current velocity of this stuntDouble.
int linearAxis()
Returns the linear axis of the rigidbody, atom and directional atom will always return -1.
RealType getMass()
Returns the mass of this stuntDouble.
virtual Mat3x3d getI()=0
Returns the inertia tensor of this stuntDouble.
bool isLinear()
Tests the if this stuntDouble is a linear rigidbody.
Vector3d getPos()
Returns the current position of this stuntDouble.
Vector3d getJ()
Returns the current angular momentum of this stuntDouble (body -fixed).
bool isDirectional()
Tests if this stuntDouble is a directional one.
Real lengthSquare()
Returns the squared length of this vector.
Vector3< Real > cross(const Vector3< Real > &v1, const Vector3< Real > &v2)
Returns the cross product of two Vectors.
Real dot(const DynamicVector< Real > &v1, const DynamicVector< Real > &v2)
Returns the dot product of two DynamicVectors.