ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/rnemd/RNEMD.cpp
Revision: 2027
Committed: Wed Oct 22 14:29:20 2014 UTC (10 years, 6 months ago) by gezelter
File size: 77418 byte(s)
Log Message:
Fixed a closing brace bug

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

Properties

Name Value
svn:keywords Author Id Revision Date