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