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