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