]> git.donarmstrong.com Git - mothur.git/blob - filterseqscommand.cpp
paralellized filter.seqs for windows
[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                         vector<unsigned long long> positions = savedPositions[s];
424                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
425                         //vector<unsigned long long> positions = m->divideFile(fastafileNames[s], processors);
426                         
427                         for (int i = 0; i < (positions.size()-1); i++) {
428                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
429                         }       
430                         
431                                 if(processors == 1){
432                                         int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
433                                         numSeqs += numFastaSeqs;
434                                 }else{
435                                         int numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta); 
436                                         numSeqs += numFastaSeqs;
437                                 }
438                                 
439                                 if (m->control_pressed) {  return 1; }
440                 #else
441             if(processors == 1){
442                 lines.push_back(new linePair(0, 1000));
443                                 int numFastaSeqs = driverRunFilter(filter, filteredFasta, fastafileNames[s], lines[0]);
444                                 numSeqs += numFastaSeqs;
445             }else {
446                 int numFastaSeqs = positions.size()-1;
447                 //positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs); 
448                 
449                 //figure out how many sequences you have to process
450                 int numSeqsPerProcessor = numFastaSeqs / processors;
451                 for (int i = 0; i < processors; i++) {
452                     int startIndex =  i * numSeqsPerProcessor;
453                     if(i == (processors - 1)){  numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;   }
454                     lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
455                 }
456                 
457                 numFastaSeqs = createProcessesRunFilter(filter, fastafileNames[s], filteredFasta); 
458                 numSeqs += numFastaSeqs;
459             }
460
461                                 if (m->control_pressed) {  return 1; }
462                 #endif
463 #endif
464                         outputNames.push_back(filteredFasta); outputTypes["fasta"].push_back(filteredFasta);
465                 }
466
467                 return 0;
468         }
469         catch(exception& e) {
470                 m->errorOut(e, "FilterSeqsCommand", "filterSequences");
471                 exit(1);
472         }
473 }
474 #ifdef USE_MPI
475 /**************************************************************************************/
476 int FilterSeqsCommand::driverMPIRun(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long long>& MPIPos) {        
477         try {
478                 string outputString = "";
479                 int count = 0;
480                 MPI_Status status; 
481                 
482                 for(int i=0;i<num;i++){
483                 
484                         if (m->control_pressed) { return 0; }
485                 
486                         //read next sequence
487                         int length = MPIPos[start+i+1] - MPIPos[start+i];
488                         char* buf4 = new char[length];
489                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
490                         
491                         string tempBuf = buf4;
492                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
493                         istringstream iss (tempBuf,istringstream::in);
494                         delete buf4;
495         
496                         Sequence seq(iss);  m->gobble(iss);
497                         
498                         if (seq.getName() != "") {
499                                 string align = seq.getAligned();
500                                 string filterSeq = "";
501                                         
502                                 for(int j=0;j<alignmentLength;j++){
503                                         if(filter[j] == '1'){
504                                                 filterSeq += align[j];
505                                         }
506                                 }
507                                 
508                                 count++;
509                                 outputString += ">" + seq.getName() + "\n" + filterSeq + "\n";
510                                 
511                                 if(count % 10 == 0){ //output to file 
512                                         //send results to parent
513                                         int length = outputString.length();
514                                         char* buf = new char[length];
515                                         memcpy(buf, outputString.c_str(), length);
516                                 
517                                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
518                                         outputString = "";
519                                         delete buf;
520                                 }
521
522                         }
523                         
524                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
525                 }
526                 
527                 if(outputString != ""){ //output to file 
528                         //send results to parent
529                         int length = outputString.length();
530                         char* buf = new char[length];
531                         memcpy(buf, outputString.c_str(), length);
532                         
533                         MPI_File_write_shared(outMPI, buf, length, MPI_CHAR, &status);
534                         outputString = "";
535                         delete buf;
536                 }
537                 
538                 if((num) % 100 != 0){   cout << (num) << endl;   m->mothurOutJustToLog(toString(num) + "\n");   }
539                         
540                 return 0;
541         }
542         catch(exception& e) {
543                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
544                 exit(1);
545         }
546 }
547 #endif
548 /**************************************************************************************/
549 int FilterSeqsCommand::driverRunFilter(string F, string outputFilename, string inputFilename, linePair* filePos) {      
550         try {
551                 ofstream out;
552                 m->openOutputFile(outputFilename, out);
553                 
554                 ifstream in;
555                 m->openInputFile(inputFilename, in);
556                                 
557                 in.seekg(filePos->start);
558
559                 bool done = false;
560                 int count = 0;
561         
562                 while (!done) {
563                                 
564                                 if (m->control_pressed) { in.close(); out.close(); return 0; }
565                                 
566                                 Sequence seq(in); m->gobble(in);
567                                 if (seq.getName() != "") {
568                                         string align = seq.getAligned();
569                                         string filterSeq = "";
570                                         
571                                         for(int j=0;j<alignmentLength;j++){
572                                                 if(filter[j] == '1'){
573                                                         filterSeq += align[j];
574                                                 }
575                                         }
576                                         
577                                         out << '>' << seq.getName() << endl << filterSeq << endl;
578                                 count++;
579                         }
580                         
581                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
582                                 unsigned long long pos = in.tellg();
583                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
584                         #else
585                                 if (in.eof()) { break; }
586                         #endif
587                         
588                         //report progress
589                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
590                 }
591                 //report progress
592                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
593                 
594                 
595                 out.close();
596                 in.close();
597                 
598                 return count;
599         }
600         catch(exception& e) {
601                 m->errorOut(e, "FilterSeqsCommand", "driverRunFilter");
602                 exit(1);
603         }
604 }
605 /**************************************************************************************************/
606
607 int FilterSeqsCommand::createProcessesRunFilter(string F, string filename, string filteredFastaName) {
608         try {
609         
610         int process = 1;
611                 int num = 0;
612                 processIDS.clear();
613         
614 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
615                 
616                 
617                 //loop through and create all the processes you want
618                 while (process != processors) {
619                         int pid = fork();
620                         
621                         if (pid > 0) {
622                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
623                                 process++;
624                         }else if (pid == 0){
625                                 string filteredFasta = filename + toString(getpid()) + ".temp";
626                                 num = driverRunFilter(F, filteredFasta, filename, lines[process]);
627                                 
628                                 //pass numSeqs to parent
629                                 ofstream out;
630                                 string tempFile = filename +  toString(getpid()) + ".num.temp";
631                                 m->openOutputFile(tempFile, out);
632                                 out << num << endl;
633                                 out.close();
634                                 
635                                 exit(0);
636                         }else { 
637                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
638                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
639                                 exit(0);
640                         }
641                 }
642                 
643         num = driverRunFilter(F, filteredFastaName, filename, lines[0]);
644         
645                 //force parent to wait until all the processes are done
646                 for (int i=0;i<processIDS.size();i++) { 
647                         int temp = processIDS[i];
648                         wait(&temp);
649                 }       
650                                         
651                 for (int i = 0; i < processIDS.size(); i++) {
652                         ifstream in;
653                         string tempFile =  filename + toString(processIDS[i]) + ".num.temp";
654                         m->openInputFile(tempFile, in);
655                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
656                         in.close(); m->mothurRemove(tempFile);
657             
658             m->appendFiles((filename + toString(processIDS[i]) + ".temp"), filteredFastaName);
659             m->mothurRemove((filename + toString(processIDS[i]) + ".temp"));
660                 }
661                
662 #else
663         
664         //////////////////////////////////////////////////////////////////////////////////////////////////////
665                 //Windows version shared memory, so be careful when passing variables through the filterData struct. 
666                 //Above fork() will clone, so memory is separate, but that's not the case with windows, 
667                 //Taking advantage of shared memory to allow both threads to add info to F.
668                 //////////////////////////////////////////////////////////////////////////////////////////////////////
669                 
670                 vector<filterRunData*> pDataArray; 
671                 DWORD   dwThreadIdArray[processors-1];
672                 HANDLE  hThreadArray[processors-1]; 
673                 
674                 //Create processor worker threads.
675                 for( int i=0; i<processors-1; i++){
676                         
677             string extension = "";
678                         if (i != 0) { extension = toString(i) + ".temp"; }
679             
680                         filterRunData* tempFilter = new filterRunData(filter, filename, (filteredFastaName + extension), m, lines[i]->start, lines[i]->end, alignmentLength, i);
681                         pDataArray.push_back(tempFilter);
682                         processIDS.push_back(i);
683             
684                         hThreadArray[i] = CreateThread(NULL, 0, MyRunFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);   
685                 }
686         
687         num = driverRunFilter(F, (filteredFastaName + toString(processors-1) + ".temp"), filename, lines[processors-1]);
688         
689                 //Wait until all threads have terminated.
690                 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
691                 
692                 //Close all thread handles and free memory allocations.
693                 for(int i=0; i < pDataArray.size(); i++){
694                         num += pDataArray[i]->count;
695             CloseHandle(hThreadArray[i]);
696                         delete pDataArray[i];
697                 }
698         
699         for (int i = 1; i < processors; i++) {
700             m->appendFiles((filteredFastaName + toString(i) + ".temp"), filteredFastaName);
701             m->mothurRemove((filteredFastaName + toString(i) + ".temp"));
702                 }
703 #endif  
704         
705         return num;
706         
707         }
708         catch(exception& e) {
709                 m->errorOut(e, "FilterSeqsCommand", "createProcessesRunFilter");
710                 exit(1);
711         }
712 }
713 /**************************************************************************************/
714 string FilterSeqsCommand::createFilter() {      
715         try {
716                 string filterString = "";                       
717                 Filters F;
718                 
719                 if (soft != 0)                  {  F.setSoft(soft);             }
720                 if (trump != '*')               {  F.setTrump(trump);   }
721                 
722                 F.setLength(alignmentLength);
723                 
724                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
725                         F.initialize();
726                 }
727                 
728                 if(hard.compare("") != 0)       {       F.doHard(hard);         }
729                 else                                            {       F.setFilter(string(alignmentLength, '1'));      }
730                 
731                 numSeqs = 0;
732                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
733                         for (int s = 0; s < fastafileNames.size(); s++) {
734                         
735                                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
736                         
737 #ifdef USE_MPI  
738                                 int pid, numSeqsPerProcessor, num; 
739                                 int tag = 2001;
740                                 vector<unsigned long long> MPIPos;
741                                 
742                                 MPI_Status status; 
743                                 MPI_File inMPI; 
744                                 MPI_Comm_size(MPI_COMM_WORLD, &processors);
745                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
746                                                         
747                                 //char* tempFileName = new char(fastafileNames[s].length());
748                                 //tempFileName = &(fastafileNames[s][0]);
749                                 
750                                 char tempFileName[1024];
751                                 strcpy(tempFileName, fastafileNames[s].c_str());
752                 
753                                 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
754                                 
755                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
756                                 
757                                 if (pid == 0) { //you are the root process
758                                                 MPIPos = m->setFilePosFasta(fastafileNames[s], num); //fills MPIPos, returns numSeqs
759                                                 numSeqs += num;
760                                                 
761                                                 //send file positions to all processes
762                                                 for(int i = 1; i < processors; i++) { 
763                                                         MPI_Send(&num, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
764                                                         MPI_Send(&MPIPos[0], (num+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
765                                                 }
766                                                                 
767                                                 //figure out how many sequences you have to do
768                                                 numSeqsPerProcessor = num / processors;
769                                                 int startIndex =  pid * numSeqsPerProcessor;
770                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
771                                                 
772                                 
773                                                 //do your part
774                                                 MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI, MPIPos);
775                                                 
776                                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
777                                                                                                 
778                                 }else { //i am the child process
779                                         MPI_Recv(&num, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
780                                         MPIPos.resize(num+1);
781                                         numSeqs += num;
782                                         MPI_Recv(&MPIPos[0], (num+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
783                                         
784                                         //figure out how many sequences you have to align
785                                         numSeqsPerProcessor = num / processors;
786                                         int startIndex =  pid * numSeqsPerProcessor;
787                                         if(pid == (processors - 1)){    numSeqsPerProcessor = num - pid * numSeqsPerProcessor;  }
788                                         
789                                         
790                                         //do your part
791                                         MPICreateFilter(startIndex, numSeqsPerProcessor, F, inMPI,  MPIPos);
792                                         
793                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  return 0;  }
794                                 }
795                                 
796                                 MPI_File_close(&inMPI);
797                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
798                                 
799 #else
800                                 
801                 vector<unsigned long long> positions;
802                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
803                                 positions = m->divideFile(fastafileNames[s], processors);
804                                 for (int i = 0; i < (positions.size()-1); i++) {
805                                         lines.push_back(new linePair(positions[i], positions[(i+1)]));
806                                 }       
807                                 
808                                 if(processors == 1){
809                                         int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
810                                         numSeqs += numFastaSeqs;
811                                 }else{
812                                         int numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]); 
813                                         numSeqs += numFastaSeqs;
814                                 }
815                 #else
816                 if(processors == 1){
817                     lines.push_back(new linePair(0, 1000));
818                     int numFastaSeqs = driverCreateFilter(F, fastafileNames[s], lines[0]);
819                     numSeqs += numFastaSeqs;
820                                 }else {
821                     int numFastaSeqs = 0;
822                     positions = m->setFilePosFasta(fastafileNames[s], numFastaSeqs); 
823                     
824                     //figure out how many sequences you have to process
825                     int numSeqsPerProcessor = numFastaSeqs / processors;
826                     for (int i = 0; i < processors; i++) {
827                         int startIndex =  i * numSeqsPerProcessor;
828                         if(i == (processors - 1)){      numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;   }
829                         lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
830                     }
831                     
832                     numFastaSeqs = createProcessesCreateFilter(F, fastafileNames[s]); 
833                                         numSeqs += numFastaSeqs;
834                 }
835                 #endif
836                 //save the file positions so we can reuse them in the runFilter function
837                 savedPositions[s] = positions;
838                 
839                                 if (m->control_pressed) {  return filterString; }
840 #endif
841                         
842                         }
843                 }
844
845
846 #ifdef USE_MPI  
847                 int pid;
848                 int Atag = 1; int Ttag = 2; int Ctag = 3; int Gtag = 4; int Gaptag = 5;
849                 MPI_Status status;
850                 
851                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
852                 
853                 if(trump != '*' || m->isTrue(vertical) || soft != 0){
854                         
855                         if (pid == 0) { //only one process should output the filter
856                         
857                                 vector<int> temp; temp.resize(alignmentLength+1);
858                                                                 
859                                 //get the frequencies from the child processes
860                                 for(int i = 1; i < processors; i++) { 
861                                 
862                                         for (int j = 0; j < 5; j++) {
863                                         
864                                                 MPI_Recv(&temp[0], (alignmentLength+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status); 
865                                                 int receiveTag = temp[temp.size()-1];  //child process added a int to the end to indicate what letter count this is for
866                                                 
867                                                 if (receiveTag == Atag) { //you are recieveing the A frequencies
868                                                         for (int k = 0; k < alignmentLength; k++) {             F.a[k] += temp[k];      }
869                                                 }else if (receiveTag == Ttag) { //you are recieveing the T frequencies
870                                                         for (int k = 0; k < alignmentLength; k++) {             F.t[k] += temp[k];      }
871                                                 }else if (receiveTag == Ctag) { //you are recieveing the C frequencies
872                                                         for (int k = 0; k < alignmentLength; k++) {             F.c[k] += temp[k];      }
873                                                 }else if (receiveTag == Gtag) { //you are recieveing the G frequencies
874                                                         for (int k = 0; k < alignmentLength; k++) {             F.g[k] += temp[k];      }
875                                                 }else if (receiveTag == Gaptag) { //you are recieveing the gap frequencies
876                                                         for (int k = 0; k < alignmentLength; k++) {             F.gap[k] += temp[k];    }
877                                                 }
878                                         }
879                                 } 
880                         }else{
881                         
882                                 //send my fequency counts
883                                 F.a.push_back(Atag);
884                                 int ierr = MPI_Send(&(F.a[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
885                                 F.t.push_back(Ttag);
886                                 ierr = MPI_Send (&(F.t[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
887                                 F.c.push_back(Ctag);
888                                 ierr = MPI_Send(&(F.c[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
889                                 F.g.push_back(Gtag);
890                                 ierr = MPI_Send(&(F.g[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
891                                 F.gap.push_back(Gaptag);
892                                 ierr = MPI_Send(&(F.gap[0]), (alignmentLength+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
893                         }
894                         
895                 }
896                 
897                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
898                 
899                 if (pid == 0) { //only one process should output the filter
900 #endif
901
902                 F.setNumSeqs(numSeqs);
903                 if(m->isTrue(vertical) == 1)    {       F.doVertical(); }
904                 if(soft != 0)                           {       F.doSoft();             }
905                 filterString = F.getFilter();
906                 
907 #ifdef USE_MPI
908                 //send filter string to kids
909                 //for(int i = 1; i < processors; i++) { 
910                 //      MPI_Send(&filterString[0], alignmentLength, MPI_CHAR, i, 2001, MPI_COMM_WORLD);
911                 //}
912                 MPI_Bcast(&filterString[0], alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
913         }else{
914                 //recieve filterString
915                 char* tempBuf = new char[alignmentLength];
916                 //MPI_Recv(&tempBuf[0], alignmentLength, MPI_CHAR, 0, 2001, MPI_COMM_WORLD, &status);
917                 MPI_Bcast(tempBuf, alignmentLength, MPI_CHAR, 0, MPI_COMM_WORLD);
918                 
919                 filterString = tempBuf;
920                 if (filterString.length() > alignmentLength) { filterString = filterString.substr(0, alignmentLength);  }
921                 delete tempBuf; 
922         }
923         
924         MPI_Barrier(MPI_COMM_WORLD);
925 #endif
926             
927                 return filterString;
928         }
929         catch(exception& e) {
930                 m->errorOut(e, "FilterSeqsCommand", "createFilter");
931                 exit(1);
932         }
933 }
934 /**************************************************************************************/
935 int FilterSeqsCommand::driverCreateFilter(Filters& F, string filename, linePair* filePos) {     
936         try {
937                 
938                 ifstream in;
939                 m->openInputFile(filename, in);
940                                 
941                 in.seekg(filePos->start);
942
943                 bool done = false;
944                 int count = 0;
945         
946                 while (!done) {
947                                 
948                         if (m->control_pressed) { in.close(); return 1; }
949                                         
950                         Sequence seq(in); m->gobble(in);
951                         if (seq.getName() != "") {
952                                         if (seq.getAligned().length() != alignmentLength) { m->mothurOut("Sequences are not all the same length, please correct."); m->mothurOutEndLine(); m->control_pressed = true;  }
953                                         
954                                         if(trump != '*')                        {       F.doTrump(seq);         }
955                                         if(m->isTrue(vertical) || soft != 0)    {       F.getFreqs(seq);        }
956                                         cout.flush();
957                                         count++;
958                         }
959                         
960                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
961                                 unsigned long long pos = in.tellg();
962                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
963                         #else
964                                 if (in.eof()) { break; }
965                         #endif
966                         
967                         //report progress
968                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
969                 }
970                 //report progress
971                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
972                 in.close();
973                 
974                 return count;
975         }
976         catch(exception& e) {
977                 m->errorOut(e, "FilterSeqsCommand", "driverCreateFilter");
978                 exit(1);
979         }
980 }
981 #ifdef USE_MPI
982 /**************************************************************************************/
983 int FilterSeqsCommand::MPICreateFilter(int start, int num, Filters& F, MPI_File& inMPI, vector<unsigned long long>& MPIPos) {   
984         try {
985                 
986                 MPI_Status status; 
987                 int pid;
988                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
989                 
990                 for(int i=0;i<num;i++){
991                         
992                         if (m->control_pressed) { return 0; }
993                         
994                         //read next sequence
995                         int length = MPIPos[start+i+1] - MPIPos[start+i];
996         
997                         char* buf4 = new char[length];
998                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
999                         
1000                         string tempBuf = buf4;
1001                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
1002                         istringstream iss (tempBuf,istringstream::in);
1003                         delete buf4;
1004
1005                         Sequence seq(iss);  
1006
1007                         if (seq.getAligned().length() != alignmentLength) {  cout << "Alignment length is " << alignmentLength << " and sequence " << seq.getName() << " has length " << seq.getAligned().length() << ", please correct." << endl; exit(1);  }
1008                         
1009                         if(trump != '*'){       F.doTrump(seq); }
1010                         if(m->isTrue(vertical) || soft != 0){   F.getFreqs(seq);        }
1011                         cout.flush();
1012                                                 
1013                         //report progress
1014                         if((i+1) % 100 == 0){   cout << (i+1) << endl;   m->mothurOutJustToLog(toString(i+1) + "\n");   }
1015                 }
1016                 
1017                 //report progress
1018                 if((num) % 100 != 0){   cout << num << endl; m->mothurOutJustToLog(toString(num) + "\n");       }
1019                 
1020                 return 0;
1021         }
1022         catch(exception& e) {
1023                 m->errorOut(e, "FilterSeqsCommand", "MPICreateFilter");
1024                 exit(1);
1025         }
1026 }
1027 #endif
1028 /**************************************************************************************************/
1029
1030 int FilterSeqsCommand::createProcessesCreateFilter(Filters& F, string filename) {
1031         try {
1032         int process = 1;
1033                 int num = 0;
1034                 processIDS.clear();
1035
1036 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1037                                 
1038                 //loop through and create all the processes you want
1039                 while (process != processors) {
1040                         int pid = fork();
1041                         
1042                         if (pid > 0) {
1043                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
1044                                 process++;
1045                         }else if (pid == 0){
1046                                 //reset child's filter counts to 0;
1047                                 F.a.clear(); F.a.resize(alignmentLength, 0);
1048                                 F.t.clear(); F.t.resize(alignmentLength, 0);
1049                                 F.g.clear(); F.g.resize(alignmentLength, 0);
1050                                 F.c.clear(); F.c.resize(alignmentLength, 0);
1051                                 F.gap.clear(); F.gap.resize(alignmentLength, 0);
1052                                 
1053                                 num = driverCreateFilter(F, filename, lines[process]);
1054                                 
1055                                 //write out filter counts to file
1056                                 filename += toString(getpid()) + "filterValues.temp";
1057                                 ofstream out;
1058                                 m->openOutputFile(filename, out);
1059                                 
1060                                 out << num << endl;
1061                                 out << F.getFilter() << endl;
1062                                 for (int k = 0; k < alignmentLength; k++) {             out << F.a[k] << '\t'; }  out << endl;
1063                                 for (int k = 0; k < alignmentLength; k++) {             out << F.t[k] << '\t'; }  out << endl;
1064                                 for (int k = 0; k < alignmentLength; k++) {             out << F.g[k] << '\t'; }  out << endl;
1065                                 for (int k = 0; k < alignmentLength; k++) {             out << F.c[k] << '\t'; }  out << endl;
1066                                 for (int k = 0; k < alignmentLength; k++) {             out << F.gap[k] << '\t'; }  out << endl;
1067
1068                                 //cout << F.getFilter() << endl;
1069                                 out.close();
1070                                 
1071                                 exit(0);
1072                         }else { 
1073                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
1074                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
1075                                 exit(0);
1076                         }
1077                 }
1078                 
1079                 //parent do your part
1080                 num = driverCreateFilter(F, filename, lines[0]);
1081                 
1082                 //force parent to wait until all the processes are done
1083                 for (int i=0;i<(processors-1);i++) { 
1084                         int temp = processIDS[i];
1085                         wait(&temp);
1086                 }
1087                 
1088                 //parent reads in and combines Filter info
1089                 for (int i = 0; i < processIDS.size(); i++) {
1090                         string tempFilename = filename + toString(processIDS[i]) + "filterValues.temp";
1091                         ifstream in;
1092                         m->openInputFile(tempFilename, in);
1093                         
1094                         int temp, tempNum;
1095                         string tempFilterString;
1096
1097                         in >> tempNum; m->gobble(in); num += tempNum;
1098
1099                         in >> tempFilterString;
1100                         F.mergeFilter(tempFilterString);
1101
1102                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.a[k] += temp; }           m->gobble(in);
1103                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.t[k] += temp; }           m->gobble(in);
1104                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.g[k] += temp; }           m->gobble(in);
1105                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.c[k] += temp; }           m->gobble(in);
1106                         for (int k = 0; k < alignmentLength; k++) {             in >> temp; F.gap[k] += temp; } m->gobble(in);
1107                                 
1108                         in.close();
1109                         m->mothurRemove(tempFilename);
1110                 }
1111                 
1112                 
1113 #else
1114         
1115         //////////////////////////////////////////////////////////////////////////////////////////////////////
1116                 //Windows version shared memory, so be careful when passing variables through the filterData struct. 
1117                 //Above fork() will clone, so memory is separate, but that's not the case with windows, 
1118                 //Taking advantage of shared memory to allow both threads to add info to F.
1119                 //////////////////////////////////////////////////////////////////////////////////////////////////////
1120                 
1121                 vector<filterData*> pDataArray; 
1122                 DWORD   dwThreadIdArray[processors];
1123                 HANDLE  hThreadArray[processors]; 
1124                 
1125                 //Create processor worker threads.
1126                 for( int i=0; i<processors; i++ ){
1127                         
1128                         filterData* tempFilter = new filterData(filename, m, lines[i]->start, lines[i]->end, alignmentLength, trump, vertical, soft, hard, i);
1129                         pDataArray.push_back(tempFilter);
1130                         processIDS.push_back(i);
1131             
1132                         hThreadArray[i] = CreateThread(NULL, 0, MyCreateFilterThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);   
1133                 }
1134         
1135                 //Wait until all threads have terminated.
1136                 WaitForMultipleObjects(processors, hThreadArray, TRUE, INFINITE);
1137                 
1138                 //Close all thread handles and free memory allocations.
1139                 for(int i=0; i < pDataArray.size(); i++){
1140                         num += pDataArray[i]->count;
1141             F.mergeFilter(pDataArray[i]->F.getFilter());
1142             
1143                         for (int k = 0; k < alignmentLength; k++) {      F.a[k] += pDataArray[i]->F.a[k];       }
1144                         for (int k = 0; k < alignmentLength; k++) {      F.t[k] += pDataArray[i]->F.t[k];       }
1145                         for (int k = 0; k < alignmentLength; k++) {      F.g[k] += pDataArray[i]->F.g[k];       }
1146                         for (int k = 0; k < alignmentLength; k++) {      F.c[k] += pDataArray[i]->F.c[k];       }
1147                         for (int k = 0; k < alignmentLength; k++) {      F.gap[k] += pDataArray[i]->F.gap[k];   }
1148
1149                         CloseHandle(hThreadArray[i]);
1150                         delete pDataArray[i];
1151                 }
1152                 
1153 #endif  
1154         return num;
1155         
1156         }
1157         catch(exception& e) {
1158                 m->errorOut(e, "FilterSeqsCommand", "createProcessesCreateFilter");
1159                 exit(1);
1160         }
1161 }
1162 /**************************************************************************************/