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