]> git.donarmstrong.com Git - mothur.git/blob - chimeraslayercommand.cpp
fixed bug with blastdb.cpp
[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 #include "deconvolutecommand.h"
13
14 //**********************************************************************************************************************
15 vector<string> ChimeraSlayerCommand::setParameters(){   
16         try {
17                 CommandParameter ptemplate("reference", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(ptemplate);
18                 CommandParameter pfasta("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pfasta);
19                 CommandParameter pname("name", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(pname);
20                 CommandParameter pwindow("window", "Number", "", "50", "", "", "",false,false); parameters.push_back(pwindow);
21                 CommandParameter pksize("ksize", "Number", "", "7", "", "", "",false,false); parameters.push_back(pksize);
22                 CommandParameter pmatch("match", "Number", "", "5.0", "", "", "",false,false); parameters.push_back(pmatch);
23                 CommandParameter pmismatch("mismatch", "Number", "", "-4.0", "", "", "",false,false); parameters.push_back(pmismatch);
24                 CommandParameter pminsim("minsim", "Number", "", "90", "", "", "",false,false); parameters.push_back(pminsim);
25                 CommandParameter pmincov("mincov", "Number", "", "70", "", "", "",false,false); parameters.push_back(pmincov);
26                 CommandParameter pminsnp("minsnp", "Number", "", "10", "", "", "",false,false); parameters.push_back(pminsnp);
27                 CommandParameter pminbs("minbs", "Number", "", "90", "", "", "",false,false); parameters.push_back(pminbs);
28                 CommandParameter psearch("search", "Multiple", "kmer-blast", "blast", "", "", "",false,false); parameters.push_back(psearch);
29                 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
30                 CommandParameter prealign("realign", "Boolean", "", "T", "", "", "",false,false); parameters.push_back(prealign);
31                 CommandParameter ptrim("trim", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(ptrim);
32                 CommandParameter psplit("split", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(psplit);
33                 CommandParameter pnumwanted("numwanted", "Number", "", "15", "", "", "",false,false); parameters.push_back(pnumwanted);
34                 CommandParameter piters("iters", "Number", "", "1000", "", "", "",false,false); parameters.push_back(piters);
35                 CommandParameter pdivergence("divergence", "Number", "", "1.007", "", "", "",false,false); parameters.push_back(pdivergence);
36                 CommandParameter pparents("parents", "Number", "", "3", "", "", "",false,false); parameters.push_back(pparents);
37                 CommandParameter pincrement("increment", "Number", "", "5", "", "", "",false,false); parameters.push_back(pincrement);
38                 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
39                 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
40                 
41                 vector<string> myArray;
42                 for (int i = 0; i < parameters.size(); i++) {   myArray.push_back(parameters[i].name);          }
43                 return myArray;
44         }
45         catch(exception& e) {
46                 m->errorOut(e, "ChimeraSlayerCommand", "setParameters");
47                 exit(1);
48         }
49 }
50 //**********************************************************************************************************************
51 string ChimeraSlayerCommand::getHelpString(){   
52         try {
53                 string helpString = "";
54                 helpString += "The chimera.slayer command reads a fastafile and referencefile and outputs potentially chimeric sequences.\n";
55                 helpString += "This command was modeled after the chimeraSlayer written by the Broad Institute.\n";
56                 helpString += "The chimera.slayer command parameters are fasta, name, template, processors, trim, ksize, window, match, mismatch, divergence. minsim, mincov, minbs, minsnp, parents, search, iters, increment, numwanted, and realign.\n";
57                 helpString += "The fasta parameter allows you to enter the fasta file containing your potentially chimeric sequences, and is required, unless you have a valid current fasta file. \n";
58                 helpString += "The name parameter allows you to provide a name file, if you are using template=self. \n";
59                 helpString += "You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amazon.fasta \n";
60                 helpString += "The reference parameter allows you to enter a reference file containing known non-chimeric sequences, and is required. You may also set template=self, in this case the abundant sequences will be used as potential parents. \n";
61                 helpString += "The processors parameter allows you to specify how many processors you would like to use.  The default is 1. \n";
62 #ifdef USE_MPI
63                 helpString += "When using MPI, the processors parameter is set to the number of MPI processes running. \n";
64 #endif
65                 helpString += "The trim parameter allows you to output a new fasta file containing your sequences with the chimeric ones trimmed to include only their longest piece, default=F. \n";
66                 helpString += "The split parameter allows you to check both pieces of non-chimeric sequence for chimeras, thus looking for trimeras and quadmeras. default=F. \n";
67                 helpString += "The window parameter allows you to specify the window size for searching for chimeras, default=50. \n";
68                 helpString += "The increment parameter allows you to specify how far you move each window while finding chimeric sequences, default=5.\n";
69                 helpString += "The numwanted parameter allows you to specify how many sequences you would each query sequence compared with, default=15.\n";
70                 helpString += "The ksize parameter allows you to input kmersize, default is 7, used if search is kmer. \n";
71                 helpString += "The match parameter allows you to reward matched bases in blast search, default is 5. \n";
72                 helpString += "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";
73                 helpString += "The mismatch parameter allows you to penalize mismatched bases in blast search, default is -4. \n";
74                 helpString += "The divergence parameter allows you to set a cutoff for chimera determination, default is 1.007. \n";
75                 helpString += "The iters parameter allows you to specify the number of bootstrap iters to do with the chimeraslayer method, default=1000.\n";
76                 helpString += "The minsim parameter allows you to specify a minimum similarity with the parent fragments, default=90. \n";
77                 helpString += "The mincov parameter allows you to specify minimum coverage by closest matches found in template. Default is 70, meaning 70%. \n";
78                 helpString += "The minbs parameter allows you to specify minimum bootstrap support for calling a sequence chimeric. Default is 90, meaning 90%. \n";
79                 helpString += "The minsnp parameter allows you to specify percent of SNPs to sample on each side of breakpoint for computing bootstrap support (default: 10) \n";
80                 helpString += "The search parameter allows you to specify search method for finding the closest parent. Choices are blast, and kmer, default blast. \n";
81                 helpString += "The realign parameter allows you to realign the query to the potential parents. Choices are true or false, default true.  \n";
82                 helpString += "The chimera.slayer command should be in the following format: \n";
83                 helpString += "chimera.slayer(fasta=yourFastaFile, reference=yourTemplate, search=yourSearch) \n";
84                 helpString += "Example: chimera.slayer(fasta=AD.align, reference=core_set_aligned.imputed.fasta, search=kmer) \n";
85                 helpString += "Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n";       
86                 return helpString;
87         }
88         catch(exception& e) {
89                 m->errorOut(e, "ChimeraSlayerCommand", "getHelpString");
90                 exit(1);
91         }
92 }
93 //**********************************************************************************************************************
94 ChimeraSlayerCommand::ChimeraSlayerCommand(){   
95         try {
96                 abort = true; calledHelp = true;
97                 setParameters();
98                 vector<string> tempOutNames;
99                 outputTypes["chimera"] = tempOutNames;
100                 outputTypes["accnos"] = tempOutNames;
101                 outputTypes["fasta"] = tempOutNames;
102         }
103         catch(exception& e) {
104                 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
105                 exit(1);
106         }
107 }
108 //***************************************************************************************************************
109 ChimeraSlayerCommand::ChimeraSlayerCommand(string option)  {
110         try {
111                 abort = false; calledHelp = false;   
112                 
113                 //allow user to run help
114                 if(option == "help") { help(); abort = true; calledHelp = true; }
115                 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
116                 
117                 else {
118                         vector<string> myArray = setParameters();
119                         
120                         OptionParser parser(option);
121                         map<string,string> parameters = parser.getParameters();
122                         
123                         ValidParameters validParameter("chimera.slayer");
124                         map<string,string>::iterator it;
125                         
126                         //check to make sure all parameters are valid for command
127                         for (it = parameters.begin(); it != parameters.end(); it++) { 
128                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
129                         }
130                         
131                         vector<string> tempOutNames;
132                         outputTypes["chimera"] = tempOutNames;
133                         outputTypes["accnos"] = tempOutNames;
134                         outputTypes["fasta"] = tempOutNames;
135                 
136                         //if the user changes the input directory command factory will send this info to us in the output parameter 
137                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
138                         if (inputDir == "not found"){   inputDir = "";          }
139                                                 
140                         //check for required parameters
141                         fastafile = validParameter.validFile(parameters, "fasta", false);
142                         if (fastafile == "not found") {                                 
143                                 //if there is a current fasta file, use it
144                                 string filename = m->getFastaFile(); 
145                                 if (filename != "") { fastaFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
146                                 else {  m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
147                         }else { 
148                                 m->splitAtDash(fastafile, fastaFileNames);
149                                 
150                                 //go through files and make sure they are good, if not, then disregard them
151                                 for (int i = 0; i < fastaFileNames.size(); i++) {
152                                         
153                                         bool ignore = false;
154                                         if (fastaFileNames[i] == "current") { 
155                                                 fastaFileNames[i] = m->getFastaFile(); 
156                                                 if (fastaFileNames[i] != "") {  m->mothurOut("Using " + fastaFileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
157                                                 else {  
158                                                         m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true; 
159                                                         //erase from file list
160                                                         fastaFileNames.erase(fastaFileNames.begin()+i);
161                                                         i--;
162                                                 }
163                                         }
164                                         
165                                         if (!ignore) {
166                                                 
167                                                 if (inputDir != "") {
168                                                         string path = m->hasPath(fastaFileNames[i]);
169                                                         //if the user has not given a path then, add inputdir. else leave path alone.
170                                                         if (path == "") {       fastaFileNames[i] = inputDir + fastaFileNames[i];               }
171                                                 }
172                 
173                                                 int ableToOpen;
174                                                 ifstream in;
175                                                 
176                                                 ableToOpen = m->openInputFile(fastaFileNames[i], in, "noerror");
177                                         
178                                                 //if you can't open it, try default location
179                                                 if (ableToOpen == 1) {
180                                                         if (m->getDefaultPath() != "") { //default path is set
181                                                                 string tryPath = m->getDefaultPath() + m->getSimpleName(fastaFileNames[i]);
182                                                                 m->mothurOut("Unable to open " + fastaFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
183                                                                 ifstream in2;
184                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
185                                                                 in2.close();
186                                                                 fastaFileNames[i] = tryPath;
187                                                         }
188                                                 }
189                                                 
190                                                 if (ableToOpen == 1) {
191                                                         if (m->getOutputDir() != "") { //default path is set
192                                                                 string tryPath = m->getOutputDir() + m->getSimpleName(fastaFileNames[i]);
193                                                                 m->mothurOut("Unable to open " + fastaFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
194                                                                 ifstream in2;
195                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
196                                                                 in2.close();
197                                                                 fastaFileNames[i] = tryPath;
198                                                         }
199                                                 }
200                                                 
201                                                 in.close();
202                                                 
203                                                 if (ableToOpen == 1) { 
204                                                         m->mothurOut("Unable to open " + fastaFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
205                                                         //erase from file list
206                                                         fastaFileNames.erase(fastaFileNames.begin()+i);
207                                                         i--;
208                                                 }
209                                         }
210                                 }
211                                 
212                                 //make sure there is at least one valid file left
213                                 if (fastaFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid files."); m->mothurOutEndLine(); abort = true; }
214                         }
215                         
216                         
217                         //check for required parameters
218                         bool hasName = true;
219                         namefile = validParameter.validFile(parameters, "name", false);
220                         if (namefile == "not found") { namefile = "";  hasName = false; }
221                         else { 
222                                 m->splitAtDash(namefile, nameFileNames);
223                                 
224                                 //go through files and make sure they are good, if not, then disregard them
225                                 for (int i = 0; i < nameFileNames.size(); i++) {
226                                         
227                                         bool ignore = false;
228                                         if (nameFileNames[i] == "current") { 
229                                                 nameFileNames[i] = m->getNameFile(); 
230                                                 if (nameFileNames[i] != "") {  m->mothurOut("Using " + nameFileNames[i] + " as input file for the name parameter where you had given current."); m->mothurOutEndLine(); }
231                                                 else {  
232                                                         m->mothurOut("You have no current namefile, ignoring current."); m->mothurOutEndLine(); ignore=true; 
233                                                         //erase from file list
234                                                         nameFileNames.erase(nameFileNames.begin()+i);
235                                                         i--;
236                                                 }
237                                         }
238                                         
239                                         if (!ignore) {
240                                                 
241                                                 if (inputDir != "") {
242                                                         string path = m->hasPath(nameFileNames[i]);
243                                                         //if the user has not given a path then, add inputdir. else leave path alone.
244                                                         if (path == "") {       nameFileNames[i] = inputDir + nameFileNames[i];         }
245                                                 }
246                                                 
247                                                 int ableToOpen;
248                                                 ifstream in;
249                                                 
250                                                 ableToOpen = m->openInputFile(nameFileNames[i], in, "noerror");
251                                                 
252                                                 //if you can't open it, try default location
253                                                 if (ableToOpen == 1) {
254                                                         if (m->getDefaultPath() != "") { //default path is set
255                                                                 string tryPath = m->getDefaultPath() + m->getSimpleName(nameFileNames[i]);
256                                                                 m->mothurOut("Unable to open " + nameFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
257                                                                 ifstream in2;
258                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
259                                                                 in2.close();
260                                                                 nameFileNames[i] = tryPath;
261                                                         }
262                                                 }
263                                                 
264                                                 if (ableToOpen == 1) {
265                                                         if (m->getOutputDir() != "") { //default path is set
266                                                                 string tryPath = m->getOutputDir() + m->getSimpleName(nameFileNames[i]);
267                                                                 m->mothurOut("Unable to open " + nameFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
268                                                                 ifstream in2;
269                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
270                                                                 in2.close();
271                                                                 nameFileNames[i] = tryPath;
272                                                         }
273                                                 }
274                                                 
275                                                 in.close();
276                                                 
277                                                 if (ableToOpen == 1) { 
278                                                         m->mothurOut("Unable to open " + nameFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
279                                                         //erase from file list
280                                                         nameFileNames.erase(nameFileNames.begin()+i);
281                                                         i--;
282                                                 }
283                                         }
284                                 }
285                                 
286                                 //make sure there is at least one valid file left
287                                 if (nameFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid name files."); m->mothurOutEndLine(); abort = true; }
288                         }
289                         
290                         if (hasName && (nameFileNames.size() != fastaFileNames.size())) { m->mothurOut("[ERROR]: The number of namefiles does not match the number of fastafiles, please correct."); m->mothurOutEndLine(); abort=true; }
291                         
292                         //if the user changes the output directory command factory will send this info to us in the output parameter 
293                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  outputDir = ""; }
294                         
295                         
296                         string path;
297                         it = parameters.find("reference");
298                         //user has given a template file
299                         if(it != parameters.end()){ 
300                                 if (it->second == "self") { templatefile = "self"; }
301                                 else {
302                                         path = m->hasPath(it->second);
303                                         //if the user has not given a path then, add inputdir. else leave path alone.
304                                         if (path == "") {       parameters["reference"] = inputDir + it->second;                }
305                                         
306                                         templatefile = validParameter.validFile(parameters, "reference", true);
307                                         if (templatefile == "not open") { abort = true; }
308                                         else if (templatefile == "not found") { templatefile = "";  m->mothurOut("reference is a required parameter for the chimera.slayer command, unless and namefile is given."); m->mothurOutEndLine(); abort = true;  }    
309                                 }
310                         }else if (hasName) {  templatefile = "self"; }
311                         else { templatefile = "";  m->mothurOut("reference is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true;  }   
312                         
313                         string temp = validParameter.validFile(parameters, "processors", false);        if (temp == "not found"){       temp = m->getProcessors();      }
314                         m->setProcessors(temp);
315                         convert(temp, processors);
316                         
317                         temp = validParameter.validFile(parameters, "ksize", false);                    if (temp == "not found") { temp = "7"; }
318                         convert(temp, ksize);
319                                                 
320                         temp = validParameter.validFile(parameters, "window", false);                   if (temp == "not found") { temp = "50"; }                       
321                         convert(temp, window);
322                         
323                         temp = validParameter.validFile(parameters, "match", false);                    if (temp == "not found") { temp = "5"; }
324                         convert(temp, match);
325                         
326                         temp = validParameter.validFile(parameters, "mismatch", false);                 if (temp == "not found") { temp = "-4"; }
327                         convert(temp, mismatch);
328                         
329                         temp = validParameter.validFile(parameters, "divergence", false);               if (temp == "not found") { temp = "1.007"; }
330                         convert(temp, divR);
331                         
332                         temp = validParameter.validFile(parameters, "minsim", false);                   if (temp == "not found") { temp = "90"; }
333                         convert(temp, minSimilarity);
334                         
335                         temp = validParameter.validFile(parameters, "mincov", false);                   if (temp == "not found") { temp = "70"; }
336                         convert(temp, minCoverage);
337                         
338                         temp = validParameter.validFile(parameters, "minbs", false);                    if (temp == "not found") { temp = "90"; }
339                         convert(temp, minBS);
340                         
341                         temp = validParameter.validFile(parameters, "minsnp", false);                   if (temp == "not found") { temp = "10"; }
342                         convert(temp, minSNP);
343
344                         temp = validParameter.validFile(parameters, "parents", false);                  if (temp == "not found") { temp = "3"; }
345                         convert(temp, parents); 
346                         
347                         temp = validParameter.validFile(parameters, "realign", false);                  if (temp == "not found") { temp = "t"; }
348                         realign = m->isTrue(temp); 
349                         
350                         temp = validParameter.validFile(parameters, "trim", false);                             if (temp == "not found") { temp = "f"; }
351                         trim = m->isTrue(temp); 
352                         
353                         temp = validParameter.validFile(parameters, "split", false);                    if (temp == "not found") { temp = "f"; }
354                         trimera = m->isTrue(temp); 
355                         
356                         search = validParameter.validFile(parameters, "search", false);                 if (search == "not found") { search = "blast"; }
357                         
358                         temp = validParameter.validFile(parameters, "iters", false);                    if (temp == "not found") { temp = "1000"; }             
359                         convert(temp, iters); 
360                          
361                         temp = validParameter.validFile(parameters, "increment", false);                if (temp == "not found") { temp = "5"; }
362                         convert(temp, increment);
363                         
364                         temp = validParameter.validFile(parameters, "numwanted", false);                if (temp == "not found") { temp = "15"; }               
365                         convert(temp, numwanted);
366
367                         if ((search != "blast") && (search != "kmer")) { m->mothurOut(search + " is not a valid search."); m->mothurOutEndLine(); abort = true;  }
368                 }
369         }
370         catch(exception& e) {
371                 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
372                 exit(1);
373         }
374 }
375 //***************************************************************************************************************
376
377 int ChimeraSlayerCommand::execute(){
378         try{
379                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
380                 
381                 for (int s = 0; s < fastaFileNames.size(); s++) {
382                                 
383                         m->mothurOut("Checking sequences from " + fastaFileNames[s] + " ..." ); m->mothurOutEndLine();
384                 
385                         int start = time(NULL); 
386                         
387                         if (templatefile != "self") { //you want to run slayer with a refernce template
388                                 chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, trim, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);     
389                         }else {
390                                 if (processors != 1) { m->mothurOut("When using template=self, mothur can only use 1 processor, continuing."); m->mothurOutEndLine(); processors = 1; }
391                                 string nameFile = "";
392                                 if (nameFileNames.size() != 0) { //you provided a namefile and we don't need to create one
393                                         nameFile = nameFileNames[s];
394                                 }else {
395                                         m->mothurOutEndLine(); m->mothurOut("No namesfile given, running unique.seqs command to generate one."); m->mothurOutEndLine(); m->mothurOutEndLine();
396                                         
397                                         //use unique.seqs to create new name and fastafile
398                                         string inputString = "fasta=" + fastaFileNames[s];
399                                         m->mothurOut("/******************************************/"); m->mothurOutEndLine(); 
400                                         m->mothurOut("Running command: unique.seqs(" + inputString + ")"); m->mothurOutEndLine(); 
401                                                                  
402                                         Command* uniqueCommand = new DeconvoluteCommand(inputString);
403                                         uniqueCommand->execute();
404                                         
405                                         map<string, vector<string> > filenames = uniqueCommand->getOutputFiles();
406                                         
407                                         delete uniqueCommand;
408                                         
409                                         m->mothurOut("/******************************************/"); m->mothurOutEndLine(); 
410                                         
411                                         nameFile = filenames["name"][0];
412                                         fastaFileNames[s] = filenames["fasta"][0];
413                                 }
414                                 
415                                 //sort fastafile by abundance, returns new sorted fastafile name
416                                 m->mothurOut("Sorting fastafile according to abundance..."); cout.flush(); 
417                                 map<string, int> priority = sortFastaFile(fastaFileNames[s], nameFile);
418                                 m->mothurOut("Done."); m->mothurOutEndLine();
419                                 
420                                 if (m->control_pressed) {  for (int j = 0; j < outputNames.size(); j++) {       remove(outputNames[j].c_str()); }  return 0;    }
421
422                                 chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, trim, priority, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);   
423                         }
424                                 
425                         if (outputDir == "") { outputDir = m->hasPath(fastaFileNames[s]);  }//if user entered a file with a path then preserve it                               
426                         string outputFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s])) + "slayer.chimera";
427                         string accnosFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s]))  + "slayer.accnos";
428                         string trimFastaFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s]))  + "slayer.fasta";
429                         
430                         if (m->control_pressed) { delete chimera; for (int j = 0; j < outputNames.size(); j++) {        remove(outputNames[j].c_str()); }  return 0;    }
431                         
432                         if (chimera->getUnaligned()) { 
433                                 m->mothurOut("Your template sequences are different lengths, please correct."); m->mothurOutEndLine(); 
434                                 delete chimera;
435                                 return 0; 
436                         }
437                         templateSeqsLength = chimera->getLength();
438                         
439                 #ifdef USE_MPI  
440                         int pid, numSeqsPerProcessor; 
441                                 int tag = 2001;
442                                 vector<unsigned long int> MPIPos;
443                                 
444                                 MPI_Status status; 
445                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
446                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); 
447
448                                 MPI_File inMPI;
449                                 MPI_File outMPI;
450                                 MPI_File outMPIAccnos;
451                                 MPI_File outMPIFasta;
452                                 
453                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
454                                 int inMode=MPI_MODE_RDONLY; 
455                                 
456                                 char outFilename[1024];
457                                 strcpy(outFilename, outputFileName.c_str());
458                                 
459                                 char outAccnosFilename[1024];
460                                 strcpy(outAccnosFilename, accnosFileName.c_str());
461                         
462                                 char outFastaFilename[1024];
463                                 strcpy(outFastaFilename, trimFastaFileName.c_str());
464                                 
465                                 char inFileName[1024];
466                                 strcpy(inFileName, fastaFileNames[s].c_str());
467
468                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
469                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
470                                 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
471                                 if (trim) { MPI_File_open(MPI_COMM_WORLD, outFastaFilename, outMode, MPI_INFO_NULL, &outMPIFasta); }
472
473                         if (m->control_pressed) { outputTypes.clear();  MPI_File_close(&inMPI);  MPI_File_close(&outMPI); if (trim) {  MPI_File_close(&outMPIFasta);  } MPI_File_close(&outMPIAccnos); for (int j = 0; j < outputNames.size(); j++) {   remove(outputNames[j].c_str()); }   delete chimera; return 0;  }
474                         
475                                 if (pid == 0) { //you are the root process 
476                                         m->mothurOutEndLine();
477                                         m->mothurOut("Only reporting sequence supported by " + toString(minBS) + "% of bootstrapped results.");
478                                         m->mothurOutEndLine();
479                 
480                                         string outTemp = "Name\tLeftParent\tRightParent\tDivQLAQRB\tPerIDQLAQRB\tBootStrapA\tDivQLBQRA\tPerIDQLBQRA\tBootStrapB\tFlag\tLeftWindow\tRightWindow\n";
481                                         
482                                         //print header
483                                         int length = outTemp.length();
484                                         char* buf2 = new char[length];
485                                         memcpy(buf2, outTemp.c_str(), length);
486
487                                         MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &status);
488                                         delete buf2;
489
490                                         MPIPos = m->setFilePosFasta(fastaFileNames[s], numSeqs); //fills MPIPos, returns numSeqs
491                                         
492                                         if (templatefile != "self") { //if template=self we can only use 1 processor
493                                                 //send file positions to all processes
494                                                 for(int i = 1; i < processors; i++) { 
495                                                         MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
496                                                         MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
497                                                 }
498                                         }
499                                         //figure out how many sequences you have to align
500                                         numSeqsPerProcessor = numSeqs / processors;
501                                         int startIndex =  pid * numSeqsPerProcessor;
502                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
503                                         
504                                         if (templatefile == "self") { //if template=self we can only use 1 processor
505                                                 startIndex = 0;
506                                                 numSeqsPerProcessor = numSeqs;
507                                         }
508                                         
509                                         //do your part
510                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, outMPIFasta, MPIPos);
511                                         
512                                         int numNoParents = chimera->getNumNoParents();
513                                         int temp;
514                                         for(int i = 1; i < processors; i++) { 
515                                                 MPI_Recv(&temp, 1, MPI_INT, 1, tag, MPI_COMM_WORLD, &status);
516                                                 numNoParents += temp;
517                                         }
518                                         
519                                         
520                                         if (numSeqs == numNoParents) {  m->mothurOut("[WARNING]: megablast returned 0 potential parents for all your sequences. This could be due to formatdb.exe not being setup properly, please check formatdb.log for errors."); m->mothurOutEndLine(); }
521                                         
522                                         if (m->control_pressed) { outputTypes.clear();  MPI_File_close(&inMPI);  MPI_File_close(&outMPI); if (trim) { MPI_File_close(&outMPIFasta); }  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;  }
523
524                                 }else{ //you are a child process
525                                         if (templatefile != "self") { //if template=self we can only use 1 processor
526                                                 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
527                                                 MPIPos.resize(numSeqs+1);
528                                                 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
529                                         
530                                                 //figure out how many sequences you have to align
531                                                 numSeqsPerProcessor = numSeqs / processors;
532                                                 int startIndex =  pid * numSeqsPerProcessor;
533                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
534                                         
535                                                 //do your part
536                                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, outMPIFasta, MPIPos);
537                                                 
538                                                 int numNoParents = chimera->getNumNoParents();
539                                                 MPI_Send(&numNoParents, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
540                                         
541                                                 if (m->control_pressed) { outputTypes.clear();  MPI_File_close(&inMPI);  MPI_File_close(&outMPI); if (trim) { MPI_File_close(&outMPIFasta); }  MPI_File_close(&outMPIAccnos);  for (int j = 0; j < outputNames.size(); j++) {   remove(outputNames[j].c_str()); }  delete chimera; return 0;  }
542                                 
543                                         }
544                                 }
545                                 
546                                 //close files 
547                                 MPI_File_close(&inMPI);
548                                 MPI_File_close(&outMPI);
549                                 MPI_File_close(&outMPIAccnos); 
550                                 if (trim) { MPI_File_close(&outMPIFasta); }
551                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
552                                 
553                 #else
554                         ofstream outHeader;
555                         string tempHeader = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s])) + "slayer.chimeras.tempHeader";
556                         m->openOutputFile(tempHeader, outHeader);
557                         
558                         chimera->printHeader(outHeader);
559                         outHeader.close();
560                         
561                         vector<unsigned long int> positions = m->divideFile(fastaFileNames[s], processors);
562                                 
563                         for (int i = 0; i < (positions.size()-1); i++) {
564                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
565                         }       
566
567                         //break up file
568                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
569                                 if(processors == 1){
570                                         numSeqs = driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName);
571                                         
572                                         int numNoParents = chimera->getNumNoParents();
573                                         if (numNoParents == numSeqs) { m->mothurOut("[WARNING]: megablast returned 0 potential parents for all your sequences. This could be due to formatdb.exe not being setup properly, please check formatdb.log for errors."); m->mothurOutEndLine(); }
574                                         
575                                         if (m->control_pressed) { outputTypes.clear(); if (trim) { remove(trimFastaFileName.c_str()); } remove(outputFileName.c_str()); remove(tempHeader.c_str()); remove(accnosFileName.c_str()); for (int j = 0; j < outputNames.size(); j++) {      remove(outputNames[j].c_str()); } for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear(); delete chimera; return 0; }
576                                         
577                                 }else{
578                                         processIDS.resize(0);
579                                         
580                                         numSeqs = createProcesses(outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName); 
581                                 
582                                         rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
583                                         rename((accnosFileName + toString(processIDS[0]) + ".temp").c_str(), accnosFileName.c_str());
584                                         if (trim) {  rename((trimFastaFileName + toString(processIDS[0]) + ".temp").c_str(), trimFastaFileName.c_str()); }
585                                                 
586                                         //append output files
587                                         for(int i=1;i<processors;i++){
588                                                 m->appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
589                                                 remove((outputFileName + toString(processIDS[i]) + ".temp").c_str());
590                                         }
591                                         
592                                         //append output files
593                                         for(int i=1;i<processors;i++){
594                                                 m->appendFiles((accnosFileName + toString(processIDS[i]) + ".temp"), accnosFileName);
595                                                 remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());
596                                         }
597                                         
598                                         if (trim) {
599                                                 for(int i=1;i<processors;i++){
600                                                         m->appendFiles((trimFastaFileName + toString(processIDS[i]) + ".temp"), trimFastaFileName);
601                                                         remove((trimFastaFileName + toString(processIDS[i]) + ".temp").c_str());
602                                                 }
603                                         }
604                                         
605                                         if (m->control_pressed) { outputTypes.clear(); if (trim) { remove(trimFastaFileName.c_str()); } remove(outputFileName.c_str()); remove(accnosFileName.c_str()); for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear(); delete chimera; return 0; }
606                                 }
607
608                         #else
609                                 numSeqs = driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName);
610                                 
611                                 int numNoParents = chimera->getNumNoParents();
612                                 if (numNoParents == numSeqs) { m->mothurOut("[WARNING]: megablast returned 0 potential parents for all your sequences. This could be due to formatdb.exe not being setup properly, please check formatdb.log for errors."); m->mothurOutEndLine(); }
613
614                                 
615                                 if (m->control_pressed) { outputTypes.clear(); if (trim) { remove(trimFastaFileName.c_str()); } remove(outputFileName.c_str()); remove(tempHeader.c_str()); remove(accnosFileName.c_str()); for (int j = 0; j < outputNames.size(); j++) {      remove(outputNames[j].c_str()); } for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear(); delete chimera; return 0; }
616                                 
617                         #endif
618                         
619                         m->appendFiles(outputFileName, tempHeader);
620                 
621                         remove(outputFileName.c_str());
622                         rename(tempHeader.c_str(), outputFileName.c_str());
623                         
624                 #endif
625                         delete chimera;
626                         
627                         
628                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
629                         
630                         outputNames.push_back(outputFileName); outputTypes["chimera"].push_back(outputFileName);
631                         outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName);
632                         if (trim) {  outputNames.push_back(trimFastaFileName); outputTypes["fasta"].push_back(trimFastaFileName); }
633                         
634                         m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine();
635                 }
636                 
637                 //set accnos file as new current accnosfile
638                 string current = "";
639                 itTypes = outputTypes.find("accnos");
640                 if (itTypes != outputTypes.end()) {
641                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setAccnosFile(current); }
642                 }
643                 
644                 if (trim) {
645                         itTypes = outputTypes.find("fasta");
646                         if (itTypes != outputTypes.end()) {
647                                 if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
648                         }
649                 }
650                 
651                 m->mothurOutEndLine();
652                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
653                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }       
654                 m->mothurOutEndLine();
655
656                 return 0;
657                 
658         }
659         catch(exception& e) {
660                 m->errorOut(e, "ChimeraSlayerCommand", "execute");
661                 exit(1);
662         }
663 }
664 //**********************************************************************************************************************
665
666 int ChimeraSlayerCommand::driver(linePair* filePos, string outputFName, string filename, string accnos, string fasta){
667         try {
668                 ofstream out;
669                 m->openOutputFile(outputFName, out);
670                 
671                 ofstream out2;
672                 m->openOutputFile(accnos, out2);
673                 
674                 ofstream out3;
675                 if (trim) {  m->openOutputFile(fasta, out3); }
676                 
677                 ifstream inFASTA;
678                 m->openInputFile(filename, inFASTA);
679
680                 inFASTA.seekg(filePos->start);
681
682                 bool done = false;
683                 int count = 0;
684         
685                 while (!done) {
686                 
687                         if (m->control_pressed) {       out.close(); out2.close(); if (trim) { out3.close(); } inFASTA.close(); return 1;       }
688                 
689                         Sequence* candidateSeq = new Sequence(inFASTA);  m->gobble(inFASTA);
690                         string candidateAligned = candidateSeq->getAligned();
691                         
692                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
693                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
694                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
695                                 }else{
696                                         //find chimeras
697                                         chimera->getChimeras(candidateSeq);
698                                         
699                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
700                                                 
701                                         //if you are not chimeric, then check each half
702                                         data_results wholeResults = chimera->getResults();
703                                         
704                                         //determine if we need to split
705                                         bool isChimeric = false;
706                                         
707                                         if (wholeResults.flag == "yes") {
708                                                 string chimeraFlag = "no";
709                                                 if(  (wholeResults.results[0].bsa >= minBS && wholeResults.results[0].divr_qla_qrb >= divR)
710                                                    ||
711                                                    (wholeResults.results[0].bsb >= minBS && wholeResults.results[0].divr_qlb_qra >= divR) ) { chimeraFlag = "yes"; }
712                                                 
713                                                 
714                                                 if (chimeraFlag == "yes") {     
715                                                         if ((wholeResults.results[0].bsa >= minBS) || (wholeResults.results[0].bsb >= minBS)) { isChimeric = true; }
716                                                 }
717                                         }
718                                         
719                                         if ((!isChimeric) && trimera) {
720                                                 
721                                                 //split sequence in half by bases
722                                                 string leftQuery, rightQuery;
723                                                 Sequence tempSeq(candidateSeq->getName(), candidateAligned);
724                                                 divideInHalf(tempSeq, leftQuery, rightQuery);
725                                                 
726                                                 //run chimeraSlayer on each piece
727                                                 Sequence* left = new Sequence(candidateSeq->getName(), leftQuery);
728                                                 Sequence* right = new Sequence(candidateSeq->getName(), rightQuery);
729                                                 
730                                                 //find chimeras
731                                                 chimera->getChimeras(left);
732                                                 data_results leftResults = chimera->getResults();
733                                                 
734                                                 chimera->getChimeras(right);
735                                                 data_results rightResults = chimera->getResults();
736                                                 
737                                                 //if either piece is chimeric then report
738                                                 Sequence trimmed = chimera->print(out, out2, leftResults, rightResults);
739                                                 if (trim) { trimmed.printSequence(out3);  }
740                                                 
741                                                 delete left; delete right;
742                                                 
743                                         }else { //already chimeric
744                                                 //print results
745                                                 Sequence trimmed = chimera->print(out, out2);
746                                                 if (trim) { trimmed.printSequence(out3);  }
747                                         }
748                                         
749                                         
750                                 }
751                                 count++;
752                         }
753                         
754                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
755                                 unsigned long int pos = inFASTA.tellg();
756                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
757                         #else
758                                 if (inFASTA.eof()) { break; }
759                         #endif
760                         
761                         delete candidateSeq;
762                         //report progress
763                         if((count) % 100 == 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
764                 }
765                 //report progress
766                 if((count) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
767                 
768                 out.close();
769                 out2.close();
770                 if (trim) { out3.close(); }
771                 inFASTA.close();
772                                 
773                 return count;
774         }
775         catch(exception& e) {
776                 m->errorOut(e, "ChimeraSlayerCommand", "driver");
777                 exit(1);
778         }
779 }
780 //**********************************************************************************************************************
781 #ifdef USE_MPI
782 int ChimeraSlayerCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, MPI_File& outAccMPI, MPI_File& outFastaMPI, vector<unsigned long int>& MPIPos){
783         try {                           
784                 MPI_Status status; 
785                 int pid;
786                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
787                 
788                 for(int i=0;i<num;i++){
789                         
790                         if (m->control_pressed) {       return 1;       }
791                         
792                         //read next sequence
793                         int length = MPIPos[start+i+1] - MPIPos[start+i];
794
795                         char* buf4 = new char[length];
796                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
797         
798                         string tempBuf = buf4;
799                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
800                         istringstream iss (tempBuf,istringstream::in);
801
802                         delete buf4;
803
804                         Sequence* candidateSeq = new Sequence(iss);  m->gobble(iss);
805                         string candidateAligned = candidateSeq->getAligned();
806                 
807                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
808                                 
809                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
810                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
811                                 }else{
812                 
813                                         //find chimeras
814                                         chimera->getChimeras(candidateSeq);
815                         
816                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
817                                         
818                                         //if you are not chimeric, then check each half
819                                         data_results wholeResults = chimera->getResults();
820                                         
821                                         //determine if we need to split
822                                         bool isChimeric = false;
823                                         
824                                         if (wholeResults.flag == "yes") {
825                                                 string chimeraFlag = "no";
826                                                 if(  (wholeResults.results[0].bsa >= minBS && wholeResults.results[0].divr_qla_qrb >= divR)
827                                                    ||
828                                                    (wholeResults.results[0].bsb >= minBS && wholeResults.results[0].divr_qlb_qra >= divR) ) { chimeraFlag = "yes"; }
829                                                 
830                                                 
831                                                 if (chimeraFlag == "yes") {     
832                                                         if ((wholeResults.results[0].bsa >= minBS) || (wholeResults.results[0].bsb >= minBS)) { isChimeric = true; }
833                                                 }
834                                         }
835                                         
836                                         if ((!isChimeric) && trimera) {                                                 
837                                                 //split sequence in half by bases
838                                                 string leftQuery, rightQuery;
839                                                 Sequence tempSeq(candidateSeq->getName(), candidateAligned);
840                                                 divideInHalf(tempSeq, leftQuery, rightQuery);
841                                                 
842                                                 //run chimeraSlayer on each piece
843                                                 Sequence* left = new Sequence(candidateSeq->getName(), leftQuery);
844                                                 Sequence* right = new Sequence(candidateSeq->getName(), rightQuery);
845                                                 
846                                                 //find chimeras
847                                                 chimera->getChimeras(left);
848                                                 data_results leftResults = chimera->getResults();
849                                                 
850                                                 chimera->getChimeras(right);
851                                                 data_results rightResults = chimera->getResults();
852                                                 
853                                                 //if either piece is chimeric then report
854                                                 Sequence trimmed = chimera->print(outMPI, outAccMPI, leftResults, rightResults);
855                                                 if (trim) {  
856                                                         string outputString = ">" + trimmed.getName() + "\n" + trimmed.getAligned() + "\n";
857                                                         
858                                                         //write to accnos file
859                                                         int length = outputString.length();
860                                                         char* buf2 = new char[length];
861                                                         memcpy(buf2, outputString.c_str(), length);
862                                                         
863                                                         MPI_File_write_shared(outFastaMPI, buf2, length, MPI_CHAR, &status);
864                                                         delete buf2;
865                                                 }
866                                                 
867                                                 delete left; delete right;
868                                                 
869                                         }else { 
870                                                 //print results
871                                                 Sequence trimmed = chimera->print(outMPI, outAccMPI);
872                                                 
873                                                 if (trim) {  
874                                                         string outputString = ">" + trimmed.getName() + "\n" + trimmed.getAligned() + "\n";
875                                                         
876                                                         //write to accnos file
877                                                         int length = outputString.length();
878                                                         char* buf2 = new char[length];
879                                                         memcpy(buf2, outputString.c_str(), length);
880                                                         
881                                                         MPI_File_write_shared(outFastaMPI, buf2, length, MPI_CHAR, &status);
882                                                         delete buf2;
883                                                 }
884                                         }
885                                         
886                                 }
887                         }
888                         delete candidateSeq;
889                         
890                         //report progress
891                         if((i+1) % 100 == 0){  cout << "Processing sequence: " << (i+1) << endl;        m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n");          }
892                 }
893                 //report progress
894                 if(num % 100 != 0){             cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n");  }
895                 
896                                 
897                 return 0;
898         }
899         catch(exception& e) {
900                 m->errorOut(e, "ChimeraSlayerCommand", "driverMPI");
901                 exit(1);
902         }
903 }
904 #endif
905
906 /**************************************************************************************************/
907
908 int ChimeraSlayerCommand::createProcesses(string outputFileName, string filename, string accnos, string fasta) {
909         try {
910 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
911                 int process = 0;
912                 int num = 0;
913                 
914                 //loop through and create all the processes you want
915                 while (process != processors) {
916                         int pid = fork();
917                         
918                         if (pid > 0) {
919                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
920                                 process++;
921                         }else if (pid == 0){
922                                 num = driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename, accnos + toString(getpid()) + ".temp", fasta + toString(getpid()) + ".temp");
923                                 
924                                 //pass numSeqs to parent
925                                 ofstream out;
926                                 string tempFile = outputFileName + toString(getpid()) + ".num.temp";
927                                 m->openOutputFile(tempFile, out);
928                                 out << num << '\t' << chimera->getNumNoParents() << endl;
929                                 out.close();
930                                 exit(0);
931                         }else { 
932                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
933                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
934                                 exit(0);
935                         }
936                 }
937                 
938                 //force parent to wait until all the processes are done
939                 for (int i=0;i<processors;i++) { 
940                         int temp = processIDS[i];
941                         wait(&temp);
942                 }
943                 
944                 int numNoParents = 0;
945                 for (int i = 0; i < processIDS.size(); i++) {
946                         ifstream in;
947                         string tempFile =  outputFileName + toString(processIDS[i]) + ".num.temp";
948                         m->openInputFile(tempFile, in);
949                         if (!in.eof()) { int tempNum = 0; int tempNumParents = 0; in >> tempNum >> tempNumParents; num += tempNum; numNoParents += tempNumParents; }
950                         in.close(); remove(tempFile.c_str());
951                 }
952                 
953                 if (num == numNoParents) {  m->mothurOut("[WARNING]: megablast returned 0 potential parents for all your sequences. This could be due to formatdb.exe not being setup properly, please check formatdb.log for errors."); m->mothurOutEndLine(); }
954                 
955                 return num;
956 #endif          
957         }
958         catch(exception& e) {
959                 m->errorOut(e, "ChimeraSlayerCommand", "createProcesses");
960                 exit(1);
961         }
962 }
963
964 /**************************************************************************************************/
965
966 int ChimeraSlayerCommand::divideInHalf(Sequence querySeq, string& leftQuery, string& rightQuery) {
967         try {
968                 
969                 string queryUnAligned = querySeq.getUnaligned();
970                 int numBases = int(queryUnAligned.length() * 0.5);
971                 
972                 string queryAligned = querySeq.getAligned();
973                 leftQuery = querySeq.getAligned();
974                 rightQuery = querySeq.getAligned();
975                 
976                 int baseCount = 0;
977                 int leftSpot = 0;
978                 for (int i = 0; i < queryAligned.length(); i++) {
979                         //if you are a base
980                         if (isalpha(queryAligned[i])) {         
981                                 baseCount++; 
982                         }
983                         
984                         //if you have half
985                         if (baseCount >= numBases) {  leftSpot = i; break; } //first half
986                 }
987                 
988                 //blank out right side
989                 for (int i = leftSpot; i < leftQuery.length(); i++) { leftQuery[i] = '.'; }
990                 
991                 //blank out left side
992                 for (int i = 0; i < leftSpot; i++) { rightQuery[i] = '.'; }
993                 
994                 return 0;
995                 
996         }
997         catch(exception& e) {
998                 m->errorOut(e, "ChimeraSlayerCommand", "divideInHalf");
999                 exit(1);
1000         }
1001 }
1002 /**************************************************************************************************/
1003 map<string, int> ChimeraSlayerCommand::sortFastaFile(string fastaFile, string nameFile) {
1004         try {
1005                 map<string, int> nameAbund;
1006                 
1007                 //read through fastafile and store info
1008                 map<string, string> seqs;
1009                 ifstream in;
1010                 m->openInputFile(fastaFile, in);
1011                 
1012                 while (!in.eof()) {
1013                         
1014                         if (m->control_pressed) { in.close(); return nameAbund; }
1015                         
1016                         Sequence seq(in); m->gobble(in);
1017                         seqs[seq.getName()] = seq.getAligned();
1018                 }
1019                 
1020                 in.close();
1021                 
1022                 //read namefile
1023                 vector<seqPriorityNode> nameMapCount;
1024                 int error = m->readNames(nameFile, nameMapCount, seqs);
1025                 
1026                 if (m->control_pressed) { return nameAbund; }
1027                 
1028                 if (error == 1) { m->control_pressed = true; return nameAbund; }
1029                 if (seqs.size() != nameMapCount.size()) { m->mothurOut( "The number of sequences in your fastafile does not match the number of sequences in your namefile, aborting."); m->mothurOutEndLine(); m->control_pressed = true; return nameAbund; }
1030
1031                 sort(nameMapCount.begin(), nameMapCount.end(), compareSeqPriorityNodes);
1032                 
1033                 string newFasta = fastaFile + ".temp";
1034                 ofstream out;
1035                 m->openOutputFile(newFasta, out);
1036                 
1037                 //print new file in order of
1038                 for (int i = 0; i < nameMapCount.size(); i++) {
1039                         out << ">" << nameMapCount[i].name << endl << nameMapCount[i].seq << endl;
1040                         nameAbund[nameMapCount[i].name] = nameMapCount[i].numIdentical;
1041                 }
1042                 out.close();
1043                 
1044                 rename(newFasta.c_str(), fastaFile.c_str());
1045                                 
1046                 return nameAbund;
1047                 
1048         }
1049         catch(exception& e) {
1050                 m->errorOut(e, "ChimeraSlayerCommand", "sortFastaFile");
1051                 exit(1);
1052         }
1053 }
1054
1055 /**************************************************************************************************/
1056