]> git.donarmstrong.com Git - mothur.git/blob - chimeraslayercommand.cpp
paralellized screen.seqs and added mpi code to it. fixed bug with all mpi commands...
[mothur.git] / chimeraslayercommand.cpp
1 /*
2  *  chimeraslayercommand.cpp
3  *  Mothur
4  *
5  *  Created by westcott on 3/31/10.
6  *  Copyright 2010 Schloss Lab. All rights reserved.
7  *
8  */
9
10 #include "chimeraslayercommand.h"
11 #include "chimeraslayer.h"
12
13
14 //***************************************************************************************************************
15
16 ChimeraSlayerCommand::ChimeraSlayerCommand(string option)  {
17         try {
18                 abort = false;
19                 
20                 //allow user to run help
21                 if(option == "help") { help(); abort = true; }
22                 
23                 else {
24                         //valid paramters for this command
25                         string Array[] =  {"fasta", "processors", "window", "template","numwanted", "ksize", "match","mismatch", 
26                         "divergence", "minsim","mincov","minbs", "minsnp","parents", "iters","outputdir","inputdir", "search","realign" };
27                         vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
28                         
29                         OptionParser parser(option);
30                         map<string,string> parameters = parser.getParameters();
31                         
32                         ValidParameters validParameter;
33                         map<string,string>::iterator it;
34                         
35                         //check to make sure all parameters are valid for command
36                         for (it = parameters.begin(); it != parameters.end(); it++) { 
37                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
38                         }
39                         
40                         //if the user changes the input directory command factory will send this info to us in the output parameter 
41                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
42                         if (inputDir == "not found"){   inputDir = "";          }
43                         else {
44                                 string path;
45                                 it = parameters.find("fasta");
46                                 //user has given a template file
47                                 if(it != parameters.end()){ 
48                                         path = hasPath(it->second);
49                                         //if the user has not given a path then, add inputdir. else leave path alone.
50                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
51                                 }
52                                 
53                                 it = parameters.find("template");
54                                 //user has given a template file
55                                 if(it != parameters.end()){ 
56                                         path = hasPath(it->second);
57                                         //if the user has not given a path then, add inputdir. else leave path alone.
58                                         if (path == "") {       parameters["template"] = inputDir + it->second;         }
59                                 }
60                         }
61
62                         
63                         //check for required parameters
64                         fastafile = validParameter.validFile(parameters, "fasta", true);
65                         if (fastafile == "not open") { abort = true; }
66                         else if (fastafile == "not found") { fastafile = ""; m->mothurOut("fasta is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true;  }     
67                         
68                         //if the user changes the output directory command factory will send this info to us in the output parameter 
69                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
70                                 outputDir = ""; 
71                                 outputDir += hasPath(fastafile); //if user entered a file with a path then preserve it  
72                         }
73
74                         templatefile = validParameter.validFile(parameters, "template", true);
75                         if (templatefile == "not open") { abort = true; }
76                         else if (templatefile == "not found") { templatefile = "";  m->mothurOut("template is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true;  }   
77                                                 
78                         string temp = validParameter.validFile(parameters, "processors", false);                if (temp == "not found") { temp = "1"; }
79                         convert(temp, processors);
80                         
81                         temp = validParameter.validFile(parameters, "ksize", false);                    if (temp == "not found") { temp = "7"; }
82                         convert(temp, ksize);
83                                                 
84                         temp = validParameter.validFile(parameters, "window", false);                   if (temp == "not found") { temp = "50"; }                       
85                         convert(temp, window);
86                         
87                         temp = validParameter.validFile(parameters, "match", false);                    if (temp == "not found") { temp = "5"; }
88                         convert(temp, match);
89                         
90                         temp = validParameter.validFile(parameters, "mismatch", false);                 if (temp == "not found") { temp = "-4"; }
91                         convert(temp, mismatch);
92                         
93                         temp = validParameter.validFile(parameters, "divergence", false);               if (temp == "not found") { temp = "1.007"; }
94                         convert(temp, divR);
95                         
96                         temp = validParameter.validFile(parameters, "minsim", false);                   if (temp == "not found") { temp = "90"; }
97                         convert(temp, minSimilarity);
98                         
99                         temp = validParameter.validFile(parameters, "mincov", false);                   if (temp == "not found") { temp = "70"; }
100                         convert(temp, minCoverage);
101                         
102                         temp = validParameter.validFile(parameters, "minbs", false);                    if (temp == "not found") { temp = "90"; }
103                         convert(temp, minBS);
104                         
105                         temp = validParameter.validFile(parameters, "minsnp", false);                   if (temp == "not found") { temp = "10"; }
106                         convert(temp, minSNP);
107
108                         temp = validParameter.validFile(parameters, "parents", false);                  if (temp == "not found") { temp = "3"; }
109                         convert(temp, parents); 
110                         
111                         temp = validParameter.validFile(parameters, "realign", false);                  if (temp == "not found") { temp = "f"; }
112                         realign = isTrue(temp); 
113                         
114                         search = validParameter.validFile(parameters, "search", false);                 if (search == "not found") { search = "distance"; }
115                         
116                         temp = validParameter.validFile(parameters, "iters", false);                    if (temp == "not found") { temp = "100"; }              
117                         convert(temp, iters); 
118                          
119                         temp = validParameter.validFile(parameters, "increment", false);                if (temp == "not found") { temp = "5"; }
120                         convert(temp, increment);
121                         
122                         temp = validParameter.validFile(parameters, "numwanted", false);                if (temp == "not found") { temp = "15"; }               
123                         convert(temp, numwanted);
124
125                         if ((search != "distance") && (search != "blast") && (search != "kmer")) { m->mothurOut(search + " is not a valid search."); m->mothurOutEndLine(); abort = true;  }
126                 }
127         }
128         catch(exception& e) {
129                 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
130                 exit(1);
131         }
132 }
133 //**********************************************************************************************************************
134
135 void ChimeraSlayerCommand::help(){
136         try {
137         
138                 m->mothurOut("The chimera.slayer command reads a fastafile and templatefile and outputs potentially chimeric sequences.\n");
139                 m->mothurOut("This command was modeled after the chimeraSlayer written by the Broad Institute.\n");
140                 m->mothurOut("The chimera.slayer command parameters are fasta, template, processors, ksize, window, match, mismatch, divergence. minsim, mincov, minbs, minsnp, parents, search, iters, increment and numwanted.\n"); //realign,
141                 m->mothurOut("The fasta parameter allows you to enter the fasta file containing your potentially chimeric sequences, and is required. \n");
142                 m->mothurOut("The template parameter allows you to enter a template file containing known non-chimeric sequences, and is required. \n");
143                 m->mothurOut("The processors parameter allows you to specify how many processors you would like to use.  The default is 1. \n");
144                 #ifdef USE_MPI
145                 m->mothurOut("When using MPI, the processors parameter is set to the number of MPI processes running. \n");
146                 #endif
147                 m->mothurOut("The window parameter allows you to specify the window size for searching for chimeras, default=50. \n");
148                 m->mothurOut("The increment parameter allows you to specify how far you move each window while finding chimeric sequences, default=5.\n");
149                 m->mothurOut("The numwanted parameter allows you to specify how many sequences you would each query sequence compared with, default=15.\n");
150                 m->mothurOut("The ksize parameter allows you to input kmersize, default is 7, used if search is kmer. \n");
151                 m->mothurOut("The match parameter allows you to reward matched bases in blast search, default is 5. \n");
152                 m->mothurOut("The parents parameter allows you to select the number of potential parents to investigate from the numwanted best matches after rating them, default is 3. \n");
153                 m->mothurOut("The mismatch parameter allows you to penalize mismatched bases in blast search, default is -4. \n");
154                 m->mothurOut("The divergence parameter allows you to set a cutoff for chimera determination, default is 1.007. \n");
155                 m->mothurOut("The iters parameter allows you to specify the number of bootstrap iters to do with the chimeraslayer method, default=100.\n");
156                 m->mothurOut("The minsim parameter allows you to specify a minimum similarity with the parent fragments, default=90. \n");
157                 m->mothurOut("The mincov parameter allows you to specify minimum coverage by closest matches found in template. Default is 70, meaning 70%. \n");
158                 m->mothurOut("The minbs parameter allows you to specify minimum bootstrap support for calling a sequence chimeric. Default is 90, meaning 90%. \n");
159                 m->mothurOut("The minsnp parameter allows you to specify percent of SNPs to sample on each side of breakpoint for computing bootstrap support (default: 10) \n");
160                 m->mothurOut("The search parameter allows you to specify search method for finding the closest parent. Choices are distance, blast, and kmer, default distance. \n");
161                 //m->mothurOut("The realign parameter allows you to realign the query to the potential parents. Choices are true or false, default false. Found to make results worse. \n");
162                 m->mothurOut("The chimera.slayer command should be in the following format: \n");
163                 m->mothurOut("chimera.slayer(fasta=yourFastaFile, template=yourTemplate, search=yourSearch) \n");
164                 m->mothurOut("Example: chimera.slayer(fasta=AD.align, template=core_set_aligned.imputed.fasta, search=kmer) \n");
165                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");     
166         }
167         catch(exception& e) {
168                 m->errorOut(e, "ChimeraSlayerCommand", "help");
169                 exit(1);
170         }
171 }
172
173 //***************************************************************************************************************
174
175 ChimeraSlayerCommand::~ChimeraSlayerCommand(){  /*      do nothing      */      }
176
177 //***************************************************************************************************************
178
179 int ChimeraSlayerCommand::execute(){
180         try{
181                 
182                 if (abort == true) { return 0; }
183                 
184                 int start = time(NULL); 
185                 
186                 chimera = new ChimeraSlayer(fastafile, templatefile, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);   
187                                                 
188                 string outputFileName = outputDir + getRootName(getSimpleName(fastafile)) + "slayer.chimeras";
189                 string accnosFileName = outputDir + getRootName(getSimpleName(fastafile))  + "slayer.accnos";
190                 bool hasAccnos = true;
191                 
192                 if (m->control_pressed) { delete chimera;       return 0;       }
193                 
194                 if (chimera->getUnaligned()) { 
195                         m->mothurOut("Your template sequences are different lengths, please correct."); m->mothurOutEndLine(); 
196                         delete chimera;
197                         return 0; 
198                 }
199                 templateSeqsLength = chimera->getLength();
200                 
201         #ifdef USE_MPI  
202                 int pid, end, numSeqsPerProcessor; 
203                         int tag = 2001;
204                         vector<long> MPIPos;
205                         MPIWroteAccnos = false;
206                         
207                         MPI_Status status; 
208                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
209                         MPI_Comm_size(MPI_COMM_WORLD, &processors); 
210
211                         MPI_File inMPI;
212                         MPI_File outMPI;
213                         MPI_File outMPIAccnos;
214                         
215                         int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
216                         int inMode=MPI_MODE_RDONLY; 
217                         
218                         //char* outFilename = new char[outputFileName.length()];
219                         //memcpy(outFilename, outputFileName.c_str(), outputFileName.length());
220                         
221                         char outFilename[1024];
222                         strcpy(outFilename, outputFileName.c_str());
223                         
224                         //char* outAccnosFilename = new char[accnosFileName.length()];
225                         //memcpy(outAccnosFilename, accnosFileName.c_str(), accnosFileName.length());
226                         
227                         char outAccnosFilename[1024];
228                         strcpy(outAccnosFilename, accnosFileName.c_str());
229
230                         //char* inFileName = new char[fastafile.length()];
231                         //memcpy(inFileName, fastafile.c_str(), fastafile.length());
232                         
233                         char inFileName[1024];
234                         strcpy(inFileName, fastafile.c_str());
235
236                         MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
237                         MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
238                         MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
239                         
240                         //delete inFileName;
241                         //delete outFilename;
242                         //delete outAccnosFilename;
243
244                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   MPI_File_close(&outMPIAccnos);  delete chimera; return 0;  }
245                 
246                         if (pid == 0) { //you are the root process 
247                                 m->mothurOutEndLine();
248                                 m->mothurOut("Only reporting sequence supported by " + toString(minBS) + "% of bootstrapped results.");
249                                 m->mothurOutEndLine();
250         
251                                 string outTemp = "Name\tLeftParent\tRightParent\tDivQLAQRB\tPerIDQLAQRB\tBootStrapA\tDivQLBQRA\tPerIDQLBQRA\tBootStrapB\tFlag\tLeftWindow\tRightWindow\n";
252                                 
253                                 //print header
254                                 int length = outTemp.length();
255                                 char* buf2 = new char[length];
256                                 memcpy(buf2, outTemp.c_str(), length);
257
258                                 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &status);
259                                 delete buf2;
260
261                                 MPIPos = setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
262                                 
263                                 //send file positions to all processes
264                                 MPI_Bcast(&numSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD);  //send numSeqs
265                                 MPI_Bcast(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //send file pos        
266                                 
267                                 //figure out how many sequences you have to align
268                                 numSeqsPerProcessor = numSeqs / processors;
269                                 int startIndex =  pid * numSeqsPerProcessor;
270                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
271                                 
272                         
273                                 //align your part
274                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
275                                 
276                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   MPI_File_close(&outMPIAccnos);  remove(outputFileName.c_str());  remove(accnosFileName.c_str());  delete chimera; return 0;  }
277
278                                 for (int i = 1; i < processors; i++) {
279                                         bool tempResult;
280                                         MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
281                                         if (tempResult != 0) { MPIWroteAccnos = true; }
282                                 }
283                         }else{ //you are a child process
284                                 MPI_Bcast(&numSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD); //get numSeqs
285                                 MPIPos.resize(numSeqs+1);
286                                 MPI_Bcast(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //get file positions
287                                 
288                                 //figure out how many sequences you have to align
289                                 numSeqsPerProcessor = numSeqs / processors;
290                                 int startIndex =  pid * numSeqsPerProcessor;
291                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
292                                 
293                                 
294                                 //align your part
295                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
296                                 
297                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   MPI_File_close(&outMPIAccnos);  delete chimera; return 0;  }
298
299                                 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD); 
300                         }
301                         
302                         //close files 
303                         MPI_File_close(&inMPI);
304                         MPI_File_close(&outMPI);
305                         MPI_File_close(&outMPIAccnos);
306                         
307                         //delete accnos file if blank
308                         if (pid == 0) {
309                                 if (!MPIWroteAccnos) { 
310                                         //MPI_Info info;
311                                         //MPI_File_delete(outAccnosFilename, info);
312                                         hasAccnos = false;      
313                                         remove(accnosFileName.c_str()); 
314                                 }
315                         }
316                 
317         #else
318                 ofstream outHeader;
319                 string tempHeader = outputDir + getRootName(getSimpleName(fastafile)) + "slayer.chimeras.tempHeader";
320                 openOutputFile(tempHeader, outHeader);
321                 
322                 chimera->printHeader(outHeader);
323                 outHeader.close();
324                 
325                 //break up file
326                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
327                         if(processors == 1){
328                                 ifstream inFASTA;
329                                 openInputFile(fastafile, inFASTA);
330                                 numSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
331                                 inFASTA.close();
332                                 
333                                 lines.push_back(new linePair(0, numSeqs));
334                                 
335                                 driver(lines[0], outputFileName, fastafile, accnosFileName);
336                                 
337                                 if (m->control_pressed) { 
338                                         remove(outputFileName.c_str()); 
339                                         remove(tempHeader.c_str()); 
340                                         remove(accnosFileName.c_str());
341                                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
342                                         delete chimera;
343                                         return 0;
344                                 }
345                                 
346                                 //delete accnos file if its blank 
347                                 if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
348                                                                 
349                         }else{
350                                 vector<int> positions;
351                                 processIDS.resize(0);
352                                 
353                                 ifstream inFASTA;
354                                 openInputFile(fastafile, inFASTA);
355                                 
356                                 string input;
357                                 while(!inFASTA.eof()){
358                                         input = getline(inFASTA);
359                                         if (input.length() != 0) {
360                                                 if(input[0] == '>'){    long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);  }
361                                         }
362                                 }
363                                 inFASTA.close();
364                                 
365                                 numSeqs = positions.size();
366                                 
367                                 int numSeqsPerProcessor = numSeqs / processors;
368                                 
369                                 for (int i = 0; i < processors; i++) {
370                                         long int startPos = positions[ i * numSeqsPerProcessor ];
371                                         if(i == processors - 1){
372                                                 numSeqsPerProcessor = numSeqs - i * numSeqsPerProcessor;
373                                         }
374                                         lines.push_back(new linePair(startPos, numSeqsPerProcessor));
375                                 }
376                                 
377                                 
378                                 createProcesses(outputFileName, fastafile, accnosFileName); 
379                         
380                                 rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
381                                         
382                                 //append output files
383                                 for(int i=1;i<processors;i++){
384                                         appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
385                                         remove((outputFileName + toString(processIDS[i]) + ".temp").c_str());
386                                 }
387                                 
388                                 vector<string> nonBlankAccnosFiles;
389                                 //delete blank accnos files generated with multiple processes
390                                 for(int i=0;i<processors;i++){  
391                                         if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
392                                                 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
393                                         }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());  }
394                                 }
395                                 
396                                 //append accnos files
397                                 if (nonBlankAccnosFiles.size() != 0) { 
398                                         rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
399                                         
400                                         for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
401                                                 appendFiles(nonBlankAccnosFiles[h], accnosFileName);
402                                                 remove(nonBlankAccnosFiles[h].c_str());
403                                         }
404                                 }else{ hasAccnos = false;  }
405                                 
406                                 if (m->control_pressed) { 
407                                         remove(outputFileName.c_str()); 
408                                         remove(accnosFileName.c_str());
409                                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
410                                         delete chimera;
411                                         return 0;
412                                 }
413
414                         }
415
416                 #else
417                         ifstream inFASTA;
418                         openInputFile(fastafile, inFASTA);
419                         numSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
420                         inFASTA.close();
421                         lines.push_back(new linePair(0, numSeqs));
422                         
423                         driver(lines[0], outputFileName, fastafile, accnosFileName);
424                         
425                         if (m->control_pressed) { 
426                                         remove(outputFileName.c_str()); 
427                                         remove(tempHeader.c_str()); 
428                                         remove(accnosFileName.c_str());
429                                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
430                                         delete chimera;
431                                         return 0;
432                         }
433                         
434                         //delete accnos file if its blank 
435                         if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
436                 #endif
437                 
438                 appendFiles(outputFileName, tempHeader);
439         
440                 remove(outputFileName.c_str());
441                 rename(tempHeader.c_str(), outputFileName.c_str());
442                 
443         #endif
444                 delete chimera;
445                 
446                 m->mothurOutEndLine();
447                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
448                 m->mothurOut(outputFileName); m->mothurOutEndLine();    
449                 if (hasAccnos) {  m->mothurOut(accnosFileName); m->mothurOutEndLine();  }
450                 m->mothurOutEndLine();
451
452                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
453                 
454                 m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine();
455                 
456                 return 0;
457                 
458         }
459         catch(exception& e) {
460                 m->errorOut(e, "ChimeraSlayerCommand", "execute");
461                 exit(1);
462         }
463 }
464 //**********************************************************************************************************************
465
466 int ChimeraSlayerCommand::driver(linePair* line, string outputFName, string filename, string accnos){
467         try {
468                 ofstream out;
469                 openOutputFile(outputFName, out);
470                 
471                 ofstream out2;
472                 openOutputFile(accnos, out2);
473                 
474                 ifstream inFASTA;
475                 openInputFile(filename, inFASTA);
476
477                 inFASTA.seekg(line->start);
478                 
479                 for(int i=0;i<line->numSeqs;i++){
480                 
481                         if (m->control_pressed) {       return 1;       }
482                 
483                         Sequence* candidateSeq = new Sequence(inFASTA);  gobble(inFASTA);
484                                 
485                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
486                                 
487                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
488                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
489                                 }else{
490                                         //find chimeras
491                                         chimera->getChimeras(candidateSeq);
492                                         
493                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
494                 
495                                         //print results
496                                         chimera->print(out, out2);
497                                 }
498                         }
499                         delete candidateSeq;
500                         
501                         //report progress
502                         if((i+1) % 100 == 0){   m->mothurOut("Processing sequence: " + toString(i+1)); m->mothurOutEndLine();           }
503                 }
504                 //report progress
505                 if((line->numSeqs) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(line->numSeqs)); m->mothurOutEndLine();         }
506                 
507                 out.close();
508                 out2.close();
509                 inFASTA.close();
510                                 
511                 return 0;
512         }
513         catch(exception& e) {
514                 m->errorOut(e, "ChimeraSlayerCommand", "driver");
515                 exit(1);
516         }
517 }
518 //**********************************************************************************************************************
519 #ifdef USE_MPI
520 int ChimeraSlayerCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, MPI_File& outAccMPI, vector<long>& MPIPos){
521         try {                           
522                 MPI_Status status; 
523                 int pid;
524                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
525                 
526                 for(int i=0;i<num;i++){
527                         
528                         if (m->control_pressed) {       return 1;       }
529                         
530                         //read next sequence
531                         int length = MPIPos[start+i+1] - MPIPos[start+i];
532
533                         char* buf4 = new char[length];
534                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
535         
536                         string tempBuf = buf4;
537                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
538                         istringstream iss (tempBuf,istringstream::in);
539
540                         delete buf4;
541
542                         Sequence* candidateSeq = new Sequence(iss);  gobble(iss);
543                 
544                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
545                                 
546                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
547                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
548                                 }else{
549                 
550                                         //find chimeras
551                                         chimera->getChimeras(candidateSeq);
552                         
553                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
554                 //cout << "about to print" << endl;
555                                         //print results
556                                         bool isChimeric = chimera->print(outMPI, outAccMPI);
557                                         if (isChimeric) { MPIWroteAccnos = true;  }
558         
559                                 }
560                         }
561                         delete candidateSeq;
562                         
563                         //report progress
564                         if((i+1) % 100 == 0){  cout << "Processing sequence: " << (i+1) << endl;        m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n");          }
565                 }
566                 //report progress
567                 if(num % 100 != 0){             cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n");  }
568                 
569                                 
570                 return 0;
571         }
572         catch(exception& e) {
573                 m->errorOut(e, "ChimeraSlayerCommand", "driverMPI");
574                 exit(1);
575         }
576 }
577 #endif
578
579 /**************************************************************************************************/
580
581 int ChimeraSlayerCommand::createProcesses(string outputFileName, string filename, string accnos) {
582         try {
583 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
584                 int process = 0;
585                 //              processIDS.resize(0);
586                 
587                 //loop through and create all the processes you want
588                 while (process != processors) {
589                         int pid = fork();
590                         
591                         if (pid > 0) {
592                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
593                                 process++;
594                         }else if (pid == 0){
595                                 driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename, accnos + toString(getpid()) + ".temp");
596                                 exit(0);
597                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
598                 }
599                 
600                 //force parent to wait until all the processes are done
601                 for (int i=0;i<processors;i++) { 
602                         int temp = processIDS[i];
603                         wait(&temp);
604                 }
605                 
606                 return 0;
607 #endif          
608         }
609         catch(exception& e) {
610                 m->errorOut(e, "ChimeraSlayerCommand", "createProcesses");
611                 exit(1);
612         }
613 }
614
615 /**************************************************************************************************/
616
617