ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/rnemd/RNEMD.cpp
Revision: 1774
Committed: Wed Aug 8 16:03:02 2012 UTC (12 years, 8 months ago) by gezelter
File size: 56094 byte(s)
Log Message:
Cleaning up RNEMD log, adding XY area to Snapshot

File Contents

# Content
1 /*
2 * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 *
4 * The University of Notre Dame grants you ("Licensee") a
5 * non-exclusive, royalty free, license to use, modify and
6 * redistribute this software in source and binary code form, provided
7 * that the following conditions are met:
8 *
9 * 1. Redistributions of source code must retain the above copyright
10 * notice, this list of conditions and the following disclaimer.
11 *
12 * 2. Redistributions in binary form must reproduce the above copyright
13 * notice, this list of conditions and the following disclaimer in the
14 * documentation and/or other materials provided with the
15 * distribution.
16 *
17 * This software is provided "AS IS," without a warranty of any
18 * kind. All express or implied conditions, representations and
19 * warranties, including any implied warranty of merchantability,
20 * fitness for a particular purpose or non-infringement, are hereby
21 * excluded. The University of Notre Dame and its licensors shall not
22 * be liable for any damages suffered by licensee as a result of
23 * using, modifying or distributing the software or its
24 * derivatives. In no event will the University of Notre Dame or its
25 * licensors be liable for any lost revenue, profit or data, or for
26 * direct, indirect, special, consequential, incidental or punitive
27 * damages, however caused and regardless of the theory of liability,
28 * arising out of the use of or inability to use software, even if the
29 * University of Notre Dame has been advised of the possibility of
30 * such damages.
31 *
32 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33 * research, please cite the appropriate papers when you publish your
34 * work. Good starting points are:
35 *
36 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 * [4] Vardeman & Gezelter, in progress (2009).
40 */
41
42 #include <cmath>
43 #include "rnemd/RNEMD.hpp"
44 #include "math/Vector3.hpp"
45 #include "math/Vector.hpp"
46 #include "math/SquareMatrix3.hpp"
47 #include "math/Polynomial.hpp"
48 #include "primitives/Molecule.hpp"
49 #include "primitives/StuntDouble.hpp"
50 #include "utils/PhysicalConstants.hpp"
51 #include "utils/Tuple.hpp"
52 #ifdef IS_MPI
53 #include <mpi.h>
54 #endif
55
56 #define HONKING_LARGE_VALUE 1.0e10
57
58 using namespace std;
59 namespace OpenMD {
60
61 RNEMD::RNEMD(SimInfo* info) : info_(info), evaluator_(info), seleMan_(info),
62 usePeriodicBoundaryConditions_(info->getSimParams()->getUsePeriodicBoundaryConditions()) {
63
64 trialCount_ = 0;
65 failTrialCount_ = 0;
66 failRootCount_ = 0;
67
68 int seedValue;
69 Globals * simParams = info->getSimParams();
70 RNEMDParameters* rnemdParams = simParams->getRNEMDParameters();
71
72 stringToMethod_["Swap"] = rnemdSwap;
73 stringToMethod_["NIVS"] = rnemdNIVS;
74 stringToMethod_["VSS"] = rnemdVSS;
75
76 stringToFluxType_["KE"] = rnemdKE;
77 stringToFluxType_["Px"] = rnemdPx;
78 stringToFluxType_["Py"] = rnemdPy;
79 stringToFluxType_["Pz"] = rnemdPz;
80 stringToFluxType_["KE+Px"] = rnemdKePx;
81 stringToFluxType_["KE+Py"] = rnemdKePy;
82 stringToFluxType_["KE+Pvector"] = rnemdKePvector;
83
84 runTime_ = simParams->getRunTime();
85 statusTime_ = simParams->getStatusTime();
86
87 rnemdObjectSelection_ = rnemdParams->getObjectSelection();
88 evaluator_.loadScriptString(rnemdObjectSelection_);
89 seleMan_.setSelectionSet(evaluator_.evaluate());
90
91 const string methStr = rnemdParams->getMethod();
92 bool hasFluxType = rnemdParams->haveFluxType();
93
94 string fluxStr;
95 if (hasFluxType) {
96 fluxStr = rnemdParams->getFluxType();
97 } else {
98 sprintf(painCave.errMsg,
99 "RNEMD: No fluxType was set in the md file. This parameter,\n"
100 "\twhich must be one of the following values:\n"
101 "\tKE, Px, Py, Pz, KE+Px, KE+Py, KE+Pvector, must be set to\n"
102 "\tuse RNEMD\n");
103 painCave.isFatal = 1;
104 painCave.severity = OPENMD_ERROR;
105 simError();
106 }
107
108 bool hasKineticFlux = rnemdParams->haveKineticFlux();
109 bool hasMomentumFlux = rnemdParams->haveMomentumFlux();
110 bool hasMomentumFluxVector = rnemdParams->haveMomentumFluxVector();
111 bool hasSlabWidth = rnemdParams->haveSlabWidth();
112 bool hasSlabACenter = rnemdParams->haveSlabACenter();
113 bool hasSlabBCenter = rnemdParams->haveSlabBCenter();
114 bool hasOutputFileName = rnemdParams->haveOutputFileName();
115 bool hasOutputFields = rnemdParams->haveOutputFields();
116
117 map<string, RNEMDMethod>::iterator i;
118 i = stringToMethod_.find(methStr);
119 if (i != stringToMethod_.end())
120 rnemdMethod_ = i->second;
121 else {
122 sprintf(painCave.errMsg,
123 "RNEMD: The current method,\n"
124 "\t\t%s is not one of the recognized\n"
125 "\texchange methods: Swap, NIVS, or VSS\n",
126 methStr.c_str());
127 painCave.isFatal = 1;
128 painCave.severity = OPENMD_ERROR;
129 simError();
130 }
131
132 map<string, RNEMDFluxType>::iterator j;
133 j = stringToFluxType_.find(fluxStr);
134 if (j != stringToFluxType_.end())
135 rnemdFluxType_ = j->second;
136 else {
137 sprintf(painCave.errMsg,
138 "RNEMD: The current fluxType,\n"
139 "\t\t%s\n"
140 "\tis not one of the recognized flux types.\n",
141 fluxStr.c_str());
142 painCave.isFatal = 1;
143 painCave.severity = OPENMD_ERROR;
144 simError();
145 }
146
147 bool methodFluxMismatch = false;
148 bool hasCorrectFlux = false;
149 switch(rnemdMethod_) {
150 case rnemdSwap:
151 switch (rnemdFluxType_) {
152 case rnemdKE:
153 hasCorrectFlux = hasKineticFlux;
154 break;
155 case rnemdPx:
156 case rnemdPy:
157 case rnemdPz:
158 hasCorrectFlux = hasMomentumFlux;
159 break;
160 default :
161 methodFluxMismatch = true;
162 break;
163 }
164 break;
165 case rnemdNIVS:
166 switch (rnemdFluxType_) {
167 case rnemdKE:
168 case rnemdRotKE:
169 case rnemdFullKE:
170 hasCorrectFlux = hasKineticFlux;
171 break;
172 case rnemdPx:
173 case rnemdPy:
174 case rnemdPz:
175 hasCorrectFlux = hasMomentumFlux;
176 break;
177 case rnemdKePx:
178 case rnemdKePy:
179 hasCorrectFlux = hasMomentumFlux && hasKineticFlux;
180 break;
181 default:
182 methodFluxMismatch = true;
183 break;
184 }
185 break;
186 case rnemdVSS:
187 switch (rnemdFluxType_) {
188 case rnemdKE:
189 case rnemdRotKE:
190 case rnemdFullKE:
191 hasCorrectFlux = hasKineticFlux;
192 break;
193 case rnemdPx:
194 case rnemdPy:
195 case rnemdPz:
196 hasCorrectFlux = hasMomentumFlux;
197 break;
198 case rnemdPvector:
199 hasCorrectFlux = hasMomentumFluxVector;
200 case rnemdKePx:
201 case rnemdKePy:
202 hasCorrectFlux = hasMomentumFlux && hasKineticFlux;
203 break;
204 case rnemdKePvector:
205 hasCorrectFlux = hasMomentumFluxVector && hasKineticFlux;
206 break;
207 default:
208 methodFluxMismatch = true;
209 break;
210 }
211 default:
212 break;
213 }
214
215 if (methodFluxMismatch) {
216 sprintf(painCave.errMsg,
217 "RNEMD: The current method,\n"
218 "\t\t%s\n"
219 "\tcannot be used with the current flux type, %s\n",
220 methStr.c_str(), fluxStr.c_str());
221 painCave.isFatal = 1;
222 painCave.severity = OPENMD_ERROR;
223 simError();
224 }
225 if (!hasCorrectFlux) {
226 sprintf(painCave.errMsg,
227 "RNEMD: The current method,\n"
228 "\t%s, and flux type %s\n"
229 "\tdid not have the correct flux value specified. Options\n"
230 "\tinclude: kineticFlux, momentumFlux, and momentumFluxVector\n",
231 methStr.c_str(), fluxStr.c_str());
232 painCave.isFatal = 1;
233 painCave.severity = OPENMD_ERROR;
234 simError();
235 }
236
237 if (hasKineticFlux) {
238 kineticFlux_ = rnemdParams->getKineticFlux();
239 } else {
240 kineticFlux_ = 0.0;
241 }
242 if (hasMomentumFluxVector) {
243 momentumFluxVector_ = rnemdParams->getMomentumFluxVector();
244 } else {
245 momentumFluxVector_ = V3Zero;
246 if (hasMomentumFlux) {
247 RealType momentumFlux = rnemdParams->getMomentumFlux();
248 switch (rnemdFluxType_) {
249 case rnemdPx:
250 momentumFluxVector_.x() = momentumFlux;
251 break;
252 case rnemdPy:
253 momentumFluxVector_.y() = momentumFlux;
254 break;
255 case rnemdPz:
256 momentumFluxVector_.z() = momentumFlux;
257 break;
258 case rnemdKePx:
259 momentumFluxVector_.x() = momentumFlux;
260 break;
261 case rnemdKePy:
262 momentumFluxVector_.y() = momentumFlux;
263 break;
264 default:
265 break;
266 }
267 }
268 }
269
270 // do some sanity checking
271
272 int selectionCount = seleMan_.getSelectionCount();
273 int nIntegrable = info->getNGlobalIntegrableObjects();
274
275 if (selectionCount > nIntegrable) {
276 sprintf(painCave.errMsg,
277 "RNEMD: The current objectSelection,\n"
278 "\t\t%s\n"
279 "\thas resulted in %d selected objects. However,\n"
280 "\tthe total number of integrable objects in the system\n"
281 "\tis only %d. This is almost certainly not what you want\n"
282 "\tto do. A likely cause of this is forgetting the _RB_0\n"
283 "\tselector in the selection script!\n",
284 rnemdObjectSelection_.c_str(),
285 selectionCount, nIntegrable);
286 painCave.isFatal = 0;
287 painCave.severity = OPENMD_WARNING;
288 simError();
289 }
290
291 areaAccumulator_ = new Accumulator();
292
293 nBins_ = rnemdParams->getOutputBins();
294
295 data_.resize(RNEMD::ENDINDEX);
296 OutputData z;
297 z.units = "Angstroms";
298 z.title = "Z";
299 z.dataType = "RealType";
300 z.accumulator.reserve(nBins_);
301 for (unsigned int i = 0; i < nBins_; i++)
302 z.accumulator.push_back( new Accumulator() );
303 data_[Z] = z;
304 outputMap_["Z"] = Z;
305
306 OutputData temperature;
307 temperature.units = "K";
308 temperature.title = "Temperature";
309 temperature.dataType = "RealType";
310 temperature.accumulator.reserve(nBins_);
311 for (unsigned int i = 0; i < nBins_; i++)
312 temperature.accumulator.push_back( new Accumulator() );
313 data_[TEMPERATURE] = temperature;
314 outputMap_["TEMPERATURE"] = TEMPERATURE;
315
316 OutputData velocity;
317 velocity.units = "amu/fs";
318 velocity.title = "Velocity";
319 velocity.dataType = "Vector3d";
320 velocity.accumulator.reserve(nBins_);
321 for (unsigned int i = 0; i < nBins_; i++)
322 velocity.accumulator.push_back( new VectorAccumulator() );
323 data_[VELOCITY] = velocity;
324 outputMap_["VELOCITY"] = VELOCITY;
325
326 OutputData density;
327 density.units = "g cm^-3";
328 density.title = "Density";
329 density.dataType = "RealType";
330 density.accumulator.reserve(nBins_);
331 for (unsigned int i = 0; i < nBins_; i++)
332 density.accumulator.push_back( new Accumulator() );
333 data_[DENSITY] = density;
334 outputMap_["DENSITY"] = DENSITY;
335
336 if (hasOutputFields) {
337 parseOutputFileFormat(rnemdParams->getOutputFields());
338 } else {
339 outputMask_.set(Z);
340 switch (rnemdFluxType_) {
341 case rnemdKE:
342 case rnemdRotKE:
343 case rnemdFullKE:
344 outputMask_.set(TEMPERATURE);
345 break;
346 case rnemdPx:
347 case rnemdPy:
348 outputMask_.set(VELOCITY);
349 break;
350 case rnemdPz:
351 case rnemdPvector:
352 outputMask_.set(VELOCITY);
353 outputMask_.set(DENSITY);
354 break;
355 case rnemdKePx:
356 case rnemdKePy:
357 outputMask_.set(TEMPERATURE);
358 outputMask_.set(VELOCITY);
359 break;
360 case rnemdKePvector:
361 outputMask_.set(TEMPERATURE);
362 outputMask_.set(VELOCITY);
363 outputMask_.set(DENSITY);
364 break;
365 default:
366 break;
367 }
368 }
369
370 if (hasOutputFileName) {
371 rnemdFileName_ = rnemdParams->getOutputFileName();
372 } else {
373 rnemdFileName_ = getPrefix(info->getFinalConfigFileName()) + ".rnemd";
374 }
375
376 exchangeTime_ = rnemdParams->getExchangeTime();
377
378 Snapshot* currentSnap_ = info->getSnapshotManager()->getCurrentSnapshot();
379 Mat3x3d hmat = currentSnap_->getHmat();
380
381 // Target exchange quantities (in each exchange) = 2 Lx Ly dt flux
382 // Lx, Ly = box dimensions in x & y
383 // dt = exchange time interval
384 // flux = target flux
385
386 RealType area = currentSnap_->getXYarea();
387 kineticTarget_ = 2.0 * kineticFlux_ * exchangeTime_ * area;
388 momentumTarget_ = 2.0 * momentumFluxVector_ * exchangeTime_ * area;
389
390 // total exchange sums are zeroed out at the beginning:
391
392 kineticExchange_ = 0.0;
393 momentumExchange_ = V3Zero;
394
395 if (hasSlabWidth)
396 slabWidth_ = rnemdParams->getSlabWidth();
397 else
398 slabWidth_ = hmat(2,2) / 10.0;
399
400 if (hasSlabACenter)
401 slabACenter_ = rnemdParams->getSlabACenter();
402 else
403 slabACenter_ = 0.0;
404
405 if (hasSlabBCenter)
406 slabBCenter_ = rnemdParams->getSlabBCenter();
407 else
408 slabBCenter_ = hmat(2,2) / 2.0;
409
410 }
411
412 RNEMD::~RNEMD() {
413
414 #ifdef IS_MPI
415 if (worldRank == 0) {
416 #endif
417
418 writeOutputFile();
419
420 rnemdFile_.close();
421
422 #ifdef IS_MPI
423 }
424 #endif
425 }
426
427 bool RNEMD::inSlabA(Vector3d pos) {
428 return (abs(pos.z() - slabACenter_) < 0.5*slabWidth_);
429 }
430 bool RNEMD::inSlabB(Vector3d pos) {
431 return (abs(pos.z() - slabBCenter_) < 0.5*slabWidth_);
432 }
433
434 void RNEMD::doSwap() {
435
436 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
437 Mat3x3d hmat = currentSnap_->getHmat();
438
439 seleMan_.setSelectionSet(evaluator_.evaluate());
440
441 int selei;
442 StuntDouble* sd;
443 int idx;
444
445 RealType min_val;
446 bool min_found = false;
447 StuntDouble* min_sd;
448
449 RealType max_val;
450 bool max_found = false;
451 StuntDouble* max_sd;
452
453 for (sd = seleMan_.beginSelected(selei); sd != NULL;
454 sd = seleMan_.nextSelected(selei)) {
455
456 idx = sd->getLocalIndex();
457
458 Vector3d pos = sd->getPos();
459
460 // wrap the stuntdouble's position back into the box:
461
462 if (usePeriodicBoundaryConditions_)
463 currentSnap_->wrapVector(pos);
464 bool inA = inSlabA(pos);
465 bool inB = inSlabB(pos);
466
467 if (inA || inB) {
468
469 RealType mass = sd->getMass();
470 Vector3d vel = sd->getVel();
471 RealType value;
472
473 switch(rnemdFluxType_) {
474 case rnemdKE :
475
476 value = mass * vel.lengthSquare();
477
478 if (sd->isDirectional()) {
479 Vector3d angMom = sd->getJ();
480 Mat3x3d I = sd->getI();
481
482 if (sd->isLinear()) {
483 int i = sd->linearAxis();
484 int j = (i + 1) % 3;
485 int k = (i + 2) % 3;
486 value += angMom[j] * angMom[j] / I(j, j) +
487 angMom[k] * angMom[k] / I(k, k);
488 } else {
489 value += angMom[0]*angMom[0]/I(0, 0)
490 + angMom[1]*angMom[1]/I(1, 1)
491 + angMom[2]*angMom[2]/I(2, 2);
492 }
493 } //angular momenta exchange enabled
494 //energyConvert temporarily disabled
495 //make kineticExchange_ comparable between swap & scale
496 //value = value * 0.5 / PhysicalConstants::energyConvert;
497 value *= 0.5;
498 break;
499 case rnemdPx :
500 value = mass * vel[0];
501 break;
502 case rnemdPy :
503 value = mass * vel[1];
504 break;
505 case rnemdPz :
506 value = mass * vel[2];
507 break;
508 default :
509 break;
510 }
511
512 if (inA == 0) {
513 if (!min_found) {
514 min_val = value;
515 min_sd = sd;
516 min_found = true;
517 } else {
518 if (min_val > value) {
519 min_val = value;
520 min_sd = sd;
521 }
522 }
523 } else {
524 if (!max_found) {
525 max_val = value;
526 max_sd = sd;
527 max_found = true;
528 } else {
529 if (max_val < value) {
530 max_val = value;
531 max_sd = sd;
532 }
533 }
534 }
535 }
536 }
537
538 #ifdef IS_MPI
539 int nProc, worldRank;
540
541 nProc = MPI::COMM_WORLD.Get_size();
542 worldRank = MPI::COMM_WORLD.Get_rank();
543
544 bool my_min_found = min_found;
545 bool my_max_found = max_found;
546
547 // Even if we didn't find a minimum, did someone else?
548 MPI::COMM_WORLD.Allreduce(&my_min_found, &min_found, 1, MPI::BOOL, MPI::LOR);
549 // Even if we didn't find a maximum, did someone else?
550 MPI::COMM_WORLD.Allreduce(&my_max_found, &max_found, 1, MPI::BOOL, MPI::LOR);
551 #endif
552
553 if (max_found && min_found) {
554
555 #ifdef IS_MPI
556 struct {
557 RealType val;
558 int rank;
559 } max_vals, min_vals;
560
561 if (my_min_found) {
562 min_vals.val = min_val;
563 } else {
564 min_vals.val = HONKING_LARGE_VALUE;
565 }
566 min_vals.rank = worldRank;
567
568 // Who had the minimum?
569 MPI::COMM_WORLD.Allreduce(&min_vals, &min_vals,
570 1, MPI::REALTYPE_INT, MPI::MINLOC);
571 min_val = min_vals.val;
572
573 if (my_max_found) {
574 max_vals.val = max_val;
575 } else {
576 max_vals.val = -HONKING_LARGE_VALUE;
577 }
578 max_vals.rank = worldRank;
579
580 // Who had the maximum?
581 MPI::COMM_WORLD.Allreduce(&max_vals, &max_vals,
582 1, MPI::REALTYPE_INT, MPI::MAXLOC);
583 max_val = max_vals.val;
584 #endif
585
586 if (min_val < max_val) {
587
588 #ifdef IS_MPI
589 if (max_vals.rank == worldRank && min_vals.rank == worldRank) {
590 // I have both maximum and minimum, so proceed like a single
591 // processor version:
592 #endif
593
594 Vector3d min_vel = min_sd->getVel();
595 Vector3d max_vel = max_sd->getVel();
596 RealType temp_vel;
597
598 switch(rnemdFluxType_) {
599 case rnemdKE :
600 min_sd->setVel(max_vel);
601 max_sd->setVel(min_vel);
602 if (min_sd->isDirectional() && max_sd->isDirectional()) {
603 Vector3d min_angMom = min_sd->getJ();
604 Vector3d max_angMom = max_sd->getJ();
605 min_sd->setJ(max_angMom);
606 max_sd->setJ(min_angMom);
607 }//angular momenta exchange enabled
608 //assumes same rigid body identity
609 break;
610 case rnemdPx :
611 temp_vel = min_vel.x();
612 min_vel.x() = max_vel.x();
613 max_vel.x() = temp_vel;
614 min_sd->setVel(min_vel);
615 max_sd->setVel(max_vel);
616 break;
617 case rnemdPy :
618 temp_vel = min_vel.y();
619 min_vel.y() = max_vel.y();
620 max_vel.y() = temp_vel;
621 min_sd->setVel(min_vel);
622 max_sd->setVel(max_vel);
623 break;
624 case rnemdPz :
625 temp_vel = min_vel.z();
626 min_vel.z() = max_vel.z();
627 max_vel.z() = temp_vel;
628 min_sd->setVel(min_vel);
629 max_sd->setVel(max_vel);
630 break;
631 default :
632 break;
633 }
634
635 #ifdef IS_MPI
636 // the rest of the cases only apply in parallel simulations:
637 } else if (max_vals.rank == worldRank) {
638 // I had the max, but not the minimum
639
640 Vector3d min_vel;
641 Vector3d max_vel = max_sd->getVel();
642 MPI::Status status;
643
644 // point-to-point swap of the velocity vector
645 MPI::COMM_WORLD.Sendrecv(max_vel.getArrayPointer(), 3, MPI::REALTYPE,
646 min_vals.rank, 0,
647 min_vel.getArrayPointer(), 3, MPI::REALTYPE,
648 min_vals.rank, 0, status);
649
650 switch(rnemdFluxType_) {
651 case rnemdKE :
652 max_sd->setVel(min_vel);
653 //angular momenta exchange enabled
654 if (max_sd->isDirectional()) {
655 Vector3d min_angMom;
656 Vector3d max_angMom = max_sd->getJ();
657
658 // point-to-point swap of the angular momentum vector
659 MPI::COMM_WORLD.Sendrecv(max_angMom.getArrayPointer(), 3,
660 MPI::REALTYPE, min_vals.rank, 1,
661 min_angMom.getArrayPointer(), 3,
662 MPI::REALTYPE, min_vals.rank, 1,
663 status);
664
665 max_sd->setJ(min_angMom);
666 }
667 break;
668 case rnemdPx :
669 max_vel.x() = min_vel.x();
670 max_sd->setVel(max_vel);
671 break;
672 case rnemdPy :
673 max_vel.y() = min_vel.y();
674 max_sd->setVel(max_vel);
675 break;
676 case rnemdPz :
677 max_vel.z() = min_vel.z();
678 max_sd->setVel(max_vel);
679 break;
680 default :
681 break;
682 }
683 } else if (min_vals.rank == worldRank) {
684 // I had the minimum but not the maximum:
685
686 Vector3d max_vel;
687 Vector3d min_vel = min_sd->getVel();
688 MPI::Status status;
689
690 // point-to-point swap of the velocity vector
691 MPI::COMM_WORLD.Sendrecv(min_vel.getArrayPointer(), 3, MPI::REALTYPE,
692 max_vals.rank, 0,
693 max_vel.getArrayPointer(), 3, MPI::REALTYPE,
694 max_vals.rank, 0, status);
695
696 switch(rnemdFluxType_) {
697 case rnemdKE :
698 min_sd->setVel(max_vel);
699 //angular momenta exchange enabled
700 if (min_sd->isDirectional()) {
701 Vector3d min_angMom = min_sd->getJ();
702 Vector3d max_angMom;
703
704 // point-to-point swap of the angular momentum vector
705 MPI::COMM_WORLD.Sendrecv(min_angMom.getArrayPointer(), 3,
706 MPI::REALTYPE, max_vals.rank, 1,
707 max_angMom.getArrayPointer(), 3,
708 MPI::REALTYPE, max_vals.rank, 1,
709 status);
710
711 min_sd->setJ(max_angMom);
712 }
713 break;
714 case rnemdPx :
715 min_vel.x() = max_vel.x();
716 min_sd->setVel(min_vel);
717 break;
718 case rnemdPy :
719 min_vel.y() = max_vel.y();
720 min_sd->setVel(min_vel);
721 break;
722 case rnemdPz :
723 min_vel.z() = max_vel.z();
724 min_sd->setVel(min_vel);
725 break;
726 default :
727 break;
728 }
729 }
730 #endif
731
732 switch(rnemdFluxType_) {
733 case rnemdKE:
734 cerr << "KE\n";
735 kineticExchange_ += max_val - min_val;
736 break;
737 case rnemdPx:
738 momentumExchange_.x() += max_val - min_val;
739 break;
740 case rnemdPy:
741 momentumExchange_.y() += max_val - min_val;
742 break;
743 case rnemdPz:
744 momentumExchange_.z() += max_val - min_val;
745 break;
746 default:
747 cerr << "default\n";
748 break;
749 }
750 } else {
751 sprintf(painCave.errMsg,
752 "RNEMD::doSwap exchange NOT performed because min_val > max_val\n");
753 painCave.isFatal = 0;
754 painCave.severity = OPENMD_INFO;
755 simError();
756 failTrialCount_++;
757 }
758 } else {
759 sprintf(painCave.errMsg,
760 "RNEMD::doSwap exchange NOT performed because selected object\n"
761 "\twas not present in at least one of the two slabs.\n");
762 painCave.isFatal = 0;
763 painCave.severity = OPENMD_INFO;
764 simError();
765 failTrialCount_++;
766 }
767 }
768
769 void RNEMD::doNIVS() {
770
771 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
772 Mat3x3d hmat = currentSnap_->getHmat();
773
774 seleMan_.setSelectionSet(evaluator_.evaluate());
775
776 int selei;
777 StuntDouble* sd;
778 int idx;
779
780 vector<StuntDouble*> hotBin, coldBin;
781
782 RealType Phx = 0.0;
783 RealType Phy = 0.0;
784 RealType Phz = 0.0;
785 RealType Khx = 0.0;
786 RealType Khy = 0.0;
787 RealType Khz = 0.0;
788 RealType Khw = 0.0;
789 RealType Pcx = 0.0;
790 RealType Pcy = 0.0;
791 RealType Pcz = 0.0;
792 RealType Kcx = 0.0;
793 RealType Kcy = 0.0;
794 RealType Kcz = 0.0;
795 RealType Kcw = 0.0;
796
797 for (sd = seleMan_.beginSelected(selei); sd != NULL;
798 sd = seleMan_.nextSelected(selei)) {
799
800 idx = sd->getLocalIndex();
801
802 Vector3d pos = sd->getPos();
803
804 // wrap the stuntdouble's position back into the box:
805
806 if (usePeriodicBoundaryConditions_)
807 currentSnap_->wrapVector(pos);
808
809 // which bin is this stuntdouble in?
810 bool inA = inSlabA(pos);
811 bool inB = inSlabB(pos);
812
813 if (inA || inB) {
814
815 RealType mass = sd->getMass();
816 Vector3d vel = sd->getVel();
817
818 if (inA) {
819 hotBin.push_back(sd);
820 Phx += mass * vel.x();
821 Phy += mass * vel.y();
822 Phz += mass * vel.z();
823 Khx += mass * vel.x() * vel.x();
824 Khy += mass * vel.y() * vel.y();
825 Khz += mass * vel.z() * vel.z();
826 if (sd->isDirectional()) {
827 Vector3d angMom = sd->getJ();
828 Mat3x3d I = sd->getI();
829 if (sd->isLinear()) {
830 int i = sd->linearAxis();
831 int j = (i + 1) % 3;
832 int k = (i + 2) % 3;
833 Khw += angMom[j] * angMom[j] / I(j, j) +
834 angMom[k] * angMom[k] / I(k, k);
835 } else {
836 Khw += angMom[0]*angMom[0]/I(0, 0)
837 + angMom[1]*angMom[1]/I(1, 1)
838 + angMom[2]*angMom[2]/I(2, 2);
839 }
840 }
841 } else {
842 coldBin.push_back(sd);
843 Pcx += mass * vel.x();
844 Pcy += mass * vel.y();
845 Pcz += mass * vel.z();
846 Kcx += mass * vel.x() * vel.x();
847 Kcy += mass * vel.y() * vel.y();
848 Kcz += mass * vel.z() * vel.z();
849 if (sd->isDirectional()) {
850 Vector3d angMom = sd->getJ();
851 Mat3x3d I = sd->getI();
852 if (sd->isLinear()) {
853 int i = sd->linearAxis();
854 int j = (i + 1) % 3;
855 int k = (i + 2) % 3;
856 Kcw += angMom[j] * angMom[j] / I(j, j) +
857 angMom[k] * angMom[k] / I(k, k);
858 } else {
859 Kcw += angMom[0]*angMom[0]/I(0, 0)
860 + angMom[1]*angMom[1]/I(1, 1)
861 + angMom[2]*angMom[2]/I(2, 2);
862 }
863 }
864 }
865 }
866 }
867
868 Khx *= 0.5;
869 Khy *= 0.5;
870 Khz *= 0.5;
871 Khw *= 0.5;
872 Kcx *= 0.5;
873 Kcy *= 0.5;
874 Kcz *= 0.5;
875 Kcw *= 0.5;
876
877 #ifdef IS_MPI
878 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phx, 1, MPI::REALTYPE, MPI::SUM);
879 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phy, 1, MPI::REALTYPE, MPI::SUM);
880 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phz, 1, MPI::REALTYPE, MPI::SUM);
881 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcx, 1, MPI::REALTYPE, MPI::SUM);
882 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcy, 1, MPI::REALTYPE, MPI::SUM);
883 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcz, 1, MPI::REALTYPE, MPI::SUM);
884
885 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khx, 1, MPI::REALTYPE, MPI::SUM);
886 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khy, 1, MPI::REALTYPE, MPI::SUM);
887 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khz, 1, MPI::REALTYPE, MPI::SUM);
888 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khw, 1, MPI::REALTYPE, MPI::SUM);
889
890 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcx, 1, MPI::REALTYPE, MPI::SUM);
891 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcy, 1, MPI::REALTYPE, MPI::SUM);
892 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcz, 1, MPI::REALTYPE, MPI::SUM);
893 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcw, 1, MPI::REALTYPE, MPI::SUM);
894 #endif
895
896 //solve coldBin coeff's first
897 RealType px = Pcx / Phx;
898 RealType py = Pcy / Phy;
899 RealType pz = Pcz / Phz;
900 RealType c, x, y, z;
901 bool successfulScale = false;
902 if ((rnemdFluxType_ == rnemdFullKE) ||
903 (rnemdFluxType_ == rnemdRotKE)) {
904 //may need sanity check Khw & Kcw > 0
905
906 if (rnemdFluxType_ == rnemdFullKE) {
907 c = 1.0 - kineticTarget_ / (Kcx + Kcy + Kcz + Kcw);
908 } else {
909 c = 1.0 - kineticTarget_ / Kcw;
910 }
911
912 if ((c > 0.81) && (c < 1.21)) {//restrict scaling coefficients
913 c = sqrt(c);
914 //std::cerr << "cold slab scaling coefficient: " << c << endl;
915 //now convert to hotBin coefficient
916 RealType w = 0.0;
917 if (rnemdFluxType_ == rnemdFullKE) {
918 x = 1.0 + px * (1.0 - c);
919 y = 1.0 + py * (1.0 - c);
920 z = 1.0 + pz * (1.0 - c);
921 /* more complicated way
922 w = 1.0 + (Kcw - Kcw * c * c - (c * c * (Kcx + Kcy + Kcz
923 + Khx * px * px + Khy * py * py + Khz * pz * pz)
924 - 2.0 * c * (Khx * px * (1.0 + px) + Khy * py * (1.0 + py)
925 + Khz * pz * (1.0 + pz)) + Khx * px * (2.0 + px)
926 + Khy * py * (2.0 + py) + Khz * pz * (2.0 + pz)
927 - Kcx - Kcy - Kcz)) / Khw; the following is simpler
928 */
929 if ((fabs(x - 1.0) < 0.1) && (fabs(y - 1.0) < 0.1) &&
930 (fabs(z - 1.0) < 0.1)) {
931 w = 1.0 + (kineticTarget_
932 + Khx * (1.0 - x * x) + Khy * (1.0 - y * y)
933 + Khz * (1.0 - z * z)) / Khw;
934 }//no need to calculate w if x, y or z is out of range
935 } else {
936 w = 1.0 + kineticTarget_ / Khw;
937 }
938 if ((w > 0.81) && (w < 1.21)) {//restrict scaling coefficients
939 //if w is in the right range, so should be x, y, z.
940 vector<StuntDouble*>::iterator sdi;
941 Vector3d vel;
942 for (sdi = coldBin.begin(); sdi != coldBin.end(); sdi++) {
943 if (rnemdFluxType_ == rnemdFullKE) {
944 vel = (*sdi)->getVel() * c;
945 (*sdi)->setVel(vel);
946 }
947 if ((*sdi)->isDirectional()) {
948 Vector3d angMom = (*sdi)->getJ() * c;
949 (*sdi)->setJ(angMom);
950 }
951 }
952 w = sqrt(w);
953 // std::cerr << "xh= " << x << "\tyh= " << y << "\tzh= " << z
954 // << "\twh= " << w << endl;
955 for (sdi = hotBin.begin(); sdi != hotBin.end(); sdi++) {
956 if (rnemdFluxType_ == rnemdFullKE) {
957 vel = (*sdi)->getVel();
958 vel.x() *= x;
959 vel.y() *= y;
960 vel.z() *= z;
961 (*sdi)->setVel(vel);
962 }
963 if ((*sdi)->isDirectional()) {
964 Vector3d angMom = (*sdi)->getJ() * w;
965 (*sdi)->setJ(angMom);
966 }
967 }
968 successfulScale = true;
969 kineticExchange_ += kineticTarget_;
970 }
971 }
972 } else {
973 RealType a000, a110, c0, a001, a111, b01, b11, c1;
974 switch(rnemdFluxType_) {
975 case rnemdKE :
976 /* used hotBin coeff's & only scale x & y dimensions
977 RealType px = Phx / Pcx;
978 RealType py = Phy / Pcy;
979 a110 = Khy;
980 c0 = - Khx - Khy - kineticTarget_;
981 a000 = Khx;
982 a111 = Kcy * py * py;
983 b11 = -2.0 * Kcy * py * (1.0 + py);
984 c1 = Kcy * py * (2.0 + py) + Kcx * px * ( 2.0 + px) + kineticTarget_;
985 b01 = -2.0 * Kcx * px * (1.0 + px);
986 a001 = Kcx * px * px;
987 */
988 //scale all three dimensions, let c_x = c_y
989 a000 = Kcx + Kcy;
990 a110 = Kcz;
991 c0 = kineticTarget_ - Kcx - Kcy - Kcz;
992 a001 = Khx * px * px + Khy * py * py;
993 a111 = Khz * pz * pz;
994 b01 = -2.0 * (Khx * px * (1.0 + px) + Khy * py * (1.0 + py));
995 b11 = -2.0 * Khz * pz * (1.0 + pz);
996 c1 = Khx * px * (2.0 + px) + Khy * py * (2.0 + py)
997 + Khz * pz * (2.0 + pz) - kineticTarget_;
998 break;
999 case rnemdPx :
1000 c = 1 - momentumTarget_.x() / Pcx;
1001 a000 = Kcy;
1002 a110 = Kcz;
1003 c0 = Kcx * c * c - Kcx - Kcy - Kcz;
1004 a001 = py * py * Khy;
1005 a111 = pz * pz * Khz;
1006 b01 = -2.0 * Khy * py * (1.0 + py);
1007 b11 = -2.0 * Khz * pz * (1.0 + pz);
1008 c1 = Khy * py * (2.0 + py) + Khz * pz * (2.0 + pz)
1009 + Khx * (fastpow(c * px - px - 1.0, 2) - 1.0);
1010 break;
1011 case rnemdPy :
1012 c = 1 - momentumTarget_.y() / Pcy;
1013 a000 = Kcx;
1014 a110 = Kcz;
1015 c0 = Kcy * c * c - Kcx - Kcy - Kcz;
1016 a001 = px * px * Khx;
1017 a111 = pz * pz * Khz;
1018 b01 = -2.0 * Khx * px * (1.0 + px);
1019 b11 = -2.0 * Khz * pz * (1.0 + pz);
1020 c1 = Khx * px * (2.0 + px) + Khz * pz * (2.0 + pz)
1021 + Khy * (fastpow(c * py - py - 1.0, 2) - 1.0);
1022 break;
1023 case rnemdPz ://we don't really do this, do we?
1024 c = 1 - momentumTarget_.z() / Pcz;
1025 a000 = Kcx;
1026 a110 = Kcy;
1027 c0 = Kcz * c * c - Kcx - Kcy - Kcz;
1028 a001 = px * px * Khx;
1029 a111 = py * py * Khy;
1030 b01 = -2.0 * Khx * px * (1.0 + px);
1031 b11 = -2.0 * Khy * py * (1.0 + py);
1032 c1 = Khx * px * (2.0 + px) + Khy * py * (2.0 + py)
1033 + Khz * (fastpow(c * pz - pz - 1.0, 2) - 1.0);
1034 break;
1035 default :
1036 break;
1037 }
1038
1039 RealType v1 = a000 * a111 - a001 * a110;
1040 RealType v2 = a000 * b01;
1041 RealType v3 = a000 * b11;
1042 RealType v4 = a000 * c1 - a001 * c0;
1043 RealType v8 = a110 * b01;
1044 RealType v10 = - b01 * c0;
1045
1046 RealType u0 = v2 * v10 - v4 * v4;
1047 RealType u1 = -2.0 * v3 * v4;
1048 RealType u2 = -v2 * v8 - v3 * v3 - 2.0 * v1 * v4;
1049 RealType u3 = -2.0 * v1 * v3;
1050 RealType u4 = - v1 * v1;
1051 //rescale coefficients
1052 RealType maxAbs = fabs(u0);
1053 if (maxAbs < fabs(u1)) maxAbs = fabs(u1);
1054 if (maxAbs < fabs(u2)) maxAbs = fabs(u2);
1055 if (maxAbs < fabs(u3)) maxAbs = fabs(u3);
1056 if (maxAbs < fabs(u4)) maxAbs = fabs(u4);
1057 u0 /= maxAbs;
1058 u1 /= maxAbs;
1059 u2 /= maxAbs;
1060 u3 /= maxAbs;
1061 u4 /= maxAbs;
1062 //max_element(start, end) is also available.
1063 Polynomial<RealType> poly; //same as DoublePolynomial poly;
1064 poly.setCoefficient(4, u4);
1065 poly.setCoefficient(3, u3);
1066 poly.setCoefficient(2, u2);
1067 poly.setCoefficient(1, u1);
1068 poly.setCoefficient(0, u0);
1069 vector<RealType> realRoots = poly.FindRealRoots();
1070
1071 vector<RealType>::iterator ri;
1072 RealType r1, r2, alpha0;
1073 vector<pair<RealType,RealType> > rps;
1074 for (ri = realRoots.begin(); ri !=realRoots.end(); ri++) {
1075 r2 = *ri;
1076 //check if FindRealRoots() give the right answer
1077 if ( fabs(u0 + r2 * (u1 + r2 * (u2 + r2 * (u3 + r2 * u4)))) > 1e-6 ) {
1078 sprintf(painCave.errMsg,
1079 "RNEMD Warning: polynomial solve seems to have an error!");
1080 painCave.isFatal = 0;
1081 simError();
1082 failRootCount_++;
1083 }
1084 //might not be useful w/o rescaling coefficients
1085 alpha0 = -c0 - a110 * r2 * r2;
1086 if (alpha0 >= 0.0) {
1087 r1 = sqrt(alpha0 / a000);
1088 if (fabs(c1 + r1 * (b01 + r1 * a001) + r2 * (b11 + r2 * a111))
1089 < 1e-6)
1090 { rps.push_back(make_pair(r1, r2)); }
1091 if (r1 > 1e-6) { //r1 non-negative
1092 r1 = -r1;
1093 if (fabs(c1 + r1 * (b01 + r1 * a001) + r2 * (b11 + r2 * a111))
1094 < 1e-6)
1095 { rps.push_back(make_pair(r1, r2)); }
1096 }
1097 }
1098 }
1099 // Consider combining together the solving pair part w/ the searching
1100 // best solution part so that we don't need the pairs vector
1101 if (!rps.empty()) {
1102 RealType smallestDiff = HONKING_LARGE_VALUE;
1103 RealType diff;
1104 pair<RealType,RealType> bestPair = make_pair(1.0, 1.0);
1105 vector<pair<RealType,RealType> >::iterator rpi;
1106 for (rpi = rps.begin(); rpi != rps.end(); rpi++) {
1107 r1 = (*rpi).first;
1108 r2 = (*rpi).second;
1109 switch(rnemdFluxType_) {
1110 case rnemdKE :
1111 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1112 + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcx, 2)
1113 + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcy, 2);
1114 break;
1115 case rnemdPx :
1116 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1117 + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcy, 2);
1118 break;
1119 case rnemdPy :
1120 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1121 + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcx, 2);
1122 break;
1123 case rnemdPz :
1124 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1125 + fastpow(r1 * r1 / r2 / r2 - Kcy/Kcx, 2);
1126 default :
1127 break;
1128 }
1129 if (diff < smallestDiff) {
1130 smallestDiff = diff;
1131 bestPair = *rpi;
1132 }
1133 }
1134 #ifdef IS_MPI
1135 if (worldRank == 0) {
1136 #endif
1137 // sprintf(painCave.errMsg,
1138 // "RNEMD: roots r1= %lf\tr2 = %lf\n",
1139 // bestPair.first, bestPair.second);
1140 // painCave.isFatal = 0;
1141 // painCave.severity = OPENMD_INFO;
1142 // simError();
1143 #ifdef IS_MPI
1144 }
1145 #endif
1146
1147 switch(rnemdFluxType_) {
1148 case rnemdKE :
1149 x = bestPair.first;
1150 y = bestPair.first;
1151 z = bestPair.second;
1152 break;
1153 case rnemdPx :
1154 x = c;
1155 y = bestPair.first;
1156 z = bestPair.second;
1157 break;
1158 case rnemdPy :
1159 x = bestPair.first;
1160 y = c;
1161 z = bestPair.second;
1162 break;
1163 case rnemdPz :
1164 x = bestPair.first;
1165 y = bestPair.second;
1166 z = c;
1167 break;
1168 default :
1169 break;
1170 }
1171 vector<StuntDouble*>::iterator sdi;
1172 Vector3d vel;
1173 for (sdi = coldBin.begin(); sdi != coldBin.end(); sdi++) {
1174 vel = (*sdi)->getVel();
1175 vel.x() *= x;
1176 vel.y() *= y;
1177 vel.z() *= z;
1178 (*sdi)->setVel(vel);
1179 }
1180 //convert to hotBin coefficient
1181 x = 1.0 + px * (1.0 - x);
1182 y = 1.0 + py * (1.0 - y);
1183 z = 1.0 + pz * (1.0 - z);
1184 for (sdi = hotBin.begin(); sdi != hotBin.end(); sdi++) {
1185 vel = (*sdi)->getVel();
1186 vel.x() *= x;
1187 vel.y() *= y;
1188 vel.z() *= z;
1189 (*sdi)->setVel(vel);
1190 }
1191 successfulScale = true;
1192 switch(rnemdFluxType_) {
1193 case rnemdKE :
1194 kineticExchange_ += kineticTarget_;
1195 break;
1196 case rnemdPx :
1197 case rnemdPy :
1198 case rnemdPz :
1199 momentumExchange_ += momentumTarget_;
1200 break;
1201 default :
1202 break;
1203 }
1204 }
1205 }
1206 if (successfulScale != true) {
1207 sprintf(painCave.errMsg,
1208 "RNEMD::doNIVS exchange NOT performed - roots that solve\n"
1209 "\tthe constraint equations may not exist or there may be\n"
1210 "\tno selected objects in one or both slabs.\n");
1211 painCave.isFatal = 0;
1212 painCave.severity = OPENMD_INFO;
1213 simError();
1214 failTrialCount_++;
1215 }
1216 }
1217
1218 void RNEMD::doVSS() {
1219
1220 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1221 RealType time = currentSnap_->getTime();
1222 Mat3x3d hmat = currentSnap_->getHmat();
1223
1224 seleMan_.setSelectionSet(evaluator_.evaluate());
1225
1226 int selei;
1227 StuntDouble* sd;
1228 int idx;
1229
1230 vector<StuntDouble*> hotBin, coldBin;
1231
1232 Vector3d Ph(V3Zero);
1233 RealType Mh = 0.0;
1234 RealType Kh = 0.0;
1235 Vector3d Pc(V3Zero);
1236 RealType Mc = 0.0;
1237 RealType Kc = 0.0;
1238
1239
1240 for (sd = seleMan_.beginSelected(selei); sd != NULL;
1241 sd = seleMan_.nextSelected(selei)) {
1242
1243 idx = sd->getLocalIndex();
1244
1245 Vector3d pos = sd->getPos();
1246
1247 // wrap the stuntdouble's position back into the box:
1248
1249 if (usePeriodicBoundaryConditions_)
1250 currentSnap_->wrapVector(pos);
1251
1252 // which bin is this stuntdouble in?
1253 bool inA = inSlabA(pos);
1254 bool inB = inSlabB(pos);
1255
1256 if (inA || inB) {
1257
1258 RealType mass = sd->getMass();
1259 Vector3d vel = sd->getVel();
1260
1261 if (inA) {
1262 hotBin.push_back(sd);
1263 //std::cerr << "before, velocity = " << vel << endl;
1264 Ph += mass * vel;
1265 //std::cerr << "after, velocity = " << vel << endl;
1266 Mh += mass;
1267 Kh += mass * vel.lengthSquare();
1268 if (rnemdFluxType_ == rnemdFullKE) {
1269 if (sd->isDirectional()) {
1270 Vector3d angMom = sd->getJ();
1271 Mat3x3d I = sd->getI();
1272 if (sd->isLinear()) {
1273 int i = sd->linearAxis();
1274 int j = (i + 1) % 3;
1275 int k = (i + 2) % 3;
1276 Kh += angMom[j] * angMom[j] / I(j, j) +
1277 angMom[k] * angMom[k] / I(k, k);
1278 } else {
1279 Kh += angMom[0] * angMom[0] / I(0, 0) +
1280 angMom[1] * angMom[1] / I(1, 1) +
1281 angMom[2] * angMom[2] / I(2, 2);
1282 }
1283 }
1284 }
1285 } else { //midBin_
1286 coldBin.push_back(sd);
1287 Pc += mass * vel;
1288 Mc += mass;
1289 Kc += mass * vel.lengthSquare();
1290 if (rnemdFluxType_ == rnemdFullKE) {
1291 if (sd->isDirectional()) {
1292 Vector3d angMom = sd->getJ();
1293 Mat3x3d I = sd->getI();
1294 if (sd->isLinear()) {
1295 int i = sd->linearAxis();
1296 int j = (i + 1) % 3;
1297 int k = (i + 2) % 3;
1298 Kc += angMom[j] * angMom[j] / I(j, j) +
1299 angMom[k] * angMom[k] / I(k, k);
1300 } else {
1301 Kc += angMom[0] * angMom[0] / I(0, 0) +
1302 angMom[1] * angMom[1] / I(1, 1) +
1303 angMom[2] * angMom[2] / I(2, 2);
1304 }
1305 }
1306 }
1307 }
1308 }
1309 }
1310
1311 Kh *= 0.5;
1312 Kc *= 0.5;
1313
1314 // std::cerr << "Mh= " << Mh << "\tKh= " << Kh << "\tMc= " << Mc
1315 // << "\tKc= " << Kc << endl;
1316 // std::cerr << "Ph= " << Ph << "\tPc= " << Pc << endl;
1317
1318 #ifdef IS_MPI
1319 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Ph[0], 3, MPI::REALTYPE, MPI::SUM);
1320 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pc[0], 3, MPI::REALTYPE, MPI::SUM);
1321 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Mh, 1, MPI::REALTYPE, MPI::SUM);
1322 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kh, 1, MPI::REALTYPE, MPI::SUM);
1323 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Mc, 1, MPI::REALTYPE, MPI::SUM);
1324 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kc, 1, MPI::REALTYPE, MPI::SUM);
1325 #endif
1326
1327 bool successfulExchange = false;
1328 if ((Mh > 0.0) && (Mc > 0.0)) {//both slabs are not empty
1329 Vector3d vc = Pc / Mc;
1330 Vector3d ac = -momentumTarget_ / Mc + vc;
1331 Vector3d acrec = -momentumTarget_ / Mc;
1332 RealType cNumerator = Kc - kineticTarget_ - 0.5 * Mc * ac.lengthSquare();
1333 if (cNumerator > 0.0) {
1334 RealType cDenominator = Kc - 0.5 * Mc * vc.lengthSquare();
1335 if (cDenominator > 0.0) {
1336 RealType c = sqrt(cNumerator / cDenominator);
1337 if ((c > 0.9) && (c < 1.1)) {//restrict scaling coefficients
1338 Vector3d vh = Ph / Mh;
1339 Vector3d ah = momentumTarget_ / Mh + vh;
1340 Vector3d ahrec = momentumTarget_ / Mh;
1341 RealType hNumerator = Kh + kineticTarget_
1342 - 0.5 * Mh * ah.lengthSquare();
1343 if (hNumerator > 0.0) {
1344 RealType hDenominator = Kh - 0.5 * Mh * vh.lengthSquare();
1345 if (hDenominator > 0.0) {
1346 RealType h = sqrt(hNumerator / hDenominator);
1347 if ((h > 0.9) && (h < 1.1)) {
1348 // std::cerr << "cold slab scaling coefficient: " << c << "\n";
1349 // std::cerr << "hot slab scaling coefficient: " << h << "\n";
1350 vector<StuntDouble*>::iterator sdi;
1351 Vector3d vel;
1352 for (sdi = coldBin.begin(); sdi != coldBin.end(); sdi++) {
1353 //vel = (*sdi)->getVel();
1354 vel = ((*sdi)->getVel() - vc) * c + ac;
1355 (*sdi)->setVel(vel);
1356 if (rnemdFluxType_ == rnemdFullKE) {
1357 if ((*sdi)->isDirectional()) {
1358 Vector3d angMom = (*sdi)->getJ() * c;
1359 (*sdi)->setJ(angMom);
1360 }
1361 }
1362 }
1363 for (sdi = hotBin.begin(); sdi != hotBin.end(); sdi++) {
1364 //vel = (*sdi)->getVel();
1365 vel = ((*sdi)->getVel() - vh) * h + ah;
1366 (*sdi)->setVel(vel);
1367 if (rnemdFluxType_ == rnemdFullKE) {
1368 if ((*sdi)->isDirectional()) {
1369 Vector3d angMom = (*sdi)->getJ() * h;
1370 (*sdi)->setJ(angMom);
1371 }
1372 }
1373 }
1374 successfulExchange = true;
1375 kineticExchange_ += kineticTarget_;
1376 momentumExchange_ += momentumTarget_;
1377 }
1378 }
1379 }
1380 }
1381 }
1382 }
1383 }
1384 if (successfulExchange != true) {
1385 sprintf(painCave.errMsg,
1386 "RNEMD::doVSS exchange NOT performed - roots that solve\n"
1387 "\tthe constraint equations may not exist or there may be\n"
1388 "\tno selected objects in one or both slabs.\n");
1389 painCave.isFatal = 0;
1390 painCave.severity = OPENMD_INFO;
1391 simError();
1392 failTrialCount_++;
1393 }
1394 }
1395
1396 void RNEMD::doRNEMD() {
1397
1398 trialCount_++;
1399 switch(rnemdMethod_) {
1400 case rnemdSwap:
1401 doSwap();
1402 break;
1403 case rnemdNIVS:
1404 doNIVS();
1405 break;
1406 case rnemdVSS:
1407 doVSS();
1408 break;
1409 case rnemdUnkownMethod:
1410 default :
1411 break;
1412 }
1413 }
1414
1415 void RNEMD::collectData() {
1416
1417 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1418 Mat3x3d hmat = currentSnap_->getHmat();
1419
1420 areaAccumulator_->add(currentSnap_->getXYarea());
1421
1422 seleMan_.setSelectionSet(evaluator_.evaluate());
1423
1424 int selei;
1425 StuntDouble* sd;
1426 int idx;
1427
1428 vector<RealType> binMass(nBins_, 0.0);
1429 vector<RealType> binPx(nBins_, 0.0);
1430 vector<RealType> binPy(nBins_, 0.0);
1431 vector<RealType> binPz(nBins_, 0.0);
1432 vector<RealType> binKE(nBins_, 0.0);
1433 vector<int> binDOF(nBins_, 0);
1434 vector<int> binCount(nBins_, 0);
1435
1436 // alternative approach, track all molecules instead of only those
1437 // selected for scaling/swapping:
1438 /*
1439 SimInfo::MoleculeIterator miter;
1440 vector<StuntDouble*>::iterator iiter;
1441 Molecule* mol;
1442 StuntDouble* sd;
1443 for (mol = info_->beginMolecule(miter); mol != NULL;
1444 mol = info_->nextMolecule(miter))
1445 sd is essentially sd
1446 for (sd = mol->beginIntegrableObject(iiter);
1447 sd != NULL;
1448 sd = mol->nextIntegrableObject(iiter))
1449 */
1450 for (sd = seleMan_.beginSelected(selei); sd != NULL;
1451 sd = seleMan_.nextSelected(selei)) {
1452
1453 idx = sd->getLocalIndex();
1454
1455 Vector3d pos = sd->getPos();
1456
1457 // wrap the stuntdouble's position back into the box:
1458
1459 if (usePeriodicBoundaryConditions_)
1460 currentSnap_->wrapVector(pos);
1461
1462
1463 // which bin is this stuntdouble in?
1464 // wrapped positions are in the range [-0.5*hmat(2,2), +0.5*hmat(2,2)]
1465 // Shift molecules by half a box to have bins start at 0
1466 // The modulo operator is used to wrap the case when we are
1467 // beyond the end of the bins back to the beginning.
1468 int binNo = int(nBins_ * (pos.z() / hmat(2,2) + 0.5)) % nBins_;
1469
1470 RealType mass = sd->getMass();
1471 Vector3d vel = sd->getVel();
1472
1473 binCount[binNo]++;
1474 binMass[binNo] += mass;
1475 binPx[binNo] += mass*vel.x();
1476 binPy[binNo] += mass*vel.y();
1477 binPz[binNo] += mass*vel.z();
1478 binKE[binNo] += 0.5 * (mass * vel.lengthSquare());
1479 binDOF[binNo] += 3;
1480
1481 if (sd->isDirectional()) {
1482 Vector3d angMom = sd->getJ();
1483 Mat3x3d I = sd->getI();
1484 if (sd->isLinear()) {
1485 int i = sd->linearAxis();
1486 int j = (i + 1) % 3;
1487 int k = (i + 2) % 3;
1488 binKE[binNo] += 0.5 * (angMom[j] * angMom[j] / I(j, j) +
1489 angMom[k] * angMom[k] / I(k, k));
1490 binDOF[binNo] += 2;
1491 } else {
1492 binKE[binNo] += 0.5 * (angMom[0] * angMom[0] / I(0, 0) +
1493 angMom[1] * angMom[1] / I(1, 1) +
1494 angMom[2] * angMom[2] / I(2, 2));
1495 binDOF[binNo] += 3;
1496 }
1497 }
1498 }
1499
1500
1501 #ifdef IS_MPI
1502 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binCount[0],
1503 nBins_, MPI::INT, MPI::SUM);
1504 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binMass[0],
1505 nBins_, MPI::REALTYPE, MPI::SUM);
1506 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPx[0],
1507 nBins_, MPI::REALTYPE, MPI::SUM);
1508 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPy[0],
1509 nBins_, MPI::REALTYPE, MPI::SUM);
1510 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPz[0],
1511 nBins_, MPI::REALTYPE, MPI::SUM);
1512 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binKE[0],
1513 nBins_, MPI::REALTYPE, MPI::SUM);
1514 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binDOF[0],
1515 nBins_, MPI::INT, MPI::SUM);
1516 #endif
1517
1518 Vector3d vel;
1519 RealType den;
1520 RealType temp;
1521 RealType z;
1522 for (int i = 0; i < nBins_; i++) {
1523 z = (((RealType)i + 0.5) / (RealType)nBins_) * hmat(2,2);
1524 vel.x() = binPx[i] / binMass[i];
1525 vel.y() = binPy[i] / binMass[i];
1526 vel.z() = binPz[i] / binMass[i];
1527 den = binCount[i] * nBins_ / currentSnap_->getVolume();
1528 temp = 2.0 * binKE[i] / (binDOF[i] * PhysicalConstants::kb *
1529 PhysicalConstants::energyConvert);
1530
1531 for (unsigned int j = 0; j < outputMask_.size(); ++j) {
1532 if(outputMask_[j]) {
1533 switch(j) {
1534 case Z:
1535 (data_[j].accumulator[i])->add(z);
1536 break;
1537 case TEMPERATURE:
1538 data_[j].accumulator[i]->add(temp);
1539 break;
1540 case VELOCITY:
1541 dynamic_cast<VectorAccumulator *>(data_[j].accumulator[i])->add(vel);
1542 break;
1543 case DENSITY:
1544 data_[j].accumulator[i]->add(den);
1545 break;
1546 }
1547 }
1548 }
1549 }
1550 }
1551
1552 void RNEMD::getStarted() {
1553 collectData();
1554 writeOutputFile();
1555 }
1556
1557 void RNEMD::parseOutputFileFormat(const std::string& format) {
1558 StringTokenizer tokenizer(format, " ,;|\t\n\r");
1559
1560 while(tokenizer.hasMoreTokens()) {
1561 std::string token(tokenizer.nextToken());
1562 toUpper(token);
1563 OutputMapType::iterator i = outputMap_.find(token);
1564 if (i != outputMap_.end()) {
1565 outputMask_.set(i->second);
1566 } else {
1567 sprintf( painCave.errMsg,
1568 "RNEMD::parseOutputFileFormat: %s is not a recognized\n"
1569 "\toutputFileFormat keyword.\n", token.c_str() );
1570 painCave.isFatal = 0;
1571 painCave.severity = OPENMD_ERROR;
1572 simError();
1573 }
1574 }
1575 }
1576
1577 void RNEMD::writeOutputFile() {
1578
1579 #ifdef IS_MPI
1580 // If we're the root node, should we print out the results
1581 int worldRank = MPI::COMM_WORLD.Get_rank();
1582 if (worldRank == 0) {
1583 #endif
1584 rnemdFile_.open(rnemdFileName_.c_str(), std::ios::out | std::ios::trunc );
1585
1586 if( !rnemdFile_ ){
1587 sprintf( painCave.errMsg,
1588 "Could not open \"%s\" for RNEMD output.\n",
1589 rnemdFileName_.c_str());
1590 painCave.isFatal = 1;
1591 simError();
1592 }
1593
1594 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1595
1596 RealType time = currentSnap_->getTime();
1597 RealType avgArea;
1598 areaAccumulator_->getAverage(avgArea);
1599 RealType Jz = kineticExchange_ / (2.0 * time * avgArea);
1600 Vector3d JzP = momentumExchange_ / (2.0 * time * avgArea);
1601
1602 rnemdFile_ << "#######################################################\n";
1603 rnemdFile_ << "# RNEMD {\n";
1604
1605 map<string, RNEMDMethod>::iterator mi;
1606 for(mi = stringToMethod_.begin(); mi != stringToMethod_.end(); ++mi) {
1607 if ( (*mi).second == rnemdMethod_)
1608 rnemdFile_ << "# exchangeMethod = \"" << (*mi).first << "\";\n";
1609 }
1610 map<string, RNEMDFluxType>::iterator fi;
1611 for(fi = stringToFluxType_.begin(); fi != stringToFluxType_.end(); ++fi) {
1612 if ( (*fi).second == rnemdFluxType_)
1613 rnemdFile_ << "# fluxType = \"" << (*fi).first << "\";\n";
1614 }
1615
1616 rnemdFile_ << "# exchangeTime = " << exchangeTime_ << "; fs\n";
1617
1618 rnemdFile_ << "# objectSelection = \""
1619 << rnemdObjectSelection_ << "\";\n";
1620 rnemdFile_ << "# slabWidth = " << slabWidth_ << ";\n";
1621 rnemdFile_ << "# slabAcenter = " << slabACenter_ << ";\n";
1622 rnemdFile_ << "# slabBcenter = " << slabBCenter_ << ";\n";
1623 rnemdFile_ << "# }\n";
1624 rnemdFile_ << "#######################################################\n";
1625 rnemdFile_ << "# RNEMD report:\n";
1626 rnemdFile_ << "# running time = " << time << " fs\n";
1627 rnemdFile_ << "# target flux:\n";
1628 rnemdFile_ << "# kinetic = " << kineticFlux_ << "\n";
1629 rnemdFile_ << "# momentum = " << momentumFluxVector_ << "\n";
1630 rnemdFile_ << "# target one-time exchanges:\n";
1631 rnemdFile_ << "# kinetic = " << kineticTarget_ << "\n";
1632 rnemdFile_ << "# momentum = " << momentumTarget_ << "\n";
1633 rnemdFile_ << "# actual exchange totals:\n";
1634 rnemdFile_ << "# kinetic = " << kineticExchange_ << "\n";
1635 rnemdFile_ << "# momentum = " << momentumExchange_ << "\n";
1636 rnemdFile_ << "# actual flux:\n";
1637 rnemdFile_ << "# kinetic = " << Jz << "\n";
1638 rnemdFile_ << "# momentum = " << JzP << "\n";
1639 rnemdFile_ << "# exchange statistics:\n";
1640 rnemdFile_ << "# attempted = " << trialCount_ << "\n";
1641 rnemdFile_ << "# failed = " << failTrialCount_ << "\n";
1642 if (rnemdMethod_ == rnemdNIVS) {
1643 rnemdFile_ << "# NIVS root-check errors = "
1644 << failRootCount_ << "\n";
1645 }
1646 rnemdFile_ << "#######################################################\n";
1647
1648
1649
1650 //write title
1651 rnemdFile_ << "#";
1652 for (unsigned int i = 0; i < outputMask_.size(); ++i) {
1653 if (outputMask_[i]) {
1654 rnemdFile_ << "\t" << data_[i].title <<
1655 "(" << data_[i].units << ")";
1656 }
1657 }
1658 rnemdFile_ << std::endl;
1659
1660 rnemdFile_.precision(8);
1661
1662 for (unsigned int j = 0; j < nBins_; j++) {
1663
1664 for (unsigned int i = 0; i < outputMask_.size(); ++i) {
1665 if (outputMask_[i]) {
1666 if (data_[i].dataType == "RealType")
1667 writeReal(i,j);
1668 else if (data_[i].dataType == "Vector3d")
1669 writeVector(i,j);
1670 else {
1671 sprintf( painCave.errMsg,
1672 "RNEMD found an unknown data type for: %s ",
1673 data_[i].title.c_str());
1674 painCave.isFatal = 1;
1675 simError();
1676 }
1677 }
1678 }
1679 rnemdFile_ << std::endl;
1680
1681 }
1682
1683 rnemdFile_ << "#######################################################\n";
1684 rnemdFile_ << "# Standard Deviations in those quantities follow:\n";
1685 rnemdFile_ << "#######################################################\n";
1686
1687
1688 for (unsigned int j = 0; j < nBins_; j++) {
1689 rnemdFile_ << "#";
1690 for (unsigned int i = 0; i < outputMask_.size(); ++i) {
1691 if (outputMask_[i]) {
1692 if (data_[i].dataType == "RealType")
1693 writeRealStdDev(i,j);
1694 else if (data_[i].dataType == "Vector3d")
1695 writeVectorStdDev(i,j);
1696 else {
1697 sprintf( painCave.errMsg,
1698 "RNEMD found an unknown data type for: %s ",
1699 data_[i].title.c_str());
1700 painCave.isFatal = 1;
1701 simError();
1702 }
1703 }
1704 }
1705 rnemdFile_ << std::endl;
1706
1707 }
1708
1709 rnemdFile_.flush();
1710 rnemdFile_.close();
1711
1712 #ifdef IS_MPI
1713 }
1714 #endif
1715
1716 }
1717
1718 void RNEMD::writeReal(int index, unsigned int bin) {
1719 assert(index >=0 && index < ENDINDEX);
1720 assert(bin < nBins_);
1721 RealType s;
1722
1723 data_[index].accumulator[bin]->getAverage(s);
1724
1725 if (! isinf(s) && ! isnan(s)) {
1726 rnemdFile_ << "\t" << s;
1727 } else{
1728 sprintf( painCave.errMsg,
1729 "RNEMD detected a numerical error writing: %s for bin %d",
1730 data_[index].title.c_str(), bin);
1731 painCave.isFatal = 1;
1732 simError();
1733 }
1734 }
1735
1736 void RNEMD::writeVector(int index, unsigned int bin) {
1737 assert(index >=0 && index < ENDINDEX);
1738 assert(bin < nBins_);
1739 Vector3d s;
1740 dynamic_cast<VectorAccumulator*>(data_[index].accumulator[bin])->getAverage(s);
1741 if (isinf(s[0]) || isnan(s[0]) ||
1742 isinf(s[1]) || isnan(s[1]) ||
1743 isinf(s[2]) || isnan(s[2]) ) {
1744 sprintf( painCave.errMsg,
1745 "RNEMD detected a numerical error writing: %s for bin %d",
1746 data_[index].title.c_str(), bin);
1747 painCave.isFatal = 1;
1748 simError();
1749 } else {
1750 rnemdFile_ << "\t" << s[0] << "\t" << s[1] << "\t" << s[2];
1751 }
1752 }
1753
1754 void RNEMD::writeRealStdDev(int index, unsigned int bin) {
1755 assert(index >=0 && index < ENDINDEX);
1756 assert(bin < nBins_);
1757 RealType s;
1758
1759 data_[index].accumulator[bin]->getStdDev(s);
1760
1761 if (! isinf(s) && ! isnan(s)) {
1762 rnemdFile_ << "\t" << s;
1763 } else{
1764 sprintf( painCave.errMsg,
1765 "RNEMD detected a numerical error writing: %s std. dev. for bin %d",
1766 data_[index].title.c_str(), bin);
1767 painCave.isFatal = 1;
1768 simError();
1769 }
1770 }
1771
1772 void RNEMD::writeVectorStdDev(int index, unsigned int bin) {
1773 assert(index >=0 && index < ENDINDEX);
1774 assert(bin < nBins_);
1775 Vector3d s;
1776 dynamic_cast<VectorAccumulator*>(data_[index].accumulator[bin])->getStdDev(s);
1777 if (isinf(s[0]) || isnan(s[0]) ||
1778 isinf(s[1]) || isnan(s[1]) ||
1779 isinf(s[2]) || isnan(s[2]) ) {
1780 sprintf( painCave.errMsg,
1781 "RNEMD detected a numerical error writing: %s std. dev. for bin %d",
1782 data_[index].title.c_str(), bin);
1783 painCave.isFatal = 1;
1784 simError();
1785 } else {
1786 rnemdFile_ << "\t" << s[0] << "\t" << s[1] << "\t" << s[2];
1787 }
1788 }
1789 }
1790

Properties

Name Value
svn:keywords Author Id Revision Date