ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/rnemd/RNEMD.cpp
Revision: 1803
Committed: Wed Oct 3 14:20:07 2012 UTC (12 years, 7 months ago) by gezelter
File size: 56472 byte(s)
Log Message:
Merging trunk changes back to development branch

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

Properties

Name Value
svn:keywords Author Id Revision Date