]> git.donarmstrong.com Git - mothur.git/blob - filterseqscommand.cpp
93e18e82f374c9fa5bad6a90f6fc2ff56e20246a
[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 vector<string> FilterSeqsCommand::setParameters(){      
16         try {
17                 CommandParameter pfasta("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pfasta);
18                 CommandParameter phard("hard", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(phard);
19                 CommandParameter ptrump("trump", "String", "", "*", "", "", "",false,false); parameters.push_back(ptrump);
20                 CommandParameter psoft("soft", "Number", "", "0", "", "", "",false,false); parameters.push_back(psoft);
21                 CommandParameter pvertical("vertical", "Boolean", "", "T", "", "", "",false,false); parameters.push_back(pvertical);
22                 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
23                 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
24                 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
25                 
26                 vector<string> myArray;
27                 for (int i = 0; i < parameters.size(); i++) {   myArray.push_back(parameters[i].name);          }
28                 return myArray;
29         }
30         catch(exception& e) {
31                 m->errorOut(e, "FilterSeqsCommand", "setParameters");
32                 exit(1);
33         }
34 }
35 //**********************************************************************************************************************
36 string FilterSeqsCommand::getHelpString(){      
37         try {
38                 string helpString = "";
39                 helpString += "The filter.seqs command reads a file containing sequences and creates a .filter and .filter.fasta file.\n";
40                 helpString += "The filter.seqs command parameters are fasta, trump, soft, hard, processors and vertical. \n";
41                 helpString += "The fasta parameter is required, unless you have a valid current fasta file. You may enter several fasta files to build the filter from and filter, by separating their names with -'s.\n";
42                 helpString += "For example: fasta=abrecovery.fasta-amazon.fasta \n";
43                 helpString += "The trump option will remove a column if the trump character is found at that position in any sequence of the alignment. Default=*, meaning no trump. \n";
44                 helpString += "A soft mask removes any column where the dominant base (i.e. A, T, G, C, or U) does not occur in at least a designated percentage of sequences. Default=0.\n";
45                 helpString += "The hard parameter allows you to enter a file containing the filter you want to use.\n";
46                 helpString += "The vertical parameter removes columns where all sequences contain a gap character. The default is T.\n";
47                 helpString += "The processors parameter allows you to specify the number of processors to use. The default is 1.\n";
48                 helpString += "The filter.seqs command should be in the following format: \n";
49                 helpString += "filter.seqs(fasta=yourFastaFile, trump=yourTrump) \n";
50                 helpString += "Example filter.seqs(fasta=abrecovery.fasta, trump=.).\n";
51                 helpString += "Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n";
52                 return helpString;
53         }
54         catch(exception& e) {
55                 m->errorOut(e, "FilterSeqsCommand", "getHelpString");
56                 exit(1);
57         }
58 }
59 //**********************************************************************************************************************
60 FilterSeqsCommand::FilterSeqsCommand(){ 
61         try {
62                 abort = true; calledHelp = true; 
63                 setParameters();
64                 vector<string> tempOutNames;
65                 outputTypes["fasta"] = tempOutNames;
66                 outputTypes["filter"] = tempOutNames;
67         }
68         catch(exception& e) {
69                 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
70                 exit(1);
71         }
72 }
73 /**************************************************************************************/
74 FilterSeqsCommand::FilterSeqsCommand(string option)  {
75         try {
76                 abort = false; calledHelp = false;   
77                 filterFileName = "";
78                 
79                 //allow user to run help
80                 if(option == "help") { help(); abort = true; calledHelp = true; }
81                 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
82                 
83                 else {
84                         vector<string> myArray = setParameters();
85                         
86                         OptionParser parser(option);
87                         map<string,string> parameters = parser.getParameters();
88                         
89                         ValidParameters validParameter("filter.seqs");
90                         map<string,string>::iterator it;
91                         
92                         //check to make sure all parameters are valid for command
93                         for (it = parameters.begin(); it != parameters.end(); it++) { 
94                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
95                         }
96                         
97                         //initialize outputTypes
98                         vector<string> tempOutNames;
99                         outputTypes["fasta"] = tempOutNames;
100                         outputTypes["filter"] = tempOutNames;
101                 
102                         //if the user changes the input directory command factory will send this info to us in the output parameter 
103                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
104                         if (inputDir == "not found"){   inputDir = "";          }
105                         else {
106                                 string path;
107                                 it = parameters.find("fasta");
108                                 //user has given a template file
109                                 if(it != parameters.end()){ 
110                                         path = m->hasPath(it->second);
111                                         //if the user has not given a path then, add inputdir. else leave path alone.
112                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
113                                 }
114                                 
115                                 it = parameters.find("hard");
116                                 //user has given a template file
117                                 if(it != parameters.end()){ 
118                                         path = m->hasPath(it->second);
119                                         //if the user has not given a path then, add inputdir. else leave path alone.
120                                         if (path == "") {       parameters["hard"] = inputDir + it->second;             }
121                                 }
122                         }
123                         
124                         //check for required parameters
125                         fasta = validParameter.validFile(parameters, "fasta", false);
126                         if (fasta == "not found") {                             
127                                 fasta = m->getFastaFile(); 
128                                 if (fasta != "") { fastafileNames.push_back(fasta);  m->mothurOut("Using " + fasta + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
129                                 else {  m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
130                         }
131                         else { 
132                                 m->splitAtDash(fasta, fastafileNames);
133                                 
134                                 //go through files and make sure they are good, if not, then disregard them
135                                 for (int i = 0; i < fastafileNames.size(); i++) {
136                                         
137                                         bool ignore = false;
138                                         if (fastafileNames[i] == "current") { 
139                                                 fastafileNames[i] = m->getFastaFile(); 
140                                                 if (fastafileNames[i] != "") {  m->mothurOut("Using " + fastafileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
141                                                 else {  
142                                                         m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true; 
143                                                         //erase from file list
144                                                         fastafileNames.erase(fastafileNames.begin()+i);
145                                                         i--;
146                                                 }
147                                         }
148                                         
149                                         if (!ignore) {
150                                                 if (inputDir != "") {
151                                                         string path = m->hasPath(fastafileNames[i]);
152                                                         //if the user has not given a path then, add inputdir. else leave path alone.
153                                                         if (path == "") {       fastafileNames[i] = inputDir + fastafileNames[i];               }
154                                                 }
155
156                                                 ifstream in;
157                                                 int ableToOpen = m->openInputFile(fastafileNames[i], in, "noerror");
158                                         
159                                                 //if you can't open it, try default location
160                                                 if (ableToOpen == 1) {
161                                                         if (m->getDefaultPath() != "") { //default path is set
162                                                                 string tryPath = m->getDefaultPath() + m->getSimpleName(fastafileNames[i]);
163                                                                 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
164                                                                 ifstream in2;
165                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
166                                                                 in2.close();
167                                                                 fastafileNames[i] = tryPath;
168                                                         }
169                                                 }
170                                                 
171                                                 //if you can't open it, try default location
172                                                 if (ableToOpen == 1) {
173                                                         if (m->getOutputDir() != "") { //default path is set
174                                                                 string tryPath = m->getOutputDir() + m->getSimpleName(fastafileNames[i]);
175                                                                 m->mothurOut("Unable to open " + fastafileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
176                                                                 ifstream in2;
177                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
178                                                                 in2.close();
179                                                                 fastafileNames[i] = tryPath;
180                                                         }
181                                                 }
182                                                 
183                                                 in.close();
184                                                 
185                                                 if (ableToOpen == 1) { 
186                                                         m->mothurOut("Unable to open " + fastafileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
187                                                         //erase from file list
188                                                         fastafileNames.erase(fastafileNames.begin()+i);
189                                                         i--;
190                                                 }else{  
191                                                         string simpleName = m->getSimpleName(fastafileNames[i]);
192                                                         filterFileName += simpleName.substr(0, simpleName.find_first_of('.'));
193                                                         m->setFastaFile(fastafileNames[i]);
194                                                 }
195                                                 in.close();
196                                         }
197                                 }
198                                 
199                                 //make sure there is at least one valid file left
200                                 if (fastafileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
201                         }
202                         
203                         if (!abort) {
204                                 //if the user changes the output directory command factory will send this info to us in the output parameter 
205                                 outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
206                                         outputDir = ""; 
207                                         outputDir += m->hasPath(fastafileNames[0]); //if user entered a file with a path then preserve it       
208                                 }
209                         }
210                         //check for optional parameter and set defaults
211                         // ...at some point should added some additional type checking...
212                         
213                         string temp;
214                         hard = validParameter.validFile(parameters, "hard", true);                              if (hard == "not found") { hard = ""; }
215                         else if (hard == "not open") { hard = ""; abort = true; }       
216
217                         temp = validParameter.validFile(parameters, "trump", false);                    if (temp == "not found") { temp = "*"; }
218                         trump = temp[0];
219                         
220                         temp = validParameter.validFile(parameters, "soft", false);                             if (temp == "not found") { soft = 0; }
221                         else {  soft = (float)atoi(temp.c_str()) / 100.0;  }
222                         
223                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = m->getProcessors();      }
224                         m->setProcessors(temp);
225                         m->mothurConvert(temp, processors); 
226                         
227                         vertical = validParameter.validFile(parameters, "vertical", false);             
228                         if (vertical == "not found") { 
229                                 if ((hard == "") && (trump == '*') && (soft == 0)) { vertical = "T"; } //you have not given a hard file or set the trump char.
230                                 else { vertical = "F";  }
231                         }
232                         
233                         numSeqs = 0;
234                 }
235                 
236         }
237         catch(exception& e) {
238                 m->errorOut(e, "FilterSeqsCommand", "FilterSeqsCommand");
239                 exit(1);
240         }
241 }
242 /**************************************************************************************/
243
244 int FilterSeqsCommand::execute() {      
245         try {
246         
247                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
248                 
249                 ifstream inFASTA;
250                 m->openInputFile(fastafileNames[0], inFASTA);
251                 
252                 Sequence testSeq(inFASTA);
253                 alignmentLength = testSeq.getAlignLength();
254                 inFASTA.close();
255                 
256                 ////////////create filter/////////////////
257                 m->mothurOut("Creating Filter... "); m->mothurOutEndLine();
258                 
259                 filter = createFilter();
260                 
261                 m->mothurOutEndLine();  m->mothurOutEndLine();
262                 
263                 if (m->control_pressed) { outputTypes.clear(); return 0; }
264                 
265                 #ifdef USE_MPI
266                         int pid;
267                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
268                                         
269                         if (pid == 0) { //only one process should output the filter
270                 #endif
271                 
272                 ofstream outFilter;
273                 
274                 //prevent giantic file name
275                 string filterFile;
276                 if (fastafileNames.size() > 3) { filterFile = outputDir + "merge.filter"; }
277                 else {  filterFile = outputDir + filterFileName + ".filter";  }
278                 
279                 m->openOutputFile(filterFile, outFilter);
280                 outFilter << filter << endl;
281                 outFilter.close();
282                 outputNames.push_back(filterFile); outputTypes["filter"].push_back(filterFile);
283                 
284                 #ifdef USE_MPI
285                         }
286                 #endif
287                 
288                 ////////////run filter/////////////////
289                 
290                 m->mothurOut("Running Filter... "); m->mothurOutEndLine();
291                 
292                 filterSequences();
293                 
294                 m->mothurOutEndLine();  m->mothurOutEndLine();
295                                         
296                 int filteredLength = 0;
297                 for(int i=0;i<alignmentLength;i++){
298                         if(filter[i] == '1'){   filteredLength++;       }
299                 }
300                 
301                 if (m->control_pressed) {  outputTypes.clear(); for(int i = 0; i < outputNames.size(); i++) { m->mothurRemove(outputNames[i]); }  return 0; }
302
303                 
304                 m->mothurOutEndLine();
305                 m->mothurOut("Length of filtered alignment: " + toString(filteredLength)); m->mothurOutEndLine();
306                 m->mothurOut("Number of columns removed: " + toString((alignmentLength-filteredLength))); m->mothurOutEndLine();
307                 m->mothurOut("Length of the original alignment: " + toString(alignmentLength)); m->mothurOutEndLine();
308                 m->mothurOut("Number of sequences used to construct filter: " + toString(numSeqs)); m->mothurOutEndLine();
309                 
310                 //set fasta file as new current fastafile
311                 string current = "";
312                 itTypes = outputTypes.find("fasta");
313                 if (itTypes != outputTypes.end()) {
314                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
315                 }
316                 
317                 m->mothurOutEndLine();
318                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
319                 for(int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();      }
320                 m->mothurOutEndLine();
321                 
322                 return 0;
323                 
324         }
325         catch(exception& e) {
326                 m->errorOut(e, "FilterSeqsCommand", "execute");
327                 exit(1);
328         }
329 }
330 /**************************************************************************************/
331 int FilterSeqsCommand::filterSequences() {      
332         try {
333                 
334                 numSeqs = 0;
335                 
336                 for (int s = 0; s < fastafileNames.size(); s++) {
337                         
338                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
339                                 
340                                 string filteredFasta = outputDir + m->getRootName(m->getSimpleName(fastafileNames[s])) + "filter.fasta";
341 #ifdef USE_MPI  
342                                 int pid, numSeqsPerProcessor, num; 
343                                 int tag = 2001;
344                                 vector<unsigned long long>MPIPos;
345                                                 
346                                 MPI_Status status; 
347                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
348                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
349                                 
350                                 MPI_File outMPI;
351                                 MPI_File inMPI;
352                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
353                                 int inMode=MPI_MODE_RDONLY; 
354                                 
355                                 char outFilename[1024];
356                                 strcpy(outFilename, filteredFasta.c_str());
357                         
358                                 char inFileName[1024];
359                                 strcpy(inFileName, fastafileNames[s].c_str());
360                                 
361                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
362                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
363
364                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
365
366                                 if (pid == 0) { //you are the root process 
367                                         
368                                         MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
369                                         numSeqs += num;
370                                         
371                                         //send file positions to all processes
372                                         for(int i = 1; i < processors; i++) { 
373                                                 MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
374                                                 MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
375                                         }
376                                         
377                                         //figure out how many sequences you have to do
378                                         numSeqsPerProcessor = num / processors;
379                                         int startIndex =  pid * numSeqsPerProcessor;
380                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
381                                         
382                                 
383                                         //do your part
384                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
385                                         
386                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
387                                         
388                                         //wait on chidren
389                                         for(int i = 1; i < processors; i++) { 
390                                                 char buf[5];
391                                                 MPI_Recv(buf, 5, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status); 
392                                         }
393                                         
394                                 }else { //you are a child process
395                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
396                                         MPIPos.resize(num+1);
397                                         numSeqs += num;
398                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
399                                         
400                                         //figure out how many sequences you have to align
401                                         numSeqsPerProcessor = num / processors;
402                                         int startIndex =  pid * numSeqsPerProcessor;
403                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
404                                         
405                                         
406                                         //align your part
407                                         driverMPIRun(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);           
408                                         
409                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  return 0;  }
410                                         
411                                         char buf[5];
412                                         strcpy(buf, "done"); 
413                                         
414                                         //tell parent you are done.
415                                         MPI_Send(buf, 5, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
416                                 }
417                                 
418                                 MPI_File_close(&outMPI);
419                                 MPI_File_close(&inMPI);
420                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
421                                 
422 #else
423             
424             vector<unsigned long long> positions;
425             if (savedPositions.size() != 0) { positions = savedPositions[s]; }
426             else {
427 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
428                                 positions = m->divideFile(fastafileNames[s], processors);
429 #else
430                 if(processors != 1){
431                     int numFastaSeqs = 0;
432                     positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs); 
433                     if (positions.size() < processors) { processors = positions.size(); }
434                 }
435 #endif
436             }
437                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
438                         //vector<unsigned long long> positions = m->divideFile(fastafileNames[s], processors);
439                         
440                         for (int i = 0; i < (positions.size()-1); i++) {
441                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
442                         }       
443                         
444                                 if(processors == 1){
445                                         int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
446                                         numSeqs += numFastaSeqs;
447                                 }else{
448                                         int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta); 
449                                         numSeqs += numFastaSeqs;
450                                 }
451                                 
452                                 if (m->control_pressed) {  return 1; }
453                 #else
454             if(processors == 1){
455                 lines.push_back(new linePair(0, 1000));
456                                 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
457                                 numSeqs += numFastaSeqs;
458             }else {
459                 int numFastaSeqs = positions.size()-1;
460                 //positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs); 
461                 
462                 //figure out how many sequences you have to process
463                 int numSeqsPerProcessor = numFastaSeqs / processors;
464                 for (int i = 0; i < processors; i++) {
465                     int startIndex =  i * numSeqsPerProcessor;
466                     if(i == (processors - 1)){  numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;   }
467                     lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
468                 }
469                 
470                 numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta); 
471                 numSeqs += numFastaSeqs;
472             }
473
474                                 if (m->control_pressed) {  return 1; }
475                 #endif
476 #endif
477                         outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
478                 }
479
480                 return 0;
481         }
482         catch(exception& e) {
483                 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
484                 exit(1);
485         }
486 }
487 #ifdef USE_MPI
488 /**************************************************************************************/
489 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long long>& MPIPos) {        
490         try {
491                 string outputString = "";
492                 int count = 0;
493                 MPI_Status status; 
494                 
495                 for(int i=0;i<num;i++){
496                 
497                         if (m->control_pressed) { return 0; }
498                 
499                         //read next sequence
500                         int length = MPIPos[start+i+1] - MPIPos[start+i];
501                         char* buf4 = new char[length];
502                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
503                         
504                         string tempBuf = buf4;
505                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
506                         istringstream iss (tempBuf,istringstream::in);
507                         delete buf4;
508         
509                         Sequence seq(iss);  m->gobble(iss);
510                         
511                         if (seq.getName() != "") {
512                                 string align = seq.getAligned();
513                                 string filterSeq = "";
514                                         
515                                 for(int j=0;j<alignmentLength;j++){
516                                         if(filter[j] == '1'){
517                                                 filterSeq += align[j];
518                                         }
519                                 }
520                                 
521                                 count++;
522                                 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
523                                 
524                                 if(count % 10 == 0){ //output to file 
525                                         //send results to parent
526                                         int length = outputString.length();
527                                         char* buf = new char[length];
528                                         memcpy(buf, outputString.c_str(), length);
529                                 
530                                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
531                                         outputString = "";
532                                         delete buf;
533                                 }
534
535                         }
536                         
537                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
538                 }
539                 
540                 if(outputString != ""){ //output to file 
541                         //send results to parent
542                         int length = outputString.length();
543                         char* buf = new char[length];
544                         memcpy(buf, outputString.c_str(), length);
545                         
546                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
547                         outputString = "";
548                         delete buf;
549                 }
550                 
551                 if((num) % 100 != 0){   cout << (num) << endl;   m->mothurOutJustToLog(toString(num) + "\n");   }
552                         
553                 return 0;
554         }
555         catch(exception& e) {
556                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
557                 exit(1);
558         }
559 }
560 #endif
561 /**************************************************************************************/
562 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {      
563         try {
564                 ofstream out;
565                 m->openOutputFile(outputFilename, out);
566                 
567                 ifstream in;
568                 m->openInputFile(inputFilename, in);
569                                 
570                 in.seekg(filePos->start);
571
572                 bool done = false;
573                 int count = 0;
574         
575                 while (!done) {
576                                 
577                                 if (m->control_pressed) { in.close(); out.close(); return 0; }
578                                 
579                                 Sequence seq(in); m->gobble(in);
580                                 if (seq.getName() != "") {
581                                         string align = seq.getAligned();
582                                         string filterSeq = "";
583                                         
584                                         for(int j=0;j<alignmentLength;j++){
585                                                 if(filter[j] == '1'){
586                                                         filterSeq += align[j];
587                                                 }
588                                         }
589                                         
590                                         out << '>' << seq.getName() << endl << filterSeq << endl;
591                                 count++;
592                         }
593                         
594                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
595                                 unsigned long long pos = in.tellg();
596                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
597                         #else
598                                 if (in.eof()) { break; }
599                         #endif
600                         
601                         //report progress
602                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
603                 }
604                 //report progress
605                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
606                 
607                 
608                 out.close();
609                 in.close();
610                 
611                 return count;
612         }
613         catch(exception& e) {
614                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
615                 exit(1);
616         }
617 }
618 /**************************************************************************************************/
619
620 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename, string filteredFastaName) {
621         try {
622         
623         int process = 1;
624                 int num = 0;
625                 processIDS.clear();
626         
627 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
628                 
629                 
630                 //loop through and create all the processes you want
631                 while (process != processors) {
632                         int pid = fork();
633                         
634                         if (pid > 0) {
635                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
636                                 process++;
637                         }else if (pid == 0){
638                                 string filteredFasta = filename + toString(getpid()) + ".temp";
639                                 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
640                                 
641                                 //pass numSeqs to parent
642                                 ofstream out;
643                                 string tempFile = filename +  toString(getpid()) + ".num.temp";
644                                 m->openOutputFile(tempFile, out);
645                                 out << num << endl;
646                                 out.close();
647                                 
648                                 exit(0);
649                         }else { 
650                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
651                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
652                                 exit(0);
653                         }
654                 }
655                 
656         num = driverRunFilter(F, filteredFastaName, filename, lines[0]);
657         
658                 //force parent to wait until all the processes are done
659                 for (int i=0;i<processIDS.size();i++) { 
660                         int temp = processIDS[i];
661                         wait(&temp);
662                 }       
663                                         
664                 for (int i = 0; i < processIDS.size(); i++) {
665                         ifstream in;
666                         string tempFile =  filename + toString(processIDS[i]) + ".num.temp";
667                         m->openInputFile(tempFile, in);
668                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
669                         in.close(); m->mothurRemove(tempFile);
670             
671             m->appendFiles((filename + toString(processIDS[i]) + ".temp"), filteredFastaName);
672             m->mothurRemove((filename + toString(processIDS[i]) + ".temp"));
673                 }
674                
675 #else
676         
677         //////////////////////////////////////////////////////////////////////////////////////////////////////
678                 //Windows version shared memory, so be careful when passing variables through the filterData struct. 
679                 //Above fork() will clone, so memory is separate, but that's not the case with windows, 
680                 //Taking advantage of shared memory to allow both threads to add info to F.
681                 //////////////////////////////////////////////////////////////////////////////////////////////////////
682                 
683                 vector<filterRunData*> pDataArray; 
684                 DWORD   dwThreadIdArray[processors-1];
685                 HANDLE  hThreadArray[processors-1]; 
686                 
687                 //Create processor worker threads.
688                 for( int i=0; i<processors-1; i++){
689                         
690             string extension = "";
691                         if (i != 0) { extension = toString(i) + ".temp"; }
692             
693                         filterRunData* tempFilter = new filterRunData(filter, filename, (filteredFastaName + extension), m, lines[i]->start, lines[i]->end, alignmentLength, i);
694                         pDataArray.push_back(tempFilter);
695                         processIDS.push_back(i);
696             
697                         hThreadArray[i] = CreateThread(NULL, 0, MyRunFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);   
698                 }
699         
700         num = driverRunFilter(F, (filteredFastaName + toString(processors-1) + ".temp"), filename, lines[processors-1]);
701         
702                 //Wait until all threads have terminated.
703                 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
704                 
705                 //Close all thread handles and free memory allocations.
706                 for(int i=0; i < pDataArray.size(); i++){
707                         num += pDataArray[i]->count;
708             CloseHandle(hThreadArray[i]);
709                         delete pDataArray[i];
710                 }
711         
712         for (int i = 1; i < processors; i++) {
713             m->appendFiles((filteredFastaName + toString(i) + ".temp"), filteredFastaName);
714             m->mothurRemove((filteredFastaName + toString(i) + ".temp"));
715                 }
716 #endif  
717         
718         return num;
719         
720         }
721         catch(exception& e) {
722                 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
723                 exit(1);
724         }
725 }
726 /**************************************************************************************/
727 string FilterSeqsCommand::createFilter() {      
728         try {
729                 string filterString = "";                       
730                 Filters F;
731                 
732                 if (soft != 0)                  {  F.setSoft(soft);             }
733                 if (trump != '*')               {  F.setTrump(trump);   }
734                 
735                 F.setLength(alignmentLength);
736                 
737                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
738                         F.initialize();
739                 }
740                 
741                 if(hard.compare("") != 0)       {       F.doHard(hard);         }
742                 else                                            {       F.setFilter(string(alignmentLength, '1'));      }
743                 
744                 numSeqs = 0;
745                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
746                         for (int s = 0; s < fastafileNames.size(); s++) {
747                         
748                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
749                         
750 #ifdef USE_MPI  
751                                 int pid, numSeqsPerProcessor, num; 
752                                 int tag = 2001;
753                                 vector<unsigned long long> MPIPos;
754                                 
755                                 MPI_Status status; 
756                                 MPI_File inMPI; 
757                                 MPI_Comm_size(MPI_COMM_WORLD, &processors);
758                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
759                                                         
760                                 //char* tempFileName = new char(fastafileNames[s].length());
761                                 //tempFileName = &(fastafileNames[s][0]);
762                                 
763                                 char tempFileName[1024];
764                                 strcpy(tempFileName, fastafileNames[s].c_str());
765                 
766                                 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
767                                 
768                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
769                                 
770                                 if (pid == 0) { //you are the root process
771                                                 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
772                                                 numSeqs += num;
773                                                 
774                                                 //send file positions to all processes
775                                                 for(int i = 1; i < processors; i++) { 
776                                                         MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
777                                                         MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
778                                                 }
779                                                                 
780                                                 //figure out how many sequences you have to do
781                                                 numSeqsPerProcessor = num / processors;
782                                                 int startIndex =  pid * numSeqsPerProcessor;
783                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
784                                                 
785                                 
786                                                 //do your part
787                                                 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
788                                                 
789                                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
790                                                                                                 
791                                 }else { //i am the child process
792                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
793                                         MPIPos.resize(num+1);
794                                         numSeqs += num;
795                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
796                                         
797                                         //figure out how many sequences you have to align
798                                         numSeqsPerProcessor = num / processors;
799                                         int startIndex =  pid * numSeqsPerProcessor;
800                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
801                                         
802                                         
803                                         //do your part
804                                         MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI,  MPIPos);
805                                         
806                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
807                                 }
808                                 
809                                 MPI_File_close(&inMPI);
810                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
811                                 
812 #else
813                                 
814                 vector<unsigned long long> positions;
815                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
816                                 positions = m->divideFile(fastafileNames[s], processors);
817                                 for (int i = 0; i < (positions.size()-1); i++) {
818                                         lines.push_back(new linePair(positions[i], positions[(i+1)]));
819                                 }       
820                                 
821                                 if(processors == 1){
822                                         int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
823                                         numSeqs += numFastaSeqs;
824                                 }else{
825                                         int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]); 
826                                         numSeqs += numFastaSeqs;
827                                 }
828                 #else
829                 if(processors == 1){
830                     lines.push_back(new linePair(0, 1000));
831                     int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
832                     numSeqs += numFastaSeqs;
833                                 }else {
834                     int numFastaSeqs = 0;
835                     positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs); 
836                     if (positions.size() < processors) { processors = positions.size(); }
837                     
838                     //figure out how many sequences you have to process
839                     int numSeqsPerProcessor = numFastaSeqs / processors;
840                     for (int i = 0; i < processors; i++) {
841                         int startIndex =  i * numSeqsPerProcessor;
842                         if(i == (processors - 1)){      numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;   }
843                         lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
844                     }
845                     
846                     numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]); 
847                                         numSeqs += numFastaSeqs;
848                 }
849                 #endif
850                 //save the file positions so we can reuse them in the runFilter function
851                 savedPositions[s] = positions;
852                 
853                                 if (m->control_pressed) {  return filterString; }
854 #endif
855                         
856                         }
857                 }
858
859
860 #ifdef USE_MPI  
861                 int pid;
862                 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
863                 MPI_Status status;
864                 
865                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
866                 
867                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
868                         
869                         if (pid == 0) { //only one process should output the filter
870                         
871                                 vector<int> temp; temp.resize(alignmentLength+1);
872                                                                 
873                                 //get the frequencies from the child processes
874                                 for(int i = 1; i < processors; i++) { 
875                                 
876                                         for (int j = 0; j < 5; j++) {
877                                         
878                                                 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status); 
879                                                 int receiveTag = temp[temp.size()-1];  //child process added a int to the end to indicate what letter count this is for
880                                                 
881                                                 if (receiveTag == Atag) { //you are recieveing the A frequencies
882                                                         for (int k = 0; k < alignmentLength; k++) {             F.a[k] += temp[k];      }
883                                                 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
884                                                         for (int k = 0; k < alignmentLength; k++) {             F.t[k] += temp[k];      }
885                                                 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
886                                                         for (int k = 0; k < alignmentLength; k++) {             F.c[k] += temp[k];      }
887                                                 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
888                                                         for (int k = 0; k < alignmentLength; k++) {             F.g[k] += temp[k];      }
889                                                 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
890                                                         for (int k = 0; k < alignmentLength; k++) {             F.gap[k] += temp[k];    }
891                                                 }
892                                         }
893                                 } 
894                         }else{
895                         
896                                 //send my fequency counts
897                                 F.a.push_back(Atag);
898                                 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
899                                 F.t.push_back(Ttag);
900                                 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
901                                 F.c.push_back(Ctag);
902                                 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
903                                 F.g.push_back(Gtag);
904                                 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
905                                 F.gap.push_back(Gaptag);
906                                 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
907                         }
908                         
909                 }
910                 
911                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
912                 
913                 if (pid == 0) { //only one process should output the filter
914 #endif
915
916                 F.setNumSeqs(numSeqs);
917                 if(m->isTrue(vertical) == 1)    {       F.doVertical(); }
918                 if(soft != 0)                           {       F.doSoft();             }
919                 filterString = F.getFilter();
920                 
921 #ifdef USE_MPI
922                 //send filter string to kids
923                 //for(int i = 1; i < processors; i++) { 
924                 //      MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
925                 //}
926                 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
927         }else{
928                 //recieve filterString
929                 char* tempBuf = new char[alignmentLength];
930                 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
931                 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
932                 
933                 filterString = tempBuf;
934                 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength);  }
935                 delete tempBuf; 
936         }
937         
938         MPI_Barrier(MPI_COMM_WORLD);
939 #endif
940             
941                 return filterString;
942         }
943         catch(exception& e) {
944                 m->errorOut(e, "FilterSeqsCommand", "createFilter");
945                 exit(1);
946         }
947 }
948 /**************************************************************************************/
949 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {     
950         try {
951                 
952                 ifstream in;
953                 m->openInputFile(filename, in);
954                                 
955                 in.seekg(filePos->start);
956
957                 bool done = false;
958                 int count = 0;
959         
960                 while (!done) {
961                                 
962                         if (m->control_pressed) { in.close(); return 1; }
963                                         
964                         Sequence seq(in); m->gobble(in);
965                         if (seq.getName() != "") {
966                                         if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true;  }
967                                         
968                                         if(trump != '*')                        {       F.doTrump(seq);         }
969                                         if(m->isTrue(vertical) || soft != 0)    {       F.getFreqs(seq);        }
970                                         cout.flush();
971                                         count++;
972                         }
973                         
974                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
975                                 unsigned long long pos = in.tellg();
976                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
977                         #else
978                                 if (in.eof()) { break; }
979                         #endif
980                         
981                         //report progress
982                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
983                 }
984                 //report progress
985                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
986                 in.close();
987                 
988                 return count;
989         }
990         catch(exception& e) {
991                 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
992                 exit(1);
993         }
994 }
995 #ifdef USE_MPI
996 /**************************************************************************************/
997 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long long>& MPIPos) {   
998         try {
999                 
1000                 MPI_Status status; 
1001                 int pid;
1002                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
1003                 
1004                 for(int i=0;i<num;i++){
1005                         
1006                         if (m->control_pressed) { return 0; }
1007                         
1008                         //read next sequence
1009                         int length = MPIPos[start+i+1] - MPIPos[start+i];
1010         
1011                         char* buf4 = new char[length];
1012                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
1013                         
1014                         string tempBuf = buf4;
1015                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
1016                         istringstream iss (tempBuf,istringstream::in);
1017                         delete buf4;
1018
1019                         Sequence seq(iss);  
1020
1021                         if (seq.getAligned().length() != alignmentLength) {  cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1);  }
1022                         
1023                         if(trump != '*'){       F.doTrump(seq); }
1024                         if(m->isTrue(vertical) || soft != 0){   F.getFreqs(seq);        }
1025                         cout.flush();
1026                                                 
1027                         //report progress
1028                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
1029                 }
1030                 
1031                 //report progress
1032                 if((num) % 100 != 0){   cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n");       }
1033                 
1034                 return 0;
1035         }
1036         catch(exception& e) {
1037                 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
1038                 exit(1);
1039         }
1040 }
1041 #endif
1042 /**************************************************************************************************/
1043
1044 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
1045         try {
1046         int process = 1;
1047                 int num = 0;
1048                 processIDS.clear();
1049
1050 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
1051                                 
1052                 //loop through and create all the processes you want
1053                 while (process != processors) {
1054                         int pid = fork();
1055                         
1056                         if (pid > 0) {
1057                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
1058                                 process++;
1059                         }else if (pid == 0){
1060                                 //reset child's filter counts to 0;
1061                                 F.a.clear(); F.a.resize(alignmentLength, 0);
1062                                 F.t.clear(); F.t.resize(alignmentLength, 0);
1063                                 F.g.clear(); F.g.resize(alignmentLength, 0);
1064                                 F.c.clear(); F.c.resize(alignmentLength, 0);
1065                                 F.gap.clear(); F.gap.resize(alignmentLength, 0);
1066                                 
1067                                 num = driverCreateFilter(F, filename, lines[process]);
1068                                 
1069                                 //write out filter counts to file
1070                                 filename += toString(getpid()) + "filterValues.temp";
1071                                 ofstream out;
1072                                 m->openOutputFile(filename, out);
1073                                 
1074                                 out << num << endl;
1075                                 out << F.getFilter() << endl;
1076                                 for (int k = 0; k < alignmentLength; k++) {             out << F.a[k] << '\t'; }  out << endl;
1077                                 for (int k = 0; k < alignmentLength; k++) {             out << F.t[k] << '\t'; }  out << endl;
1078                                 for (int k = 0; k < alignmentLength; k++) {             out << F.g[k] << '\t'; }  out << endl;
1079                                 for (int k = 0; k < alignmentLength; k++) {             out << F.c[k] << '\t'; }  out << endl;
1080                                 for (int k = 0; k < alignmentLength; k++) {             out << F.gap[k] << '\t'; }  out << endl;
1081
1082                                 //cout << F.getFilter() << endl;
1083                                 out.close();
1084                                 
1085                                 exit(0);
1086                         }else { 
1087                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
1088                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
1089                                 exit(0);
1090                         }
1091                 }
1092                 
1093                 //parent do your part
1094                 num = driverCreateFilter(F, filename, lines[0]);
1095                 
1096                 //force parent to wait until all the processes are done
1097                 for (int i=0;i<(processors-1);i++) { 
1098                         int temp = processIDS[i];
1099                         wait(&temp);
1100                 }
1101                 
1102                 //parent reads in and combines Filter info
1103                 for (int i = 0; i < processIDS.size(); i++) {
1104                         string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
1105                         ifstream in;
1106                         m->openInputFile(tempFilename, in);
1107                         
1108                         int temp, tempNum;
1109                         string tempFilterString;
1110
1111                         in >> tempNum; m->gobble(in); num += tempNum;
1112
1113                         in >> tempFilterString;
1114                         F.mergeFilter(tempFilterString);
1115
1116                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.a[k] += temp; }           m->gobble(in);
1117                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.t[k] += temp; }           m->gobble(in);
1118                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.g[k] += temp; }           m->gobble(in);
1119                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.c[k] += temp; }           m->gobble(in);
1120                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.gap[k] += temp; } m->gobble(in);
1121                                 
1122                         in.close();
1123                         m->mothurRemove(tempFilename);
1124                 }
1125                 
1126                 
1127 #else
1128         
1129         //////////////////////////////////////////////////////////////////////////////////////////////////////
1130                 //Windows version shared memory, so be careful when passing variables through the filterData struct. 
1131                 //Above fork() will clone, so memory is separate, but that's not the case with windows, 
1132                 //Taking advantage of shared memory to allow both threads to add info to F.
1133                 //////////////////////////////////////////////////////////////////////////////////////////////////////
1134                 
1135                 vector<filterData*> pDataArray; 
1136                 DWORD   dwThreadIdArray[processors];
1137                 HANDLE  hThreadArray[processors]; 
1138                 
1139                 //Create processor worker threads.
1140                 for( int i=0; i<processors; i++ ){
1141                         
1142                         filterData* tempFilter = new filterData(filename, m, lines[i]->start, lines[i]->end, alignmentLength, trump, vertical, soft, hard, i);
1143                         pDataArray.push_back(tempFilter);
1144                         processIDS.push_back(i);
1145             
1146                         hThreadArray[i] = CreateThread(NULL, 0, MyCreateFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);   
1147                 }
1148         
1149                 //Wait until all threads have terminated.
1150                 WaitForMultipleObjects(processors, hThreadArray, TRUE, INFINITE);
1151                 
1152                 //Close all thread handles and free memory allocations.
1153                 for(int i=0; i < pDataArray.size(); i++){
1154                         num += pDataArray[i]->count;
1155             F.mergeFilter(pDataArray[i]->F.getFilter());
1156             
1157                         for (int k = 0; k < alignmentLength; k++) {      F.a[k] += pDataArray[i]->F.a[k];       }
1158                         for (int k = 0; k < alignmentLength; k++) {      F.t[k] += pDataArray[i]->F.t[k];       }
1159                         for (int k = 0; k < alignmentLength; k++) {      F.g[k] += pDataArray[i]->F.g[k];       }
1160                         for (int k = 0; k < alignmentLength; k++) {      F.c[k] += pDataArray[i]->F.c[k];       }
1161                         for (int k = 0; k < alignmentLength; k++) {      F.gap[k] += pDataArray[i]->F.gap[k];   }
1162
1163                         CloseHandle(hThreadArray[i]);
1164                         delete pDataArray[i];
1165                 }
1166                 
1167 #endif  
1168         return num;
1169         
1170         }
1171         catch(exception& e) {
1172                 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1173                 exit(1);
1174         }
1175 }
1176 /**************************************************************************************/