ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpWriter.cpp
(Generate patch)

Comparing trunk/OOPSE/libmdtools/DumpWriter.cpp (file contents):
Revision 912 by gezelter, Thu Jan 8 18:59:36 2004 UTC vs.
Revision 927 by tim, Mon Jan 12 22:54:42 2004 UTC

# Line 72 | Line 72 | void DumpWriter::writeDump( double currentTime ){
72    char writeLine[BUFFERSIZE];
73  
74    int i;
75 +
76   #ifdef IS_MPI
77 +  
78 +  int *potatoes;
79 +  int myPotato;
80 +
81 +  int nProc;
82    int j, which_node, done, which_atom, local_index;
83 <  double atomTransData[6];
84 <  double atomOrientData[7];
83 >  double atomData6[6];
84 >  double atomData13[13];
85    int isDirectional;
86    char* atomTypeString;
87    char MPIatomTypeString[MINIBUFFERSIZE];
88 <  int me;
83 <  int atomTypeTag;
84 <  int atomIsDirectionalTag;
85 <  int atomTransDataTag;
86 <  int atomOrientDataTag;
88 >
89   #else //is_mpi
90    int nAtoms = entry_plug->n_atoms;
91   #endif //is_mpi
# Line 158 | Line 160 | void DumpWriter::writeDump( double currentTime ){
160  
161   #else // is_mpi
162  
163 <  // first thing first, suspend fatalities.
164 <  painCave.isEventLoop = 1;
163 >  /* code to find maximum tag value */
164 >  
165 >  int *tagub, flag, MAXTAG;
166 >  MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, &tagub, &flag);
167 >  if (flag) {
168 >    MAXTAG = *tagub;
169 >  } else {
170 >    MAXTAG = 32767;
171 >  }  
172  
164  int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
173    int haveError;
174  
175    MPI_Status istatus;
# Line 170 | Line 178 | void DumpWriter::writeDump( double currentTime ){
178    // write out header and node 0's coordinates
179  
180    if( worldRank == 0 ){
181 +
182 +    // Node 0 needs a list of the magic potatoes for each processor;
183 +
184 +    nProc = mpiSim->getNumberProcessors();
185 +    potatoes = new int[nProc];
186 +
187 +    for (i = 0; i < nProc; i++)
188 +      potatoes[i] = 0;
189 +    
190      outFile << mpiSim->getTotAtoms() << "\n";
191  
192      outFile << currentTime << ";\t"
# Line 188 | Line 205 | void DumpWriter::writeDump( double currentTime ){
205      outFile << entry_plug->the_integrator->getAdditionalParameters();
206      outFile << endl;
207      outFile.flush();
208 +
209      for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
210 +      
211        // Get the Node number which has this atom;
212 <
212 >      
213        which_node = AtomToProcMap[i];
214 <
214 >      
215        if (which_node != 0) {
197        
198        atomTypeTag          = 4*i;
199        atomIsDirectionalTag = 4*i + 1;
200        atomTransDataTag     = 4*i + 2;
201        atomOrientDataTag    = 4*i + 3;
216  
217 +        if (potatoes[which_node] + 3 >= MAXTAG) {
218 +          // The potato was going to exceed the maximum value,
219 +          // so wrap this processor potato back to 0:        
220 +
221 +          potatoes[which_node] = 0;          
222 +          MPI_Send(0, 1, MPI_INT, which_node, 0, MPI_COMM_WORLD);
223 +          
224 +        }
225 +
226 +        myPotato = potatoes[which_node];        
227 +        
228          MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
229 <                 atomTypeTag, MPI_COMM_WORLD, &istatus);
229 >                 myPotato, MPI_COMM_WORLD, &istatus);
230          
231 <        strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
231 >        //strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
232          
233          // Null terminate the atomTypeString just in case:
234  
235 <        atomTypeString[strlen(atomTypeString) - 1] = '\0';
235 >        //atomTypeString[strlen(atomTypeString) - 1] = '\0';
236 >        atomTypeString = MPIatomTypeString;
237 >        
238 >        myPotato++;
239  
240          MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
241 <                 atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
242 <        
243 <        MPI_Recv(atomTransData, 6, MPI_DOUBLE, which_node,
216 <                 atomTransDataTag, MPI_COMM_WORLD, &istatus);
241 >                 myPotato, MPI_COMM_WORLD, &istatus);
242 >              
243 >        myPotato++;
244  
245 <        if (isDirectional) {
246 <
247 <          MPI_Recv(atomOrientData, 7, MPI_DOUBLE, which_node,
248 <                   atomOrientDataTag, MPI_COMM_WORLD, &istatus);
249 <
245 >        if (isDirectional) {          
246 >          MPI_Recv(atomData13, 13, MPI_DOUBLE, which_node,
247 >                   myPotato, MPI_COMM_WORLD, &istatus);
248 >        } else {
249 >          MPI_Recv(atomData6, 6, MPI_DOUBLE, which_node,
250 >                   myPotato, MPI_COMM_WORLD, &istatus);          
251          }
252 +        
253 +        myPotato++;
254 +        potatoes[which_node] = myPotato;
255  
256        } else {
257          
258          haveError = 0;
259          which_atom = i;
260          local_index=-1;
261 <
261 >        
262          for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
263            if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
264          }
265 <
265 >        
266          if (local_index != -1) {
267 <
267 >          
268            atomTypeString = atoms[local_index]->getType();
269  
270            atoms[local_index]->getPos(pos);
271 <          atoms[local_index]->getVel(vel);
271 >          atoms[local_index]->getVel(vel);          
272  
273 <          atomTransData[0] = pos[0];
274 <          atomTransData[1] = pos[1];
275 <          atomTransData[2] = pos[2];
273 >          atomData6[0] = pos[0];
274 >          atomData6[1] = pos[1];
275 >          atomData6[2] = pos[2];
276  
277 <          atomTransData[3] = vel[0];
278 <          atomTransData[4] = vel[1];
279 <          atomTransData[5] = vel[2];
277 >          atomData6[3] = vel[0];
278 >          atomData6[4] = vel[1];
279 >          atomData6[5] = vel[2];
280            
281            isDirectional = 0;
282  
# Line 255 | Line 286 | void DumpWriter::writeDump( double currentTime ){
286              
287              dAtom = (DirectionalAtom *)atoms[local_index];
288              dAtom->getQ( q );
258            
259            atomOrientData[0] = q[0];
260            atomOrientData[1] = q[1];
261            atomOrientData[2] = q[2];
262            atomOrientData[3] = q[3];
289  
290 <            atomOrientData[4] = dAtom->getJx();
291 <            atomOrientData[5] = dAtom->getJy();
292 <            atomOrientData[6] = dAtom->getJz();
290 >            for (int j = 0; j < 6 ; j++)
291 >              atomData13[j] = atomData6[j];            
292 >            
293 >            atomData13[6] = q[0];
294 >            atomData13[7] = q[1];
295 >            atomData13[8] = q[2];
296 >            atomData13[9] = q[3];
297 >            
298 >            atomData13[10] = dAtom->getJx();
299 >            atomData13[11] = dAtom->getJy();
300 >            atomData13[12] = dAtom->getJz();
301            }
302 <
302 >          
303          } else {
304            sprintf(painCave.errMsg,
305                    "Atom %d not found on processor %d\n",
# Line 273 | Line 307 | void DumpWriter::writeDump( double currentTime ){
307            haveError= 1;
308            simError();
309          }
310 <
310 >        
311          if(haveError) DieDieDie();
278                              
279        // If we've survived to here, format the line:
312          
281        sprintf( tempBuffer,
282                 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
283                 atomTypeString,
284                 atomTransData[0],
285                 atomTransData[1],
286                 atomTransData[2],
287                 atomTransData[3],
288                 atomTransData[4],
289                 atomTransData[5]);
290
291        strcpy( writeLine, tempBuffer );
292
293        if (isDirectional) {
294
295          sprintf( tempBuffer,
296                   "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
297                   atomOrientData[0],
298                   atomOrientData[1],
299                   atomOrientData[2],
300                   atomOrientData[3],
301                   atomOrientData[4],
302                   atomOrientData[5],
303                   atomOrientData[6]);
304          strcat( writeLine, tempBuffer );
305
306        } else {
307          strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
308        }
309
310        outFile << writeLine;
311        outFile.flush();
313        }
314 +      // If we've survived to here, format the line:
315 +      
316 +      if (!isDirectional) {
317 +        
318 +        sprintf( tempBuffer,
319 +                 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
320 +                 atomTypeString,
321 +                 atomData6[0],
322 +                 atomData6[1],
323 +                 atomData6[2],
324 +                 atomData6[3],
325 +                 atomData6[4],
326 +                 atomData6[5]);
327 +        
328 +        strcpy( writeLine, tempBuffer );
329 +        strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
330 +        
331 +      } else {
332 +        
333 +        sprintf( tempBuffer,
334 +                 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
335 +                 atomTypeString,
336 +                 atomData13[0],
337 +                 atomData13[1],
338 +                 atomData13[2],
339 +                 atomData13[3],
340 +                 atomData13[4],
341 +                 atomData13[5],
342 +                 atomData13[6],
343 +                 atomData13[7],
344 +                 atomData13[8],
345 +                 atomData13[9],
346 +                 atomData13[10],
347 +                 atomData13[11],
348 +                 atomData13[12]);
349 +        
350 +        strcpy( writeLine, tempBuffer );
351 +        
352 +      }
353 +      
354 +      outFile << writeLine;
355 +      outFile.flush();
356      }
357 +    
358  
359      outFile.flush();
360      sprintf( checkPointMsg,
361               "Sucessfully took a dump.\n");
362      MPIcheckPoint();        
363 <    
363 >    delete[] potatoes;
364    } else {
365  
366      // worldRank != 0, so I'm a remote node.  
367 +
368 +    // Set my magic potato to 0:
369 +
370 +    myPotato = 0;
371      
372      for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
373        
# Line 327 | Line 375 | void DumpWriter::writeDump( double currentTime ){
375        
376        if (AtomToProcMap[i] == worldRank) {
377  
378 <        local_index=-1;
378 >        if (myPotato + 3 >= MAXTAG) {
379 >
380 >          // The potato was going to exceed the maximum value,
381 >          // so wrap this processor potato back to 0 (and block until
382 >          // node 0 says we can go:
383 >
384 >          MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, &istatus);
385 >          
386 >        }
387 >        which_atom = i;
388 >        local_index=-1;
389          for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
390            if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
391          }
# Line 338 | Line 396 | void DumpWriter::writeDump( double currentTime ){
396            atoms[local_index]->getPos(pos);
397            atoms[local_index]->getVel(vel);
398  
399 <          atomTransData[0] = pos[0];
400 <          atomTransData[1] = pos[1];
401 <          atomTransData[2] = pos[2];
399 >          atomData6[0] = pos[0];
400 >          atomData6[1] = pos[1];
401 >          atomData6[2] = pos[2];
402  
403 <          atomTransData[3] = vel[0];
404 <          atomTransData[4] = vel[1];
405 <          atomTransData[5] = vel[2];
403 >          atomData6[3] = vel[0];
404 >          atomData6[4] = vel[1];
405 >          atomData6[5] = vel[2];
406            
407            isDirectional = 0;
408  
# Line 355 | Line 413 | void DumpWriter::writeDump( double currentTime ){
413              dAtom = (DirectionalAtom *)atoms[local_index];
414              dAtom->getQ( q );
415              
416 <            atomOrientData[0] = q[0];
417 <            atomOrientData[1] = q[1];
418 <            atomOrientData[2] = q[2];
419 <            atomOrientData[3] = q[3];
416 >            for (int j = 0; j < 6 ; j++)
417 >              atomData13[j] = atomData6[j];
418 >            
419 >            atomData13[6] = q[0];
420 >            atomData13[7] = q[1];
421 >            atomData13[8] = q[2];
422 >            atomData13[9] = q[3];
423  
424 <            atomOrientData[4] = dAtom->getJx();
425 <            atomOrientData[5] = dAtom->getJy();
426 <            atomOrientData[6] = dAtom->getJz();
424 >            atomData13[10] = dAtom->getJx();
425 >            atomData13[11] = dAtom->getJy();
426 >            atomData13[12] = dAtom->getJz();
427            }
428  
429          } else {
# Line 373 | Line 434 | void DumpWriter::writeDump( double currentTime ){
434            simError();
435          }
436  
376        // I've survived this far, so send off the data!
377
378        atomTypeTag          = 4*i;
379        atomIsDirectionalTag = 4*i + 1;
380        atomTransDataTag     = 4*i + 2;
381        atomOrientDataTag    = 4*i + 3;
382
383
437          strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
438  
439          // null terminate the string before sending (just in case):
440          MPIatomTypeString[MINIBUFFERSIZE-1] = '\0';
441  
442          MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
443 <                 atomTypeTag, MPI_COMM_WORLD);
443 >                 myPotato, MPI_COMM_WORLD);
444          
445 +        myPotato++;
446 +
447          MPI_Send(&isDirectional, 1, MPI_INT, 0,
448 <                 atomIsDirectionalTag, MPI_COMM_WORLD);
448 >                 myPotato, MPI_COMM_WORLD);
449          
450 <        MPI_Send(atomTransData, 6, MPI_DOUBLE, 0,
451 <                 atomTransDataTag, MPI_COMM_WORLD);
397 <
450 >        myPotato++;
451 >        
452          if (isDirectional) {
453  
454 <          MPI_Send(atomOrientData, 7, MPI_DOUBLE, 0,
455 <                   atomOrientDataTag, MPI_COMM_WORLD);
454 >          MPI_Send(atomData13, 13, MPI_DOUBLE, 0,
455 >                   myPotato, MPI_COMM_WORLD);
456            
457 +        } else {
458 +
459 +          MPI_Send(atomData6, 6, MPI_DOUBLE, 0,
460 +                   myPotato, MPI_COMM_WORLD);
461          }
462 <      
462 >
463 >        myPotato++;      
464        }
465      }
466  
# Line 409 | Line 468 | void DumpWriter::writeDump( double currentTime ){
468               "Sucessfully took a dump.\n");
469      MPIcheckPoint();        
470      
471 <  }
471 >  }
472    
414  painCave.isEventLoop = 0;
415
473   #endif // is_mpi
474   }
475  
# Line 431 | Line 488 | void DumpWriter::writeFinal(double finalTime){
488    Atom** atoms = entry_plug->atoms;
489    int i;
490   #ifdef IS_MPI
491 +  
492 +  int *potatoes;
493 +  int myPotato;
494 +
495 +  int nProc;
496    int j, which_node, done, which_atom, local_index;
497 <  double atomTransData[6];
498 <  double atomOrientData[7];
497 >  double atomData6[6];
498 >  double atomData13[13];
499    int isDirectional;
500    char* atomTypeString;
501    char MPIatomTypeString[MINIBUFFERSIZE];
502 <  int atomTypeTag;
441 <  int atomIsDirectionalTag;
442 <  int atomTransDataTag;
443 <  int atomOrientDataTag;
502 >
503   #else //is_mpi
504    int nAtoms = entry_plug->n_atoms;
505   #endif //is_mpi
# Line 536 | Line 595 | void DumpWriter::writeFinal(double finalTime){
595  
596   #else // is_mpi
597  
598 <  // first thing first, suspend fatalities.
599 <  painCave.isEventLoop = 1;
598 >  /* code to find maximum tag value */
599 >  int *tagub, flag, MAXTAG;
600 >  MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, &tagub, &flag);
601 >  if (flag) {
602 >    MAXTAG = *tagub;
603 >  } else {
604 >    MAXTAG = 32767;
605 >  }  
606  
542  int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
607    int haveError;
608  
609    MPI_Status istatus;
# Line 548 | Line 612 | void DumpWriter::writeFinal(double finalTime){
612    // write out header and node 0's coordinates
613  
614    if( worldRank == 0 ){
615 +
616 +    // Node 0 needs a list of the magic potatoes for each processor;
617 +
618 +    nProc = mpiSim->getNumberProcessors();
619 +    potatoes = new int[nProc];
620 +
621 +    for (i = 0; i < nProc; i++)
622 +      potatoes[i] = 0;
623 +    
624      finalOut << mpiSim->getTotAtoms() << "\n";
625  
626      finalOut << finalTime << ";\t"
# Line 566 | Line 639 | void DumpWriter::writeFinal(double finalTime){
639      finalOut << entry_plug->the_integrator->getAdditionalParameters();
640      finalOut << endl;
641      finalOut.flush();
642 +
643      for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
644 +      
645        // Get the Node number which has this atom;
646 <
646 >      
647        which_node = AtomToProcMap[i];
648 <
648 >      
649        if (which_node != 0) {
575        
576        atomTypeTag          = 4*i;
577        atomIsDirectionalTag = 4*i + 1;
578        atomTransDataTag     = 4*i + 2;
579        atomOrientDataTag    = 4*i + 3;
650  
651 +        if (potatoes[which_node] + 3 >= MAXTAG) {
652 +          // The potato was going to exceed the maximum value,
653 +          // so wrap this processor potato back to 0:        
654 +
655 +          potatoes[which_node] = 0;          
656 +          MPI_Send(0, 1, MPI_INT, which_node, 0, MPI_COMM_WORLD);
657 +          
658 +        }
659 +
660 +        myPotato = potatoes[which_node];        
661 +        
662          MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
663 <                 atomTypeTag, MPI_COMM_WORLD, &istatus);
663 >                 myPotato, MPI_COMM_WORLD, &istatus);
664          
665 <        strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
665 >        atomTypeString = MPIatomTypeString;
666 >                
667 >        myPotato++;
668  
669          MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
670 <                 atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
671 <        
672 <        MPI_Recv(atomTransData, 6, MPI_DOUBLE, which_node,
590 <                 atomTransDataTag, MPI_COMM_WORLD, &istatus);
670 >                 myPotato, MPI_COMM_WORLD, &istatus);
671 >              
672 >        myPotato++;
673  
674 <        if (isDirectional) {
675 <
676 <          MPI_Recv(atomOrientData, 7, MPI_DOUBLE, which_node,
677 <                   atomOrientDataTag, MPI_COMM_WORLD, &istatus);
678 <
674 >        if (isDirectional) {          
675 >          MPI_Recv(atomData13, 13, MPI_DOUBLE, which_node,
676 >                   myPotato, MPI_COMM_WORLD, &istatus);
677 >        } else {
678 >          MPI_Recv(atomData6, 6, MPI_DOUBLE, which_node,
679 >                   myPotato, MPI_COMM_WORLD, &istatus);          
680          }
681 +        
682 +        myPotato++;
683 +        potatoes[which_node] = myPotato;
684  
685        } else {
686          
687          haveError = 0;
688          which_atom = i;
689          local_index=-1;
690 <
690 >        
691          for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
692            if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
693          }
694 <
694 >        
695          if (local_index != -1) {
696 <
696 >          
697            atomTypeString = atoms[local_index]->getType();
698  
699            atoms[local_index]->getPos(pos);
700 <          atoms[local_index]->getVel(vel);
700 >          atoms[local_index]->getVel(vel);          
701  
702 <          atomTransData[0] = pos[0];
703 <          atomTransData[1] = pos[1];
704 <          atomTransData[2] = pos[2];
702 >          atomData6[0] = pos[0];
703 >          atomData6[1] = pos[1];
704 >          atomData6[2] = pos[2];
705  
706 <          atomTransData[3] = vel[0];
707 <          atomTransData[4] = vel[1];
708 <          atomTransData[5] = vel[2];
706 >          atomData6[3] = vel[0];
707 >          atomData6[4] = vel[1];
708 >          atomData6[5] = vel[2];
709            
710            isDirectional = 0;
711  
# Line 629 | Line 715 | void DumpWriter::writeFinal(double finalTime){
715              
716              dAtom = (DirectionalAtom *)atoms[local_index];
717              dAtom->getQ( q );
632            
633            atomOrientData[0] = q[0];
634            atomOrientData[1] = q[1];
635            atomOrientData[2] = q[2];
636            atomOrientData[3] = q[3];
718  
719 <            atomOrientData[4] = dAtom->getJx();
720 <            atomOrientData[5] = dAtom->getJy();
721 <            atomOrientData[6] = dAtom->getJz();
719 >            for (int j = 0; j < 6 ; j++)
720 >              atomData13[j] = atomData6[j];            
721 >            
722 >            atomData13[6] = q[0];
723 >            atomData13[7] = q[1];
724 >            atomData13[8] = q[2];
725 >            atomData13[9] = q[3];
726 >            
727 >            atomData13[10] = dAtom->getJx();
728 >            atomData13[11] = dAtom->getJy();
729 >            atomData13[12] = dAtom->getJz();
730            }
731 <
731 >          
732          } else {
733            sprintf(painCave.errMsg,
734                    "Atom %d not found on processor %d\n",
# Line 647 | Line 736 | void DumpWriter::writeFinal(double finalTime){
736            haveError= 1;
737            simError();
738          }
739 <
739 >        
740          if(haveError) DieDieDie();
652                              
653        // If we've survived to here, format the line:
741          
742 <        sprintf( tempBuffer,
656 <                 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
657 <                 atomTypeString,
658 <                 atomTransData[0],
659 <                 atomTransData[1],
660 <                 atomTransData[2],
661 <                 atomTransData[3],
662 <                 atomTransData[4],
663 <                 atomTransData[5]);
742 >      }
743  
665        strcpy( writeLine, tempBuffer );
744  
745 <        if (isDirectional) {
746 <
747 <          sprintf( tempBuffer,
748 <                   "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
749 <                   atomOrientData[0],
750 <                   atomOrientData[1],
751 <                   atomOrientData[2],
752 <                   atomOrientData[3],
753 <                   atomOrientData[4],
754 <                   atomOrientData[5],
755 <                   atomOrientData[6]);
756 <          strcat( writeLine, tempBuffer );
757 <
758 <        } else {
759 <          strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
760 <        }
761 <
762 <        finalOut << writeLine;
763 <        finalOut.flush();
745 >      // If we've survived to here, format the line:
746 >      
747 >      if (!isDirectional) {
748 >        
749 >        sprintf( tempBuffer,
750 >                 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
751 >                 atomTypeString,
752 >                 atomData6[0],
753 >                 atomData6[1],
754 >                 atomData6[2],
755 >                 atomData6[3],
756 >                 atomData6[4],
757 >                 atomData6[5]);
758 >        
759 >        strcpy( writeLine, tempBuffer );
760 >        strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
761 >        
762 >      } else {
763 >        
764 >        sprintf( tempBuffer,
765 >                 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
766 >                 atomTypeString,
767 >                 atomData13[0],
768 >                 atomData13[1],
769 >                 atomData13[2],
770 >                 atomData13[3],
771 >                 atomData13[4],
772 >                 atomData13[5],
773 >                 atomData13[6],
774 >                 atomData13[7],
775 >                 atomData13[8],
776 >                 atomData13[9],
777 >                 atomData13[10],
778 >                 atomData13[11],
779 >                 atomData13[12]);
780 >        
781 >        strcpy( writeLine, tempBuffer );
782 >        
783        }
784 +        
785 +      finalOut << writeLine;
786 +      finalOut.flush();
787      }
788 <
788 >  
789      finalOut.flush();
790      sprintf( checkPointMsg,
791               "Sucessfully took a dump.\n");
792 +    delete[] potatoes;
793 +    
794      MPIcheckPoint();        
795      
796    } else {
797  
798      // worldRank != 0, so I'm a remote node.  
799 +
800 +    // Set my magic potato to 0:
801 +
802 +    myPotato = 0;
803      
804      for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
805        
# Line 701 | Line 807 | void DumpWriter::writeFinal(double finalTime){
807        
808        if (AtomToProcMap[i] == worldRank) {
809  
810 +        if (myPotato + 3 >= MAXTAG) {
811 +
812 +          // The potato was going to exceed the maximum value,
813 +          // so wrap this processor potato back to 0 (and block until
814 +          // node 0 says we can go:
815 +
816 +          MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, &istatus);
817 +          
818 +        }
819 +        which_atom = i;  
820          local_index=-1;
821          for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
822            if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
# Line 712 | Line 828 | void DumpWriter::writeFinal(double finalTime){
828            atoms[local_index]->getPos(pos);
829            atoms[local_index]->getVel(vel);
830  
831 <          atomTransData[0] = pos[0];
832 <          atomTransData[1] = pos[1];
833 <          atomTransData[2] = pos[2];
831 >          atomData6[0] = pos[0];
832 >          atomData6[1] = pos[1];
833 >          atomData6[2] = pos[2];
834  
835 <          atomTransData[3] = vel[0];
836 <          atomTransData[4] = vel[1];
837 <          atomTransData[5] = vel[2];
835 >          atomData6[3] = vel[0];
836 >          atomData6[4] = vel[1];
837 >          atomData6[5] = vel[2];
838            
839            isDirectional = 0;
840  
# Line 729 | Line 845 | void DumpWriter::writeFinal(double finalTime){
845              dAtom = (DirectionalAtom *)atoms[local_index];
846              dAtom->getQ( q );
847              
848 <            atomOrientData[0] = q[0];
849 <            atomOrientData[1] = q[1];
850 <            atomOrientData[2] = q[2];
851 <            atomOrientData[3] = q[3];
848 >            for (int j = 0; j < 6 ; j++)
849 >              atomData13[j] = atomData6[j];
850 >            
851 >            atomData13[6] = q[0];
852 >            atomData13[7] = q[1];
853 >            atomData13[8] = q[2];
854 >            atomData13[9] = q[3];
855  
856 <            atomOrientData[4] = dAtom->getJx();
857 <            atomOrientData[5] = dAtom->getJy();
858 <            atomOrientData[6] = dAtom->getJz();
856 >            atomData13[10] = dAtom->getJx();
857 >            atomData13[11] = dAtom->getJy();
858 >            atomData13[12] = dAtom->getJz();
859            }
860  
861          } else {
# Line 747 | Line 866 | void DumpWriter::writeFinal(double finalTime){
866            simError();
867          }
868  
750        // I've survived this far, so send off the data!
751
752        atomTypeTag          = 4*i;
753        atomIsDirectionalTag = 4*i + 1;
754        atomTransDataTag     = 4*i + 2;
755        atomOrientDataTag    = 4*i + 3;
756
869          strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
870  
871 +        // null terminate the string before sending (just in case):
872 +        MPIatomTypeString[MINIBUFFERSIZE-1] = '\0';
873 +
874          MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
875 <                 atomTypeTag, MPI_COMM_WORLD);
875 >                 myPotato, MPI_COMM_WORLD);
876          
877 +        myPotato++;
878 +
879          MPI_Send(&isDirectional, 1, MPI_INT, 0,
880 <                 atomIsDirectionalTag, MPI_COMM_WORLD);
880 >                 myPotato, MPI_COMM_WORLD);
881          
882 <        MPI_Send(atomTransData, 6, MPI_DOUBLE, 0,
883 <                 atomTransDataTag, MPI_COMM_WORLD);
767 <
882 >        myPotato++;
883 >        
884          if (isDirectional) {
885  
886 <          MPI_Send(atomOrientData, 7, MPI_DOUBLE, 0,
887 <                   atomOrientDataTag, MPI_COMM_WORLD);
886 >          MPI_Send(atomData13, 13, MPI_DOUBLE, 0,
887 >                   myPotato, MPI_COMM_WORLD);
888            
889 +        } else {
890 +
891 +          MPI_Send(atomData6, 6, MPI_DOUBLE, 0,
892 +                   myPotato, MPI_COMM_WORLD);
893          }
894 <      
894 >
895 >        myPotato++;      
896        }
897      }
898  
899      sprintf( checkPointMsg,
900 <             "Sucessfully wrote final file.\n");
900 >             "Sucessfully took a dump.\n");
901      MPIcheckPoint();        
902      
903 <  }
903 >  }
904    
784  painCave.isEventLoop = 0;
785
905    if( worldRank == 0 ) finalOut.close();
906   #endif // is_mpi
907   }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines