ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/DumpReader.cpp
Revision: 1752
Committed: Sun Jun 10 14:05:02 2012 UTC (12 years, 10 months ago) by gezelter
File size: 21658 byte(s)
Log Message:
Added SiteData parsing and writing.

File Contents

# User Rev Content
1 gezelter 1390 /*
2     * Copyright (c) 2009 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     * 1. Redistributions of source code must retain the above copyright
10     * notice, this list of conditions and the following disclaimer.
11     *
12     * 2. Redistributions in binary form must reproduce the above copyright
13     * 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     *
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 1665 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 1390 */
42 chrisfen 721
43     #define _LARGEFILE_SOURCE64
44     #define _FILE_OFFSET_BITS 64
45    
46     #include <sys/types.h>
47     #include <sys/stat.h>
48    
49     #include <iostream>
50     #include <math.h>
51    
52     #include <stdio.h>
53     #include <stdlib.h>
54     #include <string.h>
55    
56     #include "io/DumpReader.hpp"
57     #include "primitives/Molecule.hpp"
58     #include "utils/simError.h"
59     #include "utils/MemoryUtils.hpp"
60     #include "utils/StringTokenizer.hpp"
61    
62     #ifdef IS_MPI
63    
64     #include <mpi.h>
65     #define TAKE_THIS_TAG_CHAR 0
66     #define TAKE_THIS_TAG_INT 1
67    
68     #endif // is_mpi
69    
70    
71 gezelter 1390 namespace OpenMD {
72 chrisfen 721
73     DumpReader::DumpReader(SimInfo* info, const std::string& filename)
74 gezelter 1104 : info_(info), filename_(filename), isScanned_(false), nframes_(0), needCOMprops_(false) {
75 chrisfen 996
76 chrisfen 721 #ifdef IS_MPI
77 chrisfen 996
78     if (worldRank == 0) {
79 chrisfen 721 #endif
80 chrisfen 996
81 chrisfen 721 inFile_ = new std::ifstream(filename_.c_str());
82 chrisfen 996
83     if (inFile_->fail()) {
84     sprintf(painCave.errMsg,
85     "DumpReader: Cannot open file: %s\n",
86     filename_.c_str());
87     painCave.isFatal = 1;
88     simError();
89     }
90    
91 chrisfen 721 #ifdef IS_MPI
92 chrisfen 996
93     }
94    
95     strcpy(checkPointMsg, "Dump file opened for reading successfully.");
96 gezelter 1241 errorCheckPoint();
97 chrisfen 996
98 chrisfen 721 #endif
99 chrisfen 996
100     return;
101     }
102    
103 chrisfen 721 DumpReader::~DumpReader() {
104 chrisfen 996
105 chrisfen 721 #ifdef IS_MPI
106 chrisfen 996
107 chrisfen 721 if (worldRank == 0) {
108     #endif
109 chrisfen 996
110 chrisfen 721 delete inFile_;
111 chrisfen 996
112 chrisfen 721 #ifdef IS_MPI
113 chrisfen 996
114 chrisfen 721 }
115 chrisfen 996
116 chrisfen 721 strcpy(checkPointMsg, "Dump file closed successfully.");
117 gezelter 1241 errorCheckPoint();
118 chrisfen 996
119 chrisfen 721 #endif
120 chrisfen 996
121 chrisfen 721 return;
122     }
123 chrisfen 996
124 chrisfen 721 int DumpReader::getNFrames(void) {
125    
126     if (!isScanned_)
127     scanFile();
128    
129     return nframes_;
130     }
131    
132     void DumpReader::scanFile(void) {
133 tim 1024 int lineNo = 0;
134     std::streampos prevPos;
135 chrisfen 721 std::streampos currPos;
136 tim 1024
137 chrisfen 721 #ifdef IS_MPI
138 tim 1024
139 chrisfen 721 if (worldRank == 0) {
140     #endif // is_mpi
141 tim 1024
142     currPos = inFile_->tellg();
143     prevPos = currPos;
144     bool foundOpenSnapshotTag = false;
145     bool foundClosedSnapshotTag = false;
146 gezelter 1752 bool foundOpenSiteDataTag = false;
147 tim 1024 while(inFile_->getline(buffer, bufferSize)) {
148     ++lineNo;
149    
150     std::string line = buffer;
151     currPos = inFile_->tellg();
152     if (line.find("<Snapshot>")!= std::string::npos) {
153     if (foundOpenSnapshotTag) {
154 chrisfen 721 sprintf(painCave.errMsg,
155 tim 1024 "DumpReader:<Snapshot> is multiply nested at line %d in %s \n", lineNo,
156     filename_.c_str());
157 chrisfen 721 painCave.isFatal = 1;
158 tim 1024 simError();
159     }
160     foundOpenSnapshotTag = true;
161     foundClosedSnapshotTag = false;
162     framePos_.push_back(prevPos);
163    
164     } else if (line.find("</Snapshot>") != std::string::npos){
165     if (!foundOpenSnapshotTag) {
166     sprintf(painCave.errMsg,
167     "DumpReader:</Snapshot> appears before <Snapshot> at line %d in %s \n", lineNo,
168     filename_.c_str());
169     painCave.isFatal = 1;
170 chrisfen 721 simError();
171 tim 1024 }
172    
173     if (foundClosedSnapshotTag) {
174     sprintf(painCave.errMsg,
175     "DumpReader:</Snapshot> appears multiply nested at line %d in %s \n", lineNo,
176     filename_.c_str());
177     painCave.isFatal = 1;
178     simError();
179     }
180     foundClosedSnapshotTag = true;
181     foundOpenSnapshotTag = false;
182     }
183     prevPos = currPos;
184     }
185    
186     // only found <Snapshot> for the last frame means the file is corrupted, we should discard
187     // it and give a warning message
188     if (foundOpenSnapshotTag) {
189     sprintf(painCave.errMsg,
190     "DumpReader: last frame in %s is invalid\n", filename_.c_str());
191     painCave.isFatal = 0;
192     simError();
193     framePos_.pop_back();
194     }
195    
196 chrisfen 721 nframes_ = framePos_.size();
197 tim 1024
198     if (nframes_ == 0) {
199     sprintf(painCave.errMsg,
200     "DumpReader: %s does not contain a valid frame\n", filename_.c_str());
201     painCave.isFatal = 1;
202     simError();
203     }
204 chrisfen 721 #ifdef IS_MPI
205     }
206    
207     MPI_Bcast(&nframes_, 1, MPI_INT, 0, MPI_COMM_WORLD);
208 tim 1024
209 chrisfen 721 #endif // is_mpi
210 tim 1024
211 chrisfen 721 isScanned_ = true;
212     }
213    
214     void DumpReader::readFrame(int whichFrame) {
215     if (!isScanned_)
216     scanFile();
217    
218     int storageLayout = info_->getSnapshotManager()->getStorageLayout();
219    
220     if (storageLayout & DataStorage::dslPosition) {
221     needPos_ = true;
222     } else {
223     needPos_ = false;
224     }
225    
226     if (storageLayout & DataStorage::dslVelocity) {
227     needVel_ = true;
228     } else {
229     needVel_ = false;
230     }
231    
232     if (storageLayout & DataStorage::dslAmat || storageLayout & DataStorage::dslElectroFrame) {
233     needQuaternion_ = true;
234     } else {
235     needQuaternion_ = false;
236     }
237    
238     if (storageLayout & DataStorage::dslAngularMomentum) {
239     needAngMom_ = true;
240     } else {
241     needAngMom_ = false;
242     }
243    
244     readSet(whichFrame);
245 gezelter 1104
246     if (needCOMprops_) {
247     Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
248     Vector3d com;
249     Vector3d comvel;
250     Vector3d comw;
251 chuckv 1112 if (needPos_ && needVel_){
252     info_->getComAll(com, comvel);
253     comw = info_->getAngularMomentum();
254     }else{
255     com = info_->getCom();
256     comvel = 0.0;
257     comw = 0.0;
258     }
259 gezelter 1104 s->setCOMprops(com, comvel, comw);
260     }
261    
262 chrisfen 721 }
263    
264 tim 1024 void DumpReader::readSet(int whichFrame) {
265     std::string line;
266    
267 chrisfen 721 #ifndef IS_MPI
268     inFile_->clear();
269     inFile_->seekg(framePos_[whichFrame]);
270 tim 1024
271     std::istream& inputStream = *inFile_;
272    
273     #else
274     int masterNode = 0;
275     std::stringstream sstream;
276     if (worldRank == masterNode) {
277     std::string sendBuffer;
278    
279     inFile_->clear();
280     inFile_->seekg(framePos_[whichFrame]);
281    
282     while (inFile_->getline(buffer, bufferSize)) {
283    
284     line = buffer;
285     sendBuffer += line;
286     sendBuffer += '\n';
287     if (line.find("</Snapshot>") != std::string::npos) {
288     break;
289     }
290     }
291    
292     int sendBufferSize = sendBuffer.size();
293     MPI_Bcast(&sendBufferSize, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
294     MPI_Bcast((void *)sendBuffer.c_str(), sendBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
295    
296     sstream.str(sendBuffer);
297     } else {
298     int sendBufferSize;
299     MPI_Bcast(&sendBufferSize, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
300     char * recvBuffer = new char[sendBufferSize+1];
301 gezelter 1313 assert(recvBuffer);
302     recvBuffer[sendBufferSize] = '\0';
303 tim 1024 MPI_Bcast(recvBuffer, sendBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
304     sstream.str(recvBuffer);
305 gezelter 1313 delete [] recvBuffer;
306 tim 1024 }
307    
308     std::istream& inputStream = sstream;
309     #endif
310    
311     inputStream.getline(buffer, bufferSize);
312    
313     line = buffer;
314     if (line.find("<Snapshot>") == std::string::npos) {
315 chrisfen 721 sprintf(painCave.errMsg,
316 tim 1024 "DumpReader Error: can not find <Snapshot>\n");
317 chrisfen 721 painCave.isFatal = 1;
318     simError();
319     }
320 tim 1024
321     //read frameData
322     readFrameProperties(inputStream);
323    
324     //read StuntDoubles
325     readStuntDoubles(inputStream);
326    
327     inputStream.getline(buffer, bufferSize);
328     line = buffer;
329 gezelter 1752
330     if (line.find("<SiteData>") != std::string::npos) {
331     //read SiteData
332     readSiteData(inputStream);
333     } else {
334     if (line.find("</Snapshot>") == std::string::npos) {
335     sprintf(painCave.errMsg,
336     "DumpReader Error: can not find </Snapshot>\n");
337     painCave.isFatal = 1;
338     simError();
339     }
340     }
341 chrisfen 721 }
342    
343 tim 1024 void DumpReader::parseDumpLine(const std::string& line) {
344    
345    
346 chrisfen 721 StringTokenizer tokenizer(line);
347     int nTokens;
348    
349     nTokens = tokenizer.countTokens();
350    
351 gezelter 1450 if (nTokens < 2) {
352 chrisfen 721 sprintf(painCave.errMsg,
353 tim 1024 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
354 chrisfen 721 painCave.isFatal = 1;
355     simError();
356     }
357 tim 1024
358     int index = tokenizer.nextTokenAsInt();
359    
360     StuntDouble* integrableObject = info_->getIOIndexToIntegrableObject(index);
361    
362     if (integrableObject == NULL) {
363     return;
364     }
365     std::string type = tokenizer.nextToken();
366     int size = type.size();
367 gezelter 1037
368 gezelter 1450 size_t found;
369    
370     if (needPos_) {
371     found = type.find("p");
372     if (found == std::string::npos) {
373     sprintf(painCave.errMsg,
374     "DumpReader Error: StuntDouble %d has no Position\n"
375     "\tField (\"p\") specified.\n%s\n", index,
376     line.c_str());
377     painCave.isFatal = 1;
378     simError();
379     }
380     }
381    
382     if (integrableObject->isDirectional()) {
383     if (needQuaternion_) {
384     found = type.find("q");
385     if (found == std::string::npos) {
386     sprintf(painCave.errMsg,
387     "DumpReader Error: Directional StuntDouble %d has no\n"
388     "\tQuaternion Field (\"q\") specified.\n%s\n", index,
389     line.c_str());
390     painCave.isFatal = 1;
391     simError();
392     }
393     }
394     }
395    
396 tim 1024 for(int i = 0; i < size; ++i) {
397     switch(type[i]) {
398    
399     case 'p': {
400     Vector3d pos;
401     pos[0] = tokenizer.nextTokenAsDouble();
402     pos[1] = tokenizer.nextTokenAsDouble();
403     pos[2] = tokenizer.nextTokenAsDouble();
404     if (needPos_) {
405     integrableObject->setPos(pos);
406     }
407     break;
408     }
409     case 'v' : {
410     Vector3d vel;
411     vel[0] = tokenizer.nextTokenAsDouble();
412     vel[1] = tokenizer.nextTokenAsDouble();
413     vel[2] = tokenizer.nextTokenAsDouble();
414     if (needVel_) {
415     integrableObject->setVel(vel);
416     }
417     break;
418     }
419    
420     case 'q' : {
421     Quat4d q;
422     if (integrableObject->isDirectional()) {
423    
424     q[0] = tokenizer.nextTokenAsDouble();
425     q[1] = tokenizer.nextTokenAsDouble();
426     q[2] = tokenizer.nextTokenAsDouble();
427     q[3] = tokenizer.nextTokenAsDouble();
428    
429     RealType qlen = q.length();
430 gezelter 1390 if (qlen < OpenMD::epsilon) { //check quaternion is not equal to 0
431 tim 1024
432     sprintf(painCave.errMsg,
433     "DumpReader Error: initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2) ~ 0\n");
434     painCave.isFatal = 1;
435     simError();
436    
437     }
438    
439     q.normalize();
440     if (needQuaternion_) {
441     integrableObject->setQ(q);
442     }
443     }
444     break;
445     }
446     case 'j' : {
447     Vector3d ji;
448     if (integrableObject->isDirectional()) {
449     ji[0] = tokenizer.nextTokenAsDouble();
450     ji[1] = tokenizer.nextTokenAsDouble();
451     ji[2] = tokenizer.nextTokenAsDouble();
452     if (needAngMom_) {
453     integrableObject->setJ(ji);
454     }
455     }
456 gezelter 1037 break;
457 tim 1024 }
458     case 'f': {
459    
460     Vector3d force;
461     force[0] = tokenizer.nextTokenAsDouble();
462     force[1] = tokenizer.nextTokenAsDouble();
463     force[2] = tokenizer.nextTokenAsDouble();
464     integrableObject->setFrc(force);
465     break;
466     }
467     case 't' : {
468    
469     Vector3d torque;
470     torque[0] = tokenizer.nextTokenAsDouble();
471     torque[1] = tokenizer.nextTokenAsDouble();
472     torque[2] = tokenizer.nextTokenAsDouble();
473     integrableObject->setTrq(torque);
474     break;
475     }
476 gezelter 1629 case 'u' : {
477    
478     RealType particlePot;
479     particlePot = tokenizer.nextTokenAsDouble();
480     integrableObject->setParticlePot(particlePot);
481     break;
482     }
483 gezelter 1714 case 'c' : {
484    
485     RealType flucQPos;
486     flucQPos = tokenizer.nextTokenAsDouble();
487     integrableObject->setFlucQPos(flucQPos);
488     break;
489     }
490     case 'w' : {
491    
492     RealType flucQVel;
493     flucQVel = tokenizer.nextTokenAsDouble();
494     integrableObject->setFlucQVel(flucQVel);
495     break;
496     }
497     case 'g' : {
498    
499     RealType flucQFrc;
500     flucQFrc = tokenizer.nextTokenAsDouble();
501     integrableObject->setFlucQFrc(flucQFrc);
502     break;
503     }
504     case 'e' : {
505    
506     Vector3d eField;
507     eField[0] = tokenizer.nextTokenAsDouble();
508     eField[1] = tokenizer.nextTokenAsDouble();
509     eField[2] = tokenizer.nextTokenAsDouble();
510     integrableObject->setElectricField(eField);
511     break;
512     }
513 tim 1024 default: {
514     sprintf(painCave.errMsg,
515     "DumpReader Error: %s is an unrecognized type\n", type.c_str());
516     painCave.isFatal = 1;
517     simError();
518     break;
519     }
520    
521     }
522     }
523 chrisfen 721
524 tim 1024 }
525    
526    
527 gezelter 1752 void DumpReader::parseSiteLine(const std::string& line) {
528    
529     StringTokenizer tokenizer(line);
530     int nTokens;
531    
532     nTokens = tokenizer.countTokens();
533    
534     if (nTokens < 2) {
535     sprintf(painCave.errMsg,
536     "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
537     painCave.isFatal = 1;
538     simError();
539     }
540    
541     /**
542     * The first token is the global integrable object index.
543     */
544    
545     int index = tokenizer.nextTokenAsInt();
546     StuntDouble* integrableObject = info_->getIOIndexToIntegrableObject(index);
547     if (integrableObject == NULL) {
548     return;
549     }
550     StuntDouble* sd = integrableObject;
551    
552     /**
553     * Test to see if the next token is an integer or not. If not,
554     * we've got data on the integrable object itself. If there is an
555     * integer, we're parsing data for a site on a rigid body.
556     */
557    
558     std::string indexTest = tokenizer.peekNextToken();
559     std::istringstream i(indexTest);
560     int siteIndex;
561     if (i >> siteIndex) {
562     // chew up this token and parse as an int:
563     siteIndex = tokenizer.nextTokenAsInt();
564     RigidBody* rb = static_cast<RigidBody*>(integrableObject);
565     sd = rb->getAtoms()[siteIndex];
566     }
567    
568     /**
569     * The next token contains information on what follows.
570     */
571     std::string type = tokenizer.nextToken();
572     int size = type.size();
573    
574     for(int i = 0; i < size; ++i) {
575     switch(type[i]) {
576    
577     case 'u' : {
578    
579     RealType particlePot;
580     particlePot = tokenizer.nextTokenAsDouble();
581     sd->setParticlePot(particlePot);
582     break;
583     }
584     case 'c' : {
585    
586     RealType flucQPos;
587     flucQPos = tokenizer.nextTokenAsDouble();
588     sd->setFlucQPos(flucQPos);
589     break;
590     }
591     case 'w' : {
592    
593     RealType flucQVel;
594     flucQVel = tokenizer.nextTokenAsDouble();
595     sd->setFlucQVel(flucQVel);
596     break;
597     }
598     case 'g' : {
599    
600     RealType flucQFrc;
601     flucQFrc = tokenizer.nextTokenAsDouble();
602     sd->setFlucQFrc(flucQFrc);
603     break;
604     }
605     case 'e' : {
606    
607     Vector3d eField;
608     eField[0] = tokenizer.nextTokenAsDouble();
609     eField[1] = tokenizer.nextTokenAsDouble();
610     eField[2] = tokenizer.nextTokenAsDouble();
611     sd->setElectricField(eField);
612     break;
613     }
614     default: {
615     sprintf(painCave.errMsg,
616     "DumpReader Error: %s is an unrecognized type\n", type.c_str());
617     painCave.isFatal = 1;
618     simError();
619     break;
620     }
621     }
622     }
623     }
624    
625    
626 tim 1024 void DumpReader::readStuntDoubles(std::istream& inputStream) {
627 gezelter 1752
628 tim 1024 inputStream.getline(buffer, bufferSize);
629     std::string line(buffer);
630    
631     if (line.find("<StuntDoubles>") == std::string::npos) {
632 chrisfen 721 sprintf(painCave.errMsg,
633 tim 1024 "DumpReader Error: Missing <StuntDoubles>\n");
634 chrisfen 721 painCave.isFatal = 1;
635 tim 1024 simError();
636     }
637    
638     while(inputStream.getline(buffer, bufferSize)) {
639     line = buffer;
640    
641     if(line.find("</StuntDoubles>") != std::string::npos) {
642     break;
643     }
644    
645     parseDumpLine(line);
646     }
647    
648     }
649    
650 gezelter 1752 void DumpReader::readSiteData(std::istream& inputStream) {
651    
652     inputStream.getline(buffer, bufferSize);
653     std::string line(buffer);
654    
655     if (line.find("<SiteData>") == std::string::npos) {
656     // site data isn't required for a simulation, so skip
657     return;
658     }
659    
660     while(inputStream.getline(buffer, bufferSize)) {
661     line = buffer;
662    
663     if(line.find("</SiteData>") != std::string::npos) {
664     break;
665     }
666    
667     parseSiteLine(line);
668     }
669    
670     }
671    
672 tim 1024 void DumpReader::readFrameProperties(std::istream& inputStream) {
673    
674     Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
675     inputStream.getline(buffer, bufferSize);
676     std::string line(buffer);
677    
678     if (line.find("<FrameData>") == std::string::npos) {
679     sprintf(painCave.errMsg,
680     "DumpReader Error: Missing <FrameData>\n");
681     painCave.isFatal = 1;
682     simError();
683     }
684    
685     while(inputStream.getline(buffer, bufferSize)) {
686     line = buffer;
687    
688     if(line.find("</FrameData>") != std::string::npos) {
689     break;
690     }
691    
692     StringTokenizer tokenizer(line, " ;\t\n\r{}:,");
693     if (!tokenizer.hasMoreTokens()) {
694 chrisfen 721 sprintf(painCave.errMsg,
695 tim 1024 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
696 chrisfen 721 painCave.isFatal = 1;
697 tim 1024 simError();
698     }
699    
700     std::string propertyName = tokenizer.nextToken();
701     if (propertyName == "Time") {
702     RealType currTime = tokenizer.nextTokenAsDouble();
703     s->setTime(currTime);
704     } else if (propertyName == "Hmat"){
705     Mat3x3d hmat;
706     hmat(0, 0) = tokenizer.nextTokenAsDouble();
707     hmat(0, 1) = tokenizer.nextTokenAsDouble();
708     hmat(0, 2) = tokenizer.nextTokenAsDouble();
709     hmat(1, 0) = tokenizer.nextTokenAsDouble();
710     hmat(1, 1) = tokenizer.nextTokenAsDouble();
711     hmat(1, 2) = tokenizer.nextTokenAsDouble();
712     hmat(2, 0) = tokenizer.nextTokenAsDouble();
713     hmat(2, 1) = tokenizer.nextTokenAsDouble();
714     hmat(2, 2) = tokenizer.nextTokenAsDouble();
715     s->setHmat(hmat);
716     } else if (propertyName == "Thermostat") {
717     RealType chi = tokenizer.nextTokenAsDouble();
718     RealType integralOfChiDt = tokenizer.nextTokenAsDouble();
719     s->setChi(chi);
720     s->setIntegralOfChiDt(integralOfChiDt);
721     } else if (propertyName == "Barostat") {
722     Mat3x3d eta;
723     eta(0, 0) = tokenizer.nextTokenAsDouble();
724     eta(0, 1) = tokenizer.nextTokenAsDouble();
725     eta(0, 2) = tokenizer.nextTokenAsDouble();
726     eta(1, 0) = tokenizer.nextTokenAsDouble();
727     eta(1, 1) = tokenizer.nextTokenAsDouble();
728     eta(1, 2) = tokenizer.nextTokenAsDouble();
729     eta(2, 0) = tokenizer.nextTokenAsDouble();
730     eta(2, 1) = tokenizer.nextTokenAsDouble();
731     eta(2, 2) = tokenizer.nextTokenAsDouble();
732     s->setEta(eta);
733     } else {
734     sprintf(painCave.errMsg,
735     "DumpReader Error: %s is an invalid property in <FrameData>\n", propertyName.c_str());
736     painCave.isFatal = 0;
737     simError();
738     }
739    
740     }
741    
742     }
743    
744 chrisfen 721
745 gezelter 1390 }//end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date