5 * Created by Pat Schloss on 5/30/09.
6 * Copyright 2009 Patrick D. Schloss. All rights reserved.
10 #include "seqsummarycommand.h"
11 #include "sequence.hpp"
13 //**********************************************************************************************************************
14 vector<string> SeqSummaryCommand::getValidParameters(){
16 string Array[] = {"fasta","name","processors","outputdir","inputdir"};
17 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
21 m->errorOut(e, "SeqSummaryCommand", "getValidParameters");
25 //**********************************************************************************************************************
26 SeqSummaryCommand::SeqSummaryCommand(){
28 abort = true; calledHelp = true;
29 vector<string> tempOutNames;
30 outputTypes["summary"] = tempOutNames;
33 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
37 //**********************************************************************************************************************
38 vector<string> SeqSummaryCommand::getRequiredParameters(){
40 string Array[] = {"fasta"};
41 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
45 m->errorOut(e, "SeqSummaryCommand", "getRequiredParameters");
49 //**********************************************************************************************************************
50 vector<string> SeqSummaryCommand::getRequiredFiles(){
52 vector<string> myArray;
56 m->errorOut(e, "SeqSummaryCommand", "getRequiredFiles");
60 //***************************************************************************************************************
62 SeqSummaryCommand::SeqSummaryCommand(string option) {
64 abort = false; calledHelp = false;
66 //allow user to run help
67 if(option == "help") { help(); abort = true; calledHelp = true; }
70 //valid paramters for this command
71 string Array[] = {"fasta","name","processors","outputdir","inputdir"};
72 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
74 OptionParser parser(option);
75 map<string,string> parameters = parser.getParameters();
77 ValidParameters validParameter("summary.seqs");
78 map<string,string>::iterator it;
80 //check to make sure all parameters are valid for command
81 for (it = parameters.begin(); it != parameters.end(); it++) {
82 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
85 //if the user changes the input directory command factory will send this info to us in the output parameter
86 string inputDir = validParameter.validFile(parameters, "inputdir", false);
87 if (inputDir == "not found"){ inputDir = ""; }
90 it = parameters.find("fasta");
91 //user has given a template file
92 if(it != parameters.end()){
93 path = m->hasPath(it->second);
94 //if the user has not given a path then, add inputdir. else leave path alone.
95 if (path == "") { parameters["fasta"] = inputDir + it->second; }
98 it = parameters.find("name");
99 //user has given a template file
100 if(it != parameters.end()){
101 path = m->hasPath(it->second);
102 //if the user has not given a path then, add inputdir. else leave path alone.
103 if (path == "") { parameters["name"] = inputDir + it->second; }
107 //initialize outputTypes
108 vector<string> tempOutNames;
109 outputTypes["summary"] = tempOutNames;
111 //check for required parameters
112 fastafile = validParameter.validFile(parameters, "fasta", true);
113 if (fastafile == "not open") { abort = true; }
114 else if (fastafile == "not found") { fastafile = ""; m->mothurOut("fasta is a required parameter for the summary.seqs command."); m->mothurOutEndLine(); abort = true; }
116 namefile = validParameter.validFile(parameters, "name", true);
117 if (namefile == "not open") { namefile = ""; abort = true; }
118 else if (namefile == "not found") { namefile = ""; }
120 //if the user changes the output directory command factory will send this info to us in the output parameter
121 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
123 outputDir += m->hasPath(fastafile); //if user entered a file with a path then preserve it
126 string temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
127 convert(temp, processors);
132 catch(exception& e) {
133 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
137 //**********************************************************************************************************************
139 void SeqSummaryCommand::help(){
141 m->mothurOut("The summary.seqs command reads a fastafile and summarizes the sequences.\n");
142 m->mothurOut("The summary.seqs command parameters are fasta, name and processors, fasta is required.\n");
143 m->mothurOut("The name parameter allows you to enter a name file associated with your fasta file. \n");
144 m->mothurOut("The summary.seqs command should be in the following format: \n");
145 m->mothurOut("summary.seqs(fasta=yourFastaFile, processors=2) \n");
146 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");
148 catch(exception& e) {
149 m->errorOut(e, "SeqSummaryCommand", "help");
154 //***************************************************************************************************************
156 SeqSummaryCommand::~SeqSummaryCommand(){ /* do nothing */ }
158 //***************************************************************************************************************
160 int SeqSummaryCommand::execute(){
163 if (abort == true) { if (calledHelp) { return 0; } return 2; }
165 string summaryFile = outputDir + m->getSimpleName(fastafile) + ".summary";
169 vector<int> startPosition;
170 vector<int> endPosition;
171 vector<int> seqLength;
172 vector<int> ambigBases;
173 vector<int> longHomoPolymer;
175 if (namefile != "") { nameMap = m->readNames(namefile); }
177 if (m->control_pressed) { return 0; }
180 int pid, numSeqsPerProcessor;
182 int startTag = 1; int endTag = 2; int lengthTag = 3; int baseTag = 4; int lhomoTag = 5;
183 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
184 vector<unsigned long int> MPIPos;
187 MPI_Status statusOut;
190 MPI_Comm_size(MPI_COMM_WORLD, &processors);
191 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
193 char tempFileName[1024];
194 strcpy(tempFileName, fastafile.c_str());
196 char sumFileName[1024];
197 strcpy(sumFileName, summaryFile.c_str());
199 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
200 MPI_File_open(MPI_COMM_WORLD, sumFileName, outMode, MPI_INFO_NULL, &outMPI);
202 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
204 if (pid == 0) { //you are the root process
206 string outputString = "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs\n";
207 int length = outputString.length();
208 char* buf2 = new char[length];
209 memcpy(buf2, outputString.c_str(), length);
211 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
214 MPIPos = m->setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
216 for(int i = 1; i < processors; i++) {
217 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
218 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
221 //figure out how many sequences you have to do
222 numSeqsPerProcessor = numSeqs / processors;
223 int startIndex = pid * numSeqsPerProcessor;
224 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
227 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
229 }else { //i am the child process
231 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
232 MPIPos.resize(numSeqs+1);
233 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
235 //figure out how many sequences you have to align
236 numSeqsPerProcessor = numSeqs / processors;
237 int startIndex = pid * numSeqsPerProcessor;
238 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
241 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
244 MPI_File_close(&inMPI);
245 MPI_File_close(&outMPI);
246 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
249 //get the info from the child processes
250 for(int i = 1; i < processors; i++) {
252 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
254 vector<int> temp; temp.resize(size+1);
256 for(int j = 0; j < 5; j++) {
258 MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
259 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what count this is for
261 if (receiveTag == startTag) {
262 for (int k = 0; k < size; k++) { startPosition.push_back(temp[k]); }
263 }else if (receiveTag == endTag) {
264 for (int k = 0; k < size; k++) { endPosition.push_back(temp[k]); }
265 }else if (receiveTag == lengthTag) {
266 for (int k = 0; k < size; k++) { seqLength.push_back(temp[k]); }
267 }else if (receiveTag == baseTag) {
268 for (int k = 0; k < size; k++) { ambigBases.push_back(temp[k]); }
269 }else if (receiveTag == lhomoTag) {
270 for (int k = 0; k < size; k++) { longHomoPolymer.push_back(temp[k]); }
278 int size = startPosition.size();
279 MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
281 startPosition.push_back(startTag);
282 int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
283 endPosition.push_back(endTag);
284 ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
285 seqLength.push_back(lengthTag);
286 ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
287 ambigBases.push_back(baseTag);
288 ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
289 longHomoPolymer.push_back(lhomoTag);
290 ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
293 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
295 vector<unsigned long int> positions = m->divideFile(fastafile, processors);
297 for (int i = 0; i < (positions.size()-1); i++) {
298 lines.push_back(new linePair(positions[i], positions[(i+1)]));
301 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
303 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
305 numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile);
307 rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
309 for(int i=1;i<processors;i++){
310 m->appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
311 remove((summaryFile + toString(processIDS[i]) + ".temp").c_str());
315 if (m->control_pressed) { return 0; }
317 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
318 if (m->control_pressed) { return 0; }
326 sort(startPosition.begin(), startPosition.end());
327 sort(endPosition.begin(), endPosition.end());
328 sort(seqLength.begin(), seqLength.end());
329 sort(ambigBases.begin(), ambigBases.end());
330 sort(longHomoPolymer.begin(), longHomoPolymer.end());
331 int size = startPosition.size();
333 int ptile0_25 = int(size * 0.025);
334 int ptile25 = int(size * 0.250);
335 int ptile50 = int(size * 0.500);
336 int ptile75 = int(size * 0.750);
337 int ptile97_5 = int(size * 0.975);
338 int ptile100 = size - 1;
340 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
341 if (startPosition[0] == -1) { startPosition[0] = 0; }
342 if (endPosition[0] == -1) { endPosition[0] = 0; }
344 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
346 m->mothurOutEndLine();
347 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
348 m->mothurOut("Minimum:\t" + toString(startPosition[0]) + "\t" + toString(endPosition[0]) + "\t" + toString(seqLength[0]) + "\t" + toString(ambigBases[0]) + "\t" + toString(longHomoPolymer[0])); m->mothurOutEndLine();
349 m->mothurOut("2.5%-tile:\t" + toString(startPosition[ptile0_25]) + "\t" + toString(endPosition[ptile0_25]) + "\t" + toString(seqLength[ptile0_25]) + "\t" + toString(ambigBases[ptile0_25]) + "\t"+ toString(longHomoPolymer[ptile0_25])); m->mothurOutEndLine();
350 m->mothurOut("25%-tile:\t" + toString(startPosition[ptile25]) + "\t" + toString(endPosition[ptile25]) + "\t" + toString(seqLength[ptile25]) + "\t" + toString(ambigBases[ptile25]) + "\t" + toString(longHomoPolymer[ptile25])); m->mothurOutEndLine();
351 m->mothurOut("Median: \t" + toString(startPosition[ptile50]) + "\t" + toString(endPosition[ptile50]) + "\t" + toString(seqLength[ptile50]) + "\t" + toString(ambigBases[ptile50]) + "\t" + toString(longHomoPolymer[ptile50])); m->mothurOutEndLine();
352 m->mothurOut("75%-tile:\t" + toString(startPosition[ptile75]) + "\t" + toString(endPosition[ptile75]) + "\t" + toString(seqLength[ptile75]) + "\t" + toString(ambigBases[ptile75]) + "\t" + toString(longHomoPolymer[ptile75])); m->mothurOutEndLine();
353 m->mothurOut("97.5%-tile:\t" + toString(startPosition[ptile97_5]) + "\t" + toString(endPosition[ptile97_5]) + "\t" + toString(seqLength[ptile97_5]) + "\t" + toString(ambigBases[ptile97_5]) + "\t" + toString(longHomoPolymer[ptile97_5])); m->mothurOutEndLine();
354 m->mothurOut("Maximum:\t" + toString(startPosition[ptile100]) + "\t" + toString(endPosition[ptile100]) + "\t" + toString(seqLength[ptile100]) + "\t" + toString(ambigBases[ptile100]) + "\t" + toString(longHomoPolymer[ptile100])); m->mothurOutEndLine();
355 if (namefile == "") { m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); }
356 else { m->mothurOut("# of unique seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); m->mothurOut("total # of seqs:\t" + toString(startPosition.size())); m->mothurOutEndLine(); }
358 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
360 m->mothurOutEndLine();
361 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
362 m->mothurOut(summaryFile); m->mothurOutEndLine(); outputNames.push_back(summaryFile); outputTypes["summary"].push_back(summaryFile);
363 m->mothurOutEndLine();
371 catch(exception& e) {
372 m->errorOut(e, "SeqSummaryCommand", "execute");
376 /**************************************************************************************/
377 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 m->openOutputFile(sumFile, outSummary);
383 //print header if you are process 0
384 if (filePos->start == 0) {
385 outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs" << endl;
389 m->openInputFile(filename, in);
391 in.seekg(filePos->start);
398 if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
400 Sequence current(in); m->gobble(in);
402 if (current.getName() != "") {
405 if (namefile != "") {
406 //make sure this sequence is in the namefile, else error
407 map<string, int>::iterator it = nameMap.find(current.getName());
409 if (it == nameMap.end()) { m->mothurOut("[ERROR]: " + current.getName() + " is not in your namefile, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
410 else { num = it->second; }
413 //for each sequence this sequence represents
414 for (int i = 0; i < num; i++) {
415 startPosition.push_back(current.getStartPos());
416 endPosition.push_back(current.getEndPos());
417 seqLength.push_back(current.getNumBases());
418 ambigBases.push_back(current.getAmbigBases());
419 longHomoPolymer.push_back(current.getLongHomoPolymer());
423 outSummary << current.getName() << '\t';
424 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
425 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
426 outSummary << current.getLongHomoPolymer() << '\t' << num << endl;
429 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
430 unsigned long int pos = in.tellg();
431 if ((pos == -1) || (pos >= filePos->end)) { break; }
433 if (in.eof()) { break; }
437 //if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
440 //if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
446 catch(exception& e) {
447 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
452 /**************************************************************************************/
453 int SeqSummaryCommand::MPICreateSummary(int start, int num, vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, MPI_File& inMPI, MPI_File& outMPI, vector<unsigned long int>& MPIPos) {
458 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
460 for(int i=0;i<num;i++){
462 if (m->control_pressed) { return 0; }
465 int length = MPIPos[start+i+1] - MPIPos[start+i];
467 char* buf4 = new char[length];
468 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
470 string tempBuf = buf4;
471 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
472 istringstream iss (tempBuf,istringstream::in);
475 Sequence current(iss);
477 if (current.getName() != "") {
480 if (namefile != "") {
481 //make sure this sequence is in the namefile, else error
482 map<string, int>::iterator it = nameMap.find(current.getName());
484 if (it == nameMap.end()) { cout << "[ERROR]: " << current.getName() << " is not in your namefile, please correct." << endl; m->control_pressed = true; }
485 else { num = it->second; }
488 //for each sequence this sequence represents
489 for (int i = 0; i < num; i++) {
490 startPosition.push_back(current.getStartPos());
491 endPosition.push_back(current.getEndPos());
492 seqLength.push_back(current.getNumBases());
493 ambigBases.push_back(current.getAmbigBases());
494 longHomoPolymer.push_back(current.getLongHomoPolymer());
497 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
498 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\t" + toString(num) + "\n";
501 length = outputString.length();
502 char* buf3 = new char[length];
503 memcpy(buf3, outputString.c_str(), length);
505 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
512 catch(exception& e) {
513 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
518 /**************************************************************************************************/
519 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
521 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
526 //loop through and create all the processes you want
527 while (process != processors) {
531 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
534 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
536 //pass numSeqs to parent
538 string tempFile = fastafile + toString(getpid()) + ".num.temp";
539 m->openOutputFile(tempFile, out);
542 out << startPosition.size() << endl;
543 for (int k = 0; k < startPosition.size(); k++) { out << startPosition[k] << '\t'; } out << endl;
544 for (int k = 0; k < endPosition.size(); k++) { out << endPosition[k] << '\t'; } out << endl;
545 for (int k = 0; k < seqLength.size(); k++) { out << seqLength[k] << '\t'; } out << endl;
546 for (int k = 0; k < ambigBases.size(); k++) { out << ambigBases[k] << '\t'; } out << endl;
547 for (int k = 0; k < longHomoPolymer.size(); k++) { out << longHomoPolymer[k] << '\t'; } out << endl;
553 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
554 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
559 //force parent to wait until all the processes are done
560 for (int i=0;i<processors;i++) {
561 int temp = processIDS[i];
565 //parent reads in and combine Filter info
566 for (int i = 0; i < processIDS.size(); i++) {
567 string tempFilename = fastafile + toString(processIDS[i]) + ".num.temp";
569 m->openInputFile(tempFilename, in);
572 in >> tempNum; m->gobble(in); num += tempNum;
573 in >> tempNum; m->gobble(in);
574 for (int k = 0; k < tempNum; k++) { in >> temp; startPosition.push_back(temp); } m->gobble(in);
575 for (int k = 0; k < tempNum; k++) { in >> temp; endPosition.push_back(temp); } m->gobble(in);
576 for (int k = 0; k < tempNum; k++) { in >> temp; seqLength.push_back(temp); } m->gobble(in);
577 for (int k = 0; k < tempNum; k++) { in >> temp; ambigBases.push_back(temp); } m->gobble(in);
578 for (int k = 0; k < tempNum; k++) { in >> temp; longHomoPolymer.push_back(temp); } m->gobble(in);
581 remove(tempFilename.c_str());
587 catch(exception& e) {
588 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
592 /**********************************************************************************************************************/