ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/rnemd/RNEMD.cpp
Revision: 1946
Committed: Tue Nov 12 02:18:35 2013 UTC (11 years, 5 months ago) by gezelter
File size: 76600 byte(s)
Log Message:
fixed a bug in the temperature profile calculation

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

Properties

Name Value
svn:keywords Author Id Revision Date