ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/rnemd/RNEMD.cpp
Revision: 1867
Committed: Mon Apr 29 17:53:48 2013 UTC (12 years ago) by gezelter
File size: 74787 byte(s)
Log Message:
Attempt to fix memory leak in Multiple hull calls

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

Properties

Name Value
svn:keywords Author Id Revision Date