70 |
|
namespace oopse { |
71 |
|
|
72 |
|
DumpReader::DumpReader(SimInfo* info, const std::string& filename) |
73 |
< |
: info_(info), filename_(filename), isScanned_(false), nframes_(0) { |
74 |
< |
|
73 |
> |
: info_(info), filename_(filename), isScanned_(false), nframes_(0), needCOMprops_(false) { |
74 |
> |
|
75 |
|
#ifdef IS_MPI |
76 |
< |
|
77 |
< |
if (worldRank == 0) { |
76 |
> |
|
77 |
> |
if (worldRank == 0) { |
78 |
|
#endif |
79 |
< |
|
79 |
> |
|
80 |
|
inFile_ = new std::ifstream(filename_.c_str()); |
81 |
< |
|
82 |
< |
if (inFile_->fail()) { |
83 |
< |
sprintf(painCave.errMsg, "DumpReader: Cannot open file: %s\n", filename_.c_str()); |
84 |
< |
painCave.isFatal = 1; |
85 |
< |
simError(); |
86 |
< |
} |
87 |
< |
|
88 |
< |
#ifdef IS_MPI |
89 |
< |
|
81 |
> |
|
82 |
> |
if (inFile_->fail()) { |
83 |
> |
sprintf(painCave.errMsg, |
84 |
> |
"DumpReader: Cannot open file: %s\n", |
85 |
> |
filename_.c_str()); |
86 |
> |
painCave.isFatal = 1; |
87 |
> |
simError(); |
88 |
|
} |
89 |
< |
|
90 |
< |
strcpy(checkPointMsg, "Dump file opened for reading successfully."); |
91 |
< |
MPIcheckPoint(); |
94 |
< |
|
95 |
< |
#endif |
96 |
< |
|
97 |
< |
return; |
89 |
> |
|
90 |
> |
#ifdef IS_MPI |
91 |
> |
|
92 |
|
} |
93 |
< |
|
93 |
> |
|
94 |
> |
strcpy(checkPointMsg, "Dump file opened for reading successfully."); |
95 |
> |
MPIcheckPoint(); |
96 |
> |
|
97 |
> |
#endif |
98 |
> |
|
99 |
> |
return; |
100 |
> |
} |
101 |
> |
|
102 |
|
DumpReader::~DumpReader() { |
103 |
< |
|
103 |
> |
|
104 |
|
#ifdef IS_MPI |
105 |
< |
|
105 |
> |
|
106 |
|
if (worldRank == 0) { |
107 |
|
#endif |
108 |
< |
|
108 |
> |
|
109 |
|
delete inFile_; |
110 |
< |
|
110 |
> |
|
111 |
|
#ifdef IS_MPI |
112 |
< |
|
112 |
> |
|
113 |
|
} |
114 |
< |
|
114 |
> |
|
115 |
|
strcpy(checkPointMsg, "Dump file closed successfully."); |
116 |
|
MPIcheckPoint(); |
117 |
< |
|
117 |
> |
|
118 |
|
#endif |
119 |
< |
|
119 |
> |
|
120 |
|
return; |
121 |
|
} |
122 |
< |
|
122 |
> |
|
123 |
|
int DumpReader::getNFrames(void) { |
124 |
|
|
125 |
|
if (!isScanned_) |
129 |
|
} |
130 |
|
|
131 |
|
void DumpReader::scanFile(void) { |
132 |
< |
int i, j; |
133 |
< |
int lineNum = 0; |
132 |
< |
char readBuffer[maxBufferSize]; |
132 |
> |
int lineNo = 0; |
133 |
> |
std::streampos prevPos; |
134 |
|
std::streampos currPos; |
135 |
< |
|
135 |
> |
|
136 |
|
#ifdef IS_MPI |
137 |
< |
|
137 |
> |
|
138 |
|
if (worldRank == 0) { |
139 |
|
#endif // is_mpi |
140 |
< |
|
141 |
< |
inFile_->seekg (0, std::ios::beg); |
142 |
< |
|
143 |
< |
|
144 |
< |
currPos = inFile_->tellg(); |
145 |
< |
inFile_->getline(readBuffer, sizeof(readBuffer)); |
146 |
< |
lineNum++; |
147 |
< |
|
148 |
< |
if (inFile_->eof()) { |
149 |
< |
sprintf(painCave.errMsg, |
150 |
< |
"DumpReader Error: File \"%s\" ended unexpectedly at line %d\n", |
151 |
< |
filename_.c_str(), |
151 |
< |
lineNum); |
152 |
< |
painCave.isFatal = 1; |
153 |
< |
simError(); |
154 |
< |
} |
155 |
< |
|
156 |
< |
while (!inFile_->eof()) { |
157 |
< |
framePos_.push_back(currPos); |
158 |
< |
|
159 |
< |
i = atoi(readBuffer); |
160 |
< |
|
161 |
< |
inFile_->getline(readBuffer, sizeof(readBuffer)); |
162 |
< |
lineNum++; |
163 |
< |
|
164 |
< |
if (inFile_->eof()) { |
165 |
< |
sprintf(painCave.errMsg, |
166 |
< |
"DumpReader Error: File \"%s\" ended unexpectedly at line %d\n", |
167 |
< |
filename_.c_str(), |
168 |
< |
lineNum); |
169 |
< |
painCave.isFatal = 1; |
170 |
< |
simError(); |
171 |
< |
} |
172 |
< |
|
173 |
< |
for(j = 0; j < i; j++) { |
174 |
< |
inFile_->getline(readBuffer, sizeof(readBuffer)); |
175 |
< |
lineNum++; |
176 |
< |
|
177 |
< |
if (inFile_->eof()) { |
140 |
> |
|
141 |
> |
currPos = inFile_->tellg(); |
142 |
> |
prevPos = currPos; |
143 |
> |
bool foundOpenSnapshotTag = false; |
144 |
> |
bool foundClosedSnapshotTag = 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 Error: File \"%s\" ended unexpectedly at line %d," |
154 |
< |
" with atom %d\n", filename_.c_str(), |
181 |
< |
lineNum, |
182 |
< |
j); |
183 |
< |
|
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 |
< |
|
172 |
< |
currPos = inFile_->tellg(); |
173 |
< |
inFile_->getline(readBuffer, sizeof(readBuffer)); |
174 |
< |
lineNum++; |
175 |
< |
} |
176 |
< |
|
177 |
< |
inFile_->seekg (0, std::ios::beg); |
178 |
< |
|
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_Bcast(&nframes_, 1, MPI_INT, 0, MPI_COMM_WORLD); |
206 |
< |
|
202 |
< |
strcpy(checkPointMsg, "Successfully scanned DumpFile\n"); |
203 |
< |
MPIcheckPoint(); |
204 |
< |
|
206 |
> |
|
207 |
|
#endif // is_mpi |
208 |
< |
|
208 |
> |
|
209 |
|
isScanned_ = true; |
210 |
|
} |
211 |
|
|
240 |
|
} |
241 |
|
|
242 |
|
readSet(whichFrame); |
243 |
+ |
|
244 |
+ |
if (needCOMprops_) { |
245 |
+ |
Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot(); |
246 |
+ |
Vector3d com; |
247 |
+ |
Vector3d comvel; |
248 |
+ |
Vector3d comw; |
249 |
+ |
if (needPos_ && needVel_){ |
250 |
+ |
info_->getComAll(com, comvel); |
251 |
+ |
comw = info_->getAngularMomentum(); |
252 |
+ |
}else{ |
253 |
+ |
com = info_->getCom(); |
254 |
+ |
comvel = 0.0; |
255 |
+ |
comw = 0.0; |
256 |
+ |
} |
257 |
+ |
s->setCOMprops(com, comvel, comw); |
258 |
+ |
} |
259 |
+ |
|
260 |
|
} |
261 |
|
|
262 |
< |
void DumpReader::readSet(int whichFrame) { |
263 |
< |
int i; |
264 |
< |
int nTotObjs; // the number of atoms |
246 |
< |
char read_buffer[maxBufferSize]; //the line buffer for reading |
247 |
< |
char * eof_test; // ptr to see when we reach the end of the file |
248 |
< |
|
249 |
< |
Molecule* mol; |
250 |
< |
StuntDouble* integrableObject; |
251 |
< |
SimInfo::MoleculeIterator mi; |
252 |
< |
Molecule::IntegrableObjectIterator ii; |
253 |
< |
|
262 |
> |
void DumpReader::readSet(int whichFrame) { |
263 |
> |
std::string line; |
264 |
> |
|
265 |
|
#ifndef IS_MPI |
266 |
+ |
|
267 |
|
inFile_->clear(); |
268 |
|
inFile_->seekg(framePos_[whichFrame]); |
269 |
< |
|
270 |
< |
if (!inFile_->getline(read_buffer, sizeof(read_buffer))) { |
269 |
> |
|
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 |
> |
MPI_Bcast(recvBuffer, sendBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD); |
301 |
> |
sstream.str(recvBuffer); |
302 |
> |
} |
303 |
> |
|
304 |
> |
std::istream& inputStream = sstream; |
305 |
> |
#endif |
306 |
> |
|
307 |
> |
inputStream.getline(buffer, bufferSize); |
308 |
> |
|
309 |
> |
line = buffer; |
310 |
> |
if (line.find("<Snapshot>") == std::string::npos) { |
311 |
|
sprintf(painCave.errMsg, |
312 |
< |
"DumpReader error: error reading 1st line of \"%s\"\n", |
261 |
< |
filename_.c_str()); |
312 |
> |
"DumpReader Error: can not find <Snapshot>\n"); |
313 |
|
painCave.isFatal = 1; |
314 |
|
simError(); |
315 |
|
} |
316 |
< |
|
317 |
< |
nTotObjs = atoi(read_buffer); |
318 |
< |
|
319 |
< |
if (nTotObjs != info_->getNGlobalIntegrableObjects()) { |
316 |
> |
|
317 |
> |
//read frameData |
318 |
> |
readFrameProperties(inputStream); |
319 |
> |
|
320 |
> |
//read StuntDoubles |
321 |
> |
readStuntDoubles(inputStream); |
322 |
> |
|
323 |
> |
inputStream.getline(buffer, bufferSize); |
324 |
> |
line = buffer; |
325 |
> |
if (line.find("</Snapshot>") == std::string::npos) { |
326 |
|
sprintf(painCave.errMsg, |
327 |
< |
"DumpReader error. %s nIntegrable, %d, " |
271 |
< |
"does not match the meta-data file's nIntegrable, %d.\n", |
272 |
< |
filename_.c_str(), |
273 |
< |
nTotObjs, |
274 |
< |
info_->getNGlobalIntegrableObjects()); |
275 |
< |
|
327 |
> |
"DumpReader Error: can not find </Snapshot>\n"); |
328 |
|
painCave.isFatal = 1; |
329 |
|
simError(); |
330 |
< |
} |
331 |
< |
|
280 |
< |
//read the box mat from the comment line |
281 |
< |
|
282 |
< |
|
283 |
< |
if (!inFile_->getline(read_buffer, sizeof(read_buffer))) { |
284 |
< |
sprintf(painCave.errMsg, "DumpReader Error: error in reading commment in %s\n", |
285 |
< |
filename_.c_str()); |
286 |
< |
painCave.isFatal = 1; |
287 |
< |
simError(); |
288 |
< |
} |
289 |
< |
|
290 |
< |
parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot()); |
291 |
< |
|
292 |
< |
//parse dump lines |
293 |
< |
|
294 |
< |
i = 0; |
295 |
< |
for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) { |
296 |
< |
|
297 |
< |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
298 |
< |
integrableObject = mol->nextIntegrableObject(ii)) { |
299 |
< |
|
300 |
< |
|
301 |
< |
|
302 |
< |
if (!inFile_->getline(read_buffer, sizeof(read_buffer))) { |
303 |
< |
sprintf(painCave.errMsg, |
304 |
< |
"DumpReader Error: error in reading file %s\n" |
305 |
< |
"natoms = %d; index = %d\n" |
306 |
< |
"error reading the line from the file.\n", |
307 |
< |
filename_.c_str(), |
308 |
< |
nTotObjs, |
309 |
< |
i); |
310 |
< |
|
311 |
< |
painCave.isFatal = 1; |
312 |
< |
simError(); |
313 |
< |
} |
314 |
< |
|
315 |
< |
parseDumpLine(read_buffer, integrableObject); |
316 |
< |
i++; |
317 |
< |
} |
318 |
< |
} |
319 |
< |
|
320 |
< |
// MPI Section of code.......... |
321 |
< |
|
322 |
< |
#else //IS_MPI |
323 |
< |
|
324 |
< |
// first thing first, suspend fatalities. |
325 |
< |
int masterNode = 0; |
326 |
< |
int nCurObj; |
327 |
< |
painCave.isEventLoop = 1; |
328 |
< |
|
329 |
< |
int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone |
330 |
< |
int haveError; |
331 |
< |
|
332 |
< |
MPI_Status istatus; |
333 |
< |
int nitems; |
334 |
< |
|
335 |
< |
nTotObjs = info_->getNGlobalIntegrableObjects(); |
336 |
< |
haveError = 0; |
337 |
< |
|
338 |
< |
if (worldRank == masterNode) { |
339 |
< |
inFile_->clear(); |
340 |
< |
inFile_->seekg(framePos_[whichFrame]); |
341 |
< |
|
342 |
< |
if (!inFile_->getline(read_buffer, sizeof(read_buffer))) { |
343 |
< |
sprintf(painCave.errMsg, "DumpReader Error: Error reading 1st line of %s \n ", |
344 |
< |
filename_.c_str()); |
345 |
< |
painCave.isFatal = 1; |
346 |
< |
simError(); |
347 |
< |
} |
348 |
< |
|
349 |
< |
nitems = atoi(read_buffer); |
350 |
< |
|
351 |
< |
// Check to see that the number of integrable objects in the |
352 |
< |
// intial configuration file is the same as derived from the |
353 |
< |
// meta-data file. |
354 |
< |
|
355 |
< |
if (nTotObjs != nitems) { |
356 |
< |
sprintf(painCave.errMsg, |
357 |
< |
"DumpReader Error. %s nIntegrable, %d, " |
358 |
< |
"does not match the meta-data file's nIntegrable, %d.\n", |
359 |
< |
filename_.c_str(), |
360 |
< |
nTotObjs, |
361 |
< |
info_->getNGlobalIntegrableObjects()); |
362 |
< |
|
363 |
< |
painCave.isFatal = 1; |
364 |
< |
simError(); |
365 |
< |
} |
366 |
< |
|
367 |
< |
//read the boxMat from the comment line |
368 |
< |
|
369 |
< |
|
370 |
< |
|
371 |
< |
if (!inFile_->getline(read_buffer, sizeof(read_buffer))) { |
372 |
< |
sprintf(painCave.errMsg, "DumpReader Error: error in reading commment in %s\n", |
373 |
< |
filename_.c_str()); |
374 |
< |
painCave.isFatal = 1; |
375 |
< |
simError(); |
376 |
< |
} |
377 |
< |
|
378 |
< |
//Every single processor will parse the comment line by itself |
379 |
< |
//By using this way, we might lose some efficiency, but if we want to add |
380 |
< |
//more parameters into comment line, we only need to modify function |
381 |
< |
//parseCommentLine |
382 |
< |
|
383 |
< |
MPI_Bcast(read_buffer, maxBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD); |
384 |
< |
parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot()); |
385 |
< |
|
386 |
< |
for(i = 0; i < info_->getNGlobalMolecules(); i++) { |
387 |
< |
int which_node = info_->getMolToProc(i); |
388 |
< |
|
389 |
< |
if (which_node == masterNode) { |
390 |
< |
//molecules belong to master node |
391 |
< |
|
392 |
< |
mol = info_->getMoleculeByGlobalIndex(i); |
393 |
< |
|
394 |
< |
if (mol == NULL) { |
395 |
< |
sprintf(painCave.errMsg, "DumpReader Error: Molecule not found on node %d!", worldRank); |
396 |
< |
painCave.isFatal = 1; |
397 |
< |
simError(); |
398 |
< |
} |
399 |
< |
|
400 |
< |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
401 |
< |
integrableObject = mol->nextIntegrableObject(ii)){ |
402 |
< |
|
403 |
< |
|
404 |
< |
|
405 |
< |
if (!inFile_->getline(read_buffer, sizeof(read_buffer))) { |
406 |
< |
sprintf(painCave.errMsg, |
407 |
< |
"DumpReader Error: error in reading file %s\n" |
408 |
< |
"natoms = %d; index = %d\n" |
409 |
< |
"error reading the line from the file.\n", |
410 |
< |
filename_.c_str(), |
411 |
< |
nTotObjs, |
412 |
< |
i); |
413 |
< |
|
414 |
< |
painCave.isFatal = 1; |
415 |
< |
simError(); |
416 |
< |
} |
417 |
< |
|
418 |
< |
parseDumpLine(read_buffer, integrableObject); |
419 |
< |
} |
420 |
< |
} else { |
421 |
< |
//molecule belongs to slave nodes |
422 |
< |
|
423 |
< |
MPI_Recv(&nCurObj, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT, |
424 |
< |
MPI_COMM_WORLD, &istatus); |
425 |
< |
|
426 |
< |
for(int j = 0; j < nCurObj; j++) { |
427 |
< |
|
428 |
< |
|
429 |
< |
if (!inFile_->getline(read_buffer, sizeof(read_buffer))) { |
430 |
< |
sprintf(painCave.errMsg, |
431 |
< |
"DumpReader Error: error in reading file %s\n" |
432 |
< |
"natoms = %d; index = %d\n" |
433 |
< |
"error reading the line from the file.\n", |
434 |
< |
filename_.c_str(), |
435 |
< |
nTotObjs, |
436 |
< |
i); |
437 |
< |
|
438 |
< |
painCave.isFatal = 1; |
439 |
< |
simError(); |
440 |
< |
} |
441 |
< |
|
442 |
< |
MPI_Send(read_buffer, maxBufferSize, MPI_CHAR, which_node, |
443 |
< |
TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD); |
444 |
< |
} |
445 |
< |
} |
446 |
< |
} |
447 |
< |
} else { |
448 |
< |
//actions taken at slave nodes |
449 |
< |
MPI_Bcast(read_buffer, maxBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD); |
450 |
< |
|
451 |
< |
/**@todo*/ |
452 |
< |
parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot()); |
453 |
< |
|
454 |
< |
for(i = 0; i < info_->getNGlobalMolecules(); i++) { |
455 |
< |
int which_node = info_->getMolToProc(i); |
456 |
< |
|
457 |
< |
if (which_node == worldRank) { |
458 |
< |
//molecule with global index i belongs to this processor |
459 |
< |
|
460 |
< |
mol = info_->getMoleculeByGlobalIndex(i); |
461 |
< |
if (mol == NULL) { |
462 |
< |
sprintf(painCave.errMsg, "DumpReader Error: Molecule not found on node %d!", worldRank); |
463 |
< |
painCave.isFatal = 1; |
464 |
< |
simError(); |
465 |
< |
} |
466 |
< |
|
467 |
< |
nCurObj = mol->getNIntegrableObjects(); |
468 |
< |
|
469 |
< |
MPI_Send(&nCurObj, 1, MPI_INT, masterNode, TAKE_THIS_TAG_INT, |
470 |
< |
MPI_COMM_WORLD); |
471 |
< |
|
472 |
< |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
473 |
< |
integrableObject = mol->nextIntegrableObject(ii)){ |
474 |
< |
|
475 |
< |
MPI_Recv(read_buffer, maxBufferSize, MPI_CHAR, masterNode, |
476 |
< |
TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD, &istatus); |
477 |
< |
|
478 |
< |
parseDumpLine(read_buffer, integrableObject); |
479 |
< |
} |
480 |
< |
|
481 |
< |
} |
482 |
< |
|
483 |
< |
} |
484 |
< |
|
485 |
< |
} |
486 |
< |
|
487 |
< |
#endif |
488 |
< |
|
330 |
> |
} |
331 |
> |
|
332 |
|
} |
333 |
|
|
334 |
< |
void DumpReader::parseDumpLine(char *line, StuntDouble *integrableObject) { |
335 |
< |
|
336 |
< |
Vector3d pos; // position place holders |
494 |
< |
Vector3d vel; // velocity placeholders |
495 |
< |
Quat4d q; // the quaternions |
496 |
< |
Vector3d ji; // angular velocity placeholders; |
334 |
> |
void DumpReader::parseDumpLine(const std::string& line) { |
335 |
> |
|
336 |
> |
|
337 |
|
StringTokenizer tokenizer(line); |
338 |
|
int nTokens; |
339 |
|
|
340 |
|
nTokens = tokenizer.countTokens(); |
341 |
|
|
342 |
< |
if (nTokens < 14) { |
342 |
> |
if (nTokens < 2) { |
343 |
|
sprintf(painCave.errMsg, |
344 |
< |
"DumpReader Error: Not enough Tokens.\n%s\n", line); |
344 |
> |
"DumpReader Error: Not enough Tokens.\n%s\n", line.c_str()); |
345 |
|
painCave.isFatal = 1; |
346 |
|
simError(); |
347 |
|
} |
348 |
+ |
|
349 |
+ |
int index = tokenizer.nextTokenAsInt(); |
350 |
+ |
|
351 |
+ |
StuntDouble* integrableObject = info_->getIOIndexToIntegrableObject(index); |
352 |
+ |
|
353 |
+ |
if (integrableObject == NULL) { |
354 |
+ |
return; |
355 |
+ |
} |
356 |
+ |
std::string type = tokenizer.nextToken(); |
357 |
+ |
int size = type.size(); |
358 |
+ |
|
359 |
+ |
for(int i = 0; i < size; ++i) { |
360 |
+ |
switch(type[i]) { |
361 |
+ |
|
362 |
+ |
case 'p': { |
363 |
+ |
Vector3d pos; |
364 |
+ |
pos[0] = tokenizer.nextTokenAsDouble(); |
365 |
+ |
pos[1] = tokenizer.nextTokenAsDouble(); |
366 |
+ |
pos[2] = tokenizer.nextTokenAsDouble(); |
367 |
+ |
if (needPos_) { |
368 |
+ |
integrableObject->setPos(pos); |
369 |
+ |
} |
370 |
+ |
break; |
371 |
+ |
} |
372 |
+ |
case 'v' : { |
373 |
+ |
Vector3d vel; |
374 |
+ |
vel[0] = tokenizer.nextTokenAsDouble(); |
375 |
+ |
vel[1] = tokenizer.nextTokenAsDouble(); |
376 |
+ |
vel[2] = tokenizer.nextTokenAsDouble(); |
377 |
+ |
if (needVel_) { |
378 |
+ |
integrableObject->setVel(vel); |
379 |
+ |
} |
380 |
+ |
break; |
381 |
+ |
} |
382 |
+ |
|
383 |
+ |
case 'q' : { |
384 |
+ |
Quat4d q; |
385 |
+ |
if (integrableObject->isDirectional()) { |
386 |
+ |
|
387 |
+ |
q[0] = tokenizer.nextTokenAsDouble(); |
388 |
+ |
q[1] = tokenizer.nextTokenAsDouble(); |
389 |
+ |
q[2] = tokenizer.nextTokenAsDouble(); |
390 |
+ |
q[3] = tokenizer.nextTokenAsDouble(); |
391 |
+ |
|
392 |
+ |
RealType qlen = q.length(); |
393 |
+ |
if (qlen < oopse::epsilon) { //check quaternion is not equal to 0 |
394 |
+ |
|
395 |
+ |
sprintf(painCave.errMsg, |
396 |
+ |
"DumpReader Error: initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2) ~ 0\n"); |
397 |
+ |
painCave.isFatal = 1; |
398 |
+ |
simError(); |
399 |
+ |
|
400 |
+ |
} |
401 |
+ |
|
402 |
+ |
q.normalize(); |
403 |
+ |
if (needQuaternion_) { |
404 |
+ |
integrableObject->setQ(q); |
405 |
+ |
} |
406 |
+ |
} |
407 |
+ |
break; |
408 |
+ |
} |
409 |
+ |
case 'j' : { |
410 |
+ |
Vector3d ji; |
411 |
+ |
if (integrableObject->isDirectional()) { |
412 |
+ |
ji[0] = tokenizer.nextTokenAsDouble(); |
413 |
+ |
ji[1] = tokenizer.nextTokenAsDouble(); |
414 |
+ |
ji[2] = tokenizer.nextTokenAsDouble(); |
415 |
+ |
if (needAngMom_) { |
416 |
+ |
integrableObject->setJ(ji); |
417 |
+ |
} |
418 |
+ |
} |
419 |
+ |
break; |
420 |
+ |
} |
421 |
+ |
case 'f': { |
422 |
+ |
|
423 |
+ |
Vector3d force; |
424 |
+ |
force[0] = tokenizer.nextTokenAsDouble(); |
425 |
+ |
force[1] = tokenizer.nextTokenAsDouble(); |
426 |
+ |
force[2] = tokenizer.nextTokenAsDouble(); |
427 |
+ |
integrableObject->setFrc(force); |
428 |
+ |
break; |
429 |
+ |
} |
430 |
+ |
case 't' : { |
431 |
+ |
|
432 |
+ |
Vector3d torque; |
433 |
+ |
torque[0] = tokenizer.nextTokenAsDouble(); |
434 |
+ |
torque[1] = tokenizer.nextTokenAsDouble(); |
435 |
+ |
torque[2] = tokenizer.nextTokenAsDouble(); |
436 |
+ |
integrableObject->setTrq(torque); |
437 |
+ |
break; |
438 |
+ |
} |
439 |
+ |
default: { |
440 |
+ |
sprintf(painCave.errMsg, |
441 |
+ |
"DumpReader Error: %s is an unrecognized type\n", type.c_str()); |
442 |
+ |
painCave.isFatal = 1; |
443 |
+ |
simError(); |
444 |
+ |
break; |
445 |
+ |
} |
446 |
+ |
|
447 |
+ |
} |
448 |
+ |
} |
449 |
|
|
450 |
< |
std::string name = tokenizer.nextToken(); |
451 |
< |
|
452 |
< |
if (name != integrableObject->getType()) { |
453 |
< |
|
450 |
> |
} |
451 |
> |
|
452 |
> |
|
453 |
> |
void DumpReader::readStuntDoubles(std::istream& inputStream) { |
454 |
> |
|
455 |
> |
inputStream.getline(buffer, bufferSize); |
456 |
> |
std::string line(buffer); |
457 |
> |
|
458 |
> |
if (line.find("<StuntDoubles>") == std::string::npos) { |
459 |
|
sprintf(painCave.errMsg, |
460 |
< |
"DumpReader Error: Atom type [%s] in %s does not match Atom Type [%s] in .md file.\n", |
515 |
< |
name.c_str(), filename_.c_str(), integrableObject->getType().c_str()); |
460 |
> |
"DumpReader Error: Missing <StuntDoubles>\n"); |
461 |
|
painCave.isFatal = 1; |
462 |
< |
simError(); |
463 |
< |
} |
464 |
< |
|
465 |
< |
pos[0] = tokenizer.nextTokenAsDouble(); |
466 |
< |
pos[1] = tokenizer.nextTokenAsDouble(); |
467 |
< |
pos[2] = tokenizer.nextTokenAsDouble(); |
468 |
< |
if (needPos_) { |
469 |
< |
integrableObject->setPos(pos); |
470 |
< |
} |
471 |
< |
|
472 |
< |
vel[0] = tokenizer.nextTokenAsDouble(); |
473 |
< |
vel[1] = tokenizer.nextTokenAsDouble(); |
474 |
< |
vel[2] = tokenizer.nextTokenAsDouble(); |
475 |
< |
if (needVel_) { |
476 |
< |
integrableObject->setVel(vel); |
477 |
< |
} |
478 |
< |
|
479 |
< |
if (integrableObject->isDirectional()) { |
480 |
< |
|
481 |
< |
q[0] = tokenizer.nextTokenAsDouble(); |
482 |
< |
q[1] = tokenizer.nextTokenAsDouble(); |
483 |
< |
q[2] = tokenizer.nextTokenAsDouble(); |
484 |
< |
q[3] = tokenizer.nextTokenAsDouble(); |
485 |
< |
|
486 |
< |
RealType qlen = q.length(); |
487 |
< |
if (qlen < oopse::epsilon) { //check quaternion is not equal to 0 |
488 |
< |
|
462 |
> |
simError(); |
463 |
> |
} |
464 |
> |
|
465 |
> |
while(inputStream.getline(buffer, bufferSize)) { |
466 |
> |
line = buffer; |
467 |
> |
|
468 |
> |
if(line.find("</StuntDoubles>") != std::string::npos) { |
469 |
> |
break; |
470 |
> |
} |
471 |
> |
|
472 |
> |
parseDumpLine(line); |
473 |
> |
} |
474 |
> |
|
475 |
> |
} |
476 |
> |
|
477 |
> |
void DumpReader::readFrameProperties(std::istream& inputStream) { |
478 |
> |
|
479 |
> |
Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot(); |
480 |
> |
inputStream.getline(buffer, bufferSize); |
481 |
> |
std::string line(buffer); |
482 |
> |
|
483 |
> |
if (line.find("<FrameData>") == std::string::npos) { |
484 |
> |
sprintf(painCave.errMsg, |
485 |
> |
"DumpReader Error: Missing <FrameData>\n"); |
486 |
> |
painCave.isFatal = 1; |
487 |
> |
simError(); |
488 |
> |
} |
489 |
> |
|
490 |
> |
while(inputStream.getline(buffer, bufferSize)) { |
491 |
> |
line = buffer; |
492 |
> |
|
493 |
> |
if(line.find("</FrameData>") != std::string::npos) { |
494 |
> |
break; |
495 |
> |
} |
496 |
> |
|
497 |
> |
StringTokenizer tokenizer(line, " ;\t\n\r{}:,"); |
498 |
> |
if (!tokenizer.hasMoreTokens()) { |
499 |
|
sprintf(painCave.errMsg, |
500 |
< |
"DumpReader Error: initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2) ~ 0\n"); |
500 |
> |
"DumpReader Error: Not enough Tokens.\n%s\n", line.c_str()); |
501 |
|
painCave.isFatal = 1; |
502 |
< |
simError(); |
503 |
< |
|
504 |
< |
} |
505 |
< |
|
506 |
< |
q.normalize(); |
507 |
< |
if (needQuaternion_) { |
508 |
< |
integrableObject->setQ(q); |
509 |
< |
} |
510 |
< |
|
511 |
< |
ji[0] = tokenizer.nextTokenAsDouble(); |
512 |
< |
ji[1] = tokenizer.nextTokenAsDouble(); |
513 |
< |
ji[2] = tokenizer.nextTokenAsDouble(); |
514 |
< |
if (needAngMom_) { |
515 |
< |
integrableObject->setJ(ji); |
516 |
< |
} |
517 |
< |
} |
518 |
< |
|
519 |
< |
} |
502 |
> |
simError(); |
503 |
> |
} |
504 |
> |
|
505 |
> |
std::string propertyName = tokenizer.nextToken(); |
506 |
> |
if (propertyName == "Time") { |
507 |
> |
RealType currTime = tokenizer.nextTokenAsDouble(); |
508 |
> |
s->setTime(currTime); |
509 |
> |
} else if (propertyName == "Hmat"){ |
510 |
> |
Mat3x3d hmat; |
511 |
> |
hmat(0, 0) = tokenizer.nextTokenAsDouble(); |
512 |
> |
hmat(0, 1) = tokenizer.nextTokenAsDouble(); |
513 |
> |
hmat(0, 2) = tokenizer.nextTokenAsDouble(); |
514 |
> |
hmat(1, 0) = tokenizer.nextTokenAsDouble(); |
515 |
> |
hmat(1, 1) = tokenizer.nextTokenAsDouble(); |
516 |
> |
hmat(1, 2) = tokenizer.nextTokenAsDouble(); |
517 |
> |
hmat(2, 0) = tokenizer.nextTokenAsDouble(); |
518 |
> |
hmat(2, 1) = tokenizer.nextTokenAsDouble(); |
519 |
> |
hmat(2, 2) = tokenizer.nextTokenAsDouble(); |
520 |
> |
s->setHmat(hmat); |
521 |
> |
} else if (propertyName == "Thermostat") { |
522 |
> |
RealType chi = tokenizer.nextTokenAsDouble(); |
523 |
> |
RealType integralOfChiDt = tokenizer.nextTokenAsDouble(); |
524 |
> |
s->setChi(chi); |
525 |
> |
s->setIntegralOfChiDt(integralOfChiDt); |
526 |
> |
} else if (propertyName == "Barostat") { |
527 |
> |
Mat3x3d eta; |
528 |
> |
eta(0, 0) = tokenizer.nextTokenAsDouble(); |
529 |
> |
eta(0, 1) = tokenizer.nextTokenAsDouble(); |
530 |
> |
eta(0, 2) = tokenizer.nextTokenAsDouble(); |
531 |
> |
eta(1, 0) = tokenizer.nextTokenAsDouble(); |
532 |
> |
eta(1, 1) = tokenizer.nextTokenAsDouble(); |
533 |
> |
eta(1, 2) = tokenizer.nextTokenAsDouble(); |
534 |
> |
eta(2, 0) = tokenizer.nextTokenAsDouble(); |
535 |
> |
eta(2, 1) = tokenizer.nextTokenAsDouble(); |
536 |
> |
eta(2, 2) = tokenizer.nextTokenAsDouble(); |
537 |
> |
s->setEta(eta); |
538 |
> |
} else { |
539 |
> |
sprintf(painCave.errMsg, |
540 |
> |
"DumpReader Error: %s is an invalid property in <FrameData>\n", propertyName.c_str()); |
541 |
> |
painCave.isFatal = 0; |
542 |
> |
simError(); |
543 |
> |
} |
544 |
> |
|
545 |
> |
} |
546 |
> |
|
547 |
> |
} |
548 |
> |
|
549 |
|
|
566 |
– |
|
567 |
– |
void DumpReader::parseCommentLine(char* line, Snapshot* s) { |
568 |
– |
RealType currTime; |
569 |
– |
Mat3x3d hmat; |
570 |
– |
RealType chi; |
571 |
– |
RealType integralOfChiDt; |
572 |
– |
Mat3x3d eta; |
573 |
– |
|
574 |
– |
StringTokenizer tokenizer(line); |
575 |
– |
int nTokens; |
576 |
– |
|
577 |
– |
nTokens = tokenizer.countTokens(); |
578 |
– |
|
579 |
– |
//comment line should at least contain 10 tokens: current time(1 token) and h-matrix(9 tokens) |
580 |
– |
if (nTokens < 10) { |
581 |
– |
sprintf(painCave.errMsg, |
582 |
– |
"DumpReader Error: Not enough tokens in comment line: %s", line); |
583 |
– |
painCave.isFatal = 1; |
584 |
– |
simError(); |
585 |
– |
} |
586 |
– |
|
587 |
– |
//read current time |
588 |
– |
currTime = tokenizer.nextTokenAsDouble(); |
589 |
– |
s->setTime(currTime); |
590 |
– |
|
591 |
– |
//read h-matrix |
592 |
– |
hmat(0, 0) = tokenizer.nextTokenAsDouble(); |
593 |
– |
hmat(0, 1) = tokenizer.nextTokenAsDouble(); |
594 |
– |
hmat(0, 2) = tokenizer.nextTokenAsDouble(); |
595 |
– |
hmat(1, 0) = tokenizer.nextTokenAsDouble(); |
596 |
– |
hmat(1, 1) = tokenizer.nextTokenAsDouble(); |
597 |
– |
hmat(1, 2) = tokenizer.nextTokenAsDouble(); |
598 |
– |
hmat(2, 0) = tokenizer.nextTokenAsDouble(); |
599 |
– |
hmat(2, 1) = tokenizer.nextTokenAsDouble(); |
600 |
– |
hmat(2, 2) = tokenizer.nextTokenAsDouble(); |
601 |
– |
s->setHmat(hmat); |
602 |
– |
|
603 |
– |
//read chi and integralOfChidt, they should apprear in pair |
604 |
– |
if (tokenizer.countTokens() >= 2) { |
605 |
– |
chi = tokenizer.nextTokenAsDouble(); |
606 |
– |
integralOfChiDt = tokenizer.nextTokenAsDouble(); |
607 |
– |
|
608 |
– |
s->setChi(chi); |
609 |
– |
s->setIntegralOfChiDt(integralOfChiDt); |
610 |
– |
} |
611 |
– |
|
612 |
– |
//read eta (eta is 3x3 matrix) |
613 |
– |
if (tokenizer.countTokens() >= 9) { |
614 |
– |
eta(0, 0) = tokenizer.nextTokenAsDouble(); |
615 |
– |
eta(0, 1) = tokenizer.nextTokenAsDouble(); |
616 |
– |
eta(0, 2) = tokenizer.nextTokenAsDouble(); |
617 |
– |
eta(1, 0) = tokenizer.nextTokenAsDouble(); |
618 |
– |
eta(1, 1) = tokenizer.nextTokenAsDouble(); |
619 |
– |
eta(1, 2) = tokenizer.nextTokenAsDouble(); |
620 |
– |
eta(2, 0) = tokenizer.nextTokenAsDouble(); |
621 |
– |
eta(2, 1) = tokenizer.nextTokenAsDouble(); |
622 |
– |
eta(2, 2) = tokenizer.nextTokenAsDouble(); |
623 |
– |
|
624 |
– |
s->setEta(eta); |
625 |
– |
} |
626 |
– |
|
627 |
– |
|
628 |
– |
} |
629 |
– |
|
550 |
|
}//end namespace oopse |