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

Comparing trunk/OOPSE/libmdtools/InitializeFromFile.cpp (file contents):
Revision 1113 by tim, Thu Apr 15 16:18:26 2004 UTC vs.
Revision 1198 by tim, Thu May 27 00:48:12 2004 UTC

# Line 16 | Line 16
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
19 > #define TAKE_THIS_TAG_CHAR 3134
20 > #define TAKE_THIS_TAG_INT 3135
21  
22   namespace initFile{
23    void nodeZeroError( void );
# Line 175 | Line 175 | void InitializeFromFile :: readInit( SimInfo* the_simn
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  
# Line 187 | Line 189 | void InitializeFromFile :: readInit( SimInfo* the_simn
189        simError();
190      }
191  
192 <    nTotObjs = atoi( read_buffer );
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 != simnfo->getTotIntegrableObjects()){
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",
# Line 217 | Line 219 | void InitializeFromFile :: readInit( SimInfo* the_simn
219      //parseCommentLine
220  
221      MPI_Bcast(read_buffer, BUFFERSIZE, MPI_CHAR, 0, MPI_COMM_WORLD);
222 +
223  
224      parseErr = parseCommentLine( read_buffer, simnfo);
225  
# Line 255 | Line 258 | void InitializeFromFile :: readInit( SimInfo* the_simn
258            
259            if(haveError) nodeZeroError();
260  
261 <          parseDumpLine(read_buffer, integrableObjects[i]);
261 >          parseDumpLine(read_buffer, integrableObjects[j]);
262            
263         }
264  
# Line 264 | Line 267 | void InitializeFromFile :: readInit( SimInfo* the_simn
267        else{
268        //molecule belongs to slave nodes
269  
270 <        MPI_Recv(&nCurObj, 1, MPI_INT, 0,
270 >        MPI_Recv(&nCurObj, 1, MPI_INT, which_node,
271                 TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus);
272 <      
273 <       for(j=0; j < integrableObjects.size(); j++){
272 >
273 >       for(j=0; j < nCurObj; j++){
274          
275            eof_test = fgets(read_buffer, sizeof(read_buffer), c_in_file);
276            if(eof_test == NULL){
# Line 284 | Line 287 | void InitializeFromFile :: readInit( SimInfo* the_simn
287  
288              MPI_Send(read_buffer, BUFFERSIZE, MPI_CHAR, which_node,
289                        TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD);
290 <          
290 >
291         }
292  
293        }
# Line 294 | Line 297 | void InitializeFromFile :: readInit( SimInfo* the_simn
297    }
298    else{
299    //actions taken at slave nodes
300 +
301 +    MPI_Bcast(read_buffer, BUFFERSIZE, MPI_CHAR, 0, MPI_COMM_WORLD);
302 +
303 +    parseErr = parseCommentLine( read_buffer, simnfo);
304 +
305 +    if( parseErr != NULL ){
306 +      strcpy( painCave.errMsg, parseErr );
307 +      haveError = 1;
308 +      simError();
309 +    }
310 +  
311      for (i=0 ; i < mpiSim->getTotNmol(); i++) {
312        which_node = MolToProcMap[i];
313        
# Line 312 | Line 326 | void InitializeFromFile :: readInit( SimInfo* the_simn
326  
327          nCurObj = integrableObjects.size();
328          
329 <        MPI_Recv(&nCurObj, 1, MPI_INT, 0,
330 <                        TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus);
329 >        MPI_Send(&nCurObj, 1, MPI_INT, 0,
330 >                        TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
331  
332          for(j = 0; j < integrableObjects.size(); j++){
333  
# Line 617 | Line 631 | char* InitializeFromFile::parseCommentLine(char* readL
631      
632      //push eta into SimInfo::properties which can be
633      //retrieved by integrator later
634 <    //entry_plug->setBoxM( theBoxMat3 );
634 >    
635      DoubleArrayData* etaValue = new DoubleArrayData();
636      etaValue->setID(ETAVALUE_ID);
637      etaValue->setData(eta, 9);
# Line 635 | Line 649 | void initFile::nodeZeroError( void ){
649    int j, myStatus;
650  
651    myStatus = 0;
652 <  for (j = 0; j < mpiSim->getNumberProcessors(); j++) {
652 >  for (j = 0; j < mpiSim->getNprocessors(); j++) {
653      MPI_Send( &myStatus, 1, MPI_INT, j,
654                TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
655    }

Diff Legend

Removed lines
+ Added lines
< Changed lines
> Changed lines