]> git.donarmstrong.com Git - mothur.git/blob - filterseqscommand.cpp
working on chimera change to add trim feature, fixed bug in print of distance file...
[mothur.git] / filterseqscommand.cpp
1 /*
2  *  filterseqscommand.cpp
3  *  Mothur
4  *
5  *  Created by Thomas Ryabin on 5/4/09.
6  *  Copyright 2009 Schloss Lab UMASS Amherst. All rights reserved.
7  *
8  */
9
10 #include "filterseqscommand.h"
11 #include "sequence.hpp"
12
13 //**********************************************************************************************************************
14 vector<string> FilterSeqsCommand::getValidParameters(){ 
15         try {
16                 string Array[] =  {"fasta", "trump", "soft", "hard", "vertical", "outputdir","inputdir", "processors"};
17                 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
18                 return myArray;
19         }
20         catch(exception& e) {
21                 m->errorOut(e, "FilterSeqsCommand", "getValidParameters");
22                 exit(1);
23         }
24 }
25 //**********************************************************************************************************************
26 FilterSeqsCommand::FilterSeqsCommand(){ 
27         try {
28                 abort = true;
29                 //initialize outputTypes
30                 vector<string> tempOutNames;
31                 outputTypes["fasta"] = tempOutNames;
32                 outputTypes["filter"] = tempOutNames;
33         }
34         catch(exception& e) {
35                 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
36                 exit(1);
37         }
38 }
39 //**********************************************************************************************************************
40 vector<string> FilterSeqsCommand::getRequiredParameters(){      
41         try {
42                 string Array[] =  {"fasta"};
43                 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
44                 return myArray;
45         }
46         catch(exception& e) {
47                 m->errorOut(e, "FilterSeqsCommand", "getRequiredParameters");
48                 exit(1);
49         }
50 }
51 //**********************************************************************************************************************
52 vector<string> FilterSeqsCommand::getRequiredFiles(){   
53         try {
54                 vector<string> myArray;
55                 return myArray;
56         }
57         catch(exception& e) {
58                 m->errorOut(e, "FilterSeqsCommand", "getRequiredFiles");
59                 exit(1);
60         }
61 }
62 /**************************************************************************************/
63 FilterSeqsCommand::FilterSeqsCommand(string option)  {
64         try {
65                 abort = false;
66                 filterFileName = "";
67                 
68                 //allow user to run help
69                 if(option == "help") { help(); abort = true; }
70                 
71                 else {
72                         //valid paramters for this command
73                         string Array[] =  {"fasta", "trump", "soft", "hard", "vertical", "outputdir","inputdir", "processors"};
74                         vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
75                         
76                         OptionParser parser(option);
77                         map<string,string> parameters = parser.getParameters();
78                         
79                         ValidParameters validParameter("filter.seqs");
80                         map<string,string>::iterator it;
81                         
82                         //check to make sure all parameters are valid for command
83                         for (it = parameters.begin(); it != parameters.end(); it++) { 
84                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
85                         }
86                         
87                         //initialize outputTypes
88                         vector<string> tempOutNames;
89                         outputTypes["fasta"] = tempOutNames;
90                         outputTypes["filter"] = tempOutNames;
91                 
92                         //if the user changes the input directory command factory will send this info to us in the output parameter 
93                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
94                         if (inputDir == "not found"){   inputDir = "";          }
95                         else {
96                                 string path;
97                                 it = parameters.find("fasta");
98                                 //user has given a template file
99                                 if(it != parameters.end()){ 
100                                         path = m->hasPath(it->second);
101                                         //if the user has not given a path then, add inputdir. else leave path alone.
102                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
103                                 }
104                                 
105                                 it = parameters.find("hard");
106                                 //user has given a template file
107                                 if(it != parameters.end()){ 
108                                         path = m->hasPath(it->second);
109                                         //if the user has not given a path then, add inputdir. else leave path alone.
110                                         if (path == "") {       parameters["hard"] = inputDir + it->second;             }
111                                 }
112                         }
113                         
114                         //check for required parameters
115                         fasta = validParameter.validFile(parameters, "fasta", false);
116                         if (fasta == "not found") { m->mothurOut("fasta is a required parameter for the filter.seqs command."); m->mothurOutEndLine(); abort = true;  }
117                         else { 
118                                 m->splitAtDash(fasta, fastafileNames);
119                                 
120                                 //go through files and make sure they are good, if not, then disregard them
121                                 for (int i = 0; i < fastafileNames.size(); i++) {
122                                         if (inputDir != "") {
123                                                 string path = m->hasPath(fastafileNames[i]);
124                                                 //if the user has not given a path then, add inputdir. else leave path alone.
125                                                 if (path == "") {       fastafileNames[i] = inputDir + fastafileNames[i];               }
126                                         }
127
128                                         ifstream in;
129                                         int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
130                                 
131                                         //if you can't open it, try default location
132                                         if (ableToOpen == 1) {
133                                                 if (m->getDefaultPath() != "") { //default path is set
134                                                         string tryPath = m->getDefaultPath() + m->getSimpleName(fastafileNames[i]);
135                                                         m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
136                                                         ifstream in2;
137                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
138                                                         in2.close();
139                                                         fastafileNames[i] = tryPath;
140                                                 }
141                                         }
142                                         
143                                         //if you can't open it, try default location
144                                         if (ableToOpen == 1) {
145                                                 if (m->getOutputDir() != "") { //default path is set
146                                                         string tryPath = m->getOutputDir() + m->getSimpleName(fastafileNames[i]);
147                                                         m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
148                                                         ifstream in2;
149                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
150                                                         in2.close();
151                                                         fastafileNames[i] = tryPath;
152                                                 }
153                                         }
154                                         
155                                         in.close();
156                                         
157                                         if (ableToOpen == 1) { 
158                                                 m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
159                                                 //erase from file list
160                                                 fastafileNames.erase(fastafileNames.begin()+i);
161                                                 i--;
162                                         }else{  
163                                                 string simpleName = m->getSimpleName(fastafileNames[i]);
164                                                 filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
165                                         }
166                                         in.close();
167                                 }
168                                 
169                                 //make sure there is at least one valid file left
170                                 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
171                         }
172                         
173                         if (!abort) {
174                                 //if the user changes the output directory command factory will send this info to us in the output parameter 
175                                 outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
176                                         outputDir = ""; 
177                                         outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it       
178                                 }
179                         }
180                         //check for optional parameter and set defaults
181                         // ...at some point should added some additional type checking...
182                         
183                         string temp;
184                         hard = validParameter.validFile(parameters, "hard", true);                              if (hard == "not found") { hard = ""; }
185                         else if (hard == "not open") { hard = ""; abort = true; }       
186
187                         temp = validParameter.validFile(parameters, "trump", false);                    if (temp == "not found") { temp = "*"; }
188                         trump = temp[0];
189                         
190                         temp = validParameter.validFile(parameters, "soft", false);                             if (temp == "not found") { soft = 0; }
191                         else {  soft = (float)atoi(temp.c_str()) / 100.0;  }
192                         
193                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = "1";                             }
194                         convert(temp, processors); 
195                         
196                         vertical = validParameter.validFile(parameters, "vertical", false);             
197                         if (vertical == "not found") { 
198                                 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
199                                 else { vertical = "F";  }
200                         }
201                         
202                         numSeqs = 0;
203                 }
204                 
205         }
206         catch(exception& e) {
207                 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
208                 exit(1);
209         }
210 }
211
212 //**********************************************************************************************************************
213
214 void FilterSeqsCommand::help(){
215         try {
216                                 
217                 m->mothurOut("The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n");
218                 m->mothurOut("The filter.seqs command parameters are fasta, trump, soft, hard, processors and vertical. \n");
219                 m->mothurOut("The fasta parameter is required. You may enter several fasta files to build the filter from and filter, by separating their names with -'s.\n");
220                 m->mothurOut("For example: fasta=abrecovery.fasta-amazon.fasta \n");
221                 m->mothurOut("The trump parameter .... The default is ...\n");
222                 m->mothurOut("The soft parameter .... The default is ....\n");
223                 m->mothurOut("The hard parameter allows you to enter a file containing the filter you want to use.\n");
224                 m->mothurOut("The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n");
225                 m->mothurOut("The processors parameter allows you to specify the number of processors to use. The default is 1.\n");
226                 m->mothurOut("The filter.seqs command should be in the following format: \n");
227                 m->mothurOut("filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n");
228                 m->mothurOut("Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n");
229                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n\n");
230                 
231         }
232         catch(exception& e) {
233                 m->errorOut(e, "FilterSeqsCommand", "help");
234                 exit(1);
235         }
236 }
237
238 /**************************************************************************************/
239
240 int FilterSeqsCommand::execute() {      
241         try {
242         
243                 if (abort == true) { return 0; }
244                 
245                 ifstream inFASTA;
246                 m->openInputFile(fastafileNames[0], inFASTA);
247                 
248                 Sequence testSeq(inFASTA);
249                 alignmentLength = testSeq.getAlignLength();
250                 inFASTA.close();
251                 
252                 ////////////create filter/////////////////
253                 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
254                 
255                 filter = createFilter();
256                 
257                 m->mothurOutEndLine();  m->mothurOutEndLine();
258                 
259                 if (m->control_pressed) { outputTypes.clear(); return 0; }
260                 
261                 #ifdef USE_MPI
262                         int pid;
263                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
264                                         
265                         if (pid == 0) { //only one process should output the filter
266                 #endif
267                 
268                 ofstream outFilter;
269                 
270                 //prevent giantic file name
271                 string filterFile;
272                 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
273                 else {  filterFile = outputDir + filterFileName + ".filter";  }
274                 
275                 m->openOutputFile(filterFile, outFilter);
276                 outFilter << filter << endl;
277                 outFilter.close();
278                 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
279                 
280                 #ifdef USE_MPI
281                         }
282                 #endif
283                 
284                 ////////////run filter/////////////////
285                 
286                 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
287                 
288                 filterSequences();
289                 
290                 m->mothurOutEndLine();  m->mothurOutEndLine();
291                                         
292                 int filteredLength = 0;
293                 for(int i=0;i<alignmentLength;i++){
294                         if(filter[i] == '1'){   filteredLength++;       }
295                 }
296                 
297                 if (m->control_pressed) {  outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }  return 0; }
298
299                 
300                 m->mothurOutEndLine();
301                 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
302                 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
303                 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
304                 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
305                 
306                 
307                 m->mothurOutEndLine();
308                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
309                 for(int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();      }
310                 m->mothurOutEndLine();
311                 
312                 return 0;
313                 
314         }
315         catch(exception& e) {
316                 m->errorOut(e, "FilterSeqsCommand", "execute");
317                 exit(1);
318         }
319 }
320 /**************************************************************************************/
321 int FilterSeqsCommand::filterSequences() {      
322         try {
323                 
324                 numSeqs = 0;
325                 
326                 for (int s = 0; s < fastafileNames.size(); s++) {
327                         
328                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
329                                 
330                                 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
331 #ifdef USE_MPI  
332                                 int pid, numSeqsPerProcessor, num; 
333                                 int tag = 2001;
334                                 vector<unsigned long int>MPIPos;
335                                                 
336                                 MPI_Status status; 
337                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
338                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
339                                 
340                                 MPI_File outMPI;
341                                 MPI_File inMPI;
342                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
343                                 int inMode=MPI_MODE_RDONLY; 
344                                 
345                                 char outFilename[1024];
346                                 strcpy(outFilename, filteredFasta.c_str());
347                         
348                                 char inFileName[1024];
349                                 strcpy(inFileName, fastafileNames[s].c_str());
350                                 
351                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
352                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
353
354                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
355
356                                 if (pid == 0) { //you are the root process 
357                                         
358                                         MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
359                                         numSeqs += num;
360                                         
361                                         //send file positions to all processes
362                                         for(int i = 1; i < processors; i++) { 
363                                                 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
364                                                 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
365                                         }
366                                         
367                                         //figure out how many sequences you have to do
368                                         numSeqsPerProcessor = num / processors;
369                                         int startIndex =  pid * numSeqsPerProcessor;
370                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
371                                         
372                                 
373                                         //do your part
374                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
375                                         
376                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
377                                         
378                                         //wait on chidren
379                                         for(int i = 1; i < processors; i++) { 
380                                                 char buf[5];
381                                                 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); 
382                                         }
383                                         
384                                 }else { //you are a child process
385                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
386                                         MPIPos.resize(num+1);
387                                         numSeqs += num;
388                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
389                                         
390                                         //figure out how many sequences you have to align
391                                         numSeqsPerProcessor = num / processors;
392                                         int startIndex =  pid * numSeqsPerProcessor;
393                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
394                                         
395                                         
396                                         //align your part
397                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);           
398                                         
399                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
400                                         
401                                         char buf[5];
402                                         strcpy(buf, "done"); 
403                                         
404                                         //tell parent you are done.
405                                         MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
406                                 }
407                                 
408                                 MPI_File_close(&outMPI);
409                                 MPI_File_close(&inMPI);
410                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
411                                 
412 #else
413                         vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
414                                 
415                         for (int i = 0; i < (positions.size()-1); i++) {
416                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
417                         }       
418                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
419                                 if(processors == 1){
420                                         int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
421                                         numSeqs += numFastaSeqs;
422                                 }else{
423                                         int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]); 
424                                         numSeqs += numFastaSeqs;
425                                 
426                                         rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
427                                 
428                                         //append fasta files
429                                         for(int i=1;i<processors;i++){
430                                                 m->appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
431                                                 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
432                                         }
433                                 }
434                                 
435                                 if (m->control_pressed) {  return 1; }
436                 #else
437                                 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
438                                 numSeqs += numFastaSeqs;
439
440                                 if (m->control_pressed) {  return 1; }
441                 #endif
442 #endif
443                         outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
444                 }
445
446                 return 0;
447         }
448         catch(exception& e) {
449                 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
450                 exit(1);
451         }
452 }
453 #ifdef USE_MPI
454 /**************************************************************************************/
455 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) { 
456         try {
457                 string outputString = "";
458                 int count = 0;
459                 MPI_Status status; 
460                 
461                 for(int i=0;i<num;i++){
462                 
463                         if (m->control_pressed) { return 0; }
464                 
465                         //read next sequence
466                         int length = MPIPos[start+i+1] - MPIPos[start+i];
467                         char* buf4 = new char[length];
468                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
469                         
470                         string tempBuf = buf4;
471                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
472                         istringstream iss (tempBuf,istringstream::in);
473                         delete buf4;
474         
475                         Sequence seq(iss);  m->gobble(iss);
476                         
477                         if (seq.getName() != "") {
478                                 string align = seq.getAligned();
479                                 string filterSeq = "";
480                                         
481                                 for(int j=0;j<alignmentLength;j++){
482                                         if(filter[j] == '1'){
483                                                 filterSeq += align[j];
484                                         }
485                                 }
486                                 
487                                 count++;
488                                 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
489                                 
490                                 if(count % 10 == 0){ //output to file 
491                                         //send results to parent
492                                         int length = outputString.length();
493                                         char* buf = new char[length];
494                                         memcpy(buf, outputString.c_str(), length);
495                                 
496                                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
497                                         outputString = "";
498                                         delete buf;
499                                 }
500
501                         }
502                         
503                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
504                 }
505                 
506                 if(outputString != ""){ //output to file 
507                         //send results to parent
508                         int length = outputString.length();
509                         char* buf = new char[length];
510                         memcpy(buf, outputString.c_str(), length);
511                         
512                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
513                         outputString = "";
514                         delete buf;
515                 }
516                 
517                 if((num) % 100 != 0){   cout << (num) << endl;   m->mothurOutJustToLog(toString(num) + "\n");   }
518                         
519                 return 0;
520         }
521         catch(exception& e) {
522                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
523                 exit(1);
524         }
525 }
526 #endif
527 /**************************************************************************************/
528 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {      
529         try {
530                 ofstream out;
531                 m->openOutputFile(outputFilename, out);
532                 
533                 ifstream in;
534                 m->openInputFile(inputFilename, in);
535                                 
536                 in.seekg(filePos->start);
537
538                 bool done = false;
539                 int count = 0;
540         
541                 while (!done) {
542                                 
543                                 if (m->control_pressed) { in.close(); out.close(); return 0; }
544                                 
545                                 Sequence seq(in); m->gobble(in);
546                                 if (seq.getName() != "") {
547                                         string align = seq.getAligned();
548                                         string filterSeq = "";
549                                         
550                                         for(int j=0;j<alignmentLength;j++){
551                                                 if(filter[j] == '1'){
552                                                         filterSeq += align[j];
553                                                 }
554                                         }
555                                         
556                                         out << '>' << seq.getName() << endl << filterSeq << endl;
557                                 count++;
558                         }
559                         
560                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
561                                 unsigned long int pos = in.tellg();
562                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
563                         #else
564                                 if (in.eof()) { break; }
565                         #endif
566                         
567                         //report progress
568                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
569                 }
570                 //report progress
571                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
572                 
573                 
574                 out.close();
575                 in.close();
576                 
577                 return count;
578         }
579         catch(exception& e) {
580                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
581                 exit(1);
582         }
583 }
584 /**************************************************************************************************/
585
586 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
587         try {
588 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
589                 int process = 0;
590                 int num = 0;
591                 processIDS.clear();
592                 
593                 //loop through and create all the processes you want
594                 while (process != processors) {
595                         int pid = fork();
596                         
597                         if (pid > 0) {
598                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
599                                 process++;
600                         }else if (pid == 0){
601                                 string filteredFasta = filename + toString(getpid()) + ".temp";
602                                 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
603                                 
604                                 //pass numSeqs to parent
605                                 ofstream out;
606                                 string tempFile = filename +  toString(getpid()) + ".num.temp";
607                                 m->openOutputFile(tempFile, out);
608                                 out << num << endl;
609                                 out.close();
610                                 
611                                 exit(0);
612                         }else { 
613                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
614                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
615                                 exit(0);
616                         }
617                 }
618                 
619                 //force parent to wait until all the processes are done
620                 for (int i=0;i<processors;i++) { 
621                         int temp = processIDS[i];
622                         wait(&temp);
623                 }       
624                                         
625                 for (int i = 0; i < processIDS.size(); i++) {
626                         ifstream in;
627                         string tempFile =  filename + toString(processIDS[i]) + ".num.temp";
628                         m->openInputFile(tempFile, in);
629                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
630                         in.close(); remove(tempFile.c_str());
631                 }
632
633                 
634                 return num;
635 #endif          
636         }
637         catch(exception& e) {
638                 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
639                 exit(1);
640         }
641 }
642 /**************************************************************************************/
643 string FilterSeqsCommand::createFilter() {      
644         try {
645                 string filterString = "";                       
646                 Filters F;
647                 
648                 if (soft != 0)                  {  F.setSoft(soft);             }
649                 if (trump != '*')               {  F.setTrump(trump);   }
650                 
651                 F.setLength(alignmentLength);
652                 
653                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
654                         F.initialize();
655                 }
656                 
657                 if(hard.compare("") != 0)       {       F.doHard(hard);         }
658                 else                                            {       F.setFilter(string(alignmentLength, '1'));      }
659                 
660                 numSeqs = 0;
661                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
662                         for (int s = 0; s < fastafileNames.size(); s++) {
663                         
664                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
665                         
666 #ifdef USE_MPI  
667                                 int pid, numSeqsPerProcessor, num; 
668                                 int tag = 2001;
669                                 vector<unsigned long int> MPIPos;
670                                 
671                                 MPI_Status status; 
672                                 MPI_File inMPI; 
673                                 MPI_Comm_size(MPI_COMM_WORLD, &processors);
674                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
675                                                         
676                                 //char* tempFileName = new char(fastafileNames[s].length());
677                                 //tempFileName = &(fastafileNames[s][0]);
678                                 
679                                 char tempFileName[1024];
680                                 strcpy(tempFileName, fastafileNames[s].c_str());
681                 
682                                 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
683                                 
684                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
685                                 
686                                 if (pid == 0) { //you are the root process
687                                                 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
688                                                 numSeqs += num;
689                                                 
690                                                 //send file positions to all processes
691                                                 for(int i = 1; i < processors; i++) { 
692                                                         MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
693                                                         MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
694                                                 }
695                                                                 
696                                                 //figure out how many sequences you have to do
697                                                 numSeqsPerProcessor = num / processors;
698                                                 int startIndex =  pid * numSeqsPerProcessor;
699                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
700                                                 
701                                 
702                                                 //do your part
703                                                 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
704                                                 
705                                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
706                                                                                                 
707                                 }else { //i am the child process
708                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
709                                         MPIPos.resize(num+1);
710                                         numSeqs += num;
711                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
712                                         
713                                         //figure out how many sequences you have to align
714                                         numSeqsPerProcessor = num / processors;
715                                         int startIndex =  pid * numSeqsPerProcessor;
716                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
717                                         
718                                         
719                                         //do your part
720                                         MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI,  MPIPos);
721                                         
722                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
723                                 }
724                                 
725                                 MPI_File_close(&inMPI);
726                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
727                                 
728 #else
729                 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
730                                 
731                 for (int i = 0; i < (positions.size()-1); i++) {
732                         lines.push_back(new linePair(positions[i], positions[(i+1)]));
733                 }       
734                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
735                                 if(processors == 1){
736                                         int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
737                                         numSeqs += numFastaSeqs;
738                                 }else{
739                                         int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]); 
740                                         numSeqs += numFastaSeqs;
741                                 }
742                                 
743                                 if (m->control_pressed) {  return filterString; }
744                 #else
745                                 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
746                                 numSeqs += numFastaSeqs;
747                                 if (m->control_pressed) {  return filterString; }
748                 #endif
749 #endif
750                         
751                         }
752                 }
753
754
755 #ifdef USE_MPI  
756                 int pid;
757                 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
758                 MPI_Status status;
759                 
760                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
761                 
762                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
763                         
764                         if (pid == 0) { //only one process should output the filter
765                         
766                                 vector<int> temp; temp.resize(alignmentLength+1);
767                                                                 
768                                 //get the frequencies from the child processes
769                                 for(int i = 1; i < processors; i++) { 
770                                 
771                                         for (int j = 0; j < 5; j++) {
772                                         
773                                                 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status); 
774                                                 int receiveTag = temp[temp.size()-1];  //child process added a int to the end to indicate what letter count this is for
775                                                 
776                                                 if (receiveTag == Atag) { //you are recieveing the A frequencies
777                                                         for (int k = 0; k < alignmentLength; k++) {             F.a[k] += temp[k];      }
778                                                 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
779                                                         for (int k = 0; k < alignmentLength; k++) {             F.t[k] += temp[k];      }
780                                                 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
781                                                         for (int k = 0; k < alignmentLength; k++) {             F.c[k] += temp[k];      }
782                                                 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
783                                                         for (int k = 0; k < alignmentLength; k++) {             F.g[k] += temp[k];      }
784                                                 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
785                                                         for (int k = 0; k < alignmentLength; k++) {             F.gap[k] += temp[k];    }
786                                                 }
787                                         }
788                                 } 
789                         }else{
790                         
791                                 //send my fequency counts
792                                 F.a.push_back(Atag);
793                                 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
794                                 F.t.push_back(Ttag);
795                                 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
796                                 F.c.push_back(Ctag);
797                                 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
798                                 F.g.push_back(Gtag);
799                                 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
800                                 F.gap.push_back(Gaptag);
801                                 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
802                         }
803                         
804                 }
805                 
806                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
807                 
808                 if (pid == 0) { //only one process should output the filter
809 #endif
810
811                 F.setNumSeqs(numSeqs);
812                 if(m->isTrue(vertical) == 1)    {       F.doVertical(); }
813                 if(soft != 0)                           {       F.doSoft();             }
814                 filterString = F.getFilter();
815                 
816 #ifdef USE_MPI
817                 //send filter string to kids
818                 //for(int i = 1; i < processors; i++) { 
819                 //      MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
820                 //}
821                 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
822         }else{
823                 //recieve filterString
824                 char* tempBuf = new char[alignmentLength];
825                 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
826                 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
827                 
828                 filterString = tempBuf;
829                 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength);  }
830                 delete tempBuf; 
831         }
832         
833         MPI_Barrier(MPI_COMM_WORLD);
834 #endif
835                                 
836                 return filterString;
837         }
838         catch(exception& e) {
839                 m->errorOut(e, "FilterSeqsCommand", "createFilter");
840                 exit(1);
841         }
842 }
843 /**************************************************************************************/
844 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {     
845         try {
846                 
847                 ifstream in;
848                 m->openInputFile(filename, in);
849                                 
850                 in.seekg(filePos->start);
851
852                 bool done = false;
853                 int count = 0;
854         
855                 while (!done) {
856                                 
857                         if (m->control_pressed) { in.close(); return 1; }
858                                         
859                         Sequence seq(in); m->gobble(in);
860                         if (seq.getName() != "") {
861                                         if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true;  }
862                                         
863                                         if(trump != '*')                        {       F.doTrump(seq);         }
864                                         if(m->isTrue(vertical) || soft != 0)    {       F.getFreqs(seq);        }
865                                         cout.flush();
866                                         count++;
867                         }
868                         
869                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
870                                 unsigned long int pos = in.tellg();
871                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
872                         #else
873                                 if (in.eof()) { break; }
874                         #endif
875                         
876                         //report progress
877                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
878                 }
879                 //report progress
880                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
881                 in.close();
882                 
883                 return count;
884         }
885         catch(exception& e) {
886                 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
887                 exit(1);
888         }
889 }
890 #ifdef USE_MPI
891 /**************************************************************************************/
892 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {    
893         try {
894                 
895                 MPI_Status status; 
896                 int pid;
897                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
898                 
899                 for(int i=0;i<num;i++){
900                         
901                         if (m->control_pressed) { return 0; }
902                         
903                         //read next sequence
904                         int length = MPIPos[start+i+1] - MPIPos[start+i];
905         
906                         char* buf4 = new char[length];
907                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
908                         
909                         string tempBuf = buf4;
910                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
911                         istringstream iss (tempBuf,istringstream::in);
912                         delete buf4;
913
914                         Sequence seq(iss);  
915
916                         if (seq.getAligned().length() != alignmentLength) {  cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1);  }
917                         
918                         if(trump != '*'){       F.doTrump(seq); }
919                         if(m->isTrue(vertical) || soft != 0){   F.getFreqs(seq);        }
920                         cout.flush();
921                                                 
922                         //report progress
923                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
924                 }
925                 
926                 //report progress
927                 if((num) % 100 != 0){   cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n");       }
928                 
929                 return 0;
930         }
931         catch(exception& e) {
932                 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
933                 exit(1);
934         }
935 }
936 #endif
937 /**************************************************************************************************/
938
939 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
940         try {
941 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
942                 int process = 1;
943                 int num = 0;
944                 processIDS.clear();
945                 
946                 //loop through and create all the processes you want
947                 while (process != processors) {
948                         int pid = fork();
949                         
950                         if (pid > 0) {
951                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
952                                 process++;
953                         }else if (pid == 0){
954                                 //reset child's filter counts to 0;
955                                 F.a.clear(); F.a.resize(alignmentLength, 0);
956                                 F.t.clear(); F.t.resize(alignmentLength, 0);
957                                 F.g.clear(); F.g.resize(alignmentLength, 0);
958                                 F.c.clear(); F.c.resize(alignmentLength, 0);
959                                 F.gap.clear(); F.gap.resize(alignmentLength, 0);
960                                 
961                                 num = driverCreateFilter(F, filename, lines[process]);
962                                 
963                                 //write out filter counts to file
964                                 filename += toString(getpid()) + "filterValues.temp";
965                                 ofstream out;
966                                 m->openOutputFile(filename, out);
967                                 
968                                 out << num << endl;
969                                 out << F.getFilter() << endl;
970                                 for (int k = 0; k < alignmentLength; k++) {             out << F.a[k] << '\t'; }  out << endl;
971                                 for (int k = 0; k < alignmentLength; k++) {             out << F.t[k] << '\t'; }  out << endl;
972                                 for (int k = 0; k < alignmentLength; k++) {             out << F.g[k] << '\t'; }  out << endl;
973                                 for (int k = 0; k < alignmentLength; k++) {             out << F.c[k] << '\t'; }  out << endl;
974                                 for (int k = 0; k < alignmentLength; k++) {             out << F.gap[k] << '\t'; }  out << endl;
975
976                                 //cout << F.getFilter() << endl;
977                                 out.close();
978                                 
979                                 exit(0);
980                         }else { 
981                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
982                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
983                                 exit(0);
984                         }
985                 }
986                 
987                 //parent do your part
988                 num = driverCreateFilter(F, filename, lines[0]);
989                 
990                 //force parent to wait until all the processes are done
991                 for (int i=0;i<(processors-1);i++) { 
992                         int temp = processIDS[i];
993                         wait(&temp);
994                 }
995                 
996                 //parent reads in and combines Filter info
997                 for (int i = 0; i < processIDS.size(); i++) {
998                         string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
999                         ifstream in;
1000                         m->openInputFile(tempFilename, in);
1001                         
1002                         int temp, tempNum;
1003                         string tempFilterString;
1004
1005                         in >> tempNum; m->gobble(in); num += tempNum;
1006
1007                         in >> tempFilterString;
1008                         F.mergeFilter(tempFilterString);
1009
1010                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.a[k] += temp; }           m->gobble(in);
1011                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.t[k] += temp; }           m->gobble(in);
1012                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.g[k] += temp; }           m->gobble(in);
1013                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.c[k] += temp; }           m->gobble(in);
1014                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.gap[k] += temp; } m->gobble(in);
1015                                 
1016                         in.close();
1017                         remove(tempFilename.c_str());
1018                 }
1019                 
1020                 return num;
1021 #endif          
1022         }
1023         catch(exception& e) {
1024                 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1025                 exit(1);
1026         }
1027 }
1028 /**************************************************************************************/