ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/io/DumpReader.cpp
Revision: 2300
Committed: Thu Sep 15 19:17:04 2005 UTC (18 years, 10 months ago) by tim
File size: 19198 byte(s)
Log Message:
Fix a bug in DumpReader in case readFrame is called without calling getNFrames

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