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