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