ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/io/DumpReader.cpp
Revision: 1957
Committed: Tue Jan 25 17:45:23 2005 UTC (19 years, 5 months ago) by tim
File size: 18330 byte(s)
Log Message:
(1) complete section parser's error message
(2) add GhostTorsion
(3) accumulate inertial tensor from the directional atoms before calculate rigidbody's inertial tensor

File Contents

# User Rev Content
1 gezelter 1930 /*
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 1490 #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 1930 #include "io/DumpReader.hpp"
56     #include "primitives/Molecule.hpp"
57 tim 1492 #include "utils/simError.h"
58 gezelter 1930 #include "utils/MemoryUtils.hpp"
59     #include "utils/StringTokenizer.hpp"
60 gezelter 1490
61     #ifdef IS_MPI
62 gezelter 1930
63 gezelter 1490 #include <mpi.h>
64     #define TAKE_THIS_TAG_CHAR 0
65     #define TAKE_THIS_TAG_INT 1
66 gezelter 1930
67 gezelter 1490 #endif // is_mpi
68    
69    
70 gezelter 1930 namespace oopse {
71 gezelter 1490
72 gezelter 1930 DumpReader::DumpReader(SimInfo* info, const std::string& filename)
73     : info_(info), filename_(filename), isScanned_(false), nframes_(0) {
74 gezelter 1490
75     #ifdef IS_MPI
76    
77 gezelter 1930 if (worldRank == 0) {
78 gezelter 1490 #endif
79    
80 gezelter 1930 inFile_ = fopen(filename_.c_str(), "r");
81    
82     if (inFile_ == NULL) {
83 tim 1957 sprintf(painCave.errMsg, "DumpReader: Cannot open file: %s\n", filename_.c_str());
84 gezelter 1930 painCave.isFatal = 1;
85     simError();
86     }
87    
88 gezelter 1490 #ifdef IS_MPI
89    
90 gezelter 1930 }
91 gezelter 1490
92 gezelter 1930 strcpy(checkPointMsg, "Dump file opened for reading successfully.");
93     MPIcheckPoint();
94    
95 gezelter 1490 #endif
96    
97 gezelter 1930 return;
98 gezelter 1490 }
99    
100 gezelter 1930 DumpReader::~DumpReader() {
101 gezelter 1490
102 gezelter 1930 #ifdef IS_MPI
103 gezelter 1490
104 gezelter 1930 if (worldRank == 0) {
105     #endif
106 gezelter 1490
107 gezelter 1930 int error;
108     error = fclose(inFile_);
109 gezelter 1490
110 gezelter 1930 if (error) {
111     sprintf(painCave.errMsg, "Error closing %s\n", filename_.c_str());
112     painCave.isFatal = 1;
113     simError();
114     }
115    
116     MemoryUtils::deleteVectorOfPointer(framePos_);
117    
118 gezelter 1490 #ifdef IS_MPI
119 gezelter 1930
120 gezelter 1490 }
121    
122 gezelter 1930 strcpy(checkPointMsg, "Dump file closed successfully.");
123     MPIcheckPoint();
124 gezelter 1490
125 gezelter 1930 #endif
126 gezelter 1490
127 gezelter 1930 return;
128 gezelter 1490 }
129    
130 gezelter 1930 int DumpReader::getNFrames(void) {
131 gezelter 1490
132 gezelter 1930 if (!isScanned_)
133     scanFile();
134 gezelter 1490
135 gezelter 1930 return nframes_;
136 gezelter 1490 }
137    
138 gezelter 1930 void DumpReader::scanFile(void) {
139     int i, j;
140     int lineNum = 0;
141     char readBuffer[maxBufferSize];
142     fpos_t * currPos;
143 gezelter 1490
144     #ifdef IS_MPI
145    
146 gezelter 1930 if (worldRank == 0) {
147     #endif // is_mpi
148 gezelter 1490
149 gezelter 1930 rewind(inFile_);
150 gezelter 1490
151 gezelter 1930 currPos = new fpos_t;
152     fgetpos(inFile_, currPos);
153     fgets(readBuffer, sizeof(readBuffer), inFile_);
154     lineNum++;
155 gezelter 1490
156 gezelter 1930 if (feof(inFile_)) {
157     sprintf(painCave.errMsg,
158     "File \"%s\" ended unexpectedly at line %d\n",
159     filename_.c_str(),
160     lineNum);
161     painCave.isFatal = 1;
162     simError();
163     }
164 gezelter 1490
165 gezelter 1930 while (!feof(inFile_)) {
166     framePos_.push_back(currPos);
167 gezelter 1490
168 gezelter 1930 i = atoi(readBuffer);
169 gezelter 1490
170 gezelter 1930 fgets(readBuffer, sizeof(readBuffer), inFile_);
171     lineNum++;
172 gezelter 1490
173 gezelter 1930 if (feof(inFile_)) {
174     sprintf(painCave.errMsg,
175     "File \"%s\" ended unexpectedly at line %d\n",
176     filename_.c_str(),
177     lineNum);
178     painCave.isFatal = 1;
179     simError();
180     }
181 gezelter 1490
182 gezelter 1930 for(j = 0; j < i; j++) {
183     fgets(readBuffer, sizeof(readBuffer), inFile_);
184     lineNum++;
185 gezelter 1490
186 gezelter 1930 if (feof(inFile_)) {
187     sprintf(painCave.errMsg,
188     "File \"%s\" ended unexpectedly at line %d,"
189     " with atom %d\n", filename_.c_str(),
190     lineNum,
191     j);
192 gezelter 1490
193 gezelter 1930 painCave.isFatal = 1;
194     simError();
195     }
196     }
197 gezelter 1490
198 gezelter 1930 currPos = new fpos_t;
199     fgetpos(inFile_, currPos);
200     fgets(readBuffer, sizeof(readBuffer), inFile_);
201     lineNum++;
202     }
203 gezelter 1490
204 gezelter 1930 delete currPos;
205     rewind(inFile_);
206    
207     nframes_ = framePos_.size();
208     #ifdef IS_MPI
209 gezelter 1490 }
210    
211 gezelter 1930 MPI_Bcast(&nframes_, 1, MPI_INT, 0, MPI_COMM_WORLD);
212 gezelter 1490
213 gezelter 1930 strcpy(checkPointMsg, "Successfully scanned DumpFile\n");
214     MPIcheckPoint();
215 gezelter 1490
216 gezelter 1930 #endif // is_mpi
217 gezelter 1490
218 gezelter 1930 isScanned_ = true;
219     }
220 gezelter 1490
221 gezelter 1930 void DumpReader::readFrame(int whichFrame) {
222     readSet(whichFrame);
223     }
224 gezelter 1490
225 gezelter 1930 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 1490
231 gezelter 1930 Molecule* mol;
232     StuntDouble* integrableObject;
233     SimInfo::MoleculeIterator mi;
234     Molecule::IntegrableObjectIterator ii;
235 gezelter 1490
236 gezelter 1930 #ifndef IS_MPI
237 gezelter 1490
238 gezelter 1930 fsetpos(inFile_, framePos_[whichFrame]);
239     eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
240 gezelter 1490
241 gezelter 1930 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 1490 }
248    
249 gezelter 1930 nTotObjs = atoi(read_buffer);
250 gezelter 1490
251 gezelter 1930 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 1957 info_->getNGlobalIntegrableObjects());
258 gezelter 1490
259 gezelter 1930 painCave.isFatal = 1;
260     simError();
261 gezelter 1490 }
262    
263 gezelter 1930 //read the box mat from the comment line
264 gezelter 1490
265 gezelter 1930 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
266 gezelter 1490
267 gezelter 1930 if (eof_test == NULL) {
268     sprintf(painCave.errMsg, "error in reading commment in %s\n",
269     filename_.c_str());
270     painCave.isFatal = 1;
271 gezelter 1490 simError();
272 gezelter 1930 }
273 gezelter 1490
274 gezelter 1930 parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
275 gezelter 1490
276 gezelter 1930 //parse dump lines
277 gezelter 1490
278 gezelter 1930 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
279 gezelter 1490
280 gezelter 1930 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
281     integrableObject = mol->nextIntegrableObject(ii)) {
282 gezelter 1490
283 gezelter 1930 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
284 gezelter 1490
285 gezelter 1930 if (eof_test == NULL) {
286     sprintf(painCave.errMsg,
287     "error in reading file %s\n"
288     "natoms = %d; index = %d\n"
289     "error reading the line from the file.\n",
290     filename_.c_str(),
291     nTotObjs,
292     i);
293 gezelter 1490
294 gezelter 1930 painCave.isFatal = 1;
295     simError();
296     }
297    
298     parseDumpLine(read_buffer, integrableObject);
299    
300     }
301 gezelter 1490 }
302    
303 gezelter 1930 // MPI Section of code..........
304 gezelter 1490
305 gezelter 1930 #else //IS_MPI
306 gezelter 1490
307 gezelter 1930 // first thing first, suspend fatalities.
308     int masterNode = 0;
309     int nCurObj;
310     painCave.isEventLoop = 1;
311 gezelter 1490
312 gezelter 1930 int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
313     int haveError;
314 gezelter 1490
315 gezelter 1930 MPI_Status istatus;
316     int nitems;
317 gezelter 1490
318 gezelter 1930 nTotObjs = info_->getNGlobalIntegrableObjects();
319     haveError = 0;
320 gezelter 1490
321 gezelter 1930 if (worldRank == masterNode) {
322     fsetpos(inFile_, framePos_[whichFrame]);
323 gezelter 1490
324 gezelter 1930 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
325 gezelter 1490
326 gezelter 1930 if (eof_test == NULL) {
327     sprintf(painCave.errMsg, "Error reading 1st line of %s \n ",
328     filename_.c_str());
329     painCave.isFatal = 1;
330     simError();
331 gezelter 1490 }
332    
333 gezelter 1930 nitems = atoi(read_buffer);
334 gezelter 1490
335 gezelter 1930 // 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 1490
339 gezelter 1930 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 1490
347 gezelter 1930 painCave.isFatal = 1;
348     simError();
349     }
350 gezelter 1490
351 gezelter 1930 //read the boxMat from the comment line
352 gezelter 1490
353 gezelter 1930 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
354 gezelter 1490
355 gezelter 1930 if (eof_test == NULL) {
356     sprintf(painCave.errMsg, "error in reading commment in %s\n",
357     filename_.c_str());
358     painCave.isFatal = 1;
359     simError();
360     }
361 gezelter 1490
362 gezelter 1930 //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 1490
367 gezelter 1930 MPI_Bcast(read_buffer, maxBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
368     parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
369 gezelter 1490
370 gezelter 1930 for(i = 0; i < info_->getNGlobalMolecules(); i++) {
371     int which_node = info_->getMolToProc(i);
372 gezelter 1490
373 gezelter 1930 if (which_node == masterNode) {
374     //molecules belong to master node
375 gezelter 1490
376 gezelter 1930 mol = info_->getMoleculeByGlobalIndex(i);
377 gezelter 1490
378 gezelter 1930 if (mol == NULL) {
379     sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
380     painCave.isFatal = 1;
381     simError();
382     }
383 gezelter 1490
384 gezelter 1930 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 1490
389 gezelter 1930 if (eof_test == NULL) {
390     sprintf(painCave.errMsg,
391     "error in reading file %s\n"
392     "natoms = %d; index = %d\n"
393     "error reading the line from the file.\n",
394     filename_.c_str(),
395     nTotObjs,
396     i);
397 gezelter 1490
398 gezelter 1930 painCave.isFatal = 1;
399     simError();
400     }
401 gezelter 1490
402 gezelter 1930 parseDumpLine(read_buffer, integrableObject);
403     }
404     } else {
405     //molecule belongs to slave nodes
406 gezelter 1490
407 gezelter 1930 MPI_Recv(&nCurObj, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT,
408     MPI_COMM_WORLD, &istatus);
409 gezelter 1490
410 gezelter 1930 for(int j = 0; j < nCurObj; j++) {
411     eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
412 gezelter 1490
413 gezelter 1930 if (eof_test == NULL) {
414     sprintf(painCave.errMsg,
415     "error in reading file %s\n"
416     "natoms = %d; index = %d\n"
417     "error reading the line from the file.\n",
418     filename_.c_str(),
419     nTotObjs,
420     i);
421 gezelter 1490
422 gezelter 1930 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 1490
435 gezelter 1930 /**@todo*/
436     parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
437 gezelter 1490
438 gezelter 1930 for(i = 0; i < info_->getNGlobalMolecules(); i++) {
439     int which_node = info_->getMolToProc(i);
440 gezelter 1490
441 gezelter 1930 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     sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
447     painCave.isFatal = 1;
448     simError();
449     }
450    
451     nCurObj = mol->getNIntegrableObjects();
452 gezelter 1490
453 gezelter 1930 MPI_Send(&nCurObj, 1, MPI_INT, masterNode, TAKE_THIS_TAG_INT,
454     MPI_COMM_WORLD);
455 gezelter 1490
456 gezelter 1930 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 1490
462 gezelter 1930 parseDumpLine(read_buffer, integrableObject);
463     }
464    
465     }
466    
467     }
468    
469 gezelter 1490 }
470    
471 gezelter 1930 #endif
472 gezelter 1490
473 gezelter 1930 }
474 gezelter 1490
475 gezelter 1930 void DumpReader::parseDumpLine(char *line, StuntDouble *integrableObject) {
476 gezelter 1490
477 gezelter 1930 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 1490
486 gezelter 1930 if (nTokens < 14) {
487     sprintf(painCave.errMsg,
488     "Not enough Tokens.\n");
489     painCave.isFatal = 1;
490     simError();
491 gezelter 1490 }
492    
493 gezelter 1930 std::string name = tokenizer.nextToken();
494    
495     if (name != integrableObject->getType()) {
496    
497 gezelter 1490 }
498    
499 gezelter 1930 pos[0] = tokenizer.nextTokenAsDouble();
500     pos[1] = tokenizer.nextTokenAsDouble();
501     pos[2] = tokenizer.nextTokenAsDouble();
502     integrableObject->setPos(pos);
503    
504     vel[0] = tokenizer.nextTokenAsDouble();
505     vel[1] = tokenizer.nextTokenAsDouble();
506     vel[2] = tokenizer.nextTokenAsDouble();
507     integrableObject->setVel(vel);
508 gezelter 1490
509 gezelter 1930 if (integrableObject->isDirectional()) {
510    
511     q[0] = tokenizer.nextTokenAsDouble();
512     q[1] = tokenizer.nextTokenAsDouble();
513     q[2] = tokenizer.nextTokenAsDouble();
514     q[3] = tokenizer.nextTokenAsDouble();
515 gezelter 1490
516 gezelter 1930 double qlen = q.length();
517     if (qlen < oopse::epsilon) { //check quaternion is not equal to 0
518    
519     sprintf(painCave.errMsg,
520     "initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2 ~ 0).\n");
521     painCave.isFatal = 1;
522     simError();
523    
524     }
525 gezelter 1490
526 gezelter 1930 q.normalize();
527    
528     integrableObject->setQ(q);
529    
530     ji[0] = tokenizer.nextTokenAsDouble();
531     ji[1] = tokenizer.nextTokenAsDouble();
532     ji[2] = tokenizer.nextTokenAsDouble();
533     integrableObject->setJ(ji);
534 gezelter 1490 }
535    
536     }
537    
538    
539 gezelter 1930 void DumpReader::parseCommentLine(char* line, Snapshot* s) {
540     double currTime;
541     Mat3x3d hmat;
542     double chi;
543     double integralOfChiDt;
544     Mat3x3d eta;
545 gezelter 1490
546 gezelter 1930 StringTokenizer tokenizer(line);
547     int nTokens;
548 gezelter 1490
549 gezelter 1930 nTokens = tokenizer.countTokens();
550 gezelter 1490
551 gezelter 1930 //comment line should at least contain 10 tokens: current time(1 token) and h-matrix(9 tokens)
552     if (nTokens < 10) {
553     sprintf(painCave.errMsg,
554     "Not enough tokens in comment line: %s", line);
555     painCave.isFatal = 1;
556     simError();
557 gezelter 1490 }
558    
559 gezelter 1930 //read current time
560     currTime = tokenizer.nextTokenAsDouble();
561     s->setTime(currTime);
562    
563     //read h-matrix
564     hmat(0, 0) = tokenizer.nextTokenAsDouble();
565     hmat(0, 1) = tokenizer.nextTokenAsDouble();
566     hmat(0, 2) = tokenizer.nextTokenAsDouble();
567     hmat(1, 0) = tokenizer.nextTokenAsDouble();
568     hmat(1, 1) = tokenizer.nextTokenAsDouble();
569     hmat(1, 2) = tokenizer.nextTokenAsDouble();
570     hmat(2, 0) = tokenizer.nextTokenAsDouble();
571     hmat(2, 1) = tokenizer.nextTokenAsDouble();
572     hmat(2, 2) = tokenizer.nextTokenAsDouble();
573     s->setHmat(hmat);
574    
575     //read chi and integrablOfChidt, they should apprear in pair
576     if (tokenizer.countTokens() >= 2) {
577     chi = tokenizer.nextTokenAsDouble();
578     integralOfChiDt = tokenizer.nextTokenAsDouble();
579 gezelter 1490
580 gezelter 1930 s->setChi(chi);
581     s->setIntegralOfChiDt(integralOfChiDt);
582 gezelter 1490 }
583    
584 gezelter 1930 //read eta (eta is 3x3 matrix)
585     if (tokenizer.countTokens() >= 9) {
586     eta(0, 0) = tokenizer.nextTokenAsDouble();
587     eta(0, 1) = tokenizer.nextTokenAsDouble();
588     eta(0, 2) = tokenizer.nextTokenAsDouble();
589     eta(1, 0) = tokenizer.nextTokenAsDouble();
590     eta(1, 1) = tokenizer.nextTokenAsDouble();
591     eta(1, 2) = tokenizer.nextTokenAsDouble();
592     eta(2, 0) = tokenizer.nextTokenAsDouble();
593     eta(2, 1) = tokenizer.nextTokenAsDouble();
594     eta(2, 2) = tokenizer.nextTokenAsDouble();
595    
596     s->setEta(eta);
597 gezelter 1490 }
598 gezelter 1930
599 gezelter 1490
600     }
601    
602 gezelter 1930 }//end namespace oopse