ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/DumpReader.cpp
Revision: 1798
Committed: Thu Sep 13 14:10:11 2012 UTC (12 years, 7 months ago) by gezelter
File size: 21443 byte(s)
Log Message:
Merged trunk changes into the development branch

File Contents

# Content
1 /*
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] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40 * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 */
42
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 #include "brains/Thermo.hpp"
62
63 #ifdef IS_MPI
64 #include <mpi.h>
65 #endif
66
67
68 namespace OpenMD {
69
70 DumpReader::DumpReader(SimInfo* info, const std::string& filename)
71 : info_(info), filename_(filename), isScanned_(false), nframes_(0), needCOMprops_(false) {
72
73 #ifdef IS_MPI
74
75 if (worldRank == 0) {
76 #endif
77
78 inFile_ = new std::ifstream(filename_.c_str(),
79 ifstream::in | ifstream::binary);
80
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 #ifdef IS_MPI
90
91 }
92
93 strcpy(checkPointMsg, "Dump file opened for reading successfully.");
94 errorCheckPoint();
95
96 #endif
97
98 return;
99 }
100
101 DumpReader::~DumpReader() {
102
103 #ifdef IS_MPI
104
105 if (worldRank == 0) {
106 #endif
107
108 delete inFile_;
109
110 #ifdef IS_MPI
111
112 }
113
114 strcpy(checkPointMsg, "Dump file closed successfully.");
115 errorCheckPoint();
116
117 #endif
118
119 return;
120 }
121
122 int DumpReader::getNFrames(void) {
123
124 if (!isScanned_)
125 scanFile();
126
127 return nframes_;
128 }
129
130 void DumpReader::scanFile(void) {
131 int lineNo = 0;
132 std::streampos prevPos;
133 std::streampos currPos;
134
135 #ifdef IS_MPI
136
137 if (worldRank == 0) {
138 #endif // is_mpi
139
140 currPos = inFile_->tellg();
141 prevPos = currPos;
142 bool foundOpenSnapshotTag = false;
143 bool foundClosedSnapshotTag = false;
144 bool foundOpenSiteDataTag = 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 sprintf(painCave.errMsg,
153 "DumpReader:<Snapshot> is multiply nested at line %d in %s \n", lineNo,
154 filename_.c_str());
155 painCave.isFatal = 1;
156 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 simError();
169 }
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 nframes_ = framePos_.size();
195
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 #ifdef IS_MPI
203 }
204
205 MPI::COMM_WORLD.Bcast(&nframes_, 1, MPI::INT, 0);
206
207 #endif // is_mpi
208
209 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 ||
231 storageLayout & DataStorage::dslDipole ||
232 storageLayout & DataStorage::dslQuadrupole) {
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
246 if (needCOMprops_) {
247 Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
248 Thermo thermo(info_);
249 Vector3d com;
250
251 if (needPos_ && needVel_) {
252 Vector3d comvel;
253 Vector3d comw;
254 thermo.getComAll(com, comvel);
255 comw = thermo.getAngularMomentum();
256 } else {
257 com = thermo.getCom();
258 }
259 }
260 }
261
262 void DumpReader::readSet(int whichFrame) {
263 std::string line;
264
265 #ifndef IS_MPI
266 inFile_->clear();
267 inFile_->seekg(framePos_[whichFrame]);
268
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::COMM_WORLD.Bcast(&sendBufferSize, 1, MPI::INT, masterNode);
292 MPI::COMM_WORLD.Bcast((void *)sendBuffer.c_str(), sendBufferSize,
293 MPI::CHAR, masterNode);
294
295 sstream.str(sendBuffer);
296 } else {
297 int sendBufferSize;
298 MPI::COMM_WORLD.Bcast(&sendBufferSize, 1, MPI::INT, masterNode);
299 char * recvBuffer = new char[sendBufferSize+1];
300 assert(recvBuffer);
301 recvBuffer[sendBufferSize] = '\0';
302 MPI::COMM_WORLD.Bcast(recvBuffer, sendBufferSize, MPI::CHAR, masterNode);
303 sstream.str(recvBuffer);
304 delete [] recvBuffer;
305 }
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 sprintf(painCave.errMsg,
315 "DumpReader Error: can not find <Snapshot>\n");
316 painCave.isFatal = 1;
317 simError();
318 }
319
320 //read frameData
321 readFrameProperties(inputStream);
322
323 //read StuntDoubles
324 readStuntDoubles(inputStream);
325
326 inputStream.getline(buffer, bufferSize);
327 line = buffer;
328
329 if (line.find("<SiteData>") != std::string::npos) {
330 //read SiteData
331 readSiteData(inputStream);
332 } else {
333 if (line.find("</Snapshot>") == std::string::npos) {
334 sprintf(painCave.errMsg,
335 "DumpReader Error: can not find </Snapshot>\n");
336 painCave.isFatal = 1;
337 simError();
338 }
339 }
340 }
341
342 void DumpReader::parseDumpLine(const std::string& line) {
343
344
345 StringTokenizer tokenizer(line);
346 int nTokens;
347
348 nTokens = tokenizer.countTokens();
349
350 if (nTokens < 2) {
351 sprintf(painCave.errMsg,
352 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
353 painCave.isFatal = 1;
354 simError();
355 }
356
357 int index = tokenizer.nextTokenAsInt();
358
359 StuntDouble* sd = info_->getIOIndexToIntegrableObject(index);
360
361 if (sd == NULL) {
362 return;
363 }
364 std::string type = tokenizer.nextToken();
365 int size = type.size();
366
367 size_t found;
368
369 if (needPos_) {
370 found = type.find("p");
371 if (found == std::string::npos) {
372 sprintf(painCave.errMsg,
373 "DumpReader Error: StuntDouble %d has no Position\n"
374 "\tField (\"p\") specified.\n%s\n", index,
375 line.c_str());
376 painCave.isFatal = 1;
377 simError();
378 }
379 }
380
381 if (sd->isDirectional()) {
382 if (needQuaternion_) {
383 found = type.find("q");
384 if (found == std::string::npos) {
385 sprintf(painCave.errMsg,
386 "DumpReader Error: Directional StuntDouble %d has no\n"
387 "\tQuaternion Field (\"q\") specified.\n%s\n", index,
388 line.c_str());
389 painCave.isFatal = 1;
390 simError();
391 }
392 }
393 }
394
395 for(int i = 0; i < size; ++i) {
396 switch(type[i]) {
397
398 case 'p': {
399 Vector3d pos;
400 pos[0] = tokenizer.nextTokenAsDouble();
401 pos[1] = tokenizer.nextTokenAsDouble();
402 pos[2] = tokenizer.nextTokenAsDouble();
403 if (needPos_) {
404 sd->setPos(pos);
405 }
406 break;
407 }
408 case 'v' : {
409 Vector3d vel;
410 vel[0] = tokenizer.nextTokenAsDouble();
411 vel[1] = tokenizer.nextTokenAsDouble();
412 vel[2] = tokenizer.nextTokenAsDouble();
413 if (needVel_) {
414 sd->setVel(vel);
415 }
416 break;
417 }
418
419 case 'q' : {
420 Quat4d q;
421 if (sd->isDirectional()) {
422
423 q[0] = tokenizer.nextTokenAsDouble();
424 q[1] = tokenizer.nextTokenAsDouble();
425 q[2] = tokenizer.nextTokenAsDouble();
426 q[3] = tokenizer.nextTokenAsDouble();
427
428 RealType qlen = q.length();
429 if (qlen < OpenMD::epsilon) { //check quaternion is not equal to 0
430
431 sprintf(painCave.errMsg,
432 "DumpReader Error: initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2) ~ 0\n");
433 painCave.isFatal = 1;
434 simError();
435
436 }
437
438 q.normalize();
439 if (needQuaternion_) {
440 sd->setQ(q);
441 }
442 }
443 break;
444 }
445 case 'j' : {
446 Vector3d ji;
447 if (sd->isDirectional()) {
448 ji[0] = tokenizer.nextTokenAsDouble();
449 ji[1] = tokenizer.nextTokenAsDouble();
450 ji[2] = tokenizer.nextTokenAsDouble();
451 if (needAngMom_) {
452 sd->setJ(ji);
453 }
454 }
455 break;
456 }
457 case 'f': {
458
459 Vector3d force;
460 force[0] = tokenizer.nextTokenAsDouble();
461 force[1] = tokenizer.nextTokenAsDouble();
462 force[2] = tokenizer.nextTokenAsDouble();
463 sd->setFrc(force);
464 break;
465 }
466 case 't' : {
467
468 Vector3d torque;
469 torque[0] = tokenizer.nextTokenAsDouble();
470 torque[1] = tokenizer.nextTokenAsDouble();
471 torque[2] = tokenizer.nextTokenAsDouble();
472 sd->setTrq(torque);
473 break;
474 }
475 case 'u' : {
476
477 RealType particlePot;
478 particlePot = tokenizer.nextTokenAsDouble();
479 sd->setParticlePot(particlePot);
480 break;
481 }
482 case 'c' : {
483
484 RealType flucQPos;
485 flucQPos = tokenizer.nextTokenAsDouble();
486 sd->setFlucQPos(flucQPos);
487 break;
488 }
489 case 'w' : {
490
491 RealType flucQVel;
492 flucQVel = tokenizer.nextTokenAsDouble();
493 sd->setFlucQVel(flucQVel);
494 break;
495 }
496 case 'g' : {
497
498 RealType flucQFrc;
499 flucQFrc = tokenizer.nextTokenAsDouble();
500 sd->setFlucQFrc(flucQFrc);
501 break;
502 }
503 case 'e' : {
504
505 Vector3d eField;
506 eField[0] = tokenizer.nextTokenAsDouble();
507 eField[1] = tokenizer.nextTokenAsDouble();
508 eField[2] = tokenizer.nextTokenAsDouble();
509 sd->setElectricField(eField);
510 break;
511 }
512 default: {
513 sprintf(painCave.errMsg,
514 "DumpReader Error: %s is an unrecognized type\n", type.c_str());
515 painCave.isFatal = 1;
516 simError();
517 break;
518 }
519
520 }
521 }
522
523 }
524
525
526 void DumpReader::parseSiteLine(const std::string& line) {
527
528 StringTokenizer tokenizer(line);
529 int nTokens;
530
531 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 RigidBody* rb = static_cast<RigidBody*>(sd);
563 sd = rb->getAtoms()[siteIndex];
564 }
565
566 /**
567 * The next token contains information on what follows.
568 */
569 std::string type = tokenizer.nextToken();
570 int size = type.size();
571
572 for(int i = 0; i < size; ++i) {
573 switch(type[i]) {
574
575 case 'u' : {
576
577 RealType particlePot;
578 particlePot = tokenizer.nextTokenAsDouble();
579 sd->setParticlePot(particlePot);
580 break;
581 }
582 case 'c' : {
583
584 RealType flucQPos;
585 flucQPos = tokenizer.nextTokenAsDouble();
586 sd->setFlucQPos(flucQPos);
587 break;
588 }
589 case 'w' : {
590
591 RealType flucQVel;
592 flucQVel = tokenizer.nextTokenAsDouble();
593 sd->setFlucQVel(flucQVel);
594 break;
595 }
596 case 'g' : {
597
598 RealType flucQFrc;
599 flucQFrc = tokenizer.nextTokenAsDouble();
600 sd->setFlucQFrc(flucQFrc);
601 break;
602 }
603 case 'e' : {
604
605 Vector3d eField;
606 eField[0] = tokenizer.nextTokenAsDouble();
607 eField[1] = tokenizer.nextTokenAsDouble();
608 eField[2] = tokenizer.nextTokenAsDouble();
609 sd->setElectricField(eField);
610 break;
611 }
612 default: {
613 sprintf(painCave.errMsg,
614 "DumpReader Error: %s is an unrecognized type\n", type.c_str());
615 painCave.isFatal = 1;
616 simError();
617 break;
618 }
619 }
620 }
621 }
622
623
624 void DumpReader::readStuntDoubles(std::istream& inputStream) {
625
626 inputStream.getline(buffer, bufferSize);
627 std::string line(buffer);
628
629 if (line.find("<StuntDoubles>") == std::string::npos) {
630 sprintf(painCave.errMsg,
631 "DumpReader Error: Missing <StuntDoubles>\n");
632 painCave.isFatal = 1;
633 simError();
634 }
635
636 while(inputStream.getline(buffer, bufferSize)) {
637 line = buffer;
638
639 if(line.find("</StuntDoubles>") != std::string::npos) {
640 break;
641 }
642
643 parseDumpLine(line);
644 }
645
646 }
647
648 void DumpReader::readSiteData(std::istream& inputStream) {
649
650 inputStream.getline(buffer, bufferSize);
651 std::string line(buffer);
652
653 if (line.find("<SiteData>") == std::string::npos) {
654 // site data isn't required for a simulation, so skip
655 return;
656 }
657
658 while(inputStream.getline(buffer, bufferSize)) {
659 line = buffer;
660
661 if(line.find("</SiteData>") != std::string::npos) {
662 break;
663 }
664
665 parseSiteLine(line);
666 }
667
668 }
669
670 void DumpReader::readFrameProperties(std::istream& inputStream) {
671
672 Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
673 inputStream.getline(buffer, bufferSize);
674 std::string line(buffer);
675
676 if (line.find("<FrameData>") == std::string::npos) {
677 sprintf(painCave.errMsg,
678 "DumpReader Error: Missing <FrameData>\n");
679 painCave.isFatal = 1;
680 simError();
681 }
682
683 while(inputStream.getline(buffer, bufferSize)) {
684 line = buffer;
685
686 if(line.find("</FrameData>") != std::string::npos) {
687 break;
688 }
689
690 StringTokenizer tokenizer(line, " ;\t\n\r{}:,");
691 if (!tokenizer.hasMoreTokens()) {
692 sprintf(painCave.errMsg,
693 "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
694 painCave.isFatal = 1;
695 simError();
696 }
697
698 std::string propertyName = tokenizer.nextToken();
699 if (propertyName == "Time") {
700 RealType currTime = tokenizer.nextTokenAsDouble();
701 s->setTime(currTime);
702 } else if (propertyName == "Hmat"){
703 Mat3x3d hmat;
704 hmat(0, 0) = tokenizer.nextTokenAsDouble();
705 hmat(0, 1) = tokenizer.nextTokenAsDouble();
706 hmat(0, 2) = tokenizer.nextTokenAsDouble();
707 hmat(1, 0) = tokenizer.nextTokenAsDouble();
708 hmat(1, 1) = tokenizer.nextTokenAsDouble();
709 hmat(1, 2) = tokenizer.nextTokenAsDouble();
710 hmat(2, 0) = tokenizer.nextTokenAsDouble();
711 hmat(2, 1) = tokenizer.nextTokenAsDouble();
712 hmat(2, 2) = tokenizer.nextTokenAsDouble();
713 s->setHmat(hmat);
714 } else if (propertyName == "Thermostat") {
715 pair<RealType, RealType> thermostat;
716 thermostat.first = tokenizer.nextTokenAsDouble();
717 thermostat.second = tokenizer.nextTokenAsDouble();
718 s->setThermostat(thermostat);
719 } else if (propertyName == "Barostat") {
720 Mat3x3d eta;
721 eta(0, 0) = tokenizer.nextTokenAsDouble();
722 eta(0, 1) = tokenizer.nextTokenAsDouble();
723 eta(0, 2) = tokenizer.nextTokenAsDouble();
724 eta(1, 0) = tokenizer.nextTokenAsDouble();
725 eta(1, 1) = tokenizer.nextTokenAsDouble();
726 eta(1, 2) = tokenizer.nextTokenAsDouble();
727 eta(2, 0) = tokenizer.nextTokenAsDouble();
728 eta(2, 1) = tokenizer.nextTokenAsDouble();
729 eta(2, 2) = tokenizer.nextTokenAsDouble();
730 s->setBarostat(eta);
731 } else {
732 sprintf(painCave.errMsg,
733 "DumpReader Error: %s is an invalid property in <FrameData>\n", propertyName.c_str());
734 painCave.isFatal = 0;
735 simError();
736 }
737
738 }
739
740 }
741
742
743 }//end namespace OpenMD

Properties

Name Value
svn:keywords Author Id Revision Date