ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/DumpReader.cpp
Revision: 1787
Committed: Wed Aug 29 18:13:11 2012 UTC (12 years, 8 months ago) by gezelter
File size: 21369 byte(s)
Log Message:
Massive multipole rewrite

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

Properties

Name Value
svn:keywords Author Id Revision Date