1 |
#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 |
#include "io/DumpReader.hpp" |
15 |
#include "primitives/Molecule.hpp" |
16 |
#include "utils/simError.h" |
17 |
#include "utils/MemoryUtils.hpp" |
18 |
#include "utils/StringTokenizer.hpp" |
19 |
|
20 |
#ifdef IS_MPI |
21 |
|
22 |
#include <mpi.h> |
23 |
#define TAKE_THIS_TAG_CHAR 0 |
24 |
#define TAKE_THIS_TAG_INT 1 |
25 |
|
26 |
#endif // is_mpi |
27 |
|
28 |
|
29 |
namespace oopse { |
30 |
|
31 |
DumpReader::DumpReader(SimInfo* info, const std::string& filename) |
32 |
: info_(info), filename_(filename), isScanned_(false), nframes_(0) { |
33 |
|
34 |
#ifdef IS_MPI |
35 |
|
36 |
if (worldRank == 0) { |
37 |
#endif |
38 |
|
39 |
inFile_ = fopen(filename_.c_str(), "r"); |
40 |
|
41 |
if (inFile_ == NULL) { |
42 |
sprintf(painCave.errMsg, "Cannot open file: %s\n", filename_.c_str()); |
43 |
painCave.isFatal = 1; |
44 |
simError(); |
45 |
} |
46 |
|
47 |
#ifdef IS_MPI |
48 |
|
49 |
} |
50 |
|
51 |
strcpy(checkPointMsg, "Dump file opened for reading successfully."); |
52 |
MPIcheckPoint(); |
53 |
|
54 |
#endif |
55 |
|
56 |
return; |
57 |
} |
58 |
|
59 |
DumpReader::~DumpReader() { |
60 |
|
61 |
#ifdef IS_MPI |
62 |
|
63 |
if (worldRank == 0) { |
64 |
#endif |
65 |
|
66 |
int error; |
67 |
error = fclose(inFile_); |
68 |
|
69 |
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 |
#ifdef IS_MPI |
78 |
|
79 |
} |
80 |
|
81 |
strcpy(checkPointMsg, "Dump file closed successfully."); |
82 |
MPIcheckPoint(); |
83 |
|
84 |
#endif |
85 |
|
86 |
return; |
87 |
} |
88 |
|
89 |
int DumpReader::getNFrames(void) { |
90 |
#ifdef IS_MPI |
91 |
|
92 |
if (worldRank == 0) { |
93 |
#endif |
94 |
|
95 |
if (!isScanned_) |
96 |
scanFile(); |
97 |
|
98 |
return nframes_; |
99 |
#ifdef IS_MPI |
100 |
} |
101 |
#endif |
102 |
} |
103 |
|
104 |
void DumpReader::scanFile(void) { |
105 |
int i, j; |
106 |
int lineNum = 0; |
107 |
char readBuffer[maxBufferSize]; |
108 |
fpos_t * currPos; |
109 |
|
110 |
#ifdef IS_MPI |
111 |
|
112 |
if (worldRank == 0) { |
113 |
#endif // is_mpi |
114 |
|
115 |
rewind(inFile_); |
116 |
|
117 |
currPos = new fpos_t; |
118 |
fgetpos(inFile_, currPos); |
119 |
fgets(readBuffer, sizeof(readBuffer), inFile_); |
120 |
lineNum++; |
121 |
|
122 |
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 |
|
131 |
while (!feof(inFile_)) { |
132 |
framePos_.push_back(currPos); |
133 |
|
134 |
i = atoi(readBuffer); |
135 |
|
136 |
fgets(readBuffer, sizeof(readBuffer), inFile_); |
137 |
lineNum++; |
138 |
|
139 |
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 |
|
148 |
for(j = 0; j < i; j++) { |
149 |
fgets(readBuffer, sizeof(readBuffer), inFile_); |
150 |
lineNum++; |
151 |
|
152 |
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 |
|
159 |
painCave.isFatal = 1; |
160 |
simError(); |
161 |
} |
162 |
} |
163 |
|
164 |
currPos = new fpos_t; |
165 |
fgetpos(inFile_, currPos); |
166 |
fgets(readBuffer, sizeof(readBuffer), inFile_); |
167 |
lineNum++; |
168 |
} |
169 |
|
170 |
delete currPos; |
171 |
rewind(inFile_); |
172 |
|
173 |
isScanned_ = true; |
174 |
|
175 |
nframes_ = framePos_.size(); |
176 |
#ifdef IS_MPI |
177 |
MPI_Bcast(&nframes_, 1, MPI_INT, 0, MPI_COMM_WORLD); |
178 |
|
179 |
} |
180 |
|
181 |
strcpy(checkPointMsg, "Successfully scanned DumpFile\n"); |
182 |
MPIcheckPoint(); |
183 |
|
184 |
#endif // is_mpi |
185 |
|
186 |
} |
187 |
|
188 |
void DumpReader::readFrame(int whichFrame) { |
189 |
readSet(whichFrame); |
190 |
} |
191 |
|
192 |
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 |
|
198 |
Molecule* mol; |
199 |
StuntDouble* integrableObject; |
200 |
SimInfo::MoleculeIterator mi; |
201 |
Molecule::IntegrableObjectIterator ii; |
202 |
|
203 |
#ifndef IS_MPI |
204 |
|
205 |
fsetpos(inFile_, framePos_[whichFrame]); |
206 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_); |
207 |
|
208 |
if (eof_test == NULL) { |
209 |
sprintf(painCave.errMsg, |
210 |
"DumpReader error: error reading 1st line of \"%s\"\n", |
211 |
filename_.c_str()); |
212 |
painCave.isFatal = 1; |
213 |
simError(); |
214 |
} |
215 |
|
216 |
nTotObjs = atoi(read_buffer); |
217 |
|
218 |
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 |
|
226 |
painCave.isFatal = 1; |
227 |
simError(); |
228 |
} |
229 |
|
230 |
//read the box mat from the comment line |
231 |
|
232 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_); |
233 |
|
234 |
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 |
} |
240 |
|
241 |
parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot()); |
242 |
|
243 |
//parse dump lines |
244 |
|
245 |
for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) { |
246 |
|
247 |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
248 |
integrableObject = mol->nextIntegrableObject(ii)) { |
249 |
|
250 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_); |
251 |
|
252 |
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 |
|
261 |
painCave.isFatal = 1; |
262 |
simError(); |
263 |
} |
264 |
|
265 |
parseDumpLine(read_buffer, integrableObject); |
266 |
|
267 |
} |
268 |
} |
269 |
|
270 |
// MPI Section of code.......... |
271 |
|
272 |
#else //IS_MPI |
273 |
|
274 |
// first thing first, suspend fatalities. |
275 |
int masterNode = 0; |
276 |
int nCurObj; |
277 |
painCave.isEventLoop = 1; |
278 |
|
279 |
int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone |
280 |
int haveError; |
281 |
|
282 |
MPI_Status istatus; |
283 |
int nitems; |
284 |
|
285 |
nTotObjs = info_->getNGlobalIntegrableObjects(); |
286 |
haveError = 0; |
287 |
|
288 |
if (worldRank == masterNode) { |
289 |
fsetpos(inFile_, framePos_[whichFrame]); |
290 |
|
291 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_); |
292 |
|
293 |
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 |
|
300 |
nitems = atoi(read_buffer); |
301 |
|
302 |
// 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 |
|
306 |
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 |
info_->getNGlobalIntegrableObjects()); |
313 |
|
314 |
painCave.isFatal = 1; |
315 |
simError(); |
316 |
} |
317 |
|
318 |
//read the boxMat from the comment line |
319 |
|
320 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_); |
321 |
|
322 |
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 |
} |
328 |
|
329 |
//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 |
|
334 |
MPI_Bcast(read_buffer, maxBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD); |
335 |
parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot()); |
336 |
|
337 |
for(i = 0; i < info_->getNGlobalMolecules(); i++) { |
338 |
int which_node = info_->getMolToProc(i); |
339 |
|
340 |
if (which_node == masterNode) { |
341 |
//molecules belong to master node |
342 |
|
343 |
mol = info_->getMoleculeByGlobalIndex(i); |
344 |
|
345 |
if (mol == NULL) { |
346 |
sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank); |
347 |
painCave.isFatal = 1; |
348 |
simError(); |
349 |
} |
350 |
|
351 |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
352 |
integrableObject = mol->nextIntegrableObject(ii)){ |
353 |
|
354 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_); |
355 |
|
356 |
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 |
|
365 |
painCave.isFatal = 1; |
366 |
simError(); |
367 |
} |
368 |
|
369 |
parseDumpLine(read_buffer, integrableObject); |
370 |
} |
371 |
} else { |
372 |
//molecule belongs to slave nodes |
373 |
|
374 |
MPI_Recv(&nCurObj, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT, |
375 |
MPI_COMM_WORLD, &istatus); |
376 |
|
377 |
for(int j = 0; j < nCurObj; j++) { |
378 |
eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_); |
379 |
|
380 |
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 |
|
389 |
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 |
|
402 |
/**@todo*/ |
403 |
parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot()); |
404 |
|
405 |
for(i = 0; i < info_->getNGlobalMolecules(); i++) { |
406 |
int which_node = info_->getMolToProc(i); |
407 |
|
408 |
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 |
sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank); |
414 |
painCave.isFatal = 1; |
415 |
simError(); |
416 |
} |
417 |
|
418 |
nCurObj = mol->getNIntegrableObjects(); |
419 |
|
420 |
MPI_Send(&nCurObj, 1, MPI_INT, masterNode, TAKE_THIS_TAG_INT, |
421 |
MPI_COMM_WORLD); |
422 |
|
423 |
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 |
|
429 |
parseDumpLine(read_buffer, integrableObject); |
430 |
} |
431 |
|
432 |
} |
433 |
|
434 |
} |
435 |
|
436 |
} |
437 |
|
438 |
#endif |
439 |
|
440 |
} |
441 |
|
442 |
void DumpReader::parseDumpLine(char *line, StuntDouble *integrableObject) { |
443 |
|
444 |
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 |
|
453 |
if (nTokens < 14) { |
454 |
sprintf(painCave.errMsg, |
455 |
"Not enough Tokens.\n"); |
456 |
painCave.isFatal = 1; |
457 |
simError(); |
458 |
} |
459 |
|
460 |
std::string name = tokenizer.nextToken(); |
461 |
|
462 |
if (name != integrableObject->getType()) { |
463 |
|
464 |
} |
465 |
|
466 |
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 |
|
476 |
if (integrableObject->isDirectional()) { |
477 |
|
478 |
q[0] = tokenizer.nextTokenAsFloat(); |
479 |
q[1] = tokenizer.nextTokenAsFloat(); |
480 |
q[2] = tokenizer.nextTokenAsFloat(); |
481 |
q[3] = tokenizer.nextTokenAsFloat(); |
482 |
|
483 |
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 |
} |
492 |
|
493 |
q.normalize(); |
494 |
|
495 |
integrableObject->setQ(q); |
496 |
|
497 |
ji[0] = tokenizer.nextTokenAsFloat(); |
498 |
ji[1] = tokenizer.nextTokenAsFloat(); |
499 |
ji[2] = tokenizer.nextTokenAsFloat(); |
500 |
integrableObject->setJ(ji); |
501 |
} |
502 |
|
503 |
} |
504 |
|
505 |
|
506 |
void DumpReader::parseCommentLine(char* line, Snapshot* s) { |
507 |
double currTime; |
508 |
Mat3x3d hmat; |
509 |
double chi; |
510 |
double integralOfChiDt; |
511 |
Mat3x3d eta; |
512 |
|
513 |
StringTokenizer tokenizer(line); |
514 |
int nTokens; |
515 |
|
516 |
nTokens = tokenizer.countTokens(); |
517 |
|
518 |
//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 |
} |
525 |
|
526 |
//read current time |
527 |
currTime = tokenizer.nextTokenAsFloat(); |
528 |
s->setTime(currTime); |
529 |
|
530 |
//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 |
|
542 |
//read chi and integrablOfChidt, they should apprear in pair |
543 |
if (tokenizer.countTokens() >= 2) { |
544 |
chi = tokenizer.nextTokenAsFloat(); |
545 |
integralOfChiDt = tokenizer.nextTokenAsFloat(); |
546 |
|
547 |
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 |
|
563 |
s->setEta(eta); |
564 |
} |
565 |
|
566 |
|
567 |
} |
568 |
|
569 |
}//end namespace oopse |