]> git.donarmstrong.com Git - mothur.git/blob - chimeraslayercommand.cpp
4b309c0124a2038cb49db3980345e10e535ce84d
[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, filter, mask, 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 mask parameter allows you to specify a file containing one sequence you wish to use as a mask for the your sequences. \n");
148                 m->mothurOut("The window parameter allows you to specify the window size for searching for chimeras, default=50. \n");
149                 m->mothurOut("The increment parameter allows you to specify how far you move each window while finding chimeric sequences, default=5.\n");
150                 m->mothurOut("The numwanted parameter allows you to specify how many sequences you would each query sequence compared with, default=15.\n");
151                 m->mothurOut("The ksize parameter allows you to input kmersize, default is 7, used if search is kmer. \n");
152                 m->mothurOut("The match parameter allows you to reward matched bases in blast search, default is 5. \n");
153                 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");
154                 m->mothurOut("The mismatch parameter allows you to penalize mismatched bases in blast search, default is -4. \n");
155                 m->mothurOut("The divergence parameter allows you to set a cutoff for chimera determination, default is 1.007. \n");
156                 m->mothurOut("The iters parameter allows you to specify the number of bootstrap iters to do with the chimeraslayer method, default=100.\n");
157                 m->mothurOut("The minsim parameter allows you to specify a minimum similarity with the parent fragments, default=90. \n");
158                 m->mothurOut("The mincov parameter allows you to specify minimum coverage by closest matches found in template. Default is 70, meaning 70%. \n");
159                 m->mothurOut("The minbs parameter allows you to specify minimum bootstrap support for calling a sequence chimeric. Default is 90, meaning 90%. \n");
160                 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");
161                 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");
162                 //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");
163                 m->mothurOut("NOT ALL PARAMETERS ARE USED BY ALL METHODS. Please look below for method specifics.\n\n");
164                 m->mothurOut("The chimera.slayer command should be in the following format: \n");
165                 m->mothurOut("chimera.slayer(fasta=yourFastaFile, template=yourTemplate, search=yourSearch) \n");
166                 m->mothurOut("Example: chimera.slayer(fasta=AD.align, template=core_set_aligned.imputed.fasta, search=kmer) \n");
167                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");     
168         }
169         catch(exception& e) {
170                 m->errorOut(e, "ChimeraSlayerCommand", "help");
171                 exit(1);
172         }
173 }
174
175 //***************************************************************************************************************
176
177 ChimeraSlayerCommand::~ChimeraSlayerCommand(){  /*      do nothing      */      }
178
179 //***************************************************************************************************************
180
181 int ChimeraSlayerCommand::execute(){
182         try{
183                 
184                 if (abort == true) { return 0; }
185                 
186                 int start = time(NULL); 
187                 
188                 chimera = new ChimeraSlayer(fastafile, templatefile, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);   
189                                                 
190                 string outputFileName = outputDir + getRootName(getSimpleName(fastafile)) + "slayer.chimeras";
191                 string accnosFileName = outputDir + getRootName(getSimpleName(fastafile))  + "slayer.accnos";
192                 bool hasAccnos = true;
193                 
194                 if (m->control_pressed) { delete chimera;       return 0;       }
195                 
196                 if (chimera->getUnaligned()) { 
197                         m->mothurOut("Your template sequences are different lengths, please correct."); m->mothurOutEndLine(); 
198                         delete chimera;
199                         return 0; 
200                 }
201                 templateSeqsLength = chimera->getLength();
202                 
203         #ifdef USE_MPI  
204                 int pid, end, numSeqsPerProcessor; 
205                         int tag = 2001;
206                         vector<long> MPIPos;
207                         MPIWroteAccnos = false;
208                         
209                         MPI_Status status; 
210                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
211                         MPI_Comm_size(MPI_COMM_WORLD, &processors); 
212
213                         MPI_File inMPI;
214                         MPI_File outMPI;
215                         MPI_File outMPIAccnos;
216                         
217                         int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
218                         int inMode=MPI_MODE_RDONLY; 
219                         
220                         char* outFilename = new char[outputFileName.length()];
221                         memcpy(outFilename, outputFileName.c_str(), outputFileName.length());
222                         
223                         char* outAccnosFilename = new char[accnosFileName.length()];
224                         memcpy(outAccnosFilename, accnosFileName.c_str(), accnosFileName.length());
225
226                         char* inFileName = new char[fastafile.length()];
227                         memcpy(inFileName, fastafile.c_str(), fastafile.length());
228
229                         MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
230                         MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
231                         MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
232                         
233                         delete inFileName;
234                         delete outFilename;
235                         delete outAccnosFilename;
236
237                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   MPI_File_close(&outMPIAccnos);  delete chimera; return 0;  }
238                 
239                         if (pid == 0) { //you are the root process 
240                                 m->mothurOutEndLine();
241                                 m->mothurOut("Only reporting sequence supported by " + toString(minBS) + "% of bootstrapped results.");
242                                 m->mothurOutEndLine();
243         
244                                 string outTemp = "Name\tLeftParent\tRightParent\tDivQLAQRB\tPerIDQLAQRB\tBootStrapA\tDivQLBQRA\tPerIDQLBQRA\tBootStrapB\tFlag\tLeftWindow\tRightWindow\n";
245                                 
246                                 //print header
247                                 int length = outTemp.length();
248                                 char* buf2 = new char[length];
249                                 memcpy(buf2, outTemp.c_str(), length);
250
251                                 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &status);
252                                 delete buf2;
253
254                                 MPIPos = setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
255                                 
256                                 //send file positions to all processes
257                                 MPI_Bcast(&numSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD);  //send numSeqs
258                                 MPI_Bcast(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //send file pos        
259                                 
260                                 //figure out how many sequences you have to align
261                                 numSeqsPerProcessor = numSeqs / processors;
262                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
263                                 int startIndex =  pid * numSeqsPerProcessor;
264                         
265                                 //align your part
266                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
267                                 
268                                 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;  }
269
270                                 for (int i = 1; i < processors; i++) {
271                                         bool tempResult;
272                                         MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
273                                         if (tempResult != 0) { MPIWroteAccnos = true; }
274                                 }
275                         }else{ //you are a child process
276                                 MPI_Bcast(&numSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD); //get numSeqs
277                                 MPIPos.resize(numSeqs+1);
278                                 MPI_Bcast(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //get file positions
279                                 
280                                 //figure out how many sequences you have to align
281                                 numSeqsPerProcessor = numSeqs / processors;
282                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
283                                 int startIndex =  pid * numSeqsPerProcessor;
284                                 
285                                 //align your part
286                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
287                                 
288                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   MPI_File_close(&outMPIAccnos);  delete chimera; return 0;  }
289
290                                 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD); 
291                         }
292                         
293                         //close files 
294                         MPI_File_close(&inMPI);
295                         MPI_File_close(&outMPI);
296                         MPI_File_close(&outMPIAccnos);
297                         
298                         //delete accnos file if blank
299                         if (pid == 0) {
300                                 if (!MPIWroteAccnos) { 
301                                         //MPI_Info info;
302                                         //MPI_File_delete(outAccnosFilename, info);
303                                         hasAccnos = false;      
304                                         remove(accnosFileName.c_str()); 
305                                 }
306                         }
307                 
308         #else
309                 ofstream outHeader;
310                 string tempHeader = outputDir + getRootName(getSimpleName(fastafile)) + "slayer.chimeras.tempHeader";
311                 openOutputFile(tempHeader, outHeader);
312                 
313                 chimera->printHeader(outHeader);
314                 outHeader.close();
315                 
316                 //break up file
317                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
318                         if(processors == 1){
319                                 ifstream inFASTA;
320                                 openInputFile(fastafile, inFASTA);
321                                 numSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
322                                 inFASTA.close();
323                                 
324                                 lines.push_back(new linePair(0, numSeqs));
325                                 
326                                 driver(lines[0], outputFileName, fastafile, accnosFileName);
327                                 
328                                 if (m->control_pressed) { 
329                                         remove(outputFileName.c_str()); 
330                                         remove(tempHeader.c_str()); 
331                                         remove(accnosFileName.c_str());
332                                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
333                                         delete chimera;
334                                         return 0;
335                                 }
336                                 
337                                 //delete accnos file if its blank 
338                                 if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
339                                                                 
340                         }else{
341                                 vector<int> positions;
342                                 processIDS.resize(0);
343                                 
344                                 ifstream inFASTA;
345                                 openInputFile(fastafile, inFASTA);
346                                 
347                                 string input;
348                                 while(!inFASTA.eof()){
349                                         input = getline(inFASTA);
350                                         if (input.length() != 0) {
351                                                 if(input[0] == '>'){    long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);  }
352                                         }
353                                 }
354                                 inFASTA.close();
355                                 
356                                 numSeqs = positions.size();
357                                 
358                                 int numSeqsPerProcessor = numSeqs / processors;
359                                 
360                                 for (int i = 0; i < processors; i++) {
361                                         long int startPos = positions[ i * numSeqsPerProcessor ];
362                                         if(i == processors - 1){
363                                                 numSeqsPerProcessor = numSeqs - i * numSeqsPerProcessor;
364                                         }
365                                         lines.push_back(new linePair(startPos, numSeqsPerProcessor));
366                                 }
367                                 
368                                 
369                                 createProcesses(outputFileName, fastafile, accnosFileName); 
370                         
371                                 rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
372                                         
373                                 //append output files
374                                 for(int i=1;i<processors;i++){
375                                         appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
376                                         remove((outputFileName + toString(processIDS[i]) + ".temp").c_str());
377                                 }
378                                 
379                                 vector<string> nonBlankAccnosFiles;
380                                 //delete blank accnos files generated with multiple processes
381                                 for(int i=0;i<processors;i++){  
382                                         if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
383                                                 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
384                                         }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());  }
385                                 }
386                                 
387                                 //append accnos files
388                                 if (nonBlankAccnosFiles.size() != 0) { 
389                                         rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
390                                         
391                                         for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
392                                                 appendFiles(nonBlankAccnosFiles[h], accnosFileName);
393                                                 remove(nonBlankAccnosFiles[h].c_str());
394                                         }
395                                 }else{ hasAccnos = false;  }
396                                 
397                                 if (m->control_pressed) { 
398                                         remove(outputFileName.c_str()); 
399                                         remove(accnosFileName.c_str());
400                                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
401                                         delete chimera;
402                                         return 0;
403                                 }
404
405                         }
406
407                 #else
408                         ifstream inFASTA;
409                         openInputFile(fastafile, inFASTA);
410                         numSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
411                         inFASTA.close();
412                         lines.push_back(new linePair(0, numSeqs));
413                         
414                         driver(lines[0], outputFileName, fastafile, accnosFileName);
415                         
416                         if (m->control_pressed) { 
417                                         remove(outputFileName.c_str()); 
418                                         remove(tempHeader.c_str()); 
419                                         remove(accnosFileName.c_str());
420                                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
421                                         delete chimera;
422                                         return 0;
423                         }
424                         
425                         //delete accnos file if its blank 
426                         if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
427                 #endif
428                 
429                 appendFiles(outputFileName, tempHeader);
430         
431                 remove(outputFileName.c_str());
432                 rename(tempHeader.c_str(), outputFileName.c_str());
433                 
434         #endif
435                 delete chimera;
436                 
437                 m->mothurOutEndLine();
438                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
439                 m->mothurOut(outputFileName); m->mothurOutEndLine();    
440                 if (hasAccnos) {  m->mothurOut(accnosFileName); m->mothurOutEndLine();  }
441                 m->mothurOutEndLine();
442
443                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
444                 
445                 m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine();
446                 
447                 return 0;
448                 
449         }
450         catch(exception& e) {
451                 m->errorOut(e, "ChimeraSlayerCommand", "execute");
452                 exit(1);
453         }
454 }
455 //**********************************************************************************************************************
456
457 int ChimeraSlayerCommand::driver(linePair* line, string outputFName, string filename, string accnos){
458         try {
459                 ofstream out;
460                 openOutputFile(outputFName, out);
461                 
462                 ofstream out2;
463                 openOutputFile(accnos, out2);
464                 
465                 ifstream inFASTA;
466                 openInputFile(filename, inFASTA);
467
468                 inFASTA.seekg(line->start);
469                 
470                 for(int i=0;i<line->numSeqs;i++){
471                 
472                         if (m->control_pressed) {       return 1;       }
473                 
474                         Sequence* candidateSeq = new Sequence(inFASTA);  gobble(inFASTA);
475                                 
476                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
477                                 
478                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
479                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
480                                 }else{
481                                         //find chimeras
482                                         chimera->getChimeras(candidateSeq);
483                                         
484                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
485                 
486                                         //print results
487                                         chimera->print(out, out2);
488                                 }
489                         }
490                         delete candidateSeq;
491                         
492                         //report progress
493                         if((i+1) % 100 == 0){   m->mothurOut("Processing sequence: " + toString(i+1)); m->mothurOutEndLine();           }
494                 }
495                 //report progress
496                 if((line->numSeqs) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(line->numSeqs)); m->mothurOutEndLine();         }
497                 
498                 out.close();
499                 out2.close();
500                 inFASTA.close();
501                                 
502                 return 0;
503         }
504         catch(exception& e) {
505                 m->errorOut(e, "ChimeraSlayerCommand", "driver");
506                 exit(1);
507         }
508 }
509 //**********************************************************************************************************************
510 #ifdef USE_MPI
511 int ChimeraSlayerCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, MPI_File& outAccMPI, vector<long>& MPIPos){
512         try {                           
513                 MPI_Status status; 
514                 int pid;
515                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
516                 
517                 for(int i=0;i<num;i++){
518                         
519                         if (m->control_pressed) {       return 1;       }
520                         
521                         //read next sequence
522                         int length = MPIPos[start+i+1] - MPIPos[start+i];
523
524                         char* buf4 = new char[length];
525                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
526         
527                         string tempBuf = buf4;
528                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
529                         istringstream iss (tempBuf,istringstream::in);
530
531                         delete buf4;
532
533                         Sequence* candidateSeq = new Sequence(iss);  gobble(iss);
534                 
535                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
536                                 
537                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
538                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
539                                 }else{
540                 
541                                         //find chimeras
542                                         chimera->getChimeras(candidateSeq);
543                         
544                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
545                 //cout << "about to print" << endl;
546                                         //print results
547                                         bool isChimeric = chimera->print(outMPI, outAccMPI);
548                                         if (isChimeric) { MPIWroteAccnos = true;  }
549         
550                                 }
551                         }
552                         delete candidateSeq;
553                         
554                         //report progress
555                         if((i+1) % 100 == 0){  cout << "Processing sequence: " << (i+1) << endl;        m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n");          }
556                 }
557                 //report progress
558                 if(num % 100 != 0){             cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n");  }
559                 
560                                 
561                 return 0;
562         }
563         catch(exception& e) {
564                 m->errorOut(e, "ChimeraSlayerCommand", "driverMPI");
565                 exit(1);
566         }
567 }
568 #endif
569
570 /**************************************************************************************************/
571
572 int ChimeraSlayerCommand::createProcesses(string outputFileName, string filename, string accnos) {
573         try {
574 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
575                 int process = 0;
576                 //              processIDS.resize(0);
577                 
578                 //loop through and create all the processes you want
579                 while (process != processors) {
580                         int pid = fork();
581                         
582                         if (pid > 0) {
583                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
584                                 process++;
585                         }else if (pid == 0){
586                                 driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename, accnos + toString(getpid()) + ".temp");
587                                 exit(0);
588                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
589                 }
590                 
591                 //force parent to wait until all the processes are done
592                 for (int i=0;i<processors;i++) { 
593                         int temp = processIDS[i];
594                         wait(&temp);
595                 }
596                 
597                 return 0;
598 #endif          
599         }
600         catch(exception& e) {
601                 m->errorOut(e, "ChimeraSlayerCommand", "createProcesses");
602                 exit(1);
603         }
604 }
605
606 /**************************************************************************************************/
607
608