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