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