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