OpenMD 3.1
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
VSS.cpp
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45#include "rnemd/VSS.hpp"
46
47#include <algorithm>
48#include <cmath>
49#include <map>
50#include <set>
51#include <sstream>
52#include <string>
53#include <vector>
54
55#ifdef IS_MPI
56#include <mpi.h>
57#endif
58
60#include "brains/Thermo.hpp"
61#include "io/Globals.hpp"
62#include "math/ConvexHull.hpp"
63#include "math/Polynomial.hpp"
65#include "math/Vector.hpp"
66#include "math/Vector3.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"
74
75namespace OpenMD::RNEMD {
76
77 VSSMethod::VSSMethod(SimInfo* info, ForceManager* forceMan) :
78 RNEMD {info, forceMan} {
79 rnemdMethodLabel_ = "VSS";
80
81 RNEMDParameters* rnemdParams = info->getSimParams()->getRNEMDParameters();
82
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();
89
90 bool methodFluxMismatch = false;
91 bool hasCorrectFlux = false;
92
93 switch (rnemdFluxType_) {
94 case rnemdKE:
95 case rnemdRotKE:
96 case rnemdFullKE:
97 hasCorrectFlux = hasKineticFlux;
98 break;
99 case rnemdPx:
100 case rnemdPy:
101 case rnemdPz:
102 hasCorrectFlux = hasMomentumFlux;
103 break;
104 case rnemdLx:
105 case rnemdLy:
106 case rnemdLz:
107 hasCorrectFlux = hasAngularMomentumFlux;
108 break;
109 case rnemdPvector:
110 hasCorrectFlux = hasMomentumFluxVector;
111 break;
112 case rnemdLvector:
113 hasCorrectFlux = hasAngularMomentumFluxVector;
114 break;
115 case rnemdKePx:
116 case rnemdKePy:
117 hasCorrectFlux = hasMomentumFlux && hasKineticFlux;
118 break;
119 case rnemdKeLx:
120 case rnemdKeLy:
121 case rnemdKeLz:
122 hasCorrectFlux = hasAngularMomentumFlux && hasKineticFlux;
123 break;
124 case rnemdKePvector:
125 hasCorrectFlux = hasMomentumFluxVector && hasKineticFlux;
126 break;
127 case rnemdKeLvector:
128 hasCorrectFlux = hasAngularMomentumFluxVector && hasKineticFlux;
129 break;
130 default:
131 methodFluxMismatch = true;
132 break;
133 }
134
135 if (methodFluxMismatch) {
136 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
137 "RNEMD: The current method,\n"
138 "\t\tVSS\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;
143 simError();
144 }
145
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;
155 simError();
156 }
157
158 if (hasKineticFlux) {
159 setKineticFlux(rnemdParams->getKineticFlux());
160 } else {
161 setKineticFlux(0.0);
162 }
163
164 if (hasMomentumFluxVector) {
165 setMomentumFluxVector(rnemdParams->getMomentumFluxVector());
166 } else {
167 std::vector<RealType> momentumFluxVector(3);
168
169 if (hasMomentumFlux) {
170 RealType momentumFlux = rnemdParams->getMomentumFlux();
171
172 switch (rnemdFluxType_) {
173 case rnemdPx:
174 case rnemdKePx:
175 momentumFluxVector[0] = momentumFlux;
176 break;
177 case rnemdPy:
178 case rnemdKePy:
179 momentumFluxVector[1] = momentumFlux;
180 break;
181 case rnemdPz:
182 momentumFluxVector[2] = momentumFlux;
183 break;
184 default:
185 break;
186 }
187 }
188
189 setMomentumFluxVector(momentumFluxVector);
190 }
191
192 if (hasAngularMomentumFluxVector) {
193 setAngularMomentumFluxVector(rnemdParams->getAngularMomentumFluxVector());
194 } else {
195 std::vector<RealType> angularMomentumFluxVector(3);
196
197 if (hasAngularMomentumFlux) {
198 RealType angularMomentumFlux = rnemdParams->getAngularMomentumFlux();
199
200 switch (rnemdFluxType_) {
201 case rnemdLx:
202 case rnemdKeLx:
203 angularMomentumFluxVector[0] = angularMomentumFlux;
204 break;
205 case rnemdLy:
206 case rnemdKeLy:
207 angularMomentumFluxVector[1] = angularMomentumFlux;
208 break;
209 case rnemdLz:
210 case rnemdKeLz:
211 angularMomentumFluxVector[2] = angularMomentumFlux;
212 default:
213 break;
214 }
215 }
216
217 setAngularMomentumFluxVector(angularMomentumFluxVector);
218 }
219 }
220
221 void VSSMethod::doRNEMDImpl(SelectionManager& smanA,
222 SelectionManager& smanB) {
223 if (!doRNEMD_) return;
224 int selei;
225 int selej;
226
227 StuntDouble* sd;
228
229 vector<StuntDouble*> hotBin, coldBin;
230
231 Vector3d Ph(V3Zero);
232 Vector3d Lh(V3Zero);
233 RealType Mh = 0.0;
234 Mat3x3d Ih(0.0);
235 RealType Kh = 0.0;
236 Vector3d Pc(V3Zero);
237 Vector3d Lc(V3Zero);
238 RealType Mc = 0.0;
239 Mat3x3d Ic(0.0);
240 RealType Kc = 0.0;
241
242 // Constraints can be on only the linear or angular momentum, but
243 // not both. Usually, the user will specify which they want, but
244 // in case they don't, the use of periodic boundaries should make
245 // the choice for us.
246 bool doLinearPart = false;
247 bool doAngularPart = false;
248
249 switch (rnemdFluxType_) {
250 case rnemdPx:
251 case rnemdPy:
252 case rnemdPz:
253 case rnemdPvector:
254 case rnemdKePx:
255 case rnemdKePy:
256 case rnemdKePvector:
257 doLinearPart = true;
258 break;
259 case rnemdLx:
260 case rnemdLy:
261 case rnemdLz:
262 case rnemdLvector:
263 case rnemdKeLx:
264 case rnemdKeLy:
265 case rnemdKeLz:
266 case rnemdKeLvector:
267 doAngularPart = true;
268 break;
269 case rnemdKE:
270 case rnemdRotKE:
271 case rnemdFullKE:
272 default:
273 if (usePeriodicBoundaryConditions_)
274 doLinearPart = true;
275 else
276 doAngularPart = true;
277 break;
278 }
279
280 for (sd = smanA.beginSelected(selei); sd != NULL;
281 sd = smanA.nextSelected(selei)) {
282 Vector3d pos = sd->getPos();
283
284 // wrap the stuntdouble's position back into the box:
285 if (usePeriodicBoundaryConditions_) currentSnap_->wrapVector(pos);
286
287 RealType mass = sd->getMass();
288 Vector3d vel = sd->getVel();
289 Vector3d rPos = sd->getPos() - coordinateOrigin_;
290 RealType r2;
291
292 hotBin.push_back(sd);
293 Ph += mass * vel;
294 Mh += mass;
295 Kh += mass * vel.lengthSquare();
296 Lh += mass * cross(rPos, vel);
297 Ih -= outProduct(rPos, rPos) * mass;
298 r2 = rPos.lengthSquare();
299 Ih(0, 0) += mass * r2;
300 Ih(1, 1) += mass * r2;
301 Ih(2, 2) += mass * r2;
302
303 if (rnemdFluxType_ == rnemdFullKE) {
304 if (sd->isDirectional()) {
305 Vector3d angMom = sd->getJ();
306 Mat3x3d I = sd->getI();
307 if (sd->isLinear()) {
308 int i = sd->linearAxis();
309 int j = (i + 1) % 3;
310 int k = (i + 2) % 3;
311 Kh += angMom[j] * angMom[j] / I(j, j) +
312 angMom[k] * angMom[k] / I(k, k);
313 } else {
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);
317 }
318 }
319 }
320 }
321
322 for (sd = smanB.beginSelected(selej); sd != NULL;
323 sd = smanB.nextSelected(selej)) {
324 Vector3d pos = sd->getPos();
325
326 // wrap the stuntdouble's position back into the box:
327 if (usePeriodicBoundaryConditions_) currentSnap_->wrapVector(pos);
328
329 RealType mass = sd->getMass();
330 Vector3d vel = sd->getVel();
331 Vector3d rPos = sd->getPos() - coordinateOrigin_;
332 RealType r2;
333
334 coldBin.push_back(sd);
335 Pc += mass * vel;
336 Mc += mass;
337 Kc += mass * vel.lengthSquare();
338 Lc += mass * cross(rPos, vel);
339 Ic -= outProduct(rPos, rPos) * mass;
340 r2 = rPos.lengthSquare();
341 Ic(0, 0) += mass * r2;
342 Ic(1, 1) += mass * r2;
343 Ic(2, 2) += mass * r2;
344
345 if (rnemdFluxType_ == rnemdFullKE) {
346 if (sd->isDirectional()) {
347 Vector3d angMom = sd->getJ();
348 Mat3x3d I = sd->getI();
349 if (sd->isLinear()) {
350 int i = sd->linearAxis();
351 int j = (i + 1) % 3;
352 int k = (i + 2) % 3;
353 Kc += angMom[j] * angMom[j] / I(j, j) +
354 angMom[k] * angMom[k] / I(k, k);
355 } else {
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);
359 }
360 }
361 }
362 }
363
364 Kh *= 0.5;
365 Kc *= 0.5;
366
367#ifdef IS_MPI
368 MPI_Allreduce(MPI_IN_PLACE, &Ph[0], 3, MPI_REALTYPE, MPI_SUM,
369 MPI_COMM_WORLD);
370 MPI_Allreduce(MPI_IN_PLACE, &Pc[0], 3, MPI_REALTYPE, MPI_SUM,
371 MPI_COMM_WORLD);
372 MPI_Allreduce(MPI_IN_PLACE, &Lh[0], 3, MPI_REALTYPE, MPI_SUM,
373 MPI_COMM_WORLD);
374 MPI_Allreduce(MPI_IN_PLACE, &Lc[0], 3, MPI_REALTYPE, MPI_SUM,
375 MPI_COMM_WORLD);
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,
381 MPI_COMM_WORLD);
382 MPI_Allreduce(MPI_IN_PLACE, Ic.getArrayPointer(), 9, MPI_REALTYPE, MPI_SUM,
383 MPI_COMM_WORLD);
384#endif
385
386 Vector3d ac, acrec, bc, bcrec;
387 Vector3d ah, ahrec, bh, bhrec;
388
389 bool successfulExchange = false;
390 if ((Mh > 0.0) && (Mc > 0.0)) { // both slabs are not empty
391
392 Vector3d vc = Pc / Mc;
393 ac = -momentumTarget_ / Mc + vc;
394 acrec = -momentumTarget_ / Mc;
395
396 // We now need the inverse of the inertia tensor to calculate the
397 // angular velocity of the cold slab;
398 Mat3x3d Ici = Ic.inverse();
399 Vector3d omegac = Ici * Lc;
400 bc = -(Ici * angularMomentumTarget_) + omegac;
401 bcrec = bc - omegac;
402
403 RealType cNumerator = Kc - kineticTarget_;
404 if (doLinearPart) cNumerator -= 0.5 * Mc * ac.lengthSquare();
405
406 if (doAngularPart) cNumerator -= 0.5 * (dot(bc, Ic * bc));
407
408 RealType cDenominator = Kc;
409
410 if (doLinearPart) cDenominator -= 0.5 * Mc * vc.lengthSquare();
411
412 if (doAngularPart) cDenominator -= 0.5 * (dot(omegac, Ic * omegac));
413
414 if (cNumerator / cDenominator > 0.0) {
415 RealType c = sqrt(cNumerator / cDenominator);
416
417 if ((c > 0.9) && (c < 1.1)) { // restrict scaling coefficients
418
419 Vector3d vh = Ph / Mh;
420 ah = momentumTarget_ / Mh + vh;
421 ahrec = momentumTarget_ / Mh;
422
423 // We now need the inverse of the inertia tensor to
424 // calculate the angular velocity of the hot slab;
425 Mat3x3d Ihi = Ih.inverse();
426 Vector3d omegah = Ihi * Lh;
427 bh = (Ihi * angularMomentumTarget_) + omegah;
428 bhrec = bh - omegah;
429
430 RealType hNumerator = Kh + kineticTarget_;
431 if (doLinearPart) hNumerator -= 0.5 * Mh * ah.lengthSquare();
432
433 if (doAngularPart) hNumerator -= 0.5 * (dot(bh, Ih * bh));
434
435 RealType hDenominator = Kh;
436 if (doLinearPart) hDenominator -= 0.5 * Mh * vh.lengthSquare();
437 if (doAngularPart) hDenominator -= 0.5 * (dot(omegah, Ih * omegah));
438
439 if (hNumerator / hDenominator > 0.0) {
440 RealType h = sqrt(hNumerator / hDenominator);
441
442 if ((h > 0.9) && (h < 1.1)) {
443 vector<StuntDouble*>::iterator sdi;
444 Vector3d vel;
445 Vector3d rPos;
446
447 for (sdi = coldBin.begin(); sdi != coldBin.end(); ++sdi) {
448 if (doLinearPart) vel = ((*sdi)->getVel() - vc) * c + ac;
449 if (doAngularPart) {
450 rPos = (*sdi)->getPos() - coordinateOrigin_;
451 vel = ((*sdi)->getVel() - cross(omegac, rPos)) * c +
452 cross(bc, rPos);
453 }
454
455 (*sdi)->setVel(vel);
456
457 if (rnemdFluxType_ == rnemdFullKE) {
458 if ((*sdi)->isDirectional()) {
459 Vector3d angMom = (*sdi)->getJ() * c;
460 (*sdi)->setJ(angMom);
461 }
462 }
463 }
464
465 for (sdi = hotBin.begin(); sdi != hotBin.end(); ++sdi) {
466 if (doLinearPart) vel = ((*sdi)->getVel() - vh) * h + ah;
467 if (doAngularPart) {
468 rPos = (*sdi)->getPos() - coordinateOrigin_;
469 vel = ((*sdi)->getVel() - cross(omegah, rPos)) * h +
470 cross(bh, rPos);
471 }
472
473 (*sdi)->setVel(vel);
474
475 if (rnemdFluxType_ == rnemdFullKE) {
476 if ((*sdi)->isDirectional()) {
477 Vector3d angMom = (*sdi)->getJ() * h;
478 (*sdi)->setJ(angMom);
479 }
480 }
481 }
482
483 successfulExchange = true;
484 kineticExchange_ += kineticTarget_;
485 momentumExchange_ += momentumTarget_;
486 angularMomentumExchange_ += angularMomentumTarget_;
487 }
488 }
489 }
490 }
491 }
492
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;
500 simError();
501 failTrialCount_++;
502 }
503 }
504} // namespace OpenMD::RNEMD
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...
Definition SimInfo.hpp:93
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.
Definition Vector.hpp:399
Vector3< Real > cross(const Vector3< Real > &v1, const Vector3< Real > &v2)
Returns the cross product of two Vectors.
Definition Vector3.hpp:136
Real dot(const DynamicVector< Real > &v1, const DynamicVector< Real > &v2)
Returns the dot product of two DynamicVectors.