ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/rnemd/RNEMD.cpp
Revision: 1875
Committed: Fri May 17 14:41:42 2013 UTC (11 years, 11 months ago) by gezelter
File size: 75074 byte(s)
Log Message:
Fixed a bunch of stylistic and performance issues discovered via cppcheck.

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date