ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/rnemd/RNEMD.cpp
Revision: 1776
Committed: Thu Aug 9 15:52:59 2012 UTC (12 years, 8 months ago) by gezelter
File size: 56503 byte(s)
Log Message:
Fixes to mdParser to handle vector assignments, fixes for VelocityVerletIntegrator deleting rnemd_ when it doesn't exist yet.

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date