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