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