ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/DumpReader.cpp
Revision: 1665
Committed: Tue Nov 22 20:38:56 2011 UTC (13 years, 5 months ago) by gezelter
File size: 17441 byte(s)
Log Message:
updated copyright notices

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     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     MPI_Bcast(&nframes_, 1, MPI_INT, 0, MPI_COMM_WORLD);
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     if (storageLayout & DataStorage::dslAmat || storageLayout & DataStorage::dslElectroFrame) {
232     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     Vector3d com;
248     Vector3d comvel;
249     Vector3d comw;
250 chuckv 1112 if (needPos_ && needVel_){
251     info_->getComAll(com, comvel);
252     comw = info_->getAngularMomentum();
253     }else{
254     com = info_->getCom();
255     comvel = 0.0;
256     comw = 0.0;
257     }
258 gezelter 1104 s->setCOMprops(com, comvel, comw);
259     }
260    
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     MPI_Bcast(&sendBufferSize, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
293     MPI_Bcast((void *)sendBuffer.c_str(), sendBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
294    
295     sstream.str(sendBuffer);
296     } else {
297     int sendBufferSize;
298     MPI_Bcast(&sendBufferSize, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
299     char * recvBuffer = new char[sendBufferSize+1];
300 gezelter 1313 assert(recvBuffer);
301     recvBuffer[sendBufferSize] = '\0';
302 tim 1024 MPI_Bcast(recvBuffer, sendBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
303     sstream.str(recvBuffer);
304 gezelter 1313 delete [] recvBuffer;
305 tim 1024 }
306    
307     std::istream& inputStream = sstream;
308     #endif
309    
310     inputStream.getline(buffer, bufferSize);
311    
312     line = buffer;
313     if (line.find("<Snapshot>") == std::string::npos) {
314 chrisfen 721 sprintf(painCave.errMsg,
315 tim 1024 "DumpReader Error: can not find <Snapshot>\n");
316 chrisfen 721 painCave.isFatal = 1;
317     simError();
318     }
319 tim 1024
320     //read frameData
321     readFrameProperties(inputStream);
322    
323     //read StuntDoubles
324     readStuntDoubles(inputStream);
325    
326     inputStream.getline(buffer, bufferSize);
327     line = buffer;
328     if (line.find("</Snapshot>") == std::string::npos) {
329 chrisfen 721 sprintf(painCave.errMsg,
330 tim 1024 "DumpReader Error: can not find </Snapshot>\n");
331 chrisfen 721 painCave.isFatal = 1;
332     simError();
333 tim 1024 }
334    
335 chrisfen 721 }
336    
337 tim 1024 void DumpReader::parseDumpLine(const std::string& line) {
338    
339    
340 chrisfen 721 StringTokenizer tokenizer(line);
341     int nTokens;
342    
343     nTokens = tokenizer.countTokens();
344    
345 gezelter 1450 if (nTokens < 2) {
346 chrisfen 721 sprintf(painCave.errMsg,
347 tim 1024 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
348 chrisfen 721 painCave.isFatal = 1;
349     simError();
350     }
351 tim 1024
352     int index = tokenizer.nextTokenAsInt();
353    
354     StuntDouble* integrableObject = info_->getIOIndexToIntegrableObject(index);
355    
356     if (integrableObject == NULL) {
357     return;
358     }
359     std::string type = tokenizer.nextToken();
360     int size = type.size();
361 gezelter 1037
362 gezelter 1450 size_t found;
363    
364     if (needPos_) {
365     found = type.find("p");
366     if (found == std::string::npos) {
367     sprintf(painCave.errMsg,
368     "DumpReader Error: StuntDouble %d has no Position\n"
369     "\tField (\"p\") specified.\n%s\n", index,
370     line.c_str());
371     painCave.isFatal = 1;
372     simError();
373     }
374     }
375    
376     if (integrableObject->isDirectional()) {
377     if (needQuaternion_) {
378     found = type.find("q");
379     if (found == std::string::npos) {
380     sprintf(painCave.errMsg,
381     "DumpReader Error: Directional StuntDouble %d has no\n"
382     "\tQuaternion Field (\"q\") specified.\n%s\n", index,
383     line.c_str());
384     painCave.isFatal = 1;
385     simError();
386     }
387     }
388     }
389    
390 tim 1024 for(int i = 0; i < size; ++i) {
391     switch(type[i]) {
392    
393     case 'p': {
394     Vector3d pos;
395     pos[0] = tokenizer.nextTokenAsDouble();
396     pos[1] = tokenizer.nextTokenAsDouble();
397     pos[2] = tokenizer.nextTokenAsDouble();
398     if (needPos_) {
399     integrableObject->setPos(pos);
400     }
401     break;
402     }
403     case 'v' : {
404     Vector3d vel;
405     vel[0] = tokenizer.nextTokenAsDouble();
406     vel[1] = tokenizer.nextTokenAsDouble();
407     vel[2] = tokenizer.nextTokenAsDouble();
408     if (needVel_) {
409     integrableObject->setVel(vel);
410     }
411     break;
412     }
413    
414     case 'q' : {
415     Quat4d q;
416     if (integrableObject->isDirectional()) {
417    
418     q[0] = tokenizer.nextTokenAsDouble();
419     q[1] = tokenizer.nextTokenAsDouble();
420     q[2] = tokenizer.nextTokenAsDouble();
421     q[3] = tokenizer.nextTokenAsDouble();
422    
423     RealType qlen = q.length();
424 gezelter 1390 if (qlen < OpenMD::epsilon) { //check quaternion is not equal to 0
425 tim 1024
426     sprintf(painCave.errMsg,
427     "DumpReader Error: initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2) ~ 0\n");
428     painCave.isFatal = 1;
429     simError();
430    
431     }
432    
433     q.normalize();
434     if (needQuaternion_) {
435     integrableObject->setQ(q);
436     }
437     }
438     break;
439     }
440     case 'j' : {
441     Vector3d ji;
442     if (integrableObject->isDirectional()) {
443     ji[0] = tokenizer.nextTokenAsDouble();
444     ji[1] = tokenizer.nextTokenAsDouble();
445     ji[2] = tokenizer.nextTokenAsDouble();
446     if (needAngMom_) {
447     integrableObject->setJ(ji);
448     }
449     }
450 gezelter 1037 break;
451 tim 1024 }
452     case 'f': {
453    
454     Vector3d force;
455     force[0] = tokenizer.nextTokenAsDouble();
456     force[1] = tokenizer.nextTokenAsDouble();
457     force[2] = tokenizer.nextTokenAsDouble();
458     integrableObject->setFrc(force);
459     break;
460     }
461     case 't' : {
462    
463     Vector3d torque;
464     torque[0] = tokenizer.nextTokenAsDouble();
465     torque[1] = tokenizer.nextTokenAsDouble();
466     torque[2] = tokenizer.nextTokenAsDouble();
467     integrableObject->setTrq(torque);
468     break;
469     }
470 gezelter 1629 case 'u' : {
471    
472     RealType particlePot;
473     particlePot = tokenizer.nextTokenAsDouble();
474     integrableObject->setParticlePot(particlePot);
475     break;
476     }
477 tim 1024 default: {
478     sprintf(painCave.errMsg,
479     "DumpReader Error: %s is an unrecognized type\n", type.c_str());
480     painCave.isFatal = 1;
481     simError();
482     break;
483     }
484    
485     }
486     }
487 chrisfen 721
488 tim 1024 }
489    
490    
491     void DumpReader::readStuntDoubles(std::istream& inputStream) {
492    
493     inputStream.getline(buffer, bufferSize);
494     std::string line(buffer);
495    
496     if (line.find("<StuntDoubles>") == std::string::npos) {
497 chrisfen 721 sprintf(painCave.errMsg,
498 tim 1024 "DumpReader Error: Missing <StuntDoubles>\n");
499 chrisfen 721 painCave.isFatal = 1;
500 tim 1024 simError();
501     }
502    
503     while(inputStream.getline(buffer, bufferSize)) {
504     line = buffer;
505    
506     if(line.find("</StuntDoubles>") != std::string::npos) {
507     break;
508     }
509    
510     parseDumpLine(line);
511     }
512    
513     }
514    
515     void DumpReader::readFrameProperties(std::istream& inputStream) {
516    
517     Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
518     inputStream.getline(buffer, bufferSize);
519     std::string line(buffer);
520    
521     if (line.find("<FrameData>") == std::string::npos) {
522     sprintf(painCave.errMsg,
523     "DumpReader Error: Missing <FrameData>\n");
524     painCave.isFatal = 1;
525     simError();
526     }
527    
528     while(inputStream.getline(buffer, bufferSize)) {
529     line = buffer;
530    
531     if(line.find("</FrameData>") != std::string::npos) {
532     break;
533     }
534    
535     StringTokenizer tokenizer(line, " ;\t\n\r{}:,");
536     if (!tokenizer.hasMoreTokens()) {
537 chrisfen 721 sprintf(painCave.errMsg,
538 tim 1024 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
539 chrisfen 721 painCave.isFatal = 1;
540 tim 1024 simError();
541     }
542    
543     std::string propertyName = tokenizer.nextToken();
544     if (propertyName == "Time") {
545     RealType currTime = tokenizer.nextTokenAsDouble();
546     s->setTime(currTime);
547     } else if (propertyName == "Hmat"){
548     Mat3x3d hmat;
549     hmat(0, 0) = tokenizer.nextTokenAsDouble();
550     hmat(0, 1) = tokenizer.nextTokenAsDouble();
551     hmat(0, 2) = tokenizer.nextTokenAsDouble();
552     hmat(1, 0) = tokenizer.nextTokenAsDouble();
553     hmat(1, 1) = tokenizer.nextTokenAsDouble();
554     hmat(1, 2) = tokenizer.nextTokenAsDouble();
555     hmat(2, 0) = tokenizer.nextTokenAsDouble();
556     hmat(2, 1) = tokenizer.nextTokenAsDouble();
557     hmat(2, 2) = tokenizer.nextTokenAsDouble();
558     s->setHmat(hmat);
559     } else if (propertyName == "Thermostat") {
560     RealType chi = tokenizer.nextTokenAsDouble();
561     RealType integralOfChiDt = tokenizer.nextTokenAsDouble();
562     s->setChi(chi);
563     s->setIntegralOfChiDt(integralOfChiDt);
564     } else if (propertyName == "Barostat") {
565     Mat3x3d eta;
566     eta(0, 0) = tokenizer.nextTokenAsDouble();
567     eta(0, 1) = tokenizer.nextTokenAsDouble();
568     eta(0, 2) = tokenizer.nextTokenAsDouble();
569     eta(1, 0) = tokenizer.nextTokenAsDouble();
570     eta(1, 1) = tokenizer.nextTokenAsDouble();
571     eta(1, 2) = tokenizer.nextTokenAsDouble();
572     eta(2, 0) = tokenizer.nextTokenAsDouble();
573     eta(2, 1) = tokenizer.nextTokenAsDouble();
574     eta(2, 2) = tokenizer.nextTokenAsDouble();
575     s->setEta(eta);
576     } else {
577     sprintf(painCave.errMsg,
578     "DumpReader Error: %s is an invalid property in <FrameData>\n", propertyName.c_str());
579     painCave.isFatal = 0;
580     simError();
581     }
582    
583     }
584    
585     }
586    
587 chrisfen 721
588 gezelter 1390 }//end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date