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