ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-3.0/src/io/DumpReader.cpp
Revision: 1885
Committed: Tue Dec 14 21:57:47 2004 UTC (19 years, 9 months ago) by tim
File size: 16293 byte(s)
Log Message:
fix a bug in DumpReader. Master nodes does not brocast the total number of frames

File Contents

# User Rev Content
1 gezelter 1490 #define _LARGEFILE_SOURCE64
2     #define _FILE_OFFSET_BITS 64
3    
4     #include <sys/types.h>
5     #include <sys/stat.h>
6    
7     #include <iostream>
8     #include <math.h>
9    
10     #include <stdio.h>
11     #include <stdlib.h>
12     #include <string.h>
13    
14 tim 1739 #include "io/DumpReader.hpp"
15 tim 1808 #include "primitives/Molecule.hpp"
16 tim 1492 #include "utils/simError.h"
17 tim 1784 #include "utils/MemoryUtils.hpp"
18 tim 1808 #include "utils/StringTokenizer.hpp"
19    
20 tim 1739 #ifdef IS_MPI
21 gezelter 1490
22     #include <mpi.h>
23     #define TAKE_THIS_TAG_CHAR 0
24     #define TAKE_THIS_TAG_INT 1
25 tim 1739
26 gezelter 1490 #endif // is_mpi
27    
28    
29 tim 1739 namespace oopse {
30 gezelter 1490
31 tim 1739 DumpReader::DumpReader(SimInfo* info, const std::string& filename)
32 tim 1885 : info_(info), filename_(filename), isScanned_(false), nframes_(0) {
33 gezelter 1490
34     #ifdef IS_MPI
35 tim 1739
36     if (worldRank == 0) {
37 gezelter 1490 #endif
38    
39 tim 1739 inFile_ = fopen(filename_.c_str(), "r");
40    
41     if (inFile_ == NULL) {
42 tim 1784 sprintf(painCave.errMsg, "Cannot open file: %s\n", filename_.c_str());
43 tim 1739 painCave.isFatal = 1;
44     simError();
45     }
46    
47 gezelter 1490 #ifdef IS_MPI
48 tim 1739
49     }
50    
51     strcpy(checkPointMsg, "Dump file opened for reading successfully.");
52     MPIcheckPoint();
53    
54 gezelter 1490 #endif
55 tim 1739
56     return;
57 gezelter 1490 }
58    
59 tim 1739 DumpReader::~DumpReader() {
60    
61 gezelter 1490 #ifdef IS_MPI
62 tim 1739
63     if (worldRank == 0) {
64 gezelter 1490 #endif
65    
66 tim 1739 int error;
67     error = fclose(inFile_);
68 gezelter 1490
69 tim 1739 if (error) {
70     sprintf(painCave.errMsg, "Error closing %s\n", filename_.c_str());
71     painCave.isFatal = 1;
72     simError();
73     }
74    
75     MemoryUtils::deleteVectorOfPointer(framePos_);
76    
77 gezelter 1490 #ifdef IS_MPI
78 tim 1739
79     }
80    
81     strcpy(checkPointMsg, "Dump file closed successfully.");
82     MPIcheckPoint();
83    
84 gezelter 1490 #endif
85    
86 tim 1739 return;
87 gezelter 1490 }
88    
89 tim 1739 int DumpReader::getNFrames(void) {
90 tim 1885 #ifdef IS_MPI
91 gezelter 1490
92 tim 1885 if (worldRank == 0) {
93     #endif
94    
95     if (!isScanned_)
96     scanFile();
97    
98     return nframes_;
99     #ifdef IS_MPI
100     }
101     #endif
102 gezelter 1490 }
103    
104 tim 1739 void DumpReader::scanFile(void) {
105 gezelter 1490 int i, j;
106     int lineNum = 0;
107 tim 1739 char readBuffer[maxBufferSize];
108     fpos_t * currPos;
109 gezelter 1490
110     #ifdef IS_MPI
111    
112 tim 1739 if (worldRank == 0) {
113 gezelter 1490 #endif // is_mpi
114    
115 tim 1739 rewind(inFile_);
116 gezelter 1490
117 tim 1739 currPos = new fpos_t;
118     fgetpos(inFile_, currPos);
119     fgets(readBuffer, sizeof(readBuffer), inFile_);
120     lineNum++;
121 gezelter 1490
122 tim 1739 if (feof(inFile_)) {
123     sprintf(painCave.errMsg,
124     "File \"%s\" ended unexpectedly at line %d\n",
125     filename_.c_str(),
126     lineNum);
127     painCave.isFatal = 1;
128     simError();
129     }
130 gezelter 1490
131 tim 1739 while (!feof(inFile_)) {
132     framePos_.push_back(currPos);
133 gezelter 1490
134 tim 1739 i = atoi(readBuffer);
135 gezelter 1490
136 tim 1739 fgets(readBuffer, sizeof(readBuffer), inFile_);
137     lineNum++;
138 gezelter 1490
139 tim 1739 if (feof(inFile_)) {
140     sprintf(painCave.errMsg,
141     "File \"%s\" ended unexpectedly at line %d\n",
142     filename_.c_str(),
143     lineNum);
144     painCave.isFatal = 1;
145     simError();
146     }
147 gezelter 1490
148 tim 1739 for(j = 0; j < i; j++) {
149     fgets(readBuffer, sizeof(readBuffer), inFile_);
150     lineNum++;
151 gezelter 1490
152 tim 1739 if (feof(inFile_)) {
153     sprintf(painCave.errMsg,
154     "File \"%s\" ended unexpectedly at line %d,"
155     " with atom %d\n", filename_.c_str(),
156     lineNum,
157     j);
158 gezelter 1490
159 tim 1739 painCave.isFatal = 1;
160     simError();
161     }
162     }
163 gezelter 1490
164 tim 1739 currPos = new fpos_t;
165     fgetpos(inFile_, currPos);
166     fgets(readBuffer, sizeof(readBuffer), inFile_);
167     lineNum++;
168     }
169 gezelter 1490
170 tim 1739 delete currPos;
171     rewind(inFile_);
172 gezelter 1490
173 tim 1784 isScanned_ = true;
174 gezelter 1490
175 tim 1885 nframes_ = framePos_.size();
176 tim 1739 #ifdef IS_MPI
177 tim 1885 MPI_Bcast(&nframes_, 1, MPI_INT, 0, MPI_COMM_WORLD);
178 gezelter 1490
179 tim 1739 }
180 gezelter 1490
181 tim 1739 strcpy(checkPointMsg, "Successfully scanned DumpFile\n");
182     MPIcheckPoint();
183 gezelter 1490
184 tim 1739 #endif // is_mpi
185 gezelter 1490
186 tim 1739 }
187 gezelter 1490
188 tim 1739 void DumpReader::readFrame(int whichFrame) {
189     readSet(whichFrame);
190     }
191 gezelter 1490
192 tim 1739 void DumpReader::readSet(int whichFrame) {
193     int i;
194     int nTotObjs; // the number of atoms
195     char read_buffer[maxBufferSize]; //the line buffer for reading
196     char * eof_test; // ptr to see when we reach the end of the file
197 gezelter 1490
198 tim 1739 Molecule* mol;
199     StuntDouble* integrableObject;
200 tim 1808 SimInfo::MoleculeIterator mi;
201     Molecule::IntegrableObjectIterator ii;
202 gezelter 1490
203 tim 1739 #ifndef IS_MPI
204 gezelter 1490
205 tim 1739 fsetpos(inFile_, framePos_[whichFrame]);
206     eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
207 gezelter 1490
208 tim 1739 if (eof_test == NULL) {
209 gezelter 1490 sprintf(painCave.errMsg,
210 tim 1739 "DumpReader error: error reading 1st line of \"%s\"\n",
211     filename_.c_str());
212 gezelter 1490 painCave.isFatal = 1;
213     simError();
214     }
215    
216 tim 1739 nTotObjs = atoi(read_buffer);
217 gezelter 1490
218 tim 1739 if (nTotObjs != info_->getNGlobalIntegrableObjects()) {
219     sprintf(painCave.errMsg,
220     "DumpReader error. %s nIntegrable, %d, "
221     "does not match the meta-data file's nIntegrable, %d.\n",
222     filename_.c_str(),
223     nTotObjs,
224     info_->getNIntegrableObjects());
225 gezelter 1490
226 tim 1739 painCave.isFatal = 1;
227     simError();
228     }
229 gezelter 1490
230 tim 1739 //read the box mat from the comment line
231 gezelter 1490
232 tim 1739 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
233 gezelter 1490
234 tim 1739 if (eof_test == NULL) {
235     sprintf(painCave.errMsg, "error in reading commment in %s\n",
236     filename_.c_str());
237     painCave.isFatal = 1;
238     simError();
239 gezelter 1490 }
240    
241 tim 1808 parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
242 gezelter 1490
243 tim 1739 //parse dump lines
244 gezelter 1490
245 tim 1739 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
246 gezelter 1490
247 tim 1808 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
248     integrableObject = mol->nextIntegrableObject(ii)) {
249 gezelter 1490
250 tim 1739 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
251 gezelter 1490
252 tim 1739 if (eof_test == NULL) {
253     sprintf(painCave.errMsg,
254     "error in reading file %s\n"
255     "natoms = %d; index = %d\n"
256     "error reading the line from the file.\n",
257     filename_.c_str(),
258     nTotObjs,
259     i);
260 gezelter 1490
261 tim 1739 painCave.isFatal = 1;
262     simError();
263     }
264 gezelter 1490
265 tim 1808 parseDumpLine(read_buffer, integrableObject);
266 tim 1739
267 tim 1808 }
268 gezelter 1490 }
269    
270 tim 1739 // MPI Section of code..........
271 gezelter 1490
272 tim 1739 #else //IS_MPI
273 gezelter 1490
274 tim 1739 // first thing first, suspend fatalities.
275     int masterNode = 0;
276     int nCurObj;
277     painCave.isEventLoop = 1;
278 gezelter 1490
279 tim 1739 int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
280     int haveError;
281 gezelter 1490
282 tim 1739 MPI_Status istatus;
283     int nitems;
284 gezelter 1490
285 tim 1739 nTotObjs = info_->getNGlobalIntegrableObjects();
286     haveError = 0;
287 gezelter 1490
288 tim 1739 if (worldRank == masterNode) {
289     fsetpos(inFile_, framePos_[whichFrame]);
290 gezelter 1490
291 tim 1739 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
292 gezelter 1490
293 tim 1739 if (eof_test == NULL) {
294     sprintf(painCave.errMsg, "Error reading 1st line of %s \n ",
295     filename_.c_str());
296     painCave.isFatal = 1;
297     simError();
298     }
299 gezelter 1490
300 tim 1739 nitems = atoi(read_buffer);
301 gezelter 1490
302 tim 1739 // Check to see that the number of integrable objects in the
303     // intial configuration file is the same as derived from the
304     // meta-data file.
305 gezelter 1490
306 tim 1739 if (nTotObjs != nitems) {
307     sprintf(painCave.errMsg,
308     "DumpReader Error. %s nIntegrable, %d, "
309     "does not match the meta-data file's nIntegrable, %d.\n",
310     filename_.c_str(),
311     nTotObjs,
312 tim 1883 info_->getNGlobalIntegrableObjects());
313 gezelter 1490
314 tim 1739 painCave.isFatal = 1;
315     simError();
316 gezelter 1490 }
317    
318 tim 1739 //read the boxMat from the comment line
319 gezelter 1490
320 tim 1739 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
321 gezelter 1490
322 tim 1739 if (eof_test == NULL) {
323     sprintf(painCave.errMsg, "error in reading commment in %s\n",
324     filename_.c_str());
325     painCave.isFatal = 1;
326     simError();
327 gezelter 1490 }
328    
329 tim 1739 //Every single processor will parse the comment line by itself
330     //By using this way, we might lose some efficiency, but if we want to add
331     //more parameters into comment line, we only need to modify function
332     //parseCommentLine
333 gezelter 1490
334 tim 1739 MPI_Bcast(read_buffer, maxBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
335 tim 1808 parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
336 gezelter 1490
337 tim 1739 for(i = 0; i < info_->getNGlobalMolecules(); i++) {
338 tim 1883 int which_node = info_->getMolToProc(i);
339 gezelter 1490
340 tim 1739 if (which_node == masterNode) {
341     //molecules belong to master node
342 gezelter 1490
343 tim 1739 mol = info_->getMoleculeByGlobalIndex(i);
344 gezelter 1490
345 tim 1739 if (mol == NULL) {
346 tim 1884 sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
347 tim 1739 painCave.isFatal = 1;
348     simError();
349     }
350 gezelter 1490
351 tim 1739 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
352     integrableObject = mol->nextIntegrableObject(ii)){
353    
354     eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
355 gezelter 1490
356 tim 1739 if (eof_test == NULL) {
357     sprintf(painCave.errMsg,
358     "error in reading file %s\n"
359     "natoms = %d; index = %d\n"
360     "error reading the line from the file.\n",
361     filename_.c_str(),
362     nTotObjs,
363     i);
364 gezelter 1490
365 tim 1739 painCave.isFatal = 1;
366     simError();
367     }
368 gezelter 1490
369 tim 1739 parseDumpLine(read_buffer, integrableObject);
370     }
371     } else {
372     //molecule belongs to slave nodes
373 gezelter 1490
374 tim 1739 MPI_Recv(&nCurObj, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT,
375     MPI_COMM_WORLD, &istatus);
376 gezelter 1490
377 tim 1739 for(int j = 0; j < nCurObj; j++) {
378     eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
379 gezelter 1490
380 tim 1739 if (eof_test == NULL) {
381     sprintf(painCave.errMsg,
382     "error in reading file %s\n"
383     "natoms = %d; index = %d\n"
384     "error reading the line from the file.\n",
385     filename_.c_str(),
386     nTotObjs,
387     i);
388 gezelter 1490
389 tim 1739 painCave.isFatal = 1;
390     simError();
391     }
392    
393     MPI_Send(read_buffer, maxBufferSize, MPI_CHAR, which_node,
394     TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD);
395     }
396     }
397     }
398     } else {
399     //actions taken at slave nodes
400     MPI_Bcast(read_buffer, maxBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
401 gezelter 1490
402 tim 1739 /**@todo*/
403 tim 1808 parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
404 gezelter 1490
405 tim 1739 for(i = 0; i < info_->getNGlobalMolecules(); i++) {
406 tim 1883 int which_node = info_->getMolToProc(i);
407 gezelter 1490
408 tim 1739 if (which_node == worldRank) {
409     //molecule with global index i belongs to this processor
410    
411     mol = info_->getMoleculeByGlobalIndex(i);
412     if (mol == NULL) {
413 tim 1884 sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
414 tim 1739 painCave.isFatal = 1;
415     simError();
416     }
417    
418     nCurObj = mol->getNIntegrableObjects();
419 gezelter 1490
420 tim 1739 MPI_Send(&nCurObj, 1, MPI_INT, masterNode, TAKE_THIS_TAG_INT,
421     MPI_COMM_WORLD);
422 gezelter 1490
423 tim 1739 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
424     integrableObject = mol->nextIntegrableObject(ii)){
425    
426     MPI_Recv(read_buffer, maxBufferSize, MPI_CHAR, masterNode,
427     TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD, &istatus);
428 gezelter 1490
429 tim 1739 parseDumpLine(read_buffer, integrableObject);
430     }
431    
432     }
433    
434     }
435    
436     }
437 gezelter 1490
438 tim 1739 #endif
439 gezelter 1490
440 tim 1739 }
441 gezelter 1490
442 tim 1739 void DumpReader::parseDumpLine(char *line, StuntDouble *integrableObject) {
443 gezelter 1490
444 tim 1739 Vector3d pos; // position place holders
445     Vector3d vel; // velocity placeholders
446     Quat4d q; // the quaternions
447     Vector3d ji; // angular velocity placeholders;
448     StringTokenizer tokenizer(line);
449     int nTokens;
450    
451     nTokens = tokenizer.countTokens();
452 gezelter 1490
453 tim 1808 if (nTokens < 14) {
454     sprintf(painCave.errMsg,
455     "Not enough Tokens.\n");
456     painCave.isFatal = 1;
457     simError();
458 gezelter 1490 }
459    
460 tim 1739 std::string name = tokenizer.nextToken();
461 gezelter 1490
462 tim 1739 if (name != integrableObject->getType()) {
463    
464 gezelter 1490 }
465    
466 tim 1739 pos[0] = tokenizer.nextTokenAsFloat();
467     pos[1] = tokenizer.nextTokenAsFloat();
468     pos[2] = tokenizer.nextTokenAsFloat();
469     integrableObject->setPos(pos);
470    
471     vel[0] = tokenizer.nextTokenAsFloat();
472     vel[1] = tokenizer.nextTokenAsFloat();
473     vel[2] = tokenizer.nextTokenAsFloat();
474     integrableObject->setVel(vel);
475 gezelter 1490
476 tim 1739 if (integrableObject->isDirectional()) {
477 tim 1808
478 tim 1739 q[0] = tokenizer.nextTokenAsFloat();
479     q[1] = tokenizer.nextTokenAsFloat();
480     q[2] = tokenizer.nextTokenAsFloat();
481     q[3] = tokenizer.nextTokenAsFloat();
482 gezelter 1490
483 tim 1739 double qlen = q.length();
484     if (qlen < oopse::epsilon) { //check quaternion is not equal to 0
485    
486     sprintf(painCave.errMsg,
487     "initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2 ~ 0).\n");
488     painCave.isFatal = 1;
489     simError();
490    
491 tim 1864 }
492    
493     q.normalize();
494    
495 tim 1739 integrableObject->setQ(q);
496    
497     ji[0] = tokenizer.nextTokenAsFloat();
498     ji[1] = tokenizer.nextTokenAsFloat();
499     ji[2] = tokenizer.nextTokenAsFloat();
500     integrableObject->setJ(ji);
501 gezelter 1490 }
502    
503     }
504    
505    
506 tim 1808 void DumpReader::parseCommentLine(char* line, Snapshot* s) {
507 tim 1739 double currTime;
508     Mat3x3d hmat;
509     double chi;
510     double integralOfChiDt;
511     Mat3x3d eta;
512 gezelter 1490
513 tim 1739 StringTokenizer tokenizer(line);
514     int nTokens;
515 gezelter 1490
516 tim 1739 nTokens = tokenizer.countTokens();
517 gezelter 1490
518 tim 1739 //comment line should at least contain 10 tokens: current time(1 token) and h-matrix(9 tokens)
519     if (nTokens < 10) {
520     sprintf(painCave.errMsg,
521     "Not enough tokens in comment line: %s", line);
522     painCave.isFatal = 1;
523     simError();
524 gezelter 1490 }
525    
526 tim 1739 //read current time
527     currTime = tokenizer.nextTokenAsFloat();
528 tim 1808 s->setTime(currTime);
529 gezelter 1490
530 tim 1739 //read h-matrix
531     hmat(0, 0) = tokenizer.nextTokenAsFloat();
532     hmat(0, 1) = tokenizer.nextTokenAsFloat();
533     hmat(0, 2) = tokenizer.nextTokenAsFloat();
534     hmat(1, 0) = tokenizer.nextTokenAsFloat();
535     hmat(1, 1) = tokenizer.nextTokenAsFloat();
536     hmat(1, 2) = tokenizer.nextTokenAsFloat();
537     hmat(2, 0) = tokenizer.nextTokenAsFloat();
538     hmat(2, 1) = tokenizer.nextTokenAsFloat();
539     hmat(2, 2) = tokenizer.nextTokenAsFloat();
540     s->setHmat(hmat);
541 gezelter 1490
542 tim 1808 //read chi and integrablOfChidt, they should apprear in pair
543     if (tokenizer.countTokens() >= 2) {
544     chi = tokenizer.nextTokenAsFloat();
545     integralOfChiDt = tokenizer.nextTokenAsFloat();
546 gezelter 1490
547 tim 1808 s->setChi(chi);
548     s->setIntegralOfChiDt(integralOfChiDt);
549     }
550    
551     //read eta (eta is 3x3 matrix)
552     if (tokenizer.countTokens() >= 9) {
553     eta(0, 0) = tokenizer.nextTokenAsFloat();
554     eta(0, 1) = tokenizer.nextTokenAsFloat();
555     eta(0, 2) = tokenizer.nextTokenAsFloat();
556     eta(1, 0) = tokenizer.nextTokenAsFloat();
557     eta(1, 1) = tokenizer.nextTokenAsFloat();
558     eta(1, 2) = tokenizer.nextTokenAsFloat();
559     eta(2, 0) = tokenizer.nextTokenAsFloat();
560     eta(2, 1) = tokenizer.nextTokenAsFloat();
561     eta(2, 2) = tokenizer.nextTokenAsFloat();
562 gezelter 1490
563 tim 1808 s->setEta(eta);
564 tim 1739 }
565 gezelter 1490
566 tim 1739
567 gezelter 1490 }
568    
569 tim 1739 }//end namespace oopse