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