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