]> git.donarmstrong.com Git - mothur.git/blob - screenseqscommand.cpp
numSeqs passing fix
[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
15 ScreenSeqsCommand::ScreenSeqsCommand(string option)  {
16         try {
17                 abort = false;
18                 
19                 //allow user to run help
20                 if(option == "help") { help(); abort = true; }
21                 
22                 else {
23                         //valid paramters for this command
24                         string AlignArray[] =  {"fasta", "start", "end", "maxambig", "maxhomop", "minlength", "maxlength",
25                                                                         "name", "group", "alignreport","processors","outputdir","inputdir"};
26                         vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
27                         
28                         OptionParser parser(option);
29                         map<string,string> parameters = parser.getParameters();
30                         
31                         ValidParameters validParameter("screen.seqs");
32                         map<string,string>::iterator it;
33                         
34                         //check to make sure all parameters are valid for command
35                         for (it = parameters.begin(); it != parameters.end(); it++) { 
36                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
37                         }
38                         
39                         //if the user changes the input directory command factory will send this info to us in the output parameter 
40                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
41                         if (inputDir == "not found"){   inputDir = "";          }
42                         else {
43                                 string path;
44                                 it = parameters.find("fasta");
45                                 //user has given a template file
46                                 if(it != parameters.end()){ 
47                                         path = hasPath(it->second);
48                                         //if the user has not given a path then, add inputdir. else leave path alone.
49                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
50                                 }
51                                 
52                                 it = parameters.find("group");
53                                 //user has given a template file
54                                 if(it != parameters.end()){ 
55                                         path = hasPath(it->second);
56                                         //if the user has not given a path then, add inputdir. else leave path alone.
57                                         if (path == "") {       parameters["group"] = inputDir + it->second;            }
58                                 }
59                                 
60                                 it = parameters.find("name");
61                                 //user has given a template file
62                                 if(it != parameters.end()){ 
63                                         path = hasPath(it->second);
64                                         //if the user has not given a path then, add inputdir. else leave path alone.
65                                         if (path == "") {       parameters["name"] = inputDir + it->second;             }
66                                 }
67                                 
68                                 it = parameters.find("alignreport");
69                                 //user has given a template file
70                                 if(it != parameters.end()){ 
71                                         path = hasPath(it->second);
72                                         //if the user has not given a path then, add inputdir. else leave path alone.
73                                         if (path == "") {       parameters["alignreport"] = inputDir + it->second;              }
74                                 }
75                         }
76
77                         //check for required parameters
78                         fastafile = validParameter.validFile(parameters, "fasta", true);
79                         if (fastafile == "not found") { m->mothurOut("fasta is a required parameter for the screen.seqs command."); m->mothurOutEndLine(); abort = true; }
80                         else if (fastafile == "not open") { abort = true; }     
81         
82                         groupfile = validParameter.validFile(parameters, "group", true);
83                         if (groupfile == "not open") { abort = true; }  
84                         else if (groupfile == "not found") { groupfile = ""; }
85                         
86                         namefile = validParameter.validFile(parameters, "name", true);
87                         if (namefile == "not open") { abort = true; }
88                         else if (namefile == "not found") { namefile = ""; }    
89
90                         alignreport = validParameter.validFile(parameters, "alignreport", true);
91                         if (alignreport == "not open") { abort = true; }
92                         else if (alignreport == "not found") { alignreport = ""; }      
93                         
94                         //if the user changes the output directory command factory will send this info to us in the output parameter 
95                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
96                                 outputDir = ""; 
97                                 outputDir += hasPath(fastafile); //if user entered a file with a path then preserve it  
98                         }
99
100                         //check for optional parameter and set defaults
101                         // ...at some point should added some additional type checking...
102                         string temp;
103                         temp = validParameter.validFile(parameters, "start", false);            if (temp == "not found") { temp = "-1"; }
104                         convert(temp, startPos); 
105                 
106                         temp = validParameter.validFile(parameters, "end", false);                      if (temp == "not found") { temp = "-1"; }
107                         convert(temp, endPos);  
108
109                         temp = validParameter.validFile(parameters, "maxambig", false);         if (temp == "not found") { temp = "-1"; }
110                         convert(temp, maxAmbig);  
111
112                         temp = validParameter.validFile(parameters, "maxhomop", false);         if (temp == "not found") { temp = "-1"; }
113                         convert(temp, maxHomoP);  
114
115                         temp = validParameter.validFile(parameters, "minlength", false);        if (temp == "not found") { temp = "-1"; }
116                         convert(temp, minLength); 
117                         
118                         temp = validParameter.validFile(parameters, "maxlength", false);        if (temp == "not found") { temp = "-1"; }
119                         convert(temp, maxLength); 
120                         
121                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = "1";                             }
122                         convert(temp, processors); 
123
124                 }
125
126         }
127         catch(exception& e) {
128                 m->errorOut(e, "ScreenSeqsCommand", "ScreenSeqsCommand");
129                 exit(1);
130         }
131 }
132 //**********************************************************************************************************************
133
134 void ScreenSeqsCommand::help(){
135         try {
136                 m->mothurOut("The screen.seqs command reads a fastafile and creates .....\n");
137                 m->mothurOut("The screen.seqs command parameters are fasta, start, end, maxambig, maxhomop, minlength, maxlength, name, group and processors.\n");
138                 m->mothurOut("The fasta parameter is required.\n");
139                 m->mothurOut("The start parameter .... The default is -1.\n");
140                 m->mothurOut("The end parameter .... The default is -1.\n");
141                 m->mothurOut("The maxambig parameter .... The default is -1.\n");
142                 m->mothurOut("The maxhomop parameter .... The default is -1.\n");
143                 m->mothurOut("The minlength parameter .... The default is -1.\n");
144                 m->mothurOut("The maxlength parameter .... The default is -1.\n");
145                 m->mothurOut("The processors parameter allows you to specify the number of processors to use while running the command. The default is 1.\n");
146                 m->mothurOut("The name parameter allows you to provide a namesfile, and the group parameter allows you to provide a groupfile.\n");
147                 m->mothurOut("The screen.seqs command should be in the following format: \n");
148                 m->mothurOut("screen.seqs(fasta=yourFastaFile, name=youNameFile, group=yourGroupFIle, start=yourStart, end=yourEnd, maxambig=yourMaxambig,  \n");
149                 m->mothurOut("maxhomop=yourMaxhomop, minlength=youMinlength, maxlength=yourMaxlength)  \n");    
150                 m->mothurOut("Example screen.seqs(fasta=abrecovery.fasta, name=abrecovery.names, group=abrecovery.groups, start=..., end=..., maxambig=..., maxhomop=..., minlength=..., maxlength=...).\n");
151                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFasta).\n\n");
152
153         }
154         catch(exception& e) {
155                 m->errorOut(e, "ScreenSeqsCommand", "help");
156                 exit(1);
157         }
158 }
159
160 //***************************************************************************************************************
161
162 ScreenSeqsCommand::~ScreenSeqsCommand(){        /*      do nothing      */      }
163
164 //***************************************************************************************************************
165
166 int ScreenSeqsCommand::execute(){
167         try{
168                 
169                 if (abort == true) { return 0; }
170                                 
171                 string goodSeqFile = outputDir + getRootName(getSimpleName(fastafile)) + "good" + getExtension(fastafile);
172                 string badSeqFile =  outputDir + getRootName(getSimpleName(fastafile)) + "bad" + getExtension(fastafile);
173                 string badAccnosFile =  outputDir + getRootName(getSimpleName(fastafile)) + "bad.accnos";
174                 
175                 int numFastaSeqs = 0;
176                 set<string> badSeqNames;
177                 int start = time(NULL);
178                 
179 #ifdef USE_MPI  
180                         int pid, end, numSeqsPerProcessor; 
181                         int tag = 2001;
182                         vector<unsigned long int> MPIPos;
183                         
184                         MPI_Status status; 
185                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
186                         MPI_Comm_size(MPI_COMM_WORLD, &processors); 
187
188                         MPI_File inMPI;
189                         MPI_File outMPIGood;
190                         MPI_File outMPIBad;
191                         MPI_File outMPIBadAccnos;
192                         
193                         int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
194                         int inMode=MPI_MODE_RDONLY; 
195                         
196                         char outGoodFilename[1024];
197                         strcpy(outGoodFilename, goodSeqFile.c_str());
198
199                         char outBadFilename[1024];
200                         strcpy(outBadFilename, badSeqFile.c_str());
201                         
202                         char outBadAccnosFilename[1024];
203                         strcpy(outBadAccnosFilename, badAccnosFile.c_str());
204
205                         char inFileName[1024];
206                         strcpy(inFileName, fastafile.c_str());
207                         
208                         MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
209                         MPI_File_open(MPI_COMM_WORLD, outGoodFilename, outMode, MPI_INFO_NULL, &outMPIGood);
210                         MPI_File_open(MPI_COMM_WORLD, outBadFilename, outMode, MPI_INFO_NULL, &outMPIBad);
211                         MPI_File_open(MPI_COMM_WORLD, outBadAccnosFilename, outMode, MPI_INFO_NULL, &outMPIBadAccnos);
212                         
213                         if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIGood);  MPI_File_close(&outMPIBad); MPI_File_close(&outMPIBadAccnos); return 0; }
214                         
215                         if (pid == 0) { //you are the root process 
216                                 
217                                 MPIPos = setFilePosFasta(fastafile, numFastaSeqs); //fills MPIPos, returns numSeqs
218                                 
219                                 //send file positions to all processes
220                                 for(int i = 1; i < processors; i++) { 
221                                         MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
222                                         MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
223                                 }
224                                 
225                                 //figure out how many sequences you have to align
226                                 numSeqsPerProcessor = numFastaSeqs / processors;
227                                 int startIndex =  pid * numSeqsPerProcessor;
228                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
229                                 
230                                 //align your part
231                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIGood, outMPIBad, outMPIBadAccnos, MPIPos, badSeqNames);
232
233                                 if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIGood);  MPI_File_close(&outMPIBadAccnos); MPI_File_close(&outMPIBad);  return 0; }
234
235                                 for (int i = 1; i < processors; i++) {
236                                 
237                                         //get bad lists
238                                         int badSize;
239                                         MPI_Recv(&badSize, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
240                                         /*for (int j = 0; j < badSize; j++) {
241                                                 int length;
242                                                 MPI_Recv(&length, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);  //recv the length of the name
243                                                 char* buf2 = new char[length];                                                                          //make space to recieve it
244                                                 MPI_Recv(buf2, length, MPI_CHAR, i, tag, MPI_COMM_WORLD, &status);  //get name
245                                                 
246                                                 string tempBuf = buf2;
247                                                 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
248                                                 delete buf2;
249                                                 
250                                                 badSeqNames.insert(tempBuf);
251                                         }*/
252                                 }
253                         }else{ //you are a child process
254                                 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
255                                 MPIPos.resize(numFastaSeqs+1);
256                                 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
257
258                                 //figure out how many sequences you have to align
259                                 numSeqsPerProcessor = numFastaSeqs / processors;
260                                 int startIndex =  pid * numSeqsPerProcessor;
261                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
262                                 
263                                 //align your part
264                                 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIGood, outMPIBad, outMPIBadAccnos, MPIPos, badSeqNames);
265
266                                 if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIGood);  MPI_File_close(&outMPIBad); MPI_File_close(&outMPIBadAccnos); return 0; }
267                                 
268                                 //send bad list 
269                                 int badSize = badSeqNames.size();
270                                 MPI_Send(&badSize, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
271                                 
272                                 /*
273                                 set<string>::iterator it;
274                                 for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {
275                                         string name = *it;
276                                         int length = name.length();
277                                         char* buf2 = new char[length];
278                                         memcpy(buf2, name.c_str(), length);
279                                         
280                                         MPI_Send(&length, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
281                                         MPI_Send(buf2, length, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
282                                 }*/
283                         }
284                         
285                         //close files 
286                         MPI_File_close(&inMPI);
287                         MPI_File_close(&outMPIGood);
288                         MPI_File_close(&outMPIBad);
289                         MPI_File_close(&outMPIBadAccnos);
290                         MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
291                                         
292 #else
293                         vector<unsigned long int> positions = divideFile(fastafile, processors);
294                                 
295                         for (int i = 0; i < (positions.size()-1); i++) {
296                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
297                         }       
298                                                 
299         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
300                         if(processors == 1){
301                                 numFastaSeqs = driver(lines[0], goodSeqFile, badSeqFile, badAccnosFile, fastafile, badSeqNames);
302                                 
303                                 if (m->control_pressed) { remove(goodSeqFile.c_str()); remove(badSeqFile.c_str()); return 0; }
304                                 
305                         }else{
306                                 processIDS.resize(0);
307                                 
308                                 numFastaSeqs = createProcesses(goodSeqFile, badSeqFile, badAccnosFile, fastafile, badSeqNames); 
309                                 
310                                 rename((goodSeqFile + toString(processIDS[0]) + ".temp").c_str(), goodSeqFile.c_str());
311                                 rename((badSeqFile + toString(processIDS[0]) + ".temp").c_str(), badSeqFile.c_str());
312                                 rename((badAccnosFile + toString(processIDS[0]) + ".temp").c_str(), badAccnosFile.c_str());
313                                 
314                                 //append alignment and report files
315                                 for(int i=1;i<processors;i++){
316                                         appendFiles((goodSeqFile + toString(processIDS[i]) + ".temp"), goodSeqFile);
317                                         remove((goodSeqFile + toString(processIDS[i]) + ".temp").c_str());
318                                         
319                                         appendFiles((badSeqFile + toString(processIDS[i]) + ".temp"), badSeqFile);
320                                         remove((badSeqFile + toString(processIDS[i]) + ".temp").c_str());
321                                         
322                                         appendFiles((badAccnosFile + toString(processIDS[i]) + ".temp"), badAccnosFile);
323                                         remove((badAccnosFile + toString(processIDS[i]) + ".temp").c_str());
324                                 }
325                                 
326                                 if (m->control_pressed) { remove(goodSeqFile.c_str()); remove(badSeqFile.c_str()); return 0; }
327                                 
328                                 //read badSeqs in because root process doesnt know what other "bad" seqs the children found
329                                 ifstream inBad;
330                                 int ableToOpen = openInputFile(badAccnosFile, inBad, "no error");
331                                 
332                                 if (ableToOpen == 0) {
333                                         badSeqNames.clear();
334                                         string tempName;
335                                         while (!inBad.eof()) {
336                                                 inBad >> tempName; gobble(inBad);
337                                                 badSeqNames.insert(tempName);
338                                         }
339                                         inBad.close();
340                                 }
341                         }
342         #else
343                         numFastaSeqs = driver(lines[0], goodSeqFile, badSeqFile, badAccnosFile, fastafile, badSeqNames);
344                         
345                         if (m->control_pressed) { remove(goodSeqFile.c_str()); remove(badSeqFile.c_str()); return 0; }
346                         
347         #endif
348
349 #endif          
350
351                 #ifdef USE_MPI
352                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
353                                         
354                         if (pid == 0) { //only one process should fix files
355                         
356                                 //read accnos file with all names in it, process 0 just has its names
357                                 MPI_File inMPIAccnos;
358                                 MPI_Offset size;
359                         
360                                 char inFileName[1024];
361                                 strcpy(inFileName, badAccnosFile.c_str());
362                         
363                                 MPI_File_open(MPI_COMM_SELF, inFileName, inMode, MPI_INFO_NULL, &inMPIAccnos);  //comm, filename, mode, info, filepointer
364                                 MPI_File_get_size(inMPIAccnos, &size);
365                         
366                                 char* buffer = new char[size];
367                                 MPI_File_read(inMPIAccnos, buffer, size, MPI_CHAR, &status);
368                         
369                                 string tempBuf = buffer;
370                                 if (tempBuf.length() > size) { tempBuf = tempBuf.substr(0, size);  }
371                                 istringstream iss (tempBuf,istringstream::in);
372
373                                 delete buffer;
374                                 MPI_File_close(&inMPIAccnos);
375                                 
376                                 badSeqNames.clear();
377                                 string tempName;
378                                 while (!iss.eof()) {
379                                         iss >> tempName; gobble(iss);
380                                         badSeqNames.insert(tempName);
381                                 }
382                 #endif
383                                                                                                                                                                         
384                 if(namefile != "" && groupfile != "")   {       
385                         screenNameGroupFile(badSeqNames);       
386                         if (m->control_pressed) {  remove(goodSeqFile.c_str()); remove(badSeqFile.c_str()); return 0; }
387                 }else if(namefile != "")        {       
388                         screenNameGroupFile(badSeqNames);
389                         if (m->control_pressed) {  remove(goodSeqFile.c_str()); remove(badSeqFile.c_str()); return 0; } 
390                 }else if(groupfile != "")                               {       screenGroupFile(badSeqNames);           }       // this screens just the group
391                 
392                 if (m->control_pressed) { remove(goodSeqFile.c_str()); remove(badSeqFile.c_str()); return 0; }
393
394                 if(alignreport != "")                                   {       screenAlignReport(badSeqNames);         }
395                 
396                 if (m->control_pressed) { remove(goodSeqFile.c_str()); remove(badSeqFile.c_str()); return 0; }
397                 
398                 #ifdef USE_MPI
399                         }
400                 #endif
401
402                 m->mothurOutEndLine();
403                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
404                 m->mothurOut(goodSeqFile); m->mothurOutEndLine();       
405                 m->mothurOut(badSeqFile); m->mothurOutEndLine();        
406                 m->mothurOut(badAccnosFile); m->mothurOutEndLine();     
407                 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
408                 m->mothurOutEndLine();
409                 m->mothurOutEndLine();
410
411                 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to screen " + toString(numFastaSeqs) + " sequences.");
412                 m->mothurOutEndLine();
413
414                 return 0;
415         }
416         catch(exception& e) {
417                 m->errorOut(e, "ScreenSeqsCommand", "execute");
418                 exit(1);
419         }
420 }
421
422 //***************************************************************************************************************
423
424 int ScreenSeqsCommand::screenNameGroupFile(set<string> badSeqNames){
425         try {
426                 ifstream inputNames;
427                 openInputFile(namefile, inputNames);
428                 set<string> badSeqGroups;
429                 string seqName, seqList, group;
430                 set<string>::iterator it;
431
432                 string goodNameFile = outputDir + getRootName(getSimpleName(namefile)) + "good" + getExtension(namefile);
433                 string badNameFile = outputDir + getRootName(getSimpleName(namefile)) + "bad" + getExtension(namefile);
434                 
435                 outputNames.push_back(goodNameFile);  outputNames.push_back(badNameFile);
436                 
437                 ofstream goodNameOut;   openOutputFile(goodNameFile, goodNameOut);
438                 ofstream badNameOut;    openOutputFile(badNameFile, badNameOut);                
439                 
440                 while(!inputNames.eof()){
441                         if (m->control_pressed) { goodNameOut.close(); badNameOut.close(); inputNames.close(); remove(goodNameFile.c_str()); remove(badNameFile.c_str()); return 0; }
442
443                         inputNames >> seqName >> seqList;
444                         it = badSeqNames.find(seqName);
445                         
446                         if(it != badSeqNames.end()){
447                                 badSeqNames.erase(it);
448                                 badNameOut << seqName << '\t' << seqList << endl;
449                                 if(namefile != ""){
450                                         int start = 0;
451                                         for(int i=0;i<seqList.length();i++){
452                                                 if(seqList[i] == ','){
453                                                         badSeqGroups.insert(seqList.substr(start,i-start));
454                                                         start = i+1;
455                                                 }                                       
456                                         }
457                                         badSeqGroups.insert(seqList.substr(start,seqList.length()-start));
458                                 }
459                         }
460                         else{
461                                 goodNameOut << seqName << '\t' << seqList << endl;
462                         }
463                         gobble(inputNames);
464                 }
465                 inputNames.close();
466                 goodNameOut.close();
467                 badNameOut.close();
468                 
469                 //we were unable to remove some of the bad sequences
470                 if (badSeqNames.size() != 0) {
471                         for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {  
472                                 m->mothurOut("Your namefile does not include the sequence " + *it + " please correct."); 
473                                 m->mothurOutEndLine();
474                         }
475                 }
476
477                 if(groupfile != ""){
478                         
479                         ifstream inputGroups;
480                         openInputFile(groupfile, inputGroups);
481
482                         string goodGroupFile = outputDir + getRootName(getSimpleName(groupfile)) + "good" + getExtension(groupfile);
483                         string badGroupFile = outputDir + getRootName(getSimpleName(groupfile)) + "bad" + getExtension(groupfile);
484                         
485                         outputNames.push_back(goodGroupFile);  outputNames.push_back(badGroupFile);
486                         
487                         ofstream goodGroupOut;  openOutputFile(goodGroupFile, goodGroupOut);
488                         ofstream badGroupOut;   openOutputFile(badGroupFile, badGroupOut);              
489                         
490                         while(!inputGroups.eof()){
491                                 if (m->control_pressed) { goodGroupOut.close(); badGroupOut.close(); inputGroups.close(); remove(goodNameFile.c_str()); remove(badNameFile.c_str()); remove(goodGroupFile.c_str()); remove(badGroupFile.c_str()); return 0; }
492
493                                 inputGroups >> seqName >> group;
494
495                                 it = badSeqGroups.find(seqName);
496                                 
497                                 if(it != badSeqGroups.end()){
498                                         badSeqGroups.erase(it);
499                                         badGroupOut << seqName << '\t' << group << endl;
500                                 }
501                                 else{
502                                         goodGroupOut << seqName << '\t' << group << endl;
503                                 }
504                                 gobble(inputGroups);
505                         }
506                         inputGroups.close();
507                         goodGroupOut.close();
508                         badGroupOut.close();
509                         
510                         //we were unable to remove some of the bad sequences
511                         if (badSeqGroups.size() != 0) {
512                                 for (it = badSeqGroups.begin(); it != badSeqGroups.end(); it++) {  
513                                         m->mothurOut("Your namefile does not include the sequence " + *it + " please correct."); 
514                                         m->mothurOutEndLine();
515                                 }
516                         }
517                 }
518                         
519                 return 0;
520         
521         }
522         catch(exception& e) {
523                 m->errorOut(e, "ScreenSeqsCommand", "screenNameGroupFile");
524                 exit(1);
525         }
526 }
527
528 //***************************************************************************************************************
529
530 int ScreenSeqsCommand::screenGroupFile(set<string> badSeqNames){
531         try {
532                 ifstream inputGroups;
533                 openInputFile(groupfile, inputGroups);
534                 string seqName, group;
535                 set<string>::iterator it;
536                 
537                 string goodGroupFile = outputDir + getRootName(getSimpleName(groupfile)) + "good" + getExtension(groupfile);
538                 string badGroupFile = outputDir + getRootName(getSimpleName(groupfile)) + "bad" + getExtension(groupfile);
539                 
540                 outputNames.push_back(goodGroupFile);  outputNames.push_back(badGroupFile);
541                 
542                 ofstream goodGroupOut;  openOutputFile(goodGroupFile, goodGroupOut);
543                 ofstream badGroupOut;   openOutputFile(badGroupFile, badGroupOut);              
544                 
545                 while(!inputGroups.eof()){
546                         if (m->control_pressed) { goodGroupOut.close(); badGroupOut.close(); inputGroups.close(); remove(goodGroupFile.c_str()); remove(badGroupFile.c_str()); return 0; }
547
548                         inputGroups >> seqName >> group;
549                         it = badSeqNames.find(seqName);
550                         
551                         if(it != badSeqNames.end()){
552                                 badSeqNames.erase(it);
553                                 badGroupOut << seqName << '\t' << group << endl;
554                         }
555                         else{
556                                 goodGroupOut << seqName << '\t' << group << endl;
557                         }
558                         gobble(inputGroups);
559                 }
560                 
561                 if (m->control_pressed) { goodGroupOut.close(); badGroupOut.close(); inputGroups.close(); remove(goodGroupFile.c_str()); remove(badGroupFile.c_str()); return 0; }
562
563                 //we were unable to remove some of the bad sequences
564                 if (badSeqNames.size() != 0) {
565                         for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {  
566                                 m->mothurOut("Your groupfile does not include the sequence " + *it + " please correct."); 
567                                 m->mothurOutEndLine();
568                         }
569                 }
570                 
571                 inputGroups.close();
572                 goodGroupOut.close();
573                 badGroupOut.close();
574                 
575                 if (m->control_pressed) { remove(goodGroupFile.c_str()); remove(badGroupFile.c_str());  }
576
577                 
578                 return 0;
579         
580         }
581         catch(exception& e) {
582                 m->errorOut(e, "ScreenSeqsCommand", "screenGroupFile");
583                 exit(1);
584         }
585 }
586
587 //***************************************************************************************************************
588
589 int ScreenSeqsCommand::screenAlignReport(set<string> badSeqNames){
590         try {
591                 ifstream inputAlignReport;
592                 openInputFile(alignreport, inputAlignReport);
593                 string seqName, group;
594                 set<string>::iterator it;
595                 
596                 string goodAlignReportFile = outputDir + getRootName(getSimpleName(alignreport)) + "good" + getExtension(alignreport);
597                 string badAlignReportFile = outputDir + getRootName(getSimpleName(alignreport)) + "bad" + getExtension(alignreport);
598                 
599                 outputNames.push_back(goodAlignReportFile);  outputNames.push_back(badAlignReportFile);
600                 
601                 ofstream goodAlignReportOut;    openOutputFile(goodAlignReportFile, goodAlignReportOut);
602                 ofstream badAlignReportOut;             openOutputFile(badAlignReportFile, badAlignReportOut);          
603
604                 while (!inputAlignReport.eof()) {               //      need to copy header
605                         char c = inputAlignReport.get();
606                         goodAlignReportOut << c;
607                         badAlignReportOut << c;
608                         if (c == 10 || c == 13){        break;  }       
609                 }
610
611                 while(!inputAlignReport.eof()){
612                         if (m->control_pressed) { goodAlignReportOut.close(); badAlignReportOut.close(); inputAlignReport.close(); remove(goodAlignReportFile.c_str()); remove(badAlignReportFile.c_str()); return 0; }
613
614                         inputAlignReport >> seqName;
615                         it = badSeqNames.find(seqName);
616                         string line;            
617                         while (!inputAlignReport.eof()) {               //      need to copy header
618                                 char c = inputAlignReport.get();
619                                 line += c;
620                                 if (c == 10 || c == 13){        break;  }       
621                         }
622                         
623                         if(it != badSeqNames.end()){
624                                 badSeqNames.erase(it);
625                                 badAlignReportOut << seqName << '\t' << line;
626                         }
627                         else{
628                                 goodAlignReportOut << seqName << '\t' << line;
629                         }
630                         gobble(inputAlignReport);
631                 }
632                 
633                 if (m->control_pressed) { goodAlignReportOut.close(); badAlignReportOut.close(); inputAlignReport.close(); remove(goodAlignReportFile.c_str()); remove(badAlignReportFile.c_str()); return 0; }
634
635                 //we were unable to remove some of the bad sequences
636                 if (badSeqNames.size() != 0) {
637                         for (it = badSeqNames.begin(); it != badSeqNames.end(); it++) {  
638                                 m->mothurOut("Your file does not include the sequence " + *it + " please correct."); 
639                                 m->mothurOutEndLine();
640                         }
641                 }
642
643                 inputAlignReport.close();
644                 goodAlignReportOut.close();
645                 badAlignReportOut.close();
646                                 
647                 if (m->control_pressed) {  remove(goodAlignReportFile.c_str()); remove(badAlignReportFile.c_str()); return 0; }
648                 
649                 return 0;
650         
651         }
652         catch(exception& e) {
653                 m->errorOut(e, "ScreenSeqsCommand", "screenAlignReport");
654                 exit(1);
655         }
656         
657 }
658 //**********************************************************************************************************************
659
660 int ScreenSeqsCommand::driver(linePair* filePos, string goodFName, string badFName, string badAccnosFName, string filename, set<string>& badSeqNames){
661         try {
662                 ofstream goodFile;
663                 openOutputFile(goodFName, goodFile);
664                 
665                 ofstream badFile;
666                 openOutputFile(badFName, badFile);
667                 
668                 ofstream badAccnosFile;
669                 openOutputFile(badAccnosFName, badAccnosFile);
670                 
671                 ifstream inFASTA;
672                 openInputFile(filename, inFASTA);
673
674                 inFASTA.seekg(filePos->start);
675
676                 bool done = false;
677                 int count = 0;
678         
679                 while (!done) {
680                 
681                         if (m->control_pressed) {  return 0; }
682                         
683                         Sequence currSeq(inFASTA); gobble(inFASTA);
684                         if (currSeq.getName() != "") {
685                                 bool goodSeq = 1;               //      innocent until proven guilty
686                                 if(goodSeq == 1 && startPos != -1 && startPos < currSeq.getStartPos())                  {       goodSeq = 0;    }
687                                 if(goodSeq == 1 && endPos != -1 && endPos > currSeq.getEndPos())                                {       goodSeq = 0;    }
688                                 if(goodSeq == 1 && maxAmbig != -1 && maxAmbig < currSeq.getAmbigBases())                {       goodSeq = 0;    }
689                                 if(goodSeq == 1 && maxHomoP != -1 && maxHomoP < currSeq.getLongHomoPolymer())   {       goodSeq = 0;    }
690                                 if(goodSeq == 1 && minLength != -1 && minLength > currSeq.getNumBases())                {       goodSeq = 0;    }
691                                 if(goodSeq == 1 && maxLength != -1 && maxLength < currSeq.getNumBases())                {       goodSeq = 0;    }
692                                 
693                                 if(goodSeq == 1){
694                                         currSeq.printSequence(goodFile);        
695                                 }
696                                 else{
697                                         currSeq.printSequence(badFile); 
698                                         badAccnosFile << currSeq.getName() << endl;
699                                         badSeqNames.insert(currSeq.getName());
700                                 }
701                         count++;
702                         }
703                         
704                         unsigned long int pos = inFASTA.tellg();
705                         if ((pos == -1) || (pos >= filePos->end)) { break; }
706                         
707                         //report progress
708                         if((count) % 100 == 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
709                 }
710                 //report progress
711                 if((count) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
712                 
713                         
714                 goodFile.close();
715                 inFASTA.close();
716                 badFile.close();
717                 badAccnosFile.close();
718                 
719                 return count;
720         }
721         catch(exception& e) {
722                 m->errorOut(e, "ScreenSeqsCommand", "driver");
723                 exit(1);
724         }
725 }
726 //**********************************************************************************************************************
727 #ifdef USE_MPI
728 int ScreenSeqsCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& goodFile, MPI_File& badFile, MPI_File& badAccnosFile, vector<unsigned long int>& MPIPos, set<string>& badSeqNames){
729         try {
730                 string outputString = "";
731                 MPI_Status statusGood; 
732                 MPI_Status statusBad; 
733                 MPI_Status statusBadAccnos; 
734                 MPI_Status status; 
735                 int pid;
736                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
737
738                 for(int i=0;i<num;i++){
739                 
740                         if (m->control_pressed) {  return 0; }
741                         
742                         //read next sequence
743                         int length = MPIPos[start+i+1] - MPIPos[start+i];
744
745                         char* buf4 = new char[length];
746                         memcpy(buf4, outputString.c_str(), length);
747
748                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
749                         
750                         string tempBuf = buf4;  delete buf4;
751                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
752                         istringstream iss (tempBuf,istringstream::in);
753                         
754                         Sequence currSeq(iss);                  
755                         
756                         //process seq
757                         if (currSeq.getName() != "") {
758                                 bool goodSeq = 1;               //      innocent until proven guilty
759                                 if(goodSeq == 1 && startPos != -1 && startPos < currSeq.getStartPos())                  {       goodSeq = 0;    }
760                                 if(goodSeq == 1 && endPos != -1 && endPos > currSeq.getEndPos())                                {       goodSeq = 0;    }
761                                 if(goodSeq == 1 && maxAmbig != -1 && maxAmbig < currSeq.getAmbigBases())                {       goodSeq = 0;    }
762                                 if(goodSeq == 1 && maxHomoP != -1 && maxHomoP < currSeq.getLongHomoPolymer())   {       goodSeq = 0;    }
763                                 if(goodSeq == 1 && minLength != -1 && minLength > currSeq.getNumBases())                {       goodSeq = 0;    }
764                                 if(goodSeq == 1 && maxLength != -1 && maxLength < currSeq.getNumBases())                {       goodSeq = 0;    }
765                                 
766                                 if(goodSeq == 1){
767                                         outputString =  ">" + currSeq.getName() + "\n" + currSeq.getAligned() + "\n";
768                                 
769                                         //print good seq
770                                         length = outputString.length();
771                                         char* buf2 = new char[length];
772                                         memcpy(buf2, outputString.c_str(), length);
773                                         
774                                         MPI_File_write_shared(goodFile, buf2, length, MPI_CHAR, &statusGood);
775                                         delete buf2;
776                                 }
777                                 else{
778                                         outputString =  ">" + currSeq.getName() + "\n" + currSeq.getAligned() + "\n";
779                                 
780                                         //print bad seq to fasta
781                                         length = outputString.length();
782                                         char* buf2 = new char[length];
783                                         memcpy(buf2, outputString.c_str(), length);
784                                         
785                                         MPI_File_write_shared(badFile, buf2, length, MPI_CHAR, &statusBad);
786                                         delete buf2;
787
788                                         badSeqNames.insert(currSeq.getName());
789                                         
790                                         //write to bad accnos file
791                                         outputString = currSeq.getName() + "\n";
792                                 
793                                         length = outputString.length();
794                                         char* buf3 = new char[length];
795                                         memcpy(buf3, outputString.c_str(), length);
796                                         
797                                         MPI_File_write_shared(badAccnosFile, buf3, length, MPI_CHAR, &statusBadAccnos);
798                                         delete buf3;
799                                 }
800                         }
801                 }
802                                 
803                 return 1;
804         }
805         catch(exception& e) {
806                 m->errorOut(e, "ScreenSeqsCommand", "driverMPI");
807                 exit(1);
808         }
809 }
810 #endif
811 /**************************************************************************************************/
812
813 int ScreenSeqsCommand::createProcesses(string goodFileName, string badFileName, string badAccnos, string filename, set<string>& badSeqNames) {
814         try {
815 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
816                 int process = 0;
817                 int num = 0;
818                 
819                 //loop through and create all the processes you want
820                 while (process != processors) {
821                         int pid = fork();
822                         
823                         if (pid > 0) {
824                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
825                                 process++;
826                         }else if (pid == 0){
827                                 num = driver(lines[process], goodFileName + toString(getpid()) + ".temp", badFileName + toString(getpid()) + ".temp", badAccnos + toString(getpid()) + ".temp", filename, badSeqNames);
828                                 
829                                 //pass numSeqs to parent
830                                 ofstream out;
831                                 string tempFile = filename + toString(getpid()) + ".num.temp";
832                                 openOutputFile(tempFile, out);
833                                 out << num << endl;
834                                 out.close();
835                                 
836                                 exit(0);
837                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
838                 }
839                 
840                 //force parent to wait until all the processes are done
841                 for (int i=0;i<processors;i++) { 
842                         int temp = processIDS[i];
843                         wait(&temp);
844                 }
845                 
846                 for (int i = 0; i < processIDS.size(); i++) {
847                         ifstream in;
848                         string tempFile =  filename + toString(processIDS[i]) + ".num.temp";
849                         openInputFile(tempFile, in);
850                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
851                         in.close(); remove(tempFile.c_str());
852                 }
853                 
854                 return num;
855 #endif          
856         }
857         catch(exception& e) {
858                 m->errorOut(e, "ScreenSeqsCommand", "createProcesses");
859                 exit(1);
860         }
861 }
862
863 //***************************************************************************************************************
864
865