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