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::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");
49 //**********************************************************************************************************************
50 SeqSummaryCommand::SeqSummaryCommand(){
52 abort = true; calledHelp = true;
54 vector<string> tempOutNames;
55 outputTypes["summary"] = tempOutNames;
58 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
62 //***************************************************************************************************************
64 SeqSummaryCommand::SeqSummaryCommand(string option) {
66 abort = false; calledHelp = false;
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;}
73 vector<string> myArray = setParameters();
75 OptionParser parser(option);
76 map<string,string> parameters = parser.getParameters();
78 ValidParameters validParameter("summary.seqs");
79 map<string,string>::iterator it;
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; }
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 = ""; }
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; }
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; }
108 //initialize outputTypes
109 vector<string> tempOutNames;
110 outputTypes["summary"] = tempOutNames;
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; }
121 namefile = validParameter.validFile(parameters, "name", true);
122 if (namefile == "not open") { namefile = ""; abort = true; }
123 else if (namefile == "not found") { namefile = ""; }
125 //if the user changes the output directory command factory will send this info to us in the output parameter
126 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
128 outputDir += m->hasPath(fastafile); //if user entered a file with a path then preserve it
131 string temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
132 m->setProcessors(temp);
133 convert(temp, processors);
138 catch(exception& e) {
139 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
143 //***************************************************************************************************************
145 int SeqSummaryCommand::execute(){
148 if (abort == true) { if (calledHelp) { return 0; } return 2; }
150 //set current fasta to fastafile
151 m->setFastaFile(fastafile);
153 string summaryFile = outputDir + m->getSimpleName(fastafile) + ".summary";
157 vector<int> startPosition;
158 vector<int> endPosition;
159 vector<int> seqLength;
160 vector<int> ambigBases;
161 vector<int> longHomoPolymer;
163 if (namefile != "") { nameMap = m->readNames(namefile); }
165 if (m->control_pressed) { return 0; }
168 int pid, numSeqsPerProcessor;
170 int startTag = 1; int endTag = 2; int lengthTag = 3; int baseTag = 4; int lhomoTag = 5;
171 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
172 vector<unsigned long int> MPIPos;
175 MPI_Status statusOut;
178 MPI_Comm_size(MPI_COMM_WORLD, &processors);
179 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
181 char tempFileName[1024];
182 strcpy(tempFileName, fastafile.c_str());
184 char sumFileName[1024];
185 strcpy(sumFileName, summaryFile.c_str());
187 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
188 MPI_File_open(MPI_COMM_WORLD, sumFileName, outMode, MPI_INFO_NULL, &outMPI);
190 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
192 if (pid == 0) { //you are the root process
194 string outputString = "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs\n";
195 int length = outputString.length();
196 char* buf2 = new char[length];
197 memcpy(buf2, outputString.c_str(), length);
199 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
202 MPIPos = m->setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
204 for(int i = 1; i < processors; i++) {
205 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
206 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
209 //figure out how many sequences you have to do
210 numSeqsPerProcessor = numSeqs / processors;
211 int startIndex = pid * numSeqsPerProcessor;
212 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
215 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
217 }else { //i am the child process
219 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
220 MPIPos.resize(numSeqs+1);
221 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
223 //figure out how many sequences you have to align
224 numSeqsPerProcessor = numSeqs / processors;
225 int startIndex = pid * numSeqsPerProcessor;
226 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
229 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
232 MPI_File_close(&inMPI);
233 MPI_File_close(&outMPI);
234 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
237 //get the info from the child processes
238 for(int i = 1; i < processors; i++) {
240 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
242 vector<int> temp; temp.resize(size+1);
244 for(int j = 0; j < 5; j++) {
246 MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
247 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what count this is for
249 if (receiveTag == startTag) {
250 for (int k = 0; k < size; k++) { startPosition.push_back(temp[k]); }
251 }else if (receiveTag == endTag) {
252 for (int k = 0; k < size; k++) { endPosition.push_back(temp[k]); }
253 }else if (receiveTag == lengthTag) {
254 for (int k = 0; k < size; k++) { seqLength.push_back(temp[k]); }
255 }else if (receiveTag == baseTag) {
256 for (int k = 0; k < size; k++) { ambigBases.push_back(temp[k]); }
257 }else if (receiveTag == lhomoTag) {
258 for (int k = 0; k < size; k++) { longHomoPolymer.push_back(temp[k]); }
266 int size = startPosition.size();
267 MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
269 startPosition.push_back(startTag);
270 int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
271 endPosition.push_back(endTag);
272 ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
273 seqLength.push_back(lengthTag);
274 ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
275 ambigBases.push_back(baseTag);
276 ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
277 longHomoPolymer.push_back(lhomoTag);
278 ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
281 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
283 vector<unsigned long int> positions = m->divideFile(fastafile, processors);
285 for (int i = 0; i < (positions.size()-1); i++) {
286 lines.push_back(new linePair(positions[i], positions[(i+1)]));
289 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
291 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
293 numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile);
295 rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
297 for(int i=1;i<processors;i++){
298 m->appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
299 remove((summaryFile + toString(processIDS[i]) + ".temp").c_str());
303 if (m->control_pressed) { return 0; }
305 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
306 if (m->control_pressed) { return 0; }
314 sort(startPosition.begin(), startPosition.end());
315 sort(endPosition.begin(), endPosition.end());
316 sort(seqLength.begin(), seqLength.end());
317 sort(ambigBases.begin(), ambigBases.end());
318 sort(longHomoPolymer.begin(), longHomoPolymer.end());
319 int size = startPosition.size();
321 int ptile0_25 = int(size * 0.025);
322 int ptile25 = int(size * 0.250);
323 int ptile50 = int(size * 0.500);
324 int ptile75 = int(size * 0.750);
325 int ptile97_5 = int(size * 0.975);
326 int ptile100 = size - 1;
328 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
329 if (startPosition[0] == -1) { startPosition[0] = 0; }
330 if (endPosition[0] == -1) { endPosition[0] = 0; }
332 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
334 m->mothurOutEndLine();
335 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
336 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();
337 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();
338 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();
339 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();
340 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();
341 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();
342 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();
343 if (namefile == "") { m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); }
344 else { m->mothurOut("# of unique seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); m->mothurOut("total # of seqs:\t" + toString(startPosition.size())); m->mothurOutEndLine(); }
346 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
348 m->mothurOutEndLine();
349 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
350 m->mothurOut(summaryFile); m->mothurOutEndLine(); outputNames.push_back(summaryFile); outputTypes["summary"].push_back(summaryFile);
351 m->mothurOutEndLine();
359 catch(exception& e) {
360 m->errorOut(e, "SeqSummaryCommand", "execute");
364 /**************************************************************************************/
365 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* filePos) {
369 m->openOutputFile(sumFile, outSummary);
371 //print header if you are process 0
372 if (filePos->start == 0) {
373 outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs" << endl;
377 m->openInputFile(filename, in);
379 in.seekg(filePos->start);
386 if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
388 Sequence current(in); m->gobble(in);
390 if (current.getName() != "") {
393 if (namefile != "") {
394 //make sure this sequence is in the namefile, else error
395 map<string, int>::iterator it = nameMap.find(current.getName());
397 if (it == nameMap.end()) { m->mothurOut("[ERROR]: " + current.getName() + " is not in your namefile, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
398 else { num = it->second; }
401 //for each sequence this sequence represents
402 for (int i = 0; i < num; i++) {
403 startPosition.push_back(current.getStartPos());
404 endPosition.push_back(current.getEndPos());
405 seqLength.push_back(current.getNumBases());
406 ambigBases.push_back(current.getAmbigBases());
407 longHomoPolymer.push_back(current.getLongHomoPolymer());
411 outSummary << current.getName() << '\t';
412 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
413 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
414 outSummary << current.getLongHomoPolymer() << '\t' << num << endl;
417 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
418 unsigned long int pos = in.tellg();
419 if ((pos == -1) || (pos >= filePos->end)) { break; }
421 if (in.eof()) { break; }
425 //if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
428 //if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
434 catch(exception& e) {
435 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
440 /**************************************************************************************/
441 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) {
446 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
448 for(int i=0;i<num;i++){
450 if (m->control_pressed) { return 0; }
453 int length = MPIPos[start+i+1] - MPIPos[start+i];
455 char* buf4 = new char[length];
456 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
458 string tempBuf = buf4;
459 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
460 istringstream iss (tempBuf,istringstream::in);
463 Sequence current(iss);
465 if (current.getName() != "") {
468 if (namefile != "") {
469 //make sure this sequence is in the namefile, else error
470 map<string, int>::iterator it = nameMap.find(current.getName());
472 if (it == nameMap.end()) { cout << "[ERROR]: " << current.getName() << " is not in your namefile, please correct." << endl; m->control_pressed = true; }
473 else { num = it->second; }
476 //for each sequence this sequence represents
477 for (int i = 0; i < num; i++) {
478 startPosition.push_back(current.getStartPos());
479 endPosition.push_back(current.getEndPos());
480 seqLength.push_back(current.getNumBases());
481 ambigBases.push_back(current.getAmbigBases());
482 longHomoPolymer.push_back(current.getLongHomoPolymer());
485 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
486 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\t" + toString(num) + "\n";
489 length = outputString.length();
490 char* buf3 = new char[length];
491 memcpy(buf3, outputString.c_str(), length);
493 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
500 catch(exception& e) {
501 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
506 /**************************************************************************************************/
507 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
509 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
514 //loop through and create all the processes you want
515 while (process != processors) {
519 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
522 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
524 //pass numSeqs to parent
526 string tempFile = fastafile + toString(getpid()) + ".num.temp";
527 m->openOutputFile(tempFile, out);
530 out << startPosition.size() << endl;
531 for (int k = 0; k < startPosition.size(); k++) { out << startPosition[k] << '\t'; } out << endl;
532 for (int k = 0; k < endPosition.size(); k++) { out << endPosition[k] << '\t'; } out << endl;
533 for (int k = 0; k < seqLength.size(); k++) { out << seqLength[k] << '\t'; } out << endl;
534 for (int k = 0; k < ambigBases.size(); k++) { out << ambigBases[k] << '\t'; } out << endl;
535 for (int k = 0; k < longHomoPolymer.size(); k++) { out << longHomoPolymer[k] << '\t'; } out << endl;
541 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
542 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
547 //force parent to wait until all the processes are done
548 for (int i=0;i<processors;i++) {
549 int temp = processIDS[i];
553 //parent reads in and combine Filter info
554 for (int i = 0; i < processIDS.size(); i++) {
555 string tempFilename = fastafile + toString(processIDS[i]) + ".num.temp";
557 m->openInputFile(tempFilename, in);
560 in >> tempNum; m->gobble(in); num += tempNum;
561 in >> tempNum; m->gobble(in);
562 for (int k = 0; k < tempNum; k++) { in >> temp; startPosition.push_back(temp); } m->gobble(in);
563 for (int k = 0; k < tempNum; k++) { in >> temp; endPosition.push_back(temp); } m->gobble(in);
564 for (int k = 0; k < tempNum; k++) { in >> temp; seqLength.push_back(temp); } m->gobble(in);
565 for (int k = 0; k < tempNum; k++) { in >> temp; ambigBases.push_back(temp); } m->gobble(in);
566 for (int k = 0; k < tempNum; k++) { in >> temp; longHomoPolymer.push_back(temp); } m->gobble(in);
569 remove(tempFilename.c_str());
575 catch(exception& e) {
576 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
580 /**********************************************************************************************************************/