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