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 907 by gezelter, Thu Jan 8 17:40:56 2004 UTC vs.
Revision 910 by gezelter, Thu Jan 8 18:05:37 2004 UTC

# Line 78 | Line 78 | void DumpWriter::writeDump( double currentTime ){
78    double atomOrientData[7];
79    int isDirectional;
80    char* atomTypeString;
81 +  char MPIatomTypeString[MINIBUFFERSIZE];
82    int me;
83    int atomTypeTag;
84    int atomIsDirectionalTag;
# Line 199 | Line 200 | void DumpWriter::writeDump( double currentTime ){
200          atomTransDataTag     = 4*i + 2;
201          atomOrientDataTag    = 4*i + 3;
202  
203 <        MPI_Recv(atomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
203 >        MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
204                   atomTypeTag, MPI_COMM_WORLD, &istatus);
205          
206 +        strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
207 +
208          MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
209                   atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
210          
# Line 373 | Line 376 | void DumpWriter::writeDump( double currentTime ){
376          atomTransDataTag     = 4*i + 2;
377          atomOrientDataTag    = 4*i + 3;
378  
379 <        MPI_Send(atomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
379 >
380 >        strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
381 >
382 >        MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
383                   atomTypeTag, MPI_COMM_WORLD);
384          
385          MPI_Send(&isDirectional, 1, MPI_INT, 0,
# Line 423 | Line 429 | void DumpWriter::writeFinal(double finalTime){
429    double atomOrientData[7];
430    int isDirectional;
431    char* atomTypeString;
432 +  char MPIatomTypeString[MINIBUFFERSIZE];
433    int atomTypeTag;
434    int atomIsDirectionalTag;
435    int atomTransDataTag;
# Line 564 | Line 571 | void DumpWriter::writeFinal(double finalTime){
571          atomTransDataTag     = 4*i + 2;
572          atomOrientDataTag    = 4*i + 3;
573  
574 <        MPI_Recv(atomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
574 >        MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
575                   atomTypeTag, MPI_COMM_WORLD, &istatus);
576          
577 +        strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
578 +
579          MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
580                   atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
581          
# Line 738 | Line 747 | void DumpWriter::writeFinal(double finalTime){
747          atomTransDataTag     = 4*i + 2;
748          atomOrientDataTag    = 4*i + 3;
749  
750 <        MPI_Send(atomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
750 >        strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
751 >
752 >        MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
753                   atomTypeTag, MPI_COMM_WORLD);
754          
755          MPI_Send(&isDirectional, 1, MPI_INT, 0,

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines