ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/InitializeFromFile.cpp
Revision: 1129
Committed: Thu Apr 22 03:29:30 2004 UTC (20 years, 2 months ago) by tim
File size: 15998 byte(s)
Log Message:
fixed two bugs in MPI version of InitfromFile and one unmatch MPI command in DumpWriter

File Contents

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