ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpReader.cpp
Revision: 1097
Committed: Mon Apr 12 20:32:20 2004 UTC (20 years, 2 months ago) by gezelter
File size: 20455 byte(s)
Log Message:
Changes for RigidBody dynamics (Somewhat extensive)

File Contents

# Content
1 #define _LARGEFILE_SOURCE64
2 #define _FILE_OFFSET_BITS 64
3
4 #include <sys/types.h>
5 #include <sys/stat.h>
6
7 #include <iostream>
8 #include <math.h>
9
10 #include <stdio.h>
11 #include <stdlib.h>
12 #include <string.h>
13
14
15 #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 #endif // is_mpi
24
25 namespace dumpRead{
26
27 #ifdef IS_MPI
28 void nodeZeroError( void );
29 void anonymousNodeDie( void );
30 #endif // is_mpi
31
32 }
33
34 using namespace dumpRead;
35
36
37 DumpReader :: DumpReader(const char *in_name ){
38
39 isScanned = false;
40 headFP = new FilePos;
41
42 #ifdef IS_MPI
43 if (worldRank == 0) {
44 #endif
45
46 inFile = fopen(in_name, "r");
47 if(inFile == NULL){
48 sprintf(painCave.errMsg,
49 "Cannot open file: %s\n", in_name);
50 painCave.isFatal = 1;
51 simError();
52 }
53
54 strcpy( inName, in_name);
55 #ifdef IS_MPI
56 }
57 strcpy( checkPointMsg, "Dump file opened for reading successfully." );
58 MPIcheckPoint();
59 #endif
60 return;
61 }
62
63 DumpReader :: ~DumpReader( ){
64 #ifdef IS_MPI
65 if (worldRank == 0) {
66 #endif
67 int error;
68 error = fclose( inFile );
69 if( error ){
70 sprintf( painCave.errMsg,
71 "Error closing %s\n", inName );
72 simError();
73 }
74 #ifdef IS_MPI
75 }
76 strcpy( checkPointMsg, "Dump file closed successfully." );
77 MPIcheckPoint();
78 #endif
79
80 return;
81 }
82
83 int DumpReader::getNframes( void ){
84
85 if( !isScanned ) scanFile();
86 return nFrames;
87 }
88
89 void DumpReader::scanFile( void ){
90
91 int vectorSize;
92 int i, j, k;
93 int lineNum = 0;
94 char readBuffer[2000];
95 char* foo;
96 fpos_t *currPos;
97 double time;
98
99 FilePos *currFP;
100
101
102 #ifdef IS_MPI
103 if( worldRank == 0 ){
104 #endif // is_mpi
105
106 rewind( inFile );
107
108 currPos = new fpos_t;
109 fgetpos( inFile, currPos );
110 fgets( readBuffer, sizeof( readBuffer ), inFile );
111 lineNum++;
112 if( feof( inFile ) ){
113 sprintf( painCave.errMsg,
114 "File \"%s\" ended unexpectedly at line %d\n",
115 inName,
116 lineNum );
117 painCave.isFatal = 1;
118 simError();
119 }
120
121 nFrames = 0;
122 while( !feof( inFile ) ){
123
124 headFP->add( currPos );
125 nFrames++;
126
127 i = atoi(readBuffer);
128
129 fgets( readBuffer, sizeof( readBuffer ), inFile );
130 lineNum++;
131 if( feof( inFile ) ){
132 sprintf( painCave.errMsg,
133 "File \"%s\" ended unexpectedly at line %d\n",
134 inName,
135 lineNum );
136 painCave.isFatal = 1;
137 simError();
138 }
139
140 // if(outTime){
141 // foo = strtok( readBuffer, " ,;\t" );
142 // time = atof( foo );
143 // }
144
145 for(j=0; j<i; j++){
146
147 fgets( readBuffer, sizeof( readBuffer ), inFile );
148 lineNum++;
149 if( feof( inFile ) ){
150 sprintf( painCave.errMsg,
151 "File \"%s\" ended unexpectedly at line %d,"
152 " with atom %d\n",
153 inName,
154 lineNum,
155 j );
156 painCave.isFatal = 1;
157 simError();
158 }
159
160 }
161
162 currPos = new fpos_t;
163 fgetpos( inFile, currPos );
164 fgets( readBuffer, sizeof( readBuffer ), inFile );
165 lineNum++;
166 }
167
168 delete currPos;
169 rewind( inFile );
170
171 frameStart = new FilePos*[nFrames];
172 currFP = headFP;
173 for(i=0; i<nFrames; i++){
174
175 currFP = currFP->getNext();
176 if( currFP == NULL ){
177 sprintf( painCave.errMsg,
178 "DumpReader error: scanFile FilePos mismatch at "
179 "nFrames = %d\n",
180 i );
181 painCave.isFatal = 1;
182 simError();
183 }
184
185 frameStart[i] = currFP;
186 }
187
188 isScanned = true;
189
190 #ifdef IS_MPI
191 }
192 strcpy( checkPointMsg, "Successfully scanned DumpFile\n" );
193 MPIcheckPoint();
194 #endif // is_mpi
195 }
196
197 void DumpReader :: readFrame( SimInfo* the_simnfo, int whichFrame){
198
199 simnfo = the_simnfo;
200
201 this->readSet( whichFrame );
202 }
203
204
205
206 void DumpReader :: readSet( int whichFrame ){
207
208 int i, j, done, which_node, which_atom; // loop counter
209
210 const int BUFFERSIZE = 2000; // size of the read buffer
211 int n_atoms; // the number of atoms
212 char read_buffer[BUFFERSIZE]; //the line buffer for reading
213 #ifdef IS_MPI
214 char send_buffer[BUFFERSIZE];
215 #endif
216
217 char *eof_test; // ptr to see when we reach the end of the file
218 char *parseErr;
219 int procIndex;
220 double boxMat[9];
221 double theBoxMat3[3][3];
222 double time;
223
224 fpos_t *framePos;
225 framePos = frameStart[whichFrame]->getPos();
226
227 #ifndef IS_MPI
228
229 fsetpos(inFile, framePos);
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 fsetpos(inFile, framePos);
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 ji[3]; // 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 atoms[atomIndex]->setType(foo);
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 ji[0] = 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 ji[1] = 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 ji[2] = 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->setJ( ji );
671 }
672
673 // add the positions and velocities to the atom
674
675 atoms[atomIndex]->setPos( pos );
676 atoms[atomIndex]->setVel( vel );
677
678 return NULL;
679 }
680
681
682 char* DumpReader::parseCommentLine(char* readLine, double &time,
683 double boxMat[9]){
684
685 char *foo; // the pointer to the current string token
686 int j;
687 double chi, integralOfChidt;
688 double eta[9];
689
690 // set the string tokenizer
691
692 foo = strtok(readLine, " ,;\t");
693 if(foo == NULL){
694 sprintf( painCave.errMsg,
695 "error in reading time from %s\n",
696 inName );
697 return strdup( painCave.errMsg );
698 }
699 time = atof( foo );
700
701 // get the Hx vector
702
703 foo = strtok(NULL, " ,;\t");
704 if(foo == NULL){
705 sprintf( painCave.errMsg,
706 "error in reading Hx[0] from %s\n",
707 inName );
708 return strdup( painCave.errMsg );
709 }
710 boxMat[0] = atof( foo );
711
712 foo = strtok(NULL, " ,;\t");
713 if(foo == NULL){
714 sprintf( painCave.errMsg,
715 "error in reading Hx[1] from %s\n",
716 inName );
717 return strdup( painCave.errMsg );
718 }
719 boxMat[1] = atof( foo );
720
721 foo = strtok(NULL, " ,;\t");
722 if(foo == NULL){
723 sprintf( painCave.errMsg,
724 "error in reading Hx[2] from %s\n",
725 inName );
726 return strdup( painCave.errMsg );
727 }
728 boxMat[2] = atof( foo );
729
730 // get the Hy vector
731
732 foo = strtok(NULL, " ,;\t");
733 if(foo == NULL){
734 sprintf( painCave.errMsg,
735 "error in reading Hy[0] from %s\n",
736 inName );
737 return strdup( painCave.errMsg );
738 }
739 boxMat[3] = atof( foo );
740
741 foo = strtok(NULL, " ,;\t");
742 if(foo == NULL){
743 sprintf( painCave.errMsg,
744 "error in reading Hy[1] from %s\n",
745 inName );
746 return strdup( painCave.errMsg );
747 }
748 boxMat[4] = atof( foo );
749
750 foo = strtok(NULL, " ,;\t");
751 if(foo == NULL){
752 sprintf( painCave.errMsg,
753 "error in reading Hy[2] from %s\n",
754 inName );
755 return strdup( painCave.errMsg );
756 }
757 boxMat[5] = atof( foo );
758
759 // get the Hz vector
760
761 foo = strtok(NULL, " ,;\t");
762 if(foo == NULL){
763 sprintf( painCave.errMsg,
764 "error in reading Hz[0] from %s\n",
765 inName );
766 return strdup( painCave.errMsg );
767 }
768 boxMat[6] = atof( foo );
769
770 foo = strtok(NULL, " ,;\t");
771 if(foo == NULL){
772 sprintf( painCave.errMsg,
773 "error in reading Hz[1] from %s\n",
774 inName );
775 return strdup( painCave.errMsg );
776 }
777 boxMat[7] = atof( foo );
778
779 foo = strtok(NULL, " ,;\t");
780 if(foo == NULL){
781 sprintf( painCave.errMsg,
782 "error in reading Hz[2] from %s\n",
783 inName );
784 return strdup( painCave.errMsg );
785 }
786 boxMat[8] = atof( foo );
787
788 return NULL;
789
790 //get chi and integralOfChidt, they should appear by pair
791 foo = strtok(NULL, " ,;\t\n");
792 if(foo != NULL){
793 chi = atof(foo);
794
795 foo = strtok(NULL, " ,;\t\n");
796 if(foo == NULL){
797 sprintf( painCave.errMsg,
798 "chi and integralOfChidt should appear by pair in %s\n", inName );
799 return strdup( painCave.errMsg );
800 }
801 integralOfChidt = atof( foo );
802
803 //push chi and integralOfChidt into SimInfo::properties which can be
804 //retrieved by integrator later
805 DoubleData* chiValue = new DoubleData();
806 chiValue->setID(CHIVALUE_ID);
807 chiValue->setData(chi);
808 simnfo->addProperty(chiValue);
809
810 DoubleData* integralOfChidtValue = new DoubleData();
811 integralOfChidtValue->setID(INTEGRALOFCHIDT_ID);
812 integralOfChidtValue->setData(integralOfChidt);
813 simnfo->addProperty(integralOfChidtValue);
814
815 }
816 else
817 return NULL;
818
819 //get eta
820 for(int i = 0 ; i < 9; i++){
821 foo = strtok(NULL, " ,;\t");
822 if(foo == NULL){
823 sprintf( painCave.errMsg,
824 "error in reading eta[%d] from %s\n", i, inName );
825 return strdup( painCave.errMsg );
826 }
827 eta[i] = atof( foo );
828 }
829
830 //push eta into SimInfo::properties which can be
831 //retrieved by integrator later
832 //simnfo->setBoxM( theBoxMat3 );
833 DoubleArrayData* etaValue = new DoubleArrayData();
834 etaValue->setID(ETAVALUE_ID);
835 etaValue->setData(eta, 9);
836 simnfo->addProperty(etaValue);
837
838
839 return NULL;
840
841
842
843 }
844
845
846 #ifdef IS_MPI
847
848 // a couple of functions to let us escape the read loop
849
850 void dumpRead::nodeZeroError( void ){
851 int j, myStatus;
852
853 myStatus = 0;
854 for (j = 0; j < mpiSim->getNumberProcessors(); j++) {
855 MPI_Send( &myStatus, 1, MPI_INT, j,
856 TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
857 }
858
859
860 MPI_Finalize();
861 exit (0);
862
863 }
864
865 void dumpRead::anonymousNodeDie( void ){
866
867 MPI_Finalize();
868 exit (0);
869 }
870
871 #endif //is_mpi