5 * Created by Pat Schloss on 5/30/09.
6 * Copyright 2009 Patrick D. Schloss. All rights reserved.
10 #include "seqsummarycommand.h"
13 //**********************************************************************************************************************
14 vector<string> SeqSummaryCommand::setParameters(){
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);
22 vector<string> myArray;
23 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
27 m->errorOut(e, "SeqSummaryCommand", "setParameters");
31 //**********************************************************************************************************************
32 string SeqSummaryCommand::getHelpString(){
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";
44 m->errorOut(e, "SeqSummaryCommand", "getHelpString");
48 //**********************************************************************************************************************
49 string SeqSummaryCommand::getOutputFileNameTag(string type, string inputName=""){
51 string outputFileName = "";
52 map<string, vector<string> >::iterator it;
54 //is this a type this command creates
55 it = outputTypes.find(type);
56 if (it == outputTypes.end()) { m->mothurOut("[ERROR]: this command doesn't create a " + type + " output file.\n"); }
58 if (type == "summary") { outputFileName = "summary"; }
59 else { m->mothurOut("[ERROR]: No definition for type " + type + " output file tag.\n"); m->control_pressed = true; }
61 return outputFileName;
64 m->errorOut(e, "SeqSummaryCommand", "getOutputFileNameTag");
69 //**********************************************************************************************************************
70 SeqSummaryCommand::SeqSummaryCommand(){
72 abort = true; calledHelp = true;
74 vector<string> tempOutNames;
75 outputTypes["summary"] = tempOutNames;
78 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
82 //***************************************************************************************************************
84 SeqSummaryCommand::SeqSummaryCommand(string option) {
86 abort = false; calledHelp = false;
88 //allow user to run help
89 if(option == "help") { help(); abort = true; calledHelp = true; }
90 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
93 vector<string> myArray = setParameters();
95 OptionParser parser(option);
96 map<string,string> parameters = parser.getParameters();
98 ValidParameters validParameter("summary.seqs");
99 map<string,string>::iterator it;
101 //check to make sure all parameters are valid for command
102 for (it = parameters.begin(); it != parameters.end(); it++) {
103 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
106 //if the user changes the input directory command factory will send this info to us in the output parameter
107 string inputDir = validParameter.validFile(parameters, "inputdir", false);
108 if (inputDir == "not found"){ inputDir = ""; }
111 it = parameters.find("fasta");
112 //user has given a template file
113 if(it != parameters.end()){
114 path = m->hasPath(it->second);
115 //if the user has not given a path then, add inputdir. else leave path alone.
116 if (path == "") { parameters["fasta"] = inputDir + it->second; }
119 it = parameters.find("name");
120 //user has given a template file
121 if(it != parameters.end()){
122 path = m->hasPath(it->second);
123 //if the user has not given a path then, add inputdir. else leave path alone.
124 if (path == "") { parameters["name"] = inputDir + it->second; }
128 //initialize outputTypes
129 vector<string> tempOutNames;
130 outputTypes["summary"] = tempOutNames;
132 //check for required parameters
133 fastafile = validParameter.validFile(parameters, "fasta", true);
134 if (fastafile == "not open") { abort = true; }
135 else if (fastafile == "not found") {
136 fastafile = m->getFastaFile();
137 if (fastafile != "") { m->mothurOut("Using " + fastafile + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
138 else { m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
139 }else { m->setFastaFile(fastafile); }
141 namefile = validParameter.validFile(parameters, "name", true);
142 if (namefile == "not open") { namefile = ""; abort = true; }
143 else if (namefile == "not found") { namefile = ""; }
144 else { m->setNameFile(namefile); }
146 //if the user changes the output directory command factory will send this info to us in the output parameter
147 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
149 outputDir += m->hasPath(fastafile); //if user entered a file with a path then preserve it
152 string temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
153 m->setProcessors(temp);
154 m->mothurConvert(temp, processors);
156 if (namefile == "") {
157 vector<string> files; files.push_back(fastafile);
158 parser.getNameFile(files);
163 catch(exception& e) {
164 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
168 //***************************************************************************************************************
170 int SeqSummaryCommand::execute(){
173 if (abort == true) { if (calledHelp) { return 0; } return 2; }
175 //set current fasta to fastafile
176 m->setFastaFile(fastafile);
178 string summaryFile = outputDir + m->getRootName(m->getSimpleName(fastafile)) + getOutputFileNameTag("summary");
182 vector<int> startPosition;
183 vector<int> endPosition;
184 vector<int> seqLength;
185 vector<int> ambigBases;
186 vector<int> longHomoPolymer;
188 if (namefile != "") { nameMap = m->readNames(namefile); }
190 if (m->control_pressed) { return 0; }
193 int pid, numSeqsPerProcessor;
195 int startTag = 1; int endTag = 2; int lengthTag = 3; int baseTag = 4; int lhomoTag = 5;
196 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
197 vector<unsigned long long> MPIPos;
200 MPI_Status statusOut;
203 MPI_Comm_size(MPI_COMM_WORLD, &processors);
204 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
206 char tempFileName[1024];
207 strcpy(tempFileName, fastafile.c_str());
209 char sumFileName[1024];
210 strcpy(sumFileName, summaryFile.c_str());
212 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
213 MPI_File_open(MPI_COMM_WORLD, sumFileName, outMode, MPI_INFO_NULL, &outMPI);
215 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
217 if (pid == 0) { //you are the root process
219 string outputString = "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs\n";
220 int length = outputString.length();
221 char* buf2 = new char[length];
222 memcpy(buf2, outputString.c_str(), length);
224 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
227 MPIPos = m->setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
229 for(int i = 1; i < processors; i++) {
230 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
231 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
234 //figure out how many sequences you have to do
235 numSeqsPerProcessor = numSeqs / processors;
236 int startIndex = pid * numSeqsPerProcessor;
237 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
240 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
242 }else { //i am the child process
244 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
245 MPIPos.resize(numSeqs+1);
246 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
248 //figure out how many sequences you have to align
249 numSeqsPerProcessor = numSeqs / processors;
250 int startIndex = pid * numSeqsPerProcessor;
251 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
254 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
257 MPI_File_close(&inMPI);
258 MPI_File_close(&outMPI);
259 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
262 //get the info from the child processes
263 for(int i = 1; i < processors; i++) {
265 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
267 vector<int> temp; temp.resize(size+1);
269 for(int j = 0; j < 5; j++) {
271 MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
272 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what count this is for
274 if (receiveTag == startTag) {
275 for (int k = 0; k < size; k++) { startPosition.push_back(temp[k]); }
276 }else if (receiveTag == endTag) {
277 for (int k = 0; k < size; k++) { endPosition.push_back(temp[k]); }
278 }else if (receiveTag == lengthTag) {
279 for (int k = 0; k < size; k++) { seqLength.push_back(temp[k]); }
280 }else if (receiveTag == baseTag) {
281 for (int k = 0; k < size; k++) { ambigBases.push_back(temp[k]); }
282 }else if (receiveTag == lhomoTag) {
283 for (int k = 0; k < size; k++) { longHomoPolymer.push_back(temp[k]); }
291 int size = startPosition.size();
292 MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
294 startPosition.push_back(startTag);
295 int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
296 endPosition.push_back(endTag);
297 ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
298 seqLength.push_back(lengthTag);
299 ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
300 ambigBases.push_back(baseTag);
301 ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
302 longHomoPolymer.push_back(lhomoTag);
303 ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
306 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
308 vector<unsigned long long> positions;
309 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
310 positions = m->divideFile(fastafile, processors);
311 for (int i = 0; i < (positions.size()-1); i++) { lines.push_back(new linePair(positions[i], positions[(i+1)])); }
313 positions = m->setFilePosFasta(fastafile, numSeqs);
314 if (positions.size() < processors) { processors = positions.size(); }
316 //figure out how many sequences you have to process
317 int numSeqsPerProcessor = numSeqs / processors;
318 for (int i = 0; i < processors; i++) {
319 int startIndex = i * numSeqsPerProcessor;
320 if(i == (processors - 1)){ numSeqsPerProcessor = numSeqs - i * numSeqsPerProcessor; }
321 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
327 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
329 numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile);
332 if (m->control_pressed) { return 0; }
339 sort(startPosition.begin(), startPosition.end());
340 sort(endPosition.begin(), endPosition.end());
341 sort(seqLength.begin(), seqLength.end());
342 sort(ambigBases.begin(), ambigBases.end());
343 sort(longHomoPolymer.begin(), longHomoPolymer.end());
344 int size = startPosition.size();
347 float meanStartPosition, meanEndPosition, meanSeqLength, meanAmbigBases, meanLongHomoPolymer;
348 meanStartPosition = 0; meanEndPosition = 0; meanSeqLength = 0; meanAmbigBases = 0; meanLongHomoPolymer = 0;
349 for (int i = 0; i < size; i++) {
350 meanStartPosition += startPosition[i];
351 meanEndPosition += endPosition[i];
352 meanSeqLength += seqLength[i];
353 meanAmbigBases += ambigBases[i];
354 meanLongHomoPolymer += longHomoPolymer[i];
356 //this is an int divide so the remainder is lost
357 meanStartPosition /= (float) size; meanEndPosition /= (float) size; meanLongHomoPolymer /= (float) size; meanSeqLength /= (float) size; meanAmbigBases /= (float) size;
359 int ptile0_25 = int(size * 0.025);
360 int ptile25 = int(size * 0.250);
361 int ptile50 = int(size * 0.500);
362 int ptile75 = int(size * 0.750);
363 int ptile97_5 = int(size * 0.975);
364 int ptile100 = size - 1;
366 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
367 if (startPosition[0] == -1) { startPosition[0] = 0; }
368 if (endPosition[0] == -1) { endPosition[0] = 0; }
370 if (m->control_pressed) { m->mothurRemove(summaryFile); return 0; }
372 m->mothurOutEndLine();
373 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer\tNumSeqs"); m->mothurOutEndLine();
374 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();
375 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();
376 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();
377 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();
378 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();
379 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();
380 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();
381 m->mothurOut("Mean:\t" + toString(meanStartPosition) + "\t" + toString(meanEndPosition) + "\t" + toString(meanSeqLength) + "\t" + toString(meanAmbigBases) + "\t" + toString(meanLongHomoPolymer)); m->mothurOutEndLine();
383 if (namefile == "") { m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); }
384 else { m->mothurOut("# of unique seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); m->mothurOut("total # of seqs:\t" + toString(startPosition.size())); m->mothurOutEndLine(); }
386 if (m->control_pressed) { m->mothurRemove(summaryFile); return 0; }
388 m->mothurOutEndLine();
389 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
390 m->mothurOut(summaryFile); m->mothurOutEndLine(); outputNames.push_back(summaryFile); outputTypes["summary"].push_back(summaryFile);
391 m->mothurOutEndLine();
399 catch(exception& e) {
400 m->errorOut(e, "SeqSummaryCommand", "execute");
404 /**************************************************************************************/
405 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* filePos) {
409 m->openOutputFile(sumFile, outSummary);
411 //print header if you are process 0
412 if (filePos->start == 0) {
413 outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs" << endl;
417 m->openInputFile(filename, in);
419 in.seekg(filePos->start);
426 if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
428 Sequence current(in); m->gobble(in);
430 if (current.getName() != "") {
433 if (namefile != "") {
434 //make sure this sequence is in the namefile, else error
435 map<string, int>::iterator it = nameMap.find(current.getName());
437 if (it == nameMap.end()) { m->mothurOut("[ERROR]: '" + current.getName() + "' is not in your namefile, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
438 else { num = it->second; }
441 //for each sequence this sequence represents
442 for (int i = 0; i < num; i++) {
443 startPosition.push_back(current.getStartPos());
444 endPosition.push_back(current.getEndPos());
445 seqLength.push_back(current.getNumBases());
446 ambigBases.push_back(current.getAmbigBases());
447 longHomoPolymer.push_back(current.getLongHomoPolymer());
451 outSummary << current.getName() << '\t';
452 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
453 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
454 outSummary << current.getLongHomoPolymer() << '\t' << num << endl;
457 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
458 unsigned long long pos = in.tellg();
459 if ((pos == -1) || (pos >= filePos->end)) { break; }
461 if (in.eof()) { break; }
465 //if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
468 //if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
474 catch(exception& e) {
475 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
480 /**************************************************************************************/
481 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) {
486 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
488 for(int i=0;i<num;i++){
490 if (m->control_pressed) { return 0; }
493 int length = MPIPos[start+i+1] - MPIPos[start+i];
495 char* buf4 = new char[length];
496 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
498 string tempBuf = buf4;
499 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
500 istringstream iss (tempBuf,istringstream::in);
503 Sequence current(iss);
505 if (current.getName() != "") {
508 if (namefile != "") {
509 //make sure this sequence is in the namefile, else error
510 map<string, int>::iterator it = nameMap.find(current.getName());
512 if (it == nameMap.end()) { cout << "[ERROR]: " << current.getName() << " is not in your namefile, please correct." << endl; m->control_pressed = true; }
513 else { num = it->second; }
516 //for each sequence this sequence represents
517 for (int i = 0; i < num; i++) {
518 startPosition.push_back(current.getStartPos());
519 endPosition.push_back(current.getEndPos());
520 seqLength.push_back(current.getNumBases());
521 ambigBases.push_back(current.getAmbigBases());
522 longHomoPolymer.push_back(current.getLongHomoPolymer());
525 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
526 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\t" + toString(num) + "\n";
529 length = outputString.length();
530 char* buf3 = new char[length];
531 memcpy(buf3, outputString.c_str(), length);
533 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
540 catch(exception& e) {
541 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
546 /**************************************************************************************************/
547 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
553 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
555 //loop through and create all the processes you want
556 while (process != processors) {
560 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
563 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
565 //pass numSeqs to parent
567 string tempFile = fastafile + toString(getpid()) + ".num.temp";
568 m->openOutputFile(tempFile, out);
571 out << startPosition.size() << endl;
572 for (int k = 0; k < startPosition.size(); k++) { out << startPosition[k] << '\t'; } out << endl;
573 for (int k = 0; k < endPosition.size(); k++) { out << endPosition[k] << '\t'; } out << endl;
574 for (int k = 0; k < seqLength.size(); k++) { out << seqLength[k] << '\t'; } out << endl;
575 for (int k = 0; k < ambigBases.size(); k++) { out << ambigBases[k] << '\t'; } out << endl;
576 for (int k = 0; k < longHomoPolymer.size(); k++) { out << longHomoPolymer[k] << '\t'; } out << endl;
582 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
583 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
589 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile, lines[0]);
591 //force parent to wait until all the processes are done
592 for (int i=0;i<processIDS.size();i++) {
593 int temp = processIDS[i];
597 //parent reads in and combine Filter info
598 for (int i = 0; i < processIDS.size(); i++) {
599 string tempFilename = fastafile + toString(processIDS[i]) + ".num.temp";
601 m->openInputFile(tempFilename, in);
604 in >> tempNum; m->gobble(in); num += tempNum;
605 in >> tempNum; m->gobble(in);
606 for (int k = 0; k < tempNum; k++) { in >> temp; startPosition.push_back(temp); } m->gobble(in);
607 for (int k = 0; k < tempNum; k++) { in >> temp; endPosition.push_back(temp); } m->gobble(in);
608 for (int k = 0; k < tempNum; k++) { in >> temp; seqLength.push_back(temp); } m->gobble(in);
609 for (int k = 0; k < tempNum; k++) { in >> temp; ambigBases.push_back(temp); } m->gobble(in);
610 for (int k = 0; k < tempNum; k++) { in >> temp; longHomoPolymer.push_back(temp); } m->gobble(in);
613 m->mothurRemove(tempFilename);
615 m->appendFiles((sumFile + toString(processIDS[i]) + ".temp"), sumFile);
616 m->mothurRemove((sumFile + toString(processIDS[i]) + ".temp"));
620 //////////////////////////////////////////////////////////////////////////////////////////////////////
621 //Windows version shared memory, so be careful when passing variables through the seqSumData struct.
622 //Above fork() will clone, so memory is separate, but that's not the case with windows,
623 //Taking advantage of shared memory to allow both threads to add info to vectors.
624 //////////////////////////////////////////////////////////////////////////////////////////////////////
626 vector<seqSumData*> pDataArray;
627 DWORD dwThreadIdArray[processors-1];
628 HANDLE hThreadArray[processors-1];
630 //Create processor worker threads.
631 for( int i=0; i<processors-1; i++ ){
633 string extension = "";
634 if (i != 0) { extension = toString(i) + ".temp"; processIDS.push_back(i); }
635 // Allocate memory for thread data.
636 seqSumData* tempSum = new seqSumData(filename, (sumFile+extension), m, lines[i]->start, lines[i]->end, namefile, nameMap);
637 pDataArray.push_back(tempSum);
639 //MySeqSumThreadFunction is in header. It must be global or static to work with the threads.
640 //default security attributes, thread function name, argument to thread function, use default creation flags, returns the thread identifier
641 hThreadArray[i] = CreateThread(NULL, 0, MySeqSumThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
645 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, (sumFile+toString(processors-1)+".temp"), lines[processors-1]);
646 processIDS.push_back(processors-1);
648 //Wait until all threads have terminated.
649 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
651 //Close all thread handles and free memory allocations.
652 for(int i=0; i < pDataArray.size(); i++){
653 num += pDataArray[i]->count;
654 for (int k = 0; k < pDataArray[i]->startPosition.size(); k++) { startPosition.push_back(pDataArray[i]->startPosition[k]); }
655 for (int k = 0; k < pDataArray[i]->endPosition.size(); k++) { endPosition.push_back(pDataArray[i]->endPosition[k]); }
656 for (int k = 0; k < pDataArray[i]->seqLength.size(); k++) { seqLength.push_back(pDataArray[i]->seqLength[k]); }
657 for (int k = 0; k < pDataArray[i]->ambigBases.size(); k++) { ambigBases.push_back(pDataArray[i]->ambigBases[k]); }
658 for (int k = 0; k < pDataArray[i]->longHomoPolymer.size(); k++) { longHomoPolymer.push_back(pDataArray[i]->longHomoPolymer[k]); }
659 CloseHandle(hThreadArray[i]);
660 delete pDataArray[i];
664 for(int i=0;i<processIDS.size();i++){
665 m->appendFiles((sumFile + toString(processIDS[i]) + ".temp"), sumFile);
666 m->mothurRemove((sumFile + toString(processIDS[i]) + ".temp"));
671 catch(exception& e) {
672 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
676 /**********************************************************************************************************************/