ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/rnemd/RNEMD.cpp
Revision: 1830
Committed: Wed Jan 9 22:02:30 2013 UTC (12 years, 3 months ago) by gezelter
File size: 56549 byte(s)
Log Message:
Windows compilation fixes

File Contents

# User Rev Content
1 gezelter 1329 /*
2     * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
4     * The University of Notre Dame grants you ("Licensee") a
5     * non-exclusive, royalty free, license to use, modify and
6     * redistribute this software in source and binary code form, provided
7     * that the following conditions are met:
8     *
9 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 1329 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 1329 * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31 gezelter 1390 *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     * [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 1813 #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 gezelter 1803
285 skuang 1341 int nIntegrable = info->getNGlobalIntegrableObjects();
286    
287     if (selectionCount > nIntegrable) {
288     sprintf(painCave.errMsg,
289 gezelter 1773 "RNEMD: The current objectSelection,\n"
290 skuang 1341 "\t\t%s\n"
291     "\thas resulted in %d selected objects. However,\n"
292     "\tthe total number of integrable objects in the system\n"
293     "\tis only %d. This is almost certainly not what you want\n"
294     "\tto do. A likely cause of this is forgetting the _RB_0\n"
295     "\tselector in the selection script!\n",
296     rnemdObjectSelection_.c_str(),
297     selectionCount, nIntegrable);
298     painCave.isFatal = 0;
299 gezelter 1629 painCave.severity = OPENMD_WARNING;
300 skuang 1341 simError();
301     }
302 skuang 1330
303 gezelter 1774 areaAccumulator_ = new Accumulator();
304    
305 gezelter 1773 nBins_ = rnemdParams->getOutputBins();
306 skuang 1330
307 gezelter 1773 data_.resize(RNEMD::ENDINDEX);
308     OutputData z;
309     z.units = "Angstroms";
310     z.title = "Z";
311     z.dataType = "RealType";
312     z.accumulator.reserve(nBins_);
313 gezelter 1813 for (int i = 0; i < nBins_; i++)
314 gezelter 1773 z.accumulator.push_back( new Accumulator() );
315     data_[Z] = z;
316     outputMap_["Z"] = Z;
317 skuang 1368
318 gezelter 1773 OutputData temperature;
319     temperature.units = "K";
320     temperature.title = "Temperature";
321     temperature.dataType = "RealType";
322     temperature.accumulator.reserve(nBins_);
323 gezelter 1813 for (int i = 0; i < nBins_; i++)
324 gezelter 1773 temperature.accumulator.push_back( new Accumulator() );
325     data_[TEMPERATURE] = temperature;
326     outputMap_["TEMPERATURE"] = TEMPERATURE;
327 gezelter 1722
328 gezelter 1773 OutputData velocity;
329 gezelter 1777 velocity.units = "angstroms/fs";
330 gezelter 1773 velocity.title = "Velocity";
331     velocity.dataType = "Vector3d";
332     velocity.accumulator.reserve(nBins_);
333 gezelter 1813 for (int i = 0; i < nBins_; i++)
334 gezelter 1773 velocity.accumulator.push_back( new VectorAccumulator() );
335     data_[VELOCITY] = velocity;
336     outputMap_["VELOCITY"] = VELOCITY;
337 skuang 1368
338 gezelter 1773 OutputData density;
339     density.units = "g cm^-3";
340     density.title = "Density";
341     density.dataType = "RealType";
342     density.accumulator.reserve(nBins_);
343 gezelter 1813 for (int i = 0; i < nBins_; i++)
344 gezelter 1773 density.accumulator.push_back( new Accumulator() );
345     data_[DENSITY] = density;
346     outputMap_["DENSITY"] = DENSITY;
347 skuang 1368
348 gezelter 1773 if (hasOutputFields) {
349     parseOutputFileFormat(rnemdParams->getOutputFields());
350     } else {
351     outputMask_.set(Z);
352     switch (rnemdFluxType_) {
353     case rnemdKE:
354     case rnemdRotKE:
355     case rnemdFullKE:
356     outputMask_.set(TEMPERATURE);
357     break;
358     case rnemdPx:
359     case rnemdPy:
360     outputMask_.set(VELOCITY);
361     break;
362     case rnemdPz:
363     case rnemdPvector:
364     outputMask_.set(VELOCITY);
365     outputMask_.set(DENSITY);
366     break;
367     case rnemdKePx:
368     case rnemdKePy:
369     outputMask_.set(TEMPERATURE);
370     outputMask_.set(VELOCITY);
371     break;
372     case rnemdKePvector:
373     outputMask_.set(TEMPERATURE);
374     outputMask_.set(VELOCITY);
375     outputMask_.set(DENSITY);
376     break;
377     default:
378     break;
379 gezelter 1629 }
380     }
381 gezelter 1773
382     if (hasOutputFileName) {
383     rnemdFileName_ = rnemdParams->getOutputFileName();
384 skuang 1368 } else {
385 gezelter 1773 rnemdFileName_ = getPrefix(info->getFinalConfigFileName()) + ".rnemd";
386     }
387 gezelter 1722
388 gezelter 1773 exchangeTime_ = rnemdParams->getExchangeTime();
389 skuang 1338
390 gezelter 1773 Snapshot* currentSnap_ = info->getSnapshotManager()->getCurrentSnapshot();
391     Mat3x3d hmat = currentSnap_->getHmat();
392    
393     // Target exchange quantities (in each exchange) = 2 Lx Ly dt flux
394     // Lx, Ly = box dimensions in x & y
395     // dt = exchange time interval
396     // flux = target flux
397 skuang 1368
398 gezelter 1774 RealType area = currentSnap_->getXYarea();
399     kineticTarget_ = 2.0 * kineticFlux_ * exchangeTime_ * area;
400     momentumTarget_ = 2.0 * momentumFluxVector_ * exchangeTime_ * area;
401 gezelter 1773
402     // total exchange sums are zeroed out at the beginning:
403    
404     kineticExchange_ = 0.0;
405     momentumExchange_ = V3Zero;
406    
407     if (hasSlabWidth)
408     slabWidth_ = rnemdParams->getSlabWidth();
409     else
410     slabWidth_ = hmat(2,2) / 10.0;
411    
412     if (hasSlabACenter)
413     slabACenter_ = rnemdParams->getSlabACenter();
414     else
415     slabACenter_ = 0.0;
416    
417     if (hasSlabBCenter)
418     slabBCenter_ = rnemdParams->getSlabBCenter();
419     else
420     slabBCenter_ = hmat(2,2) / 2.0;
421    
422 gezelter 1329 }
423 gezelter 1773
424 gezelter 1329 RNEMD::~RNEMD() {
425 gezelter 1776 if (!doRNEMD_) return;
426 skuang 1368 #ifdef IS_MPI
427     if (worldRank == 0) {
428     #endif
429 gezelter 1722
430 gezelter 1773 writeOutputFile();
431    
432     rnemdFile_.close();
433 jmarr 1728
434 skuang 1368 #ifdef IS_MPI
435     }
436     #endif
437 gezelter 1329 }
438 gezelter 1773
439     bool RNEMD::inSlabA(Vector3d pos) {
440     return (abs(pos.z() - slabACenter_) < 0.5*slabWidth_);
441     }
442     bool RNEMD::inSlabB(Vector3d pos) {
443     return (abs(pos.z() - slabBCenter_) < 0.5*slabWidth_);
444     }
445 skuang 1330
446 gezelter 1329 void RNEMD::doSwap() {
447 gezelter 1776 if (!doRNEMD_) return;
448 gezelter 1332 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
449     Mat3x3d hmat = currentSnap_->getHmat();
450    
451 gezelter 1331 seleMan_.setSelectionSet(evaluator_.evaluate());
452    
453 gezelter 1333 int selei;
454 gezelter 1331 StuntDouble* sd;
455    
456 skuang 1338 RealType min_val;
457     bool min_found = false;
458     StuntDouble* min_sd;
459    
460     RealType max_val;
461     bool max_found = false;
462     StuntDouble* max_sd;
463    
464 gezelter 1333 for (sd = seleMan_.beginSelected(selei); sd != NULL;
465     sd = seleMan_.nextSelected(selei)) {
466 gezelter 1332
467 gezelter 1331 Vector3d pos = sd->getPos();
468 gezelter 1332
469     // wrap the stuntdouble's position back into the box:
470    
471 gezelter 1331 if (usePeriodicBoundaryConditions_)
472 gezelter 1332 currentSnap_->wrapVector(pos);
473 gezelter 1773 bool inA = inSlabA(pos);
474     bool inB = inSlabB(pos);
475 gezelter 1332
476 gezelter 1773 if (inA || inB) {
477 gezelter 1332
478     RealType mass = sd->getMass();
479     Vector3d vel = sd->getVel();
480     RealType value;
481 gezelter 1773
482     switch(rnemdFluxType_) {
483     case rnemdKE :
484 gezelter 1332
485 gezelter 1722 value = mass * vel.lengthSquare();
486    
487     if (sd->isDirectional()) {
488 gezelter 1332 Vector3d angMom = sd->getJ();
489     Mat3x3d I = sd->getI();
490    
491     if (sd->isLinear()) {
492 gezelter 1722 int i = sd->linearAxis();
493     int j = (i + 1) % 3;
494     int k = (i + 2) % 3;
495     value += angMom[j] * angMom[j] / I(j, j) +
496     angMom[k] * angMom[k] / I(k, k);
497 gezelter 1332 } else {
498 gezelter 1722 value += angMom[0]*angMom[0]/I(0, 0)
499     + angMom[1]*angMom[1]/I(1, 1)
500     + angMom[2]*angMom[2]/I(2, 2);
501 gezelter 1332 }
502 gezelter 1722 } //angular momenta exchange enabled
503 skuang 1368 value *= 0.5;
504 gezelter 1332 break;
505     case rnemdPx :
506     value = mass * vel[0];
507     break;
508     case rnemdPy :
509     value = mass * vel[1];
510     break;
511     case rnemdPz :
512     value = mass * vel[2];
513     break;
514     default :
515     break;
516     }
517    
518 gezelter 1773 if (inA == 0) {
519 skuang 1338 if (!min_found) {
520     min_val = value;
521     min_sd = sd;
522     min_found = true;
523     } else {
524     if (min_val > value) {
525     min_val = value;
526     min_sd = sd;
527     }
528     }
529 gezelter 1773 } else {
530 skuang 1338 if (!max_found) {
531     max_val = value;
532     max_sd = sd;
533     max_found = true;
534     } else {
535     if (max_val < value) {
536     max_val = value;
537     max_sd = sd;
538     }
539     }
540     }
541 gezelter 1332 }
542 gezelter 1331 }
543 gezelter 1773
544 gezelter 1813 #ifdef IS_MPI
545     int worldRank = MPI::COMM_WORLD.Get_rank();
546 gezelter 1773
547 gezelter 1350 bool my_min_found = min_found;
548     bool my_max_found = max_found;
549    
550     // Even if we didn't find a minimum, did someone else?
551 gezelter 1629 MPI::COMM_WORLD.Allreduce(&my_min_found, &min_found, 1, MPI::BOOL, MPI::LOR);
552 gezelter 1350 // Even if we didn't find a maximum, did someone else?
553 gezelter 1629 MPI::COMM_WORLD.Allreduce(&my_max_found, &max_found, 1, MPI::BOOL, MPI::LOR);
554 gezelter 1722 #endif
555    
556     if (max_found && min_found) {
557    
558     #ifdef IS_MPI
559     struct {
560     RealType val;
561     int rank;
562     } max_vals, min_vals;
563 jmarr 1728
564 gezelter 1722 if (my_min_found) {
565 gezelter 1350 min_vals.val = min_val;
566 gezelter 1722 } else {
567 gezelter 1350 min_vals.val = HONKING_LARGE_VALUE;
568 gezelter 1722 }
569 gezelter 1350 min_vals.rank = worldRank;
570    
571     // Who had the minimum?
572     MPI::COMM_WORLD.Allreduce(&min_vals, &min_vals,
573     1, MPI::REALTYPE_INT, MPI::MINLOC);
574     min_val = min_vals.val;
575    
576 gezelter 1722 if (my_max_found) {
577 gezelter 1350 max_vals.val = max_val;
578 gezelter 1722 } else {
579 gezelter 1350 max_vals.val = -HONKING_LARGE_VALUE;
580 gezelter 1722 }
581 gezelter 1350 max_vals.rank = worldRank;
582    
583     // Who had the maximum?
584     MPI::COMM_WORLD.Allreduce(&max_vals, &max_vals,
585     1, MPI::REALTYPE_INT, MPI::MAXLOC);
586     max_val = max_vals.val;
587     #endif
588 gezelter 1722
589 gezelter 1629 if (min_val < max_val) {
590 gezelter 1722
591 gezelter 1350 #ifdef IS_MPI
592     if (max_vals.rank == worldRank && min_vals.rank == worldRank) {
593     // I have both maximum and minimum, so proceed like a single
594     // processor version:
595     #endif
596 gezelter 1722
597 gezelter 1350 Vector3d min_vel = min_sd->getVel();
598     Vector3d max_vel = max_sd->getVel();
599     RealType temp_vel;
600    
601 gezelter 1773 switch(rnemdFluxType_) {
602     case rnemdKE :
603 gezelter 1350 min_sd->setVel(max_vel);
604     max_sd->setVel(min_vel);
605 gezelter 1722 if (min_sd->isDirectional() && max_sd->isDirectional()) {
606 gezelter 1350 Vector3d min_angMom = min_sd->getJ();
607     Vector3d max_angMom = max_sd->getJ();
608     min_sd->setJ(max_angMom);
609     max_sd->setJ(min_angMom);
610 gezelter 1722 }//angular momenta exchange enabled
611     //assumes same rigid body identity
612 gezelter 1350 break;
613     case rnemdPx :
614     temp_vel = min_vel.x();
615     min_vel.x() = max_vel.x();
616     max_vel.x() = temp_vel;
617     min_sd->setVel(min_vel);
618     max_sd->setVel(max_vel);
619     break;
620     case rnemdPy :
621     temp_vel = min_vel.y();
622     min_vel.y() = max_vel.y();
623     max_vel.y() = temp_vel;
624     min_sd->setVel(min_vel);
625     max_sd->setVel(max_vel);
626     break;
627     case rnemdPz :
628     temp_vel = min_vel.z();
629     min_vel.z() = max_vel.z();
630     max_vel.z() = temp_vel;
631     min_sd->setVel(min_vel);
632     max_sd->setVel(max_vel);
633     break;
634     default :
635     break;
636     }
637 gezelter 1722
638 gezelter 1350 #ifdef IS_MPI
639     // the rest of the cases only apply in parallel simulations:
640     } else if (max_vals.rank == worldRank) {
641     // I had the max, but not the minimum
642    
643     Vector3d min_vel;
644     Vector3d max_vel = max_sd->getVel();
645     MPI::Status status;
646 skuang 1341
647 gezelter 1350 // point-to-point swap of the velocity vector
648     MPI::COMM_WORLD.Sendrecv(max_vel.getArrayPointer(), 3, MPI::REALTYPE,
649     min_vals.rank, 0,
650     min_vel.getArrayPointer(), 3, MPI::REALTYPE,
651     min_vals.rank, 0, status);
652    
653 gezelter 1773 switch(rnemdFluxType_) {
654     case rnemdKE :
655 gezelter 1350 max_sd->setVel(min_vel);
656 gezelter 1722 //angular momenta exchange enabled
657 gezelter 1350 if (max_sd->isDirectional()) {
658     Vector3d min_angMom;
659     Vector3d max_angMom = max_sd->getJ();
660 gezelter 1629
661 gezelter 1350 // point-to-point swap of the angular momentum vector
662     MPI::COMM_WORLD.Sendrecv(max_angMom.getArrayPointer(), 3,
663     MPI::REALTYPE, min_vals.rank, 1,
664     min_angMom.getArrayPointer(), 3,
665     MPI::REALTYPE, min_vals.rank, 1,
666     status);
667 gezelter 1629
668 gezelter 1350 max_sd->setJ(min_angMom);
669 gezelter 1722 }
670 gezelter 1350 break;
671     case rnemdPx :
672     max_vel.x() = min_vel.x();
673     max_sd->setVel(max_vel);
674     break;
675     case rnemdPy :
676     max_vel.y() = min_vel.y();
677     max_sd->setVel(max_vel);
678     break;
679     case rnemdPz :
680     max_vel.z() = min_vel.z();
681     max_sd->setVel(max_vel);
682     break;
683     default :
684     break;
685 skuang 1341 }
686 gezelter 1350 } else if (min_vals.rank == worldRank) {
687     // I had the minimum but not the maximum:
688    
689     Vector3d max_vel;
690     Vector3d min_vel = min_sd->getVel();
691     MPI::Status status;
692    
693     // point-to-point swap of the velocity vector
694     MPI::COMM_WORLD.Sendrecv(min_vel.getArrayPointer(), 3, MPI::REALTYPE,
695     max_vals.rank, 0,
696     max_vel.getArrayPointer(), 3, MPI::REALTYPE,
697     max_vals.rank, 0, status);
698    
699 gezelter 1773 switch(rnemdFluxType_) {
700     case rnemdKE :
701 gezelter 1350 min_sd->setVel(max_vel);
702 gezelter 1722 //angular momenta exchange enabled
703 gezelter 1350 if (min_sd->isDirectional()) {
704     Vector3d min_angMom = min_sd->getJ();
705     Vector3d max_angMom;
706 gezelter 1629
707 gezelter 1350 // point-to-point swap of the angular momentum vector
708     MPI::COMM_WORLD.Sendrecv(min_angMom.getArrayPointer(), 3,
709     MPI::REALTYPE, max_vals.rank, 1,
710     max_angMom.getArrayPointer(), 3,
711     MPI::REALTYPE, max_vals.rank, 1,
712     status);
713 gezelter 1629
714 gezelter 1350 min_sd->setJ(max_angMom);
715     }
716     break;
717     case rnemdPx :
718     min_vel.x() = max_vel.x();
719     min_sd->setVel(min_vel);
720     break;
721     case rnemdPy :
722     min_vel.y() = max_vel.y();
723     min_sd->setVel(min_vel);
724     break;
725     case rnemdPz :
726     min_vel.z() = max_vel.z();
727     min_sd->setVel(min_vel);
728     break;
729     default :
730     break;
731     }
732     }
733     #endif
734 gezelter 1773
735     switch(rnemdFluxType_) {
736     case rnemdKE:
737     kineticExchange_ += max_val - min_val;
738     break;
739     case rnemdPx:
740     momentumExchange_.x() += max_val - min_val;
741     break;
742     case rnemdPy:
743     momentumExchange_.y() += max_val - min_val;
744     break;
745     case rnemdPz:
746     momentumExchange_.z() += max_val - min_val;
747     break;
748     default:
749     break;
750     }
751 gezelter 1629 } else {
752     sprintf(painCave.errMsg,
753 gezelter 1773 "RNEMD::doSwap exchange NOT performed because min_val > max_val\n");
754 gezelter 1629 painCave.isFatal = 0;
755     painCave.severity = OPENMD_INFO;
756     simError();
757 skuang 1368 failTrialCount_++;
758 skuang 1338 }
759     } else {
760 gezelter 1629 sprintf(painCave.errMsg,
761 gezelter 1773 "RNEMD::doSwap exchange NOT performed because selected object\n"
762     "\twas not present in at least one of the two slabs.\n");
763 gezelter 1629 painCave.isFatal = 0;
764     painCave.severity = OPENMD_INFO;
765     simError();
766 skuang 1368 failTrialCount_++;
767 gezelter 1773 }
768 skuang 1338 }
769 gezelter 1350
770 gezelter 1773 void RNEMD::doNIVS() {
771 gezelter 1776 if (!doRNEMD_) return;
772 skuang 1338 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
773     Mat3x3d hmat = currentSnap_->getHmat();
774    
775     seleMan_.setSelectionSet(evaluator_.evaluate());
776    
777     int selei;
778     StuntDouble* sd;
779    
780 gezelter 1629 vector<StuntDouble*> hotBin, coldBin;
781 gezelter 1350
782 skuang 1368 RealType Phx = 0.0;
783     RealType Phy = 0.0;
784     RealType Phz = 0.0;
785     RealType Khx = 0.0;
786     RealType Khy = 0.0;
787     RealType Khz = 0.0;
788 gezelter 1722 RealType Khw = 0.0;
789 skuang 1368 RealType Pcx = 0.0;
790     RealType Pcy = 0.0;
791     RealType Pcz = 0.0;
792     RealType Kcx = 0.0;
793     RealType Kcy = 0.0;
794     RealType Kcz = 0.0;
795 gezelter 1722 RealType Kcw = 0.0;
796 skuang 1368
797 skuang 1338 for (sd = seleMan_.beginSelected(selei); sd != NULL;
798     sd = seleMan_.nextSelected(selei)) {
799 skuang 1368
800     Vector3d pos = sd->getPos();
801    
802     // wrap the stuntdouble's position back into the box:
803    
804     if (usePeriodicBoundaryConditions_)
805     currentSnap_->wrapVector(pos);
806    
807     // which bin is this stuntdouble in?
808 gezelter 1773 bool inA = inSlabA(pos);
809     bool inB = inSlabB(pos);
810 skuang 1368
811 gezelter 1773 if (inA || inB) {
812    
813 skuang 1368 RealType mass = sd->getMass();
814     Vector3d vel = sd->getVel();
815    
816 gezelter 1773 if (inA) {
817 skuang 1368 hotBin.push_back(sd);
818     Phx += mass * vel.x();
819     Phy += mass * vel.y();
820     Phz += mass * vel.z();
821     Khx += mass * vel.x() * vel.x();
822     Khy += mass * vel.y() * vel.y();
823     Khz += mass * vel.z() * vel.z();
824 gezelter 1722 if (sd->isDirectional()) {
825     Vector3d angMom = sd->getJ();
826     Mat3x3d I = sd->getI();
827     if (sd->isLinear()) {
828     int i = sd->linearAxis();
829     int j = (i + 1) % 3;
830     int k = (i + 2) % 3;
831     Khw += angMom[j] * angMom[j] / I(j, j) +
832     angMom[k] * angMom[k] / I(k, k);
833     } else {
834     Khw += angMom[0]*angMom[0]/I(0, 0)
835     + angMom[1]*angMom[1]/I(1, 1)
836     + angMom[2]*angMom[2]/I(2, 2);
837     }
838     }
839 gezelter 1773 } else {
840 skuang 1368 coldBin.push_back(sd);
841     Pcx += mass * vel.x();
842     Pcy += mass * vel.y();
843     Pcz += mass * vel.z();
844     Kcx += mass * vel.x() * vel.x();
845     Kcy += mass * vel.y() * vel.y();
846     Kcz += mass * vel.z() * vel.z();
847 gezelter 1722 if (sd->isDirectional()) {
848     Vector3d angMom = sd->getJ();
849     Mat3x3d I = sd->getI();
850     if (sd->isLinear()) {
851     int i = sd->linearAxis();
852     int j = (i + 1) % 3;
853     int k = (i + 2) % 3;
854     Kcw += angMom[j] * angMom[j] / I(j, j) +
855     angMom[k] * angMom[k] / I(k, k);
856     } else {
857     Kcw += angMom[0]*angMom[0]/I(0, 0)
858     + angMom[1]*angMom[1]/I(1, 1)
859     + angMom[2]*angMom[2]/I(2, 2);
860     }
861     }
862 skuang 1368 }
863     }
864     }
865 gezelter 1722
866 skuang 1368 Khx *= 0.5;
867     Khy *= 0.5;
868     Khz *= 0.5;
869 gezelter 1722 Khw *= 0.5;
870 skuang 1368 Kcx *= 0.5;
871     Kcy *= 0.5;
872     Kcz *= 0.5;
873 gezelter 1722 Kcw *= 0.5;
874 skuang 1368
875     #ifdef IS_MPI
876     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phx, 1, MPI::REALTYPE, MPI::SUM);
877     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phy, 1, MPI::REALTYPE, MPI::SUM);
878     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Phz, 1, MPI::REALTYPE, MPI::SUM);
879     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcx, 1, MPI::REALTYPE, MPI::SUM);
880     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcy, 1, MPI::REALTYPE, MPI::SUM);
881     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pcz, 1, MPI::REALTYPE, MPI::SUM);
882    
883     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khx, 1, MPI::REALTYPE, MPI::SUM);
884     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khy, 1, MPI::REALTYPE, MPI::SUM);
885     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khz, 1, MPI::REALTYPE, MPI::SUM);
886 gezelter 1722 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Khw, 1, MPI::REALTYPE, MPI::SUM);
887    
888 skuang 1368 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcx, 1, MPI::REALTYPE, MPI::SUM);
889     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcy, 1, MPI::REALTYPE, MPI::SUM);
890     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcz, 1, MPI::REALTYPE, MPI::SUM);
891 gezelter 1722 MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kcw, 1, MPI::REALTYPE, MPI::SUM);
892 skuang 1368 #endif
893    
894 gezelter 1722 //solve coldBin coeff's first
895 skuang 1368 RealType px = Pcx / Phx;
896     RealType py = Pcy / Phy;
897     RealType pz = Pcz / Phz;
898 gezelter 1722 RealType c, x, y, z;
899     bool successfulScale = false;
900 gezelter 1773 if ((rnemdFluxType_ == rnemdFullKE) ||
901     (rnemdFluxType_ == rnemdRotKE)) {
902 gezelter 1722 //may need sanity check Khw & Kcw > 0
903 skuang 1368
904 gezelter 1773 if (rnemdFluxType_ == rnemdFullKE) {
905     c = 1.0 - kineticTarget_ / (Kcx + Kcy + Kcz + Kcw);
906 gezelter 1722 } else {
907 gezelter 1773 c = 1.0 - kineticTarget_ / Kcw;
908 gezelter 1722 }
909 skuang 1368
910 gezelter 1722 if ((c > 0.81) && (c < 1.21)) {//restrict scaling coefficients
911     c = sqrt(c);
912 gezelter 1803
913 gezelter 1722 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     for (sdi = hotBin.begin(); sdi != hotBin.end(); sdi++) {
951 gezelter 1773 if (rnemdFluxType_ == rnemdFullKE) {
952 gezelter 1722 vel = (*sdi)->getVel();
953     vel.x() *= x;
954     vel.y() *= y;
955     vel.z() *= z;
956     (*sdi)->setVel(vel);
957     }
958     if ((*sdi)->isDirectional()) {
959     Vector3d angMom = (*sdi)->getJ() * w;
960     (*sdi)->setJ(angMom);
961     }
962     }
963     successfulScale = true;
964 gezelter 1773 kineticExchange_ += kineticTarget_;
965 gezelter 1722 }
966 skuang 1368 }
967 gezelter 1722 } else {
968     RealType a000, a110, c0, a001, a111, b01, b11, c1;
969 gezelter 1773 switch(rnemdFluxType_) {
970     case rnemdKE :
971 gezelter 1722 /* used hotBin coeff's & only scale x & y dimensions
972     RealType px = Phx / Pcx;
973     RealType py = Phy / Pcy;
974     a110 = Khy;
975 gezelter 1773 c0 = - Khx - Khy - kineticTarget_;
976 gezelter 1722 a000 = Khx;
977     a111 = Kcy * py * py;
978     b11 = -2.0 * Kcy * py * (1.0 + py);
979 gezelter 1773 c1 = Kcy * py * (2.0 + py) + Kcx * px * ( 2.0 + px) + kineticTarget_;
980 gezelter 1722 b01 = -2.0 * Kcx * px * (1.0 + px);
981     a001 = Kcx * px * px;
982     */
983     //scale all three dimensions, let c_x = c_y
984     a000 = Kcx + Kcy;
985     a110 = Kcz;
986 gezelter 1773 c0 = kineticTarget_ - Kcx - Kcy - Kcz;
987 gezelter 1722 a001 = Khx * px * px + Khy * py * py;
988     a111 = Khz * pz * pz;
989     b01 = -2.0 * (Khx * px * (1.0 + px) + Khy * py * (1.0 + py));
990     b11 = -2.0 * Khz * pz * (1.0 + pz);
991     c1 = Khx * px * (2.0 + px) + Khy * py * (2.0 + py)
992 gezelter 1773 + Khz * pz * (2.0 + pz) - kineticTarget_;
993 gezelter 1722 break;
994 gezelter 1773 case rnemdPx :
995     c = 1 - momentumTarget_.x() / Pcx;
996 gezelter 1722 a000 = Kcy;
997     a110 = Kcz;
998     c0 = Kcx * c * c - Kcx - Kcy - Kcz;
999     a001 = py * py * Khy;
1000     a111 = pz * pz * Khz;
1001     b01 = -2.0 * Khy * py * (1.0 + py);
1002     b11 = -2.0 * Khz * pz * (1.0 + pz);
1003     c1 = Khy * py * (2.0 + py) + Khz * pz * (2.0 + pz)
1004     + Khx * (fastpow(c * px - px - 1.0, 2) - 1.0);
1005     break;
1006 gezelter 1773 case rnemdPy :
1007     c = 1 - momentumTarget_.y() / Pcy;
1008 gezelter 1722 a000 = Kcx;
1009     a110 = Kcz;
1010     c0 = Kcy * c * c - Kcx - Kcy - Kcz;
1011     a001 = px * px * Khx;
1012     a111 = pz * pz * Khz;
1013     b01 = -2.0 * Khx * px * (1.0 + px);
1014     b11 = -2.0 * Khz * pz * (1.0 + pz);
1015     c1 = Khx * px * (2.0 + px) + Khz * pz * (2.0 + pz)
1016     + Khy * (fastpow(c * py - py - 1.0, 2) - 1.0);
1017     break;
1018 gezelter 1773 case rnemdPz ://we don't really do this, do we?
1019     c = 1 - momentumTarget_.z() / Pcz;
1020 gezelter 1722 a000 = Kcx;
1021     a110 = Kcy;
1022     c0 = Kcz * c * c - Kcx - Kcy - Kcz;
1023     a001 = px * px * Khx;
1024     a111 = py * py * Khy;
1025     b01 = -2.0 * Khx * px * (1.0 + px);
1026     b11 = -2.0 * Khy * py * (1.0 + py);
1027     c1 = Khx * px * (2.0 + px) + Khy * py * (2.0 + py)
1028     + Khz * (fastpow(c * pz - pz - 1.0, 2) - 1.0);
1029     break;
1030     default :
1031     break;
1032 skuang 1368 }
1033 gezelter 1722
1034     RealType v1 = a000 * a111 - a001 * a110;
1035     RealType v2 = a000 * b01;
1036     RealType v3 = a000 * b11;
1037     RealType v4 = a000 * c1 - a001 * c0;
1038     RealType v8 = a110 * b01;
1039     RealType v10 = - b01 * c0;
1040    
1041     RealType u0 = v2 * v10 - v4 * v4;
1042     RealType u1 = -2.0 * v3 * v4;
1043     RealType u2 = -v2 * v8 - v3 * v3 - 2.0 * v1 * v4;
1044     RealType u3 = -2.0 * v1 * v3;
1045     RealType u4 = - v1 * v1;
1046     //rescale coefficients
1047     RealType maxAbs = fabs(u0);
1048     if (maxAbs < fabs(u1)) maxAbs = fabs(u1);
1049     if (maxAbs < fabs(u2)) maxAbs = fabs(u2);
1050     if (maxAbs < fabs(u3)) maxAbs = fabs(u3);
1051     if (maxAbs < fabs(u4)) maxAbs = fabs(u4);
1052     u0 /= maxAbs;
1053     u1 /= maxAbs;
1054     u2 /= maxAbs;
1055     u3 /= maxAbs;
1056     u4 /= maxAbs;
1057     //max_element(start, end) is also available.
1058     Polynomial<RealType> poly; //same as DoublePolynomial poly;
1059     poly.setCoefficient(4, u4);
1060     poly.setCoefficient(3, u3);
1061     poly.setCoefficient(2, u2);
1062     poly.setCoefficient(1, u1);
1063     poly.setCoefficient(0, u0);
1064     vector<RealType> realRoots = poly.FindRealRoots();
1065    
1066     vector<RealType>::iterator ri;
1067     RealType r1, r2, alpha0;
1068     vector<pair<RealType,RealType> > rps;
1069     for (ri = realRoots.begin(); ri !=realRoots.end(); ri++) {
1070     r2 = *ri;
1071     //check if FindRealRoots() give the right answer
1072     if ( fabs(u0 + r2 * (u1 + r2 * (u2 + r2 * (u3 + r2 * u4)))) > 1e-6 ) {
1073     sprintf(painCave.errMsg,
1074     "RNEMD Warning: polynomial solve seems to have an error!");
1075     painCave.isFatal = 0;
1076     simError();
1077     failRootCount_++;
1078     }
1079     //might not be useful w/o rescaling coefficients
1080     alpha0 = -c0 - a110 * r2 * r2;
1081     if (alpha0 >= 0.0) {
1082     r1 = sqrt(alpha0 / a000);
1083     if (fabs(c1 + r1 * (b01 + r1 * a001) + r2 * (b11 + r2 * a111))
1084     < 1e-6)
1085     { rps.push_back(make_pair(r1, r2)); }
1086     if (r1 > 1e-6) { //r1 non-negative
1087     r1 = -r1;
1088     if (fabs(c1 + r1 * (b01 + r1 * a001) + r2 * (b11 + r2 * a111))
1089     < 1e-6)
1090     { rps.push_back(make_pair(r1, r2)); }
1091     }
1092     }
1093 skuang 1368 }
1094 gezelter 1722 // Consider combining together the solving pair part w/ the searching
1095     // best solution part so that we don't need the pairs vector
1096     if (!rps.empty()) {
1097     RealType smallestDiff = HONKING_LARGE_VALUE;
1098     RealType diff;
1099     pair<RealType,RealType> bestPair = make_pair(1.0, 1.0);
1100     vector<pair<RealType,RealType> >::iterator rpi;
1101     for (rpi = rps.begin(); rpi != rps.end(); rpi++) {
1102     r1 = (*rpi).first;
1103     r2 = (*rpi).second;
1104 gezelter 1773 switch(rnemdFluxType_) {
1105     case rnemdKE :
1106 gezelter 1722 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1107     + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcx, 2)
1108     + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcy, 2);
1109     break;
1110 gezelter 1773 case rnemdPx :
1111 gezelter 1722 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1112     + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcy, 2);
1113     break;
1114 gezelter 1773 case rnemdPy :
1115 gezelter 1722 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1116     + fastpow(r1 * r1 / r2 / r2 - Kcz/Kcx, 2);
1117     break;
1118 gezelter 1773 case rnemdPz :
1119 gezelter 1722 diff = fastpow(1.0 - r1, 2) + fastpow(1.0 - r2, 2)
1120     + fastpow(r1 * r1 / r2 / r2 - Kcy/Kcx, 2);
1121     default :
1122     break;
1123     }
1124     if (diff < smallestDiff) {
1125     smallestDiff = diff;
1126     bestPair = *rpi;
1127     }
1128     }
1129 skuang 1368 #ifdef IS_MPI
1130 gezelter 1722 if (worldRank == 0) {
1131 skuang 1368 #endif
1132 gezelter 1773 // sprintf(painCave.errMsg,
1133     // "RNEMD: roots r1= %lf\tr2 = %lf\n",
1134     // bestPair.first, bestPair.second);
1135     // painCave.isFatal = 0;
1136     // painCave.severity = OPENMD_INFO;
1137     // simError();
1138 skuang 1368 #ifdef IS_MPI
1139 gezelter 1722 }
1140 skuang 1368 #endif
1141 gezelter 1722
1142 gezelter 1773 switch(rnemdFluxType_) {
1143     case rnemdKE :
1144 gezelter 1722 x = bestPair.first;
1145     y = bestPair.first;
1146     z = bestPair.second;
1147     break;
1148 gezelter 1773 case rnemdPx :
1149 gezelter 1722 x = c;
1150     y = bestPair.first;
1151     z = bestPair.second;
1152     break;
1153 gezelter 1773 case rnemdPy :
1154 gezelter 1722 x = bestPair.first;
1155     y = c;
1156     z = bestPair.second;
1157     break;
1158 gezelter 1773 case rnemdPz :
1159 gezelter 1722 x = bestPair.first;
1160     y = bestPair.second;
1161     z = c;
1162     break;
1163     default :
1164     break;
1165     }
1166     vector<StuntDouble*>::iterator sdi;
1167     Vector3d vel;
1168     for (sdi = coldBin.begin(); sdi != coldBin.end(); sdi++) {
1169     vel = (*sdi)->getVel();
1170     vel.x() *= x;
1171     vel.y() *= y;
1172     vel.z() *= z;
1173     (*sdi)->setVel(vel);
1174     }
1175     //convert to hotBin coefficient
1176     x = 1.0 + px * (1.0 - x);
1177     y = 1.0 + py * (1.0 - y);
1178     z = 1.0 + pz * (1.0 - z);
1179     for (sdi = hotBin.begin(); sdi != hotBin.end(); sdi++) {
1180     vel = (*sdi)->getVel();
1181     vel.x() *= x;
1182     vel.y() *= y;
1183     vel.z() *= z;
1184     (*sdi)->setVel(vel);
1185     }
1186     successfulScale = true;
1187 gezelter 1773 switch(rnemdFluxType_) {
1188     case rnemdKE :
1189     kineticExchange_ += kineticTarget_;
1190     break;
1191     case rnemdPx :
1192     case rnemdPy :
1193     case rnemdPz :
1194     momentumExchange_ += momentumTarget_;
1195     break;
1196     default :
1197     break;
1198     }
1199 gezelter 1629 }
1200 gezelter 1722 }
1201     if (successfulScale != true) {
1202     sprintf(painCave.errMsg,
1203 gezelter 1773 "RNEMD::doNIVS exchange NOT performed - roots that solve\n"
1204     "\tthe constraint equations may not exist or there may be\n"
1205     "\tno selected objects in one or both slabs.\n");
1206 gezelter 1722 painCave.isFatal = 0;
1207     painCave.severity = OPENMD_INFO;
1208     simError();
1209     failTrialCount_++;
1210     }
1211     }
1212    
1213 gezelter 1773 void RNEMD::doVSS() {
1214 gezelter 1776 if (!doRNEMD_) return;
1215 gezelter 1722 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1216 jmarr 1728 RealType time = currentSnap_->getTime();
1217 gezelter 1722 Mat3x3d hmat = currentSnap_->getHmat();
1218    
1219     seleMan_.setSelectionSet(evaluator_.evaluate());
1220    
1221     int selei;
1222     StuntDouble* sd;
1223    
1224     vector<StuntDouble*> hotBin, coldBin;
1225    
1226     Vector3d Ph(V3Zero);
1227     RealType Mh = 0.0;
1228     RealType Kh = 0.0;
1229     Vector3d Pc(V3Zero);
1230     RealType Mc = 0.0;
1231     RealType Kc = 0.0;
1232 jmarr 1728
1233 gezelter 1722
1234     for (sd = seleMan_.beginSelected(selei); sd != NULL;
1235     sd = seleMan_.nextSelected(selei)) {
1236    
1237     Vector3d pos = sd->getPos();
1238    
1239     // wrap the stuntdouble's position back into the box:
1240    
1241     if (usePeriodicBoundaryConditions_)
1242     currentSnap_->wrapVector(pos);
1243    
1244     // which bin is this stuntdouble in?
1245 gezelter 1773 bool inA = inSlabA(pos);
1246     bool inB = inSlabB(pos);
1247    
1248     if (inA || inB) {
1249 gezelter 1722
1250     RealType mass = sd->getMass();
1251     Vector3d vel = sd->getVel();
1252    
1253 gezelter 1773 if (inA) {
1254 gezelter 1722 hotBin.push_back(sd);
1255     Ph += mass * vel;
1256     Mh += mass;
1257     Kh += mass * vel.lengthSquare();
1258 gezelter 1773 if (rnemdFluxType_ == rnemdFullKE) {
1259 gezelter 1722 if (sd->isDirectional()) {
1260     Vector3d angMom = sd->getJ();
1261     Mat3x3d I = sd->getI();
1262     if (sd->isLinear()) {
1263     int i = sd->linearAxis();
1264     int j = (i + 1) % 3;
1265     int k = (i + 2) % 3;
1266     Kh += angMom[j] * angMom[j] / I(j, j) +
1267     angMom[k] * angMom[k] / I(k, k);
1268     } else {
1269     Kh += angMom[0] * angMom[0] / I(0, 0) +
1270     angMom[1] * angMom[1] / I(1, 1) +
1271     angMom[2] * angMom[2] / I(2, 2);
1272     }
1273     }
1274     }
1275     } else { //midBin_
1276     coldBin.push_back(sd);
1277     Pc += mass * vel;
1278     Mc += mass;
1279     Kc += mass * vel.lengthSquare();
1280 gezelter 1773 if (rnemdFluxType_ == rnemdFullKE) {
1281 gezelter 1722 if (sd->isDirectional()) {
1282     Vector3d angMom = sd->getJ();
1283     Mat3x3d I = sd->getI();
1284     if (sd->isLinear()) {
1285     int i = sd->linearAxis();
1286     int j = (i + 1) % 3;
1287     int k = (i + 2) % 3;
1288     Kc += angMom[j] * angMom[j] / I(j, j) +
1289     angMom[k] * angMom[k] / I(k, k);
1290     } else {
1291     Kc += angMom[0] * angMom[0] / I(0, 0) +
1292     angMom[1] * angMom[1] / I(1, 1) +
1293     angMom[2] * angMom[2] / I(2, 2);
1294     }
1295     }
1296     }
1297     }
1298 skuang 1368 }
1299 gezelter 1722 }
1300    
1301     Kh *= 0.5;
1302     Kc *= 0.5;
1303 jmarr 1728
1304 gezelter 1722 #ifdef IS_MPI
1305     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Ph[0], 3, MPI::REALTYPE, MPI::SUM);
1306     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Pc[0], 3, MPI::REALTYPE, MPI::SUM);
1307     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Mh, 1, MPI::REALTYPE, MPI::SUM);
1308     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kh, 1, MPI::REALTYPE, MPI::SUM);
1309     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Mc, 1, MPI::REALTYPE, MPI::SUM);
1310     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &Kc, 1, MPI::REALTYPE, MPI::SUM);
1311     #endif
1312    
1313     bool successfulExchange = false;
1314     if ((Mh > 0.0) && (Mc > 0.0)) {//both slabs are not empty
1315     Vector3d vc = Pc / Mc;
1316 gezelter 1773 Vector3d ac = -momentumTarget_ / Mc + vc;
1317     Vector3d acrec = -momentumTarget_ / Mc;
1318     RealType cNumerator = Kc - kineticTarget_ - 0.5 * Mc * ac.lengthSquare();
1319 gezelter 1722 if (cNumerator > 0.0) {
1320     RealType cDenominator = Kc - 0.5 * Mc * vc.lengthSquare();
1321     if (cDenominator > 0.0) {
1322     RealType c = sqrt(cNumerator / cDenominator);
1323     if ((c > 0.9) && (c < 1.1)) {//restrict scaling coefficients
1324     Vector3d vh = Ph / Mh;
1325 gezelter 1773 Vector3d ah = momentumTarget_ / Mh + vh;
1326     Vector3d ahrec = momentumTarget_ / Mh;
1327     RealType hNumerator = Kh + kineticTarget_
1328 gezelter 1722 - 0.5 * Mh * ah.lengthSquare();
1329     if (hNumerator > 0.0) {
1330     RealType hDenominator = Kh - 0.5 * Mh * vh.lengthSquare();
1331     if (hDenominator > 0.0) {
1332     RealType h = sqrt(hNumerator / hDenominator);
1333     if ((h > 0.9) && (h < 1.1)) {
1334 gezelter 1803
1335 gezelter 1722 vector<StuntDouble*>::iterator sdi;
1336     Vector3d vel;
1337     for (sdi = coldBin.begin(); sdi != coldBin.end(); sdi++) {
1338     //vel = (*sdi)->getVel();
1339     vel = ((*sdi)->getVel() - vc) * c + ac;
1340     (*sdi)->setVel(vel);
1341 gezelter 1773 if (rnemdFluxType_ == rnemdFullKE) {
1342 gezelter 1722 if ((*sdi)->isDirectional()) {
1343     Vector3d angMom = (*sdi)->getJ() * c;
1344     (*sdi)->setJ(angMom);
1345     }
1346     }
1347     }
1348     for (sdi = hotBin.begin(); sdi != hotBin.end(); sdi++) {
1349     //vel = (*sdi)->getVel();
1350     vel = ((*sdi)->getVel() - vh) * h + ah;
1351     (*sdi)->setVel(vel);
1352 gezelter 1773 if (rnemdFluxType_ == rnemdFullKE) {
1353 gezelter 1722 if ((*sdi)->isDirectional()) {
1354     Vector3d angMom = (*sdi)->getJ() * h;
1355     (*sdi)->setJ(angMom);
1356     }
1357     }
1358     }
1359     successfulExchange = true;
1360 gezelter 1773 kineticExchange_ += kineticTarget_;
1361     momentumExchange_ += momentumTarget_;
1362 gezelter 1722 }
1363     }
1364     }
1365     }
1366     }
1367 skuang 1368 }
1368 gezelter 1722 }
1369     if (successfulExchange != true) {
1370 gezelter 1773 sprintf(painCave.errMsg,
1371     "RNEMD::doVSS exchange NOT performed - roots that solve\n"
1372     "\tthe constraint equations may not exist or there may be\n"
1373     "\tno selected objects in one or both slabs.\n");
1374     painCave.isFatal = 0;
1375     painCave.severity = OPENMD_INFO;
1376     simError();
1377 skuang 1368 failTrialCount_++;
1378     }
1379     }
1380    
1381     void RNEMD::doRNEMD() {
1382 gezelter 1776 if (!doRNEMD_) return;
1383 gezelter 1773 trialCount_++;
1384     switch(rnemdMethod_) {
1385     case rnemdSwap:
1386 skuang 1368 doSwap();
1387     break;
1388 gezelter 1773 case rnemdNIVS:
1389     doNIVS();
1390 gezelter 1722 break;
1391 gezelter 1773 case rnemdVSS:
1392     doVSS();
1393     break;
1394     case rnemdUnkownMethod:
1395 skuang 1368 default :
1396     break;
1397     }
1398     }
1399    
1400     void RNEMD::collectData() {
1401 gezelter 1776 if (!doRNEMD_) return;
1402 skuang 1368 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1403     Mat3x3d hmat = currentSnap_->getHmat();
1404    
1405 gezelter 1774 areaAccumulator_->add(currentSnap_->getXYarea());
1406    
1407 skuang 1368 seleMan_.setSelectionSet(evaluator_.evaluate());
1408    
1409 gezelter 1803 int selei(0);
1410 skuang 1368 StuntDouble* sd;
1411    
1412 gezelter 1773 vector<RealType> binMass(nBins_, 0.0);
1413     vector<RealType> binPx(nBins_, 0.0);
1414     vector<RealType> binPy(nBins_, 0.0);
1415     vector<RealType> binPz(nBins_, 0.0);
1416     vector<RealType> binKE(nBins_, 0.0);
1417     vector<int> binDOF(nBins_, 0);
1418     vector<int> binCount(nBins_, 0);
1419 jmarr 1728
1420 gezelter 1629 // alternative approach, track all molecules instead of only those
1421     // selected for scaling/swapping:
1422     /*
1423     SimInfo::MoleculeIterator miter;
1424     vector<StuntDouble*>::iterator iiter;
1425     Molecule* mol;
1426 gezelter 1769 StuntDouble* sd;
1427 gezelter 1629 for (mol = info_->beginMolecule(miter); mol != NULL;
1428 jmarr 1728 mol = info_->nextMolecule(miter))
1429 gezelter 1769 sd is essentially sd
1430     for (sd = mol->beginIntegrableObject(iiter);
1431     sd != NULL;
1432     sd = mol->nextIntegrableObject(iiter))
1433 gezelter 1629 */
1434 gezelter 1803
1435 skuang 1368 for (sd = seleMan_.beginSelected(selei); sd != NULL;
1436 gezelter 1813 sd = seleMan_.nextSelected(selei)) {
1437 gezelter 1803
1438 skuang 1338 Vector3d pos = sd->getPos();
1439    
1440     // wrap the stuntdouble's position back into the box:
1441    
1442     if (usePeriodicBoundaryConditions_)
1443     currentSnap_->wrapVector(pos);
1444 gezelter 1773
1445 gezelter 1774
1446 skuang 1338 // which bin is this stuntdouble in?
1447     // wrapped positions are in the range [-0.5*hmat(2,2), +0.5*hmat(2,2)]
1448 gezelter 1773 // Shift molecules by half a box to have bins start at 0
1449     // The modulo operator is used to wrap the case when we are
1450     // beyond the end of the bins back to the beginning.
1451     int binNo = int(nBins_ * (pos.z() / hmat(2,2) + 0.5)) % nBins_;
1452 gezelter 1803
1453 skuang 1338 RealType mass = sd->getMass();
1454     Vector3d vel = sd->getVel();
1455    
1456 gezelter 1773 binCount[binNo]++;
1457     binMass[binNo] += mass;
1458     binPx[binNo] += mass*vel.x();
1459     binPy[binNo] += mass*vel.y();
1460     binPz[binNo] += mass*vel.z();
1461     binKE[binNo] += 0.5 * (mass * vel.lengthSquare());
1462     binDOF[binNo] += 3;
1463    
1464     if (sd->isDirectional()) {
1465     Vector3d angMom = sd->getJ();
1466     Mat3x3d I = sd->getI();
1467     if (sd->isLinear()) {
1468     int i = sd->linearAxis();
1469     int j = (i + 1) % 3;
1470     int k = (i + 2) % 3;
1471     binKE[binNo] += 0.5 * (angMom[j] * angMom[j] / I(j, j) +
1472     angMom[k] * angMom[k] / I(k, k));
1473     binDOF[binNo] += 2;
1474     } else {
1475     binKE[binNo] += 0.5 * (angMom[0] * angMom[0] / I(0, 0) +
1476     angMom[1] * angMom[1] / I(1, 1) +
1477     angMom[2] * angMom[2] / I(2, 2));
1478     binDOF[binNo] += 3;
1479     }
1480 gezelter 1722 }
1481 gezelter 1773 }
1482    
1483     #ifdef IS_MPI
1484     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binCount[0],
1485     nBins_, MPI::INT, MPI::SUM);
1486     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binMass[0],
1487     nBins_, MPI::REALTYPE, MPI::SUM);
1488     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPx[0],
1489     nBins_, MPI::REALTYPE, MPI::SUM);
1490     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPy[0],
1491     nBins_, MPI::REALTYPE, MPI::SUM);
1492     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binPz[0],
1493     nBins_, MPI::REALTYPE, MPI::SUM);
1494     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binKE[0],
1495     nBins_, MPI::REALTYPE, MPI::SUM);
1496     MPI::COMM_WORLD.Allreduce(MPI::IN_PLACE, &binDOF[0],
1497     nBins_, MPI::INT, MPI::SUM);
1498     #endif
1499    
1500     Vector3d vel;
1501     RealType den;
1502     RealType temp;
1503     RealType z;
1504     for (int i = 0; i < nBins_; i++) {
1505     z = (((RealType)i + 0.5) / (RealType)nBins_) * hmat(2,2);
1506     vel.x() = binPx[i] / binMass[i];
1507     vel.y() = binPy[i] / binMass[i];
1508     vel.z() = binPz[i] / binMass[i];
1509 gezelter 1777
1510     den = binMass[i] * nBins_ * PhysicalConstants::densityConvert
1511     / currentSnap_->getVolume() ;
1512    
1513 gezelter 1812 if (binCount[i] > 0) {
1514     // only add values if there are things to add
1515     temp = 2.0 * binKE[i] / (binDOF[i] * PhysicalConstants::kb *
1516     PhysicalConstants::energyConvert);
1517    
1518     for (unsigned int j = 0; j < outputMask_.size(); ++j) {
1519     if(outputMask_[j]) {
1520     switch(j) {
1521     case Z:
1522     dynamic_cast<Accumulator *>(data_[j].accumulator[i])->add(z);
1523     break;
1524     case TEMPERATURE:
1525     dynamic_cast<Accumulator *>(data_[j].accumulator[i])->add(temp);
1526     break;
1527     case VELOCITY:
1528     dynamic_cast<VectorAccumulator *>(data_[j].accumulator[i])->add(vel);
1529     break;
1530     case DENSITY:
1531     dynamic_cast<Accumulator *>(data_[j].accumulator[i])->add(den);
1532     break;
1533     }
1534 gezelter 1773 }
1535     }
1536 gezelter 1629 }
1537 skuang 1338 }
1538 skuang 1368 }
1539    
1540     void RNEMD::getStarted() {
1541 gezelter 1776 if (!doRNEMD_) return;
1542 gezelter 1629 collectData();
1543 gezelter 1773 writeOutputFile();
1544 skuang 1368 }
1545    
1546 gezelter 1773 void RNEMD::parseOutputFileFormat(const std::string& format) {
1547 gezelter 1776 if (!doRNEMD_) return;
1548 gezelter 1773 StringTokenizer tokenizer(format, " ,;|\t\n\r");
1549    
1550     while(tokenizer.hasMoreTokens()) {
1551     std::string token(tokenizer.nextToken());
1552     toUpper(token);
1553     OutputMapType::iterator i = outputMap_.find(token);
1554     if (i != outputMap_.end()) {
1555     outputMask_.set(i->second);
1556     } else {
1557     sprintf( painCave.errMsg,
1558     "RNEMD::parseOutputFileFormat: %s is not a recognized\n"
1559     "\toutputFileFormat keyword.\n", token.c_str() );
1560     painCave.isFatal = 0;
1561     painCave.severity = OPENMD_ERROR;
1562     simError();
1563     }
1564     }
1565     }
1566    
1567     void RNEMD::writeOutputFile() {
1568 gezelter 1776 if (!doRNEMD_) return;
1569 gezelter 1773
1570 gezelter 1350 #ifdef IS_MPI
1571     // If we're the root node, should we print out the results
1572     int worldRank = MPI::COMM_WORLD.Get_rank();
1573     if (worldRank == 0) {
1574     #endif
1575 gezelter 1773 rnemdFile_.open(rnemdFileName_.c_str(), std::ios::out | std::ios::trunc );
1576    
1577     if( !rnemdFile_ ){
1578     sprintf( painCave.errMsg,
1579     "Could not open \"%s\" for RNEMD output.\n",
1580     rnemdFileName_.c_str());
1581     painCave.isFatal = 1;
1582     simError();
1583     }
1584 gezelter 1722
1585 gezelter 1773 Snapshot* currentSnap_ = info_->getSnapshotManager()->getCurrentSnapshot();
1586    
1587     RealType time = currentSnap_->getTime();
1588 gezelter 1774 RealType avgArea;
1589     areaAccumulator_->getAverage(avgArea);
1590 gezelter 1777 RealType Jz = kineticExchange_ / (2.0 * time * avgArea)
1591     / PhysicalConstants::energyConvert;
1592 gezelter 1774 Vector3d JzP = momentumExchange_ / (2.0 * time * avgArea);
1593    
1594 gezelter 1773 rnemdFile_ << "#######################################################\n";
1595     rnemdFile_ << "# RNEMD {\n";
1596    
1597     map<string, RNEMDMethod>::iterator mi;
1598     for(mi = stringToMethod_.begin(); mi != stringToMethod_.end(); ++mi) {
1599     if ( (*mi).second == rnemdMethod_)
1600 gezelter 1774 rnemdFile_ << "# exchangeMethod = \"" << (*mi).first << "\";\n";
1601 skuang 1368 }
1602 gezelter 1773 map<string, RNEMDFluxType>::iterator fi;
1603     for(fi = stringToFluxType_.begin(); fi != stringToFluxType_.end(); ++fi) {
1604     if ( (*fi).second == rnemdFluxType_)
1605 gezelter 1774 rnemdFile_ << "# fluxType = \"" << (*fi).first << "\";\n";
1606 gezelter 1722 }
1607 gezelter 1773
1608 gezelter 1775 rnemdFile_ << "# exchangeTime = " << exchangeTime_ << ";\n";
1609 gezelter 1773
1610     rnemdFile_ << "# objectSelection = \""
1611 gezelter 1774 << rnemdObjectSelection_ << "\";\n";
1612     rnemdFile_ << "# slabWidth = " << slabWidth_ << ";\n";
1613     rnemdFile_ << "# slabAcenter = " << slabACenter_ << ";\n";
1614     rnemdFile_ << "# slabBcenter = " << slabBCenter_ << ";\n";
1615 gezelter 1773 rnemdFile_ << "# }\n";
1616     rnemdFile_ << "#######################################################\n";
1617 gezelter 1774 rnemdFile_ << "# RNEMD report:\n";
1618     rnemdFile_ << "# running time = " << time << " fs\n";
1619     rnemdFile_ << "# target flux:\n";
1620 gezelter 1777 rnemdFile_ << "# kinetic = "
1621     << kineticFlux_ / PhysicalConstants::energyConvert
1622     << " (kcal/mol/A^2/fs)\n";
1623     rnemdFile_ << "# momentum = " << momentumFluxVector_
1624     << " (amu/A/fs^2)\n";
1625 gezelter 1774 rnemdFile_ << "# target one-time exchanges:\n";
1626 gezelter 1777 rnemdFile_ << "# kinetic = "
1627     << kineticTarget_ / PhysicalConstants::energyConvert
1628     << " (kcal/mol)\n";
1629     rnemdFile_ << "# momentum = " << momentumTarget_
1630     << " (amu*A/fs)\n";
1631 gezelter 1774 rnemdFile_ << "# actual exchange totals:\n";
1632 gezelter 1777 rnemdFile_ << "# kinetic = "
1633     << kineticExchange_ / PhysicalConstants::energyConvert
1634     << " (kcal/mol)\n";
1635     rnemdFile_ << "# momentum = " << momentumExchange_
1636     << " (amu*A/fs)\n";
1637 gezelter 1774 rnemdFile_ << "# actual flux:\n";
1638 gezelter 1777 rnemdFile_ << "# kinetic = " << Jz
1639     << " (kcal/mol/A^2/fs)\n";
1640     rnemdFile_ << "# momentum = " << JzP
1641     << " (amu/A/fs^2)\n";
1642 gezelter 1774 rnemdFile_ << "# exchange statistics:\n";
1643     rnemdFile_ << "# attempted = " << trialCount_ << "\n";
1644     rnemdFile_ << "# failed = " << failTrialCount_ << "\n";
1645 gezelter 1773 if (rnemdMethod_ == rnemdNIVS) {
1646 gezelter 1774 rnemdFile_ << "# NIVS root-check errors = "
1647     << failRootCount_ << "\n";
1648 gezelter 1722 }
1649 gezelter 1773 rnemdFile_ << "#######################################################\n";
1650    
1651    
1652    
1653     //write title
1654     rnemdFile_ << "#";
1655     for (unsigned int i = 0; i < outputMask_.size(); ++i) {
1656     if (outputMask_[i]) {
1657     rnemdFile_ << "\t" << data_[i].title <<
1658     "(" << data_[i].units << ")";
1659 gezelter 1777 // add some extra tabs for column alignment
1660     if (data_[i].dataType == "Vector3d") rnemdFile_ << "\t\t";
1661 skuang 1368 }
1662 gezelter 1773 }
1663     rnemdFile_ << std::endl;
1664    
1665     rnemdFile_.precision(8);
1666    
1667 gezelter 1813 for (int j = 0; j < nBins_; j++) {
1668 gezelter 1773
1669     for (unsigned int i = 0; i < outputMask_.size(); ++i) {
1670     if (outputMask_[i]) {
1671     if (data_[i].dataType == "RealType")
1672     writeReal(i,j);
1673     else if (data_[i].dataType == "Vector3d")
1674     writeVector(i,j);
1675     else {
1676     sprintf( painCave.errMsg,
1677     "RNEMD found an unknown data type for: %s ",
1678     data_[i].title.c_str());
1679     painCave.isFatal = 1;
1680     simError();
1681     }
1682     }
1683 skuang 1368 }
1684 gezelter 1773 rnemdFile_ << std::endl;
1685    
1686     }
1687 gezelter 1774
1688     rnemdFile_ << "#######################################################\n";
1689     rnemdFile_ << "# Standard Deviations in those quantities follow:\n";
1690     rnemdFile_ << "#######################################################\n";
1691    
1692    
1693 gezelter 1813 for (int j = 0; j < nBins_; j++) {
1694 gezelter 1774 rnemdFile_ << "#";
1695     for (unsigned int i = 0; i < outputMask_.size(); ++i) {
1696     if (outputMask_[i]) {
1697     if (data_[i].dataType == "RealType")
1698     writeRealStdDev(i,j);
1699     else if (data_[i].dataType == "Vector3d")
1700     writeVectorStdDev(i,j);
1701     else {
1702     sprintf( painCave.errMsg,
1703     "RNEMD found an unknown data type for: %s ",
1704     data_[i].title.c_str());
1705     painCave.isFatal = 1;
1706     simError();
1707     }
1708     }
1709     }
1710     rnemdFile_ << std::endl;
1711    
1712     }
1713 gezelter 1773
1714     rnemdFile_.flush();
1715     rnemdFile_.close();
1716    
1717 gezelter 1350 #ifdef IS_MPI
1718 gezelter 1396 }
1719 gezelter 1350 #endif
1720 jmarr 1728
1721 gezelter 1334 }
1722 gezelter 1773
1723     void RNEMD::writeReal(int index, unsigned int bin) {
1724 gezelter 1776 if (!doRNEMD_) return;
1725 gezelter 1773 assert(index >=0 && index < ENDINDEX);
1726 gezelter 1830 assert(int(bin) < nBins_);
1727 gezelter 1773 RealType s;
1728    
1729 gezelter 1813 dynamic_cast<Accumulator *>(data_[index].accumulator[bin])->getAverage(s);
1730 gezelter 1773
1731     if (! isinf(s) && ! isnan(s)) {
1732     rnemdFile_ << "\t" << s;
1733     } else{
1734     sprintf( painCave.errMsg,
1735     "RNEMD detected a numerical error writing: %s for bin %d",
1736     data_[index].title.c_str(), bin);
1737     painCave.isFatal = 1;
1738     simError();
1739     }
1740     }
1741    
1742     void RNEMD::writeVector(int index, unsigned int bin) {
1743 gezelter 1776 if (!doRNEMD_) return;
1744 gezelter 1773 assert(index >=0 && index < ENDINDEX);
1745 gezelter 1830 assert(int(bin) < nBins_);
1746 gezelter 1773 Vector3d s;
1747     dynamic_cast<VectorAccumulator*>(data_[index].accumulator[bin])->getAverage(s);
1748     if (isinf(s[0]) || isnan(s[0]) ||
1749     isinf(s[1]) || isnan(s[1]) ||
1750     isinf(s[2]) || isnan(s[2]) ) {
1751     sprintf( painCave.errMsg,
1752     "RNEMD detected a numerical error writing: %s for bin %d",
1753     data_[index].title.c_str(), bin);
1754     painCave.isFatal = 1;
1755     simError();
1756     } else {
1757     rnemdFile_ << "\t" << s[0] << "\t" << s[1] << "\t" << s[2];
1758     }
1759     }
1760 gezelter 1774
1761     void RNEMD::writeRealStdDev(int index, unsigned int bin) {
1762 gezelter 1776 if (!doRNEMD_) return;
1763 gezelter 1774 assert(index >=0 && index < ENDINDEX);
1764 gezelter 1830 assert(int(bin) < nBins_);
1765 gezelter 1774 RealType s;
1766    
1767 gezelter 1813 dynamic_cast<Accumulator *>(data_[index].accumulator[bin])->getStdDev(s);
1768 gezelter 1774
1769     if (! isinf(s) && ! isnan(s)) {
1770     rnemdFile_ << "\t" << s;
1771     } else{
1772     sprintf( painCave.errMsg,
1773     "RNEMD detected a numerical error writing: %s std. dev. for bin %d",
1774     data_[index].title.c_str(), bin);
1775     painCave.isFatal = 1;
1776     simError();
1777     }
1778     }
1779    
1780     void RNEMD::writeVectorStdDev(int index, unsigned int bin) {
1781 gezelter 1776 if (!doRNEMD_) return;
1782 gezelter 1774 assert(index >=0 && index < ENDINDEX);
1783 gezelter 1830 assert(int(bin) < nBins_);
1784 gezelter 1774 Vector3d s;
1785     dynamic_cast<VectorAccumulator*>(data_[index].accumulator[bin])->getStdDev(s);
1786     if (isinf(s[0]) || isnan(s[0]) ||
1787     isinf(s[1]) || isnan(s[1]) ||
1788     isinf(s[2]) || isnan(s[2]) ) {
1789     sprintf( painCave.errMsg,
1790     "RNEMD detected a numerical error writing: %s std. dev. for bin %d",
1791     data_[index].title.c_str(), bin);
1792     painCave.isFatal = 1;
1793     simError();
1794     } else {
1795     rnemdFile_ << "\t" << s[0] << "\t" << s[1] << "\t" << s[2];
1796     }
1797     }
1798 skuang 1338 }
1799 gezelter 1722

Properties

Name Value
svn:keywords Author Id Revision Date