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