]> git.donarmstrong.com Git - mothur.git/blob - filterseqscommand.cpp
changed how we count sequences in a fastafile to allow for '>' in sequence names
[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()];
260                                 //memcpy(outFilename, filteredFasta.c_str(), filteredFasta.length());
261                                 
262                                 char outFilename[1024];
263                                 strcpy(outFilename, filteredFasta.c_str());
264
265                                 //char* inFileName = new char[fastafileNames[s].length()];
266                                 //memcpy(inFileName, fastafileNames[s].c_str(), fastafileNames[s].length());
267                                 
268                                 char inFileName[1024];
269                                 strcpy(inFileName, fastafileNames[s].c_str());
270                                 
271                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
272                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
273                                 
274                                 //delete inFileName;
275                                 //delete outFilename;
276
277                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
278
279                                 if (pid == 0) { //you are the root process 
280                                         
281                                         MPIPos = setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
282                                         numSeqs += num;
283                                         
284                                         //send file positions to all processes
285                                         for(int i = 1; i < processors; i++) { 
286                                                 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
287                                                 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
288                                         }
289                                         
290                                         //figure out how many sequences you have to do
291                                         numSeqsPerProcessor = num / processors;
292                                         int startIndex =  pid * numSeqsPerProcessor;
293                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
294                                         
295                                 
296                                         //do your part
297                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
298                                         
299                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
300                                         
301                                         //wait on chidren
302                                         for(int i = 1; i < processors; i++) { 
303                                                 char buf[4];
304                                                 MPI_Recv(buf, 4, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); 
305                                         }
306                                         
307                                 }else { //you are a child process
308                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
309                                         MPIPos.resize(num+1);
310                                         numSeqs += num;
311                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
312                                         
313                                         //figure out how many sequences you have to align
314                                         numSeqsPerProcessor = num / processors;
315                                         int startIndex =  pid * numSeqsPerProcessor;
316                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
317                                         
318                                         
319                                         //align your part
320                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);           
321                                         
322                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
323                                         
324                                         char buf[4];
325                                         strcpy(buf, "done"); 
326                                         
327                                         //tell parent you are done.
328                                         MPI_Send(buf, 4, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
329                                 }
330                                 
331                                 MPI_File_close(&outMPI);
332                                 MPI_File_close(&inMPI);
333                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
334                                 
335 #else
336                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
337                                 if(processors == 1){
338                                         ifstream inFASTA;
339                                         int numFastaSeqs;
340                                         openInputFile(fastafileNames[s], inFASTA);
341                                         getNumSeqs(inFASTA, numFastaSeqs);
342                                         inFASTA.close();
343                                         
344                                         lines.push_back(new linePair(0, numFastaSeqs));
345                                         
346                                         numSeqs += numFastaSeqs;
347                                         
348                                         driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
349                                 }else{
350                                         setLines(fastafileNames[s]);
351                                         createProcessesRunFilter(filter, fastafileNames[s]); 
352                                 
353                                         rename((fastafileNames[s] + toString(processIDS[0]) + ".temp").c_str(), filteredFasta.c_str());
354                                 
355                                         //append fasta files
356                                         for(int i=1;i<processors;i++){
357                                                 appendFiles((fastafileNames[s] + toString(processIDS[i]) + ".temp"), filteredFasta);
358                                                 remove((fastafileNames[s] + toString(processIDS[i]) + ".temp").c_str());
359                                         }
360                                 }
361                                 
362                                 if (m->control_pressed) {  return 1; }
363                 #else
364                                 ifstream inFASTA;
365                                 int numFastaSeqs;
366                                 openInputFile(fastafileNames[s], inFASTA);
367                                 getNumSeqs(inFASTA, numFastaSeqs);
368                                 inFASTA.close();
369                                         
370                                 lines.push_back(new linePair(0, numFastaSeqs));
371                                 
372                                 numSeqs += numFastaSeqs;
373                                 
374                                 driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
375
376                                 if (m->control_pressed) {  return 1; }
377                 #endif
378 #endif
379                         outputNames.push_back(filteredFasta);
380                 }
381
382                 return 0;
383         }
384         catch(exception& e) {
385                 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
386                 exit(1);
387         }
388 }
389 #ifdef USE_MPI
390 /**************************************************************************************/
391 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<long>& MPIPos) {      
392         try {
393                 string outputString = "";
394                 int count = 0;
395                 MPI_Status status; 
396                 
397                 for(int i=0;i<num;i++){
398                 
399                         if (m->control_pressed) { return 0; }
400                 
401                         //read next sequence
402                         int length = MPIPos[start+i+1] - MPIPos[start+i];
403                         char* buf4 = new char[length];
404                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
405                         
406                         string tempBuf = buf4;
407                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
408                         istringstream iss (tempBuf,istringstream::in);
409                         delete buf4;
410         
411                         Sequence seq(iss);  gobble(iss);
412                         
413                         if (seq.getName() != "") {
414                                 string align = seq.getAligned();
415                                 string filterSeq = "";
416                                         
417                                 for(int j=0;j<alignmentLength;j++){
418                                         if(filter[j] == '1'){
419                                                 filterSeq += align[j];
420                                         }
421                                 }
422                                 
423                                 count++;
424                                 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
425                                 
426                                 if(count % 10 == 0){ //output to file 
427                                         //send results to parent
428                                         int length = outputString.length();
429                                         char* buf = new char[length];
430                                         memcpy(buf, outputString.c_str(), length);
431                                 
432                                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
433                                         outputString = "";
434                                         delete buf;
435                                 }
436
437                         }
438                         
439                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
440                 }
441                 
442                 if(outputString != ""){ //output to file 
443                         //send results to parent
444                         int length = outputString.length();
445                         char* buf = new char[length];
446                         memcpy(buf, outputString.c_str(), length);
447                         
448                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
449                         outputString = "";
450                         delete buf;
451                 }
452                 
453                 if((num) % 100 != 0){   cout << (num) << endl;   m->mothurOutJustToLog(toString(num) + "\n");   }
454                         
455                 return 0;
456         }
457         catch(exception& e) {
458                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
459                 exit(1);
460         }
461 }
462 #endif
463 /**************************************************************************************/
464 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* line) { 
465         try {
466                 ofstream out;
467                 openOutputFile(outputFilename, out);
468                 
469                 ifstream in;
470                 openInputFile(inputFilename, in);
471                                 
472                 in.seekg(line->start);
473                 
474                 for(int i=0;i<line->num;i++){
475                                 
476                                 if (m->control_pressed) { in.close(); out.close(); return 0; }
477                                 
478                                 Sequence seq(in);
479                                 if (seq.getName() != "") {
480                                         string align = seq.getAligned();
481                                         string filterSeq = "";
482                                         
483                                         for(int j=0;j<alignmentLength;j++){
484                                                 if(filter[j] == '1'){
485                                                         filterSeq += align[j];
486                                                 }
487                                         }
488                                         
489                                         out << '>' << seq.getName() << endl << filterSeq << endl;
490                                 }
491                                 gobble(in);
492                 }
493                 out.close();
494                 in.close();
495                 
496                 return 0;
497         }
498         catch(exception& e) {
499                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
500                 exit(1);
501         }
502 }
503 /**************************************************************************************************/
504
505 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename) {
506         try {
507 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
508                 int process = 0;
509                 int exitCommand = 1;
510                 processIDS.clear();
511                 
512                 //loop through and create all the processes you want
513                 while (process != processors) {
514                         int pid = fork();
515                         
516                         if (pid > 0) {
517                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
518                                 process++;
519                         }else if (pid == 0){
520                                 string filteredFasta = filename + toString(getpid()) + ".temp";
521                                 driverRunFilter(F, filteredFasta, filename, lines[process]);
522                                 exit(0);
523                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
524                 }
525                 
526                 //force parent to wait until all the processes are done
527                 for (int i=0;i<processors;i++) { 
528                         int temp = processIDS[i];
529                         wait(&temp);
530                 }
531                 
532                 return exitCommand;
533 #endif          
534         }
535         catch(exception& e) {
536                 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
537                 exit(1);
538         }
539 }
540 /**************************************************************************************/
541 string FilterSeqsCommand::createFilter() {      
542         try {
543                 string filterString = "";                       
544                 Filters F;
545                 
546                 if (soft != 0)                  {  F.setSoft(soft);             }
547                 if (trump != '*')               {  F.setTrump(trump);   }
548                 
549                 F.setLength(alignmentLength);
550                 
551                 if(trump != '*' || isTrue(vertical) || soft != 0){
552                         F.initialize();
553                 }
554                 
555                 if(hard.compare("") != 0)       {       F.doHard(hard);         }
556                 else                                            {       F.setFilter(string(alignmentLength, '1'));      }
557                 
558                 numSeqs = 0;
559                 if(trump != '*' || isTrue(vertical) || soft != 0){
560                         for (int s = 0; s < fastafileNames.size(); s++) {
561                         
562                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
563                         
564 #ifdef USE_MPI  
565                                 int pid, numSeqsPerProcessor, num; 
566                                 int tag = 2001;
567                                 vector<long> MPIPos;
568                                 
569                                 MPI_Status status; 
570                                 MPI_File inMPI; 
571                                 MPI_Comm_size(MPI_COMM_WORLD, &processors);
572                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
573                                                         
574                                 //char* tempFileName = new char(fastafileNames[s].length());
575                                 //tempFileName = &(fastafileNames[s][0]);
576                                 
577                                 char tempFileName[1024];
578                                 strcpy(tempFileName, fastafileNames[s].c_str());
579                 
580                                 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
581                                 
582                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
583                                 
584                                 if (pid == 0) { //you are the root process
585                                                 MPIPos = setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
586                                                 numSeqs += num;
587                                                 
588                                                 //send file positions to all processes
589                                                 for(int i = 1; i < processors; i++) { 
590                                                         MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
591                                                         MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
592                                                 }
593                                                                 
594                                                 //figure out how many sequences you have to do
595                                                 numSeqsPerProcessor = num / processors;
596                                                 int startIndex =  pid * numSeqsPerProcessor;
597                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
598                                                 
599                                 
600                                                 //do your part
601                                                 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
602                                                 
603                                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
604                                                                                                 
605                                 }else { //i am the child process
606                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
607                                         MPIPos.resize(num+1);
608                                         numSeqs += num;
609                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
610                                         
611                                         //figure out how many sequences you have to align
612                                         numSeqsPerProcessor = num / processors;
613                                         int startIndex =  pid * numSeqsPerProcessor;
614                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
615                                         
616                                         
617                                         //do your part
618                                         MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI,  MPIPos);
619                                         
620                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
621                                 }
622                                 
623                                 MPI_File_close(&inMPI);
624                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
625                                 
626 #else
627                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
628                                 if(processors == 1){
629                                         ifstream inFASTA;
630                                         int numFastaSeqs;
631                                         openInputFile(fastafileNames[s], inFASTA);
632                                         getNumSeqs(inFASTA, numFastaSeqs);
633                                         inFASTA.close();
634                                         
635                                         numSeqs += numFastaSeqs;
636                                         
637                                         lines.push_back(new linePair(0, numFastaSeqs));
638                                         
639                                         driverCreateFilter(F, fastafileNames[s], lines[0]);
640                                 }else{
641                                         setLines(fastafileNames[s]);                                    
642                                         createProcessesCreateFilter(F, fastafileNames[s]); 
643                                 }
644                                 
645                                 if (m->control_pressed) {  return filterString; }
646                 #else
647                                 ifstream inFASTA;
648                                 int numFastaSeqs;
649                                 openInputFile(fastafileNames[s], inFASTA);
650                                 getNumSeqs(inFASTA, numFastaSeqs);
651                                 inFASTA.close();
652                                 
653                                 numSeqs += numFastaSeqs;
654                                 
655                                 lines.push_back(new linePair(0, numFastaSeqs));
656                                 
657                                 driverCreateFilter(F, fastafileNames[s], lines[0]);
658                                 if (m->control_pressed) {  return filterString; }
659                 #endif
660 #endif
661                         
662                         }
663                 }
664
665
666 #ifdef USE_MPI  
667                 int pid;
668                 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
669                 MPI_Status status;
670                 
671                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
672                 
673                 if(trump != '*' || isTrue(vertical) || soft != 0){
674                         
675                         if (pid == 0) { //only one process should output the filter
676                         
677                                 vector<int> temp; temp.resize(alignmentLength+1);
678                                                                 
679                                 //get the frequencies from the child processes
680                                 for(int i = 1; i < processors; i++) { 
681                                 
682                                         for (int j = 0; j < 5; j++) {
683                                         
684                                                 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status); 
685                                                 int receiveTag = temp[temp.size()-1];  //child process added a int to the end to indicate what letter count this is for
686                                                 
687                                                 if (receiveTag == Atag) { //you are recieveing the A frequencies
688                                                         for (int k = 0; k < alignmentLength; k++) {             F.a[k] += temp[k];      }
689                                                 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
690                                                         for (int k = 0; k < alignmentLength; k++) {             F.t[k] += temp[k];      }
691                                                 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
692                                                         for (int k = 0; k < alignmentLength; k++) {             F.c[k] += temp[k];      }
693                                                 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
694                                                         for (int k = 0; k < alignmentLength; k++) {             F.g[k] += temp[k];      }
695                                                 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
696                                                         for (int k = 0; k < alignmentLength; k++) {             F.gap[k] += temp[k];    }
697                                                 }
698                                         }
699                                 } 
700                         }else{
701                         
702                                 //send my fequency counts
703                                 F.a.push_back(Atag);
704                                 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
705                                 F.t.push_back(Ttag);
706                                 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
707                                 F.c.push_back(Ctag);
708                                 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
709                                 F.g.push_back(Gtag);
710                                 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
711                                 F.gap.push_back(Gaptag);
712                                 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
713                         }
714                         
715                 }
716                 
717                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
718                 
719                 if (pid == 0) { //only one process should output the filter
720 #endif
721                 F.setNumSeqs(numSeqs);
722                                 
723                 if(isTrue(vertical) == 1)       {       F.doVertical(); }
724                 if(soft != 0)                           {       F.doSoft();             }
725                         
726                 filterString = F.getFilter();
727                 
728 #ifdef USE_MPI
729                 //send filter string to kids
730                 //for(int i = 1; i < processors; i++) { 
731                 //      MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
732                 //}
733                 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
734         }else{
735                 //recieve filterString
736                 char* tempBuf = new char[alignmentLength];
737                 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
738                 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
739                 
740                 filterString = tempBuf;
741                 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength);  }
742                 delete tempBuf; 
743         }
744         
745         MPI_Barrier(MPI_COMM_WORLD);
746 #endif
747                                 
748                 return filterString;
749         }
750         catch(exception& e) {
751                 m->errorOut(e, "FilterSeqsCommand", "createFilter");
752                 exit(1);
753         }
754 }
755 /**************************************************************************************/
756 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* line) {        
757         try {
758                 
759                 ifstream in;
760                 openInputFile(filename, in);
761                                 
762                 in.seekg(line->start);
763                 
764                 for(int i=0;i<line->num;i++){
765                                 
766                         if (m->control_pressed) { in.close(); return 1; }
767                                         
768                         Sequence seq(in);
769                         if (seq.getName() != "") {
770                                         if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true;  }
771                                         
772                                         if(trump != '*'){       F.doTrump(seq); }
773                                         if(isTrue(vertical) || soft != 0){      F.getFreqs(seq);        }
774                                         cout.flush();
775                         }
776                         
777                         //report progress
778                         if((i+1) % 100 == 0){   m->mothurOut(toString(i+1)); m->mothurOutEndLine();             }
779                 }
780                 
781                 //report progress
782                 if((line->num) % 100 != 0){     m->mothurOut(toString(line->num)); m->mothurOutEndLine();               }
783                 
784                 in.close();
785                 
786                 return 0;
787         }
788         catch(exception& e) {
789                 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
790                 exit(1);
791         }
792 }
793 #ifdef USE_MPI
794 /**************************************************************************************/
795 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<long>& MPIPos) { 
796         try {
797                 
798                 MPI_Status status; 
799                 int pid;
800                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
801                 
802                 for(int i=0;i<num;i++){
803                         
804                         if (m->control_pressed) { return 0; }
805                         
806                         //read next sequence
807                         int length = MPIPos[start+i+1] - MPIPos[start+i];
808         
809                         char* buf4 = new char[length];
810                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
811                         
812                         string tempBuf = buf4;
813                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
814                         istringstream iss (tempBuf,istringstream::in);
815                         delete buf4;
816
817                         Sequence seq(iss);  
818
819                         if (seq.getAligned().length() != alignmentLength) {  cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1);  }
820                         
821                         if(trump != '*'){       F.doTrump(seq); }
822                         if(isTrue(vertical) || soft != 0){      F.getFreqs(seq);        }
823                         cout.flush();
824                                                 
825                         //report progress
826                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
827                 }
828                 
829                 //report progress
830                 if((num) % 100 != 0){   cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n");       }
831                 
832                 return 0;
833         }
834         catch(exception& e) {
835                 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
836                 exit(1);
837         }
838 }
839 #endif
840 /**************************************************************************************************/
841
842 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
843         try {
844 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
845                 int process = 0;
846                 int exitCommand = 1;
847                 processIDS.clear();
848                 
849                 //loop through and create all the processes you want
850                 while (process != processors) {
851                         int pid = vfork();
852                         
853                         if (pid > 0) {
854                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
855                                 process++;
856                         }else if (pid == 0){
857                                 driverCreateFilter(F, filename, lines[process]);
858                                 exit(0);
859                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
860                 }
861                 
862                 //force parent to wait until all the processes are done
863                 for (int i=0;i<processors;i++) { 
864                         int temp = processIDS[i];
865                         wait(&temp);
866                 }
867                 
868                 return exitCommand;
869 #endif          
870         }
871         catch(exception& e) {
872                 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
873                 exit(1);
874         }
875 }
876 /**************************************************************************************************/
877
878 int FilterSeqsCommand::setLines(string filename) {
879         try {
880                 
881                 vector<long int> positions;
882                 bufferSizes.clear();
883                 
884                 ifstream inFASTA;
885                 openInputFile(filename, inFASTA);
886                         
887                 string input;
888                 while(!inFASTA.eof()){  
889                         input = getline(inFASTA);
890
891                         if (input.length() != 0) {
892                                 if(input[0] == '>'){ long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);     }
893                         }
894                 }
895                 inFASTA.close();
896                 
897                 int numFastaSeqs = positions.size();
898         
899                 FILE * pFile;
900                 long size;
901                 
902                 //get num bytes in file
903                 pFile = fopen (filename.c_str(),"rb");
904                 if (pFile==NULL) perror ("Error opening file");
905                 else{
906                         fseek (pFile, 0, SEEK_END);
907                         size=ftell (pFile);
908                         fclose (pFile);
909                 }
910         
911                 numSeqs += numFastaSeqs;
912                 
913                 int numSeqsPerProcessor = numFastaSeqs / processors;
914                 
915                 for (int i = 0; i < processors; i++) {
916
917                         long int startPos = positions[ i * numSeqsPerProcessor ];
918                         if(i == processors - 1){
919                                 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
920                                 bufferSizes.push_back(size - startPos);
921                         }else{  
922                                 long int myEnd = positions[ (i+1) * numSeqsPerProcessor ];
923                                 bufferSizes.push_back(myEnd-startPos);
924                         }
925                         lines.push_back(new linePair(startPos, numSeqsPerProcessor));
926                 }
927                 
928                 return 0;
929         }
930         catch(exception& e) {
931                 m->errorOut(e, "FilterSeqsCommand", "setLines");
932                 exit(1);
933         }
934 }
935 /**************************************************************************************/