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