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