]> git.donarmstrong.com Git - mothur.git/blob - filterseqscommand.cpp
addc5f338c9202811a7d6e341167c8a097e03ea0
[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                 string filterFile = outputDir + filterFileName + ".filter";
202                 m->openOutputFile(filterFile, outFilter);
203                 outFilter << filter << endl;
204                 outFilter.close();
205                 outputNames.push_back(filterFile);
206                 
207                 #ifdef USE_MPI
208                         }
209                 #endif
210                 
211                 ////////////run filter/////////////////
212                 
213                 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
214                 
215                 filterSequences();
216                 
217                 m->mothurOutEndLine();  m->mothurOutEndLine();
218                                         
219                 int filteredLength = 0;
220                 for(int i=0;i<alignmentLength;i++){
221                         if(filter[i] == '1'){   filteredLength++;       }
222                 }
223                 
224                 if (m->control_pressed) {  for(int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); }  return 0; }
225
226                 
227                 m->mothurOutEndLine();
228                 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
229                 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
230                 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
231                 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
232                 
233                 
234                 m->mothurOutEndLine();
235                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
236                 for(int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();      }
237                 m->mothurOutEndLine();
238                 
239                 return 0;
240                 
241         }
242         catch(exception& e) {
243                 m->errorOut(e, "FilterSeqsCommand", "execute");
244                 exit(1);
245         }
246 }
247 /**************************************************************************************/
248 int FilterSeqsCommand::filterSequences() {      
249         try {
250                 
251                 numSeqs = 0;
252                 
253                 for (int s = 0; s < fastafileNames.size(); s++) {
254                         
255                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
256                                 
257                                 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
258 #ifdef USE_MPI  
259                                 int pid, start, end, numSeqsPerProcessor, num; 
260                                 int tag = 2001;
261                                 vector<unsigned long int>MPIPos;
262                                                 
263                                 MPI_Status status; 
264                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
265                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
266                                 
267                                 MPI_File outMPI;
268                                 MPI_File tempMPI;
269                                 MPI_File inMPI;
270                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
271                                 int inMode=MPI_MODE_RDONLY; 
272                                 
273                                 char outFilename[1024];
274                                 strcpy(outFilename, filteredFasta.c_str());
275                         
276                                 char inFileName[1024];
277                                 strcpy(inFileName, fastafileNames[s].c_str());
278                                 
279                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
280                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
281
282                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
283
284                                 if (pid == 0) { //you are the root process 
285                                         
286                                         MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
287                                         numSeqs += num;
288                                         
289                                         //send file positions to all processes
290                                         for(int i = 1; i < processors; i++) { 
291                                                 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
292                                                 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
293                                         }
294                                         
295                                         //figure out how many sequences you have to do
296                                         numSeqsPerProcessor = num / processors;
297                                         int startIndex =  pid * numSeqsPerProcessor;
298                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
299                                         
300                                 
301                                         //do your part
302                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
303                                         
304                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
305                                         
306                                         //wait on chidren
307                                         for(int i = 1; i < processors; i++) { 
308                                                 char buf[4];
309                                                 MPI_Recv(buf, 4, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); 
310                                         }
311                                         
312                                 }else { //you are a child process
313                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
314                                         MPIPos.resize(num+1);
315                                         numSeqs += num;
316                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
317                                         
318                                         //figure out how many sequences you have to align
319                                         numSeqsPerProcessor = num / processors;
320                                         int startIndex =  pid * numSeqsPerProcessor;
321                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
322                                         
323                                         
324                                         //align your part
325                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);           
326                                         
327                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
328                                         
329                                         char buf[4];
330                                         strcpy(buf, "done"); 
331                                         
332                                         //tell parent you are done.
333                                         MPI_Send(buf, 4, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
334                                 }
335                                 
336                                 MPI_File_close(&outMPI);
337                                 MPI_File_close(&inMPI);
338                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
339                                 
340 #else
341                         vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
342                                 
343                         for (int i = 0; i < (positions.size()-1); i++) {
344                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
345                         }       
346                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
347                                 if(processors == 1){
348                                         int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
349                                         numSeqs += numFastaSeqs;
350                                 }else{
351                                         int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s]); 
352                                         numSeqs += numFastaSeqs;
353                                 
354                                         rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
355                                 
356                                         //append fasta files
357                                         for(int i=1;i<processors;i++){
358                                                 m->appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
359                                                 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
360                                         }
361                                 }
362                                 
363                                 if (m->control_pressed) {  return 1; }
364                 #else
365                                 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
366                                 numSeqs += numFastaSeqs;
367
368                                 if (m->control_pressed) {  return 1; }
369                 #endif
370 #endif
371                         outputNames.push_back(filteredFasta);
372                 }
373
374                 return 0;
375         }
376         catch(exception& e) {
377                 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
378                 exit(1);
379         }
380 }
381 #ifdef USE_MPI
382 /**************************************************************************************/
383 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) { 
384         try {
385                 string outputString = "";
386                 int count = 0;
387                 MPI_Status status; 
388                 
389                 for(int i=0;i<num;i++){
390                 
391                         if (m->control_pressed) { return 0; }
392                 
393                         //read next sequence
394                         int length = MPIPos[start+i+1] - MPIPos[start+i];
395                         char* buf4 = new char[length];
396                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
397                         
398                         string tempBuf = buf4;
399                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
400                         istringstream iss (tempBuf,istringstream::in);
401                         delete buf4;
402         
403                         Sequence seq(iss);  m->gobble(iss);
404                         
405                         if (seq.getName() != "") {
406                                 string align = seq.getAligned();
407                                 string filterSeq = "";
408                                         
409                                 for(int j=0;j<alignmentLength;j++){
410                                         if(filter[j] == '1'){
411                                                 filterSeq += align[j];
412                                         }
413                                 }
414                                 
415                                 count++;
416                                 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
417                                 
418                                 if(count % 10 == 0){ //output to file 
419                                         //send results to parent
420                                         int length = outputString.length();
421                                         char* buf = new char[length];
422                                         memcpy(buf, outputString.c_str(), length);
423                                 
424                                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
425                                         outputString = "";
426                                         delete buf;
427                                 }
428
429                         }
430                         
431                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
432                 }
433                 
434                 if(outputString != ""){ //output to file 
435                         //send results to parent
436                         int length = outputString.length();
437                         char* buf = new char[length];
438                         memcpy(buf, outputString.c_str(), length);
439                         
440                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
441                         outputString = "";
442                         delete buf;
443                 }
444                 
445                 if((num) % 100 != 0){   cout << (num) << endl;   m->mothurOutJustToLog(toString(num) + "\n");   }
446                         
447                 return 0;
448         }
449         catch(exception& e) {
450                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
451                 exit(1);
452         }
453 }
454 #endif
455 /**************************************************************************************/
456 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {      
457         try {
458                 ofstream out;
459                 m->openOutputFile(outputFilename, out);
460                 
461                 ifstream in;
462                 m->openInputFile(inputFilename, in);
463                                 
464                 in.seekg(filePos->start);
465
466                 bool done = false;
467                 int count = 0;
468         
469                 while (!done) {
470                                 
471                                 if (m->control_pressed) { in.close(); out.close(); return 0; }
472                                 
473                                 Sequence seq(in); m->gobble(in);
474                                 if (seq.getName() != "") {
475                                         string align = seq.getAligned();
476                                         string filterSeq = "";
477                                         
478                                         for(int j=0;j<alignmentLength;j++){
479                                                 if(filter[j] == '1'){
480                                                         filterSeq += align[j];
481                                                 }
482                                         }
483                                         
484                                         out << '>' << seq.getName() << endl << filterSeq << endl;
485                                 count++;
486                         }
487                         
488                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
489                                 unsigned long int pos = in.tellg();
490                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
491                         #else
492                                 if (in.eof()) { break; }
493                         #endif
494                         
495                         //report progress
496                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
497                 }
498                 //report progress
499                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
500                 
501                 
502                 out.close();
503                 in.close();
504                 
505                 return count;
506         }
507         catch(exception& e) {
508                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
509                 exit(1);
510         }
511 }
512 /**************************************************************************************************/
513
514 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
515         try {
516 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
517                 int process = 0;
518                 int num = 0;
519                 processIDS.clear();
520                 
521                 //loop through and create all the processes you want
522                 while (process != processors) {
523                         int pid = fork();
524                         
525                         if (pid > 0) {
526                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
527                                 process++;
528                         }else if (pid == 0){
529                                 string filteredFasta = filename + toString(getpid()) + ".temp";
530                                 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
531                                 
532                                 //pass numSeqs to parent
533                                 ofstream out;
534                                 string tempFile = filename +  toString(getpid()) + ".num.temp";
535                                 m->openOutputFile(tempFile, out);
536                                 out << num << endl;
537                                 out.close();
538                                 
539                                 exit(0);
540                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
541                 }
542                 
543                 //force parent to wait until all the processes are done
544                 for (int i=0;i<processors;i++) { 
545                         int temp = processIDS[i];
546                         wait(&temp);
547                 }       
548                                         
549                 for (int i = 0; i < processIDS.size(); i++) {
550                         ifstream in;
551                         string tempFile =  filename + toString(processIDS[i]) + ".num.temp";
552                         m->openInputFile(tempFile, in);
553                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
554                         in.close(); remove(tempFile.c_str());
555                 }
556
557                 
558                 return num;
559 #endif          
560         }
561         catch(exception& e) {
562                 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
563                 exit(1);
564         }
565 }
566 /**************************************************************************************/
567 string FilterSeqsCommand::createFilter() {      
568         try {
569                 string filterString = "";                       
570                 Filters F;
571                 
572                 if (soft != 0)                  {  F.setSoft(soft);             }
573                 if (trump != '*')               {  F.setTrump(trump);   }
574                 
575                 F.setLength(alignmentLength);
576                 
577                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
578                         F.initialize();
579                 }
580                 
581                 if(hard.compare("") != 0)       {       F.doHard(hard);         }
582                 else                                            {       F.setFilter(string(alignmentLength, '1'));      }
583                 
584                 numSeqs = 0;
585                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
586                         for (int s = 0; s < fastafileNames.size(); s++) {
587                         
588                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
589                         
590 #ifdef USE_MPI  
591                                 int pid, numSeqsPerProcessor, num; 
592                                 int tag = 2001;
593                                 vector<unsigned long int> MPIPos;
594                                 
595                                 MPI_Status status; 
596                                 MPI_File inMPI; 
597                                 MPI_Comm_size(MPI_COMM_WORLD, &processors);
598                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
599                                                         
600                                 //char* tempFileName = new char(fastafileNames[s].length());
601                                 //tempFileName = &(fastafileNames[s][0]);
602                                 
603                                 char tempFileName[1024];
604                                 strcpy(tempFileName, fastafileNames[s].c_str());
605                 
606                                 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
607                                 
608                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
609                                 
610                                 if (pid == 0) { //you are the root process
611                                                 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
612                                                 numSeqs += num;
613                                                 
614                                                 //send file positions to all processes
615                                                 for(int i = 1; i < processors; i++) { 
616                                                         MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
617                                                         MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
618                                                 }
619                                                                 
620                                                 //figure out how many sequences you have to do
621                                                 numSeqsPerProcessor = num / processors;
622                                                 int startIndex =  pid * numSeqsPerProcessor;
623                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
624                                                 
625                                 
626                                                 //do your part
627                                                 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
628                                                 
629                                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
630                                                                                                 
631                                 }else { //i am the child process
632                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
633                                         MPIPos.resize(num+1);
634                                         numSeqs += num;
635                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
636                                         
637                                         //figure out how many sequences you have to align
638                                         numSeqsPerProcessor = num / processors;
639                                         int startIndex =  pid * numSeqsPerProcessor;
640                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
641                                         
642                                         
643                                         //do your part
644                                         MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI,  MPIPos);
645                                         
646                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
647                                 }
648                                 
649                                 MPI_File_close(&inMPI);
650                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
651                                 
652 #else
653                 vector<unsigned long int> positions = m->divideFile(fastafileNames[s], processors);
654                                 
655                 for (int i = 0; i < (positions.size()-1); i++) {
656                         lines.push_back(new linePair(positions[i], positions[(i+1)]));
657                 }       
658                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
659                                 if(processors == 1){
660                                         int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
661                                         numSeqs += numFastaSeqs;
662                                 }else{
663                                         int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]); 
664                                         numSeqs += numFastaSeqs;
665                                 }
666                                 
667                                 if (m->control_pressed) {  return filterString; }
668                 #else
669                                 int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
670                                 numSeqs += numFastaSeqs;
671                                 if (m->control_pressed) {  return filterString; }
672                 #endif
673 #endif
674                         
675                         }
676                 }
677
678
679 #ifdef USE_MPI  
680                 int pid;
681                 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
682                 MPI_Status status;
683                 
684                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
685                 
686                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
687                         
688                         if (pid == 0) { //only one process should output the filter
689                         
690                                 vector<int> temp; temp.resize(alignmentLength+1);
691                                                                 
692                                 //get the frequencies from the child processes
693                                 for(int i = 1; i < processors; i++) { 
694                                 
695                                         for (int j = 0; j < 5; j++) {
696                                         
697                                                 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status); 
698                                                 int receiveTag = temp[temp.size()-1];  //child process added a int to the end to indicate what letter count this is for
699                                                 
700                                                 if (receiveTag == Atag) { //you are recieveing the A frequencies
701                                                         for (int k = 0; k < alignmentLength; k++) {             F.a[k] += temp[k];      }
702                                                 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
703                                                         for (int k = 0; k < alignmentLength; k++) {             F.t[k] += temp[k];      }
704                                                 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
705                                                         for (int k = 0; k < alignmentLength; k++) {             F.c[k] += temp[k];      }
706                                                 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
707                                                         for (int k = 0; k < alignmentLength; k++) {             F.g[k] += temp[k];      }
708                                                 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
709                                                         for (int k = 0; k < alignmentLength; k++) {             F.gap[k] += temp[k];    }
710                                                 }
711                                         }
712                                 } 
713                         }else{
714                         
715                                 //send my fequency counts
716                                 F.a.push_back(Atag);
717                                 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
718                                 F.t.push_back(Ttag);
719                                 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
720                                 F.c.push_back(Ctag);
721                                 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
722                                 F.g.push_back(Gtag);
723                                 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
724                                 F.gap.push_back(Gaptag);
725                                 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
726                         }
727                         
728                 }
729                 
730                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
731                 
732                 if (pid == 0) { //only one process should output the filter
733 #endif
734
735                 F.setNumSeqs(numSeqs);
736                 if(m->isTrue(vertical) == 1)    {       F.doVertical(); }
737                 if(soft != 0)                           {       F.doSoft();             }
738                 filterString = F.getFilter();
739                 
740 #ifdef USE_MPI
741                 //send filter string to kids
742                 //for(int i = 1; i < processors; i++) { 
743                 //      MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
744                 //}
745                 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
746         }else{
747                 //recieve filterString
748                 char* tempBuf = new char[alignmentLength];
749                 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
750                 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
751                 
752                 filterString = tempBuf;
753                 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength);  }
754                 delete tempBuf; 
755         }
756         
757         MPI_Barrier(MPI_COMM_WORLD);
758 #endif
759                                 
760                 return filterString;
761         }
762         catch(exception& e) {
763                 m->errorOut(e, "FilterSeqsCommand", "createFilter");
764                 exit(1);
765         }
766 }
767 /**************************************************************************************/
768 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {     
769         try {
770                 
771                 ifstream in;
772                 m->openInputFile(filename, in);
773                                 
774                 in.seekg(filePos->start);
775
776                 bool done = false;
777                 int count = 0;
778         
779                 while (!done) {
780                                 
781                         if (m->control_pressed) { in.close(); return 1; }
782                                         
783                         Sequence seq(in); m->gobble(in);
784                         if (seq.getName() != "") {
785                                         if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true;  }
786                                         
787                                         if(trump != '*')                        {       F.doTrump(seq);         }
788                                         if(m->isTrue(vertical) || soft != 0)    {       F.getFreqs(seq);        }
789                                         cout.flush();
790                                         count++;
791                         }
792                         
793                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
794                                 unsigned long int pos = in.tellg();
795                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
796                         #else
797                                 if (in.eof()) { break; }
798                         #endif
799                         
800                         //report progress
801                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
802                 }
803                 //report progress
804                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
805                 in.close();
806                 
807                 return count;
808         }
809         catch(exception& e) {
810                 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
811                 exit(1);
812         }
813 }
814 #ifdef USE_MPI
815 /**************************************************************************************/
816 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long int>& MPIPos) {    
817         try {
818                 
819                 MPI_Status status; 
820                 int pid;
821                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
822                 
823                 for(int i=0;i<num;i++){
824                         
825                         if (m->control_pressed) { return 0; }
826                         
827                         //read next sequence
828                         int length = MPIPos[start+i+1] - MPIPos[start+i];
829         
830                         char* buf4 = new char[length];
831                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
832                         
833                         string tempBuf = buf4;
834                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
835                         istringstream iss (tempBuf,istringstream::in);
836                         delete buf4;
837
838                         Sequence seq(iss);  
839
840                         if (seq.getAligned().length() != alignmentLength) {  cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1);  }
841                         
842                         if(trump != '*'){       F.doTrump(seq); }
843                         if(m->isTrue(vertical) || soft != 0){   F.getFreqs(seq);        }
844                         cout.flush();
845                                                 
846                         //report progress
847                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
848                 }
849                 
850                 //report progress
851                 if((num) % 100 != 0){   cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n");       }
852                 
853                 return 0;
854         }
855         catch(exception& e) {
856                 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
857                 exit(1);
858         }
859 }
860 #endif
861 /**************************************************************************************************/
862
863 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
864         try {
865 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
866                 int process = 0;
867                 int num = 0;
868                 processIDS.clear();
869                 
870                 //loop through and create all the processes you want
871                 while (process != processors) {
872                         int pid = fork();
873                         
874                         if (pid > 0) {
875                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
876                                 process++;
877                         }else if (pid == 0){
878                                 //reset child's filter counts to 0;
879                                 F.a.clear(); F.a.resize(alignmentLength, 0);
880                                 F.t.clear(); F.t.resize(alignmentLength, 0);
881                                 F.g.clear(); F.g.resize(alignmentLength, 0);
882                                 F.c.clear(); F.c.resize(alignmentLength, 0);
883                                 F.gap.clear(); F.gap.resize(alignmentLength, 0);
884                                 
885                                 num = driverCreateFilter(F, filename, lines[process]);
886                                 
887                                 //write out filter counts to file
888                                 filename += toString(getpid()) + "filterValues.temp";
889                                 ofstream out;
890                                 m->openOutputFile(filename, out);
891                                 
892                                 out << num << endl;
893                                 out << F.getFilter() << endl;
894                                 for (int k = 0; k < alignmentLength; k++) {             out << F.a[k] << '\t'; }  out << endl;
895                                 for (int k = 0; k < alignmentLength; k++) {             out << F.t[k] << '\t'; }  out << endl;
896                                 for (int k = 0; k < alignmentLength; k++) {             out << F.g[k] << '\t'; }  out << endl;
897                                 for (int k = 0; k < alignmentLength; k++) {             out << F.c[k] << '\t'; }  out << endl;
898                                 for (int k = 0; k < alignmentLength; k++) {             out << F.gap[k] << '\t'; }  out << endl;
899
900                                 //cout << F.getFilter() << endl;
901                                 out.close();
902                                 
903                                 exit(0);
904                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
905                 }
906                 
907                 //force parent to wait until all the processes are done
908                 for (int i=0;i<processors;i++) { 
909                         int temp = processIDS[i];
910                         wait(&temp);
911                 }
912                 
913                 //parent reads in and combine Filter info
914                 for (int i = 0; i < processIDS.size(); i++) {
915                         string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
916                         ifstream in;
917                         m->openInputFile(tempFilename, in);
918                         
919                         int temp, tempNum;
920                         string tempFilterString;
921
922                         in >> tempNum; m->gobble(in); num += tempNum;
923
924                         in >> tempFilterString;
925                         F.mergeFilter(tempFilterString);
926
927                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.a[k] += temp; }           m->gobble(in);
928                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.t[k] += temp; }           m->gobble(in);
929                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.g[k] += temp; }           m->gobble(in);
930                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.c[k] += temp; }           m->gobble(in);
931                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.gap[k] += temp; } m->gobble(in);
932                                 
933                         in.close();
934                         remove(tempFilename.c_str());
935                 }
936                 
937                 return num;
938 #endif          
939         }
940         catch(exception& e) {
941                 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
942                 exit(1);
943         }
944 }
945 /**************************************************************************************/