]> git.donarmstrong.com Git - mothur.git/blobdiff - chimeracheckcommand.cpp
added cluster.split command
[mothur.git] / chimeracheckcommand.cpp
index d80cce65cfe21703d9b57f94133040685b2fca76..66cf12fb8768bcc9c4995256544761e2d8baf820 100644 (file)
@@ -166,29 +166,41 @@ int ChimeraCheckCommand::execute(){
                                                
                        int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
                        int inMode=MPI_MODE_RDONLY; 
-                                                       
-                       char outFilename[outputFileName.length()];
+                       
+                       //char* outFilename = new char[outputFileName.length()];
+                       //memcpy(outFilename, outputFileName.c_str(), outputFileName.length());
+                       
+                       char outFilename[1024];
                        strcpy(outFilename, outputFileName.c_str());
+
+                       //char* inFileName = new char[fastafile.length()];
+                       //memcpy(inFileName, fastafile.c_str(), fastafile.length());
                        
-                       char inFileName[fastafile.length()];
+                       char inFileName[1024];
                        strcpy(inFileName, fastafile.c_str());
 
                        MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
                        MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
                        
+                       //delete outFilename;
+                       //delete inFileName;
+
                        if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  delete chimera; return 0;  }
                        
                        if (pid == 0) { //you are the root process 
                                MPIPos = setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
                                
                                //send file positions to all processes
-                               MPI_Bcast(&numSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD);  //send numSeqs
-                               MPI_Bcast(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //send file pos        
+                               for(int i = 1; i < processors; i++) { 
+                                       MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
+                                       MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
+                               }       
                                
                                //figure out how many sequences you have to align
                                numSeqsPerProcessor = numSeqs / processors;
-                               if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
                                int startIndex =  pid * numSeqsPerProcessor;
+                               if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
+                               
                        
                                //align your part
                                driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
@@ -201,14 +213,15 @@ int ChimeraCheckCommand::execute(){
                                        MPI_Recv(buf, 4, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); 
                                }
                        }else{ //you are a child process
-                               MPI_Bcast(&numSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD); //get numSeqs
+                               MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
                                MPIPos.resize(numSeqs+1);
-                               MPI_Bcast(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //get file positions
+                               MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
                                
                                //figure out how many sequences you have to align
                                numSeqsPerProcessor = numSeqs / processors;
-                               if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
                                int startIndex =  pid * numSeqsPerProcessor;
+                               if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
+                               
                                
                                //align your part
                                driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
@@ -224,6 +237,7 @@ int ChimeraCheckCommand::execute(){
                        //close files 
                        MPI_File_close(&inMPI);
                        MPI_File_close(&outMPI);
+                       MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
        #else
                
                //break up file
@@ -391,12 +405,13 @@ int ChimeraCheckCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File
                        //read next sequence
                        int length = MPIPos[start+i+1] - MPIPos[start+i];
        
-                       char buf4[length];
+                       char* buf4 = new char[length];
                        MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
                        
                        string tempBuf = buf4;
                        if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
                        istringstream iss (tempBuf,istringstream::in);
+                       delete buf4;
 
                        Sequence* candidateSeq = new Sequence(iss);  gobble(iss);