ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/InitializeFromFile.cpp
Revision: 643
Committed: Mon Jul 21 21:27:40 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 16414 byte(s)
Log Message:
some initial changes to Dumpwriter and friends to accomadate random file access

File Contents

# User Rev Content
1 mmeineke 377 #include <iostream>
2     #include <cmath>
3    
4     #include <stdio.h>
5     #include <stdlib.h>
6     #include <string.h>
7     #include <unistd.h>
8     #include <sys/types.h>
9     #include <sys/stat.h>
10    
11     #include "ReadWrite.hpp"
12     #include "simError.h"
13    
14     #ifdef IS_MPI
15 gezelter 417 #include <mpi.h>
16 mmeineke 377 #include "mpiSimulation.hpp"
17 chuckv 436 #define TAKE_THIS_TAG_CHAR 0
18     #define TAKE_THIS_TAG_INT 1
19    
20 mmeineke 440 namespace initFile{
21     void nodeZeroError( void );
22     void anonymousNodeDie( void );
23     }
24 chuckv 436
25 mmeineke 440 using namespace initFile;
26    
27 mmeineke 377 #endif // is_mpi
28    
29     InitializeFromFile :: InitializeFromFile( char *in_name ){
30     #ifdef IS_MPI
31     if (worldRank == 0) {
32     #endif
33    
34     c_in_file = fopen(in_name, "r");
35     if(c_in_file == NULL){
36     sprintf(painCave.errMsg,
37     "Cannot open file: %s\n", in_name);
38     painCave.isFatal = 1;
39     simError();
40     }
41    
42     strcpy( c_in_name, in_name);
43     #ifdef IS_MPI
44     }
45     strcpy( checkPointMsg, "Infile opened for reading successfully." );
46     MPIcheckPoint();
47     #endif
48     return;
49     }
50    
51     InitializeFromFile :: ~InitializeFromFile( ){
52     #ifdef IS_MPI
53     if (worldRank == 0) {
54     #endif
55     int error;
56     error = fclose( c_in_file );
57     if( error ){
58     sprintf( painCave.errMsg,
59     "Error closing %s\n", c_in_name );
60     simError();
61     }
62     #ifdef IS_MPI
63     }
64     strcpy( checkPointMsg, "Infile closed successfully." );
65     MPIcheckPoint();
66     #endif
67    
68     return;
69     }
70    
71    
72 mmeineke 641 void InitializeFromFile :: readInit( SimInfo* the_simnfo ){
73 mmeineke 377
74 gezelter 417 int i, j, done, which_node, which_atom; // loop counter
75 mmeineke 377
76     const int BUFFERSIZE = 2000; // size of the read buffer
77     int n_atoms; // the number of atoms
78     char read_buffer[BUFFERSIZE]; //the line buffer for reading
79     #ifdef IS_MPI
80     char send_buffer[BUFFERSIZE];
81     #endif
82    
83     char *eof_test; // ptr to see when we reach the end of the file
84     char *parseErr;
85     int procIndex;
86 mmeineke 586 double boxMat[9];
87 mmeineke 590 double theBoxMat3[3][3];
88 mmeineke 377
89 mmeineke 586 simnfo = the_simnfo;
90 mmeineke 377
91    
92     #ifndef IS_MPI
93     eof_test = fgets(read_buffer, sizeof(read_buffer), c_in_file);
94     if( eof_test == NULL ){
95     sprintf( painCave.errMsg,
96     "InitializeFromFile error: error reading 1st line of \"%s\"\n",
97     c_in_name );
98     painCave.isFatal = 1;
99     simError();
100     }
101    
102     n_atoms = atoi( read_buffer );
103    
104 mmeineke 586 Atom **atoms = simnfo->atoms;
105 mmeineke 377 DirectionalAtom* dAtom;
106    
107 mmeineke 586 if( n_atoms != simnfo->n_atoms ){
108 mmeineke 377 sprintf( painCave.errMsg,
109     "Initialize from File error. %s n_atoms, %d, "
110     "does not match the BASS file's n_atoms, %d.\n",
111 mmeineke 586 c_in_name, n_atoms, simnfo->n_atoms );
112 mmeineke 377 painCave.isFatal = 1;
113     simError();
114     }
115    
116 mmeineke 586 //read the box mat from the comment line
117 mmeineke 377
118     eof_test = fgets(read_buffer, sizeof(read_buffer), c_in_file);
119     if(eof_test == NULL){
120     sprintf( painCave.errMsg,
121     "error in reading commment in %s\n", c_in_name);
122     painCave.isFatal = 1;
123     simError();
124     }
125    
126 mmeineke 586 parseErr = parseBoxLine( read_buffer, boxMat );
127     if( parseErr != NULL ){
128     strcpy( painCave.errMsg, parseErr );
129     painCave.isFatal = 1;
130     simError();
131     }
132    
133 mmeineke 590 for(i=0;i<3;i++)
134     for(j=0;j<3;j++) theBoxMat3[i][j] = boxMat[3*j+i];
135 mmeineke 586
136 mmeineke 590 simnfo->setBoxM( theBoxMat3 );
137 mmeineke 586
138 mmeineke 590
139 mmeineke 377 for( i=0; i < n_atoms; i++){
140    
141     eof_test = fgets(read_buffer, sizeof(read_buffer), c_in_file);
142     if(eof_test == NULL){
143     sprintf(painCave.errMsg,
144     "error in reading file %s\n"
145     "natoms = %d; index = %d\n"
146     "error reading the line from the file.\n",
147     c_in_name, n_atoms, i );
148     painCave.isFatal = 1;
149     simError();
150     }
151    
152    
153     parseErr = parseDumpLine( read_buffer, i );
154     if( parseErr != NULL ){
155     strcpy( painCave.errMsg, parseErr );
156     painCave.isFatal = 1;
157     simError();
158 mmeineke 586 }
159 mmeineke 377 }
160    
161    
162     // MPI Section of code..........
163     #else //IS_MPI
164    
165 chuckv 436 // first thing first, suspend fatalities.
166     painCave.isEventLoop = 1;
167    
168     int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
169     int haveError;
170    
171 mmeineke 447 MPI_Status istatus;
172 gezelter 417 int *AtomToProcMap = mpiSim->getAtomToProcMap();
173 mmeineke 377
174 chuckv 436
175     haveError = 0;
176 mmeineke 377 if (worldRank == 0) {
177 chuckv 436
178 mmeineke 377 eof_test = fgets(read_buffer, sizeof(read_buffer), c_in_file);
179     if( eof_test == NULL ){
180     sprintf( painCave.errMsg,
181     "Error reading 1st line of %d \n ",c_in_name);
182 chuckv 436 haveError = 1;
183 mmeineke 377 simError();
184     }
185    
186     n_atoms = atoi( read_buffer );
187    
188 mmeineke 586 Atom **atoms = simnfo->atoms;
189 mmeineke 377 DirectionalAtom* dAtom;
190    
191     // Check to see that the number of atoms in the intial configuration file is the
192     // same as declared in simBass.
193    
194     if( n_atoms != mpiSim->getTotAtoms() ){
195     sprintf( painCave.errMsg,
196     "Initialize from File error. %s n_atoms, %d, "
197     "does not match the BASS file's n_atoms, %d.\n",
198 mmeineke 586 c_in_name, n_atoms, simnfo->n_atoms );
199 chuckv 436 haveError= 1;
200 mmeineke 377 simError();
201     }
202    
203 mmeineke 586 //read the boxMat from the comment line
204 mmeineke 377
205     eof_test = fgets(read_buffer, sizeof(read_buffer), c_in_file);
206     if(eof_test == NULL){
207     sprintf( painCave.errMsg,
208     "error in reading commment in %s\n", c_in_name);
209 mmeineke 586 haveError = 1;
210 mmeineke 377 simError();
211     }
212 mmeineke 586
213     parseErr = parseBoxLine( read_buffer, boxMat );
214     if( parseErr != NULL ){
215     strcpy( painCave.errMsg, parseErr );
216     haveError = 1;
217     simError();
218     }
219    
220     MPI_Bcast(boxMat, 9, MPI_DOUBLE, 0, MPI_COMM_WORLD );
221 mmeineke 377
222 chuckv 436 if(haveError) nodeZeroError();
223 gezelter 417
224     for (i=0 ; i < mpiSim->getTotAtoms(); i++) {
225    
226 mmeineke 377 eof_test = fgets(read_buffer, sizeof(read_buffer), c_in_file);
227     if(eof_test == NULL){
228     sprintf(painCave.errMsg,
229     "error in reading file %s\n"
230     "natoms = %d; index = %d\n"
231     "error reading the line from the file.\n",
232     c_in_name, n_atoms, i );
233 chuckv 436 haveError= 1;
234 mmeineke 377 simError();
235     }
236    
237 chuckv 436 if(haveError) nodeZeroError();
238    
239 gezelter 417 // Get the Node number which wants this atom:
240     which_node = AtomToProcMap[i];
241 chuckv 436 if (which_node == 0) {
242 gezelter 417 parseErr = parseDumpLine( read_buffer, i );
243     if( parseErr != NULL ){
244     strcpy( painCave.errMsg, parseErr );
245 chuckv 436 haveError = 1;
246 gezelter 417 simError();
247     }
248 chuckv 436 if(haveError) nodeZeroError();
249     }
250    
251     else {
252    
253     myStatus = 1;
254 mmeineke 447 MPI_Send(&myStatus, 1, MPI_INT, which_node,
255     TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
256     MPI_Send(read_buffer, BUFFERSIZE, MPI_CHAR, which_node,
257     TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD);
258     MPI_Send(&i, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT,
259     MPI_COMM_WORLD);
260     MPI_Recv(&myStatus, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT,
261     MPI_COMM_WORLD, &istatus);
262 chuckv 436
263     if(!myStatus) nodeZeroError();
264 mmeineke 377 }
265     }
266 chuckv 436 myStatus = -1;
267 gezelter 417 for (j = 0; j < mpiSim->getNumberProcessors(); j++) {
268 mmeineke 447 MPI_Send( &myStatus, 1, MPI_INT, j,
269     TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
270 gezelter 417 }
271 gezelter 419
272 gezelter 417 } else {
273 gezelter 419
274 mmeineke 586 MPI_Bcast(boxMat, 9, MPI_DOUBLE, 0, MPI_COMM_WORLD);
275    
276 gezelter 417 done = 0;
277     while (!done) {
278 chuckv 436
279 mmeineke 447 MPI_Recv(&myStatus, 1, MPI_INT, 0,
280     TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus);
281 chuckv 436
282     if(!myStatus) anonymousNodeDie();
283    
284     if(myStatus < 0) break;
285    
286 mmeineke 447 MPI_Recv(read_buffer, BUFFERSIZE, MPI_CHAR, 0,
287     TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD, &istatus);
288     MPI_Recv(&which_atom, 1, MPI_INT, 0,
289     TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus);
290 chuckv 436
291     myStatus = 1;
292     parseErr = parseDumpLine( read_buffer, which_atom );
293     if( parseErr != NULL ){
294     strcpy( painCave.errMsg, parseErr );
295     myStatus = 0;;
296     simError();
297 mmeineke 377 }
298 chuckv 436
299 mmeineke 447 MPI_Send( &myStatus, 1, MPI_INT, 0,
300     TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
301 chuckv 436
302 mmeineke 377 }
303     }
304 chuckv 436
305     // last thing last, enable fatalities.
306     painCave.isEventLoop = 0;
307 mmeineke 586
308 mmeineke 590 for(i=0;i<3;i++)
309     for(j=0;j<3;j++) theBoxMat3[i][j] = boxMat[3*j+i];
310    
311     simnfo->setBoxM( theBoxMat3 );
312    
313 chuckv 436
314 mmeineke 377 #endif
315     }
316    
317 gezelter 419 char* InitializeFromFile::parseDumpLine(char* readLine, int globalIndex){
318 mmeineke 377
319     char *foo; // the pointer to the current string token
320    
321     double rx, ry, rz; // position place holders
322     double vx, vy, vz; // velocity placeholders
323     double q[4]; // the quaternions
324     double jx, jy, jz; // angular velocity placeholders;
325     double qSqr, qLength; // needed to normalize the quaternion vector.
326    
327 mmeineke 586 Atom **atoms = simnfo->atoms;
328 mmeineke 377 DirectionalAtom* dAtom;
329    
330 gezelter 419 int j, n_atoms, atomIndex;
331 mmeineke 377
332     #ifdef IS_MPI
333     n_atoms = mpiSim->getTotAtoms();
334 gezelter 419 atomIndex=-1;
335     for (j=0; j < mpiSim->getMyNlocal(); j++) {
336     if (atoms[j]->getGlobalIndex() == globalIndex) atomIndex = j;
337     }
338     if (atomIndex == -1) {
339     sprintf( painCave.errMsg,
340     "Initialize from file error. Atom at index %d "
341     "in file %s does not exist on processor %d .\n",
342     globalIndex, c_in_name, mpiSim->getMyNode() );
343     return strdup( painCave.errMsg );
344     }
345 mmeineke 377 #else
346 mmeineke 586 n_atoms = simnfo->n_atoms;
347 gezelter 419 atomIndex = globalIndex;
348 mmeineke 377 #endif // is_mpi
349    
350     // set the string tokenizer
351    
352     foo = strtok(readLine, " ,;\t");
353    
354     // check the atom name to the current atom
355    
356     if( strcmp( foo, atoms[atomIndex]->getType() ) ){
357     sprintf( painCave.errMsg,
358     "Initialize from file error. Atom %s at index %d "
359     "in file %s does not"
360     " match the BASS atom %s.\n",
361     foo, atomIndex, c_in_name, atoms[atomIndex]->getType() );
362     return strdup( painCave.errMsg );
363     }
364    
365     // get the positions
366    
367     foo = strtok(NULL, " ,;\t");
368     if(foo == NULL){
369     sprintf( painCave.errMsg,
370     "error in reading postition x from %s\n"
371     "natoms = %d, index = %d\n",
372     c_in_name, n_atoms, atomIndex );
373     return strdup( painCave.errMsg );
374     }
375     rx = atof( foo );
376    
377     foo = strtok(NULL, " ,;\t");
378     if(foo == NULL){
379     sprintf( painCave.errMsg,
380     "error in reading postition y from %s\n"
381     "natoms = %d, index = %d\n",
382     c_in_name, n_atoms, atomIndex );
383     return strdup( painCave.errMsg );
384     }
385     ry = atof( foo );
386    
387     foo = strtok(NULL, " ,;\t");
388     if(foo == NULL){
389     sprintf( painCave.errMsg,
390     "error in reading postition z from %s\n"
391     "natoms = %d, index = %d\n",
392     c_in_name, n_atoms, atomIndex );
393     return strdup( painCave.errMsg );
394     }
395     rz = atof( foo );
396    
397    
398     // get the velocities
399    
400     foo = strtok(NULL, " ,;\t");
401     if(foo == NULL){
402     sprintf( painCave.errMsg,
403     "error in reading velocity x from %s\n"
404     "natoms = %d, index = %d\n",
405     c_in_name, n_atoms, atomIndex );
406     return strdup( painCave.errMsg );
407     }
408     vx = atof( foo );
409    
410     foo = strtok(NULL, " ,;\t");
411     if(foo == NULL){
412     sprintf( painCave.errMsg,
413     "error in reading velocity y from %s\n"
414     "natoms = %d, index = %d\n",
415     c_in_name, n_atoms, atomIndex );
416     return strdup( painCave.errMsg );
417     }
418     vy = atof( foo );
419    
420     foo = strtok(NULL, " ,;\t");
421     if(foo == NULL){
422     sprintf( painCave.errMsg,
423     "error in reading velocity z from %s\n"
424     "natoms = %d, index = %d\n",
425     c_in_name, n_atoms, atomIndex );
426     return strdup( painCave.errMsg );
427     }
428     vz = atof( foo );
429    
430    
431     // get the quaternions
432    
433     if( atoms[atomIndex]->isDirectional() ){
434    
435     foo = strtok(NULL, " ,;\t");
436     if(foo == NULL){
437     sprintf(painCave.errMsg,
438     "error in reading quaternion 0 from %s\n"
439     "natoms = %d, index = %d\n",
440     c_in_name, n_atoms, atomIndex );
441     return strdup( painCave.errMsg );
442     }
443     q[0] = atof( foo );
444    
445     foo = strtok(NULL, " ,;\t");
446     if(foo == NULL){
447     sprintf( painCave.errMsg,
448     "error in reading quaternion 1 from %s\n"
449     "natoms = %d, index = %d\n",
450     c_in_name, n_atoms, atomIndex );
451     return strdup( painCave.errMsg );
452     }
453     q[1] = atof( foo );
454    
455     foo = strtok(NULL, " ,;\t");
456     if(foo == NULL){
457     sprintf( painCave.errMsg,
458     "error in reading quaternion 2 from %s\n"
459     "natoms = %d, index = %d\n",
460     c_in_name, n_atoms, atomIndex );
461     return strdup( painCave.errMsg );
462     }
463     q[2] = atof( foo );
464    
465     foo = strtok(NULL, " ,;\t");
466     if(foo == NULL){
467     sprintf( painCave.errMsg,
468     "error in reading quaternion 3 from %s\n"
469     "natoms = %d, index = %d\n",
470     c_in_name, n_atoms, atomIndex );
471     return strdup( painCave.errMsg );
472     }
473     q[3] = atof( foo );
474    
475     // get the angular velocities
476    
477     foo = strtok(NULL, " ,;\t");
478     if(foo == NULL){
479     sprintf( painCave.errMsg,
480     "error in reading angular momentum jx from %s\n"
481     "natoms = %d, index = %d\n",
482     c_in_name, n_atoms, atomIndex );
483     return strdup( painCave.errMsg );
484     }
485     jx = atof( foo );
486    
487     foo = strtok(NULL, " ,;\t");
488     if(foo == NULL){
489     sprintf( painCave.errMsg,
490     "error in reading angular momentum jy from %s\n"
491     "natoms = %d, index = %d\n",
492     c_in_name, n_atoms, atomIndex );
493     return strdup( painCave.errMsg );
494     }
495     jy = atof(foo );
496    
497     foo = strtok(NULL, " ,;\t");
498     if(foo == NULL){
499     sprintf( painCave.errMsg,
500     "error in reading angular momentum jz from %s\n"
501     "natoms = %d, index = %d\n",
502     c_in_name, n_atoms, atomIndex );
503     return strdup( painCave.errMsg );
504     }
505     jz = atof( foo );
506    
507     dAtom = ( DirectionalAtom* )atoms[atomIndex];
508    
509     // check that the quaternion vector is normalized
510    
511     qSqr = (q[0] * q[0]) + (q[1] * q[1]) + (q[2] * q[2]) + (q[3] * q[3]);
512    
513     qLength = sqrt( qSqr );
514     q[0] = q[0] / qLength;
515     q[1] = q[1] / qLength;
516     q[2] = q[2] / qLength;
517     q[3] = q[3] / qLength;
518    
519     dAtom->setQ( q );
520    
521     // add the angular velocities
522    
523     dAtom->setJx( jx );
524     dAtom->setJy( jy );
525     dAtom->setJz( jz );
526     }
527    
528     // add the positions and velocities to the atom
529    
530     atoms[atomIndex]->setX( rx );
531     atoms[atomIndex]->setY( ry );
532     atoms[atomIndex]->setZ( rz );
533    
534     atoms[atomIndex]->set_vx( vx );
535     atoms[atomIndex]->set_vy( vy );
536     atoms[atomIndex]->set_vz( vz );
537    
538     return NULL;
539     }
540 chuckv 436
541    
542 mmeineke 586 char* InitializeFromFile::parseBoxLine(char* readLine, double boxMat[9]){
543    
544     char *foo; // the pointer to the current string token
545     int j;
546    
547     // set the string tokenizer
548    
549     foo = strtok(readLine, " ,;\t");
550 mmeineke 641 // set the timeToken.
551 mmeineke 586
552 mmeineke 641 if(foo == NULL){
553     sprintf( painCave.errMsg,
554     "error in reading Time from %s\n",
555     c_in_name );
556     return strdup( painCave.errMsg );
557     }
558 mmeineke 643 simnfo->setTime( atof( foo ) );
559 mmeineke 586
560     // get the Hx vector
561    
562     foo = strtok(NULL, " ,;\t");
563     if(foo == NULL){
564     sprintf( painCave.errMsg,
565     "error in reading Hx[0] from %s\n",
566     c_in_name );
567     return strdup( painCave.errMsg );
568     }
569     boxMat[0] = atof( foo );
570    
571     foo = strtok(NULL, " ,;\t");
572     if(foo == NULL){
573     sprintf( painCave.errMsg,
574     "error in reading Hx[1] from %s\n",
575     c_in_name );
576     return strdup( painCave.errMsg );
577     }
578     boxMat[1] = atof( foo );
579    
580     foo = strtok(NULL, " ,;\t");
581     if(foo == NULL){
582     sprintf( painCave.errMsg,
583     "error in reading Hx[2] from %s\n",
584     c_in_name );
585     return strdup( painCave.errMsg );
586     }
587     boxMat[2] = atof( foo );
588    
589     // get the Hy vector
590    
591     foo = strtok(NULL, " ,;\t");
592     if(foo == NULL){
593     sprintf( painCave.errMsg,
594     "error in reading Hy[0] from %s\n",
595     c_in_name );
596     return strdup( painCave.errMsg );
597     }
598     boxMat[3] = atof( foo );
599    
600     foo = strtok(NULL, " ,;\t");
601     if(foo == NULL){
602     sprintf( painCave.errMsg,
603     "error in reading Hy[1] from %s\n",
604     c_in_name );
605     return strdup( painCave.errMsg );
606     }
607     boxMat[4] = atof( foo );
608    
609     foo = strtok(NULL, " ,;\t");
610     if(foo == NULL){
611     sprintf( painCave.errMsg,
612     "error in reading Hy[2] from %s\n",
613     c_in_name );
614     return strdup( painCave.errMsg );
615     }
616     boxMat[5] = atof( foo );
617    
618     // get the Hz vector
619    
620     foo = strtok(NULL, " ,;\t");
621     if(foo == NULL){
622     sprintf( painCave.errMsg,
623     "error in reading Hz[0] from %s\n",
624     c_in_name );
625     return strdup( painCave.errMsg );
626     }
627     boxMat[6] = atof( foo );
628    
629     foo = strtok(NULL, " ,;\t");
630     if(foo == NULL){
631     sprintf( painCave.errMsg,
632     "error in reading Hz[1] from %s\n",
633     c_in_name );
634     return strdup( painCave.errMsg );
635     }
636     boxMat[7] = atof( foo );
637    
638     foo = strtok(NULL, " ,;\t");
639     if(foo == NULL){
640     sprintf( painCave.errMsg,
641     "error in reading Hz[2] from %s\n",
642     c_in_name );
643     return strdup( painCave.errMsg );
644     }
645     boxMat[8] = atof( foo );
646    
647     return NULL;
648     }
649    
650    
651 chuckv 436 #ifdef IS_MPI
652    
653     // a couple of functions to let us escape the read loop
654    
655 mmeineke 440 void initFile::nodeZeroError( void ){
656 chuckv 436 int j, myStatus;
657    
658     myStatus = 0;
659     for (j = 0; j < mpiSim->getNumberProcessors(); j++) {
660 mmeineke 447 MPI_Send( &myStatus, 1, MPI_INT, j,
661     TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
662 chuckv 436 }
663    
664    
665     MPI_Finalize();
666     exit (0);
667    
668     }
669    
670 mmeineke 440 void initFile::anonymousNodeDie( void ){
671 chuckv 436
672     MPI_Finalize();
673     exit (0);
674     }
675    
676     #endif //is_mpi