]> git.donarmstrong.com Git - mothur.git/blob - chimeraslayercommand.cpp
working on current change
[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                                                 }else {
209                                                         m->setFastaFile(fastaFileNames[i]);
210                                                 }
211                                         }
212                                 }
213                                 
214                                 //make sure there is at least one valid file left
215                                 if (fastaFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid files."); m->mothurOutEndLine(); abort = true; }
216                         }
217                         
218                         
219                         //check for required parameters
220                         bool hasName = true;
221                         namefile = validParameter.validFile(parameters, "name", false);
222                         if (namefile == "not found") { namefile = "";  hasName = false; }
223                         else { 
224                                 m->splitAtDash(namefile, nameFileNames);
225                                 
226                                 //go through files and make sure they are good, if not, then disregard them
227                                 for (int i = 0; i < nameFileNames.size(); i++) {
228                                         
229                                         bool ignore = false;
230                                         if (nameFileNames[i] == "current") { 
231                                                 nameFileNames[i] = m->getNameFile(); 
232                                                 if (nameFileNames[i] != "") {  m->mothurOut("Using " + nameFileNames[i] + " as input file for the name parameter where you had given current."); m->mothurOutEndLine(); }
233                                                 else {  
234                                                         m->mothurOut("You have no current namefile, ignoring current."); m->mothurOutEndLine(); ignore=true; 
235                                                         //erase from file list
236                                                         nameFileNames.erase(nameFileNames.begin()+i);
237                                                         i--;
238                                                 }
239                                         }
240                                         
241                                         if (!ignore) {
242                                                 
243                                                 if (inputDir != "") {
244                                                         string path = m->hasPath(nameFileNames[i]);
245                                                         //if the user has not given a path then, add inputdir. else leave path alone.
246                                                         if (path == "") {       nameFileNames[i] = inputDir + nameFileNames[i];         }
247                                                 }
248                                                 
249                                                 int ableToOpen;
250                                                 ifstream in;
251                                                 
252                                                 ableToOpen = m->openInputFile(nameFileNames[i], in, "noerror");
253                                                 
254                                                 //if you can't open it, try default location
255                                                 if (ableToOpen == 1) {
256                                                         if (m->getDefaultPath() != "") { //default path is set
257                                                                 string tryPath = m->getDefaultPath() + m->getSimpleName(nameFileNames[i]);
258                                                                 m->mothurOut("Unable to open " + nameFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
259                                                                 ifstream in2;
260                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
261                                                                 in2.close();
262                                                                 nameFileNames[i] = tryPath;
263                                                         }
264                                                 }
265                                                 
266                                                 if (ableToOpen == 1) {
267                                                         if (m->getOutputDir() != "") { //default path is set
268                                                                 string tryPath = m->getOutputDir() + m->getSimpleName(nameFileNames[i]);
269                                                                 m->mothurOut("Unable to open " + nameFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
270                                                                 ifstream in2;
271                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
272                                                                 in2.close();
273                                                                 nameFileNames[i] = tryPath;
274                                                         }
275                                                 }
276                                                 
277                                                 in.close();
278                                                 
279                                                 if (ableToOpen == 1) { 
280                                                         m->mothurOut("Unable to open " + nameFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
281                                                         //erase from file list
282                                                         nameFileNames.erase(nameFileNames.begin()+i);
283                                                         i--;
284                                                 }else {
285                                                         m->setNameFile(nameFileNames[i]);
286                                                 }
287                                         }
288                                 }
289                                 
290                                 //make sure there is at least one valid file left
291                                 if (nameFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid name files."); m->mothurOutEndLine(); abort = true; }
292                         }
293                         
294                         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; }
295                         
296                         //if the user changes the output directory command factory will send this info to us in the output parameter 
297                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  outputDir = ""; }
298                         
299                         
300                         string path;
301                         it = parameters.find("reference");
302                         //user has given a template file
303                         if(it != parameters.end()){ 
304                                 if (it->second == "self") { templatefile = "self"; }
305                                 else {
306                                         path = m->hasPath(it->second);
307                                         //if the user has not given a path then, add inputdir. else leave path alone.
308                                         if (path == "") {       parameters["reference"] = inputDir + it->second;                }
309                                         
310                                         templatefile = validParameter.validFile(parameters, "reference", true);
311                                         if (templatefile == "not open") { abort = true; }
312                                         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;  }    
313                                 }
314                         }else if (hasName) {  templatefile = "self"; }
315                         else { templatefile = "";  m->mothurOut("reference is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true;  }   
316                         
317                         string temp = validParameter.validFile(parameters, "processors", false);        if (temp == "not found"){       temp = m->getProcessors();      }
318                         m->setProcessors(temp);
319                         convert(temp, processors);
320                         
321                         temp = validParameter.validFile(parameters, "ksize", false);                    if (temp == "not found") { temp = "7"; }
322                         convert(temp, ksize);
323                                                 
324                         temp = validParameter.validFile(parameters, "window", false);                   if (temp == "not found") { temp = "50"; }                       
325                         convert(temp, window);
326                         
327                         temp = validParameter.validFile(parameters, "match", false);                    if (temp == "not found") { temp = "5"; }
328                         convert(temp, match);
329                         
330                         temp = validParameter.validFile(parameters, "mismatch", false);                 if (temp == "not found") { temp = "-4"; }
331                         convert(temp, mismatch);
332                         
333                         temp = validParameter.validFile(parameters, "divergence", false);               if (temp == "not found") { temp = "1.007"; }
334                         convert(temp, divR);
335                         
336                         temp = validParameter.validFile(parameters, "minsim", false);                   if (temp == "not found") { temp = "90"; }
337                         convert(temp, minSimilarity);
338                         
339                         temp = validParameter.validFile(parameters, "mincov", false);                   if (temp == "not found") { temp = "70"; }
340                         convert(temp, minCoverage);
341                         
342                         temp = validParameter.validFile(parameters, "minbs", false);                    if (temp == "not found") { temp = "90"; }
343                         convert(temp, minBS);
344                         
345                         temp = validParameter.validFile(parameters, "minsnp", false);                   if (temp == "not found") { temp = "10"; }
346                         convert(temp, minSNP);
347
348                         temp = validParameter.validFile(parameters, "parents", false);                  if (temp == "not found") { temp = "3"; }
349                         convert(temp, parents); 
350                         
351                         temp = validParameter.validFile(parameters, "realign", false);                  if (temp == "not found") { temp = "t"; }
352                         realign = m->isTrue(temp); 
353                         
354                         temp = validParameter.validFile(parameters, "trim", false);                             if (temp == "not found") { temp = "f"; }
355                         trim = m->isTrue(temp); 
356                         
357                         temp = validParameter.validFile(parameters, "split", false);                    if (temp == "not found") { temp = "f"; }
358                         trimera = m->isTrue(temp); 
359                         
360                         search = validParameter.validFile(parameters, "search", false);                 if (search == "not found") { search = "blast"; }
361                         
362                         temp = validParameter.validFile(parameters, "iters", false);                    if (temp == "not found") { temp = "1000"; }             
363                         convert(temp, iters); 
364                          
365                         temp = validParameter.validFile(parameters, "increment", false);                if (temp == "not found") { temp = "5"; }
366                         convert(temp, increment);
367                         
368                         temp = validParameter.validFile(parameters, "numwanted", false);                if (temp == "not found") { temp = "15"; }               
369                         convert(temp, numwanted);
370
371                         if ((search != "blast") && (search != "kmer")) { m->mothurOut(search + " is not a valid search."); m->mothurOutEndLine(); abort = true;  }
372                 }
373         }
374         catch(exception& e) {
375                 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
376                 exit(1);
377         }
378 }
379 //***************************************************************************************************************
380
381 int ChimeraSlayerCommand::execute(){
382         try{
383                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
384                 
385                 for (int s = 0; s < fastaFileNames.size(); s++) {
386                                 
387                         m->mothurOut("Checking sequences from " + fastaFileNames[s] + " ..." ); m->mothurOutEndLine();
388                 
389                         int start = time(NULL); 
390                         
391                         if (templatefile != "self") { //you want to run slayer with a refernce template
392                                 chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, trim, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);     
393                         }else {
394                                 if (processors != 1) { m->mothurOut("When using template=self, mothur can only use 1 processor, continuing."); m->mothurOutEndLine(); processors = 1; }
395                                 string nameFile = "";
396                                 if (nameFileNames.size() != 0) { //you provided a namefile and we don't need to create one
397                                         nameFile = nameFileNames[s];
398                                 }else {
399                                         m->mothurOutEndLine(); m->mothurOut("No namesfile given, running unique.seqs command to generate one."); m->mothurOutEndLine(); m->mothurOutEndLine();
400                                         
401                                         //use unique.seqs to create new name and fastafile
402                                         string inputString = "fasta=" + fastaFileNames[s];
403                                         m->mothurOut("/******************************************/"); m->mothurOutEndLine(); 
404                                         m->mothurOut("Running command: unique.seqs(" + inputString + ")"); m->mothurOutEndLine(); 
405                                                                  
406                                         Command* uniqueCommand = new DeconvoluteCommand(inputString);
407                                         uniqueCommand->execute();
408                                         
409                                         map<string, vector<string> > filenames = uniqueCommand->getOutputFiles();
410                                         
411                                         delete uniqueCommand;
412                                         
413                                         m->mothurOut("/******************************************/"); m->mothurOutEndLine(); 
414                                         
415                                         nameFile = filenames["name"][0];
416                                         fastaFileNames[s] = filenames["fasta"][0];
417                                 }
418                                 
419                                 //sort fastafile by abundance, returns new sorted fastafile name
420                                 m->mothurOut("Sorting fastafile according to abundance..."); cout.flush(); 
421                                 map<string, int> priority = sortFastaFile(fastaFileNames[s], nameFile);
422                                 m->mothurOut("Done."); m->mothurOutEndLine();
423                                 
424                                 if (m->control_pressed) {  for (int j = 0; j < outputNames.size(); j++) {       remove(outputNames[j].c_str()); }  return 0;    }
425
426                                 chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, trim, priority, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);   
427                         }
428                                 
429                         if (outputDir == "") { outputDir = m->hasPath(fastaFileNames[s]);  }//if user entered a file with a path then preserve it                               
430                         string outputFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s])) + "slayer.chimera";
431                         string accnosFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s]))  + "slayer.accnos";
432                         string trimFastaFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s]))  + "slayer.fasta";
433                         
434                         if (m->control_pressed) { delete chimera; for (int j = 0; j < outputNames.size(); j++) {        remove(outputNames[j].c_str()); }  return 0;    }
435                         
436                         if (chimera->getUnaligned()) { 
437                                 m->mothurOut("Your template sequences are different lengths, please correct."); m->mothurOutEndLine(); 
438                                 delete chimera;
439                                 return 0; 
440                         }
441                         templateSeqsLength = chimera->getLength();
442                         
443                 #ifdef USE_MPI  
444                         int pid, numSeqsPerProcessor; 
445                                 int tag = 2001;
446                                 vector<unsigned long int> MPIPos;
447                                 
448                                 MPI_Status status; 
449                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
450                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); 
451
452                                 MPI_File inMPI;
453                                 MPI_File outMPI;
454                                 MPI_File outMPIAccnos;
455                                 MPI_File outMPIFasta;
456                                 
457                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
458                                 int inMode=MPI_MODE_RDONLY; 
459                                 
460                                 char outFilename[1024];
461                                 strcpy(outFilename, outputFileName.c_str());
462                                 
463                                 char outAccnosFilename[1024];
464                                 strcpy(outAccnosFilename, accnosFileName.c_str());
465                         
466                                 char outFastaFilename[1024];
467                                 strcpy(outFastaFilename, trimFastaFileName.c_str());
468                                 
469                                 char inFileName[1024];
470                                 strcpy(inFileName, fastaFileNames[s].c_str());
471
472                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
473                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
474                                 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
475                                 if (trim) { MPI_File_open(MPI_COMM_WORLD, outFastaFilename, outMode, MPI_INFO_NULL, &outMPIFasta); }
476
477                         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;  }
478                         
479                                 if (pid == 0) { //you are the root process 
480                                         m->mothurOutEndLine();
481                                         m->mothurOut("Only reporting sequence supported by " + toString(minBS) + "% of bootstrapped results.");
482                                         m->mothurOutEndLine();
483                 
484                                         string outTemp = "Name\tLeftParent\tRightParent\tDivQLAQRB\tPerIDQLAQRB\tBootStrapA\tDivQLBQRA\tPerIDQLBQRA\tBootStrapB\tFlag\tLeftWindow\tRightWindow\n";
485                                         
486                                         //print header
487                                         int length = outTemp.length();
488                                         char* buf2 = new char[length];
489                                         memcpy(buf2, outTemp.c_str(), length);
490
491                                         MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &status);
492                                         delete buf2;
493
494                                         MPIPos = m->setFilePosFasta(fastaFileNames[s], numSeqs); //fills MPIPos, returns numSeqs
495                                         
496                                         if (templatefile != "self") { //if template=self we can only use 1 processor
497                                                 //send file positions to all processes
498                                                 for(int i = 1; i < processors; i++) { 
499                                                         MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
500                                                         MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
501                                                 }
502                                         }
503                                         //figure out how many sequences you have to align
504                                         numSeqsPerProcessor = numSeqs / processors;
505                                         int startIndex =  pid * numSeqsPerProcessor;
506                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
507                                         
508                                         if (templatefile == "self") { //if template=self we can only use 1 processor
509                                                 startIndex = 0;
510                                                 numSeqsPerProcessor = numSeqs;
511                                         }
512                                         
513                                         //do your part
514                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, outMPIFasta, MPIPos);
515                                         
516                                         int numNoParents = chimera->getNumNoParents();
517                                         int temp;
518                                         for(int i = 1; i < processors; i++) { 
519                                                 MPI_Recv(&temp, 1, MPI_INT, 1, tag, MPI_COMM_WORLD, &status);
520                                                 numNoParents += temp;
521                                         }
522                                         
523                                         
524                                         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(); }
525                                         
526                                         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;  }
527
528                                 }else{ //you are a child process
529                                         if (templatefile != "self") { //if template=self we can only use 1 processor
530                                                 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
531                                                 MPIPos.resize(numSeqs+1);
532                                                 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
533                                         
534                                                 //figure out how many sequences you have to align
535                                                 numSeqsPerProcessor = numSeqs / processors;
536                                                 int startIndex =  pid * numSeqsPerProcessor;
537                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
538                                         
539                                                 //do your part
540                                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, outMPIFasta, MPIPos);
541                                                 
542                                                 int numNoParents = chimera->getNumNoParents();
543                                                 MPI_Send(&numNoParents, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
544                                         
545                                                 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;  }
546                                 
547                                         }
548                                 }
549                                 
550                                 //close files 
551                                 MPI_File_close(&inMPI);
552                                 MPI_File_close(&outMPI);
553                                 MPI_File_close(&outMPIAccnos); 
554                                 if (trim) { MPI_File_close(&outMPIFasta); }
555                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
556                                 
557                 #else
558                         ofstream outHeader;
559                         string tempHeader = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s])) + "slayer.chimeras.tempHeader";
560                         m->openOutputFile(tempHeader, outHeader);
561                         
562                         chimera->printHeader(outHeader);
563                         outHeader.close();
564                         
565                         vector<unsigned long int> positions = m->divideFile(fastaFileNames[s], processors);
566                                 
567                         for (int i = 0; i < (positions.size()-1); i++) {
568                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
569                         }       
570
571                         //break up file
572                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
573                                 if(processors == 1){
574                                         numSeqs = driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName);
575                                         
576                                         int numNoParents = chimera->getNumNoParents();
577                                         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(); }
578                                         
579                                         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; }
580                                         
581                                 }else{
582                                         processIDS.resize(0);
583                                         
584                                         numSeqs = createProcesses(outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName); 
585                                 
586                                         rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
587                                         rename((accnosFileName + toString(processIDS[0]) + ".temp").c_str(), accnosFileName.c_str());
588                                         if (trim) {  rename((trimFastaFileName + toString(processIDS[0]) + ".temp").c_str(), trimFastaFileName.c_str()); }
589                                                 
590                                         //append output files
591                                         for(int i=1;i<processors;i++){
592                                                 m->appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
593                                                 remove((outputFileName + toString(processIDS[i]) + ".temp").c_str());
594                                         }
595                                         
596                                         //append output files
597                                         for(int i=1;i<processors;i++){
598                                                 m->appendFiles((accnosFileName + toString(processIDS[i]) + ".temp"), accnosFileName);
599                                                 remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());
600                                         }
601                                         
602                                         if (trim) {
603                                                 for(int i=1;i<processors;i++){
604                                                         m->appendFiles((trimFastaFileName + toString(processIDS[i]) + ".temp"), trimFastaFileName);
605                                                         remove((trimFastaFileName + toString(processIDS[i]) + ".temp").c_str());
606                                                 }
607                                         }
608                                         
609                                         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; }
610                                 }
611
612                         #else
613                                 numSeqs = driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName);
614                                 
615                                 int numNoParents = chimera->getNumNoParents();
616                                 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(); }
617
618                                 
619                                 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; }
620                                 
621                         #endif
622                         
623                         m->appendFiles(outputFileName, tempHeader);
624                 
625                         remove(outputFileName.c_str());
626                         rename(tempHeader.c_str(), outputFileName.c_str());
627                         
628                 #endif
629                         delete chimera;
630                         
631                         
632                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
633                         
634                         outputNames.push_back(outputFileName); outputTypes["chimera"].push_back(outputFileName);
635                         outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName);
636                         if (trim) {  outputNames.push_back(trimFastaFileName); outputTypes["fasta"].push_back(trimFastaFileName); }
637                         
638                         m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine();
639                 }
640                 
641                 //set accnos file as new current accnosfile
642                 string current = "";
643                 itTypes = outputTypes.find("accnos");
644                 if (itTypes != outputTypes.end()) {
645                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setAccnosFile(current); }
646                 }
647                 
648                 if (trim) {
649                         itTypes = outputTypes.find("fasta");
650                         if (itTypes != outputTypes.end()) {
651                                 if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
652                         }
653                 }
654                 
655                 m->mothurOutEndLine();
656                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
657                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }       
658                 m->mothurOutEndLine();
659
660                 return 0;
661                 
662         }
663         catch(exception& e) {
664                 m->errorOut(e, "ChimeraSlayerCommand", "execute");
665                 exit(1);
666         }
667 }
668 //**********************************************************************************************************************
669
670 int ChimeraSlayerCommand::driver(linePair* filePos, string outputFName, string filename, string accnos, string fasta){
671         try {
672                 ofstream out;
673                 m->openOutputFile(outputFName, out);
674                 
675                 ofstream out2;
676                 m->openOutputFile(accnos, out2);
677                 
678                 ofstream out3;
679                 if (trim) {  m->openOutputFile(fasta, out3); }
680                 
681                 ifstream inFASTA;
682                 m->openInputFile(filename, inFASTA);
683
684                 inFASTA.seekg(filePos->start);
685
686                 bool done = false;
687                 int count = 0;
688         
689                 while (!done) {
690                 
691                         if (m->control_pressed) {       out.close(); out2.close(); if (trim) { out3.close(); } inFASTA.close(); return 1;       }
692                 
693                         Sequence* candidateSeq = new Sequence(inFASTA);  m->gobble(inFASTA);
694                         string candidateAligned = candidateSeq->getAligned();
695                         
696                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
697                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
698                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
699                                 }else{
700                                         //find chimeras
701                                         chimera->getChimeras(candidateSeq);
702                                         
703                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
704                                                 
705                                         //if you are not chimeric, then check each half
706                                         data_results wholeResults = chimera->getResults();
707                                         
708                                         //determine if we need to split
709                                         bool isChimeric = false;
710                                         
711                                         if (wholeResults.flag == "yes") {
712                                                 string chimeraFlag = "no";
713                                                 if(  (wholeResults.results[0].bsa >= minBS && wholeResults.results[0].divr_qla_qrb >= divR)
714                                                    ||
715                                                    (wholeResults.results[0].bsb >= minBS && wholeResults.results[0].divr_qlb_qra >= divR) ) { chimeraFlag = "yes"; }
716                                                 
717                                                 
718                                                 if (chimeraFlag == "yes") {     
719                                                         if ((wholeResults.results[0].bsa >= minBS) || (wholeResults.results[0].bsb >= minBS)) { isChimeric = true; }
720                                                 }
721                                         }
722                                         
723                                         if ((!isChimeric) && trimera) {
724                                                 
725                                                 //split sequence in half by bases
726                                                 string leftQuery, rightQuery;
727                                                 Sequence tempSeq(candidateSeq->getName(), candidateAligned);
728                                                 divideInHalf(tempSeq, leftQuery, rightQuery);
729                                                 
730                                                 //run chimeraSlayer on each piece
731                                                 Sequence* left = new Sequence(candidateSeq->getName(), leftQuery);
732                                                 Sequence* right = new Sequence(candidateSeq->getName(), rightQuery);
733                                                 
734                                                 //find chimeras
735                                                 chimera->getChimeras(left);
736                                                 data_results leftResults = chimera->getResults();
737                                                 
738                                                 chimera->getChimeras(right);
739                                                 data_results rightResults = chimera->getResults();
740                                                 
741                                                 //if either piece is chimeric then report
742                                                 Sequence trimmed = chimera->print(out, out2, leftResults, rightResults);
743                                                 if (trim) { trimmed.printSequence(out3);  }
744                                                 
745                                                 delete left; delete right;
746                                                 
747                                         }else { //already chimeric
748                                                 //print results
749                                                 Sequence trimmed = chimera->print(out, out2);
750                                                 if (trim) { trimmed.printSequence(out3);  }
751                                         }
752                                         
753                                         
754                                 }
755                                 count++;
756                         }
757                         
758                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
759                                 unsigned long int pos = inFASTA.tellg();
760                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
761                         #else
762                                 if (inFASTA.eof()) { break; }
763                         #endif
764                         
765                         delete candidateSeq;
766                         //report progress
767                         if((count) % 100 == 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
768                 }
769                 //report progress
770                 if((count) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
771                 
772                 out.close();
773                 out2.close();
774                 if (trim) { out3.close(); }
775                 inFASTA.close();
776                                 
777                 return count;
778         }
779         catch(exception& e) {
780                 m->errorOut(e, "ChimeraSlayerCommand", "driver");
781                 exit(1);
782         }
783 }
784 //**********************************************************************************************************************
785 #ifdef USE_MPI
786 int ChimeraSlayerCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, MPI_File& outAccMPI, MPI_File& outFastaMPI, vector<unsigned long int>& MPIPos){
787         try {                           
788                 MPI_Status status; 
789                 int pid;
790                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
791                 
792                 for(int i=0;i<num;i++){
793                         
794                         if (m->control_pressed) {       return 1;       }
795                         
796                         //read next sequence
797                         int length = MPIPos[start+i+1] - MPIPos[start+i];
798
799                         char* buf4 = new char[length];
800                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
801         
802                         string tempBuf = buf4;
803                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
804                         istringstream iss (tempBuf,istringstream::in);
805
806                         delete buf4;
807
808                         Sequence* candidateSeq = new Sequence(iss);  m->gobble(iss);
809                         string candidateAligned = candidateSeq->getAligned();
810                 
811                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
812                                 
813                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
814                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
815                                 }else{
816                 
817                                         //find chimeras
818                                         chimera->getChimeras(candidateSeq);
819                         
820                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
821                                         
822                                         //if you are not chimeric, then check each half
823                                         data_results wholeResults = chimera->getResults();
824                                         
825                                         //determine if we need to split
826                                         bool isChimeric = false;
827                                         
828                                         if (wholeResults.flag == "yes") {
829                                                 string chimeraFlag = "no";
830                                                 if(  (wholeResults.results[0].bsa >= minBS && wholeResults.results[0].divr_qla_qrb >= divR)
831                                                    ||
832                                                    (wholeResults.results[0].bsb >= minBS && wholeResults.results[0].divr_qlb_qra >= divR) ) { chimeraFlag = "yes"; }
833                                                 
834                                                 
835                                                 if (chimeraFlag == "yes") {     
836                                                         if ((wholeResults.results[0].bsa >= minBS) || (wholeResults.results[0].bsb >= minBS)) { isChimeric = true; }
837                                                 }
838                                         }
839                                         
840                                         if ((!isChimeric) && trimera) {                                                 
841                                                 //split sequence in half by bases
842                                                 string leftQuery, rightQuery;
843                                                 Sequence tempSeq(candidateSeq->getName(), candidateAligned);
844                                                 divideInHalf(tempSeq, leftQuery, rightQuery);
845                                                 
846                                                 //run chimeraSlayer on each piece
847                                                 Sequence* left = new Sequence(candidateSeq->getName(), leftQuery);
848                                                 Sequence* right = new Sequence(candidateSeq->getName(), rightQuery);
849                                                 
850                                                 //find chimeras
851                                                 chimera->getChimeras(left);
852                                                 data_results leftResults = chimera->getResults();
853                                                 
854                                                 chimera->getChimeras(right);
855                                                 data_results rightResults = chimera->getResults();
856                                                 
857                                                 //if either piece is chimeric then report
858                                                 Sequence trimmed = chimera->print(outMPI, outAccMPI, leftResults, rightResults);
859                                                 if (trim) {  
860                                                         string outputString = ">" + trimmed.getName() + "\n" + trimmed.getAligned() + "\n";
861                                                         
862                                                         //write to accnos file
863                                                         int length = outputString.length();
864                                                         char* buf2 = new char[length];
865                                                         memcpy(buf2, outputString.c_str(), length);
866                                                         
867                                                         MPI_File_write_shared(outFastaMPI, buf2, length, MPI_CHAR, &status);
868                                                         delete buf2;
869                                                 }
870                                                 
871                                                 delete left; delete right;
872                                                 
873                                         }else { 
874                                                 //print results
875                                                 Sequence trimmed = chimera->print(outMPI, outAccMPI);
876                                                 
877                                                 if (trim) {  
878                                                         string outputString = ">" + trimmed.getName() + "\n" + trimmed.getAligned() + "\n";
879                                                         
880                                                         //write to accnos file
881                                                         int length = outputString.length();
882                                                         char* buf2 = new char[length];
883                                                         memcpy(buf2, outputString.c_str(), length);
884                                                         
885                                                         MPI_File_write_shared(outFastaMPI, buf2, length, MPI_CHAR, &status);
886                                                         delete buf2;
887                                                 }
888                                         }
889                                         
890                                 }
891                         }
892                         delete candidateSeq;
893                         
894                         //report progress
895                         if((i+1) % 100 == 0){  cout << "Processing sequence: " << (i+1) << endl;        m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n");          }
896                 }
897                 //report progress
898                 if(num % 100 != 0){             cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n");  }
899                 
900                                 
901                 return 0;
902         }
903         catch(exception& e) {
904                 m->errorOut(e, "ChimeraSlayerCommand", "driverMPI");
905                 exit(1);
906         }
907 }
908 #endif
909
910 /**************************************************************************************************/
911
912 int ChimeraSlayerCommand::createProcesses(string outputFileName, string filename, string accnos, string fasta) {
913         try {
914 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
915                 int process = 0;
916                 int num = 0;
917                 
918                 //loop through and create all the processes you want
919                 while (process != processors) {
920                         int pid = fork();
921                         
922                         if (pid > 0) {
923                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
924                                 process++;
925                         }else if (pid == 0){
926                                 num = driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename, accnos + toString(getpid()) + ".temp", fasta + toString(getpid()) + ".temp");
927                                 
928                                 //pass numSeqs to parent
929                                 ofstream out;
930                                 string tempFile = outputFileName + toString(getpid()) + ".num.temp";
931                                 m->openOutputFile(tempFile, out);
932                                 out << num << '\t' << chimera->getNumNoParents() << endl;
933                                 out.close();
934                                 exit(0);
935                         }else { 
936                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
937                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
938                                 exit(0);
939                         }
940                 }
941                 
942                 //force parent to wait until all the processes are done
943                 for (int i=0;i<processors;i++) { 
944                         int temp = processIDS[i];
945                         wait(&temp);
946                 }
947                 
948                 int numNoParents = 0;
949                 for (int i = 0; i < processIDS.size(); i++) {
950                         ifstream in;
951                         string tempFile =  outputFileName + toString(processIDS[i]) + ".num.temp";
952                         m->openInputFile(tempFile, in);
953                         if (!in.eof()) { int tempNum = 0; int tempNumParents = 0; in >> tempNum >> tempNumParents; num += tempNum; numNoParents += tempNumParents; }
954                         in.close(); remove(tempFile.c_str());
955                 }
956                 
957                 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(); }
958                 
959                 return num;
960 #endif          
961         }
962         catch(exception& e) {
963                 m->errorOut(e, "ChimeraSlayerCommand", "createProcesses");
964                 exit(1);
965         }
966 }
967
968 /**************************************************************************************************/
969
970 int ChimeraSlayerCommand::divideInHalf(Sequence querySeq, string& leftQuery, string& rightQuery) {
971         try {
972                 
973                 string queryUnAligned = querySeq.getUnaligned();
974                 int numBases = int(queryUnAligned.length() * 0.5);
975                 
976                 string queryAligned = querySeq.getAligned();
977                 leftQuery = querySeq.getAligned();
978                 rightQuery = querySeq.getAligned();
979                 
980                 int baseCount = 0;
981                 int leftSpot = 0;
982                 for (int i = 0; i < queryAligned.length(); i++) {
983                         //if you are a base
984                         if (isalpha(queryAligned[i])) {         
985                                 baseCount++; 
986                         }
987                         
988                         //if you have half
989                         if (baseCount >= numBases) {  leftSpot = i; break; } //first half
990                 }
991                 
992                 //blank out right side
993                 for (int i = leftSpot; i < leftQuery.length(); i++) { leftQuery[i] = '.'; }
994                 
995                 //blank out left side
996                 for (int i = 0; i < leftSpot; i++) { rightQuery[i] = '.'; }
997                 
998                 return 0;
999                 
1000         }
1001         catch(exception& e) {
1002                 m->errorOut(e, "ChimeraSlayerCommand", "divideInHalf");
1003                 exit(1);
1004         }
1005 }
1006 /**************************************************************************************************/
1007 map<string, int> ChimeraSlayerCommand::sortFastaFile(string fastaFile, string nameFile) {
1008         try {
1009                 map<string, int> nameAbund;
1010                 
1011                 //read through fastafile and store info
1012                 map<string, string> seqs;
1013                 ifstream in;
1014                 m->openInputFile(fastaFile, in);
1015                 
1016                 while (!in.eof()) {
1017                         
1018                         if (m->control_pressed) { in.close(); return nameAbund; }
1019                         
1020                         Sequence seq(in); m->gobble(in);
1021                         seqs[seq.getName()] = seq.getAligned();
1022                 }
1023                 
1024                 in.close();
1025                 
1026                 //read namefile
1027                 vector<seqPriorityNode> nameMapCount;
1028                 int error = m->readNames(nameFile, nameMapCount, seqs);
1029                 
1030                 if (m->control_pressed) { return nameAbund; }
1031                 
1032                 if (error == 1) { m->control_pressed = true; return nameAbund; }
1033                 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; }
1034
1035                 sort(nameMapCount.begin(), nameMapCount.end(), compareSeqPriorityNodes);
1036                 
1037                 string newFasta = fastaFile + ".temp";
1038                 ofstream out;
1039                 m->openOutputFile(newFasta, out);
1040                 
1041                 //print new file in order of
1042                 for (int i = 0; i < nameMapCount.size(); i++) {
1043                         out << ">" << nameMapCount[i].name << endl << nameMapCount[i].seq << endl;
1044                         nameAbund[nameMapCount[i].name] = nameMapCount[i].numIdentical;
1045                 }
1046                 out.close();
1047                 
1048                 rename(newFasta.c_str(), fastaFile.c_str());
1049                                 
1050                 return nameAbund;
1051                 
1052         }
1053         catch(exception& e) {
1054                 m->errorOut(e, "ChimeraSlayerCommand", "sortFastaFile");
1055                 exit(1);
1056         }
1057 }
1058
1059 /**************************************************************************************************/
1060