ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpReader.cpp
Revision: 1252
Committed: Mon Jun 7 14:26:33 2004 UTC (20 years, 1 month ago) by gezelter
File size: 17927 byte(s)
Log Message:
Fixes from gcc -Wall

File Contents

# User Rev Content
1 tim 1078 #define _LARGEFILE_SOURCE64
2 mmeineke 643 #define _FILE_OFFSET_BITS 64
3    
4 mmeineke 804 #include <sys/types.h>
5     #include <sys/stat.h>
6    
7 gezelter 637 #include <iostream>
8 gezelter 829 #include <math.h>
9 gezelter 637
10     #include <stdio.h>
11     #include <stdlib.h>
12     #include <string.h>
13    
14 mmeineke 804
15 gezelter 637 #include "ReadWrite.hpp"
16     #include "simError.h"
17    
18     #ifdef IS_MPI
19     #include <mpi.h>
20     #include "mpiSimulation.hpp"
21     #define TAKE_THIS_TAG_CHAR 0
22     #define TAKE_THIS_TAG_INT 1
23 mmeineke 647 #endif // is_mpi
24 gezelter 637
25 mmeineke 647
26 tim 1074 DumpReader :: DumpReader(const char *in_name ){
27 mmeineke 647
28     isScanned = false;
29    
30 gezelter 637 #ifdef IS_MPI
31     if (worldRank == 0) {
32     #endif
33    
34 mmeineke 647 inFile = fopen(in_name, "r");
35     if(inFile == NULL){
36 gezelter 637 sprintf(painCave.errMsg,
37     "Cannot open file: %s\n", in_name);
38     painCave.isFatal = 1;
39     simError();
40     }
41    
42 tim 1118 inFileName = in_name;
43 gezelter 637 #ifdef IS_MPI
44     }
45     strcpy( checkPointMsg, "Dump file opened for reading successfully." );
46     MPIcheckPoint();
47     #endif
48     return;
49     }
50    
51     DumpReader :: ~DumpReader( ){
52     #ifdef IS_MPI
53     if (worldRank == 0) {
54     #endif
55 tim 1118 vector<fpos_t*>::iterator i;
56    
57 gezelter 637 int error;
58 mmeineke 647 error = fclose( inFile );
59 gezelter 637 if( error ){
60     sprintf( painCave.errMsg,
61 tim 1118 "Error closing %s\n", inFileName.c_str());
62 gezelter 637 simError();
63     }
64 tim 1118
65     for(i = framePos.begin(); i != framePos.end(); ++i)
66     delete *i;
67     framePos.clear();
68    
69 gezelter 637 #ifdef IS_MPI
70     }
71     strcpy( checkPointMsg, "Dump file closed successfully." );
72     MPIcheckPoint();
73     #endif
74    
75     return;
76     }
77    
78 mmeineke 647 int DumpReader::getNframes( void ){
79 gezelter 637
80 tim 1118 if( !isScanned )
81     scanFile();
82     return framePos.size();
83 mmeineke 647 }
84 mmeineke 643
85 mmeineke 647 void DumpReader::scanFile( void ){
86    
87 gezelter 1252 int i, j;
88 mmeineke 647 int lineNum = 0;
89     char readBuffer[2000];
90     fpos_t *currPos;
91    
92     #ifdef IS_MPI
93     if( worldRank == 0 ){
94     #endif // is_mpi
95    
96     rewind( inFile );
97    
98     currPos = new fpos_t;
99     fgetpos( inFile, currPos );
100 mmeineke 804 fgets( readBuffer, sizeof( readBuffer ), inFile );
101 mmeineke 647 lineNum++;
102 mmeineke 804 if( feof( inFile ) ){
103 mmeineke 647 sprintf( painCave.errMsg,
104     "File \"%s\" ended unexpectedly at line %d\n",
105 tim 1118 inFileName.c_str(),
106 mmeineke 647 lineNum );
107     painCave.isFatal = 1;
108     simError();
109     }
110 tim 1118
111 mmeineke 804 while( !feof( inFile ) ){
112 tim 1119
113     framePos.push_back(currPos);
114 mmeineke 804
115 mmeineke 647 i = atoi(readBuffer);
116    
117 mmeineke 804 fgets( readBuffer, sizeof( readBuffer ), inFile );
118 mmeineke 647 lineNum++;
119 mmeineke 804 if( feof( inFile ) ){
120 mmeineke 647 sprintf( painCave.errMsg,
121     "File \"%s\" ended unexpectedly at line %d\n",
122 tim 1118 inFileName.c_str(),
123 mmeineke 647 lineNum );
124     painCave.isFatal = 1;
125     simError();
126     }
127 tim 1118
128 mmeineke 647 for(j=0; j<i; j++){
129    
130 mmeineke 804 fgets( readBuffer, sizeof( readBuffer ), inFile );
131 mmeineke 647 lineNum++;
132 mmeineke 804 if( feof( inFile ) ){
133 mmeineke 647 sprintf( painCave.errMsg,
134     "File \"%s\" ended unexpectedly at line %d,"
135     " with atom %d\n",
136 tim 1118 inFileName.c_str(),
137 mmeineke 647 lineNum,
138     j );
139     painCave.isFatal = 1;
140     simError();
141     }
142    
143     }
144    
145     currPos = new fpos_t;
146     fgetpos( inFile, currPos );
147 mmeineke 804 fgets( readBuffer, sizeof( readBuffer ), inFile );
148 mmeineke 647 lineNum++;
149     }
150    
151     delete currPos;
152     rewind( inFile );
153 tim 1118
154 mmeineke 804 isScanned = true;
155    
156 mmeineke 647 #ifdef IS_MPI
157     }
158     strcpy( checkPointMsg, "Successfully scanned DumpFile\n" );
159     MPIcheckPoint();
160     #endif // is_mpi
161     }
162    
163 mmeineke 656 void DumpReader :: readFrame( SimInfo* the_simnfo, int whichFrame){
164 mmeineke 647
165 mmeineke 656 simnfo = the_simnfo;
166 mmeineke 647
167 mmeineke 656 this->readSet( whichFrame );
168     }
169 gezelter 637
170 mmeineke 656
171    
172     void DumpReader :: readSet( int whichFrame ){
173    
174 gezelter 1252 int i;
175     unsigned int j;
176 gezelter 637
177     #ifdef IS_MPI
178 tim 1118 int done, which_node, which_atom; // loop counter
179     #endif //is_mpi
180 gezelter 637
181 tim 1118 const int BUFFERSIZE = 2000; // size of the read buffer
182     int nTotObjs; // the number of atoms
183     char read_buffer[BUFFERSIZE]; //the line buffer for reading
184    
185     char *eof_test; // ptr to see when we reach the end of the file
186 gezelter 637 char *parseErr;
187 mmeineke 804
188 tim 1118 vector<StuntDouble*> integrableObjects;
189    
190    
191 mmeineke 656 #ifndef IS_MPI
192 gezelter 637
193 tim 1118 fsetpos(inFile, framePos[whichFrame]);
194 mmeineke 647 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
195 gezelter 637 if( eof_test == NULL ){
196     sprintf( painCave.errMsg,
197 tim 1119 "DumpReader error: error reading 1st line of \"%s\"\n",
198 tim 1118 inFileName.c_str() );
199 gezelter 637 painCave.isFatal = 1;
200     simError();
201     }
202    
203 tim 1118 nTotObjs = atoi( read_buffer );
204 gezelter 637
205 tim 1118 if( nTotObjs != simnfo->getTotIntegrableObjects() ){
206 gezelter 637 sprintf( painCave.errMsg,
207 tim 1119 "DumpReader error. %s n_atoms, %d, "
208 gezelter 637 "does not match the BASS file's n_atoms, %d.\n",
209 tim 1118 inFileName.c_str(), nTotObjs, simnfo->getTotIntegrableObjects());
210 gezelter 637 painCave.isFatal = 1;
211     simError();
212     }
213 tim 1118
214     //read the box mat from the comment line
215    
216 mmeineke 647 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
217 gezelter 637 if(eof_test == NULL){
218     sprintf( painCave.errMsg,
219 tim 1118 "error in reading commment in %s\n", inFileName.c_str());
220 gezelter 637 painCave.isFatal = 1;
221     simError();
222     }
223    
224 tim 1118 parseErr = parseCommentLine( read_buffer, simnfo);
225 gezelter 637 if( parseErr != NULL ){
226     strcpy( painCave.errMsg, parseErr );
227     painCave.isFatal = 1;
228     simError();
229     }
230    
231 tim 1118 //parse dump lines
232 gezelter 637
233 tim 1118 for( i=0; i < simnfo->n_mol; i++){
234 gezelter 637
235 tim 1118 integrableObjects = (simnfo->molecules[i]).getIntegrableObjects();
236 gezelter 637
237 tim 1118 for(j = 0; j < integrableObjects.size(); j++){
238    
239     eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
240     if(eof_test == NULL){
241     sprintf(painCave.errMsg,
242     "error in reading file %s\n"
243     "natoms = %d; index = %d\n"
244     "error reading the line from the file.\n",
245     inFileName.c_str(), nTotObjs, i );
246     painCave.isFatal = 1;
247     simError();
248     }
249    
250     parseErr = parseDumpLine( read_buffer, integrableObjects[j]);
251     if( parseErr != NULL ){
252     strcpy( painCave.errMsg, parseErr );
253     painCave.isFatal = 1;
254     simError();
255     }
256 gezelter 637 }
257     }
258    
259     // MPI Section of code..........
260     #else //IS_MPI
261    
262     // first thing first, suspend fatalities.
263     painCave.isEventLoop = 1;
264    
265     int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
266     int haveError;
267 tim 1118
268 gezelter 637 MPI_Status istatus;
269 tim 1118 int *MolToProcMap = mpiSim->getMolToProcMap();
270     int localIndex;
271     int nCurObj;
272 tim 1129 int nitems;
273 gezelter 637
274 tim 1129 nTotObjs = simnfo->getTotIntegrableObjects();
275 gezelter 637 haveError = 0;
276     if (worldRank == 0) {
277 tim 1118 fsetpos(inFile, framePos[whichFrame]);
278    
279 mmeineke 647 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
280 gezelter 637 if( eof_test == NULL ){
281     sprintf( painCave.errMsg,
282 tim 1118 "Error reading 1st line of %s \n ",inFileName.c_str());
283 gezelter 637 haveError = 1;
284     simError();
285     }
286 tim 1118
287 tim 1129 nitems = atoi( read_buffer );
288 tim 1118
289     // Check to see that the number of integrable objects in the intial configuration file is the
290 gezelter 637 // same as declared in simBass.
291 tim 1118
292 tim 1129 if( nTotObjs != nitems){
293 gezelter 637 sprintf( painCave.errMsg,
294 tim 1119 "DumpReadererror. %s n_atoms, %d, "
295 gezelter 637 "does not match the BASS file's n_atoms, %d.\n",
296 tim 1118 inFileName.c_str(), nTotObjs, simnfo->getTotIntegrableObjects());
297 gezelter 637 haveError= 1;
298     simError();
299     }
300 tim 1118
301     //read the boxMat from the comment line
302    
303 mmeineke 647 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
304 gezelter 637 if(eof_test == NULL){
305     sprintf( painCave.errMsg,
306 tim 1118 "error in reading commment in %s\n", inFileName.c_str());
307 gezelter 637 haveError = 1;
308     simError();
309     }
310 tim 1118
311     //Every single processor will parse the comment line by itself
312     //By using this way, we might lose some efficiency, but if we want to add
313     //more parameters into comment line, we only need to modify function
314     //parseCommentLine
315    
316     MPI_Bcast(read_buffer, BUFFERSIZE, MPI_CHAR, 0, MPI_COMM_WORLD);
317    
318     parseErr = parseCommentLine( read_buffer, simnfo);
319    
320 gezelter 637 if( parseErr != NULL ){
321     strcpy( painCave.errMsg, parseErr );
322     haveError = 1;
323     simError();
324     }
325    
326 gezelter 1203 for (i=0 ; i < mpiSim->getNMolGlobal(); i++) {
327 tim 1118 which_node = MolToProcMap[i];
328     if(which_node == 0){
329     //molecules belong to master node
330    
331     localIndex = mpiSim->getGlobalToLocalMol(i);
332    
333     if(localIndex == -1) {
334     strcpy(painCave.errMsg, "Molecule not found on node 0!");
335     haveError = 1;
336     simError();
337 gezelter 637 }
338    
339 tim 1118 integrableObjects = (simnfo->molecules[localIndex]).getIntegrableObjects();
340     for(j=0; j < integrableObjects.size(); j++){
341    
342     eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
343     if(eof_test == NULL){
344     sprintf(painCave.errMsg,
345     "error in reading file %s\n"
346     "natoms = %d; index = %d\n"
347     "error reading the line from the file.\n",
348     inFileName.c_str(), nTotObjs, i );
349     haveError= 1;
350     simError();
351     }
352    
353     if(haveError) nodeZeroError();
354    
355 tim 1129 parseDumpLine(read_buffer, integrableObjects[j]);
356 tim 1118
357     }
358    
359    
360 gezelter 637 }
361 tim 1118 else{
362     //molecule belongs to slave nodes
363 gezelter 637
364 tim 1129 MPI_Recv(&nCurObj, 1, MPI_INT, which_node,
365 gezelter 637 TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus);
366    
367 tim 1130 for(j=0; j < nCurObj; j++){
368 tim 1118
369     eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
370     if(eof_test == NULL){
371     sprintf(painCave.errMsg,
372     "error in reading file %s\n"
373     "natoms = %d; index = %d\n"
374     "error reading the line from the file.\n",
375     inFileName.c_str(), nTotObjs, i );
376     haveError= 1;
377     simError();
378     }
379    
380     if(haveError) nodeZeroError();
381 gezelter 637
382 tim 1118 MPI_Send(read_buffer, BUFFERSIZE, MPI_CHAR, which_node,
383     TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD);
384    
385     }
386    
387 gezelter 637 }
388    
389     }
390 tim 1118
391 gezelter 637 }
392 tim 1118 else{
393     //actions taken at slave nodes
394 tim 1129 MPI_Bcast(read_buffer, BUFFERSIZE, MPI_CHAR, 0, MPI_COMM_WORLD);
395    
396     parseErr = parseCommentLine( read_buffer, simnfo);
397    
398     if( parseErr != NULL ){
399     strcpy( painCave.errMsg, parseErr );
400     haveError = 1;
401     simError();
402     }
403    
404 gezelter 1203 for (i=0 ; i < mpiSim->getNMolGlobal(); i++) {
405 tim 1118 which_node = MolToProcMap[i];
406    
407     if(which_node == worldRank){
408     //molecule with global index i belongs to this processor
409    
410     localIndex = mpiSim->getGlobalToLocalMol(i);
411 gezelter 637
412 tim 1118 if(localIndex == -1) {
413     sprintf(painCave.errMsg, "Molecule not found on node %d\n", worldRank);
414     haveError = 1;
415     simError();
416     }
417 gezelter 637
418 tim 1118 integrableObjects = (simnfo->molecules[localIndex]).getIntegrableObjects();
419 gezelter 637
420 tim 1118 nCurObj = integrableObjects.size();
421    
422 tim 1129 MPI_Send(&nCurObj, 1, MPI_INT, 0,
423     TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
424 gezelter 637
425 tim 1118 for(j = 0; j < integrableObjects.size(); j++){
426    
427     MPI_Recv(read_buffer, BUFFERSIZE, MPI_CHAR, 0,
428     TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD, &istatus);
429    
430     parseErr = parseDumpLine(read_buffer, integrableObjects[j]);
431    
432     if( parseErr != NULL ){
433     strcpy( painCave.errMsg, parseErr );
434     simError();
435     }
436    
437     }
438    
439     }
440    
441     }
442    
443     }
444    
445 gezelter 637 #endif
446     }
447    
448 tim 1118 char* DumpReader::parseDumpLine(char* readLine, StuntDouble* sd){
449 gezelter 637
450 tim 1118 char *foo; // the pointer to the current string token
451    
452     double pos[3]; // position place holders
453     double vel[3]; // velocity placeholders
454     double q[4]; // the quaternions
455     double ji[3]; // angular velocity placeholders;
456 gezelter 637 double qSqr, qLength; // needed to normalize the quaternion vector.
457    
458    
459     // set the string tokenizer
460 tim 1118
461 gezelter 637 foo = strtok(readLine, " ,;\t");
462 tim 1118
463 gezelter 637 // check the atom name to the current atom
464 tim 1118
465     if( strcmp( foo, sd->getType() ) ){
466     sprintf( painCave.errMsg,
467 tim 1119 "DumpReader error. Does not"
468 tim 1118 " match the BASS atom %s.\n",
469     sd->getType() );
470     return strdup( painCave.errMsg );
471     }
472    
473 gezelter 637 // get the positions
474    
475     foo = strtok(NULL, " ,;\t");
476     if(foo == NULL){
477     sprintf( painCave.errMsg,
478 tim 1118 "error in reading postition x from %s\n",
479     inFileName.c_str());
480 gezelter 637 return strdup( painCave.errMsg );
481     }
482 mmeineke 804 pos[0] = atof( foo );
483 tim 1118
484 gezelter 637 foo = strtok(NULL, " ,;\t");
485     if(foo == NULL){
486     sprintf( painCave.errMsg,
487 tim 1118 "error in reading postition y from %s\n",
488     inFileName.c_str());
489 gezelter 637 return strdup( painCave.errMsg );
490     }
491 mmeineke 804 pos[1] = atof( foo );
492 tim 1118
493 gezelter 637 foo = strtok(NULL, " ,;\t");
494     if(foo == NULL){
495     sprintf( painCave.errMsg,
496 tim 1118 "error in reading postition z from %s\n",
497     inFileName.c_str());
498 gezelter 637 return strdup( painCave.errMsg );
499     }
500 tim 1118 pos[2] = atof( foo );
501 gezelter 637
502    
503     // get the velocities
504    
505     foo = strtok(NULL, " ,;\t");
506     if(foo == NULL){
507     sprintf( painCave.errMsg,
508 tim 1118 "error in reading velocity x from %s\n",
509     inFileName.c_str() );
510 gezelter 637 return strdup( painCave.errMsg );
511     }
512 mmeineke 804 vel[0] = atof( foo );
513 tim 1118
514 gezelter 637 foo = strtok(NULL, " ,;\t");
515     if(foo == NULL){
516     sprintf( painCave.errMsg,
517 tim 1118 "error in reading velocity x from %s\n",
518     inFileName.c_str() );
519 gezelter 637 return strdup( painCave.errMsg );
520     }
521 mmeineke 804 vel[1] = atof( foo );
522 tim 1118
523 gezelter 637 foo = strtok(NULL, " ,;\t");
524     if(foo == NULL){
525     sprintf( painCave.errMsg,
526 tim 1118 "error in reading velocity x from %s\n",
527     inFileName.c_str() );
528 gezelter 637 return strdup( painCave.errMsg );
529     }
530 mmeineke 804 vel[2] = atof( foo );
531 tim 1118
532    
533     // add the positions and velocities to the atom
534    
535     sd->setPos( pos );
536     sd->setVel( vel );
537    
538     if (!sd->isDirectional())
539     return NULL;
540    
541 gezelter 637 // get the quaternions
542 tim 1118
543     if( sd->isDirectional() ){
544    
545 gezelter 637 foo = strtok(NULL, " ,;\t");
546     if(foo == NULL){
547 tim 1118 sprintf( painCave.errMsg,
548     "error in reading velocity x from %s\n",
549     inFileName.c_str() );
550 gezelter 637 return strdup( painCave.errMsg );
551     }
552     q[0] = atof( foo );
553 tim 1118
554 gezelter 637 foo = strtok(NULL, " ,;\t");
555     if(foo == NULL){
556     sprintf( painCave.errMsg,
557 tim 1118 "error in reading velocity x from %s\n",
558     inFileName.c_str() );
559 gezelter 637 return strdup( painCave.errMsg );
560     }
561     q[1] = atof( foo );
562 tim 1118
563 gezelter 637 foo = strtok(NULL, " ,;\t");
564     if(foo == NULL){
565     sprintf( painCave.errMsg,
566 tim 1118 "error in reading velocity x from %s\n",
567     inFileName.c_str() );
568 gezelter 637 return strdup( painCave.errMsg );
569     }
570     q[2] = atof( foo );
571 tim 1118
572 gezelter 637 foo = strtok(NULL, " ,;\t");
573     if(foo == NULL){
574     sprintf( painCave.errMsg,
575 tim 1118 "error in reading velocity x from %s\n",
576     inFileName.c_str() );
577 gezelter 637 return strdup( painCave.errMsg );
578     }
579     q[3] = atof( foo );
580 tim 1118
581 gezelter 637 // get the angular velocities
582 tim 1118
583 gezelter 637 foo = strtok(NULL, " ,;\t");
584     if(foo == NULL){
585     sprintf( painCave.errMsg,
586 tim 1118 "error in reading velocity x from %s\n",
587     inFileName.c_str() );
588 gezelter 637 return strdup( painCave.errMsg );
589     }
590 gezelter 1097 ji[0] = atof( foo );
591 tim 1118
592 gezelter 637 foo = strtok(NULL, " ,;\t");
593     if(foo == NULL){
594     sprintf( painCave.errMsg,
595 tim 1118 "error in reading velocity x from %s\n",
596     inFileName.c_str() );
597 gezelter 637 return strdup( painCave.errMsg );
598     }
599 gezelter 1097 ji[1] = atof(foo );
600 tim 1118
601 gezelter 637 foo = strtok(NULL, " ,;\t");
602     if(foo == NULL){
603     sprintf( painCave.errMsg,
604 tim 1118 "error in reading velocity x from %s\n",
605     inFileName.c_str() );
606 gezelter 637 return strdup( painCave.errMsg );
607     }
608 gezelter 1097 ji[2] = atof( foo );
609 gezelter 637
610 tim 1118
611 gezelter 637 // check that the quaternion vector is normalized
612    
613     qSqr = (q[0] * q[0]) + (q[1] * q[1]) + (q[2] * q[2]) + (q[3] * q[3]);
614 tim 1118
615 gezelter 637 qLength = sqrt( qSqr );
616     q[0] = q[0] / qLength;
617     q[1] = q[1] / qLength;
618     q[2] = q[2] / qLength;
619     q[3] = q[3] / qLength;
620    
621 tim 1118 // add quaternion and angular velocities
622    
623     sd->setQ( q );
624     sd->setJ( ji );
625 gezelter 637 }
626    
627 mmeineke 804
628 tim 1118
629 gezelter 637 return NULL;
630     }
631    
632    
633 tim 1118 char* DumpReader::parseCommentLine(char* readLine, SimInfo* entry_plug){
634 gezelter 637
635 tim 1118 double currTime;
636     double boxMat[9];
637     double theBoxMat3[3][3];
638     double chi;
639     double integralOfChidt;
640 mmeineke 847 double eta[9];
641 gezelter 637
642 tim 1118 char *foo; // the pointer to the current string token
643    
644 gezelter 637 // set the string tokenizer
645 tim 1118
646 gezelter 637 foo = strtok(readLine, " ,;\t");
647 tim 1118 // set the timeToken.
648    
649 gezelter 637 if(foo == NULL){
650     sprintf( painCave.errMsg,
651 tim 1118 "error in reading Time from %s\n",
652     inFileName.c_str() );
653 gezelter 637 return strdup( painCave.errMsg );
654     }
655    
656 tim 1118 currTime = atof( foo );
657     entry_plug->setTime( currTime );
658 gezelter 637
659 tim 1118 //get H-Matrix
660 gezelter 637
661 tim 1118 for(int i = 0 ; i < 9; i++){
662     foo = strtok(NULL, " ,;\t");
663     if(foo == NULL){
664     sprintf( painCave.errMsg,
665     "error in reading H[%d] from %s\n", i, inFileName.c_str() );
666     return strdup( painCave.errMsg );
667     }
668     boxMat[i] = atof( foo );
669 gezelter 637 }
670    
671 tim 1118 for(int i=0;i<3;i++)
672     for(int j=0;j<3;j++) theBoxMat3[i][j] = boxMat[3*j+i];
673 gezelter 637
674 tim 1118 //set H-Matrix
675     entry_plug->setBoxM( theBoxMat3 );
676 gezelter 637
677 tim 1118 //get chi and integralOfChidt, they should appear by pair
678 mmeineke 847
679 tim 1118 if( entry_plug->useInitXSstate ){
680     foo = strtok(NULL, " ,;\t\n");
681     if(foo != NULL){
682     chi = atof(foo);
683    
684     foo = strtok(NULL, " ,;\t\n");
685     if(foo == NULL){
686     sprintf( painCave.errMsg,
687     "chi and integralOfChidt should appear by pair in %s\n", inFileName.c_str() );
688     return strdup( painCave.errMsg );
689     }
690     integralOfChidt = atof( foo );
691    
692     //push chi and integralOfChidt into SimInfo::properties which can be
693     //retrieved by integrator later
694     DoubleData* chiValue = new DoubleData();
695     chiValue->setID(CHIVALUE_ID);
696     chiValue->setData(chi);
697     entry_plug->addProperty(chiValue);
698    
699     DoubleData* integralOfChidtValue = new DoubleData();
700     integralOfChidtValue->setID(INTEGRALOFCHIDT_ID);
701     integralOfChidtValue->setData(integralOfChidt);
702     entry_plug->addProperty(integralOfChidtValue);
703    
704     }
705     else
706     return NULL;
707 mmeineke 847
708 tim 1118 //get eta
709 mmeineke 847 foo = strtok(NULL, " ,;\t\n");
710 tim 1118 if(foo != NULL ){
711    
712     for(int i = 0 ; i < 9; i++){
713    
714     if(foo == NULL){
715     sprintf( painCave.errMsg,
716     "error in reading eta[%d] from %s\n", i, inFileName.c_str() );
717     return strdup( painCave.errMsg );
718     }
719     eta[i] = atof( foo );
720     foo = strtok(NULL, " ,;\t\n");
721     }
722 mmeineke 847 }
723 tim 1118 else
724     return NULL;
725 mmeineke 847
726 tim 1118 //push eta into SimInfo::properties which can be
727 mmeineke 847 //retrieved by integrator later
728 tim 1118 //entry_plug->setBoxM( theBoxMat3 );
729     DoubleArrayData* etaValue = new DoubleArrayData();
730     etaValue->setID(ETAVALUE_ID);
731     etaValue->setData(eta, 9);
732     entry_plug->addProperty(etaValue);
733 mmeineke 847 }
734 tim 1118
735 mmeineke 847 return NULL;
736 gezelter 637 }
737    
738     #ifdef IS_MPI
739 tim 1118 void DumpReader::nodeZeroError( void ){
740     int j, myStatus;
741 gezelter 637
742     myStatus = 0;
743 gezelter 1203 for (j = 0; j < mpiSim->getNProcessors(); j++) {
744 tim 1118 MPI_Send( &myStatus, 1, MPI_INT, j,
745 gezelter 637 TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
746 tim 1118 }
747 gezelter 637
748 tim 1118
749 gezelter 637 MPI_Finalize();
750     exit (0);
751 tim 1118
752 gezelter 637 }
753    
754 tim 1118 void DumpReader::anonymousNodeDie( void ){
755 gezelter 637
756     MPI_Finalize();
757     exit (0);
758     }
759 tim 1198 #endif