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