]> git.donarmstrong.com Git - mothur.git/blob - screenseqscommand.cpp
working on current change
[mothur.git] / screenseqscommand.cpp
1 /*
2  *  screenseqscommand.cpp
3  *  Mothur
4  *
5  *  Created by Pat Schloss on 6/3/09.
6  *  Copyright 2009 Patrick D. Schloss. All rights reserved.
7  *
8  */
9
10 #include "screenseqscommand.h"
11 #include "sequence.hpp"
12
13 //**********************************************************************************************************************
14 vector<string> ScreenSeqsCommand::setParameters(){      
15         try {
16                 CommandParameter pfasta("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pfasta);
17                 CommandParameter pname("name", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(pname);
18                 CommandParameter pgroup("group", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(pgroup);
19                 CommandParameter pqfile("qfile", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(pqfile);
20                 CommandParameter palignreport("alignreport", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(palignreport);
21                 CommandParameter pstart("start", "Number", "", "-1", "", "", "",false,false); parameters.push_back(pstart);
22                 CommandParameter pend("end", "Number", "", "-1", "", "", "",false,false); parameters.push_back(pend);
23                 CommandParameter pmaxambig("maxambig", "Number", "", "-1", "", "", "",false,false); parameters.push_back(pmaxambig);
24                 CommandParameter pmaxhomop("maxhomop", "Number", "", "-1", "", "", "",false,false); parameters.push_back(pmaxhomop);
25                 CommandParameter pminlength("minlength", "Number", "", "-1", "", "", "",false,false); parameters.push_back(pminlength);
26                 CommandParameter pmaxlength("maxlength", "Number", "", "-1", "", "", "",false,false); parameters.push_back(pmaxlength);
27                 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
28                 CommandParameter pcriteria("criteria", "Number", "", "90", "", "", "",false,false); parameters.push_back(pcriteria);
29                 CommandParameter poptimize("optimize", "Multiple", "none-start-end-maxambig-maxhomop-minlength-maxlength", "none", "", "", "",true,false); parameters.push_back(poptimize);
30                 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
31                 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
32                 
33                 vector<string> myArray;
34                 for (int i = 0; i < parameters.size(); i++) {   myArray.push_back(parameters[i].name);          }
35                 return myArray;
36         }
37         catch(exception& e) {
38                 m->errorOut(e, "ScreenSeqsCommand", "setParameters");
39                 exit(1);
40         }
41 }
42 //**********************************************************************************************************************
43 string ScreenSeqsCommand::getHelpString(){      
44         try {
45                 string helpString = "";
46                 helpString += "The screen.seqs command reads a fastafile and creates .....\n";
47                 helpString += "The screen.seqs command parameters are fasta, start, end, maxambig, maxhomop, minlength, maxlength, name, group, qfile, optimize, criteria and processors.\n";
48                 helpString += "The fasta parameter is required.\n";
49                 helpString += "The start parameter .... The default is -1.\n";
50                 helpString += "The end parameter .... The default is -1.\n";
51                 helpString += "The maxambig parameter allows you to set the maximum number of ambigious bases allowed. The default is -1.\n";
52                 helpString += "The maxhomop parameter allows you to set a maximum homopolymer length. \n";
53                 helpString += "The minlength parameter allows you to set and minimum sequence length. \n";
54                 helpString += "The maxlength parameter allows you to set and maximum sequence length. \n";
55                 helpString += "The processors parameter allows you to specify the number of processors to use while running the command. The default is 1.\n";
56                 helpString += "The optimize and criteria parameters allow you set the start, end, maxabig, maxhomop, minlength and maxlength parameters relative to your set of sequences .\n";
57                 helpString += "For example optimize=start-end, criteria=90, would set the start and end values to the position 90% of your sequences started and ended.\n";
58                 helpString += "The name parameter allows you to provide a namesfile, and the group parameter allows you to provide a groupfile.\n";
59                 helpString += "The screen.seqs command should be in the following format: \n";
60                 helpString += "screen.seqs(fasta=yourFastaFile, name=youNameFile, group=yourGroupFIle, start=yourStart, end=yourEnd, maxambig=yourMaxambig,  \n";
61                 helpString += "maxhomop=yourMaxhomop, minlength=youMinlength, maxlength=yourMaxlength)  \n";    
62                 helpString += "Example screen.seqs(fasta=abrecovery.fasta, name=abrecovery.names, group=abrecovery.groups, start=..., end=..., maxambig=..., maxhomop=..., minlength=..., maxlength=...).\n";
63                 helpString += "Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n";
64                 return helpString;
65         }
66         catch(exception& e) {
67                 m->errorOut(e, "ScreenSeqsCommand", "getHelpString");
68                 exit(1);
69         }
70 }
71 //**********************************************************************************************************************
72 ScreenSeqsCommand::ScreenSeqsCommand(){ 
73         try {
74                 abort = true; calledHelp = true; 
75                 setParameters();
76                 vector<string> tempOutNames;
77                 outputTypes["fasta"] = tempOutNames;
78                 outputTypes["name"] = tempOutNames;
79                 outputTypes["group"] = tempOutNames;
80                 outputTypes["alignreport"] = tempOutNames;
81                 outputTypes["accnos"] = tempOutNames;
82                 outputTypes["qfile"] = tempOutNames;
83         }
84         catch(exception& e) {
85                 m->errorOut(e, "ScreenSeqsCommand", "ScreenSeqsCommand");
86                 exit(1);
87         }
88 }
89 //***************************************************************************************************************
90
91 ScreenSeqsCommand::ScreenSeqsCommand(string option)  {
92         try {
93                 abort = false; calledHelp = false;   
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("screen.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["name"] = tempOutNames;
117                         outputTypes["group"] = tempOutNames;
118                         outputTypes["alignreport"] = tempOutNames;
119                         outputTypes["accnos"] = tempOutNames;
120                         outputTypes["qfile"] = tempOutNames;
121                         
122                         //if the user changes the input directory command factory will send this info to us in the output parameter 
123                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
124                         if (inputDir == "not found"){   inputDir = "";          }
125                         else {
126                                 string path;
127                                 it = parameters.find("fasta");
128                                 //user has given a template file
129                                 if(it != parameters.end()){ 
130                                         path = m->hasPath(it->second);
131                                         //if the user has not given a path then, add inputdir. else leave path alone.
132                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
133                                 }
134                                 
135                                 it = parameters.find("group");
136                                 //user has given a template file
137                                 if(it != parameters.end()){ 
138                                         path = m->hasPath(it->second);
139                                         //if the user has not given a path then, add inputdir. else leave path alone.
140                                         if (path == "") {       parameters["group"] = inputDir + it->second;            }
141                                 }
142                                 
143                                 it = parameters.find("name");
144                                 //user has given a template file
145                                 if(it != parameters.end()){ 
146                                         path = m->hasPath(it->second);
147                                         //if the user has not given a path then, add inputdir. else leave path alone.
148                                         if (path == "") {       parameters["name"] = inputDir + it->second;             }
149                                 }
150                                 
151                                 it = parameters.find("alignreport");
152                                 //user has given a template file
153                                 if(it != parameters.end()){ 
154                                         path = m->hasPath(it->second);
155                                         //if the user has not given a path then, add inputdir. else leave path alone.
156                                         if (path == "") {       parameters["alignreport"] = inputDir + it->second;              }
157                                 }
158                                 
159                                 it = parameters.find("qfile");
160                                 //user has given a template file
161                                 if(it != parameters.end()){ 
162                                         path = m->hasPath(it->second);
163                                         //if the user has not given a path then, add inputdir. else leave path alone.
164                                         if (path == "") {       parameters["qfile"] = inputDir + it->second;            }
165                                 }
166                                 
167                         }
168
169                         //check for required parameters
170                         fastafile = validParameter.validFile(parameters, "fasta", true);
171                         if (fastafile == "not found") {                         
172                                 fastafile = m->getFastaFile(); 
173                                 if (fastafile != "") { m->mothurOut("Using " + fastafile + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
174                                 else {  m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
175                         }
176                         else if (fastafile == "not open") { abort = true; }
177                         else { m->setFastaFile(fastafile); }
178         
179                         groupfile = validParameter.validFile(parameters, "group", true);
180                         if (groupfile == "not open") { abort = true; }  
181                         else if (groupfile == "not found") { groupfile = ""; }
182                         else { m->setGroupFile(groupfile); }
183                         
184                         qualfile = validParameter.validFile(parameters, "qfile", true);
185                         if (qualfile == "not open") { abort = true; }   
186                         else if (qualfile == "not found") { qualfile = ""; }
187                         else { m->setQualFile(qualfile); }
188                         
189                         namefile = validParameter.validFile(parameters, "name", true);
190                         if (namefile == "not open") { namefile = ""; abort = true; }
191                         else if (namefile == "not found") { namefile = ""; }    
192                         else { m->setNameFile(namefile); }
193                         
194                         alignreport = validParameter.validFile(parameters, "alignreport", true);
195                         if (alignreport == "not open") { abort = true; }
196                         else if (alignreport == "not found") { alignreport = ""; }      
197                         
198                         //if the user changes the output directory command factory will send this info to us in the output parameter 
199                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
200                                 outputDir = ""; 
201                                 outputDir += m->hasPath(fastafile); //if user entered a file with a path then preserve it       
202                         }
203
204                         //check for optional parameter and set defaults
205                         // ...at some point should added some additional type checking...
206                         string temp;
207                         temp = validParameter.validFile(parameters, "start", false);            if (temp == "not found") { temp = "-1"; }
208                         convert(temp, startPos); 
209                 
210                         temp = validParameter.validFile(parameters, "end", false);                      if (temp == "not found") { temp = "-1"; }
211                         convert(temp, endPos);  
212
213                         temp = validParameter.validFile(parameters, "maxambig", false);         if (temp == "not found") { temp = "-1"; }
214                         convert(temp, maxAmbig);  
215
216                         temp = validParameter.validFile(parameters, "maxhomop", false);         if (temp == "not found") { temp = "-1"; }
217                         convert(temp, maxHomoP);  
218
219                         temp = validParameter.validFile(parameters, "minlength", false);        if (temp == "not found") { temp = "-1"; }
220                         convert(temp, minLength); 
221                         
222                         temp = validParameter.validFile(parameters, "maxlength", false);        if (temp == "not found") { temp = "-1"; }
223                         convert(temp, maxLength); 
224                         
225                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = m->getProcessors();      }
226                         m->setProcessors(temp);
227                         convert(temp, processors);
228                         
229                         temp = validParameter.validFile(parameters, "optimize", false); //optimizing trumps the optimized values original value
230                         if (temp == "not found"){       temp = "none";          }
231                         m->splitAtDash(temp, optimize);         
232                         
233                         //check for invalid optimize options
234                         set<string> validOptimizers;
235                         validOptimizers.insert("none"); validOptimizers.insert("start"); validOptimizers.insert("end"); validOptimizers.insert("maxambig"); validOptimizers.insert("maxhomop"); validOptimizers.insert("minlength"); validOptimizers.insert("maxlength");
236                         for (int i = 0; i < optimize.size(); i++) { 
237                                 if (validOptimizers.count(optimize[i]) == 0) { 
238                                         m->mothurOut(optimize[i] + " is not a valid optimizer. Valid options are start, end, maxambig, maxhomop, minlength and maxlength."); m->mothurOutEndLine();
239                                         optimize.erase(optimize.begin()+i);
240                                         i--;
241                                 }
242                         }
243                         
244                         if (optimize.size() == 1) { if (optimize[0] == "none") { optimize.clear(); } }
245                         
246                         temp = validParameter.validFile(parameters, "criteria", false); if (temp == "not found"){       temp = "90";                            }
247                         convert(temp, criteria); 
248                 }
249
250         }
251         catch(exception& e) {
252                 m->errorOut(e, "ScreenSeqsCommand", "ScreenSeqsCommand");
253                 exit(1);
254         }
255 }
256 //***************************************************************************************************************
257
258 int ScreenSeqsCommand::execute(){
259         try{
260                 
261                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
262                 
263                 //if the user want to optimize we need to know the 90% mark
264                 vector<unsigned long int> positions;
265                 if (optimize.size() != 0) {  //get summary is paralellized so we need to divideFile, no need to do this step twice so I moved it here
266                         //use the namefile to optimize correctly
267                         if (namefile != "") { nameMap = m->readNames(namefile); }
268                         getSummary(positions); 
269                 } 
270                 else { 
271                         positions = m->divideFile(fastafile, processors);
272                         for (int i = 0; i < (positions.size()-1); i++) {
273                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
274                         }       
275                 }
276                                 
277                 string goodSeqFile = outputDir + m->getRootName(m->getSimpleName(fastafile)) + "good" + m->getExtension(fastafile);
278                 string badAccnosFile =  outputDir + m->getRootName(m->getSimpleName(fastafile)) + "bad.accnos";
279                 
280                 int numFastaSeqs = 0;
281                 set<string> badSeqNames;
282                 int start = time(NULL);
283                 
284 #ifdef USE_MPI  
285                         int pid, numSeqsPerProcessor; 
286                         int tag = 2001;
287                         vector<unsigned long int> MPIPos;
288                         
289                         MPI_Status status; 
290                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
291                         MPI_Comm_size(MPI_COMM_WORLD, &processors); 
292
293                         MPI_File inMPI;
294                         MPI_File outMPIGood;
295                         MPI_File outMPIBadAccnos;
296                         
297                         int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
298                         int inMode=MPI_MODE_RDONLY; 
299                         
300                         char outGoodFilename[1024];
301                         strcpy(outGoodFilename, goodSeqFile.c_str());
302
303                         char outBadAccnosFilename[1024];
304                         strcpy(outBadAccnosFilename, badAccnosFile.c_str());
305
306                         char inFileName[1024];
307                         strcpy(inFileName, fastafile.c_str());
308                         
309                         MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
310                         MPI_File_open(MPI_COMM_WORLD, outGoodFilename, outMode, MPI_INFO_NULL, &outMPIGood);
311                         MPI_File_open(MPI_COMM_WORLD, outBadAccnosFilename, outMode, MPI_INFO_NULL, &outMPIBadAccnos);
312                         
313                         if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIGood); MPI_File_close(&outMPIBadAccnos); return 0; }
314                         
315                         if (pid == 0) { //you are the root process 
316                                 
317                                 MPIPos = m->setFilePosFasta(fastafile, numFastaSeqs); //fills MPIPos, returns numSeqs
318                                 
319                                 //send file positions to all processes
320                                 for(int i = 1; i < processors; i++) { 
321                                         MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
322                                         MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
323                                 }
324                                 
325                                 //figure out how many sequences you have to align
326                                 numSeqsPerProcessor = numFastaSeqs / processors;
327                                 int startIndex =  pid * numSeqsPerProcessor;
328                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
329                 //      cout << pid << '\t' << numSeqsPerProcessor << '\t' <<   startIndex << endl;
330                                 //align your part
331                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIGood, outMPIBadAccnos, MPIPos, badSeqNames);
332                 //cout << pid << " done" << endl;
333                                 if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIGood);  MPI_File_close(&outMPIBadAccnos);  return 0; }
334
335                                 for (int i = 1; i < processors; i++) {
336                                 
337                                         //get bad lists
338                                         int badSize;
339                                         MPI_Recv(&badSize, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
340                                         /*for (int j = 0; j < badSize; j++) {
341                                                 int length;
342                                                 MPI_Recv(&length, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);  //recv the length of the name
343                                                 char* buf2 = new char[length];                                                                          //make space to recieve it
344                                                 MPI_Recv(buf2, length, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);  //get name
345                                                 
346                                                 string tempBuf = buf2;
347                                                 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
348                                                 delete buf2;
349                                                 
350                                                 badSeqNames.insert(tempBuf);
351                                         }*/
352                                 }
353                         }else{ //you are a child process
354                                 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
355                                 MPIPos.resize(numFastaSeqs+1);
356                                 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
357
358                                 //figure out how many sequences you have to align
359                                 numSeqsPerProcessor = numFastaSeqs / processors;
360                                 int startIndex =  pid * numSeqsPerProcessor;
361                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
362                 //cout << pid << '\t' << numSeqsPerProcessor << '\t' <<         startIndex << endl;             
363                                 //align your part
364                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIGood, outMPIBadAccnos, MPIPos, badSeqNames);
365 //cout << pid << " done" << endl;
366                                 if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIGood);  MPI_File_close(&outMPIBadAccnos); return 0; }
367                                 
368                                 //send bad list 
369                                 int badSize = badSeqNames.size();
370                                 MPI_Send(&badSize, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
371                                 
372                                 /*
373                                 set<string>::iterator it;
374                                 for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {
375                                         string name = *it;
376                                         int length = name.length();
377                                         char* buf2 = new char[length];
378                                         memcpy(buf2, name.c_str(), length);
379                                         
380                                         MPI_Send(&length, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
381                                         MPI_Send(buf2, length, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
382                                 }*/
383                         }
384                         
385                         //close files 
386                         MPI_File_close(&inMPI);
387                         MPI_File_close(&outMPIGood);
388                         MPI_File_close(&outMPIBadAccnos);
389                         MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
390                                         
391 #else
392                                                 
393         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
394                         if(processors == 1){
395                                 numFastaSeqs = driver(lines[0], goodSeqFile, badAccnosFile, fastafile, badSeqNames);
396                                 
397                                 if (m->control_pressed) { remove(goodSeqFile.c_str()); return 0; }
398                                 
399                         }else{
400                                 processIDS.resize(0);
401                                 
402                                 numFastaSeqs = createProcesses(goodSeqFile, badAccnosFile, fastafile, badSeqNames); 
403                                 
404                                 rename((goodSeqFile + toString(processIDS[0]) + ".temp").c_str(), goodSeqFile.c_str());
405                                 rename((badAccnosFile + toString(processIDS[0]) + ".temp").c_str(), badAccnosFile.c_str());
406                                 
407                                 //append alignment and report files
408                                 for(int i=1;i<processors;i++){
409                                         m->appendFiles((goodSeqFile + toString(processIDS[i]) + ".temp"), goodSeqFile);
410                                         remove((goodSeqFile + toString(processIDS[i]) + ".temp").c_str());
411                         
412                                         m->appendFiles((badAccnosFile + toString(processIDS[i]) + ".temp"), badAccnosFile);
413                                         remove((badAccnosFile + toString(processIDS[i]) + ".temp").c_str());
414                                 }
415                                 
416                                 if (m->control_pressed) { remove(goodSeqFile.c_str()); return 0; }
417                                 
418                                 //read badSeqs in because root process doesnt know what other "bad" seqs the children found
419                                 ifstream inBad;
420                                 int ableToOpen = m->openInputFile(badAccnosFile, inBad, "no error");
421                                 
422                                 if (ableToOpen == 0) {
423                                         badSeqNames.clear();
424                                         string tempName;
425                                         while (!inBad.eof()) {
426                                                 inBad >> tempName; m->gobble(inBad);
427                                                 badSeqNames.insert(tempName);
428                                         }
429                                         inBad.close();
430                                 }
431                         }
432         #else
433                         numFastaSeqs = driver(lines[0], goodSeqFile, badAccnosFile, fastafile, badSeqNames);
434                         
435                         if (m->control_pressed) { remove(goodSeqFile.c_str()); return 0; }
436                         
437         #endif
438
439 #endif          
440
441                 #ifdef USE_MPI
442                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
443                                         
444                         if (pid == 0) { //only one process should fix files
445                         
446                                 //read accnos file with all names in it, process 0 just has its names
447                                 MPI_File inMPIAccnos;
448                                 MPI_Offset size;
449                         
450                                 char inFileName[1024];
451                                 strcpy(inFileName, badAccnosFile.c_str());
452                         
453                                 MPI_File_open(MPI_COMM_SELF, inFileName, inMode, MPI_INFO_NULL, &inMPIAccnos);  //comm, filename, mode, info, filepointer
454                                 MPI_File_get_size(inMPIAccnos, &size);
455                         
456                                 char* buffer = new char[size];
457                                 MPI_File_read(inMPIAccnos, buffer, size, MPI_CHAR, &status);
458                         
459                                 string tempBuf = buffer;
460                                 if (tempBuf.length() > size) { tempBuf = tempBuf.substr(0, size);  }
461                                 istringstream iss (tempBuf,istringstream::in);
462
463                                 delete buffer;
464                                 MPI_File_close(&inMPIAccnos);
465                                 
466                                 badSeqNames.clear();
467                                 string tempName;
468                                 while (!iss.eof()) {
469                                         iss >> tempName; m->gobble(iss);
470                                         badSeqNames.insert(tempName);
471                                 }
472                 #endif
473                                                                                                                                                                         
474                 if(namefile != "" && groupfile != "")   {       
475                         screenNameGroupFile(badSeqNames);       
476                         if (m->control_pressed) {  remove(goodSeqFile.c_str()); return 0; }
477                 }else if(namefile != "")        {       
478                         screenNameGroupFile(badSeqNames);
479                         if (m->control_pressed) {  remove(goodSeqFile.c_str());  return 0; }    
480                 }else if(groupfile != "")                               {       screenGroupFile(badSeqNames);           }       // this screens just the group
481                 
482                 if (m->control_pressed) { remove(goodSeqFile.c_str());  return 0; }
483
484                 if(alignreport != "")                                   {       screenAlignReport(badSeqNames);         }
485                 if(qualfile != "")                                              {       screenQual(badSeqNames);                        }
486                 
487                 if (m->control_pressed) { remove(goodSeqFile.c_str());  return 0; }
488                 
489                 #ifdef USE_MPI
490                         }
491                 #endif
492
493                 m->mothurOutEndLine();
494                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
495                 m->mothurOut(goodSeqFile); m->mothurOutEndLine();       outputTypes["fasta"].push_back(goodSeqFile);
496                 m->mothurOut(badAccnosFile); m->mothurOutEndLine();      outputTypes["accnos"].push_back(badAccnosFile);
497                 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
498                 m->mothurOutEndLine();
499                 m->mothurOutEndLine();
500                 
501                 //set fasta file as new current fastafile
502                 string current = "";
503                 itTypes = outputTypes.find("fasta");
504                 if (itTypes != outputTypes.end()) {
505                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
506                 }
507                 
508                 itTypes = outputTypes.find("name");
509                 if (itTypes != outputTypes.end()) {
510                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setNameFile(current); }
511                 }
512                 
513                 itTypes = outputTypes.find("group");
514                 if (itTypes != outputTypes.end()) {
515                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setGroupFile(current); }
516                 }
517                 
518                 itTypes = outputTypes.find("qfile");
519                 if (itTypes != outputTypes.end()) {
520                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setQualFile(current); }
521                 }
522
523                 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to screen " + toString(numFastaSeqs) + " sequences.");
524                 m->mothurOutEndLine();
525
526                 return 0;
527         }
528         catch(exception& e) {
529                 m->errorOut(e, "ScreenSeqsCommand", "execute");
530                 exit(1);
531         }
532 }
533
534 //***************************************************************************************************************
535
536 int ScreenSeqsCommand::screenNameGroupFile(set<string> badSeqNames){
537         try {
538                 ifstream inputNames;
539                 m->openInputFile(namefile, inputNames);
540                 set<string> badSeqGroups;
541                 string seqName, seqList, group;
542                 set<string>::iterator it;
543
544                 string goodNameFile = outputDir + m->getRootName(m->getSimpleName(namefile)) + "good" + m->getExtension(namefile);
545                 outputNames.push_back(goodNameFile);  outputTypes["name"].push_back(goodNameFile);
546                 
547                 ofstream goodNameOut;   m->openOutputFile(goodNameFile, goodNameOut);
548                 
549                 while(!inputNames.eof()){
550                         if (m->control_pressed) { goodNameOut.close();  inputNames.close(); remove(goodNameFile.c_str());  return 0; }
551
552                         inputNames >> seqName >> seqList;
553                         it = badSeqNames.find(seqName);
554                                 
555                         if(it != badSeqNames.end()){
556                                 badSeqNames.erase(it);
557                                 
558                                 if(namefile != ""){
559                                         int start = 0;
560                                         for(int i=0;i<seqList.length();i++){
561                                                 if(seqList[i] == ','){
562                                                         badSeqGroups.insert(seqList.substr(start,i-start));
563                                                         start = i+1;
564                                                 }                                       
565                                         }
566                                         badSeqGroups.insert(seqList.substr(start,seqList.length()-start));
567                                 }
568                         }
569                         else{
570                                 goodNameOut << seqName << '\t' << seqList << endl;
571                         }
572                         m->gobble(inputNames);
573                 }
574                 inputNames.close();
575                 goodNameOut.close();
576         
577                 //we were unable to remove some of the bad sequences
578                 if (badSeqNames.size() != 0) {
579                         for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {  
580                                 m->mothurOut("Your namefile does not include the sequence " + *it + " please correct."); 
581                                 m->mothurOutEndLine();
582                         }
583                 }
584
585                 if(groupfile != ""){
586                         
587                         ifstream inputGroups;
588                         m->openInputFile(groupfile, inputGroups);
589
590                         string goodGroupFile = outputDir + m->getRootName(m->getSimpleName(groupfile)) + "good" + m->getExtension(groupfile);
591                         outputNames.push_back(goodGroupFile);   outputTypes["group"].push_back(goodGroupFile);
592                         
593                         ofstream goodGroupOut;  m->openOutputFile(goodGroupFile, goodGroupOut);
594                         
595                         while(!inputGroups.eof()){
596                                 if (m->control_pressed) { goodGroupOut.close(); inputGroups.close(); remove(goodNameFile.c_str());  remove(goodGroupFile.c_str()); return 0; }
597
598                                 inputGroups >> seqName >> group;
599                                 
600                                 it = badSeqGroups.find(seqName);
601                                 
602                                 if(it != badSeqGroups.end()){
603                                         badSeqGroups.erase(it);
604                                 }
605                                 else{
606                                         goodGroupOut << seqName << '\t' << group << endl;
607                                 }
608                                 m->gobble(inputGroups);
609                         }
610                         inputGroups.close();
611                         goodGroupOut.close();
612                         
613                         //we were unable to remove some of the bad sequences
614                         if (badSeqGroups.size() != 0) {
615                                 for (it = badSeqGroups.begin(); it != badSeqGroups.end(); it++) {  
616                                         m->mothurOut("Your groupfile does not include the sequence " + *it + " please correct."); 
617                                         m->mothurOutEndLine();
618                                 }
619                         }
620                 }
621                 
622                 
623                 return 0;
624         
625         }
626         catch(exception& e) {
627                 m->errorOut(e, "ScreenSeqsCommand", "screenNameGroupFile");
628                 exit(1);
629         }
630 }
631 //***************************************************************************************************************
632 int ScreenSeqsCommand::getSummary(vector<unsigned long int>& positions){
633         try {
634                 
635                 vector<int> startPosition;
636                 vector<int> endPosition;
637                 vector<int> seqLength;
638                 vector<int> ambigBases;
639                 vector<int> longHomoPolymer;
640                 
641                 vector<unsigned long int> positions = m->divideFile(fastafile, processors);
642                                 
643                 for (int i = 0; i < (positions.size()-1); i++) {
644                         lines.push_back(new linePair(positions[i], positions[(i+1)]));
645                 }       
646                 
647                 int numSeqs = 0;
648                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
649                         if(processors == 1){
650                                 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, lines[0]);
651                         }else{
652                                 numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile); 
653                         }
654                                 
655                         if (m->control_pressed) {  return 0; }
656                 #else
657                         numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, lines[0]);
658                         if (m->control_pressed) {  return 0; }
659                 #endif
660
661                 sort(startPosition.begin(), startPosition.end());
662                 sort(endPosition.begin(), endPosition.end());
663                 sort(seqLength.begin(), seqLength.end());
664                 sort(ambigBases.begin(), ambigBases.end());
665                 sort(longHomoPolymer.begin(), longHomoPolymer.end());
666                 
667                 //numSeqs is the number of unique seqs, startPosition.size() is the total number of seqs, we want to optimize using all seqs
668                 int criteriaPercentile  = int(startPosition.size() * (criteria / (float) 100));
669                 
670                 for (int i = 0; i < optimize.size(); i++) {
671                         if (optimize[i] == "start") { startPos = startPosition[criteriaPercentile]; m->mothurOut("Optimizing start to " + toString(startPos) + "."); m->mothurOutEndLine(); }
672                         else if (optimize[i] == "end") { int endcriteriaPercentile = int(endPosition.size() * ((100 - criteria) / (float) 100));  endPos = endPosition[endcriteriaPercentile]; m->mothurOut("Optimizing end to " + toString(endPos) + "."); m->mothurOutEndLine();}
673                         else if (optimize[i] == "maxambig") { maxAmbig = ambigBases[criteriaPercentile]; m->mothurOut("Optimizing maxambig to " + toString(maxAmbig) + "."); m->mothurOutEndLine(); }
674                         else if (optimize[i] == "maxhomop") { maxHomoP = longHomoPolymer[criteriaPercentile]; m->mothurOut("Optimizing maxhomop to " + toString(maxHomoP) + "."); m->mothurOutEndLine(); }
675                         else if (optimize[i] == "minlength") { int mincriteriaPercentile = int(seqLength.size() * ((100 - criteria) / (float) 100)); minLength = seqLength[mincriteriaPercentile]; m->mothurOut("Optimizing minlength to " + toString(minLength) + "."); m->mothurOutEndLine(); }
676                         else if (optimize[i] == "maxlength") { maxLength = seqLength[criteriaPercentile]; m->mothurOut("Optimizing maxlength to " + toString(maxLength) + "."); m->mothurOutEndLine(); }
677                 }
678
679                 return 0;
680         }
681         catch(exception& e) {
682                 m->errorOut(e, "ScreenSeqsCommand", "getSummary");
683                 exit(1);
684         }
685 }
686 /**************************************************************************************/
687 int ScreenSeqsCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, linePair* filePos) {   
688         try {
689                 
690                 ifstream in;
691                 m->openInputFile(filename, in);
692                                 
693                 in.seekg(filePos->start);
694
695                 bool done = false;
696                 int count = 0;
697         
698                 while (!done) {
699                                 
700                         if (m->control_pressed) { in.close(); return 1; }
701                                         
702                         Sequence current(in); m->gobble(in);
703         
704                         if (current.getName() != "") {
705                                 int num = 1;
706                                 if (namefile != "") {
707                                         //make sure this sequence is in the namefile, else error 
708                                         map<string, int>::iterator it = nameMap.find(current.getName());
709                                         
710                                         if (it == nameMap.end()) { m->mothurOut("[ERROR]: " + current.getName() + " is not in your namefile, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
711                                         else { num = it->second; }
712                                 }
713                                 
714                                 //for each sequence this sequence represents
715                                 for (int i = 0; i < num; i++) {
716                                         startPosition.push_back(current.getStartPos());
717                                         endPosition.push_back(current.getEndPos());
718                                         seqLength.push_back(current.getNumBases());
719                                         ambigBases.push_back(current.getAmbigBases());
720                                         longHomoPolymer.push_back(current.getLongHomoPolymer());
721                                 }
722                                 
723                                 count++;
724                         }
725                         
726                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
727                                 unsigned long int pos = in.tellg();
728                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
729                         #else
730                                 if (in.eof()) { break; }
731                         #endif
732                         
733                 }
734                 
735                 in.close();
736                 
737                 return count;
738         }
739         catch(exception& e) {
740                 m->errorOut(e, "ScreenSeqsCommand", "driverCreateSummary");
741                 exit(1);
742         }
743 }
744 /**************************************************************************************************/
745 int ScreenSeqsCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename) {
746         try {
747 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
748                 int process = 1;
749                 int num = 0;
750                 processIDS.clear();
751                 
752                 //loop through and create all the processes you want
753                 while (process != processors) {
754                         int pid = fork();
755                         
756                         if (pid > 0) {
757                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
758                                 process++;
759                         }else if (pid == 0){
760                                 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, lines[process]);
761                                 
762                                 //pass numSeqs to parent
763                                 ofstream out;
764                                 string tempFile = fastafile + toString(getpid()) + ".num.temp";
765                                 m->openOutputFile(tempFile, out);
766                                 
767                                 out << num << endl;
768                                 out << startPosition.size() << endl;
769                                 for (int k = 0; k < startPosition.size(); k++)          {               out << startPosition[k] << '\t'; }  out << endl;
770                                 for (int k = 0; k < endPosition.size(); k++)            {               out << endPosition[k] << '\t'; }  out << endl;
771                                 for (int k = 0; k < seqLength.size(); k++)                      {               out << seqLength[k] << '\t'; }  out << endl;
772                                 for (int k = 0; k < ambigBases.size(); k++)                     {               out << ambigBases[k] << '\t'; }  out << endl;
773                                 for (int k = 0; k < longHomoPolymer.size(); k++)        {               out << longHomoPolymer[k] << '\t'; }  out << endl;
774                                 
775                                 out.close();
776                                 
777                                 exit(0);
778                         }else { 
779                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
780                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
781                                 exit(0);
782                         }
783                 }
784                 
785                 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, lines[0]);
786                 
787                 //force parent to wait until all the processes are done
788                 for (int i=0;i<processIDS.size();i++) { 
789                         int temp = processIDS[i];
790                         wait(&temp);
791                 }
792                 
793                 //parent reads in and combine Filter info
794                 for (int i = 0; i < processIDS.size(); i++) {
795                         string tempFilename = fastafile + toString(processIDS[i]) + ".num.temp";
796                         ifstream in;
797                         m->openInputFile(tempFilename, in);
798                         
799                         int temp, tempNum;
800                         in >> tempNum; m->gobble(in); num += tempNum;
801                         in >> tempNum; m->gobble(in);
802                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; startPosition.push_back(temp);              }               m->gobble(in);
803                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; endPosition.push_back(temp);                }               m->gobble(in);
804                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; seqLength.push_back(temp);                  }               m->gobble(in);
805                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; ambigBases.push_back(temp);                 }               m->gobble(in);
806                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; longHomoPolymer.push_back(temp);    }               m->gobble(in);
807                                 
808                         in.close();
809                         remove(tempFilename.c_str());
810                 }
811                 
812                 return num;
813 #endif          
814         }
815         catch(exception& e) {
816                 m->errorOut(e, "ScreenSeqsCommand", "createProcessesCreateSummary");
817                 exit(1);
818         }
819 }
820
821 //***************************************************************************************************************
822
823 int ScreenSeqsCommand::screenGroupFile(set<string> badSeqNames){
824         try {
825                 ifstream inputGroups;
826                 m->openInputFile(groupfile, inputGroups);
827                 string seqName, group;
828                 set<string>::iterator it;
829                 
830                 string goodGroupFile = outputDir + m->getRootName(m->getSimpleName(groupfile)) + "good" + m->getExtension(groupfile);
831                 outputNames.push_back(goodGroupFile);  outputTypes["group"].push_back(goodGroupFile);
832                 ofstream goodGroupOut;  m->openOutputFile(goodGroupFile, goodGroupOut);
833                 
834                 while(!inputGroups.eof()){
835                         if (m->control_pressed) { goodGroupOut.close(); inputGroups.close(); remove(goodGroupFile.c_str()); return 0; }
836
837                         inputGroups >> seqName >> group;
838                         it = badSeqNames.find(seqName);
839                         
840                         if(it != badSeqNames.end()){
841                                 badSeqNames.erase(it);
842                         }
843                         else{
844                                 goodGroupOut << seqName << '\t' << group << endl;
845                         }
846                         m->gobble(inputGroups);
847                 }
848                 
849                 if (m->control_pressed) { goodGroupOut.close();  inputGroups.close(); remove(goodGroupFile.c_str());  return 0; }
850
851                 //we were unable to remove some of the bad sequences
852                 if (badSeqNames.size() != 0) {
853                         for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {  
854                                 m->mothurOut("Your groupfile does not include the sequence " + *it + " please correct."); 
855                                 m->mothurOutEndLine();
856                         }
857                 }
858                 
859                 inputGroups.close();
860                 goodGroupOut.close();
861                 
862                 if (m->control_pressed) { remove(goodGroupFile.c_str());   }
863                 
864                 return 0;
865         
866         }
867         catch(exception& e) {
868                 m->errorOut(e, "ScreenSeqsCommand", "screenGroupFile");
869                 exit(1);
870         }
871 }
872
873 //***************************************************************************************************************
874
875 int ScreenSeqsCommand::screenAlignReport(set<string> badSeqNames){
876         try {
877                 ifstream inputAlignReport;
878                 m->openInputFile(alignreport, inputAlignReport);
879                 string seqName, group;
880                 set<string>::iterator it;
881                 
882                 string goodAlignReportFile = outputDir + m->getRootName(m->getSimpleName(alignreport)) + "good" + m->getExtension(alignreport);
883                 outputNames.push_back(goodAlignReportFile);  outputTypes["alignreport"].push_back(goodAlignReportFile);
884                 ofstream goodAlignReportOut;    m->openOutputFile(goodAlignReportFile, goodAlignReportOut);
885
886                 while (!inputAlignReport.eof()) {               //      need to copy header
887                         char c = inputAlignReport.get();
888                         goodAlignReportOut << c;
889                         if (c == 10 || c == 13){        break;  }       
890                 }
891
892                 while(!inputAlignReport.eof()){
893                         if (m->control_pressed) { goodAlignReportOut.close(); inputAlignReport.close(); remove(goodAlignReportFile.c_str()); return 0; }
894
895                         inputAlignReport >> seqName;
896                         it = badSeqNames.find(seqName);
897                         string line;            
898                         while (!inputAlignReport.eof()) {               //      need to copy header
899                                 char c = inputAlignReport.get();
900                                 line += c;
901                                 if (c == 10 || c == 13){        break;  }       
902                         }
903                         
904                         if(it != badSeqNames.end()){
905                                 badSeqNames.erase(it);
906                         }
907                         else{
908                                 goodAlignReportOut << seqName << '\t' << line;
909                         }
910                         m->gobble(inputAlignReport);
911                 }
912                 
913                 if (m->control_pressed) { goodAlignReportOut.close();  inputAlignReport.close(); remove(goodAlignReportFile.c_str());  return 0; }
914
915                 //we were unable to remove some of the bad sequences
916                 if (badSeqNames.size() != 0) {
917                         for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {  
918                                 m->mothurOut("Your alignreport file does not include the sequence " + *it + " please correct."); 
919                                 m->mothurOutEndLine();
920                         }
921                 }
922
923                 inputAlignReport.close();
924                 goodAlignReportOut.close();
925                                 
926                 if (m->control_pressed) {  remove(goodAlignReportFile.c_str());  return 0; }
927                 
928                 return 0;
929         
930         }
931         catch(exception& e) {
932                 m->errorOut(e, "ScreenSeqsCommand", "screenAlignReport");
933                 exit(1);
934         }
935         
936 }
937 //***************************************************************************************************************
938
939 int ScreenSeqsCommand::screenQual(set<string> badSeqNames){
940         try {
941                 ifstream in;
942                 m->openInputFile(qualfile, in);
943                 set<string>::iterator it;
944                 
945                 string goodQualFile = outputDir + m->getRootName(m->getSimpleName(qualfile)) + "good" + m->getExtension(qualfile);
946                 outputNames.push_back(goodQualFile);  outputTypes["qfile"].push_back(goodQualFile);
947                 ofstream goodQual;      m->openOutputFile(goodQualFile, goodQual);
948                 
949                 while(!in.eof()){       
950                         
951                         if (m->control_pressed) { goodQual.close(); in.close(); remove(goodQualFile.c_str()); return 0; }
952
953                         string saveName = "";
954                         string name = "";
955                         string scores = "";
956                         
957                         in >> name; 
958                         
959                         if (name.length() != 0) { 
960                                 saveName = name.substr(1);
961                                 while (!in.eof())       {       
962                                         char c = in.get(); 
963                                         if (c == 10 || c == 13){        break;  }
964                                         else { name += c; }     
965                                 } 
966                                 m->gobble(in);
967                         }
968                         
969                         while(in){
970                                 char letter= in.get();
971                                 if(letter == '>'){      in.putback(letter);     break;  }
972                                 else{ scores += letter; }
973                         }
974                         
975                         m->gobble(in);
976                         
977                         it = badSeqNames.find(saveName);
978                         
979                         if(it != badSeqNames.end()){
980                                 badSeqNames.erase(it);
981                         }else{                          
982                                 goodQual << name << endl << scores;
983                         }
984                         
985                         m->gobble(in);
986                 }
987                 
988                 in.close();
989                 goodQual.close();
990                 
991                 //we were unable to remove some of the bad sequences
992                 if (badSeqNames.size() != 0) {
993                         for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {  
994                                 m->mothurOut("Your qual file does not include the sequence " + *it + " please correct."); 
995                                 m->mothurOutEndLine();
996                         }
997                 }
998                 
999                 if (m->control_pressed) {  remove(goodQualFile.c_str());  return 0; }
1000                 
1001                 return 0;
1002                 
1003         }
1004         catch(exception& e) {
1005                 m->errorOut(e, "ScreenSeqsCommand", "screenQual");
1006                 exit(1);
1007         }
1008         
1009 }
1010 //**********************************************************************************************************************
1011
1012 int ScreenSeqsCommand::driver(linePair* filePos, string goodFName, string badAccnosFName, string filename, set<string>& badSeqNames){
1013         try {
1014                 ofstream goodFile;
1015                 m->openOutputFile(goodFName, goodFile);
1016                 
1017                 ofstream badAccnosFile;
1018                 m->openOutputFile(badAccnosFName, badAccnosFile);
1019                 
1020                 ifstream inFASTA;
1021                 m->openInputFile(filename, inFASTA);
1022
1023                 inFASTA.seekg(filePos->start);
1024
1025                 bool done = false;
1026                 int count = 0;
1027         
1028                 while (!done) {
1029                 
1030                         if (m->control_pressed) {  return 0; }
1031                         
1032                         Sequence currSeq(inFASTA); m->gobble(inFASTA);
1033                         if (currSeq.getName() != "") {
1034                                 bool goodSeq = 1;               //      innocent until proven guilty
1035                                 if(goodSeq == 1 && startPos != -1 && startPos < currSeq.getStartPos())                  {       goodSeq = 0;    }
1036                                 if(goodSeq == 1 && endPos != -1 && endPos > currSeq.getEndPos())                                {       goodSeq = 0;    }
1037                                 if(goodSeq == 1 && maxAmbig != -1 && maxAmbig < currSeq.getAmbigBases())                {       goodSeq = 0;    }
1038                                 if(goodSeq == 1 && maxHomoP != -1 && maxHomoP < currSeq.getLongHomoPolymer())   {       goodSeq = 0;    }
1039                                 if(goodSeq == 1 && minLength != -1 && minLength > currSeq.getNumBases())                {       goodSeq = 0;    }
1040                                 if(goodSeq == 1 && maxLength != -1 && maxLength < currSeq.getNumBases())                {       goodSeq = 0;    }
1041                                 
1042                                 if(goodSeq == 1){
1043                                         currSeq.printSequence(goodFile);        
1044                                 }
1045                                 else{
1046                                         badAccnosFile << currSeq.getName() << endl;
1047                                         badSeqNames.insert(currSeq.getName());
1048                                 }
1049                         count++;
1050                         }
1051                         
1052                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1053                                 unsigned long int pos = inFASTA.tellg();
1054                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
1055                         #else
1056                                 if (inFASTA.eof()) { break; }
1057                         #endif
1058                         
1059                         //report progress
1060                         if((count) % 100 == 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
1061                 }
1062                 //report progress
1063                 if((count) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
1064                 
1065                         
1066                 goodFile.close();
1067                 inFASTA.close();
1068                 badAccnosFile.close();
1069                 
1070                 return count;
1071         }
1072         catch(exception& e) {
1073                 m->errorOut(e, "ScreenSeqsCommand", "driver");
1074                 exit(1);
1075         }
1076 }
1077 //**********************************************************************************************************************
1078 #ifdef USE_MPI
1079 int ScreenSeqsCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& goodFile, MPI_File& badAccnosFile, vector<unsigned long int>& MPIPos, set<string>& badSeqNames){
1080         try {
1081                 string outputString = "";
1082                 MPI_Status statusGood; 
1083                 MPI_Status statusBadAccnos; 
1084                 MPI_Status status; 
1085                 int pid;
1086                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
1087
1088                 for(int i=0;i<num;i++){
1089                 
1090                         if (m->control_pressed) {  return 0; }
1091                         
1092                         //read next sequence
1093                         int length = MPIPos[start+i+1] - MPIPos[start+i];
1094
1095                         char* buf4 = new char[length];
1096                         memcpy(buf4, outputString.c_str(), length);
1097
1098                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
1099                         
1100                         string tempBuf = buf4;  delete buf4;
1101                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
1102                         istringstream iss (tempBuf,istringstream::in);
1103                         
1104                         Sequence currSeq(iss);                  
1105                         
1106                         //process seq
1107                         if (currSeq.getName() != "") {
1108                                 bool goodSeq = 1;               //      innocent until proven guilty
1109                                 if(goodSeq == 1 && startPos != -1 && startPos < currSeq.getStartPos())                  {       goodSeq = 0;    }
1110                                 if(goodSeq == 1 && endPos != -1 && endPos > currSeq.getEndPos())                                {       goodSeq = 0;    }
1111                                 if(goodSeq == 1 && maxAmbig != -1 && maxAmbig < currSeq.getAmbigBases())                {       goodSeq = 0;    }
1112                                 if(goodSeq == 1 && maxHomoP != -1 && maxHomoP < currSeq.getLongHomoPolymer())   {       goodSeq = 0;    }
1113                                 if(goodSeq == 1 && minLength != -1 && minLength > currSeq.getNumBases())                {       goodSeq = 0;    }
1114                                 if(goodSeq == 1 && maxLength != -1 && maxLength < currSeq.getNumBases())                {       goodSeq = 0;    }
1115                                 
1116                                 if(goodSeq == 1){
1117                                         outputString =  ">" + currSeq.getName() + "\n" + currSeq.getAligned() + "\n";
1118                                 
1119                                         //print good seq
1120                                         length = outputString.length();
1121                                         char* buf2 = new char[length];
1122                                         memcpy(buf2, outputString.c_str(), length);
1123                                         
1124                                         MPI_File_write_shared(goodFile, buf2, length, MPI_CHAR, &statusGood);
1125                                         delete buf2;
1126                                 }
1127                                 else{
1128
1129                                         badSeqNames.insert(currSeq.getName());
1130                                         
1131                                         //write to bad accnos file
1132                                         outputString = currSeq.getName() + "\n";
1133                                 
1134                                         length = outputString.length();
1135                                         char* buf3 = new char[length];
1136                                         memcpy(buf3, outputString.c_str(), length);
1137                                         
1138                                         MPI_File_write_shared(badAccnosFile, buf3, length, MPI_CHAR, &statusBadAccnos);
1139                                         delete buf3;
1140                                 }
1141                         }
1142                         
1143                         //report progress
1144                         if((i) % 100 == 0){     m->mothurOut("Processing sequence: " + toString(i)); m->mothurOutEndLine();             }
1145                 }
1146                                 
1147                 return 1;
1148         }
1149         catch(exception& e) {
1150                 m->errorOut(e, "ScreenSeqsCommand", "driverMPI");
1151                 exit(1);
1152         }
1153 }
1154 #endif
1155 /**************************************************************************************************/
1156
1157 int ScreenSeqsCommand::createProcesses(string goodFileName, string badAccnos, string filename, set<string>& badSeqNames) {
1158         try {
1159 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1160                 int process = 0;
1161                 int num = 0;
1162                 
1163                 //loop through and create all the processes you want
1164                 while (process != processors) {
1165                         int pid = fork();
1166                         
1167                         if (pid > 0) {
1168                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
1169                                 process++;
1170                         }else if (pid == 0){
1171                                 num = driver(lines[process], goodFileName + toString(getpid()) + ".temp", badAccnos + toString(getpid()) + ".temp", filename, badSeqNames);
1172                                 
1173                                 //pass numSeqs to parent
1174                                 ofstream out;
1175                                 string tempFile = filename + toString(getpid()) + ".num.temp";
1176                                 m->openOutputFile(tempFile, out);
1177                                 out << num << endl;
1178                                 out.close();
1179                                 
1180                                 exit(0);
1181                         }else { 
1182                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
1183                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
1184                                 exit(0);
1185                         }
1186                 }
1187                 
1188                 //force parent to wait until all the processes are done
1189                 for (int i=0;i<processors;i++) { 
1190                         int temp = processIDS[i];
1191                         wait(&temp);
1192                 }
1193                 
1194                 for (int i = 0; i < processIDS.size(); i++) {
1195                         ifstream in;
1196                         string tempFile =  filename + toString(processIDS[i]) + ".num.temp";
1197                         m->openInputFile(tempFile, in);
1198                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
1199                         in.close(); remove(tempFile.c_str());
1200                 }
1201                 
1202                 return num;
1203 #endif          
1204         }
1205         catch(exception& e) {
1206                 m->errorOut(e, "ScreenSeqsCommand", "createProcesses");
1207                 exit(1);
1208         }
1209 }
1210
1211 //***************************************************************************************************************
1212
1213