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