]> git.donarmstrong.com Git - mothur.git/blob - chimeraslayercommand.cpp
modified chimera commands to process multiple fasta files and added checks to pintail...
[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-amzon.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                         bool hasAccnos = true;
234                         
235                         if (m->control_pressed) { delete chimera; for (int j = 0; j < outputNames.size(); j++) {        remove(outputNames[j].c_str()); }  return 0;    }
236                         
237                         if (chimera->getUnaligned()) { 
238                                 m->mothurOut("Your template sequences are different lengths, please correct."); m->mothurOutEndLine(); 
239                                 delete chimera;
240                                 return 0; 
241                         }
242                         templateSeqsLength = chimera->getLength();
243                         
244                 #ifdef USE_MPI  
245                         int pid, end, numSeqsPerProcessor; 
246                                 int tag = 2001;
247                                 vector<long> MPIPos;
248                                 MPIWroteAccnos = false;
249                                 
250                                 MPI_Status status; 
251                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
252                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); 
253
254                                 MPI_File inMPI;
255                                 MPI_File outMPI;
256                                 MPI_File outMPIAccnos;
257                                 
258                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
259                                 int inMode=MPI_MODE_RDONLY; 
260                                 
261                                 char outFilename[1024];
262                                 strcpy(outFilename, outputFileName.c_str());
263                                 
264                                 char outAccnosFilename[1024];
265                                 strcpy(outAccnosFilename, accnosFileName.c_str());
266                                 
267                                 char inFileName[1024];
268                                 strcpy(inFileName, fastaFileNames[s].c_str());
269
270                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
271                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
272                                 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
273
274                                 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;  }
275                         
276                                 if (pid == 0) { //you are the root process 
277                                         m->mothurOutEndLine();
278                                         m->mothurOut("Only reporting sequence supported by " + toString(minBS) + "% of bootstrapped results.");
279                                         m->mothurOutEndLine();
280                 
281                                         string outTemp = "Name\tLeftParent\tRightParent\tDivQLAQRB\tPerIDQLAQRB\tBootStrapA\tDivQLBQRA\tPerIDQLBQRA\tBootStrapB\tFlag\tLeftWindow\tRightWindow\n";
282                                         
283                                         //print header
284                                         int length = outTemp.length();
285                                         char* buf2 = new char[length];
286                                         memcpy(buf2, outTemp.c_str(), length);
287
288                                         MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &status);
289                                         delete buf2;
290
291                                         MPIPos = setFilePosFasta(fastaFileNames[s], numSeqs); //fills MPIPos, returns numSeqs
292                                         
293                                         //send file positions to all processes
294                                         for(int i = 1; i < processors; i++) { 
295                                                 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
296                                                 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
297                                         }
298                                         
299                                         //figure out how many sequences you have to align
300                                         numSeqsPerProcessor = numSeqs / processors;
301                                         int startIndex =  pid * numSeqsPerProcessor;
302                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
303                                 
304                                         //do your part
305                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
306                                         
307                                         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;  }
308
309                                         for (int i = 1; i < processors; i++) {
310                                                 bool tempResult;
311                                                 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
312                                                 if (tempResult != 0) { MPIWroteAccnos = true; }
313                                         }
314                                 }else{ //you are a child process
315                                         MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
316                                         MPIPos.resize(numSeqs+1);
317                                         MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
318                                         
319                                         //figure out how many sequences you have to align
320                                         numSeqsPerProcessor = numSeqs / processors;
321                                         int startIndex =  pid * numSeqsPerProcessor;
322                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
323                                         
324                                         //do your part
325                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
326                                         
327                                         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;  }
328
329                                         MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD); 
330                                 }
331                                 
332                                 //close files 
333                                 MPI_File_close(&inMPI);
334                                 MPI_File_close(&outMPI);
335                                 MPI_File_close(&outMPIAccnos);
336                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
337                                 
338                                 //delete accnos file if blank
339                                 if (pid == 0) {
340                                         if (!MPIWroteAccnos) { 
341                                                 hasAccnos = false;      
342                                                 remove(accnosFileName.c_str()); 
343                                         }
344                                 }
345                         
346                 #else
347                         ofstream outHeader;
348                         string tempHeader = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + "slayer.chimeras.tempHeader";
349                         openOutputFile(tempHeader, outHeader);
350                         
351                         chimera->printHeader(outHeader);
352                         outHeader.close();
353                         
354                         //break up file
355                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
356                                 if(processors == 1){
357                                         ifstream inFASTA;
358                                         openInputFile(fastaFileNames[s], inFASTA);
359                                         getNumSeqs(inFASTA, numSeqs);
360                                         inFASTA.close();
361                                         
362                                         lines.push_back(new linePair(0, numSeqs));
363                                         
364                                         driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName);
365                                         
366                                         if (m->control_pressed) { 
367                                                 remove(outputFileName.c_str()); 
368                                                 remove(tempHeader.c_str()); 
369                                                 remove(accnosFileName.c_str());
370                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
371                                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
372                                                 delete chimera;
373                                                 return 0;
374                                         }
375                                         
376                                         //delete accnos file if its blank 
377                                         if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
378                                                                         
379                                 }else{
380                                         vector<int> positions;
381                                         processIDS.resize(0);
382                                         
383                                         ifstream inFASTA;
384                                         openInputFile(fastaFileNames[s], inFASTA);
385                                         
386                                         string input;
387                                         while(!inFASTA.eof()){
388                                                 input = getline(inFASTA);
389                                                 if (input.length() != 0) {
390                                                         if(input[0] == '>'){    long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);  }
391                                                 }
392                                         }
393                                         inFASTA.close();
394                                         
395                                         numSeqs = positions.size();
396                                         
397                                         int numSeqsPerProcessor = numSeqs / processors;
398                                         
399                                         for (int i = 0; i < processors; i++) {
400                                                 long int startPos = positions[ i * numSeqsPerProcessor ];
401                                                 if(i == processors - 1){
402                                                         numSeqsPerProcessor = numSeqs - i * numSeqsPerProcessor;
403                                                 }
404                                                 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
405                                         }
406                                         
407                                         createProcesses(outputFileName, fastaFileNames[s], accnosFileName); 
408                                 
409                                         rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
410                                                 
411                                         //append output files
412                                         for(int i=1;i<processors;i++){
413                                                 appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
414                                                 remove((outputFileName + toString(processIDS[i]) + ".temp").c_str());
415                                         }
416                                         
417                                         vector<string> nonBlankAccnosFiles;
418                                         //delete blank accnos files generated with multiple processes
419                                         for(int i=0;i<processors;i++){  
420                                                 if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
421                                                         nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
422                                                 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());  }
423                                         }
424                                         
425                                         //append accnos files
426                                         if (nonBlankAccnosFiles.size() != 0) { 
427                                                 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
428                                                 
429                                                 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
430                                                         appendFiles(nonBlankAccnosFiles[h], accnosFileName);
431                                                         remove(nonBlankAccnosFiles[h].c_str());
432                                                 }
433                                         }else{ hasAccnos = false;  }
434                                         
435                                         if (m->control_pressed) { 
436                                                 remove(outputFileName.c_str()); 
437                                                 remove(accnosFileName.c_str());
438                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
439                                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
440                                                 delete chimera;
441                                                 return 0;
442                                         }
443
444                                 }
445
446                         #else
447                                 ifstream inFASTA;
448                                 openInputFile(fastaFileNames[s], inFASTA);
449                                 getNumSeqs(inFASTA, numSeqs);
450                                 inFASTA.close();
451                                 lines.push_back(new linePair(0, numSeqs));
452                                 
453                                 driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName);
454                                 
455                                 if (m->control_pressed) { 
456                                                 remove(outputFileName.c_str()); 
457                                                 remove(tempHeader.c_str()); 
458                                                 remove(accnosFileName.c_str());
459                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
460                                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
461                                                 delete chimera;
462                                                 return 0;
463                                 }
464                                 
465                                 //delete accnos file if its blank 
466                                 if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
467                         #endif
468                         
469                         appendFiles(outputFileName, tempHeader);
470                 
471                         remove(outputFileName.c_str());
472                         rename(tempHeader.c_str(), outputFileName.c_str());
473                         
474                 #endif
475                         delete chimera;
476                         
477                         
478                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
479                         
480                         outputNames.push_back(outputFileName);
481                         if (hasAccnos) { outputNames.push_back(accnosFileName); }
482                         
483                         m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine();
484                         }
485                 
486                 m->mothurOutEndLine();
487                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
488                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }       
489                 m->mothurOutEndLine();
490
491                 return 0;
492                 
493         }
494         catch(exception& e) {
495                 m->errorOut(e, "ChimeraSlayerCommand", "execute");
496                 exit(1);
497         }
498 }
499 //**********************************************************************************************************************
500
501 int ChimeraSlayerCommand::driver(linePair* line, string outputFName, string filename, string accnos){
502         try {
503                 ofstream out;
504                 openOutputFile(outputFName, out);
505                 
506                 ofstream out2;
507                 openOutputFile(accnos, out2);
508                 
509                 ifstream inFASTA;
510                 openInputFile(filename, inFASTA);
511
512                 inFASTA.seekg(line->start);
513                 
514                 for(int i=0;i<line->numSeqs;i++){
515                 
516                         if (m->control_pressed) {       return 1;       }
517                 
518                         Sequence* candidateSeq = new Sequence(inFASTA);  gobble(inFASTA);
519                                 
520                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
521                                 
522                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
523                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
524                                 }else{
525                                         //find chimeras
526                                         chimera->getChimeras(candidateSeq);
527                                         
528                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
529                 
530                                         //print results
531                                         chimera->print(out, out2);
532                                 }
533                         }
534                         delete candidateSeq;
535                         
536                         //report progress
537                         if((i+1) % 100 == 0){   m->mothurOut("Processing sequence: " + toString(i+1)); m->mothurOutEndLine();           }
538                 }
539                 //report progress
540                 if((line->numSeqs) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(line->numSeqs)); m->mothurOutEndLine();         }
541                 
542                 out.close();
543                 out2.close();
544                 inFASTA.close();
545                                 
546                 return 0;
547         }
548         catch(exception& e) {
549                 m->errorOut(e, "ChimeraSlayerCommand", "driver");
550                 exit(1);
551         }
552 }
553 //**********************************************************************************************************************
554 #ifdef USE_MPI
555 int ChimeraSlayerCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, MPI_File& outAccMPI, vector<long>& MPIPos){
556         try {                           
557                 MPI_Status status; 
558                 int pid;
559                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
560                 
561                 for(int i=0;i<num;i++){
562                         
563                         if (m->control_pressed) {       return 1;       }
564                         
565                         //read next sequence
566                         int length = MPIPos[start+i+1] - MPIPos[start+i];
567
568                         char* buf4 = new char[length];
569                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
570         
571                         string tempBuf = buf4;
572                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
573                         istringstream iss (tempBuf,istringstream::in);
574
575                         delete buf4;
576
577                         Sequence* candidateSeq = new Sequence(iss);  gobble(iss);
578                 
579                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
580                                 
581                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
582                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
583                                 }else{
584                 
585                                         //find chimeras
586                                         chimera->getChimeras(candidateSeq);
587                         
588                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
589                 //cout << "about to print" << endl;
590                                         //print results
591                                         bool isChimeric = chimera->print(outMPI, outAccMPI);
592                                         if (isChimeric) { MPIWroteAccnos = true;  }
593         
594                                 }
595                         }
596                         delete candidateSeq;
597                         
598                         //report progress
599                         if((i+1) % 100 == 0){  cout << "Processing sequence: " << (i+1) << endl;        m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n");          }
600                 }
601                 //report progress
602                 if(num % 100 != 0){             cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n");  }
603                 
604                                 
605                 return 0;
606         }
607         catch(exception& e) {
608                 m->errorOut(e, "ChimeraSlayerCommand", "driverMPI");
609                 exit(1);
610         }
611 }
612 #endif
613
614 /**************************************************************************************************/
615
616 int ChimeraSlayerCommand::createProcesses(string outputFileName, string filename, string accnos) {
617         try {
618 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
619                 int process = 0;
620                 //              processIDS.resize(0);
621                 
622                 //loop through and create all the processes you want
623                 while (process != processors) {
624                         int pid = fork();
625                         
626                         if (pid > 0) {
627                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
628                                 process++;
629                         }else if (pid == 0){
630                                 driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename, accnos + toString(getpid()) + ".temp");
631                                 exit(0);
632                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
633                 }
634                 
635                 //force parent to wait until all the processes are done
636                 for (int i=0;i<processors;i++) { 
637                         int temp = processIDS[i];
638                         wait(&temp);
639                 }
640                 
641                 return 0;
642 #endif          
643         }
644         catch(exception& e) {
645                 m->errorOut(e, "ChimeraSlayerCommand", "createProcesses");
646                 exit(1);
647         }
648 }
649
650 /**************************************************************************************************/
651
652