ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/rnemd/RNEMD.cpp
Revision: 1813
Committed: Fri Nov 16 21:39:55 2012 UTC (12 years, 5 months ago) by gezelter
File size: 56529 byte(s)
Log Message:
Forgot to merge some important accumulator changes

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

Properties

Name Value
svn:keywords Author Id Revision Date