]> git.donarmstrong.com Git - mothur.git/blob - chimeraslayercommand.cpp
modified calculators to use doubles, added numotu and fontsize parameters to heatmap...
[mothur.git] / chimeraslayercommand.cpp
1 /*
2  *  chimeraslayercommand.cpp
3  *  Mothur
4  *
5  *  Created by westcott on 3/31/10.
6  *  Copyright 2010 Schloss Lab. All rights reserved.
7  *
8  */
9
10 #include "chimeraslayercommand.h"
11 #include "chimeraslayer.h"
12
13
14 //***************************************************************************************************************
15
16 ChimeraSlayerCommand::ChimeraSlayerCommand(string option)  {
17         try {
18                 abort = false;
19                 
20                 //allow user to run help
21                 if(option == "help") { help(); abort = true; }
22                 
23                 else {
24                         //valid paramters for this command
25                         string Array[] =  {"fasta", "processors", "window", "template","numwanted", "ksize", "match","mismatch", 
26                         "divergence", "minsim","mincov","minbs", "minsnp","parents", "iters","outputdir","inputdir", "search","realign" };
27                         vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
28                         
29                         OptionParser parser(option);
30                         map<string,string> parameters = parser.getParameters();
31                         
32                         ValidParameters validParameter("chimera.slayer");
33                         map<string,string>::iterator it;
34                         
35                         //check to make sure all parameters are valid for command
36                         for (it = parameters.begin(); it != parameters.end(); it++) { 
37                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
38                         }
39                         
40                         //if the user changes the input directory command factory will send this info to us in the output parameter 
41                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
42                         if (inputDir == "not found"){   inputDir = "";          }
43                         else {
44                                 string path;
45                                 it = parameters.find("template");
46                                 //user has given a template file
47                                 if(it != parameters.end()){ 
48                                         path = hasPath(it->second);
49                                         //if the user has not given a path then, add inputdir. else leave path alone.
50                                         if (path == "") {       parameters["template"] = inputDir + it->second;         }
51                                 }
52                         }
53
54                         
55                         //check for required parameters
56                         fastafile = validParameter.validFile(parameters, "fasta", false);
57                         if (fastafile == "not found") { fastafile = ""; m->mothurOut("fasta is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true;  }
58                         else { 
59                                 splitAtDash(fastafile, fastaFileNames);
60                                 
61                                 //go through files and make sure they are good, if not, then disregard them
62                                 for (int i = 0; i < fastaFileNames.size(); i++) {
63                                         if (inputDir != "") {
64                                                 string path = hasPath(fastaFileNames[i]);
65                                                 //if the user has not given a path then, add inputdir. else leave path alone.
66                                                 if (path == "") {       fastaFileNames[i] = inputDir + fastaFileNames[i];               }
67                                         }
68         
69                                         int ableToOpen;
70                                         ifstream in;
71                                         
72                                         ableToOpen = openInputFile(fastaFileNames[i], in, "noerror");
73                                 
74                                         //if you can't open it, try default location
75                                         if (ableToOpen == 1) {
76                                                 if (m->getDefaultPath() != "") { //default path is set
77                                                         string tryPath = m->getDefaultPath() + getSimpleName(fastaFileNames[i]);
78                                                         m->mothurOut("Unable to open " + fastaFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
79                                                         ableToOpen = openInputFile(tryPath, in, "noerror");
80                                                         fastaFileNames[i] = tryPath;
81                                                 }
82                                         }
83                                         in.close();
84                                         
85                                         if (ableToOpen == 1) { 
86                                                 m->mothurOut("Unable to open " + fastaFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
87                                                 //erase from file list
88                                                 fastaFileNames.erase(fastaFileNames.begin()+i);
89                                                 i--;
90                                         }
91                                 }
92                                 
93                                 //make sure there is at least one valid file left
94                                 if (fastaFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
95                         }
96                         
97                         //if the user changes the output directory command factory will send this info to us in the output parameter 
98                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
99                                 outputDir = ""; 
100                                 outputDir += hasPath(fastafile); //if user entered a file with a path then preserve it  
101                         }
102
103                         templatefile = validParameter.validFile(parameters, "template", true);
104                         if (templatefile == "not open") { abort = true; }
105                         else if (templatefile == "not found") { templatefile = "";  m->mothurOut("template is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true;  }   
106                                                 
107                         string temp = validParameter.validFile(parameters, "processors", false);                if (temp == "not found") { temp = "1"; }
108                         convert(temp, processors);
109                         
110                         temp = validParameter.validFile(parameters, "ksize", false);                    if (temp == "not found") { temp = "7"; }
111                         convert(temp, ksize);
112                                                 
113                         temp = validParameter.validFile(parameters, "window", false);                   if (temp == "not found") { temp = "50"; }                       
114                         convert(temp, window);
115                         
116                         temp = validParameter.validFile(parameters, "match", false);                    if (temp == "not found") { temp = "5"; }
117                         convert(temp, match);
118                         
119                         temp = validParameter.validFile(parameters, "mismatch", false);                 if (temp == "not found") { temp = "-4"; }
120                         convert(temp, mismatch);
121                         
122                         temp = validParameter.validFile(parameters, "divergence", false);               if (temp == "not found") { temp = "1.007"; }
123                         convert(temp, divR);
124                         
125                         temp = validParameter.validFile(parameters, "minsim", false);                   if (temp == "not found") { temp = "90"; }
126                         convert(temp, minSimilarity);
127                         
128                         temp = validParameter.validFile(parameters, "mincov", false);                   if (temp == "not found") { temp = "70"; }
129                         convert(temp, minCoverage);
130                         
131                         temp = validParameter.validFile(parameters, "minbs", false);                    if (temp == "not found") { temp = "90"; }
132                         convert(temp, minBS);
133                         
134                         temp = validParameter.validFile(parameters, "minsnp", false);                   if (temp == "not found") { temp = "100"; }
135                         convert(temp, minSNP);
136
137                         temp = validParameter.validFile(parameters, "parents", false);                  if (temp == "not found") { temp = "3"; }
138                         convert(temp, parents); 
139                         
140                         temp = validParameter.validFile(parameters, "realign", false);                  if (temp == "not found") { temp = "f"; }
141                         realign = isTrue(temp); 
142                         
143                         search = validParameter.validFile(parameters, "search", false);                 if (search == "not found") { search = "distance"; }
144                         
145                         temp = validParameter.validFile(parameters, "iters", false);                    if (temp == "not found") { temp = "100"; }              
146                         convert(temp, iters); 
147                          
148                         temp = validParameter.validFile(parameters, "increment", false);                if (temp == "not found") { temp = "5"; }
149                         convert(temp, increment);
150                         
151                         temp = validParameter.validFile(parameters, "numwanted", false);                if (temp == "not found") { temp = "15"; }               
152                         convert(temp, numwanted);
153
154                         if ((search != "distance") && (search != "blast") && (search != "kmer")) { m->mothurOut(search + " is not a valid search."); m->mothurOutEndLine(); abort = true;  }
155                 }
156         }
157         catch(exception& e) {
158                 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
159                 exit(1);
160         }
161 }
162 //**********************************************************************************************************************
163
164 void ChimeraSlayerCommand::help(){
165         try {
166         
167                 m->mothurOut("The chimera.slayer command reads a fastafile and templatefile and outputs potentially chimeric sequences.\n");
168                 m->mothurOut("This command was modeled after the chimeraSlayer written by the Broad Institute.\n");
169                 m->mothurOut("The chimera.slayer command parameters are fasta, template, processors, ksize, window, match, mismatch, divergence. minsim, mincov, minbs, minsnp, parents, search, iters, increment and numwanted.\n"); //realign,
170                 m->mothurOut("The fasta parameter allows you to enter the fasta file containing your potentially chimeric sequences, and is required. \n");
171                 m->mothurOut("You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amazon.fasta \n");
172                 m->mothurOut("The template parameter allows you to enter a template file containing known non-chimeric sequences, and is required. \n");
173                 m->mothurOut("The processors parameter allows you to specify how many processors you would like to use.  The default is 1. \n");
174                 #ifdef USE_MPI
175                 m->mothurOut("When using MPI, the processors parameter is set to the number of MPI processes running. \n");
176                 #endif
177                 m->mothurOut("The window parameter allows you to specify the window size for searching for chimeras, default=50. \n");
178                 m->mothurOut("The increment parameter allows you to specify how far you move each window while finding chimeric sequences, default=5.\n");
179                 m->mothurOut("The numwanted parameter allows you to specify how many sequences you would each query sequence compared with, default=15.\n");
180                 m->mothurOut("The ksize parameter allows you to input kmersize, default is 7, used if search is kmer. \n");
181                 m->mothurOut("The match parameter allows you to reward matched bases in blast search, default is 5. \n");
182                 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");
183                 m->mothurOut("The mismatch parameter allows you to penalize mismatched bases in blast search, default is -4. \n");
184                 m->mothurOut("The divergence parameter allows you to set a cutoff for chimera determination, default is 1.007. \n");
185                 m->mothurOut("The iters parameter allows you to specify the number of bootstrap iters to do with the chimeraslayer method, default=100.\n");
186                 m->mothurOut("The minsim parameter allows you to specify a minimum similarity with the parent fragments, default=90. \n");
187                 m->mothurOut("The mincov parameter allows you to specify minimum coverage by closest matches found in template. Default is 70, meaning 70%. \n");
188                 m->mothurOut("The minbs parameter allows you to specify minimum bootstrap support for calling a sequence chimeric. Default is 90, meaning 90%. \n");
189                 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");
190                 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");
191                 m->mothurOut("The realign parameter allows you to realign the query to the potential parents. Choices are true or false, default false.  \n");
192                 m->mothurOut("The chimera.slayer command should be in the following format: \n");
193                 m->mothurOut("chimera.slayer(fasta=yourFastaFile, template=yourTemplate, search=yourSearch) \n");
194                 m->mothurOut("Example: chimera.slayer(fasta=AD.align, template=core_set_aligned.imputed.fasta, search=kmer) \n");
195                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");     
196         }
197         catch(exception& e) {
198                 m->errorOut(e, "ChimeraSlayerCommand", "help");
199                 exit(1);
200         }
201 }
202
203 //***************************************************************************************************************
204
205 ChimeraSlayerCommand::~ChimeraSlayerCommand(){  /*      do nothing      */      }
206
207 //***************************************************************************************************************
208
209 int ChimeraSlayerCommand::execute(){
210         try{
211                 
212                 if (abort == true) { return 0; }
213                 
214                 for (int s = 0; s < fastaFileNames.size(); s++) {
215                                 
216                         m->mothurOut("Checking sequences from " + fastaFileNames[s] + " ..." ); m->mothurOutEndLine();
217                 
218                         int start = time(NULL); 
219                         
220                         chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);   
221                                                         
222                         string outputFileName = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + "slayer.chimeras";
223                         string accnosFileName = outputDir + getRootName(getSimpleName(fastaFileNames[s]))  + "slayer.accnos";
224                         
225                         if (m->control_pressed) { delete chimera; for (int j = 0; j < outputNames.size(); j++) {        remove(outputNames[j].c_str()); }  return 0;    }
226                         
227                         if (chimera->getUnaligned()) { 
228                                 m->mothurOut("Your template sequences are different lengths, please correct."); m->mothurOutEndLine(); 
229                                 delete chimera;
230                                 return 0; 
231                         }
232                         templateSeqsLength = chimera->getLength();
233                         
234                 #ifdef USE_MPI  
235                         int pid, end, numSeqsPerProcessor; 
236                                 int tag = 2001;
237                                 vector<unsigned long int> MPIPos;
238                                 
239                                 MPI_Status status; 
240                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
241                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); 
242
243                                 MPI_File inMPI;
244                                 MPI_File outMPI;
245                                 MPI_File outMPIAccnos;
246                                 
247                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
248                                 int inMode=MPI_MODE_RDONLY; 
249                                 
250                                 char outFilename[1024];
251                                 strcpy(outFilename, outputFileName.c_str());
252                                 
253                                 char outAccnosFilename[1024];
254                                 strcpy(outAccnosFilename, accnosFileName.c_str());
255                                 
256                                 char inFileName[1024];
257                                 strcpy(inFileName, fastaFileNames[s].c_str());
258
259                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
260                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
261                                 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
262
263                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   MPI_File_close(&outMPIAccnos); for (int j = 0; j < outputNames.size(); j++) {    remove(outputNames[j].c_str()); }   delete chimera; return 0;  }
264                         
265                                 if (pid == 0) { //you are the root process 
266                                         m->mothurOutEndLine();
267                                         m->mothurOut("Only reporting sequence supported by " + toString(minBS) + "% of bootstrapped results.");
268                                         m->mothurOutEndLine();
269                 
270                                         string outTemp = "Name\tLeftParent\tRightParent\tDivQLAQRB\tPerIDQLAQRB\tBootStrapA\tDivQLBQRA\tPerIDQLBQRA\tBootStrapB\tFlag\tLeftWindow\tRightWindow\n";
271                                         
272                                         //print header
273                                         int length = outTemp.length();
274                                         char* buf2 = new char[length];
275                                         memcpy(buf2, outTemp.c_str(), length);
276
277                                         MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &status);
278                                         delete buf2;
279
280                                         MPIPos = setFilePosFasta(fastaFileNames[s], numSeqs); //fills MPIPos, returns numSeqs
281                                         
282                                         //send file positions to all processes
283                                         for(int i = 1; i < processors; i++) { 
284                                                 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
285                                                 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
286                                         }
287                                         
288                                         //figure out how many sequences you have to align
289                                         numSeqsPerProcessor = numSeqs / processors;
290                                         int startIndex =  pid * numSeqsPerProcessor;
291                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
292                                 
293                                         //do your part
294                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
295                                         
296                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   MPI_File_close(&outMPIAccnos);  for (int j = 0; j < outputNames.size(); j++) {   remove(outputNames[j].c_str()); }  remove(outputFileName.c_str());  remove(accnosFileName.c_str());  delete chimera; return 0;  }
297
298                                 }else{ //you are a child process
299                                         MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
300                                         MPIPos.resize(numSeqs+1);
301                                         MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
302                                         
303                                         //figure out how many sequences you have to align
304                                         numSeqsPerProcessor = numSeqs / processors;
305                                         int startIndex =  pid * numSeqsPerProcessor;
306                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
307                                         
308                                         //do your part
309                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
310                                         
311                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   MPI_File_close(&outMPIAccnos);  for (int j = 0; j < outputNames.size(); j++) {   remove(outputNames[j].c_str()); }  delete chimera; return 0;  }
312                                 }
313                                 
314                                 //close files 
315                                 MPI_File_close(&inMPI);
316                                 MPI_File_close(&outMPI);
317                                 MPI_File_close(&outMPIAccnos);
318                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
319                                 
320                 #else
321                         ofstream outHeader;
322                         string tempHeader = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + "slayer.chimeras.tempHeader";
323                         openOutputFile(tempHeader, outHeader);
324                         
325                         chimera->printHeader(outHeader);
326                         outHeader.close();
327                         
328                         //break up file
329                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
330                                 if(processors == 1){
331                                         ifstream inFASTA;
332                                         openInputFile(fastaFileNames[s], inFASTA);
333                                         getNumSeqs(inFASTA, numSeqs);
334                                         inFASTA.close();
335                                         
336                                         lines.push_back(new linePair(0, numSeqs));
337                                         
338                                         driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName);
339                                         
340                                         if (m->control_pressed) { 
341                                                 remove(outputFileName.c_str()); 
342                                                 remove(tempHeader.c_str()); 
343                                                 remove(accnosFileName.c_str());
344                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
345                                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
346                                                 delete chimera;
347                                                 return 0;
348                                         }
349                                         
350                                 }else{
351                                         vector<unsigned long int> positions;
352                                         processIDS.resize(0);
353                                         
354                                         ifstream inFASTA;
355                                         openInputFile(fastaFileNames[s], inFASTA);
356                                         
357                                         string input;
358                                         while(!inFASTA.eof()){
359                                                 input = getline(inFASTA);
360                                                 if (input.length() != 0) {
361                                                         if(input[0] == '>'){    unsigned long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
362                                                 }
363                                         }
364                                         inFASTA.close();
365                                         
366                                         numSeqs = positions.size();
367                                         
368                                         int numSeqsPerProcessor = numSeqs / processors;
369                                         
370                                         for (int i = 0; i < processors; i++) {
371                                                 unsigned long int startPos = positions[ i * numSeqsPerProcessor ];
372                                                 if(i == processors - 1){
373                                                         numSeqsPerProcessor = numSeqs - i * numSeqsPerProcessor;
374                                                 }
375                                                 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
376                                         }
377                                         
378                                         createProcesses(outputFileName, fastaFileNames[s], accnosFileName); 
379                                 
380                                         rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
381                                         rename((accnosFileName + toString(processIDS[0]) + ".temp").c_str(), accnosFileName.c_str());
382                                                 
383                                         //append output files
384                                         for(int i=1;i<processors;i++){
385                                                 appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
386                                                 remove((outputFileName + toString(processIDS[i]) + ".temp").c_str());
387                                         }
388                                         
389                                         //append output files
390                                         for(int i=1;i<processors;i++){
391                                                 appendFiles((accnosFileName + toString(processIDS[i]) + ".temp"), accnosFileName);
392                                                 remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());
393                                         }
394                                         
395                                         if (m->control_pressed) { 
396                                                 remove(outputFileName.c_str()); 
397                                                 remove(accnosFileName.c_str());
398                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
399                                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
400                                                 delete chimera;
401                                                 return 0;
402                                         }
403
404                                 }
405
406                         #else
407                                 ifstream inFASTA;
408                                 openInputFile(fastaFileNames[s], inFASTA);
409                                 getNumSeqs(inFASTA, numSeqs);
410                                 inFASTA.close();
411                                 lines.push_back(new linePair(0, numSeqs));
412                                 
413                                 driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName);
414                                 
415                                 if (m->control_pressed) { 
416                                                 remove(outputFileName.c_str()); 
417                                                 remove(tempHeader.c_str()); 
418                                                 remove(accnosFileName.c_str());
419                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
420                                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
421                                                 delete chimera;
422                                                 return 0;
423                                 }
424                                 
425                         #endif
426                         
427                         appendFiles(outputFileName, tempHeader);
428                 
429                         remove(outputFileName.c_str());
430                         rename(tempHeader.c_str(), outputFileName.c_str());
431                         
432                 #endif
433                         delete chimera;
434                         
435                         
436                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
437                         
438                         outputNames.push_back(outputFileName);
439                         outputNames.push_back(accnosFileName); 
440                         
441                         m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine();
442                 }
443                 
444                 m->mothurOutEndLine();
445                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
446                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }       
447                 m->mothurOutEndLine();
448
449                 return 0;
450                 
451         }
452         catch(exception& e) {
453                 m->errorOut(e, "ChimeraSlayerCommand", "execute");
454                 exit(1);
455         }
456 }
457 //**********************************************************************************************************************
458
459 int ChimeraSlayerCommand::driver(linePair* line, string outputFName, string filename, string accnos){
460         try {
461                 ofstream out;
462                 openOutputFile(outputFName, out);
463                 
464                 ofstream out2;
465                 openOutputFile(accnos, out2);
466                 
467                 ifstream inFASTA;
468                 openInputFile(filename, inFASTA);
469
470                 inFASTA.seekg(line->start);
471                 
472                 for(int i=0;i<line->numSeqs;i++){
473                 
474                         if (m->control_pressed) {       return 1;       }
475                 
476                         Sequence* candidateSeq = new Sequence(inFASTA);  gobble(inFASTA);
477                                 
478                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
479                                 
480                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
481                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
482                                 }else{
483                                         //find chimeras
484                                         chimera->getChimeras(candidateSeq);
485                                         
486                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
487                 
488                                         //print results
489                                         chimera->print(out, out2);
490                                 }
491                         }
492                         delete candidateSeq;
493                         
494                         //report progress
495                         if((i+1) % 100 == 0){   m->mothurOut("Processing sequence: " + toString(i+1)); m->mothurOutEndLine();           }
496                 }
497                 //report progress
498                 if((line->numSeqs) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(line->numSeqs)); m->mothurOutEndLine();         }
499                 
500                 out.close();
501                 out2.close();
502                 inFASTA.close();
503                                 
504                 return 0;
505         }
506         catch(exception& e) {
507                 m->errorOut(e, "ChimeraSlayerCommand", "driver");
508                 exit(1);
509         }
510 }
511 //**********************************************************************************************************************
512 #ifdef USE_MPI
513 int ChimeraSlayerCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, MPI_File& outAccMPI, vector<unsigned long int>& MPIPos){
514         try {                           
515                 MPI_Status status; 
516                 int pid;
517                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
518                 
519                 for(int i=0;i<num;i++){
520                         
521                         if (m->control_pressed) {       return 1;       }
522                         
523                         //read next sequence
524                         int length = MPIPos[start+i+1] - MPIPos[start+i];
525
526                         char* buf4 = new char[length];
527                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
528         
529                         string tempBuf = buf4;
530                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
531                         istringstream iss (tempBuf,istringstream::in);
532
533                         delete buf4;
534
535                         Sequence* candidateSeq = new Sequence(iss);  gobble(iss);
536                 
537                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
538                                 
539                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
540                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
541                                 }else{
542                 
543                                         //find chimeras
544                                         chimera->getChimeras(candidateSeq);
545                         
546                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
547                 //cout << "about to print" << endl;
548                                         //print results
549                                         bool isChimeric = chimera->print(outMPI, outAccMPI);
550                                 }
551                         }
552                         delete candidateSeq;
553                         
554                         //report progress
555                         if((i+1) % 100 == 0){  cout << "Processing sequence: " << (i+1) << endl;        m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n");          }
556                 }
557                 //report progress
558                 if(num % 100 != 0){             cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n");  }
559                 
560                                 
561                 return 0;
562         }
563         catch(exception& e) {
564                 m->errorOut(e, "ChimeraSlayerCommand", "driverMPI");
565                 exit(1);
566         }
567 }
568 #endif
569
570 /**************************************************************************************************/
571
572 int ChimeraSlayerCommand::createProcesses(string outputFileName, string filename, string accnos) {
573         try {
574 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
575                 int process = 0;
576                 //              processIDS.resize(0);
577                 
578                 //loop through and create all the processes you want
579                 while (process != processors) {
580                         int pid = fork();
581                         
582                         if (pid > 0) {
583                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
584                                 process++;
585                         }else if (pid == 0){
586                                 driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename, accnos + toString(getpid()) + ".temp");
587                                 exit(0);
588                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
589                 }
590                 
591                 //force parent to wait until all the processes are done
592                 for (int i=0;i<processors;i++) { 
593                         int temp = processIDS[i];
594                         wait(&temp);
595                 }
596                 
597                 return 0;
598 #endif          
599         }
600         catch(exception& e) {
601                 m->errorOut(e, "ChimeraSlayerCommand", "createProcesses");
602                 exit(1);
603         }
604 }
605
606 /**************************************************************************************************/
607
608