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