]> git.donarmstrong.com Git - mothur.git/blob - seqsummarycommand.cpp
added pipeline commands which involved change to command factory and command class...
[mothur.git] / seqsummarycommand.cpp
1 /*
2  *  seqcoordcommand.cpp
3  *  Mothur
4  *
5  *  Created by Pat Schloss on 5/30/09.
6  *  Copyright 2009 Patrick D. Schloss. All rights reserved.
7  *
8  */
9
10 #include "seqsummarycommand.h"
11 #include "sequence.hpp"
12
13 //**********************************************************************************************************************
14 vector<string> SeqSummaryCommand::getValidParameters(){ 
15         try {
16                 string Array[] =  {"fasta","processors","outputdir","inputdir"};
17                 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
18                 return myArray;
19         }
20         catch(exception& e) {
21                 m->errorOut(e, "SeqSummaryCommand", "getValidParameters");
22                 exit(1);
23         }
24 }
25 //**********************************************************************************************************************
26 SeqSummaryCommand::SeqSummaryCommand(){ 
27         try {
28                 //initialize outputTypes
29                 vector<string> tempOutNames;
30                 outputTypes["summary"] = tempOutNames;
31         }
32         catch(exception& e) {
33                 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
34                 exit(1);
35         }
36 }
37 //**********************************************************************************************************************
38 vector<string> SeqSummaryCommand::getRequiredParameters(){      
39         try {
40                 string Array[] =  {"fasta"};
41                 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
42                 return myArray;
43         }
44         catch(exception& e) {
45                 m->errorOut(e, "SeqSummaryCommand", "getRequiredParameters");
46                 exit(1);
47         }
48 }
49 //**********************************************************************************************************************
50 vector<string> SeqSummaryCommand::getRequiredFiles(){   
51         try {
52                 vector<string> myArray;
53                 return myArray;
54         }
55         catch(exception& e) {
56                 m->errorOut(e, "SeqSummaryCommand", "getRequiredFiles");
57                 exit(1);
58         }
59 }
60 //***************************************************************************************************************
61
62 SeqSummaryCommand::SeqSummaryCommand(string option)  {
63         try {
64                 abort = false;
65                 
66                 //allow user to run help
67                 if(option == "help") { help(); abort = true; }
68                 
69                 else {
70                         //valid paramters for this command
71                         string Array[] =  {"fasta","processors","outputdir","inputdir"};
72                         vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
73                         
74                         OptionParser parser(option);
75                         map<string,string> parameters = parser.getParameters();
76                         
77                         ValidParameters validParameter("summary.seqs");
78                         map<string,string>::iterator it;
79                         
80                         //check to make sure all parameters are valid for command
81                         for (it = parameters.begin(); it != parameters.end(); it++) { 
82                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
83                         }
84                         
85                         //if the user changes the input directory command factory will send this info to us in the output parameter 
86                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
87                         if (inputDir == "not found"){   inputDir = "";          }
88                         else {
89                                 string path;
90                                 it = parameters.find("fasta");
91                                 //user has given a template file
92                                 if(it != parameters.end()){ 
93                                         path = m->hasPath(it->second);
94                                         //if the user has not given a path then, add inputdir. else leave path alone.
95                                         if (path == "") {       parameters["fasta"] = inputDir + it->second;            }
96                                 }
97                         }
98                         
99                         //initialize outputTypes
100                         vector<string> tempOutNames;
101                         outputTypes["summary"] = tempOutNames;
102                         
103                         //check for required parameters
104                         fastafile = validParameter.validFile(parameters, "fasta", true);
105                         if (fastafile == "not open") { abort = true; }
106                         else if (fastafile == "not found") { fastafile = ""; m->mothurOut("fasta is a required parameter for the summary.seqs command."); m->mothurOutEndLine(); abort = true;  }       
107                         
108                         //if the user changes the output directory command factory will send this info to us in the output parameter 
109                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  
110                                 outputDir = ""; 
111                                 outputDir += m->hasPath(fastafile); //if user entered a file with a path then preserve it       
112                         }
113                         
114                         string temp = validParameter.validFile(parameters, "processors", false);        if (temp == "not found"){       temp = "1";                             }
115                         convert(temp, processors); 
116
117
118                 }
119         }
120         catch(exception& e) {
121                 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
122                 exit(1);
123         }
124 }
125 //**********************************************************************************************************************
126
127 void SeqSummaryCommand::help(){
128         try {
129                 m->mothurOut("The summary.seqs command reads a fastafile and summarizes the sequences.\n");
130                 m->mothurOut("The summary.seqs command parameters are fasta and processors, fasta is required.\n");
131                 m->mothurOut("The summary.seqs command should be in the following format: \n");
132                 m->mothurOut("summary.seqs(fasta=yourFastaFile, processors=2) \n");
133                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");     
134         }
135         catch(exception& e) {
136                 m->errorOut(e, "SeqSummaryCommand", "help");
137                 exit(1);
138         }
139 }
140
141 //***************************************************************************************************************
142
143 SeqSummaryCommand::~SeqSummaryCommand(){        /*      do nothing      */      }
144
145 //***************************************************************************************************************
146
147 int SeqSummaryCommand::execute(){
148         try{
149                 
150                 if (abort == true) { return 0; }
151                 
152                 string summaryFile = outputDir + m->getSimpleName(fastafile) + ".summary";
153                                 
154                 int numSeqs = 0;
155                 
156                 vector<int> startPosition;
157                 vector<int> endPosition;
158                 vector<int> seqLength;
159                 vector<int> ambigBases;
160                 vector<int> longHomoPolymer;
161                         
162 #ifdef USE_MPI  
163                                 int pid, numSeqsPerProcessor; 
164                                 int tag = 2001;
165                                 int startTag = 1; int endTag = 2; int lengthTag = 3; int baseTag = 4; int lhomoTag = 5;
166                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
167                                 vector<unsigned long int> MPIPos;
168                                 
169                                 MPI_Status status; 
170                                 MPI_Status statusOut;
171                                 MPI_File inMPI; 
172                                 MPI_File outMPI; 
173                                 MPI_Comm_size(MPI_COMM_WORLD, &processors);
174                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
175                                                         
176                                 char tempFileName[1024];
177                                 strcpy(tempFileName, fastafile.c_str());
178                                 
179                                 char sumFileName[1024];
180                                 strcpy(sumFileName, summaryFile.c_str());
181                 
182                                 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
183                                 MPI_File_open(MPI_COMM_WORLD, sumFileName, outMode, MPI_INFO_NULL, &outMPI);
184                                 
185                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI); return 0;  }
186                                 
187                                 if (pid == 0) { //you are the root process
188                                                 //print header
189                                                 string outputString = "seqname\tstart\tend\tnbases\tambigs\tpolymer\n"; 
190                                                 int length = outputString.length();
191                                                 char* buf2 = new char[length];
192                                                 memcpy(buf2, outputString.c_str(), length);
193                                         
194                                                 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
195                                                 delete buf2;
196                                                 
197                                                 MPIPos = m->setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
198                                         
199                                                 for(int i = 1; i < processors; i++) { 
200                                                         MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
201                                                         MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
202                                                 }
203                                                 
204                                                 //figure out how many sequences you have to do
205                                                 numSeqsPerProcessor = numSeqs / processors;
206                                                 int startIndex =  pid * numSeqsPerProcessor;
207                                                 if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
208                                                 
209                                                 //do your part
210                                                 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
211                                                 
212                                 }else { //i am the child process
213                         
214                                         MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
215                                         MPIPos.resize(numSeqs+1);
216                                         MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
217                                 
218                                         //figure out how many sequences you have to align
219                                         numSeqsPerProcessor = numSeqs / processors;
220                                         int startIndex =  pid * numSeqsPerProcessor;
221                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
222                                 
223                                         //do your part
224                                         MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
225                                 }
226                                 
227                                 MPI_File_close(&inMPI);
228                                 MPI_File_close(&outMPI);
229                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
230                                 
231                                 if (pid == 0) {
232                                         //get the info from the child processes
233                                         for(int i = 1; i < processors; i++) { 
234                                                 int size;
235                                                 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
236
237                                                 vector<int> temp; temp.resize(size+1);
238                                                 
239                                                 for(int j = 0; j < 5; j++) { 
240                                                 
241                                                         MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status); 
242                                                         int receiveTag = temp[temp.size()-1];  //child process added a int to the end to indicate what count this is for
243                                                         
244                                                         if (receiveTag == startTag) { 
245                                                                 for (int k = 0; k < size; k++) {                startPosition.push_back(temp[k]);       }
246                                                         }else if (receiveTag == endTag) { 
247                                                                 for (int k = 0; k < size; k++) {                endPosition.push_back(temp[k]); }
248                                                         }else if (receiveTag == lengthTag) { 
249                                                                 for (int k = 0; k < size; k++) {                seqLength.push_back(temp[k]);   }
250                                                         }else if (receiveTag == baseTag) { 
251                                                                 for (int k = 0; k < size; k++) {                ambigBases.push_back(temp[k]);  }
252                                                         }else if (receiveTag == lhomoTag) { 
253                                                                 for (int k = 0; k < size; k++) {                longHomoPolymer.push_back(temp[k]);     }
254                                                         }
255                                                 } 
256                                         }
257
258                                 }else{
259                                 
260                                         //send my counts
261                                         int size = startPosition.size();
262                                         MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
263                                         
264                                         startPosition.push_back(startTag);
265                                         int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
266                                         endPosition.push_back(endTag);
267                                         ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
268                                         seqLength.push_back(lengthTag);
269                                         ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
270                                         ambigBases.push_back(baseTag);
271                                         ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
272                                         longHomoPolymer.push_back(lhomoTag);
273                                         ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
274                                 }
275                                 
276                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
277 #else
278                         vector<unsigned long int> positions = m->divideFile(fastafile, processors);
279                                 
280                         for (int i = 0; i < (positions.size()-1); i++) {
281                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
282                         }       
283
284                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
285                                 if(processors == 1){
286                                         numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
287                                 }else{
288                                         numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile); 
289                                         
290                                         rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
291                                         //append files
292                                         for(int i=1;i<processors;i++){
293                                                 m->appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
294                                                 remove((summaryFile + toString(processIDS[i]) + ".temp").c_str());
295                                         }
296                                 }
297                                 
298                                 if (m->control_pressed) {  return 0; }
299                 #else
300                                 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
301                                 if (m->control_pressed) {  return 0; }
302                 #endif
303 #endif
304                         
305                 #ifdef USE_MPI
306                         if (pid == 0) { 
307                 #endif
308                 
309                 sort(startPosition.begin(), startPosition.end());
310                 sort(endPosition.begin(), endPosition.end());
311                 sort(seqLength.begin(), seqLength.end());
312                 sort(ambigBases.begin(), ambigBases.end());
313                 sort(longHomoPolymer.begin(), longHomoPolymer.end());
314                 
315                 int ptile0_25   = int(numSeqs * 0.025);
316                 int ptile25             = int(numSeqs * 0.250);
317                 int ptile50             = int(numSeqs * 0.500);
318                 int ptile75             = int(numSeqs * 0.750);
319                 int ptile97_5   = int(numSeqs * 0.975);
320                 int ptile100    = numSeqs - 1;
321                 
322                 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
323                 if (startPosition[0] == -1) {  startPosition[0] = 0;    }
324                 if (endPosition[0] == -1)       {  endPosition[0] = 0;          }
325                 
326                 if (m->control_pressed) {  remove(summaryFile.c_str()); return 0; }
327                 
328                 m->mothurOutEndLine();
329                 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
330                 m->mothurOut("Minimum:\t" + toString(startPosition[0]) + "\t" + toString(endPosition[0]) + "\t" + toString(seqLength[0]) + "\t" + toString(ambigBases[0]) + "\t" + toString(longHomoPolymer[0])); m->mothurOutEndLine();
331                 m->mothurOut("2.5%-tile:\t" + toString(startPosition[ptile0_25]) + "\t" + toString(endPosition[ptile0_25]) + "\t" + toString(seqLength[ptile0_25]) + "\t" + toString(ambigBases[ptile0_25]) + "\t"+ toString(longHomoPolymer[ptile0_25])); m->mothurOutEndLine();
332                 m->mothurOut("25%-tile:\t" + toString(startPosition[ptile25]) + "\t" + toString(endPosition[ptile25]) + "\t" + toString(seqLength[ptile25]) + "\t" + toString(ambigBases[ptile25]) + "\t" + toString(longHomoPolymer[ptile25])); m->mothurOutEndLine();
333                 m->mothurOut("Median: \t" + toString(startPosition[ptile50]) + "\t" + toString(endPosition[ptile50]) + "\t" + toString(seqLength[ptile50]) + "\t" + toString(ambigBases[ptile50]) + "\t" + toString(longHomoPolymer[ptile50])); m->mothurOutEndLine();
334                 m->mothurOut("75%-tile:\t" + toString(startPosition[ptile75]) + "\t" + toString(endPosition[ptile75]) + "\t" + toString(seqLength[ptile75]) + "\t" + toString(ambigBases[ptile75]) + "\t" + toString(longHomoPolymer[ptile75])); m->mothurOutEndLine();
335                 m->mothurOut("97.5%-tile:\t" + toString(startPosition[ptile97_5]) + "\t" + toString(endPosition[ptile97_5]) + "\t" + toString(seqLength[ptile97_5]) + "\t" + toString(ambigBases[ptile97_5]) + "\t" + toString(longHomoPolymer[ptile97_5])); m->mothurOutEndLine();
336                 m->mothurOut("Maximum:\t" + toString(startPosition[ptile100]) + "\t" + toString(endPosition[ptile100]) + "\t" + toString(seqLength[ptile100]) + "\t" + toString(ambigBases[ptile100]) + "\t" + toString(longHomoPolymer[ptile100])); m->mothurOutEndLine();
337                 m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine();
338                 
339                 if (m->control_pressed) {  remove(summaryFile.c_str()); return 0; }
340                 
341                 m->mothurOutEndLine();
342                 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
343                 m->mothurOut(summaryFile); m->mothurOutEndLine();       outputNames.push_back(summaryFile); outputTypes["summary"].push_back(summaryFile);
344                 m->mothurOutEndLine();
345                 
346                 #ifdef USE_MPI
347                         }
348                 #endif
349
350                 return 0;
351         }
352         catch(exception& e) {
353                 m->errorOut(e, "SeqSummaryCommand", "execute");
354                 exit(1);
355         }
356 }
357 /**************************************************************************************/
358 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* filePos) {   
359         try {
360                 
361                 ofstream outSummary;
362                 m->openOutputFile(sumFile, outSummary);
363                 
364                 //print header if you are process 0
365                 if (filePos->start == 0) {
366                         outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer" << endl;   
367                 }
368                                 
369                 ifstream in;
370                 m->openInputFile(filename, in);
371                                 
372                 in.seekg(filePos->start);
373
374                 bool done = false;
375                 int count = 0;
376         
377                 while (!done) {
378                                 
379                         if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
380                                         
381                         Sequence current(in); m->gobble(in);
382         
383                         if (current.getName() != "") {
384                                 startPosition.push_back(current.getStartPos());
385                                 endPosition.push_back(current.getEndPos());
386                                 seqLength.push_back(current.getNumBases());
387                                 ambigBases.push_back(current.getAmbigBases());
388                                 longHomoPolymer.push_back(current.getLongHomoPolymer());
389                                 
390                                 outSummary << current.getName() << '\t';
391                                 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
392                                 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
393                                 outSummary << current.getLongHomoPolymer() << endl;
394                                 count++;
395                         }
396                         
397                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
398                                 unsigned long int pos = in.tellg();
399                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
400                         #else
401                                 if (in.eof()) { break; }
402                         #endif
403                         
404                         //report progress
405                         //if((count) % 100 == 0){       m->mothurOut(toString(count)); m->mothurOutEndLine();           }
406                 }
407                 //report progress
408                 //if((count) % 100 != 0){       m->mothurOut(toString(count)); m->mothurOutEndLine();           }
409                 
410                 in.close();
411                 
412                 return count;
413         }
414         catch(exception& e) {
415                 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
416                 exit(1);
417         }
418 }
419 #ifdef USE_MPI
420 /**************************************************************************************/
421 int SeqSummaryCommand::MPICreateSummary(int start, int num, vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {        
422         try {
423                 
424                 int pid;
425                 MPI_Status status; 
426                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
427
428                 
429                 for(int i=0;i<num;i++){
430                         
431                         if (m->control_pressed) { return 0; }
432                         
433                         //read next sequence
434                         int length = MPIPos[start+i+1] - MPIPos[start+i];
435         
436                         char* buf4 = new char[length];
437                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
438                         
439                         string tempBuf = buf4;
440                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
441                         istringstream iss (tempBuf,istringstream::in);
442                         delete buf4;
443
444                         Sequence current(iss);  
445
446                         if (current.getName() != "") {
447                                 startPosition.push_back(current.getStartPos());
448                                 endPosition.push_back(current.getEndPos());
449                                 seqLength.push_back(current.getNumBases());
450                                 ambigBases.push_back(current.getAmbigBases());
451                                 longHomoPolymer.push_back(current.getLongHomoPolymer());
452                                 
453                                 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
454                                 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\n";
455                                 
456                                 //output to file
457                                 length = outputString.length();
458                                 char* buf3 = new char[length];
459                                 memcpy(buf3, outputString.c_str(), length);
460                                         
461                                 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
462                                 delete buf3;
463                         }       
464                 }
465                 
466                 return 0;
467         }
468         catch(exception& e) {
469                 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
470                 exit(1);
471         }
472 }
473 #endif
474 /**************************************************************************************************/
475 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
476         try {
477 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
478                 int process = 0;
479                 int num = 0;
480                 processIDS.clear();
481                 
482                 //loop through and create all the processes you want
483                 while (process != processors) {
484                         int pid = fork();
485                         
486                         if (pid > 0) {
487                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
488                                 process++;
489                         }else if (pid == 0){
490                                 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
491                                 
492                                 //pass numSeqs to parent
493                                 ofstream out;
494                                 string tempFile = fastafile + toString(getpid()) + ".num.temp";
495                                 m->openOutputFile(tempFile, out);
496                                 
497                                 out << num << endl;
498                                 for (int k = 0; k < startPosition.size(); k++)          {               out << startPosition[k] << '\t'; }  out << endl;
499                                 for (int k = 0; k < endPosition.size(); k++)            {               out << endPosition[k] << '\t'; }  out << endl;
500                                 for (int k = 0; k < seqLength.size(); k++)                      {               out << seqLength[k] << '\t'; }  out << endl;
501                                 for (int k = 0; k < ambigBases.size(); k++)                     {               out << ambigBases[k] << '\t'; }  out << endl;
502                                 for (int k = 0; k < longHomoPolymer.size(); k++)        {               out << longHomoPolymer[k] << '\t'; }  out << endl;
503                                 
504                                 out.close();
505                                 
506                                 exit(0);
507                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
508                 }
509                 
510                 //force parent to wait until all the processes are done
511                 for (int i=0;i<processors;i++) { 
512                         int temp = processIDS[i];
513                         wait(&temp);
514                 }
515                 
516                 //parent reads in and combine Filter info
517                 for (int i = 0; i < processIDS.size(); i++) {
518                         string tempFilename = fastafile + toString(processIDS[i]) + ".num.temp";
519                         ifstream in;
520                         m->openInputFile(tempFilename, in);
521                         
522                         int temp, tempNum;
523                         in >> tempNum; m->gobble(in); num += tempNum;
524                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; startPosition.push_back(temp);              }               m->gobble(in);
525                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; endPosition.push_back(temp);                }               m->gobble(in);
526                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; seqLength.push_back(temp);                  }               m->gobble(in);
527                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; ambigBases.push_back(temp);                 }               m->gobble(in);
528                         for (int k = 0; k < tempNum; k++)                       {               in >> temp; longHomoPolymer.push_back(temp);    }               m->gobble(in);
529                                 
530                         in.close();
531                         remove(tempFilename.c_str());
532                 }
533                 
534                 return num;
535 #endif          
536         }
537         catch(exception& e) {
538                 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
539                 exit(1);
540         }
541 }
542 //***************************************************************************************************************
543
544