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