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