ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/rnemd/RNEMD.cpp
Revision: 1793
Committed: Fri Aug 31 21:16:10 2012 UTC (12 years, 8 months ago) by gezelter
File size: 57017 byte(s)
Log Message:
Cleaning up some warning messages on linux

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

Properties

Name Value
svn:keywords Author Id Revision Date