ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/DumpReader.cpp
Revision: 1442
Committed: Mon May 10 17:28:26 2010 UTC (14 years, 11 months ago) by gezelter
Original Path: trunk/src/io/DumpReader.cpp
File size: 16321 byte(s)
Log Message:
Adding property set to svn entries

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     * [4] Vardeman & Gezelter, in progress (2009).
40     */
41 chrisfen 721
42     #define _LARGEFILE_SOURCE64
43     #define _FILE_OFFSET_BITS 64
44    
45     #include <sys/types.h>
46     #include <sys/stat.h>
47    
48     #include <iostream>
49     #include <math.h>
50    
51     #include <stdio.h>
52     #include <stdlib.h>
53     #include <string.h>
54    
55     #include "io/DumpReader.hpp"
56     #include "primitives/Molecule.hpp"
57     #include "utils/simError.h"
58     #include "utils/MemoryUtils.hpp"
59     #include "utils/StringTokenizer.hpp"
60    
61     #ifdef IS_MPI
62    
63     #include <mpi.h>
64     #define TAKE_THIS_TAG_CHAR 0
65     #define TAKE_THIS_TAG_INT 1
66    
67     #endif // is_mpi
68    
69    
70 gezelter 1390 namespace OpenMD {
71 chrisfen 721
72     DumpReader::DumpReader(SimInfo* info, const std::string& filename)
73 gezelter 1104 : info_(info), filename_(filename), isScanned_(false), nframes_(0), needCOMprops_(false) {
74 chrisfen 996
75 chrisfen 721 #ifdef IS_MPI
76 chrisfen 996
77     if (worldRank == 0) {
78 chrisfen 721 #endif
79 chrisfen 996
80 chrisfen 721 inFile_ = new std::ifstream(filename_.c_str());
81 chrisfen 996
82     if (inFile_->fail()) {
83     sprintf(painCave.errMsg,
84     "DumpReader: Cannot open file: %s\n",
85     filename_.c_str());
86     painCave.isFatal = 1;
87     simError();
88     }
89    
90 chrisfen 721 #ifdef IS_MPI
91 chrisfen 996
92     }
93    
94     strcpy(checkPointMsg, "Dump file opened for reading successfully.");
95 gezelter 1241 errorCheckPoint();
96 chrisfen 996
97 chrisfen 721 #endif
98 chrisfen 996
99     return;
100     }
101    
102 chrisfen 721 DumpReader::~DumpReader() {
103 chrisfen 996
104 chrisfen 721 #ifdef IS_MPI
105 chrisfen 996
106 chrisfen 721 if (worldRank == 0) {
107     #endif
108 chrisfen 996
109 chrisfen 721 delete inFile_;
110 chrisfen 996
111 chrisfen 721 #ifdef IS_MPI
112 chrisfen 996
113 chrisfen 721 }
114 chrisfen 996
115 chrisfen 721 strcpy(checkPointMsg, "Dump file closed successfully.");
116 gezelter 1241 errorCheckPoint();
117 chrisfen 996
118 chrisfen 721 #endif
119 chrisfen 996
120 chrisfen 721 return;
121     }
122 chrisfen 996
123 chrisfen 721 int DumpReader::getNFrames(void) {
124    
125     if (!isScanned_)
126     scanFile();
127    
128     return nframes_;
129     }
130    
131     void DumpReader::scanFile(void) {
132 tim 1024 int lineNo = 0;
133     std::streampos prevPos;
134 chrisfen 721 std::streampos currPos;
135 tim 1024
136 chrisfen 721 #ifdef IS_MPI
137 tim 1024
138 chrisfen 721 if (worldRank == 0) {
139     #endif // is_mpi
140 tim 1024
141     currPos = inFile_->tellg();
142     prevPos = currPos;
143     bool foundOpenSnapshotTag = false;
144     bool foundClosedSnapshotTag = false;
145     while(inFile_->getline(buffer, bufferSize)) {
146     ++lineNo;
147    
148     std::string line = buffer;
149     currPos = inFile_->tellg();
150     if (line.find("<Snapshot>")!= std::string::npos) {
151     if (foundOpenSnapshotTag) {
152 chrisfen 721 sprintf(painCave.errMsg,
153 tim 1024 "DumpReader:<Snapshot> is multiply nested at line %d in %s \n", lineNo,
154     filename_.c_str());
155 chrisfen 721 painCave.isFatal = 1;
156 tim 1024 simError();
157     }
158     foundOpenSnapshotTag = true;
159     foundClosedSnapshotTag = false;
160     framePos_.push_back(prevPos);
161    
162     } else if (line.find("</Snapshot>") != std::string::npos){
163     if (!foundOpenSnapshotTag) {
164     sprintf(painCave.errMsg,
165     "DumpReader:</Snapshot> appears before <Snapshot> at line %d in %s \n", lineNo,
166     filename_.c_str());
167     painCave.isFatal = 1;
168 chrisfen 721 simError();
169 tim 1024 }
170    
171     if (foundClosedSnapshotTag) {
172     sprintf(painCave.errMsg,
173     "DumpReader:</Snapshot> appears multiply nested at line %d in %s \n", lineNo,
174     filename_.c_str());
175     painCave.isFatal = 1;
176     simError();
177     }
178     foundClosedSnapshotTag = true;
179     foundOpenSnapshotTag = false;
180     }
181     prevPos = currPos;
182     }
183    
184     // only found <Snapshot> for the last frame means the file is corrupted, we should discard
185     // it and give a warning message
186     if (foundOpenSnapshotTag) {
187     sprintf(painCave.errMsg,
188     "DumpReader: last frame in %s is invalid\n", filename_.c_str());
189     painCave.isFatal = 0;
190     simError();
191     framePos_.pop_back();
192     }
193    
194 chrisfen 721 nframes_ = framePos_.size();
195 tim 1024
196     if (nframes_ == 0) {
197     sprintf(painCave.errMsg,
198     "DumpReader: %s does not contain a valid frame\n", filename_.c_str());
199     painCave.isFatal = 1;
200     simError();
201     }
202 chrisfen 721 #ifdef IS_MPI
203     }
204    
205     MPI_Bcast(&nframes_, 1, MPI_INT, 0, MPI_COMM_WORLD);
206 tim 1024
207 chrisfen 721 #endif // is_mpi
208 tim 1024
209 chrisfen 721 isScanned_ = true;
210     }
211    
212     void DumpReader::readFrame(int whichFrame) {
213     if (!isScanned_)
214     scanFile();
215    
216     int storageLayout = info_->getSnapshotManager()->getStorageLayout();
217    
218     if (storageLayout & DataStorage::dslPosition) {
219     needPos_ = true;
220     } else {
221     needPos_ = false;
222     }
223    
224     if (storageLayout & DataStorage::dslVelocity) {
225     needVel_ = true;
226     } else {
227     needVel_ = false;
228     }
229    
230     if (storageLayout & DataStorage::dslAmat || storageLayout & DataStorage::dslElectroFrame) {
231     needQuaternion_ = true;
232     } else {
233     needQuaternion_ = false;
234     }
235    
236     if (storageLayout & DataStorage::dslAngularMomentum) {
237     needAngMom_ = true;
238     } else {
239     needAngMom_ = false;
240     }
241    
242     readSet(whichFrame);
243 gezelter 1104
244     if (needCOMprops_) {
245     Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
246     Vector3d com;
247     Vector3d comvel;
248     Vector3d comw;
249 chuckv 1112 if (needPos_ && needVel_){
250     info_->getComAll(com, comvel);
251     comw = info_->getAngularMomentum();
252     }else{
253     com = info_->getCom();
254     comvel = 0.0;
255     comw = 0.0;
256     }
257 gezelter 1104 s->setCOMprops(com, comvel, comw);
258     }
259    
260 chrisfen 721 }
261    
262 tim 1024 void DumpReader::readSet(int whichFrame) {
263     std::string line;
264    
265 chrisfen 721 #ifndef IS_MPI
266     inFile_->clear();
267     inFile_->seekg(framePos_[whichFrame]);
268 tim 1024
269     std::istream& inputStream = *inFile_;
270    
271     #else
272     int masterNode = 0;
273     std::stringstream sstream;
274     if (worldRank == masterNode) {
275     std::string sendBuffer;
276    
277     inFile_->clear();
278     inFile_->seekg(framePos_[whichFrame]);
279    
280     while (inFile_->getline(buffer, bufferSize)) {
281    
282     line = buffer;
283     sendBuffer += line;
284     sendBuffer += '\n';
285     if (line.find("</Snapshot>") != std::string::npos) {
286     break;
287     }
288     }
289    
290     int sendBufferSize = sendBuffer.size();
291     MPI_Bcast(&sendBufferSize, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
292     MPI_Bcast((void *)sendBuffer.c_str(), sendBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
293    
294     sstream.str(sendBuffer);
295     } else {
296     int sendBufferSize;
297     MPI_Bcast(&sendBufferSize, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
298     char * recvBuffer = new char[sendBufferSize+1];
299 gezelter 1313 assert(recvBuffer);
300     recvBuffer[sendBufferSize] = '\0';
301 tim 1024 MPI_Bcast(recvBuffer, sendBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
302     sstream.str(recvBuffer);
303 gezelter 1313 delete [] recvBuffer;
304 tim 1024 }
305    
306     std::istream& inputStream = sstream;
307     #endif
308    
309     inputStream.getline(buffer, bufferSize);
310    
311     line = buffer;
312     if (line.find("<Snapshot>") == std::string::npos) {
313 chrisfen 721 sprintf(painCave.errMsg,
314 tim 1024 "DumpReader Error: can not find <Snapshot>\n");
315 chrisfen 721 painCave.isFatal = 1;
316     simError();
317     }
318 tim 1024
319     //read frameData
320     readFrameProperties(inputStream);
321    
322     //read StuntDoubles
323     readStuntDoubles(inputStream);
324    
325     inputStream.getline(buffer, bufferSize);
326     line = buffer;
327     if (line.find("</Snapshot>") == std::string::npos) {
328 chrisfen 721 sprintf(painCave.errMsg,
329 tim 1024 "DumpReader Error: can not find </Snapshot>\n");
330 chrisfen 721 painCave.isFatal = 1;
331     simError();
332 tim 1024 }
333    
334 chrisfen 721 }
335    
336 tim 1024 void DumpReader::parseDumpLine(const std::string& line) {
337    
338    
339 chrisfen 721 StringTokenizer tokenizer(line);
340     int nTokens;
341    
342     nTokens = tokenizer.countTokens();
343    
344 tim 1024 if (nTokens < 2) {
345 chrisfen 721 sprintf(painCave.errMsg,
346 tim 1024 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
347 chrisfen 721 painCave.isFatal = 1;
348     simError();
349     }
350 tim 1024
351     int index = tokenizer.nextTokenAsInt();
352    
353     StuntDouble* integrableObject = info_->getIOIndexToIntegrableObject(index);
354    
355     if (integrableObject == NULL) {
356     return;
357     }
358     std::string type = tokenizer.nextToken();
359     int size = type.size();
360 gezelter 1037
361 tim 1024 for(int i = 0; i < size; ++i) {
362     switch(type[i]) {
363    
364     case 'p': {
365     Vector3d pos;
366     pos[0] = tokenizer.nextTokenAsDouble();
367     pos[1] = tokenizer.nextTokenAsDouble();
368     pos[2] = tokenizer.nextTokenAsDouble();
369     if (needPos_) {
370     integrableObject->setPos(pos);
371     }
372     break;
373     }
374     case 'v' : {
375     Vector3d vel;
376     vel[0] = tokenizer.nextTokenAsDouble();
377     vel[1] = tokenizer.nextTokenAsDouble();
378     vel[2] = tokenizer.nextTokenAsDouble();
379     if (needVel_) {
380     integrableObject->setVel(vel);
381     }
382     break;
383     }
384    
385     case 'q' : {
386     Quat4d q;
387     if (integrableObject->isDirectional()) {
388    
389     q[0] = tokenizer.nextTokenAsDouble();
390     q[1] = tokenizer.nextTokenAsDouble();
391     q[2] = tokenizer.nextTokenAsDouble();
392     q[3] = tokenizer.nextTokenAsDouble();
393    
394     RealType qlen = q.length();
395 gezelter 1390 if (qlen < OpenMD::epsilon) { //check quaternion is not equal to 0
396 tim 1024
397     sprintf(painCave.errMsg,
398     "DumpReader Error: initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2) ~ 0\n");
399     painCave.isFatal = 1;
400     simError();
401    
402     }
403    
404     q.normalize();
405     if (needQuaternion_) {
406     integrableObject->setQ(q);
407     }
408     }
409     break;
410     }
411     case 'j' : {
412     Vector3d ji;
413     if (integrableObject->isDirectional()) {
414     ji[0] = tokenizer.nextTokenAsDouble();
415     ji[1] = tokenizer.nextTokenAsDouble();
416     ji[2] = tokenizer.nextTokenAsDouble();
417     if (needAngMom_) {
418     integrableObject->setJ(ji);
419     }
420     }
421 gezelter 1037 break;
422 tim 1024 }
423     case 'f': {
424    
425     Vector3d force;
426     force[0] = tokenizer.nextTokenAsDouble();
427     force[1] = tokenizer.nextTokenAsDouble();
428     force[2] = tokenizer.nextTokenAsDouble();
429     integrableObject->setFrc(force);
430     break;
431     }
432     case 't' : {
433    
434     Vector3d torque;
435     torque[0] = tokenizer.nextTokenAsDouble();
436     torque[1] = tokenizer.nextTokenAsDouble();
437     torque[2] = tokenizer.nextTokenAsDouble();
438     integrableObject->setTrq(torque);
439     break;
440     }
441     default: {
442     sprintf(painCave.errMsg,
443     "DumpReader Error: %s is an unrecognized type\n", type.c_str());
444     painCave.isFatal = 1;
445     simError();
446     break;
447     }
448    
449     }
450     }
451 chrisfen 721
452 tim 1024 }
453    
454    
455     void DumpReader::readStuntDoubles(std::istream& inputStream) {
456    
457     inputStream.getline(buffer, bufferSize);
458     std::string line(buffer);
459    
460     if (line.find("<StuntDoubles>") == std::string::npos) {
461 chrisfen 721 sprintf(painCave.errMsg,
462 tim 1024 "DumpReader Error: Missing <StuntDoubles>\n");
463 chrisfen 721 painCave.isFatal = 1;
464 tim 1024 simError();
465     }
466    
467     while(inputStream.getline(buffer, bufferSize)) {
468     line = buffer;
469    
470     if(line.find("</StuntDoubles>") != std::string::npos) {
471     break;
472     }
473    
474     parseDumpLine(line);
475     }
476    
477     }
478    
479     void DumpReader::readFrameProperties(std::istream& inputStream) {
480    
481     Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
482     inputStream.getline(buffer, bufferSize);
483     std::string line(buffer);
484    
485     if (line.find("<FrameData>") == std::string::npos) {
486     sprintf(painCave.errMsg,
487     "DumpReader Error: Missing <FrameData>\n");
488     painCave.isFatal = 1;
489     simError();
490     }
491    
492     while(inputStream.getline(buffer, bufferSize)) {
493     line = buffer;
494    
495     if(line.find("</FrameData>") != std::string::npos) {
496     break;
497     }
498    
499     StringTokenizer tokenizer(line, " ;\t\n\r{}:,");
500     if (!tokenizer.hasMoreTokens()) {
501 chrisfen 721 sprintf(painCave.errMsg,
502 tim 1024 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
503 chrisfen 721 painCave.isFatal = 1;
504 tim 1024 simError();
505     }
506    
507     std::string propertyName = tokenizer.nextToken();
508     if (propertyName == "Time") {
509     RealType currTime = tokenizer.nextTokenAsDouble();
510     s->setTime(currTime);
511     } else if (propertyName == "Hmat"){
512     Mat3x3d hmat;
513     hmat(0, 0) = tokenizer.nextTokenAsDouble();
514     hmat(0, 1) = tokenizer.nextTokenAsDouble();
515     hmat(0, 2) = tokenizer.nextTokenAsDouble();
516     hmat(1, 0) = tokenizer.nextTokenAsDouble();
517     hmat(1, 1) = tokenizer.nextTokenAsDouble();
518     hmat(1, 2) = tokenizer.nextTokenAsDouble();
519     hmat(2, 0) = tokenizer.nextTokenAsDouble();
520     hmat(2, 1) = tokenizer.nextTokenAsDouble();
521     hmat(2, 2) = tokenizer.nextTokenAsDouble();
522     s->setHmat(hmat);
523     } else if (propertyName == "Thermostat") {
524     RealType chi = tokenizer.nextTokenAsDouble();
525     RealType integralOfChiDt = tokenizer.nextTokenAsDouble();
526     s->setChi(chi);
527     s->setIntegralOfChiDt(integralOfChiDt);
528     } else if (propertyName == "Barostat") {
529     Mat3x3d eta;
530     eta(0, 0) = tokenizer.nextTokenAsDouble();
531     eta(0, 1) = tokenizer.nextTokenAsDouble();
532     eta(0, 2) = tokenizer.nextTokenAsDouble();
533     eta(1, 0) = tokenizer.nextTokenAsDouble();
534     eta(1, 1) = tokenizer.nextTokenAsDouble();
535     eta(1, 2) = tokenizer.nextTokenAsDouble();
536     eta(2, 0) = tokenizer.nextTokenAsDouble();
537     eta(2, 1) = tokenizer.nextTokenAsDouble();
538     eta(2, 2) = tokenizer.nextTokenAsDouble();
539     s->setEta(eta);
540     } else {
541     sprintf(painCave.errMsg,
542     "DumpReader Error: %s is an invalid property in <FrameData>\n", propertyName.c_str());
543     painCave.isFatal = 0;
544     simError();
545     }
546    
547     }
548    
549     }
550    
551 chrisfen 721
552 gezelter 1390 }//end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date