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