ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/io/DumpReader.cpp
Revision: 2204
Committed: Fri Apr 15 22:04:00 2005 UTC (19 years, 3 months ago) by gezelter
File size: 19150 byte(s)
Log Message:
xemacs has been drafted to perform our indentation services

File Contents

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