ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/rnemd/RNEMD.cpp
Revision: 1856
Committed: Tue Apr 2 21:30:34 2013 UTC (12 years ago) by gezelter
File size: 73072 byte(s)
Log Message:
Bug fixes to get new parallel decomposition methods working without periodic boundary conditions.  
Updated volume in stats file correctly.
LangevinHull fixes for low viscosity values.

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

Properties

Name Value
svn:keywords Author Id Revision Date