45#include "rnemd/NIVS.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"
75#define HONKING_LARGE_VALUE 1.0e10
77namespace OpenMD::RNEMD {
80 RNEMD {info, forceMan} {
81 rnemdMethodLabel_ =
"NIVS";
83 RNEMDParameters* rnemdParams = info->getSimParams()->getRNEMDParameters();
85 bool hasKineticFlux = rnemdParams->haveKineticFlux();
86 bool hasMomentumFlux = rnemdParams->haveMomentumFlux();
88 bool methodFluxMismatch =
false;
89 bool hasCorrectFlux =
false;
91 switch (rnemdFluxType_) {
95 hasCorrectFlux = hasKineticFlux;
100 hasCorrectFlux = hasMomentumFlux;
104 hasCorrectFlux = hasMomentumFlux && hasKineticFlux;
107 methodFluxMismatch =
true;
111 if (methodFluxMismatch) {
112 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
113 "RNEMD: The current method,\n"
115 "\tcannot be used with the current flux type, %s\n",
116 rnemdFluxTypeLabel_.c_str());
117 painCave.isFatal = 1;
118 painCave.severity = OPENMD_ERROR;
122 if (!hasCorrectFlux) {
123 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
124 "RNEMD: The current method, NIVS, and flux type, %s,\n"
125 "\tdid not have the correct flux value specified. Options\n"
126 "\tinclude: kineticFlux and momentumFlux.\n",
127 rnemdFluxTypeLabel_.c_str());
128 painCave.isFatal = 1;
129 painCave.severity = OPENMD_ERROR;
133 if (hasKineticFlux) {
134 setKineticFlux(rnemdParams->getKineticFlux());
139 if (hasMomentumFlux) {
140 RealType momentumFlux = rnemdParams->getMomentumFlux();
141 std::vector<RealType> momentumFluxVector(3);
143 switch (rnemdFluxType_) {
145 momentumFluxVector[0] = momentumFlux;
148 momentumFluxVector[1] = momentumFlux;
151 momentumFluxVector[2] = momentumFlux;
157 setMomentumFluxVector(momentumFluxVector);
163 if (!doRNEMD_)
return;
169 std::vector<StuntDouble*> hotBin, coldBin;
192 if (usePeriodicBoundaryConditions_) currentSnap_->wrapVector(pos);
197 hotBin.push_back(sd);
198 Phx += mass * vel.
x();
199 Phy += mass * vel.
y();
200 Phz += mass * vel.
z();
201 Khx += mass * vel.
x() * vel.
x();
202 Khy += mass * vel.
y() * vel.
y();
203 Khz += mass * vel.
z() * vel.
z();
212 angMom[j] * angMom[j] / I(j, j) + angMom[k] * angMom[k] / I(k, k);
214 Khw += angMom[0] * angMom[0] / I(0, 0) +
215 angMom[1] * angMom[1] / I(1, 1) +
216 angMom[2] * angMom[2] / I(2, 2);
226 if (usePeriodicBoundaryConditions_) currentSnap_->wrapVector(pos);
231 coldBin.push_back(sd);
232 Pcx += mass * vel.
x();
233 Pcy += mass * vel.
y();
234 Pcz += mass * vel.
z();
235 Kcx += mass * vel.
x() * vel.
x();
236 Kcy += mass * vel.
y() * vel.
y();
237 Kcz += mass * vel.
z() * vel.
z();
246 angMom[j] * angMom[j] / I(j, j) + angMom[k] * angMom[k] / I(k, k);
248 Kcw += angMom[0] * angMom[0] / I(0, 0) +
249 angMom[1] * angMom[1] / I(1, 1) +
250 angMom[2] * angMom[2] / I(2, 2);
265 MPI_Allreduce(MPI_IN_PLACE, &Phx, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
266 MPI_Allreduce(MPI_IN_PLACE, &Phy, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
267 MPI_Allreduce(MPI_IN_PLACE, &Phz, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
268 MPI_Allreduce(MPI_IN_PLACE, &Pcx, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
269 MPI_Allreduce(MPI_IN_PLACE, &Pcy, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
270 MPI_Allreduce(MPI_IN_PLACE, &Pcz, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
272 MPI_Allreduce(MPI_IN_PLACE, &Khx, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
273 MPI_Allreduce(MPI_IN_PLACE, &Khy, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
274 MPI_Allreduce(MPI_IN_PLACE, &Khz, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
275 MPI_Allreduce(MPI_IN_PLACE, &Khw, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
277 MPI_Allreduce(MPI_IN_PLACE, &Kcx, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
278 MPI_Allreduce(MPI_IN_PLACE, &Kcy, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
279 MPI_Allreduce(MPI_IN_PLACE, &Kcz, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
280 MPI_Allreduce(MPI_IN_PLACE, &Kcw, 1, MPI_REALTYPE, MPI_SUM, MPI_COMM_WORLD);
284 RealType px = Pcx / Phx;
285 RealType py = Pcy / Phy;
286 RealType pz = Pcz / Phz;
287 RealType c(0.0), x(0.0), y(0.0), z(0.0);
288 bool successfulScale =
false;
289 if ((rnemdFluxType_ == rnemdFullKE) || (rnemdFluxType_ == rnemdRotKE)) {
292 if (rnemdFluxType_ == rnemdFullKE) {
293 c = 1.0 - kineticTarget_ / (Kcx + Kcy + Kcz + Kcw);
295 c = 1.0 - kineticTarget_ / Kcw;
298 if ((c > 0.81) && (c < 1.21)) {
302 if (rnemdFluxType_ == rnemdFullKE) {
303 x = 1.0 + px * (1.0 - c);
304 y = 1.0 + py * (1.0 - c);
305 z = 1.0 + pz * (1.0 - c);
314 if ((std::fabs(x - 1.0) < 0.1) && (std::fabs(y - 1.0) < 0.1) &&
315 (std::fabs(z - 1.0) < 0.1)) {
316 w = 1.0 + (kineticTarget_ + Khx * (1.0 - x * x) +
317 Khy * (1.0 - y * y) + Khz * (1.0 - z * z)) /
321 w = 1.0 + kineticTarget_ / Khw;
323 if ((w > 0.81) && (w < 1.21)) {
325 std::vector<StuntDouble*>::iterator sdi;
327 for (sdi = coldBin.begin(); sdi != coldBin.end(); ++sdi) {
328 if (rnemdFluxType_ == rnemdFullKE) {
329 vel = (*sdi)->getVel() * c;
332 if ((*sdi)->isDirectional()) {
333 Vector3d angMom = (*sdi)->getJ() * c;
334 (*sdi)->setJ(angMom);
338 for (sdi = hotBin.begin(); sdi != hotBin.end(); ++sdi) {
339 if (rnemdFluxType_ == rnemdFullKE) {
340 vel = (*sdi)->getVel();
346 if ((*sdi)->isDirectional()) {
347 Vector3d angMom = (*sdi)->getJ() * w;
348 (*sdi)->setJ(angMom);
351 successfulScale =
true;
352 kineticExchange_ += kineticTarget_;
356 RealType a000(0.0), a110(0.0), c0(0.0);
357 RealType a001(0.0), a111(0.0), b01(0.0), b11(0.0), c1(0.0);
358 switch (rnemdFluxType_) {
375 c0 = kineticTarget_ - Kcx - Kcy - Kcz;
376 a001 = Khx * px * px + Khy * py * py;
377 a111 = Khz * pz * pz;
378 b01 = -2.0 * (Khx * px * (1.0 + px) + Khy * py * (1.0 + py));
379 b11 = -2.0 * Khz * pz * (1.0 + pz);
380 c1 = Khx * px * (2.0 + px) + Khy * py * (2.0 + py) +
381 Khz * pz * (2.0 + pz) - kineticTarget_;
384 c = 1 - momentumTarget_.x() / Pcx;
387 c0 = Kcx * c * c - Kcx - Kcy - Kcz;
388 a001 = py * py * Khy;
389 a111 = pz * pz * Khz;
390 b01 = -2.0 * Khy * py * (1.0 + py);
391 b11 = -2.0 * Khz * pz * (1.0 + pz);
392 c1 = Khy * py * (2.0 + py) + Khz * pz * (2.0 + pz) +
393 Khx * (fastpow(c * px - px - 1.0, 2) - 1.0);
396 c = 1 - momentumTarget_.y() / Pcy;
399 c0 = Kcy * c * c - Kcx - Kcy - Kcz;
400 a001 = px * px * Khx;
401 a111 = pz * pz * Khz;
402 b01 = -2.0 * Khx * px * (1.0 + px);
403 b11 = -2.0 * Khz * pz * (1.0 + pz);
404 c1 = Khx * px * (2.0 + px) + Khz * pz * (2.0 + pz) +
405 Khy * (fastpow(c * py - py - 1.0, 2) - 1.0);
408 c = 1 - momentumTarget_.z() / Pcz;
411 c0 = Kcz * c * c - Kcx - Kcy - Kcz;
412 a001 = px * px * Khx;
413 a111 = py * py * Khy;
414 b01 = -2.0 * Khx * px * (1.0 + px);
415 b11 = -2.0 * Khy * py * (1.0 + py);
416 c1 = Khx * px * (2.0 + px) + Khy * py * (2.0 + py) +
417 Khz * (fastpow(c * pz - pz - 1.0, 2) - 1.0);
423 RealType v1 = a000 * a111 - a001 * a110;
424 RealType v2 = a000 * b01;
425 RealType v3 = a000 * b11;
426 RealType v4 = a000 * c1 - a001 * c0;
427 RealType v8 = a110 * b01;
428 RealType v10 = -b01 * c0;
430 RealType u0 = v2 * v10 - v4 * v4;
431 RealType u1 = -2.0 * v3 * v4;
432 RealType u2 = -v2 * v8 - v3 * v3 - 2.0 * v1 * v4;
433 RealType u3 = -2.0 * v1 * v3;
434 RealType u4 = -v1 * v1;
436 RealType maxAbs = fabs(u0);
437 if (maxAbs < fabs(u1)) maxAbs = fabs(u1);
438 if (maxAbs < fabs(u2)) maxAbs = fabs(u2);
439 if (maxAbs < fabs(u3)) maxAbs = fabs(u3);
440 if (maxAbs < fabs(u4)) maxAbs = fabs(u4);
453 vector<RealType> realRoots = poly.FindRealRoots();
455 vector<RealType>::iterator ri;
456 RealType r1, r2, alpha0;
457 vector<pair<RealType, RealType>> rps;
458 for (ri = realRoots.begin(); ri != realRoots.end(); ++ri) {
461 if (fabs(u0 + r2 * (u1 + r2 * (u2 + r2 * (u3 + r2 * u4)))) > 1e-6) {
462 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
463 "RNEMD Warning: polynomial solve seems to have an error!");
464 painCave.isFatal = 0;
469 alpha0 = -c0 - a110 * r2 * r2;
471 r1 = sqrt(alpha0 / a000);
472 if (fabs(c1 + r1 * (b01 + r1 * a001) + r2 * (b11 + r2 * a111)) <
474 rps.push_back(make_pair(r1, r2));
478 if (fabs(c1 + r1 * (b01 + r1 * a001) + r2 * (b11 + r2 * a111)) <
480 rps.push_back(make_pair(r1, r2));
489 RealType smallestDiff = HONKING_LARGE_VALUE;
491 std::pair<RealType, RealType> bestPair = std::make_pair(1.0, 1.0);
492 std::vector<std::pair<RealType, RealType>>::iterator rpi;
493 for (rpi = rps.begin(); rpi != rps.end(); ++rpi) {
496 switch (rnemdFluxType_) {
498 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2) +
499 fastpow(r1 * r1 / r2 / r2 - Kcz / Kcx, 2) +
500 fastpow(r1 * r1 / r2 / r2 - Kcz / Kcy, 2);
503 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2) +
504 fastpow(r1 * r1 / r2 / r2 - Kcz / Kcy, 2);
507 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2) +
508 fastpow(r1 * r1 / r2 / r2 - Kcz / Kcx, 2);
511 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2) +
512 fastpow(r1 * r1 / r2 / r2 - Kcy / Kcx, 2);
516 if (diff < smallestDiff) {
522 if (worldRank == 0) {
534 switch (rnemdFluxType_) {
558 vector<StuntDouble*>::iterator sdi;
560 for (sdi = coldBin.begin(); sdi != coldBin.end(); ++sdi) {
561 vel = (*sdi)->getVel();
568 x = 1.0 + px * (1.0 - x);
569 y = 1.0 + py * (1.0 - y);
570 z = 1.0 + pz * (1.0 - z);
571 for (sdi = hotBin.begin(); sdi != hotBin.end(); ++sdi) {
572 vel = (*sdi)->getVel();
578 successfulScale =
true;
579 switch (rnemdFluxType_) {
581 kineticExchange_ += kineticTarget_;
586 momentumExchange_ += momentumTarget_;
593 if (successfulScale !=
true) {
594 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
595 "NIVS exchange NOT performed - roots that solve\n"
596 "\tthe constraint equations may not exist or there may be\n"
597 "\tno selected objects in one or both slabs.\n");
598 painCave.isFatal = 0;
599 painCave.severity = OPENMD_INFO;
ForceManager is responsible for calculating both the short range (bonded) interactions and long range...
A generic Polynomial class.
void setCoefficient(int exponent, const Real &coefficient)
Set the coefficent of the specified exponent, if the coefficient is already there,...
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...
"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 & z()
Returns reference of the third element of Vector3.
Real & x()
Returns reference of the first element of Vector3.
Real & y()
Returns reference of the second element of Vector3.