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