ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpReader.cpp
Revision: 847
Committed: Fri Oct 31 18:28:52 2003 UTC (20 years, 8 months ago) by mmeineke
File size: 20329 byte(s)
Log Message:
added template stuff to the Maikefile template

little changes to some printf format statements

File Contents

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