ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/rnemd/RNEMD.cpp
Revision: 2068
Committed: Thu Mar 5 15:40:58 2015 UTC (10 years, 2 months ago) by gezelter
File size: 77351 byte(s)
Log Message:
more g++ warning fixes

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

Properties

Name Value
svn:keywords Author Id Revision Date