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