ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/rnemd/RNEMD.cpp
Revision: 1938
Committed: Thu Oct 31 15:32:17 2013 UTC (11 years, 6 months ago) by gezelter
File size: 76187 byte(s)
Log Message:
Some MPI include re-ordering to work with the Intel MPI implementation.

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, 234107 (2008).
39 * [4] Vardeman & Gezelter, in progress (2009).
40 */
41 #ifdef IS_MPI
42 #include <mpi.h>
43 #endif
44
45 #include <cmath>
46 #include <sstream>
47 #include <string>
48
49 #include "rnemd/RNEMD.hpp"
50 #include "math/Vector3.hpp"
51 #include "math/Vector.hpp"
52 #include "math/SquareMatrix3.hpp"
53 #include "math/Polynomial.hpp"
54 #include "primitives/Molecule.hpp"
55 #include "primitives/StuntDouble.hpp"
56 #include "utils/PhysicalConstants.hpp"
57 #include "utils/Tuple.hpp"
58 #include "brains/Thermo.hpp"
59 #include "math/ConvexHull.hpp"
60
61 #ifdef _MSC_VER
62 #define isnan(x) _isnan((x))
63 #define isinf(x) (!_finite(x) && !_isnan(x))
64 #endif
65
66 #define HONKING_LARGE_VALUE 1.0e10
67
68 using namespace std;
69 namespace OpenMD {
70
71 RNEMD::RNEMD(SimInfo* info) : info_(info), evaluator_(info), seleMan_(info),
72 evaluatorA_(info), seleManA_(info),
73 commonA_(info), evaluatorB_(info),
74 seleManB_(info), commonB_(info),
75 hasData_(false), hasDividingArea_(false),
76 usePeriodicBoundaryConditions_(info->getSimParams()->getUsePeriodicBoundaryConditions()) {
77
78 trialCount_ = 0;
79 failTrialCount_ = 0;
80 failRootCount_ = 0;
81
82 Globals* simParams = info->getSimParams();
83 RNEMDParameters* rnemdParams = simParams->getRNEMDParameters();
84
85 doRNEMD_ = rnemdParams->getUseRNEMD();
86 if (!doRNEMD_) return;
87
88 stringToMethod_["Swap"] = rnemdSwap;
89 stringToMethod_["NIVS"] = rnemdNIVS;
90 stringToMethod_["VSS"] = rnemdVSS;
91
92 stringToFluxType_["KE"] = rnemdKE;
93 stringToFluxType_["Px"] = rnemdPx;
94 stringToFluxType_["Py"] = rnemdPy;
95 stringToFluxType_["Pz"] = rnemdPz;
96 stringToFluxType_["Pvector"] = rnemdPvector;
97 stringToFluxType_["Lx"] = rnemdLx;
98 stringToFluxType_["Ly"] = rnemdLy;
99 stringToFluxType_["Lz"] = rnemdLz;
100 stringToFluxType_["Lvector"] = rnemdLvector;
101 stringToFluxType_["KE+Px"] = rnemdKePx;
102 stringToFluxType_["KE+Py"] = rnemdKePy;
103 stringToFluxType_["KE+Pvector"] = rnemdKePvector;
104 stringToFluxType_["KE+Lx"] = rnemdKeLx;
105 stringToFluxType_["KE+Ly"] = rnemdKeLy;
106 stringToFluxType_["KE+Lz"] = rnemdKeLz;
107 stringToFluxType_["KE+Lvector"] = rnemdKeLvector;
108
109 runTime_ = simParams->getRunTime();
110 statusTime_ = simParams->getStatusTime();
111
112 const string methStr = rnemdParams->getMethod();
113 bool hasFluxType = rnemdParams->haveFluxType();
114
115 rnemdObjectSelection_ = rnemdParams->getObjectSelection();
116
117 string fluxStr;
118 if (hasFluxType) {
119 fluxStr = rnemdParams->getFluxType();
120 } else {
121 sprintf(painCave.errMsg,
122 "RNEMD: No fluxType was set in the md file. This parameter,\n"
123 "\twhich must be one of the following values:\n"
124 "\tKE, Px, Py, Pz, Pvector, Lx, Ly, Lz, Lvector,\n"
125 "\tKE+Px, KE+Py, KE+Pvector, KE+Lx, KE+Ly, KE+Lz, KE+Lvector\n"
126 "\tmust be set to use RNEMD\n");
127 painCave.isFatal = 1;
128 painCave.severity = OPENMD_ERROR;
129 simError();
130 }
131
132 bool hasKineticFlux = rnemdParams->haveKineticFlux();
133 bool hasMomentumFlux = rnemdParams->haveMomentumFlux();
134 bool hasMomentumFluxVector = rnemdParams->haveMomentumFluxVector();
135 bool hasAngularMomentumFlux = rnemdParams->haveAngularMomentumFlux();
136 bool hasAngularMomentumFluxVector = rnemdParams->haveAngularMomentumFluxVector();
137 hasSelectionA_ = rnemdParams->haveSelectionA();
138 hasSelectionB_ = rnemdParams->haveSelectionB();
139 bool hasSlabWidth = rnemdParams->haveSlabWidth();
140 bool hasSlabACenter = rnemdParams->haveSlabACenter();
141 bool hasSlabBCenter = rnemdParams->haveSlabBCenter();
142 bool hasSphereARadius = rnemdParams->haveSphereARadius();
143 hasSphereBRadius_ = rnemdParams->haveSphereBRadius();
144 bool hasCoordinateOrigin = rnemdParams->haveCoordinateOrigin();
145 bool hasOutputFileName = rnemdParams->haveOutputFileName();
146 bool hasOutputFields = rnemdParams->haveOutputFields();
147
148 map<string, RNEMDMethod>::iterator i;
149 i = stringToMethod_.find(methStr);
150 if (i != stringToMethod_.end())
151 rnemdMethod_ = i->second;
152 else {
153 sprintf(painCave.errMsg,
154 "RNEMD: The current method,\n"
155 "\t\t%s is not one of the recognized\n"
156 "\texchange methods: Swap, NIVS, or VSS\n",
157 methStr.c_str());
158 painCave.isFatal = 1;
159 painCave.severity = OPENMD_ERROR;
160 simError();
161 }
162
163 map<string, RNEMDFluxType>::iterator j;
164 j = stringToFluxType_.find(fluxStr);
165 if (j != stringToFluxType_.end())
166 rnemdFluxType_ = j->second;
167 else {
168 sprintf(painCave.errMsg,
169 "RNEMD: The current fluxType,\n"
170 "\t\t%s\n"
171 "\tis not one of the recognized flux types.\n",
172 fluxStr.c_str());
173 painCave.isFatal = 1;
174 painCave.severity = OPENMD_ERROR;
175 simError();
176 }
177
178 bool methodFluxMismatch = false;
179 bool hasCorrectFlux = false;
180 switch(rnemdMethod_) {
181 case rnemdSwap:
182 switch (rnemdFluxType_) {
183 case rnemdKE:
184 hasCorrectFlux = hasKineticFlux;
185 break;
186 case rnemdPx:
187 case rnemdPy:
188 case rnemdPz:
189 hasCorrectFlux = hasMomentumFlux;
190 break;
191 default :
192 methodFluxMismatch = true;
193 break;
194 }
195 break;
196 case rnemdNIVS:
197 switch (rnemdFluxType_) {
198 case rnemdKE:
199 case rnemdRotKE:
200 case rnemdFullKE:
201 hasCorrectFlux = hasKineticFlux;
202 break;
203 case rnemdPx:
204 case rnemdPy:
205 case rnemdPz:
206 hasCorrectFlux = hasMomentumFlux;
207 break;
208 case rnemdKePx:
209 case rnemdKePy:
210 hasCorrectFlux = hasMomentumFlux && hasKineticFlux;
211 break;
212 default:
213 methodFluxMismatch = true;
214 break;
215 }
216 break;
217 case rnemdVSS:
218 switch (rnemdFluxType_) {
219 case rnemdKE:
220 case rnemdRotKE:
221 case rnemdFullKE:
222 hasCorrectFlux = hasKineticFlux;
223 break;
224 case rnemdPx:
225 case rnemdPy:
226 case rnemdPz:
227 hasCorrectFlux = hasMomentumFlux;
228 break;
229 case rnemdLx:
230 case rnemdLy:
231 case rnemdLz:
232 hasCorrectFlux = hasAngularMomentumFlux;
233 break;
234 case rnemdPvector:
235 hasCorrectFlux = hasMomentumFluxVector;
236 break;
237 case rnemdLvector:
238 hasCorrectFlux = hasAngularMomentumFluxVector;
239 break;
240 case rnemdKePx:
241 case rnemdKePy:
242 hasCorrectFlux = hasMomentumFlux && hasKineticFlux;
243 break;
244 case rnemdKeLx:
245 case rnemdKeLy:
246 case rnemdKeLz:
247 hasCorrectFlux = hasAngularMomentumFlux && hasKineticFlux;
248 break;
249 case rnemdKePvector:
250 hasCorrectFlux = hasMomentumFluxVector && hasKineticFlux;
251 break;
252 case rnemdKeLvector:
253 hasCorrectFlux = hasAngularMomentumFluxVector && hasKineticFlux;
254 break;
255 default:
256 methodFluxMismatch = true;
257 break;
258 }
259 default:
260 break;
261 }
262
263 if (methodFluxMismatch) {
264 sprintf(painCave.errMsg,
265 "RNEMD: The current method,\n"
266 "\t\t%s\n"
267 "\tcannot be used with the current flux type, %s\n",
268 methStr.c_str(), fluxStr.c_str());
269 painCave.isFatal = 1;
270 painCave.severity = OPENMD_ERROR;
271 simError();
272 }
273 if (!hasCorrectFlux) {
274 sprintf(painCave.errMsg,
275 "RNEMD: The current method, %s, and flux type, %s,\n"
276 "\tdid not have the correct flux value specified. Options\n"
277 "\tinclude: kineticFlux, momentumFlux, angularMomentumFlux,\n"
278 "\tmomentumFluxVector, and angularMomentumFluxVector.\n",
279 methStr.c_str(), fluxStr.c_str());
280 painCave.isFatal = 1;
281 painCave.severity = OPENMD_ERROR;
282 simError();
283 }
284
285 if (hasKineticFlux) {
286 // convert the kcal / mol / Angstroms^2 / fs values in the md file
287 // into amu / fs^3:
288 kineticFlux_ = rnemdParams->getKineticFlux()
289 * PhysicalConstants::energyConvert;
290 } else {
291 kineticFlux_ = 0.0;
292 }
293 if (hasMomentumFluxVector) {
294 momentumFluxVector_ = rnemdParams->getMomentumFluxVector();
295 } else {
296 momentumFluxVector_ = V3Zero;
297 if (hasMomentumFlux) {
298 RealType momentumFlux = rnemdParams->getMomentumFlux();
299 switch (rnemdFluxType_) {
300 case rnemdPx:
301 momentumFluxVector_.x() = momentumFlux;
302 break;
303 case rnemdPy:
304 momentumFluxVector_.y() = momentumFlux;
305 break;
306 case rnemdPz:
307 momentumFluxVector_.z() = momentumFlux;
308 break;
309 case rnemdKePx:
310 momentumFluxVector_.x() = momentumFlux;
311 break;
312 case rnemdKePy:
313 momentumFluxVector_.y() = momentumFlux;
314 break;
315 default:
316 break;
317 }
318 }
319 if (hasAngularMomentumFluxVector) {
320 angularMomentumFluxVector_ = rnemdParams->getAngularMomentumFluxVector();
321 } else {
322 angularMomentumFluxVector_ = V3Zero;
323 if (hasAngularMomentumFlux) {
324 RealType angularMomentumFlux = rnemdParams->getAngularMomentumFlux();
325 switch (rnemdFluxType_) {
326 case rnemdLx:
327 angularMomentumFluxVector_.x() = angularMomentumFlux;
328 break;
329 case rnemdLy:
330 angularMomentumFluxVector_.y() = angularMomentumFlux;
331 break;
332 case rnemdLz:
333 angularMomentumFluxVector_.z() = angularMomentumFlux;
334 break;
335 case rnemdKeLx:
336 angularMomentumFluxVector_.x() = angularMomentumFlux;
337 break;
338 case rnemdKeLy:
339 angularMomentumFluxVector_.y() = angularMomentumFlux;
340 break;
341 case rnemdKeLz:
342 angularMomentumFluxVector_.z() = angularMomentumFlux;
343 break;
344 default:
345 break;
346 }
347 }
348 }
349
350 if (hasCoordinateOrigin) {
351 coordinateOrigin_ = rnemdParams->getCoordinateOrigin();
352 } else {
353 coordinateOrigin_ = V3Zero;
354 }
355
356 // do some sanity checking
357
358 int selectionCount = seleMan_.getSelectionCount();
359
360 int nIntegrable = info->getNGlobalIntegrableObjects();
361
362 if (selectionCount > nIntegrable) {
363 sprintf(painCave.errMsg,
364 "RNEMD: The current objectSelection,\n"
365 "\t\t%s\n"
366 "\thas resulted in %d selected objects. However,\n"
367 "\tthe total number of integrable objects in the system\n"
368 "\tis only %d. This is almost certainly not what you want\n"
369 "\tto do. A likely cause of this is forgetting the _RB_0\n"
370 "\tselector in the selection script!\n",
371 rnemdObjectSelection_.c_str(),
372 selectionCount, nIntegrable);
373 painCave.isFatal = 0;
374 painCave.severity = OPENMD_WARNING;
375 simError();
376 }
377
378 areaAccumulator_ = new Accumulator();
379
380 nBins_ = rnemdParams->getOutputBins();
381 binWidth_ = rnemdParams->getOutputBinWidth();
382
383 data_.resize(RNEMD::ENDINDEX);
384 OutputData z;
385 z.units = "Angstroms";
386 z.title = "Z";
387 z.dataType = "RealType";
388 z.accumulator.reserve(nBins_);
389 for (int i = 0; i < nBins_; i++)
390 z.accumulator.push_back( new Accumulator() );
391 data_[Z] = z;
392 outputMap_["Z"] = Z;
393
394 OutputData r;
395 r.units = "Angstroms";
396 r.title = "R";
397 r.dataType = "RealType";
398 r.accumulator.reserve(nBins_);
399 for (int i = 0; i < nBins_; i++)
400 r.accumulator.push_back( new Accumulator() );
401 data_[R] = r;
402 outputMap_["R"] = R;
403
404 OutputData temperature;
405 temperature.units = "K";
406 temperature.title = "Temperature";
407 temperature.dataType = "RealType";
408 temperature.accumulator.reserve(nBins_);
409 for (int i = 0; i < nBins_; i++)
410 temperature.accumulator.push_back( new Accumulator() );
411 data_[TEMPERATURE] = temperature;
412 outputMap_["TEMPERATURE"] = TEMPERATURE;
413
414 OutputData velocity;
415 velocity.units = "angstroms/fs";
416 velocity.title = "Velocity";
417 velocity.dataType = "Vector3d";
418 velocity.accumulator.reserve(nBins_);
419 for (int i = 0; i < nBins_; i++)
420 velocity.accumulator.push_back( new VectorAccumulator() );
421 data_[VELOCITY] = velocity;
422 outputMap_["VELOCITY"] = VELOCITY;
423
424 OutputData angularVelocity;
425 angularVelocity.units = "angstroms^2/fs";
426 angularVelocity.title = "AngularVelocity";
427 angularVelocity.dataType = "Vector3d";
428 angularVelocity.accumulator.reserve(nBins_);
429 for (int i = 0; i < nBins_; i++)
430 angularVelocity.accumulator.push_back( new VectorAccumulator() );
431 data_[ANGULARVELOCITY] = angularVelocity;
432 outputMap_["ANGULARVELOCITY"] = ANGULARVELOCITY;
433
434 OutputData density;
435 density.units = "g cm^-3";
436 density.title = "Density";
437 density.dataType = "RealType";
438 density.accumulator.reserve(nBins_);
439 for (int i = 0; i < nBins_; i++)
440 density.accumulator.push_back( new Accumulator() );
441 data_[DENSITY] = density;
442 outputMap_["DENSITY"] = DENSITY;
443
444 if (hasOutputFields) {
445 parseOutputFileFormat(rnemdParams->getOutputFields());
446 } else {
447 if (usePeriodicBoundaryConditions_)
448 outputMask_.set(Z);
449 else
450 outputMask_.set(R);
451 switch (rnemdFluxType_) {
452 case rnemdKE:
453 case rnemdRotKE:
454 case rnemdFullKE:
455 outputMask_.set(TEMPERATURE);
456 break;
457 case rnemdPx:
458 case rnemdPy:
459 outputMask_.set(VELOCITY);
460 break;
461 case rnemdPz:
462 case rnemdPvector:
463 outputMask_.set(VELOCITY);
464 outputMask_.set(DENSITY);
465 break;
466 case rnemdLx:
467 case rnemdLy:
468 case rnemdLz:
469 case rnemdLvector:
470 outputMask_.set(ANGULARVELOCITY);
471 break;
472 case rnemdKeLx:
473 case rnemdKeLy:
474 case rnemdKeLz:
475 case rnemdKeLvector:
476 outputMask_.set(TEMPERATURE);
477 outputMask_.set(ANGULARVELOCITY);
478 break;
479 case rnemdKePx:
480 case rnemdKePy:
481 outputMask_.set(TEMPERATURE);
482 outputMask_.set(VELOCITY);
483 break;
484 case rnemdKePvector:
485 outputMask_.set(TEMPERATURE);
486 outputMask_.set(VELOCITY);
487 outputMask_.set(DENSITY);
488 break;
489 default:
490 break;
491 }
492 }
493
494 if (hasOutputFileName) {
495 rnemdFileName_ = rnemdParams->getOutputFileName();
496 } else {
497 rnemdFileName_ = getPrefix(info->getFinalConfigFileName()) + ".rnemd";
498 }
499
500 exchangeTime_ = rnemdParams->getExchangeTime();
501
502 Snapshot* currentSnap_ = info->getSnapshotManager()->getCurrentSnapshot();
503 // total exchange sums are zeroed out at the beginning:
504
505 kineticExchange_ = 0.0;
506 momentumExchange_ = V3Zero;
507 angularMomentumExchange_ = V3Zero;
508
509 std::ostringstream selectionAstream;
510 std::ostringstream selectionBstream;
511
512 if (hasSelectionA_) {
513 selectionA_ = rnemdParams->getSelectionA();
514 } else {
515 if (usePeriodicBoundaryConditions_) {
516 Mat3x3d hmat = currentSnap_->getHmat();
517
518 if (hasSlabWidth)
519 slabWidth_ = rnemdParams->getSlabWidth();
520 else
521 slabWidth_ = hmat(2,2) / 10.0;
522
523 if (hasSlabACenter)
524 slabACenter_ = rnemdParams->getSlabACenter();
525 else
526 slabACenter_ = 0.0;
527
528 selectionAstream << "select wrappedz > "
529 << slabACenter_ - 0.5*slabWidth_
530 << " && wrappedz < "
531 << slabACenter_ + 0.5*slabWidth_;
532 selectionA_ = selectionAstream.str();
533 } else {
534 if (hasSphereARadius)
535 sphereARadius_ = rnemdParams->getSphereARadius();
536 else {
537 // use an initial guess to the size of the inner slab to be 1/10 the
538 // radius of an approximately spherical hull:
539 Thermo thermo(info);
540 RealType hVol = thermo.getHullVolume();
541 sphereARadius_ = 0.1 * pow((3.0 * hVol / (4.0 * M_PI)), 1.0/3.0);
542 }
543 selectionAstream << "select r < " << sphereARadius_;
544 selectionA_ = selectionAstream.str();
545 }
546 }
547
548 if (hasSelectionB_) {
549 selectionB_ = rnemdParams->getSelectionB();
550
551 } else {
552 if (usePeriodicBoundaryConditions_) {
553 Mat3x3d hmat = currentSnap_->getHmat();
554
555 if (hasSlabWidth)
556 slabWidth_ = rnemdParams->getSlabWidth();
557 else
558 slabWidth_ = hmat(2,2) / 10.0;
559
560 if (hasSlabBCenter)
561 slabBCenter_ = rnemdParams->getSlabBCenter();
562 else
563 slabBCenter_ = hmat(2,2) / 2.0;
564
565 selectionBstream << "select wrappedz > "
566 << slabBCenter_ - 0.5*slabWidth_
567 << " && wrappedz < "
568 << slabBCenter_ + 0.5*slabWidth_;
569 selectionB_ = selectionBstream.str();
570 } else {
571 if (hasSphereBRadius_) {
572 sphereBRadius_ = rnemdParams->getSphereBRadius();
573 selectionBstream << "select r > " << sphereBRadius_;
574 selectionB_ = selectionBstream.str();
575 } else {
576 selectionB_ = "select hull";
577 BisHull_ = true;
578 hasSelectionB_ = true;
579 }
580 }
581 }
582 }
583
584 // object evaluator:
585 evaluator_.loadScriptString(rnemdObjectSelection_);
586 seleMan_.setSelectionSet(evaluator_.evaluate());
587 evaluatorA_.loadScriptString(selectionA_);
588 evaluatorB_.loadScriptString(selectionB_);
589 seleManA_.setSelectionSet(evaluatorA_.evaluate());
590 seleManB_.setSelectionSet(evaluatorB_.evaluate());
591 commonA_ = seleManA_ & seleMan_;
592 commonB_ = seleManB_ & seleMan_;
593 }
594
595
596 RNEMD::~RNEMD() {
597 if (!doRNEMD_) return;
598 #ifdef IS_MPI
599 if (worldRank == 0) {
600 #endif
601
602 writeOutputFile();
603
604 rnemdFile_.close();
605
606 #ifdef IS_MPI
607 }
608 #endif
609
610 // delete all of the objects we created:
611 delete areaAccumulator_;
612 data_.clear();
613 }
614
615 void RNEMD::doSwap(SelectionManager& smanA, SelectionManager& smanB) {
616 if (!doRNEMD_) return;
617 int selei;
618 int selej;
619
620 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
621 Mat3x3d hmat = currentSnap_->getHmat();
622
623 StuntDouble* sd;
624
625 RealType min_val;
626 bool min_found = false;
627 StuntDouble* min_sd;
628
629 RealType max_val;
630 bool max_found = false;
631 StuntDouble* max_sd;
632
633 for (sd = seleManA_.beginSelected(selei); sd != NULL;
634 sd = seleManA_.nextSelected(selei)) {
635
636 Vector3d pos = sd->getPos();
637
638 // wrap the stuntdouble's position back into the box:
639
640 if (usePeriodicBoundaryConditions_)
641 currentSnap_->wrapVector(pos);
642
643 RealType mass = sd->getMass();
644 Vector3d vel = sd->getVel();
645 RealType value;
646
647 switch(rnemdFluxType_) {
648 case rnemdKE :
649
650 value = mass * vel.lengthSquare();
651
652 if (sd->isDirectional()) {
653 Vector3d angMom = sd->getJ();
654 Mat3x3d I = sd->getI();
655
656 if (sd->isLinear()) {
657 int i = sd->linearAxis();
658 int j = (i + 1) % 3;
659 int k = (i + 2) % 3;
660 value += angMom[j] * angMom[j] / I(j, j) +
661 angMom[k] * angMom[k] / I(k, k);
662 } else {
663 value += angMom[0]*angMom[0]/I(0, 0)
664 + angMom[1]*angMom[1]/I(1, 1)
665 + angMom[2]*angMom[2]/I(2, 2);
666 }
667 } //angular momenta exchange enabled
668 value *= 0.5;
669 break;
670 case rnemdPx :
671 value = mass * vel[0];
672 break;
673 case rnemdPy :
674 value = mass * vel[1];
675 break;
676 case rnemdPz :
677 value = mass * vel[2];
678 break;
679 default :
680 break;
681 }
682 if (!max_found) {
683 max_val = value;
684 max_sd = sd;
685 max_found = true;
686 } else {
687 if (max_val < value) {
688 max_val = value;
689 max_sd = sd;
690 }
691 }
692 }
693
694 for (sd = seleManB_.beginSelected(selej); sd != NULL;
695 sd = seleManB_.nextSelected(selej)) {
696
697 Vector3d pos = sd->getPos();
698
699 // wrap the stuntdouble's position back into the box:
700
701 if (usePeriodicBoundaryConditions_)
702 currentSnap_->wrapVector(pos);
703
704 RealType mass = sd->getMass();
705 Vector3d vel = sd->getVel();
706 RealType value;
707
708 switch(rnemdFluxType_) {
709 case rnemdKE :
710
711 value = mass * vel.lengthSquare();
712
713 if (sd->isDirectional()) {
714 Vector3d angMom = sd->getJ();
715 Mat3x3d I = sd->getI();
716
717 if (sd->isLinear()) {
718 int i = sd->linearAxis();
719 int j = (i + 1) % 3;
720 int k = (i + 2) % 3;
721 value += angMom[j] * angMom[j] / I(j, j) +
722 angMom[k] * angMom[k] / I(k, k);
723 } else {
724 value += angMom[0]*angMom[0]/I(0, 0)
725 + angMom[1]*angMom[1]/I(1, 1)
726 + angMom[2]*angMom[2]/I(2, 2);
727 }
728 } //angular momenta exchange enabled
729 value *= 0.5;
730 break;
731 case rnemdPx :
732 value = mass * vel[0];
733 break;
734 case rnemdPy :
735 value = mass * vel[1];
736 break;
737 case rnemdPz :
738 value = mass * vel[2];
739 break;
740 default :
741 break;
742 }
743
744 if (!min_found) {
745 min_val = value;
746 min_sd = sd;
747 min_found = true;
748 } else {
749 if (min_val > value) {
750 min_val = value;
751 min_sd = sd;
752 }
753 }
754 }
755
756 #ifdef IS_MPI
757 int worldRank = MPI::COMM_WORLD.Get_rank();
758
759 bool my_min_found = min_found;
760 bool my_max_found = max_found;
761
762 // Even if we didn't find a minimum, did someone else?
763 MPI::COMM_WORLD.Allreduce(&my_min_found, &min_found, 1, MPI::BOOL, MPI::LOR);
764 // Even if we didn't find a maximum, did someone else?
765 MPI::COMM_WORLD.Allreduce(&my_max_found, &max_found, 1, MPI::BOOL, MPI::LOR);
766 #endif
767
768 if (max_found && min_found) {
769
770 #ifdef IS_MPI
771 struct {
772 RealType val;
773 int rank;
774 } max_vals, min_vals;
775
776 if (my_min_found) {
777 min_vals.val = min_val;
778 } else {
779 min_vals.val = HONKING_LARGE_VALUE;
780 }
781 min_vals.rank = worldRank;
782
783 // Who had the minimum?
784 MPI::COMM_WORLD.Allreduce(&min_vals, &min_vals,
785 1, MPI::REALTYPE_INT, MPI::MINLOC);
786 min_val = min_vals.val;
787
788 if (my_max_found) {
789 max_vals.val = max_val;
790 } else {
791 max_vals.val = -HONKING_LARGE_VALUE;
792 }
793 max_vals.rank = worldRank;
794
795 // Who had the maximum?
796 MPI::COMM_WORLD.Allreduce(&max_vals, &max_vals,
797 1, MPI::REALTYPE_INT, MPI::MAXLOC);
798 max_val = max_vals.val;
799 #endif
800
801 if (min_val < max_val) {
802
803 #ifdef IS_MPI
804 if (max_vals.rank == worldRank && min_vals.rank == worldRank) {
805 // I have both maximum and minimum, so proceed like a single
806 // processor version:
807 #endif
808
809 Vector3d min_vel = min_sd->getVel();
810 Vector3d max_vel = max_sd->getVel();
811 RealType temp_vel;
812
813 switch(rnemdFluxType_) {
814 case rnemdKE :
815 min_sd->setVel(max_vel);
816 max_sd->setVel(min_vel);
817 if (min_sd->isDirectional() && max_sd->isDirectional()) {
818 Vector3d min_angMom = min_sd->getJ();
819 Vector3d max_angMom = max_sd->getJ();
820 min_sd->setJ(max_angMom);
821 max_sd->setJ(min_angMom);
822 }//angular momenta exchange enabled
823 //assumes same rigid body identity
824 break;
825 case rnemdPx :
826 temp_vel = min_vel.x();
827 min_vel.x() = max_vel.x();
828 max_vel.x() = temp_vel;
829 min_sd->setVel(min_vel);
830 max_sd->setVel(max_vel);
831 break;
832 case rnemdPy :
833 temp_vel = min_vel.y();
834 min_vel.y() = max_vel.y();
835 max_vel.y() = temp_vel;
836 min_sd->setVel(min_vel);
837 max_sd->setVel(max_vel);
838 break;
839 case rnemdPz :
840 temp_vel = min_vel.z();
841 min_vel.z() = max_vel.z();
842 max_vel.z() = temp_vel;
843 min_sd->setVel(min_vel);
844 max_sd->setVel(max_vel);
845 break;
846 default :
847 break;
848 }
849
850 #ifdef IS_MPI
851 // the rest of the cases only apply in parallel simulations:
852 } else if (max_vals.rank == worldRank) {
853 // I had the max, but not the minimum
854
855 Vector3d min_vel;
856 Vector3d max_vel = max_sd->getVel();
857 MPI::Status status;
858
859 // point-to-point swap of the velocity vector
860 MPI::COMM_WORLD.Sendrecv(max_vel.getArrayPointer(), 3, MPI::REALTYPE,
861 min_vals.rank, 0,
862 min_vel.getArrayPointer(), 3, MPI::REALTYPE,
863 min_vals.rank, 0, status);
864
865 switch(rnemdFluxType_) {
866 case rnemdKE :
867 max_sd->setVel(min_vel);
868 //angular momenta exchange enabled
869 if (max_sd->isDirectional()) {
870 Vector3d min_angMom;
871 Vector3d max_angMom = max_sd->getJ();
872
873 // point-to-point swap of the angular momentum vector
874 MPI::COMM_WORLD.Sendrecv(max_angMom.getArrayPointer(), 3,
875 MPI::REALTYPE, min_vals.rank, 1,
876 min_angMom.getArrayPointer(), 3,
877 MPI::REALTYPE, min_vals.rank, 1,
878 status);
879
880 max_sd->setJ(min_angMom);
881 }
882 break;
883 case rnemdPx :
884 max_vel.x() = min_vel.x();
885 max_sd->setVel(max_vel);
886 break;
887 case rnemdPy :
888 max_vel.y() = min_vel.y();
889 max_sd->setVel(max_vel);
890 break;
891 case rnemdPz :
892 max_vel.z() = min_vel.z();
893 max_sd->setVel(max_vel);
894 break;
895 default :
896 break;
897 }
898 } else if (min_vals.rank == worldRank) {
899 // I had the minimum but not the maximum:
900
901 Vector3d max_vel;
902 Vector3d min_vel = min_sd->getVel();
903 MPI::Status status;
904
905 // point-to-point swap of the velocity vector
906 MPI::COMM_WORLD.Sendrecv(min_vel.getArrayPointer(), 3, MPI::REALTYPE,
907 max_vals.rank, 0,
908 max_vel.getArrayPointer(), 3, MPI::REALTYPE,
909 max_vals.rank, 0, status);
910
911 switch(rnemdFluxType_) {
912 case rnemdKE :
913 min_sd->setVel(max_vel);
914 //angular momenta exchange enabled
915 if (min_sd->isDirectional()) {
916 Vector3d min_angMom = min_sd->getJ();
917 Vector3d max_angMom;
918
919 // point-to-point swap of the angular momentum vector
920 MPI::COMM_WORLD.Sendrecv(min_angMom.getArrayPointer(), 3,
921 MPI::REALTYPE, max_vals.rank, 1,
922 max_angMom.getArrayPointer(), 3,
923 MPI::REALTYPE, max_vals.rank, 1,
924 status);
925
926 min_sd->setJ(max_angMom);
927 }
928 break;
929 case rnemdPx :
930 min_vel.x() = max_vel.x();
931 min_sd->setVel(min_vel);
932 break;
933 case rnemdPy :
934 min_vel.y() = max_vel.y();
935 min_sd->setVel(min_vel);
936 break;
937 case rnemdPz :
938 min_vel.z() = max_vel.z();
939 min_sd->setVel(min_vel);
940 break;
941 default :
942 break;
943 }
944 }
945 #endif
946
947 switch(rnemdFluxType_) {
948 case rnemdKE:
949 kineticExchange_ += max_val - min_val;
950 break;
951 case rnemdPx:
952 momentumExchange_.x() += max_val - min_val;
953 break;
954 case rnemdPy:
955 momentumExchange_.y() += max_val - min_val;
956 break;
957 case rnemdPz:
958 momentumExchange_.z() += max_val - min_val;
959 break;
960 default:
961 break;
962 }
963 } else {
964 sprintf(painCave.errMsg,
965 "RNEMD::doSwap exchange NOT performed because min_val > max_val\n");
966 painCave.isFatal = 0;
967 painCave.severity = OPENMD_INFO;
968 simError();
969 failTrialCount_++;
970 }
971 } else {
972 sprintf(painCave.errMsg,
973 "RNEMD::doSwap exchange NOT performed because selected object\n"
974 "\twas not present in at least one of the two slabs.\n");
975 painCave.isFatal = 0;
976 painCave.severity = OPENMD_INFO;
977 simError();
978 failTrialCount_++;
979 }
980 }
981
982 void RNEMD::doNIVS(SelectionManager& smanA, SelectionManager& smanB) {
983 if (!doRNEMD_) return;
984 int selei;
985 int selej;
986
987 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
988 RealType time = currentSnap_->getTime();
989 Mat3x3d hmat = currentSnap_->getHmat();
990
991 StuntDouble* sd;
992
993 vector<StuntDouble*> hotBin, coldBin;
994
995 RealType Phx = 0.0;
996 RealType Phy = 0.0;
997 RealType Phz = 0.0;
998 RealType Khx = 0.0;
999 RealType Khy = 0.0;
1000 RealType Khz = 0.0;
1001 RealType Khw = 0.0;
1002 RealType Pcx = 0.0;
1003 RealType Pcy = 0.0;
1004 RealType Pcz = 0.0;
1005 RealType Kcx = 0.0;
1006 RealType Kcy = 0.0;
1007 RealType Kcz = 0.0;
1008 RealType Kcw = 0.0;
1009
1010 for (sd = smanA.beginSelected(selei); sd != NULL;
1011 sd = smanA.nextSelected(selei)) {
1012
1013 Vector3d pos = sd->getPos();
1014
1015 // wrap the stuntdouble's position back into the box:
1016
1017 if (usePeriodicBoundaryConditions_)
1018 currentSnap_->wrapVector(pos);
1019
1020
1021 RealType mass = sd->getMass();
1022 Vector3d vel = sd->getVel();
1023
1024 hotBin.push_back(sd);
1025 Phx += mass * vel.x();
1026 Phy += mass * vel.y();
1027 Phz += mass * vel.z();
1028 Khx += mass * vel.x() * vel.x();
1029 Khy += mass * vel.y() * vel.y();
1030 Khz += mass * vel.z() * vel.z();
1031 if (sd->isDirectional()) {
1032 Vector3d angMom = sd->getJ();
1033 Mat3x3d I = sd->getI();
1034 if (sd->isLinear()) {
1035 int i = sd->linearAxis();
1036 int j = (i + 1) % 3;
1037 int k = (i + 2) % 3;
1038 Khw += angMom[j] * angMom[j] / I(j, j) +
1039 angMom[k] * angMom[k] / I(k, k);
1040 } else {
1041 Khw += angMom[0]*angMom[0]/I(0, 0)
1042 + angMom[1]*angMom[1]/I(1, 1)
1043 + angMom[2]*angMom[2]/I(2, 2);
1044 }
1045 }
1046 }
1047 for (sd = smanB.beginSelected(selej); sd != NULL;
1048 sd = smanB.nextSelected(selej)) {
1049 Vector3d pos = sd->getPos();
1050
1051 // wrap the stuntdouble's position back into the box:
1052
1053 if (usePeriodicBoundaryConditions_)
1054 currentSnap_->wrapVector(pos);
1055
1056 RealType mass = sd->getMass();
1057 Vector3d vel = sd->getVel();
1058
1059 coldBin.push_back(sd);
1060 Pcx += mass * vel.x();
1061 Pcy += mass * vel.y();
1062 Pcz += mass * vel.z();
1063 Kcx += mass * vel.x() * vel.x();
1064 Kcy += mass * vel.y() * vel.y();
1065 Kcz += mass * vel.z() * vel.z();
1066 if (sd->isDirectional()) {
1067 Vector3d angMom = sd->getJ();
1068 Mat3x3d I = sd->getI();
1069 if (sd->isLinear()) {
1070 int i = sd->linearAxis();
1071 int j = (i + 1) % 3;
1072 int k = (i + 2) % 3;
1073 Kcw += angMom[j] * angMom[j] / I(j, j) +
1074 angMom[k] * angMom[k] / I(k, k);
1075 } else {
1076 Kcw += angMom[0]*angMom[0]/I(0, 0)
1077 + angMom[1]*angMom[1]/I(1, 1)
1078 + angMom[2]*angMom[2]/I(2, 2);
1079 }
1080 }
1081 }
1082
1083 Khx *= 0.5;
1084 Khy *= 0.5;
1085 Khz *= 0.5;
1086 Khw *= 0.5;
1087 Kcx *= 0.5;
1088 Kcy *= 0.5;
1089 Kcz *= 0.5;
1090 Kcw *= 0.5;
1091
1092 #ifdef IS_MPI
1093 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phx, 1, MPI::REALTYPE, MPI::SUM);
1094 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phy, 1, MPI::REALTYPE, MPI::SUM);
1095 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phz, 1, MPI::REALTYPE, MPI::SUM);
1096 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcx, 1, MPI::REALTYPE, MPI::SUM);
1097 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcy, 1, MPI::REALTYPE, MPI::SUM);
1098 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcz, 1, MPI::REALTYPE, MPI::SUM);
1099
1100 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khx, 1, MPI::REALTYPE, MPI::SUM);
1101 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khy, 1, MPI::REALTYPE, MPI::SUM);
1102 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khz, 1, MPI::REALTYPE, MPI::SUM);
1103 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khw, 1, MPI::REALTYPE, MPI::SUM);
1104
1105 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcx, 1, MPI::REALTYPE, MPI::SUM);
1106 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcy, 1, MPI::REALTYPE, MPI::SUM);
1107 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcz, 1, MPI::REALTYPE, MPI::SUM);
1108 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcw, 1, MPI::REALTYPE, MPI::SUM);
1109 #endif
1110
1111 //solve coldBin coeff's first
1112 RealType px = Pcx / Phx;
1113 RealType py = Pcy / Phy;
1114 RealType pz = Pcz / Phz;
1115 RealType c, x, y, z;
1116 bool successfulScale = false;
1117 if ((rnemdFluxType_ == rnemdFullKE) ||
1118 (rnemdFluxType_ == rnemdRotKE)) {
1119 //may need sanity check Khw & Kcw > 0
1120
1121 if (rnemdFluxType_ == rnemdFullKE) {
1122 c = 1.0 - kineticTarget_ / (Kcx + Kcy + Kcz + Kcw);
1123 } else {
1124 c = 1.0 - kineticTarget_ / Kcw;
1125 }
1126
1127 if ((c > 0.81) && (c < 1.21)) {//restrict scaling coefficients
1128 c = sqrt(c);
1129
1130 RealType w = 0.0;
1131 if (rnemdFluxType_ == rnemdFullKE) {
1132 x = 1.0 + px * (1.0 - c);
1133 y = 1.0 + py * (1.0 - c);
1134 z = 1.0 + pz * (1.0 - c);
1135 /* more complicated way
1136 w = 1.0 + (Kcw - Kcw * c * c - (c * c * (Kcx + Kcy + Kcz
1137 + Khx * px * px + Khy * py * py + Khz * pz * pz)
1138 - 2.0 * c * (Khx * px * (1.0 + px) + Khy * py * (1.0 + py)
1139 + Khz * pz * (1.0 + pz)) + Khx * px * (2.0 + px)
1140 + Khy * py * (2.0 + py) + Khz * pz * (2.0 + pz)
1141 - Kcx - Kcy - Kcz)) / Khw; the following is simpler
1142 */
1143 if ((fabs(x - 1.0) < 0.1) && (fabs(y - 1.0) < 0.1) &&
1144 (fabs(z - 1.0) < 0.1)) {
1145 w = 1.0 + (kineticTarget_
1146 + Khx * (1.0 - x * x) + Khy * (1.0 - y * y)
1147 + Khz * (1.0 - z * z)) / Khw;
1148 }//no need to calculate w if x, y or z is out of range
1149 } else {
1150 w = 1.0 + kineticTarget_ / Khw;
1151 }
1152 if ((w > 0.81) && (w < 1.21)) {//restrict scaling coefficients
1153 //if w is in the right range, so should be x, y, z.
1154 vector<StuntDouble*>::iterator sdi;
1155 Vector3d vel;
1156 for (sdi = coldBin.begin(); sdi != coldBin.end(); ++sdi) {
1157 if (rnemdFluxType_ == rnemdFullKE) {
1158 vel = (*sdi)->getVel() * c;
1159 (*sdi)->setVel(vel);
1160 }
1161 if ((*sdi)->isDirectional()) {
1162 Vector3d angMom = (*sdi)->getJ() * c;
1163 (*sdi)->setJ(angMom);
1164 }
1165 }
1166 w = sqrt(w);
1167 for (sdi = hotBin.begin(); sdi != hotBin.end(); ++sdi) {
1168 if (rnemdFluxType_ == rnemdFullKE) {
1169 vel = (*sdi)->getVel();
1170 vel.x() *= x;
1171 vel.y() *= y;
1172 vel.z() *= z;
1173 (*sdi)->setVel(vel);
1174 }
1175 if ((*sdi)->isDirectional()) {
1176 Vector3d angMom = (*sdi)->getJ() * w;
1177 (*sdi)->setJ(angMom);
1178 }
1179 }
1180 successfulScale = true;
1181 kineticExchange_ += kineticTarget_;
1182 }
1183 }
1184 } else {
1185 RealType a000, a110, c0, a001, a111, b01, b11, c1;
1186 switch(rnemdFluxType_) {
1187 case rnemdKE :
1188 /* used hotBin coeff's & only scale x & y dimensions
1189 RealType px = Phx / Pcx;
1190 RealType py = Phy / Pcy;
1191 a110 = Khy;
1192 c0 = - Khx - Khy - kineticTarget_;
1193 a000 = Khx;
1194 a111 = Kcy * py * py;
1195 b11 = -2.0 * Kcy * py * (1.0 + py);
1196 c1 = Kcy * py * (2.0 + py) + Kcx * px * ( 2.0 + px) + kineticTarget_;
1197 b01 = -2.0 * Kcx * px * (1.0 + px);
1198 a001 = Kcx * px * px;
1199 */
1200 //scale all three dimensions, let c_x = c_y
1201 a000 = Kcx + Kcy;
1202 a110 = Kcz;
1203 c0 = kineticTarget_ - Kcx - Kcy - Kcz;
1204 a001 = Khx * px * px + Khy * py * py;
1205 a111 = Khz * pz * pz;
1206 b01 = -2.0 * (Khx * px * (1.0 + px) + Khy * py * (1.0 + py));
1207 b11 = -2.0 * Khz * pz * (1.0 + pz);
1208 c1 = Khx * px * (2.0 + px) + Khy * py * (2.0 + py)
1209 + Khz * pz * (2.0 + pz) - kineticTarget_;
1210 break;
1211 case rnemdPx :
1212 c = 1 - momentumTarget_.x() / Pcx;
1213 a000 = Kcy;
1214 a110 = Kcz;
1215 c0 = Kcx * c * c - Kcx - Kcy - Kcz;
1216 a001 = py * py * Khy;
1217 a111 = pz * pz * Khz;
1218 b01 = -2.0 * Khy * py * (1.0 + py);
1219 b11 = -2.0 * Khz * pz * (1.0 + pz);
1220 c1 = Khy * py * (2.0 + py) + Khz * pz * (2.0 + pz)
1221 + Khx * (fastpow(c * px - px - 1.0, 2) - 1.0);
1222 break;
1223 case rnemdPy :
1224 c = 1 - momentumTarget_.y() / Pcy;
1225 a000 = Kcx;
1226 a110 = Kcz;
1227 c0 = Kcy * c * c - Kcx - Kcy - Kcz;
1228 a001 = px * px * Khx;
1229 a111 = pz * pz * Khz;
1230 b01 = -2.0 * Khx * px * (1.0 + px);
1231 b11 = -2.0 * Khz * pz * (1.0 + pz);
1232 c1 = Khx * px * (2.0 + px) + Khz * pz * (2.0 + pz)
1233 + Khy * (fastpow(c * py - py - 1.0, 2) - 1.0);
1234 break;
1235 case rnemdPz ://we don't really do this, do we?
1236 c = 1 - momentumTarget_.z() / Pcz;
1237 a000 = Kcx;
1238 a110 = Kcy;
1239 c0 = Kcz * c * c - Kcx - Kcy - Kcz;
1240 a001 = px * px * Khx;
1241 a111 = py * py * Khy;
1242 b01 = -2.0 * Khx * px * (1.0 + px);
1243 b11 = -2.0 * Khy * py * (1.0 + py);
1244 c1 = Khx * px * (2.0 + px) + Khy * py * (2.0 + py)
1245 + Khz * (fastpow(c * pz - pz - 1.0, 2) - 1.0);
1246 break;
1247 default :
1248 break;
1249 }
1250
1251 RealType v1 = a000 * a111 - a001 * a110;
1252 RealType v2 = a000 * b01;
1253 RealType v3 = a000 * b11;
1254 RealType v4 = a000 * c1 - a001 * c0;
1255 RealType v8 = a110 * b01;
1256 RealType v10 = - b01 * c0;
1257
1258 RealType u0 = v2 * v10 - v4 * v4;
1259 RealType u1 = -2.0 * v3 * v4;
1260 RealType u2 = -v2 * v8 - v3 * v3 - 2.0 * v1 * v4;
1261 RealType u3 = -2.0 * v1 * v3;
1262 RealType u4 = - v1 * v1;
1263 //rescale coefficients
1264 RealType maxAbs = fabs(u0);
1265 if (maxAbs < fabs(u1)) maxAbs = fabs(u1);
1266 if (maxAbs < fabs(u2)) maxAbs = fabs(u2);
1267 if (maxAbs < fabs(u3)) maxAbs = fabs(u3);
1268 if (maxAbs < fabs(u4)) maxAbs = fabs(u4);
1269 u0 /= maxAbs;
1270 u1 /= maxAbs;
1271 u2 /= maxAbs;
1272 u3 /= maxAbs;
1273 u4 /= maxAbs;
1274 //max_element(start, end) is also available.
1275 Polynomial<RealType> poly; //same as DoublePolynomial poly;
1276 poly.setCoefficient(4, u4);
1277 poly.setCoefficient(3, u3);
1278 poly.setCoefficient(2, u2);
1279 poly.setCoefficient(1, u1);
1280 poly.setCoefficient(0, u0);
1281 vector<RealType> realRoots = poly.FindRealRoots();
1282
1283 vector<RealType>::iterator ri;
1284 RealType r1, r2, alpha0;
1285 vector<pair<RealType,RealType> > rps;
1286 for (ri = realRoots.begin(); ri !=realRoots.end(); ++ri) {
1287 r2 = *ri;
1288 //check if FindRealRoots() give the right answer
1289 if ( fabs(u0 + r2 * (u1 + r2 * (u2 + r2 * (u3 + r2 * u4)))) > 1e-6 ) {
1290 sprintf(painCave.errMsg,
1291 "RNEMD Warning: polynomial solve seems to have an error!");
1292 painCave.isFatal = 0;
1293 simError();
1294 failRootCount_++;
1295 }
1296 //might not be useful w/o rescaling coefficients
1297 alpha0 = -c0 - a110 * r2 * r2;
1298 if (alpha0 >= 0.0) {
1299 r1 = sqrt(alpha0 / a000);
1300 if (fabs(c1 + r1 * (b01 + r1 * a001) + r2 * (b11 + r2 * a111))
1301 < 1e-6)
1302 { rps.push_back(make_pair(r1, r2)); }
1303 if (r1 > 1e-6) { //r1 non-negative
1304 r1 = -r1;
1305 if (fabs(c1 + r1 * (b01 + r1 * a001) + r2 * (b11 + r2 * a111))
1306 < 1e-6)
1307 { rps.push_back(make_pair(r1, r2)); }
1308 }
1309 }
1310 }
1311 // Consider combining together the solving pair part w/ the searching
1312 // best solution part so that we don't need the pairs vector
1313 if (!rps.empty()) {
1314 RealType smallestDiff = HONKING_LARGE_VALUE;
1315 RealType diff;
1316 pair<RealType,RealType> bestPair = make_pair(1.0, 1.0);
1317 vector<pair<RealType,RealType> >::iterator rpi;
1318 for (rpi = rps.begin(); rpi != rps.end(); ++rpi) {
1319 r1 = (*rpi).first;
1320 r2 = (*rpi).second;
1321 switch(rnemdFluxType_) {
1322 case rnemdKE :
1323 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1324 + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcx, 2)
1325 + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcy, 2);
1326 break;
1327 case rnemdPx :
1328 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1329 + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcy, 2);
1330 break;
1331 case rnemdPy :
1332 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1333 + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcx, 2);
1334 break;
1335 case rnemdPz :
1336 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1337 + fastpow(r1 * r1 / r2 / r2 - Kcy/Kcx, 2);
1338 default :
1339 break;
1340 }
1341 if (diff < smallestDiff) {
1342 smallestDiff = diff;
1343 bestPair = *rpi;
1344 }
1345 }
1346 #ifdef IS_MPI
1347 if (worldRank == 0) {
1348 #endif
1349 // sprintf(painCave.errMsg,
1350 // "RNEMD: roots r1= %lf\tr2 = %lf\n",
1351 // bestPair.first, bestPair.second);
1352 // painCave.isFatal = 0;
1353 // painCave.severity = OPENMD_INFO;
1354 // simError();
1355 #ifdef IS_MPI
1356 }
1357 #endif
1358
1359 switch(rnemdFluxType_) {
1360 case rnemdKE :
1361 x = bestPair.first;
1362 y = bestPair.first;
1363 z = bestPair.second;
1364 break;
1365 case rnemdPx :
1366 x = c;
1367 y = bestPair.first;
1368 z = bestPair.second;
1369 break;
1370 case rnemdPy :
1371 x = bestPair.first;
1372 y = c;
1373 z = bestPair.second;
1374 break;
1375 case rnemdPz :
1376 x = bestPair.first;
1377 y = bestPair.second;
1378 z = c;
1379 break;
1380 default :
1381 break;
1382 }
1383 vector<StuntDouble*>::iterator sdi;
1384 Vector3d vel;
1385 for (sdi = coldBin.begin(); sdi != coldBin.end(); ++sdi) {
1386 vel = (*sdi)->getVel();
1387 vel.x() *= x;
1388 vel.y() *= y;
1389 vel.z() *= z;
1390 (*sdi)->setVel(vel);
1391 }
1392 //convert to hotBin coefficient
1393 x = 1.0 + px * (1.0 - x);
1394 y = 1.0 + py * (1.0 - y);
1395 z = 1.0 + pz * (1.0 - z);
1396 for (sdi = hotBin.begin(); sdi != hotBin.end(); ++sdi) {
1397 vel = (*sdi)->getVel();
1398 vel.x() *= x;
1399 vel.y() *= y;
1400 vel.z() *= z;
1401 (*sdi)->setVel(vel);
1402 }
1403 successfulScale = true;
1404 switch(rnemdFluxType_) {
1405 case rnemdKE :
1406 kineticExchange_ += kineticTarget_;
1407 break;
1408 case rnemdPx :
1409 case rnemdPy :
1410 case rnemdPz :
1411 momentumExchange_ += momentumTarget_;
1412 break;
1413 default :
1414 break;
1415 }
1416 }
1417 }
1418 if (successfulScale != true) {
1419 sprintf(painCave.errMsg,
1420 "RNEMD::doNIVS exchange NOT performed - roots that solve\n"
1421 "\tthe constraint equations may not exist or there may be\n"
1422 "\tno selected objects in one or both slabs.\n");
1423 painCave.isFatal = 0;
1424 painCave.severity = OPENMD_INFO;
1425 simError();
1426 failTrialCount_++;
1427 }
1428 }
1429
1430 void RNEMD::doVSS(SelectionManager& smanA, SelectionManager& smanB) {
1431 if (!doRNEMD_) return;
1432 int selei;
1433 int selej;
1434
1435 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1436 RealType time = currentSnap_->getTime();
1437 Mat3x3d hmat = currentSnap_->getHmat();
1438
1439 StuntDouble* sd;
1440
1441 vector<StuntDouble*> hotBin, coldBin;
1442
1443 Vector3d Ph(V3Zero);
1444 Vector3d Lh(V3Zero);
1445 RealType Mh = 0.0;
1446 Mat3x3d Ih(0.0);
1447 RealType Kh = 0.0;
1448 Vector3d Pc(V3Zero);
1449 Vector3d Lc(V3Zero);
1450 RealType Mc = 0.0;
1451 Mat3x3d Ic(0.0);
1452 RealType Kc = 0.0;
1453
1454 // Constraints can be on only the linear or angular momentum, but
1455 // not both. Usually, the user will specify which they want, but
1456 // in case they don't, the use of periodic boundaries should make
1457 // the choice for us.
1458 bool doLinearPart = false;
1459 bool doAngularPart = false;
1460
1461 switch (rnemdFluxType_) {
1462 case rnemdPx:
1463 case rnemdPy:
1464 case rnemdPz:
1465 case rnemdPvector:
1466 case rnemdKePx:
1467 case rnemdKePy:
1468 case rnemdKePvector:
1469 doLinearPart = true;
1470 break;
1471 case rnemdLx:
1472 case rnemdLy:
1473 case rnemdLz:
1474 case rnemdLvector:
1475 case rnemdKeLx:
1476 case rnemdKeLy:
1477 case rnemdKeLz:
1478 case rnemdKeLvector:
1479 doAngularPart = true;
1480 break;
1481 case rnemdKE:
1482 case rnemdRotKE:
1483 case rnemdFullKE:
1484 default:
1485 if (usePeriodicBoundaryConditions_)
1486 doLinearPart = true;
1487 else
1488 doAngularPart = true;
1489 break;
1490 }
1491
1492 for (sd = smanA.beginSelected(selei); sd != NULL;
1493 sd = smanA.nextSelected(selei)) {
1494
1495 Vector3d pos = sd->getPos();
1496
1497 // wrap the stuntdouble's position back into the box:
1498
1499 if (usePeriodicBoundaryConditions_)
1500 currentSnap_->wrapVector(pos);
1501
1502 RealType mass = sd->getMass();
1503 Vector3d vel = sd->getVel();
1504 Vector3d rPos = sd->getPos() - coordinateOrigin_;
1505 RealType r2;
1506
1507 hotBin.push_back(sd);
1508 Ph += mass * vel;
1509 Mh += mass;
1510 Kh += mass * vel.lengthSquare();
1511 Lh += mass * cross(rPos, vel);
1512 Ih -= outProduct(rPos, rPos) * mass;
1513 r2 = rPos.lengthSquare();
1514 Ih(0, 0) += mass * r2;
1515 Ih(1, 1) += mass * r2;
1516 Ih(2, 2) += mass * r2;
1517
1518 if (rnemdFluxType_ == rnemdFullKE) {
1519 if (sd->isDirectional()) {
1520 Vector3d angMom = sd->getJ();
1521 Mat3x3d I = sd->getI();
1522 if (sd->isLinear()) {
1523 int i = sd->linearAxis();
1524 int j = (i + 1) % 3;
1525 int k = (i + 2) % 3;
1526 Kh += angMom[j] * angMom[j] / I(j, j) +
1527 angMom[k] * angMom[k] / I(k, k);
1528 } else {
1529 Kh += angMom[0] * angMom[0] / I(0, 0) +
1530 angMom[1] * angMom[1] / I(1, 1) +
1531 angMom[2] * angMom[2] / I(2, 2);
1532 }
1533 }
1534 }
1535 }
1536 for (sd = smanB.beginSelected(selej); sd != NULL;
1537 sd = smanB.nextSelected(selej)) {
1538
1539 Vector3d pos = sd->getPos();
1540
1541 // wrap the stuntdouble's position back into the box:
1542
1543 if (usePeriodicBoundaryConditions_)
1544 currentSnap_->wrapVector(pos);
1545
1546 RealType mass = sd->getMass();
1547 Vector3d vel = sd->getVel();
1548 Vector3d rPos = sd->getPos() - coordinateOrigin_;
1549 RealType r2;
1550
1551 coldBin.push_back(sd);
1552 Pc += mass * vel;
1553 Mc += mass;
1554 Kc += mass * vel.lengthSquare();
1555 Lc += mass * cross(rPos, vel);
1556 Ic -= outProduct(rPos, rPos) * mass;
1557 r2 = rPos.lengthSquare();
1558 Ic(0, 0) += mass * r2;
1559 Ic(1, 1) += mass * r2;
1560 Ic(2, 2) += mass * r2;
1561
1562 if (rnemdFluxType_ == rnemdFullKE) {
1563 if (sd->isDirectional()) {
1564 Vector3d angMom = sd->getJ();
1565 Mat3x3d I = sd->getI();
1566 if (sd->isLinear()) {
1567 int i = sd->linearAxis();
1568 int j = (i + 1) % 3;
1569 int k = (i + 2) % 3;
1570 Kc += angMom[j] * angMom[j] / I(j, j) +
1571 angMom[k] * angMom[k] / I(k, k);
1572 } else {
1573 Kc += angMom[0] * angMom[0] / I(0, 0) +
1574 angMom[1] * angMom[1] / I(1, 1) +
1575 angMom[2] * angMom[2] / I(2, 2);
1576 }
1577 }
1578 }
1579 }
1580
1581 Kh *= 0.5;
1582 Kc *= 0.5;
1583
1584 #ifdef IS_MPI
1585 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Ph[0], 3, MPI::REALTYPE, MPI::SUM);
1586 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pc[0], 3, MPI::REALTYPE, MPI::SUM);
1587 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Lh[0], 3, MPI::REALTYPE, MPI::SUM);
1588 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Lc[0], 3, MPI::REALTYPE, MPI::SUM);
1589 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Mh, 1, MPI::REALTYPE, MPI::SUM);
1590 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kh, 1, MPI::REALTYPE, MPI::SUM);
1591 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Mc, 1, MPI::REALTYPE, MPI::SUM);
1592 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kc, 1, MPI::REALTYPE, MPI::SUM);
1593 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, Ih.getArrayPointer(), 9,
1594 MPI::REALTYPE, MPI::SUM);
1595 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, Ic.getArrayPointer(), 9,
1596 MPI::REALTYPE, MPI::SUM);
1597 #endif
1598
1599
1600 Vector3d ac, acrec, bc, bcrec;
1601 Vector3d ah, ahrec, bh, bhrec;
1602
1603 bool successfulExchange = false;
1604 if ((Mh > 0.0) && (Mc > 0.0)) {//both slabs are not empty
1605 Vector3d vc = Pc / Mc;
1606 ac = -momentumTarget_ / Mc + vc;
1607 acrec = -momentumTarget_ / Mc;
1608
1609 // We now need the inverse of the inertia tensor to calculate the
1610 // angular velocity of the cold slab;
1611 Mat3x3d Ici = Ic.inverse();
1612 Vector3d omegac = Ici * Lc;
1613 bc = -(Ici * angularMomentumTarget_) + omegac;
1614 bcrec = bc - omegac;
1615
1616 RealType cNumerator = Kc - kineticTarget_;
1617 if (doLinearPart)
1618 cNumerator -= 0.5 * Mc * ac.lengthSquare();
1619
1620 if (doAngularPart)
1621 cNumerator -= 0.5 * ( dot(bc, Ic * bc));
1622
1623 if (cNumerator > 0.0) {
1624
1625 RealType cDenominator = Kc;
1626
1627 if (doLinearPart)
1628 cDenominator -= 0.5 * Mc * vc.lengthSquare();
1629
1630 if (doAngularPart)
1631 cDenominator -= 0.5*(dot(omegac, Ic * omegac));
1632
1633 if (cDenominator > 0.0) {
1634 RealType c = sqrt(cNumerator / cDenominator);
1635 if ((c > 0.9) && (c < 1.1)) {//restrict scaling coefficients
1636
1637 Vector3d vh = Ph / Mh;
1638 ah = momentumTarget_ / Mh + vh;
1639 ahrec = momentumTarget_ / Mh;
1640
1641 // We now need the inverse of the inertia tensor to
1642 // calculate the angular velocity of the hot slab;
1643 Mat3x3d Ihi = Ih.inverse();
1644 Vector3d omegah = Ihi * Lh;
1645 bh = (Ihi * angularMomentumTarget_) + omegah;
1646 bhrec = bh - omegah;
1647
1648 RealType hNumerator = Kh + kineticTarget_;
1649 if (doLinearPart)
1650 hNumerator -= 0.5 * Mh * ah.lengthSquare();
1651
1652 if (doAngularPart)
1653 hNumerator -= 0.5 * ( dot(bh, Ih * bh));
1654
1655 if (hNumerator > 0.0) {
1656
1657 RealType hDenominator = Kh;
1658 if (doLinearPart)
1659 hDenominator -= 0.5 * Mh * vh.lengthSquare();
1660 if (doAngularPart)
1661 hDenominator -= 0.5*(dot(omegah, Ih * omegah));
1662
1663 if (hDenominator > 0.0) {
1664 RealType h = sqrt(hNumerator / hDenominator);
1665 if ((h > 0.9) && (h < 1.1)) {
1666
1667 vector<StuntDouble*>::iterator sdi;
1668 Vector3d vel;
1669 Vector3d rPos;
1670
1671 for (sdi = coldBin.begin(); sdi != coldBin.end(); ++sdi) {
1672 //vel = (*sdi)->getVel();
1673 rPos = (*sdi)->getPos() - coordinateOrigin_;
1674 if (doLinearPart)
1675 vel = ((*sdi)->getVel() - vc) * c + ac;
1676 if (doAngularPart)
1677 vel = ((*sdi)->getVel() - cross(omegac, rPos)) * c + cross(bc, rPos);
1678
1679 (*sdi)->setVel(vel);
1680 if (rnemdFluxType_ == rnemdFullKE) {
1681 if ((*sdi)->isDirectional()) {
1682 Vector3d angMom = (*sdi)->getJ() * c;
1683 (*sdi)->setJ(angMom);
1684 }
1685 }
1686 }
1687 for (sdi = hotBin.begin(); sdi != hotBin.end(); ++sdi) {
1688 //vel = (*sdi)->getVel();
1689 rPos = (*sdi)->getPos() - coordinateOrigin_;
1690 if (doLinearPart)
1691 vel = ((*sdi)->getVel() - vh) * h + ah;
1692 if (doAngularPart)
1693 vel = ((*sdi)->getVel() - cross(omegah, rPos)) * h + cross(bh, rPos);
1694
1695 (*sdi)->setVel(vel);
1696 if (rnemdFluxType_ == rnemdFullKE) {
1697 if ((*sdi)->isDirectional()) {
1698 Vector3d angMom = (*sdi)->getJ() * h;
1699 (*sdi)->setJ(angMom);
1700 }
1701 }
1702 }
1703 successfulExchange = true;
1704 kineticExchange_ += kineticTarget_;
1705 momentumExchange_ += momentumTarget_;
1706 angularMomentumExchange_ += angularMomentumTarget_;
1707 }
1708 }
1709 }
1710 }
1711 }
1712 }
1713 }
1714 if (successfulExchange != true) {
1715 sprintf(painCave.errMsg,
1716 "RNEMD::doVSS exchange NOT performed - roots that solve\n"
1717 "\tthe constraint equations may not exist or there may be\n"
1718 "\tno selected objects in one or both slabs.\n");
1719 painCave.isFatal = 0;
1720 painCave.severity = OPENMD_INFO;
1721 simError();
1722 failTrialCount_++;
1723 }
1724 }
1725
1726 RealType RNEMD::getDividingArea() {
1727
1728 if (hasDividingArea_) return dividingArea_;
1729
1730 RealType areaA, areaB;
1731 Snapshot* snap = info_->getSnapshotManager()->getCurrentSnapshot();
1732
1733 if (hasSelectionA_) {
1734
1735 if (evaluatorA_.hasSurfaceArea())
1736 areaA = evaluatorA_.getSurfaceArea();
1737 else {
1738
1739 cerr << "selection A did not have surface area, recomputing\n";
1740 int isd;
1741 StuntDouble* sd;
1742 vector<StuntDouble*> aSites;
1743 seleManA_.setSelectionSet(evaluatorA_.evaluate());
1744 for (sd = seleManA_.beginSelected(isd); sd != NULL;
1745 sd = seleManA_.nextSelected(isd)) {
1746 aSites.push_back(sd);
1747 }
1748 #if defined(HAVE_QHULL)
1749 ConvexHull* surfaceMeshA = new ConvexHull();
1750 surfaceMeshA->computeHull(aSites);
1751 areaA = surfaceMeshA->getArea();
1752 delete surfaceMeshA;
1753 #else
1754 sprintf( painCave.errMsg,
1755 "RNEMD::getDividingArea : Hull calculation is not possible\n"
1756 "\twithout libqhull. Please rebuild OpenMD with qhull enabled.");
1757 painCave.severity = OPENMD_ERROR;
1758 painCave.isFatal = 1;
1759 simError();
1760 #endif
1761 }
1762
1763 } else {
1764 if (usePeriodicBoundaryConditions_) {
1765 // in periodic boundaries, the surface area is twice the x-y
1766 // area of the current box:
1767 areaA = 2.0 * snap->getXYarea();
1768 } else {
1769 // in non-periodic simulations, without explicitly setting
1770 // selections, the sphere radius sets the surface area of the
1771 // dividing surface:
1772 areaA = 4.0 * M_PI * pow(sphereARadius_, 2);
1773 }
1774 }
1775
1776 if (hasSelectionB_) {
1777 if (evaluatorB_.hasSurfaceArea())
1778 areaB = evaluatorB_.getSurfaceArea();
1779 else {
1780 cerr << "selection B did not have surface area, recomputing\n";
1781
1782 int isd;
1783 StuntDouble* sd;
1784 vector<StuntDouble*> bSites;
1785 seleManB_.setSelectionSet(evaluatorB_.evaluate());
1786 for (sd = seleManB_.beginSelected(isd); sd != NULL;
1787 sd = seleManB_.nextSelected(isd)) {
1788 bSites.push_back(sd);
1789 }
1790
1791 #if defined(HAVE_QHULL)
1792 ConvexHull* surfaceMeshB = new ConvexHull();
1793 surfaceMeshB->computeHull(bSites);
1794 areaB = surfaceMeshB->getArea();
1795 delete surfaceMeshB;
1796 #else
1797 sprintf( painCave.errMsg,
1798 "RNEMD::getDividingArea : Hull calculation is not possible\n"
1799 "\twithout libqhull. Please rebuild OpenMD with qhull enabled.");
1800 painCave.severity = OPENMD_ERROR;
1801 painCave.isFatal = 1;
1802 simError();
1803 #endif
1804 }
1805
1806 } else {
1807 if (usePeriodicBoundaryConditions_) {
1808 // in periodic boundaries, the surface area is twice the x-y
1809 // area of the current box:
1810 areaB = 2.0 * snap->getXYarea();
1811 } else {
1812 // in non-periodic simulations, without explicitly setting
1813 // selections, but if a sphereBradius has been set, just use that:
1814 areaB = 4.0 * M_PI * pow(sphereBRadius_, 2);
1815 }
1816 }
1817
1818 dividingArea_ = min(areaA, areaB);
1819 hasDividingArea_ = true;
1820 return dividingArea_;
1821 }
1822
1823 void RNEMD::doRNEMD() {
1824 if (!doRNEMD_) return;
1825 trialCount_++;
1826
1827 // object evaluator:
1828 evaluator_.loadScriptString(rnemdObjectSelection_);
1829 seleMan_.setSelectionSet(evaluator_.evaluate());
1830
1831 evaluatorA_.loadScriptString(selectionA_);
1832 evaluatorB_.loadScriptString(selectionB_);
1833
1834 seleManA_.setSelectionSet(evaluatorA_.evaluate());
1835 seleManB_.setSelectionSet(evaluatorB_.evaluate());
1836
1837 commonA_ = seleManA_ & seleMan_;
1838 commonB_ = seleManB_ & seleMan_;
1839
1840 // Target exchange quantities (in each exchange) = dividingArea * dt * flux
1841 // dt = exchange time interval
1842 // flux = target flux
1843 // dividingArea = smallest dividing surface between the two regions
1844
1845 hasDividingArea_ = false;
1846 RealType area = getDividingArea();
1847
1848 kineticTarget_ = kineticFlux_ * exchangeTime_ * area;
1849 momentumTarget_ = momentumFluxVector_ * exchangeTime_ * area;
1850 angularMomentumTarget_ = angularMomentumFluxVector_ * exchangeTime_ * area;
1851
1852 switch(rnemdMethod_) {
1853 case rnemdSwap:
1854 doSwap(commonA_, commonB_);
1855 break;
1856 case rnemdNIVS:
1857 doNIVS(commonA_, commonB_);
1858 break;
1859 case rnemdVSS:
1860 doVSS(commonA_, commonB_);
1861 break;
1862 case rnemdUnkownMethod:
1863 default :
1864 break;
1865 }
1866 }
1867
1868 void RNEMD::collectData() {
1869 if (!doRNEMD_) return;
1870 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1871
1872 // collectData can be called more frequently than the doRNEMD, so use the
1873 // computed area from the last exchange time:
1874 RealType area = getDividingArea();
1875 areaAccumulator_->add(area);
1876 Mat3x3d hmat = currentSnap_->getHmat();
1877 seleMan_.setSelectionSet(evaluator_.evaluate());
1878
1879 int selei(0);
1880 StuntDouble* sd;
1881 int binNo;
1882
1883 vector<RealType> binMass(nBins_, 0.0);
1884 vector<RealType> binPx(nBins_, 0.0);
1885 vector<RealType> binPy(nBins_, 0.0);
1886 vector<RealType> binPz(nBins_, 0.0);
1887 vector<RealType> binOmegax(nBins_, 0.0);
1888 vector<RealType> binOmegay(nBins_, 0.0);
1889 vector<RealType> binOmegaz(nBins_, 0.0);
1890 vector<RealType> binKE(nBins_, 0.0);
1891 vector<int> binDOF(nBins_, 0);
1892 vector<int> binCount(nBins_, 0);
1893
1894 // alternative approach, track all molecules instead of only those
1895 // selected for scaling/swapping:
1896 /*
1897 SimInfo::MoleculeIterator miter;
1898 vector<StuntDouble*>::iterator iiter;
1899 Molecule* mol;
1900 StuntDouble* sd;
1901 for (mol = info_->beginMolecule(miter); mol != NULL;
1902 mol = info_->nextMolecule(miter))
1903 sd is essentially sd
1904 for (sd = mol->beginIntegrableObject(iiter);
1905 sd != NULL;
1906 sd = mol->nextIntegrableObject(iiter))
1907 */
1908
1909 for (sd = seleMan_.beginSelected(selei); sd != NULL;
1910 sd = seleMan_.nextSelected(selei)) {
1911
1912 Vector3d pos = sd->getPos();
1913
1914 // wrap the stuntdouble's position back into the box:
1915
1916 if (usePeriodicBoundaryConditions_) {
1917 currentSnap_->wrapVector(pos);
1918 // which bin is this stuntdouble in?
1919 // wrapped positions are in the range [-0.5*hmat(2,2), +0.5*hmat(2,2)]
1920 // Shift molecules by half a box to have bins start at 0
1921 // The modulo operator is used to wrap the case when we are
1922 // beyond the end of the bins back to the beginning.
1923 binNo = int(nBins_ * (pos.z() / hmat(2,2) + 0.5)) % nBins_;
1924 } else {
1925 Vector3d rPos = pos - coordinateOrigin_;
1926 binNo = int(rPos.length() / binWidth_);
1927 }
1928
1929 RealType mass = sd->getMass();
1930 Vector3d vel = sd->getVel();
1931 Vector3d rPos = sd->getPos() - coordinateOrigin_;
1932 Vector3d aVel = cross(rPos, vel);
1933
1934 if (binNo >= 0 && binNo < nBins_) {
1935 binCount[binNo]++;
1936 binMass[binNo] += mass;
1937 binPx[binNo] += mass*vel.x();
1938 binPy[binNo] += mass*vel.y();
1939 binPz[binNo] += mass*vel.z();
1940 binOmegax[binNo] += aVel.x();
1941 binOmegay[binNo] += aVel.y();
1942 binOmegaz[binNo] += aVel.z();
1943 binKE[binNo] += 0.5 * (mass * vel.lengthSquare());
1944 binDOF[binNo] += 3;
1945
1946 if (sd->isDirectional()) {
1947 Vector3d angMom = sd->getJ();
1948 Mat3x3d I = sd->getI();
1949 if (sd->isLinear()) {
1950 int i = sd->linearAxis();
1951 int j = (i + 1) % 3;
1952 int k = (i + 2) % 3;
1953 binKE[binNo] += 0.5 * (angMom[j] * angMom[j] / I(j, j) +
1954 angMom[k] * angMom[k] / I(k, k));
1955 binDOF[binNo] += 2;
1956 } else {
1957 binKE[binNo] += 0.5 * (angMom[0] * angMom[0] / I(0, 0) +
1958 angMom[1] * angMom[1] / I(1, 1) +
1959 angMom[2] * angMom[2] / I(2, 2));
1960 binDOF[binNo] += 3;
1961 }
1962 }
1963 }
1964 }
1965
1966 #ifdef IS_MPI
1967 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binCount[0],
1968 nBins_, MPI::INT, MPI::SUM);
1969 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binMass[0],
1970 nBins_, MPI::REALTYPE, MPI::SUM);
1971 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPx[0],
1972 nBins_, MPI::REALTYPE, MPI::SUM);
1973 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPy[0],
1974 nBins_, MPI::REALTYPE, MPI::SUM);
1975 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPz[0],
1976 nBins_, MPI::REALTYPE, MPI::SUM);
1977 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binOmegax[0],
1978 nBins_, MPI::REALTYPE, MPI::SUM);
1979 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binOmegay[0],
1980 nBins_, MPI::REALTYPE, MPI::SUM);
1981 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binOmegaz[0],
1982 nBins_, MPI::REALTYPE, MPI::SUM);
1983 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binKE[0],
1984 nBins_, MPI::REALTYPE, MPI::SUM);
1985 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binDOF[0],
1986 nBins_, MPI::INT, MPI::SUM);
1987 #endif
1988
1989 Vector3d vel;
1990 Vector3d aVel;
1991 RealType den;
1992 RealType temp;
1993 RealType z;
1994 RealType r;
1995 for (int i = 0; i < nBins_; i++) {
1996 if (usePeriodicBoundaryConditions_) {
1997 z = (((RealType)i + 0.5) / (RealType)nBins_) * hmat(2,2);
1998 den = binMass[i] * nBins_ * PhysicalConstants::densityConvert
1999 / currentSnap_->getVolume() ;
2000 } else {
2001 r = (((RealType)i + 0.5) * binWidth_);
2002 RealType rinner = (RealType)i * binWidth_;
2003 RealType router = (RealType)(i+1) * binWidth_;
2004 den = binMass[i] * 3.0 * PhysicalConstants::densityConvert
2005 / (4.0 * M_PI * (pow(router,3) - pow(rinner,3)));
2006 }
2007 vel.x() = binPx[i] / binMass[i];
2008 vel.y() = binPy[i] / binMass[i];
2009 vel.z() = binPz[i] / binMass[i];
2010 aVel.x() = binOmegax[i] / binCount[i];
2011 aVel.y() = binOmegay[i] / binCount[i];
2012 aVel.z() = binOmegaz[i] / binCount[i];
2013
2014 if (binCount[i] > 0) {
2015 // only add values if there are things to add
2016 temp = 2.0 * binKE[i] / (binDOF[i] * PhysicalConstants::kb *
2017 PhysicalConstants::energyConvert);
2018
2019 for (unsigned int j = 0; j < outputMask_.size(); ++j) {
2020 if(outputMask_[j]) {
2021 switch(j) {
2022 case Z:
2023 dynamic_cast<Accumulator *>(data_[j].accumulator[i])->add(z);
2024 break;
2025 case R:
2026 dynamic_cast<Accumulator *>(data_[j].accumulator[i])->add(r);
2027 break;
2028 case TEMPERATURE:
2029 dynamic_cast<Accumulator *>(data_[j].accumulator[i])->add(temp);
2030 break;
2031 case VELOCITY:
2032 dynamic_cast<VectorAccumulator *>(data_[j].accumulator[i])->add(vel);
2033 break;
2034 case ANGULARVELOCITY:
2035 dynamic_cast<VectorAccumulator *>(data_[j].accumulator[i])->add(aVel);
2036 break;
2037 case DENSITY:
2038 dynamic_cast<Accumulator *>(data_[j].accumulator[i])->add(den);
2039 break;
2040 }
2041 }
2042 }
2043 }
2044 }
2045 hasData_ = true;
2046 }
2047
2048 void RNEMD::getStarted() {
2049 if (!doRNEMD_) return;
2050 hasDividingArea_ = false;
2051 collectData();
2052 writeOutputFile();
2053 }
2054
2055 void RNEMD::parseOutputFileFormat(const std::string& format) {
2056 if (!doRNEMD_) return;
2057 StringTokenizer tokenizer(format, " ,;|\t\n\r");
2058
2059 while(tokenizer.hasMoreTokens()) {
2060 std::string token(tokenizer.nextToken());
2061 toUpper(token);
2062 OutputMapType::iterator i = outputMap_.find(token);
2063 if (i != outputMap_.end()) {
2064 outputMask_.set(i->second);
2065 } else {
2066 sprintf( painCave.errMsg,
2067 "RNEMD::parseOutputFileFormat: %s is not a recognized\n"
2068 "\toutputFileFormat keyword.\n", token.c_str() );
2069 painCave.isFatal = 0;
2070 painCave.severity = OPENMD_ERROR;
2071 simError();
2072 }
2073 }
2074 }
2075
2076 void RNEMD::writeOutputFile() {
2077 if (!doRNEMD_) return;
2078 if (!hasData_) return;
2079
2080 #ifdef IS_MPI
2081 // If we're the root node, should we print out the results
2082 int worldRank = MPI::COMM_WORLD.Get_rank();
2083 if (worldRank == 0) {
2084 #endif
2085 rnemdFile_.open(rnemdFileName_.c_str(), std::ios::out | std::ios::trunc );
2086
2087 if( !rnemdFile_ ){
2088 sprintf( painCave.errMsg,
2089 "Could not open \"%s\" for RNEMD output.\n",
2090 rnemdFileName_.c_str());
2091 painCave.isFatal = 1;
2092 simError();
2093 }
2094
2095 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
2096
2097 RealType time = currentSnap_->getTime();
2098 RealType avgArea;
2099 areaAccumulator_->getAverage(avgArea);
2100
2101 RealType Jz(0.0);
2102 Vector3d JzP(V3Zero);
2103 Vector3d JzL(V3Zero);
2104 if (time >= info_->getSimParams()->getDt()) {
2105 Jz = kineticExchange_ / (time * avgArea)
2106 / PhysicalConstants::energyConvert;
2107 JzP = momentumExchange_ / (time * avgArea);
2108 JzL = angularMomentumExchange_ / (time * avgArea);
2109 }
2110
2111 rnemdFile_ << "#######################################################\n";
2112 rnemdFile_ << "# RNEMD {\n";
2113
2114 map<string, RNEMDMethod>::iterator mi;
2115 for(mi = stringToMethod_.begin(); mi != stringToMethod_.end(); ++mi) {
2116 if ( (*mi).second == rnemdMethod_)
2117 rnemdFile_ << "# exchangeMethod = \"" << (*mi).first << "\";\n";
2118 }
2119 map<string, RNEMDFluxType>::iterator fi;
2120 for(fi = stringToFluxType_.begin(); fi != stringToFluxType_.end(); ++fi) {
2121 if ( (*fi).second == rnemdFluxType_)
2122 rnemdFile_ << "# fluxType = \"" << (*fi).first << "\";\n";
2123 }
2124
2125 rnemdFile_ << "# exchangeTime = " << exchangeTime_ << ";\n";
2126
2127 rnemdFile_ << "# objectSelection = \""
2128 << rnemdObjectSelection_ << "\";\n";
2129 rnemdFile_ << "# selectionA = \"" << selectionA_ << "\";\n";
2130 rnemdFile_ << "# selectionB = \"" << selectionB_ << "\";\n";
2131 rnemdFile_ << "# }\n";
2132 rnemdFile_ << "#######################################################\n";
2133 rnemdFile_ << "# RNEMD report:\n";
2134 rnemdFile_ << "# running time = " << time << " fs\n";
2135 rnemdFile_ << "# Target flux:\n";
2136 rnemdFile_ << "# kinetic = "
2137 << kineticFlux_ / PhysicalConstants::energyConvert
2138 << " (kcal/mol/A^2/fs)\n";
2139 rnemdFile_ << "# momentum = " << momentumFluxVector_
2140 << " (amu/A/fs^2)\n";
2141 rnemdFile_ << "# angular momentum = " << angularMomentumFluxVector_
2142 << " (amu/A^2/fs^2)\n";
2143 rnemdFile_ << "# Target one-time exchanges:\n";
2144 rnemdFile_ << "# kinetic = "
2145 << kineticTarget_ / PhysicalConstants::energyConvert
2146 << " (kcal/mol)\n";
2147 rnemdFile_ << "# momentum = " << momentumTarget_
2148 << " (amu*A/fs)\n";
2149 rnemdFile_ << "# angular momentum = " << angularMomentumTarget_
2150 << " (amu*A^2/fs)\n";
2151 rnemdFile_ << "# Actual exchange totals:\n";
2152 rnemdFile_ << "# kinetic = "
2153 << kineticExchange_ / PhysicalConstants::energyConvert
2154 << " (kcal/mol)\n";
2155 rnemdFile_ << "# momentum = " << momentumExchange_
2156 << " (amu*A/fs)\n";
2157 rnemdFile_ << "# angular momentum = " << angularMomentumExchange_
2158 << " (amu*A^2/fs)\n";
2159 rnemdFile_ << "# Actual flux:\n";
2160 rnemdFile_ << "# kinetic = " << Jz
2161 << " (kcal/mol/A^2/fs)\n";
2162 rnemdFile_ << "# momentum = " << JzP
2163 << " (amu/A/fs^2)\n";
2164 rnemdFile_ << "# angular momentum = " << JzL
2165 << " (amu/A^2/fs^2)\n";
2166 rnemdFile_ << "# Exchange statistics:\n";
2167 rnemdFile_ << "# attempted = " << trialCount_ << "\n";
2168 rnemdFile_ << "# failed = " << failTrialCount_ << "\n";
2169 if (rnemdMethod_ == rnemdNIVS) {
2170 rnemdFile_ << "# NIVS root-check errors = "
2171 << failRootCount_ << "\n";
2172 }
2173 rnemdFile_ << "#######################################################\n";
2174
2175
2176
2177 //write title
2178 rnemdFile_ << "#";
2179 for (unsigned int i = 0; i < outputMask_.size(); ++i) {
2180 if (outputMask_[i]) {
2181 rnemdFile_ << "\t" << data_[i].title <<
2182 "(" << data_[i].units << ")";
2183 // add some extra tabs for column alignment
2184 if (data_[i].dataType == "Vector3d") rnemdFile_ << "\t\t";
2185 }
2186 }
2187 rnemdFile_ << std::endl;
2188
2189 rnemdFile_.precision(8);
2190
2191 for (int j = 0; j < nBins_; j++) {
2192
2193 for (unsigned int i = 0; i < outputMask_.size(); ++i) {
2194 if (outputMask_[i]) {
2195 if (data_[i].dataType == "RealType")
2196 writeReal(i,j);
2197 else if (data_[i].dataType == "Vector3d")
2198 writeVector(i,j);
2199 else {
2200 sprintf( painCave.errMsg,
2201 "RNEMD found an unknown data type for: %s ",
2202 data_[i].title.c_str());
2203 painCave.isFatal = 1;
2204 simError();
2205 }
2206 }
2207 }
2208 rnemdFile_ << std::endl;
2209
2210 }
2211
2212 rnemdFile_ << "#######################################################\n";
2213 rnemdFile_ << "# Standard Deviations in those quantities follow:\n";
2214 rnemdFile_ << "#######################################################\n";
2215
2216
2217 for (int j = 0; j < nBins_; j++) {
2218 rnemdFile_ << "#";
2219 for (unsigned int i = 0; i < outputMask_.size(); ++i) {
2220 if (outputMask_[i]) {
2221 if (data_[i].dataType == "RealType")
2222 writeRealStdDev(i,j);
2223 else if (data_[i].dataType == "Vector3d")
2224 writeVectorStdDev(i,j);
2225 else {
2226 sprintf( painCave.errMsg,
2227 "RNEMD found an unknown data type for: %s ",
2228 data_[i].title.c_str());
2229 painCave.isFatal = 1;
2230 simError();
2231 }
2232 }
2233 }
2234 rnemdFile_ << std::endl;
2235
2236 }
2237
2238 rnemdFile_.flush();
2239 rnemdFile_.close();
2240
2241 #ifdef IS_MPI
2242 }
2243 #endif
2244
2245 }
2246
2247 void RNEMD::writeReal(int index, unsigned int bin) {
2248 if (!doRNEMD_) return;
2249 assert(index >=0 && index < ENDINDEX);
2250 assert(int(bin) < nBins_);
2251 RealType s;
2252 int count;
2253
2254 count = data_[index].accumulator[bin]->count();
2255 if (count == 0) return;
2256
2257 dynamic_cast<Accumulator *>(data_[index].accumulator[bin])->getAverage(s);
2258
2259 if (! isinf(s) && ! isnan(s)) {
2260 rnemdFile_ << "\t" << s;
2261 } else{
2262 sprintf( painCave.errMsg,
2263 "RNEMD detected a numerical error writing: %s for bin %u",
2264 data_[index].title.c_str(), bin);
2265 painCave.isFatal = 1;
2266 simError();
2267 }
2268 }
2269
2270 void RNEMD::writeVector(int index, unsigned int bin) {
2271 if (!doRNEMD_) return;
2272 assert(index >=0 && index < ENDINDEX);
2273 assert(int(bin) < nBins_);
2274 Vector3d s;
2275 int count;
2276
2277 count = data_[index].accumulator[bin]->count();
2278
2279 if (count == 0) return;
2280
2281 dynamic_cast<VectorAccumulator*>(data_[index].accumulator[bin])->getAverage(s);
2282 if (isinf(s[0]) || isnan(s[0]) ||
2283 isinf(s[1]) || isnan(s[1]) ||
2284 isinf(s[2]) || isnan(s[2]) ) {
2285 sprintf( painCave.errMsg,
2286 "RNEMD detected a numerical error writing: %s for bin %u",
2287 data_[index].title.c_str(), bin);
2288 painCave.isFatal = 1;
2289 simError();
2290 } else {
2291 rnemdFile_ << "\t" << s[0] << "\t" << s[1] << "\t" << s[2];
2292 }
2293 }
2294
2295 void RNEMD::writeRealStdDev(int index, unsigned int bin) {
2296 if (!doRNEMD_) return;
2297 assert(index >=0 && index < ENDINDEX);
2298 assert(int(bin) < nBins_);
2299 RealType s;
2300 int count;
2301
2302 count = data_[index].accumulator[bin]->count();
2303 if (count == 0) return;
2304
2305 dynamic_cast<Accumulator *>(data_[index].accumulator[bin])->getStdDev(s);
2306
2307 if (! isinf(s) && ! isnan(s)) {
2308 rnemdFile_ << "\t" << s;
2309 } else{
2310 sprintf( painCave.errMsg,
2311 "RNEMD detected a numerical error writing: %s std. dev. for bin %u",
2312 data_[index].title.c_str(), bin);
2313 painCave.isFatal = 1;
2314 simError();
2315 }
2316 }
2317
2318 void RNEMD::writeVectorStdDev(int index, unsigned int bin) {
2319 if (!doRNEMD_) return;
2320 assert(index >=0 && index < ENDINDEX);
2321 assert(int(bin) < nBins_);
2322 Vector3d s;
2323 int count;
2324
2325 count = data_[index].accumulator[bin]->count();
2326 if (count == 0) return;
2327
2328 dynamic_cast<VectorAccumulator*>(data_[index].accumulator[bin])->getStdDev(s);
2329 if (isinf(s[0]) || isnan(s[0]) ||
2330 isinf(s[1]) || isnan(s[1]) ||
2331 isinf(s[2]) || isnan(s[2]) ) {
2332 sprintf( painCave.errMsg,
2333 "RNEMD detected a numerical error writing: %s std. dev. for bin %u",
2334 data_[index].title.c_str(), bin);
2335 painCave.isFatal = 1;
2336 simError();
2337 } else {
2338 rnemdFile_ << "\t" << s[0] << "\t" << s[1] << "\t" << s[2];
2339 }
2340 }
2341 }
2342

Properties

Name Value
svn:keywords Author Id Revision Date