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