ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/io/DumpReader.cpp
(Generate patch)

Comparing:
trunk/src/io/DumpReader.cpp (file contents), Revision 273 by tim, Tue Jan 25 17:45:23 2005 UTC vs.
branches/development/src/io/DumpReader.cpp (file contents), Revision 1825 by gezelter, Wed Jan 9 19:27:52 2013 UTC

# Line 1 | Line 1
1 < /*
2 < * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
1 > /*
2 > * Copyright (c) 2009 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
9 > * 1. Redistributions of source code must retain the above copyright
10   *    notice, this list of conditions and the following disclaimer.
11   *
12 < * 3. Redistributions in binary form must reproduce the above copyright
12 > * 2. Redistributions in binary form must reproduce the above copyright
13   *    notice, this list of conditions and the following disclaimer in the
14   *    documentation and/or other materials provided with the
15   *    distribution.
# Line 37 | Line 28
28   * arising out of the use of or inability to use software, even if the
29   * University of Notre Dame has been advised of the possibility of
30   * such damages.
31 + *
32 + * SUPPORT OPEN SCIENCE!  If you use OpenMD or its source code in your
33 + * research, please cite the appropriate papers when you publish your
34 + * work.  Good starting points are:
35 + *                                                                      
36 + * [1]  Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).            
37 + * [2]  Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).          
38 + * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).          
39 + * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010).
40 + * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41   */
42 +  
43 + #define _LARGEFILE_SOURCE64
44 + #define _FILE_OFFSET_BITS 64
45  
46 < #define _LARGEFILE_SOURCE64
47 < #define _FILE_OFFSET_BITS 64
46 > #include <sys/types.h>
47 > #include <sys/stat.h>
48 >
49 > #include <iostream>
50 > #include <math.h>
51 >
52 > #include <stdio.h>
53 > #include <stdlib.h>
54 > #include <string.h>
55 >
56 > #include "io/DumpReader.hpp"
57 > #include "primitives/Molecule.hpp"
58 > #include "utils/simError.h"
59 > #include "utils/MemoryUtils.hpp"
60 > #include "utils/StringTokenizer.hpp"
61 > #include "brains/Thermo.hpp"
62 >
63 > #ifdef IS_MPI
64 > #include <mpi.h>
65 > #endif
66 >
67 >
68 > namespace OpenMD {
69 >  
70 >  DumpReader::DumpReader(SimInfo* info, const std::string& filename)
71 >    : info_(info), filename_(filename), isScanned_(false), nframes_(0), needCOMprops_(false) {
72 >    
73 > #ifdef IS_MPI
74 >    
75 >    if (worldRank == 0) {
76 > #endif
77 >      
78 >      inFile_ = new std::ifstream(filename_.c_str(),  
79 >                                  ifstream::in | ifstream::binary);
80 >      
81 >      if (inFile_->fail()) {
82 >        sprintf(painCave.errMsg,
83 >                "DumpReader: Cannot open file: %s\n",
84 >                filename_.c_str());
85 >        painCave.isFatal = 1;
86 >        simError();
87 >      }
88 >      
89 > #ifdef IS_MPI
90 >      
91 >    }
92 >    
93 >    strcpy(checkPointMsg, "Dump file opened for reading successfully.");
94 >    errorCheckPoint();
95 >    
96 > #endif
97 >    
98 >    return;
99 >  }
100 >  
101 >  DumpReader::~DumpReader() {
102 >    
103 > #ifdef IS_MPI
104 >    
105 >    if (worldRank == 0) {
106 > #endif
107 >      
108 >      delete inFile_;
109 >      
110 > #ifdef IS_MPI
111 >      
112 >    }
113 >    
114 >    strcpy(checkPointMsg, "Dump file closed successfully.");
115 >    errorCheckPoint();
116 >    
117 > #endif
118 >    
119 >    return;
120 >  }
121 >  
122 >  int DumpReader::getNFrames(void) {
123 >    
124 >    if (!isScanned_)
125 >      scanFile();
126 >    
127 >    return nframes_;
128 >  }
129 >  
130 >  void DumpReader::scanFile(void) {
131 >    int lineNo = 0;
132 >    std::streampos prevPos;
133 >    std::streampos  currPos;
134 >    
135 > #ifdef IS_MPI
136 >    
137 >    if (worldRank == 0) {
138 > #endif // is_mpi
139 >      
140 >      currPos = inFile_->tellg();
141 >      prevPos = currPos;
142 >      bool foundOpenSnapshotTag = false;
143 >      bool foundClosedSnapshotTag = false;
144  
145 < #include <sys/types.h>
146 < #include <sys/stat.h>
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:<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 >          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::COMM_WORLD.Bcast(&nframes_, 1, MPI::INT, 0);
206 >    
207 > #endif // is_mpi
208 >    
209 >    isScanned_ = true;
210 >  }
211 >  
212 >  void DumpReader::readFrame(int whichFrame) {
213 >    if (!isScanned_)
214 >      scanFile();
215 >        
216 >    int storageLayout = info_->getSnapshotManager()->getStorageLayout();
217 >    
218 >    if (storageLayout & DataStorage::dslPosition) {
219 >      needPos_ = true;
220 >    } else {
221 >      needPos_ = false;
222 >    }
223 >    
224 >    if (storageLayout & DataStorage::dslVelocity) {
225 >      needVel_ = true;
226 >    } else {
227 >      needVel_ = false;
228 >    }
229 >    
230 >    if (storageLayout & DataStorage::dslAmat ||
231 >        storageLayout & DataStorage::dslDipole ||
232 >        storageLayout & DataStorage::dslQuadrupole) {
233 >      needQuaternion_ = true;
234 >    } else {
235 >      needQuaternion_ = false;
236 >    }
237 >    
238 >    if (storageLayout & DataStorage::dslAngularMomentum) {
239 >      needAngMom_ = true;
240 >    } else {
241 >      needAngMom_ = false;    
242 >    }
243 >    
244 >    readSet(whichFrame);
245  
246 < #include <iostream>
247 < #include <math.h>
246 >    if (needCOMprops_) {
247 >      Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
248 >      Thermo thermo(info_);
249 >      Vector3d com;
250  
251 < #include <stdio.h>
252 < #include <stdlib.h>
253 < #include <string.h>
251 >      if (needPos_ && needVel_) {
252 >        Vector3d comvel;
253 >        Vector3d comw;
254 >        thermo.getComAll(com, comvel);
255 >        comw = thermo.getAngularMomentum();
256 >      } else {
257 >        com = thermo.getCom();
258 >      }                    
259 >    }
260 >  }
261 >  
262 >  void DumpReader::readSet(int whichFrame) {    
263 >    std::string line;
264  
265 < #include "io/DumpReader.hpp"
266 < #include "primitives/Molecule.hpp"
267 < #include "utils/simError.h"
58 < #include "utils/MemoryUtils.hpp"
59 < #include "utils/StringTokenizer.hpp"
265 > #ifndef IS_MPI
266 >    inFile_->clear();  
267 >    inFile_->seekg(framePos_[whichFrame]);
268  
269 < #ifdef IS_MPI
269 >    std::istream& inputStream = *inFile_;    
270  
271 < #include <mpi.h>
272 < #define TAKE_THIS_TAG_CHAR 0
273 < #define TAKE_THIS_TAG_INT 1
271 > #else
272 >    int masterNode = 0;
273 >    std::stringstream sstream;
274 >    if (worldRank == masterNode) {
275 >      std::string sendBuffer;
276  
277 < #endif // is_mpi
277 >      inFile_->clear();  
278 >      inFile_->seekg(framePos_[whichFrame]);
279 >      
280 >      while (inFile_->getline(buffer, bufferSize)) {
281  
282 +        line = buffer;
283 +        sendBuffer += line;
284 +        sendBuffer += '\n';
285 +        if (line.find("</Snapshot>") != std::string::npos) {
286 +          break;
287 +        }        
288 +      }
289  
290 < namespace oopse {
290 >      int sendBufferSize = sendBuffer.size();
291 >      MPI::COMM_WORLD.Bcast(&sendBufferSize, 1, MPI::INT, masterNode);    
292 >      MPI::COMM_WORLD.Bcast((void *)sendBuffer.c_str(), sendBufferSize,
293 >                            MPI::CHAR, masterNode);    
294 >      
295 >      sstream.str(sendBuffer);
296 >    } else {
297 >      int sendBufferSize;
298 >      MPI::COMM_WORLD.Bcast(&sendBufferSize, 1, MPI::INT, masterNode);
299 >      char * recvBuffer = new char[sendBufferSize+1];
300 >      assert(recvBuffer);
301 >      recvBuffer[sendBufferSize] = '\0';
302 >      MPI::COMM_WORLD.Bcast(recvBuffer, sendBufferSize, MPI::CHAR, masterNode);
303 >      sstream.str(recvBuffer);
304 >      delete [] recvBuffer;
305 >    }      
306  
307 < DumpReader::DumpReader(SimInfo* info, const std::string& filename)
73 <                     : info_(info), filename_(filename), isScanned_(false), nframes_(0) {
74 <
75 < #ifdef IS_MPI
76 <
77 <    if (worldRank == 0) {
307 >    std::istream& inputStream = sstream;  
308   #endif
309  
310 <        inFile_ = fopen(filename_.c_str(), "r");
310 >    inputStream.getline(buffer, bufferSize);
311  
312 <        if (inFile_ == NULL) {
313 <            sprintf(painCave.errMsg, "DumpReader: Cannot open file: %s\n", filename_.c_str());
314 <            painCave.isFatal = 1;
315 <            simError();
316 <        }
317 <
318 < #ifdef IS_MPI
312 >    line = buffer;
313 >    if (line.find("<Snapshot>") == std::string::npos) {
314 >      sprintf(painCave.errMsg,
315 >              "DumpReader Error: can not find <Snapshot>\n");
316 >      painCave.isFatal = 1;
317 >      simError();
318 >    }
319 >    
320 >    //read frameData
321 >    readFrameProperties(inputStream);
322  
323 <    }
323 >    //read StuntDoubles
324 >    readStuntDoubles(inputStream);    
325  
326 <    strcpy(checkPointMsg, "Dump file opened for reading successfully.");
327 <    MPIcheckPoint();
326 >    inputStream.getline(buffer, bufferSize);
327 >    line = buffer;
328  
329 < #endif
329 >    if (line.find("<SiteData>") != std::string::npos) {
330 >      //read SiteData
331 >      readSiteData(inputStream);        
332 >    } else {
333 >      if (line.find("</Snapshot>") == std::string::npos) {
334 >        sprintf(painCave.errMsg,
335 >                "DumpReader Error: can not find </Snapshot>\n");
336 >        painCave.isFatal = 1;
337 >        simError();
338 >      }        
339 >    }
340 >  }
341 >  
342 >  void DumpReader::parseDumpLine(const std::string& line) {
343  
344 <    return;
345 < }
344 >      
345 >    StringTokenizer tokenizer(line);
346 >    int nTokens;
347 >    
348 >    nTokens = tokenizer.countTokens();
349 >    
350 >    if (nTokens < 2) {  
351 >      sprintf(painCave.errMsg,
352 >              "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
353 >      painCave.isFatal = 1;
354 >      simError();
355 >    }
356  
357 < DumpReader::~DumpReader() {
357 >    int index = tokenizer.nextTokenAsInt();
358 >
359 >    StuntDouble* sd = info_->getIOIndexToIntegrableObject(index);
360  
361 < #ifdef IS_MPI
361 >    if (sd == NULL) {
362 >      return;
363 >    }
364 >    std::string type = tokenizer.nextToken();
365 >    int size = type.size();
366  
367 <    if (worldRank == 0) {
368 < #endif
369 <
370 <        int error;
371 <        error = fclose(inFile_);
372 <
373 <        if (error) {
374 <            sprintf(painCave.errMsg, "Error closing %s\n", filename_.c_str());
375 <            painCave.isFatal = 1;            
376 <            simError();
367 >    size_t found;
368 >    
369 >    if (needPos_) {
370 >      found = type.find("p");      
371 >      if (found == std::string::npos) {
372 >        sprintf(painCave.errMsg,
373 >                "DumpReader Error: StuntDouble %d has no Position\n"
374 >                "\tField (\"p\") specified.\n%s\n", index,
375 >                line.c_str());  
376 >        painCave.isFatal = 1;
377 >        simError();
378 >      }
379 >    }
380 >    
381 >    if (sd->isDirectional()) {
382 >      if (needQuaternion_) {
383 >        found = type.find("q");      
384 >        if (found == std::string::npos) {
385 >          sprintf(painCave.errMsg,
386 >                  "DumpReader Error: Directional StuntDouble %d has no\n"
387 >                  "\tQuaternion Field (\"q\") specified.\n%s\n", index,
388 >                  line.c_str());  
389 >          painCave.isFatal = 1;
390 >          simError();
391          }
392 <
116 <        MemoryUtils::deleteVectorOfPointer(framePos_);
117 <
118 < #ifdef IS_MPI
119 <
392 >      }      
393      }
394  
395 <    strcpy(checkPointMsg, "Dump file closed successfully.");
396 <    MPIcheckPoint();
395 >    for(int i = 0; i < size; ++i) {
396 >      switch(type[i]) {
397 >        
398 >        case 'p': {
399 >            Vector3d pos;
400 >            pos[0] = tokenizer.nextTokenAsDouble();
401 >            pos[1] = tokenizer.nextTokenAsDouble();
402 >            pos[2] = tokenizer.nextTokenAsDouble();
403 >            if (needPos_) {
404 >              sd->setPos(pos);
405 >            }            
406 >            break;
407 >        }
408 >        case 'v' : {
409 >            Vector3d vel;
410 >            vel[0] = tokenizer.nextTokenAsDouble();
411 >            vel[1] = tokenizer.nextTokenAsDouble();
412 >            vel[2] = tokenizer.nextTokenAsDouble();
413 >            if (needVel_) {
414 >              sd->setVel(vel);
415 >            }
416 >            break;
417 >        }
418  
419 < #endif
419 >        case 'q' : {
420 >           Quat4d q;
421 >           if (sd->isDirectional()) {
422 >              
423 >             q[0] = tokenizer.nextTokenAsDouble();
424 >             q[1] = tokenizer.nextTokenAsDouble();
425 >             q[2] = tokenizer.nextTokenAsDouble();
426 >             q[3] = tokenizer.nextTokenAsDouble();
427 >              
428 >             RealType qlen = q.length();
429 >             if (qlen < OpenMD::epsilon) { //check quaternion is not equal to 0
430 >                
431 >               sprintf(painCave.errMsg,
432 >                       "DumpReader Error: initial quaternion error (q0^2 + q1^2 + q2^2 + q3^2) ~ 0\n");
433 >               painCave.isFatal = 1;
434 >               simError();
435 >                
436 >             }  
437 >              
438 >             q.normalize();
439 >             if (needQuaternion_) {            
440 >               sd->setQ(q);
441 >             }              
442 >           }            
443 >           break;
444 >        }  
445 >        case 'j' : {
446 >          Vector3d ji;
447 >          if (sd->isDirectional()) {
448 >             ji[0] = tokenizer.nextTokenAsDouble();
449 >             ji[1] = tokenizer.nextTokenAsDouble();
450 >             ji[2] = tokenizer.nextTokenAsDouble();
451 >             if (needAngMom_) {
452 >               sd->setJ(ji);
453 >             }
454 >          }
455 >          break;
456 >        }  
457 >        case 'f': {
458  
459 <    return;
460 < }
459 >          Vector3d force;
460 >          force[0] = tokenizer.nextTokenAsDouble();
461 >          force[1] = tokenizer.nextTokenAsDouble();
462 >          force[2] = tokenizer.nextTokenAsDouble();          
463 >          sd->setFrc(force);
464 >          break;
465 >        }
466 >        case 't' : {
467  
468 < int DumpReader::getNFrames(void) {
468 >           Vector3d torque;
469 >           torque[0] = tokenizer.nextTokenAsDouble();
470 >           torque[1] = tokenizer.nextTokenAsDouble();
471 >           torque[2] = tokenizer.nextTokenAsDouble();          
472 >           sd->setTrq(torque);          
473 >           break;
474 >        }
475 >        case 'u' : {
476  
477 <    if (!isScanned_)
478 <        scanFile();
477 >           RealType particlePot;
478 >           particlePot = tokenizer.nextTokenAsDouble();
479 >           sd->setParticlePot(particlePot);          
480 >           break;
481 >        }
482 >        case 'c' : {
483  
484 <    return nframes_;
485 < }
484 >           RealType flucQPos;
485 >           flucQPos = tokenizer.nextTokenAsDouble();
486 >           sd->setFlucQPos(flucQPos);          
487 >           break;
488 >        }
489 >        case 'w' : {
490  
491 < void DumpReader::scanFile(void) {
492 <  int i, j;
493 <  int lineNum = 0;
494 <  char readBuffer[maxBufferSize];
495 <  fpos_t * currPos;
491 >           RealType flucQVel;
492 >           flucQVel = tokenizer.nextTokenAsDouble();
493 >           sd->setFlucQVel(flucQVel);          
494 >           break;
495 >        }
496 >        case 'g' : {
497  
498 < #ifdef IS_MPI
498 >           RealType flucQFrc;
499 >           flucQFrc = tokenizer.nextTokenAsDouble();
500 >           sd->setFlucQFrc(flucQFrc);          
501 >           break;
502 >        }
503 >        case 'e' : {
504  
505 <    if (worldRank == 0) {
506 < #endif // is_mpi
507 <
508 <        rewind(inFile_);
509 <
510 <        currPos = new fpos_t;
152 <        fgetpos(inFile_, currPos);
153 <        fgets(readBuffer, sizeof(readBuffer), inFile_);
154 <        lineNum++;
155 <
156 <        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();
505 >           Vector3d eField;
506 >           eField[0] = tokenizer.nextTokenAsDouble();
507 >           eField[1] = tokenizer.nextTokenAsDouble();
508 >           eField[2] = tokenizer.nextTokenAsDouble();          
509 >           sd->setElectricField(eField);          
510 >           break;
511          }
512 +        default: {
513 +               sprintf(painCave.errMsg,
514 +                       "DumpReader Error: %s is an unrecognized type\n", type.c_str());
515 +               painCave.isFatal = 1;
516 +               simError();
517 +          break;  
518 +        }
519  
520 <        while (!feof(inFile_)) {
521 <            framePos_.push_back(currPos);
520 >      }
521 >    }
522 >    
523 >  }
524 >  
525  
526 <            i = atoi(readBuffer);
526 >  void DumpReader::parseSiteLine(const std::string& line) {
527  
528 <            fgets(readBuffer, sizeof(readBuffer), inFile_);
529 <            lineNum++;
530 <
531 <            if (feof(inFile_)) {
532 <                sprintf(painCave.errMsg,
533 <                        "File \"%s\" ended unexpectedly at line %d\n",
534 <                        filename_.c_str(),
535 <                        lineNum);
536 <                painCave.isFatal = 1;
537 <                simError();
538 <            }
181 <
182 <            for(j = 0; j < i; j++) {
183 <                fgets(readBuffer, sizeof(readBuffer), inFile_);
184 <                lineNum++;
185 <
186 <                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 <
193 <                    painCave.isFatal = 1;
194 <                    simError();
195 <                }
196 <            }
197 <
198 <            currPos = new fpos_t;
199 <            fgetpos(inFile_, currPos);
200 <            fgets(readBuffer, sizeof(readBuffer), inFile_);
201 <            lineNum++;
202 <        }
203 <
204 <        delete currPos;
205 <        rewind(inFile_);
206 <        
207 <        nframes_ = framePos_.size();
208 < #ifdef IS_MPI
209 <    }
528 >    StringTokenizer tokenizer(line);
529 >    int nTokens;
530 >    
531 >    nTokens = tokenizer.countTokens();
532 >    
533 >    if (nTokens < 2) {  
534 >      sprintf(painCave.errMsg,
535 >              "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
536 >      painCave.isFatal = 1;
537 >      simError();
538 >    }
539  
540 <    MPI_Bcast(&nframes_, 1, MPI_INT, 0, MPI_COMM_WORLD);
540 >    /**
541 >     * The first token is the global integrable object index.
542 >     */
543  
544 <    strcpy(checkPointMsg, "Successfully scanned DumpFile\n");
545 <    MPIcheckPoint();
544 >    int index = tokenizer.nextTokenAsInt();
545 >    StuntDouble* sd = info_->getIOIndexToIntegrableObject(index);
546 >    if (sd == NULL) {
547 >      return;
548 >    }
549  
550 < #endif // is_mpi
550 >    /**
551 >     * Test to see if the next token is an integer or not.  If not,
552 >     * we've got data on the integrable object itself.  If there is an
553 >     * integer, we're parsing data for a site on a rigid body.
554 >     */
555  
556 <    isScanned_ = true;
557 < }
556 >    std::string indexTest = tokenizer.peekNextToken();
557 >    std::istringstream i(indexTest);
558 >    int siteIndex;
559 >    if (i >> siteIndex) {
560 >      // chew up this token and parse as an int:
561 >      siteIndex = tokenizer.nextTokenAsInt();
562 >      RigidBody* rb = static_cast<RigidBody*>(sd);
563 >      sd = rb->getAtoms()[siteIndex];
564 >    }
565  
566 < void DumpReader::readFrame(int whichFrame) {
567 <    readSet(whichFrame);
568 < }
566 >    /**
567 >     * The next token contains information on what follows.
568 >     */
569 >    std::string type = tokenizer.nextToken();
570 >    int size = type.size();
571 >    
572 >    for(int i = 0; i < size; ++i) {
573 >      switch(type[i]) {
574 >        
575 >      case 'u' : {
576 >        
577 >        RealType particlePot;
578 >        particlePot = tokenizer.nextTokenAsDouble();
579 >        sd->setParticlePot(particlePot);
580 >        break;
581 >      }
582 >      case 'c' : {
583 >        
584 >        RealType flucQPos;
585 >        flucQPos = tokenizer.nextTokenAsDouble();
586 >        sd->setFlucQPos(flucQPos);
587 >        break;
588 >      }
589 >      case 'w' : {
590 >        
591 >        RealType flucQVel;
592 >        flucQVel = tokenizer.nextTokenAsDouble();
593 >        sd->setFlucQVel(flucQVel);
594 >        break;
595 >      }
596 >      case 'g' : {
597 >        
598 >        RealType flucQFrc;
599 >        flucQFrc = tokenizer.nextTokenAsDouble();
600 >        sd->setFlucQFrc(flucQFrc);
601 >        break;
602 >      }
603 >      case 'e' : {
604 >        
605 >        Vector3d eField;
606 >        eField[0] = tokenizer.nextTokenAsDouble();
607 >        eField[1] = tokenizer.nextTokenAsDouble();
608 >        eField[2] = tokenizer.nextTokenAsDouble();  
609 >        sd->setElectricField(eField);          
610 >        break;
611 >      }
612 >      default: {
613 >        sprintf(painCave.errMsg,
614 >                "DumpReader Error: %s is an unrecognized type\n", type.c_str());
615 >        painCave.isFatal = 1;
616 >        simError();
617 >        break;  
618 >      }
619 >      }
620 >    }    
621 >  }
622 >  
623 >  
624 >  void  DumpReader::readStuntDoubles(std::istream& inputStream) {
625 >    
626 >    inputStream.getline(buffer, bufferSize);
627 >    std::string line(buffer);
628 >    
629 >    if (line.find("<StuntDoubles>") == std::string::npos) {
630 >      sprintf(painCave.errMsg,
631 >              "DumpReader Error: Missing <StuntDoubles>\n");
632 >      painCave.isFatal = 1;
633 >      simError();
634 >    }
635  
636 < void DumpReader::readSet(int whichFrame) {
637 <  int i;
638 <  int nTotObjs;                  // the number of atoms
639 <  char read_buffer[maxBufferSize];  //the line buffer for reading
640 <  char * eof_test;               // ptr to see when we reach the end of the file
636 >    while(inputStream.getline(buffer, bufferSize)) {
637 >      line = buffer;
638 >      
639 >      if(line.find("</StuntDoubles>") != std::string::npos) {
640 >        break;
641 >      }
642  
643 <  Molecule* mol;
644 <  StuntDouble* integrableObject;
645 <  SimInfo::MoleculeIterator mi;
646 <  Molecule::IntegrableObjectIterator ii;
643 >      parseDumpLine(line);
644 >    }
645 >  
646 >  }
647  
648 < #ifndef IS_MPI
648 >  void  DumpReader::readSiteData(std::istream& inputStream) {
649  
650 <    fsetpos(inFile_, framePos_[whichFrame]);
651 <    eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
652 <
653 <    if (eof_test == NULL) {
654 <        sprintf(painCave.errMsg,
655 <                "DumpReader error: error reading 1st line of \"%s\"\n",
244 <                filename_.c_str());
245 <        painCave.isFatal = 1;
246 <        simError();
650 >    inputStream.getline(buffer, bufferSize);
651 >    std::string line(buffer);
652 >    
653 >    if (line.find("<SiteData>") == std::string::npos) {
654 >      // site data isn't required for a simulation, so skip
655 >      return;
656      }
657  
658 <    nTotObjs = atoi(read_buffer);
658 >    while(inputStream.getline(buffer, bufferSize)) {
659 >      line = buffer;
660 >      
661 >      if(line.find("</SiteData>") != std::string::npos) {
662 >        break;
663 >      }
664  
665 <    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 <                info_->getNGlobalIntegrableObjects());
258 <
259 <        painCave.isFatal = 1;
260 <        simError();
665 >      parseSiteLine(line);
666      }
667 +  
668 +  }
669  
670 <    //read the box mat from the comment line
670 >  void DumpReader::readFrameProperties(std::istream& inputStream) {
671  
672 <    eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
672 >    Snapshot* s = info_->getSnapshotManager()->getCurrentSnapshot();
673 >    inputStream.getline(buffer, bufferSize);
674 >    std::string line(buffer);
675  
676 <    if (eof_test == NULL) {
677 <        sprintf(painCave.errMsg, "error in reading commment in %s\n",
678 <                filename_.c_str());
679 <        painCave.isFatal = 1;
680 <        simError();
676 >    if (line.find("<FrameData>") == std::string::npos) {
677 >      sprintf(painCave.errMsg,
678 >              "DumpReader Error: Missing <FrameData>\n");
679 >      painCave.isFatal = 1;
680 >      simError();
681      }
682  
683 <    parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
683 >    while(inputStream.getline(buffer, bufferSize)) {
684 >      line = buffer;
685 >      
686 >      if(line.find("</FrameData>") != std::string::npos) {
687 >        break;
688 >      }
689 >      
690 >      StringTokenizer tokenizer(line, " ;\t\n\r{}:,");
691 >      if (!tokenizer.hasMoreTokens()) {
692 >        sprintf(painCave.errMsg,
693 >                "DumpReader Error: Not enough Tokens.\n%s\n", line.c_str());
694 >        painCave.isFatal = 1;
695 >        simError();      
696 >      }
697  
698 <    //parse dump lines
699 <
700 <    for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
701 <
702 <        for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
703 <            integrableObject = mol->nextIntegrableObject(ii)) {          
704 <
705 <            eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
706 <
707 <            if (eof_test == NULL) {
708 <                sprintf(painCave.errMsg,
709 <                        "error in reading file %s\n"
710 <                            "natoms  = %d; index = %d\n"
711 <                            "error reading the line from the file.\n",
712 <                        filename_.c_str(),
713 <                        nTotObjs,
714 <                        i);
715 <
716 <                painCave.isFatal = 1;
717 <                simError();
718 <            }
719 <
720 <            parseDumpLine(read_buffer, integrableObject);
721 <            
722 <            }
698 >      std::string propertyName = tokenizer.nextToken();
699 >      if (propertyName == "Time") {
700 >        RealType currTime = tokenizer.nextTokenAsDouble();
701 >        s->setTime(currTime);
702 >      } else if (propertyName == "Hmat"){
703 >        Mat3x3d hmat;
704 >        hmat(0, 0) = tokenizer.nextTokenAsDouble();
705 >        hmat(0, 1) = tokenizer.nextTokenAsDouble();
706 >        hmat(0, 2) = tokenizer.nextTokenAsDouble();
707 >        hmat(1, 0) = tokenizer.nextTokenAsDouble();
708 >        hmat(1, 1) = tokenizer.nextTokenAsDouble();
709 >        hmat(1, 2) = tokenizer.nextTokenAsDouble();
710 >        hmat(2, 0) = tokenizer.nextTokenAsDouble();
711 >        hmat(2, 1) = tokenizer.nextTokenAsDouble();
712 >        hmat(2, 2) = tokenizer.nextTokenAsDouble();
713 >        s->setHmat(hmat);      
714 >      } else if (propertyName == "Thermostat") {
715 >        pair<RealType, RealType> thermostat;
716 >        thermostat.first = tokenizer.nextTokenAsDouble();
717 >        thermostat.second = tokenizer.nextTokenAsDouble();
718 >        s->setThermostat(thermostat);
719 >     } else if (propertyName == "Barostat") {
720 >        Mat3x3d eta;
721 >        eta(0, 0) = tokenizer.nextTokenAsDouble();
722 >        eta(0, 1) = tokenizer.nextTokenAsDouble();
723 >        eta(0, 2) = tokenizer.nextTokenAsDouble();
724 >        eta(1, 0) = tokenizer.nextTokenAsDouble();
725 >        eta(1, 1) = tokenizer.nextTokenAsDouble();
726 >        eta(1, 2) = tokenizer.nextTokenAsDouble();
727 >        eta(2, 0) = tokenizer.nextTokenAsDouble();
728 >        eta(2, 1) = tokenizer.nextTokenAsDouble();
729 >        eta(2, 2) = tokenizer.nextTokenAsDouble();
730 >        s->setBarostat(eta);
731 >      } else {
732 >        sprintf(painCave.errMsg,
733 >                "DumpReader Error: %s is an invalid property in <FrameData>\n", propertyName.c_str());
734 >        painCave.isFatal = 0;
735 >        simError();        
736 >      }
737 >      
738      }
739  
740 <    // MPI Section of code..........
740 >  }
741  
742 < #else //IS_MPI
743 <
307 <    // first thing first, suspend fatalities.
308 <    int masterNode = 0;
309 <    int nCurObj;
310 <    painCave.isEventLoop = 1;
311 <
312 <    int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
313 <    int haveError;
314 <
315 <    MPI_Status istatus;
316 <    int nitems;
317 <
318 <    nTotObjs = info_->getNGlobalIntegrableObjects();
319 <    haveError = 0;
320 <
321 <    if (worldRank == masterNode) {
322 <        fsetpos(inFile_, framePos_[whichFrame]);
323 <
324 <        eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
325 <
326 <        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 <        }
332 <
333 <        nitems = atoi(read_buffer);
334 <
335 <        // 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 <
339 <        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 <
347 <            painCave.isFatal = 1;
348 <            simError();
349 <        }
350 <
351 <        //read the boxMat from the comment line
352 <
353 <        eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
354 <
355 <        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 <
362 <        //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 <
367 <        MPI_Bcast(read_buffer, maxBufferSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
368 <        parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
369 <
370 <        for(i = 0; i < info_->getNGlobalMolecules(); i++) {
371 <            int which_node = info_->getMolToProc(i);
372 <
373 <            if (which_node == masterNode) {
374 <                //molecules belong to master node
375 <
376 <                mol = info_->getMoleculeByGlobalIndex(i);
377 <
378 <                if (mol == NULL) {
379 <                    sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
380 <                        painCave.isFatal = 1;
381 <                    simError();
382 <                }
383 <
384 <                for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
385 <                       integrableObject = mol->nextIntegrableObject(ii)){
386 <                        
387 <                    eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
388 <
389 <                    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 <
398 <                        painCave.isFatal = 1;
399 <                        simError();
400 <                    }
401 <
402 <                    parseDumpLine(read_buffer, integrableObject);
403 <                }
404 <            } else {
405 <                //molecule belongs to slave nodes
406 <
407 <                MPI_Recv(&nCurObj, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT,
408 <                         MPI_COMM_WORLD, &istatus);
409 <
410 <                for(int j = 0; j < nCurObj; j++) {
411 <                    eof_test = fgets(read_buffer, sizeof(read_buffer), inFile_);
412 <
413 <                    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 <
422 <                        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 <
435 <        /**@todo*/
436 <        parseCommentLine(read_buffer, info_->getSnapshotManager()->getCurrentSnapshot());
437 <
438 <        for(i = 0; i < info_->getNGlobalMolecules(); i++) {
439 <            int which_node = info_->getMolToProc(i);
440 <
441 <            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 <
453 <                MPI_Send(&nCurObj, 1, MPI_INT, masterNode, TAKE_THIS_TAG_INT,
454 <                         MPI_COMM_WORLD);
455 <
456 <                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 <
462 <                    parseDumpLine(read_buffer, integrableObject);
463 <                }
464 <                      
465 <            }
466 <            
467 <        }
468 <        
469 <    }
470 <
471 < #endif
472 <
473 < }
474 <
475 < void DumpReader::parseDumpLine(char *line, StuntDouble *integrableObject) {
476 <
477 <    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 <
486 <    if (nTokens < 14) {
487 <            sprintf(painCave.errMsg,
488 <                    "Not enough Tokens.\n");
489 <            painCave.isFatal = 1;
490 <            simError();
491 <    }
492 <
493 <    std::string name = tokenizer.nextToken();
494 <
495 <    if (name != integrableObject->getType()) {
496 <        
497 <    }
498 <
499 <    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 <
509 <    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 <
516 <        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 <
526 <        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 <    }
535 <
536 < }
537 <
538 <
539 < void DumpReader::parseCommentLine(char* line, Snapshot* s) {
540 <    double currTime;
541 <    Mat3x3d hmat;
542 <    double chi;
543 <    double integralOfChiDt;
544 <    Mat3x3d eta;
545 <
546 <    StringTokenizer tokenizer(line);
547 <    int nTokens;
548 <
549 <    nTokens = tokenizer.countTokens();
550 <
551 <    //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 <    }
558 <
559 <    //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 <
580 <        s->setChi(chi);
581 <        s->setIntegralOfChiDt(integralOfChiDt);
582 <    }
583 <    
584 <    //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 <    }
598 <
599 <    
600 < }
601 <
602 < }//end namespace oopse
742 >  
743 > }//end namespace OpenMD

Comparing:
trunk/src/io/DumpReader.cpp (property svn:keywords), Revision 273 by tim, Tue Jan 25 17:45:23 2005 UTC vs.
branches/development/src/io/DumpReader.cpp (property svn:keywords), Revision 1825 by gezelter, Wed Jan 9 19:27:52 2013 UTC

# Line 0 | Line 1
1 + Author Id Revision Date

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines