--- trunk/OOPSE/libmdtools/InitializeFromFile.cpp 2004/04/22 14:55:17 1130 +++ trunk/OOPSE/libmdtools/InitializeFromFile.cpp 2004/05/27 18:59:17 1203 @@ -220,7 +220,6 @@ void InitializeFromFile :: readInit( SimInfo* the_simn MPI_Bcast(read_buffer, BUFFERSIZE, MPI_CHAR, 0, MPI_COMM_WORLD); - cerr << "node " << worldRank << " finished MPI_Bcast" << endl; parseErr = parseCommentLine( read_buffer, simnfo); @@ -230,7 +229,7 @@ void InitializeFromFile :: readInit( SimInfo* the_simn simError(); } - for (i=0 ; i < mpiSim->getTotNmol(); i++) { + for (i=0 ; i < mpiSim->getNMolGlobal(); i++) { which_node = MolToProcMap[i]; if(which_node == 0){ //molecules belong to master node @@ -270,7 +269,7 @@ void InitializeFromFile :: readInit( SimInfo* the_simn MPI_Recv(&nCurObj, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus); - cerr << "node " << worldRank << " finished MPI_Send" << endl; + for(j=0; j < nCurObj; j++){ eof_test = fgets(read_buffer, sizeof(read_buffer), c_in_file); @@ -288,7 +287,7 @@ void InitializeFromFile :: readInit( SimInfo* the_simn MPI_Send(read_buffer, BUFFERSIZE, MPI_CHAR, which_node, TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD); - cerr << "node " << worldRank << " finished MPI_Send" << endl; + } } @@ -301,7 +300,6 @@ void InitializeFromFile :: readInit( SimInfo* the_simn MPI_Bcast(read_buffer, BUFFERSIZE, MPI_CHAR, 0, MPI_COMM_WORLD); - cerr << "node " << worldRank << " finished MPI_Bcast" << endl; parseErr = parseCommentLine( read_buffer, simnfo); if( parseErr != NULL ){ @@ -310,7 +308,7 @@ void InitializeFromFile :: readInit( SimInfo* the_simn simError(); } - for (i=0 ; i < mpiSim->getTotNmol(); i++) { + for (i=0 ; i < mpiSim->getNMolGlobal(); i++) { which_node = MolToProcMap[i]; if(which_node == worldRank){ @@ -331,13 +329,11 @@ void InitializeFromFile :: readInit( SimInfo* the_simn MPI_Send(&nCurObj, 1, MPI_INT, 0, TAKE_THIS_TAG_INT, MPI_COMM_WORLD); - cerr << "node " << worldRank << " finished MPI_Send" << endl; for(j = 0; j < integrableObjects.size(); j++){ MPI_Recv(read_buffer, BUFFERSIZE, MPI_CHAR, 0, TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD, &istatus); - cerr << "node " << worldRank << " finished MPI_Recv" << endl; parseErr = parseDumpLine(read_buffer, integrableObjects[j]); if( parseErr != NULL ){ @@ -653,7 +649,7 @@ void initFile::nodeZeroError( void ){ int j, myStatus; myStatus = 0; - for (j = 0; j < mpiSim->getNumberProcessors(); j++) { + for (j = 0; j < mpiSim->getNProcessors(); j++) { MPI_Send( &myStatus, 1, MPI_INT, j, TAKE_THIS_TAG_INT, MPI_COMM_WORLD); }