ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/io/DumpReader.cpp
Revision: 1964
Committed: Thu Jan 16 16:00:43 2014 UTC (11 years, 3 months ago) by gezelter
File size: 21424 byte(s)
Log Message:
Fixing some leaks

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 gezelter 1879 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
39 gezelter 1782 * [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 gezelter 1938 #ifdef IS_MPI
47     #include <mpi.h>
48     #endif
49    
50 chrisfen 721 #include <sys/types.h>
51     #include <sys/stat.h>
52    
53     #include <iostream>
54     #include <math.h>
55    
56     #include <stdio.h>
57     #include <stdlib.h>
58     #include <string.h>
59    
60     #include "io/DumpReader.hpp"
61     #include "primitives/Molecule.hpp"
62     #include "utils/simError.h"
63     #include "utils/MemoryUtils.hpp"
64     #include "utils/StringTokenizer.hpp"
65 gezelter 1782 #include "brains/Thermo.hpp"
66 chrisfen 721
67    
68 gezelter 1390 namespace OpenMD {
69 chrisfen 721
70     DumpReader::DumpReader(SimInfo* info, const std::string& filename)
71 gezelter 1104 : info_(info), filename_(filename), isScanned_(false), nframes_(0), needCOMprops_(false) {
72 chrisfen 996
73 chrisfen 721 #ifdef IS_MPI
74 chrisfen 996
75     if (worldRank == 0) {
76 chrisfen 721 #endif
77 chrisfen 996
78 gezelter 1790 inFile_ = new std::ifstream(filename_.c_str(),
79     ifstream::in | ifstream::binary);
80 chrisfen 996
81     if (inFile_->fail()) {
82     sprintf(painCave.errMsg,
83     "DumpReader: Cannot open file: %s\n",
84     filename_.c_str());
85     painCave.isFatal = 1;
86     simError();
87     }
88    
89 chrisfen 721 #ifdef IS_MPI
90 chrisfen 996
91     }
92    
93     strcpy(checkPointMsg, "Dump file opened for reading successfully.");
94 gezelter 1241 errorCheckPoint();
95 chrisfen 996
96 chrisfen 721 #endif
97 chrisfen 996
98     return;
99     }
100    
101 chrisfen 721 DumpReader::~DumpReader() {
102 chrisfen 996
103 chrisfen 721 #ifdef IS_MPI
104 chrisfen 996
105 chrisfen 721 if (worldRank == 0) {
106     #endif
107 chrisfen 996
108 chrisfen 721 delete inFile_;
109 chrisfen 996
110 chrisfen 721 #ifdef IS_MPI
111 chrisfen 996
112 chrisfen 721 }
113 chrisfen 996
114 chrisfen 721 strcpy(checkPointMsg, "Dump file closed successfully.");
115 gezelter 1241 errorCheckPoint();
116 chrisfen 996
117 chrisfen 721 #endif
118 chrisfen 996
119 chrisfen 721 return;
120     }
121 chrisfen 996
122 chrisfen 721 int DumpReader::getNFrames(void) {
123    
124     if (!isScanned_)
125     scanFile();
126    
127     return nframes_;
128     }
129    
130     void DumpReader::scanFile(void) {
131 gezelter 1879
132 tim 1024 std::streampos prevPos;
133 chrisfen 721 std::streampos currPos;
134 tim 1024
135 chrisfen 721 #ifdef IS_MPI
136 tim 1024
137 chrisfen 721 if (worldRank == 0) {
138     #endif // is_mpi
139 tim 1024
140     currPos = inFile_->tellg();
141     prevPos = currPos;
142     bool foundOpenSnapshotTag = false;
143     bool foundClosedSnapshotTag = false;
144 gezelter 1793
145 gezelter 1879 int lineNo = 0;
146 tim 1024 while(inFile_->getline(buffer, bufferSize)) {
147     ++lineNo;
148    
149     std::string line = buffer;
150     currPos = inFile_->tellg();
151     if (line.find("<Snapshot>")!= std::string::npos) {
152     if (foundOpenSnapshotTag) {
153 chrisfen 721 sprintf(painCave.errMsg,
154 tim 1024 "DumpReader:<Snapshot> is multiply nested at line %d in %s \n", lineNo,
155     filename_.c_str());
156 chrisfen 721 painCave.isFatal = 1;
157 tim 1024 simError();
158     }
159     foundOpenSnapshotTag = true;
160     foundClosedSnapshotTag = false;
161     framePos_.push_back(prevPos);
162    
163     } else if (line.find("</Snapshot>") != std::string::npos){
164     if (!foundOpenSnapshotTag) {
165     sprintf(painCave.errMsg,
166     "DumpReader:</Snapshot> appears before <Snapshot> at line %d in %s \n", lineNo,
167     filename_.c_str());
168     painCave.isFatal = 1;
169 chrisfen 721 simError();
170 tim 1024 }
171    
172     if (foundClosedSnapshotTag) {
173     sprintf(painCave.errMsg,
174     "DumpReader:</Snapshot> appears multiply nested at line %d in %s \n", lineNo,
175     filename_.c_str());
176     painCave.isFatal = 1;
177     simError();
178     }
179     foundClosedSnapshotTag = true;
180     foundOpenSnapshotTag = false;
181     }
182     prevPos = currPos;
183     }
184    
185     // only found <Snapshot> for the last frame means the file is corrupted, we should discard
186     // it and give a warning message
187     if (foundOpenSnapshotTag) {
188     sprintf(painCave.errMsg,
189     "DumpReader: last frame in %s is invalid\n", filename_.c_str());
190     painCave.isFatal = 0;
191     simError();
192     framePos_.pop_back();
193     }
194    
195 chrisfen 721 nframes_ = framePos_.size();
196 tim 1024
197     if (nframes_ == 0) {
198     sprintf(painCave.errMsg,
199     "DumpReader: %s does not contain a valid frame\n", filename_.c_str());
200     painCave.isFatal = 1;
201     simError();
202     }
203 chrisfen 721 #ifdef IS_MPI
204     }
205    
206 gezelter 1796 MPI::COMM_WORLD.Bcast(&nframes_, 1, MPI::INT, 0);
207 tim 1024
208 chrisfen 721 #endif // is_mpi
209 tim 1024
210 chrisfen 721 isScanned_ = true;
211     }
212    
213     void DumpReader::readFrame(int whichFrame) {
214     if (!isScanned_)
215     scanFile();
216    
217     int storageLayout = info_->getSnapshotManager()->getStorageLayout();
218    
219     if (storageLayout & DataStorage::dslPosition) {
220     needPos_ = true;
221     } else {
222     needPos_ = false;
223     }
224    
225     if (storageLayout & DataStorage::dslVelocity) {
226     needVel_ = true;
227     } else {
228     needVel_ = false;
229     }
230    
231 gezelter 1879 if (storageLayout & DataStorage::dslAmat ||
232     storageLayout & DataStorage::dslDipole ||
233     storageLayout & DataStorage::dslQuadrupole) {
234 chrisfen 721 needQuaternion_ = true;
235     } else {
236     needQuaternion_ = false;
237     }
238    
239     if (storageLayout & DataStorage::dslAngularMomentum) {
240     needAngMom_ = true;
241     } else {
242     needAngMom_ = false;
243     }
244    
245     readSet(whichFrame);
246 gezelter 1104
247     if (needCOMprops_) {
248     Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
249 gezelter 1782 Thermo thermo(info_);
250 gezelter 1104 Vector3d com;
251 gezelter 1782
252     if (needPos_ && needVel_) {
253     Vector3d comvel;
254     Vector3d comw;
255     thermo.getComAll(com, comvel);
256     comw = thermo.getAngularMomentum();
257     } else {
258     com = thermo.getCom();
259     }
260 gezelter 1104 }
261 chrisfen 721 }
262    
263 tim 1024 void DumpReader::readSet(int whichFrame) {
264     std::string line;
265    
266 chrisfen 721 #ifndef IS_MPI
267     inFile_->clear();
268     inFile_->seekg(framePos_[whichFrame]);
269 tim 1024
270     std::istream& inputStream = *inFile_;
271    
272     #else
273     int masterNode = 0;
274     std::stringstream sstream;
275     if (worldRank == masterNode) {
276     std::string sendBuffer;
277    
278     inFile_->clear();
279     inFile_->seekg(framePos_[whichFrame]);
280    
281     while (inFile_->getline(buffer, bufferSize)) {
282    
283     line = buffer;
284     sendBuffer += line;
285     sendBuffer += '\n';
286     if (line.find("</Snapshot>") != std::string::npos) {
287     break;
288     }
289     }
290    
291     int sendBufferSize = sendBuffer.size();
292 gezelter 1796 MPI::COMM_WORLD.Bcast(&sendBufferSize, 1, MPI::INT, masterNode);
293     MPI::COMM_WORLD.Bcast((void *)sendBuffer.c_str(), sendBufferSize,
294     MPI::CHAR, masterNode);
295 tim 1024
296     sstream.str(sendBuffer);
297     } else {
298     int sendBufferSize;
299 gezelter 1796 MPI::COMM_WORLD.Bcast(&sendBufferSize, 1, MPI::INT, masterNode);
300 tim 1024 char * recvBuffer = new char[sendBufferSize+1];
301 gezelter 1313 assert(recvBuffer);
302     recvBuffer[sendBufferSize] = '\0';
303 gezelter 1796 MPI::COMM_WORLD.Bcast(recvBuffer, sendBufferSize, MPI::CHAR, masterNode);
304 tim 1024 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 1782
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 gezelter 1782 StuntDouble* sd = info_->getIOIndexToIntegrableObject(index);
361 tim 1024
362 gezelter 1782 if (sd == NULL) {
363 tim 1024 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 gezelter 1782 if (sd->isDirectional()) {
383 gezelter 1450 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 gezelter 1782 sd->setPos(pos);
406 tim 1024 }
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 gezelter 1782 sd->setVel(vel);
416 tim 1024 }
417     break;
418     }
419    
420     case 'q' : {
421     Quat4d q;
422 gezelter 1782 if (sd->isDirectional()) {
423 tim 1024
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 gezelter 1782 sd->setQ(q);
442 tim 1024 }
443     }
444     break;
445     }
446     case 'j' : {
447     Vector3d ji;
448 gezelter 1782 if (sd->isDirectional()) {
449 tim 1024 ji[0] = tokenizer.nextTokenAsDouble();
450     ji[1] = tokenizer.nextTokenAsDouble();
451     ji[2] = tokenizer.nextTokenAsDouble();
452     if (needAngMom_) {
453 gezelter 1782 sd->setJ(ji);
454 tim 1024 }
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 gezelter 1782 sd->setFrc(force);
465 tim 1024 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 gezelter 1782 sd->setTrq(torque);
474 tim 1024 break;
475     }
476 gezelter 1610 case 'u' : {
477    
478     RealType particlePot;
479     particlePot = tokenizer.nextTokenAsDouble();
480 gezelter 1782 sd->setParticlePot(particlePot);
481 gezelter 1610 break;
482     }
483 gezelter 1782 case 'c' : {
484    
485     RealType flucQPos;
486     flucQPos = tokenizer.nextTokenAsDouble();
487     sd->setFlucQPos(flucQPos);
488     break;
489     }
490     case 'w' : {
491    
492     RealType flucQVel;
493     flucQVel = tokenizer.nextTokenAsDouble();
494     sd->setFlucQVel(flucQVel);
495     break;
496     }
497     case 'g' : {
498    
499     RealType flucQFrc;
500     flucQFrc = tokenizer.nextTokenAsDouble();
501     sd->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     sd->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     }
524    
525    
526 gezelter 1782 void DumpReader::parseSiteLine(const std::string& line) {
527    
528     StringTokenizer tokenizer(line);
529     int nTokens;
530 gezelter 1879
531 gezelter 1782 nTokens = tokenizer.countTokens();
532    
533     if (nTokens < 2) {
534     sprintf(painCave.errMsg,
535     "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
536     painCave.isFatal = 1;
537     simError();
538     }
539    
540     /**
541     * The first token is the global integrable object index.
542     */
543    
544     int index = tokenizer.nextTokenAsInt();
545     StuntDouble* sd = info_->getIOIndexToIntegrableObject(index);
546     if (sd == NULL) {
547     return;
548     }
549    
550     /**
551     * Test to see if the next token is an integer or not. If not,
552     * we've got data on the integrable object itself. If there is an
553     * integer, we're parsing data for a site on a rigid body.
554     */
555    
556     std::string indexTest = tokenizer.peekNextToken();
557     std::istringstream i(indexTest);
558     int siteIndex;
559     if (i >> siteIndex) {
560     // chew up this token and parse as an int:
561     siteIndex = tokenizer.nextTokenAsInt();
562 gezelter 1897
563     if (sd->isRigidBody()) {
564     RigidBody* rb = static_cast<RigidBody*>(sd);
565     sd = rb->getAtoms()[siteIndex];
566     }
567 gezelter 1782 }
568    
569     /**
570     * The next token contains information on what follows.
571     */
572     std::string type = tokenizer.nextToken();
573     int size = type.size();
574    
575     for(int i = 0; i < size; ++i) {
576     switch(type[i]) {
577    
578     case 'u' : {
579    
580     RealType particlePot;
581     particlePot = tokenizer.nextTokenAsDouble();
582     sd->setParticlePot(particlePot);
583     break;
584     }
585     case 'c' : {
586    
587     RealType flucQPos;
588     flucQPos = tokenizer.nextTokenAsDouble();
589     sd->setFlucQPos(flucQPos);
590     break;
591     }
592     case 'w' : {
593    
594     RealType flucQVel;
595     flucQVel = tokenizer.nextTokenAsDouble();
596     sd->setFlucQVel(flucQVel);
597     break;
598     }
599     case 'g' : {
600    
601     RealType flucQFrc;
602     flucQFrc = tokenizer.nextTokenAsDouble();
603     sd->setFlucQFrc(flucQFrc);
604     break;
605     }
606     case 'e' : {
607    
608     Vector3d eField;
609     eField[0] = tokenizer.nextTokenAsDouble();
610     eField[1] = tokenizer.nextTokenAsDouble();
611     eField[2] = tokenizer.nextTokenAsDouble();
612     sd->setElectricField(eField);
613     break;
614     }
615     default: {
616     sprintf(painCave.errMsg,
617     "DumpReader Error: %s is an unrecognized type\n", type.c_str());
618     painCave.isFatal = 1;
619     simError();
620     break;
621     }
622     }
623     }
624     }
625    
626    
627 tim 1024 void DumpReader::readStuntDoubles(std::istream& inputStream) {
628 gezelter 1782
629 tim 1024 inputStream.getline(buffer, bufferSize);
630     std::string line(buffer);
631    
632     if (line.find("<StuntDoubles>") == std::string::npos) {
633 chrisfen 721 sprintf(painCave.errMsg,
634 tim 1024 "DumpReader Error: Missing <StuntDoubles>\n");
635 chrisfen 721 painCave.isFatal = 1;
636 tim 1024 simError();
637     }
638    
639     while(inputStream.getline(buffer, bufferSize)) {
640     line = buffer;
641    
642     if(line.find("</StuntDoubles>") != std::string::npos) {
643     break;
644     }
645    
646     parseDumpLine(line);
647     }
648    
649     }
650    
651 gezelter 1782 void DumpReader::readSiteData(std::istream& inputStream) {
652    
653     std::string line(buffer);
654 gezelter 1879
655     // We already found the starting <SiteData> tag or we wouldn't be
656     // here, so just start parsing until we get to the ending
657     // </SiteData> tag:
658 gezelter 1782
659     while(inputStream.getline(buffer, bufferSize)) {
660     line = buffer;
661    
662     if(line.find("</SiteData>") != std::string::npos) {
663     break;
664     }
665    
666     parseSiteLine(line);
667     }
668    
669     }
670    
671 tim 1024 void DumpReader::readFrameProperties(std::istream& inputStream) {
672    
673     Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
674     inputStream.getline(buffer, bufferSize);
675     std::string line(buffer);
676    
677     if (line.find("<FrameData>") == std::string::npos) {
678     sprintf(painCave.errMsg,
679     "DumpReader Error: Missing <FrameData>\n");
680     painCave.isFatal = 1;
681     simError();
682     }
683    
684     while(inputStream.getline(buffer, bufferSize)) {
685     line = buffer;
686    
687     if(line.find("</FrameData>") != std::string::npos) {
688     break;
689     }
690    
691     StringTokenizer tokenizer(line, " ;\t\n\r{}:,");
692     if (!tokenizer.hasMoreTokens()) {
693 chrisfen 721 sprintf(painCave.errMsg,
694 tim 1024 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
695 chrisfen 721 painCave.isFatal = 1;
696 tim 1024 simError();
697     }
698    
699     std::string propertyName = tokenizer.nextToken();
700     if (propertyName == "Time") {
701     RealType currTime = tokenizer.nextTokenAsDouble();
702     s->setTime(currTime);
703     } else if (propertyName == "Hmat"){
704     Mat3x3d hmat;
705     hmat(0, 0) = tokenizer.nextTokenAsDouble();
706     hmat(0, 1) = tokenizer.nextTokenAsDouble();
707     hmat(0, 2) = tokenizer.nextTokenAsDouble();
708     hmat(1, 0) = tokenizer.nextTokenAsDouble();
709     hmat(1, 1) = tokenizer.nextTokenAsDouble();
710     hmat(1, 2) = tokenizer.nextTokenAsDouble();
711     hmat(2, 0) = tokenizer.nextTokenAsDouble();
712     hmat(2, 1) = tokenizer.nextTokenAsDouble();
713     hmat(2, 2) = tokenizer.nextTokenAsDouble();
714     s->setHmat(hmat);
715     } else if (propertyName == "Thermostat") {
716 gezelter 1782 pair<RealType, RealType> thermostat;
717     thermostat.first = tokenizer.nextTokenAsDouble();
718     thermostat.second = tokenizer.nextTokenAsDouble();
719     s->setThermostat(thermostat);
720 tim 1024 } else if (propertyName == "Barostat") {
721     Mat3x3d eta;
722     eta(0, 0) = tokenizer.nextTokenAsDouble();
723     eta(0, 1) = tokenizer.nextTokenAsDouble();
724     eta(0, 2) = tokenizer.nextTokenAsDouble();
725     eta(1, 0) = tokenizer.nextTokenAsDouble();
726     eta(1, 1) = tokenizer.nextTokenAsDouble();
727     eta(1, 2) = tokenizer.nextTokenAsDouble();
728     eta(2, 0) = tokenizer.nextTokenAsDouble();
729     eta(2, 1) = tokenizer.nextTokenAsDouble();
730     eta(2, 2) = tokenizer.nextTokenAsDouble();
731 gezelter 1782 s->setBarostat(eta);
732 tim 1024 } else {
733     sprintf(painCave.errMsg,
734     "DumpReader Error: %s is an invalid property in <FrameData>\n", propertyName.c_str());
735     painCave.isFatal = 0;
736     simError();
737     }
738    
739     }
740    
741     }
742    
743 chrisfen 721
744 gezelter 1390 }//end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date