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 pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
19 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
21 vector<string> myArray;
22 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
26 m->errorOut(e, "SeqSummaryCommand", "setParameters");
30 //**********************************************************************************************************************
31 string SeqSummaryCommand::getHelpString(){
33 string helpString = "";
34 helpString += "The summary.seqs command reads a fastafile and summarizes the sequences.\n";
35 helpString += "The summary.seqs command parameters are fasta, name and processors, fasta is required, unless you have a valid current fasta file.\n";
36 helpString += "The name parameter allows you to enter a name file associated with your fasta file. \n";
37 helpString += "The summary.seqs command should be in the following format: \n";
38 helpString += "summary.seqs(fasta=yourFastaFile, processors=2) \n";
39 helpString += "Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n";
43 m->errorOut(e, "SeqSummaryCommand", "getHelpString");
48 //**********************************************************************************************************************
49 SeqSummaryCommand::SeqSummaryCommand(){
51 abort = true; calledHelp = true;
53 vector<string> tempOutNames;
54 outputTypes["summary"] = tempOutNames;
57 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
61 //***************************************************************************************************************
63 SeqSummaryCommand::SeqSummaryCommand(string option) {
65 abort = false; calledHelp = false;
67 //allow user to run help
68 if(option == "help") { help(); abort = true; calledHelp = true; }
69 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
72 vector<string> myArray = setParameters();
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") {
115 fastafile = m->getFastaFile();
116 if (fastafile != "") { m->mothurOut("Using " + fastafile + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
117 else { m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
120 namefile = validParameter.validFile(parameters, "name", true);
121 if (namefile == "not open") { namefile = ""; abort = true; }
122 else if (namefile == "not found") { namefile = ""; }
124 //if the user changes the output directory command factory will send this info to us in the output parameter
125 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
127 outputDir += m->hasPath(fastafile); //if user entered a file with a path then preserve it
130 string temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
131 m->setProcessors(temp);
132 convert(temp, processors);
137 catch(exception& e) {
138 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
142 //***************************************************************************************************************
144 int SeqSummaryCommand::execute(){
147 if (abort == true) { if (calledHelp) { return 0; } return 2; }
149 //set current fasta to fastafile
150 m->setFastaFile(fastafile);
152 string summaryFile = outputDir + m->getSimpleName(fastafile) + ".summary";
156 vector<int> startPosition;
157 vector<int> endPosition;
158 vector<int> seqLength;
159 vector<int> ambigBases;
160 vector<int> longHomoPolymer;
162 if (namefile != "") { nameMap = m->readNames(namefile); }
164 if (m->control_pressed) { return 0; }
167 int pid, numSeqsPerProcessor;
169 int startTag = 1; int endTag = 2; int lengthTag = 3; int baseTag = 4; int lhomoTag = 5;
170 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
171 vector<unsigned long int> MPIPos;
174 MPI_Status statusOut;
177 MPI_Comm_size(MPI_COMM_WORLD, &processors);
178 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
180 char tempFileName[1024];
181 strcpy(tempFileName, fastafile.c_str());
183 char sumFileName[1024];
184 strcpy(sumFileName, summaryFile.c_str());
186 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
187 MPI_File_open(MPI_COMM_WORLD, sumFileName, outMode, MPI_INFO_NULL, &outMPI);
189 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
191 if (pid == 0) { //you are the root process
193 string outputString = "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs\n";
194 int length = outputString.length();
195 char* buf2 = new char[length];
196 memcpy(buf2, outputString.c_str(), length);
198 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
201 MPIPos = m->setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
203 for(int i = 1; i < processors; i++) {
204 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
205 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
208 //figure out how many sequences you have to do
209 numSeqsPerProcessor = numSeqs / processors;
210 int startIndex = pid * numSeqsPerProcessor;
211 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
214 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
216 }else { //i am the child process
218 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
219 MPIPos.resize(numSeqs+1);
220 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
222 //figure out how many sequences you have to align
223 numSeqsPerProcessor = numSeqs / processors;
224 int startIndex = pid * numSeqsPerProcessor;
225 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
228 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
231 MPI_File_close(&inMPI);
232 MPI_File_close(&outMPI);
233 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
236 //get the info from the child processes
237 for(int i = 1; i < processors; i++) {
239 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
241 vector<int> temp; temp.resize(size+1);
243 for(int j = 0; j < 5; j++) {
245 MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
246 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what count this is for
248 if (receiveTag == startTag) {
249 for (int k = 0; k < size; k++) { startPosition.push_back(temp[k]); }
250 }else if (receiveTag == endTag) {
251 for (int k = 0; k < size; k++) { endPosition.push_back(temp[k]); }
252 }else if (receiveTag == lengthTag) {
253 for (int k = 0; k < size; k++) { seqLength.push_back(temp[k]); }
254 }else if (receiveTag == baseTag) {
255 for (int k = 0; k < size; k++) { ambigBases.push_back(temp[k]); }
256 }else if (receiveTag == lhomoTag) {
257 for (int k = 0; k < size; k++) { longHomoPolymer.push_back(temp[k]); }
265 int size = startPosition.size();
266 MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
268 startPosition.push_back(startTag);
269 int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
270 endPosition.push_back(endTag);
271 ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
272 seqLength.push_back(lengthTag);
273 ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
274 ambigBases.push_back(baseTag);
275 ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
276 longHomoPolymer.push_back(lhomoTag);
277 ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
280 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
282 vector<unsigned long int> positions = m->divideFile(fastafile, processors);
284 for (int i = 0; i < (positions.size()-1); i++) {
285 lines.push_back(new linePair(positions[i], positions[(i+1)]));
288 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
290 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
292 numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile);
294 rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
296 for(int i=1;i<processors;i++){
297 m->appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
298 remove((summaryFile + toString(processIDS[i]) + ".temp").c_str());
302 if (m->control_pressed) { return 0; }
304 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
305 if (m->control_pressed) { return 0; }
313 sort(startPosition.begin(), startPosition.end());
314 sort(endPosition.begin(), endPosition.end());
315 sort(seqLength.begin(), seqLength.end());
316 sort(ambigBases.begin(), ambigBases.end());
317 sort(longHomoPolymer.begin(), longHomoPolymer.end());
318 int size = startPosition.size();
320 int ptile0_25 = int(size * 0.025);
321 int ptile25 = int(size * 0.250);
322 int ptile50 = int(size * 0.500);
323 int ptile75 = int(size * 0.750);
324 int ptile97_5 = int(size * 0.975);
325 int ptile100 = size - 1;
327 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
328 if (startPosition[0] == -1) { startPosition[0] = 0; }
329 if (endPosition[0] == -1) { endPosition[0] = 0; }
331 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
333 m->mothurOutEndLine();
334 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
335 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();
336 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();
337 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();
338 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();
339 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();
340 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();
341 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();
342 if (namefile == "") { m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); }
343 else { m->mothurOut("# of unique seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); m->mothurOut("total # of seqs:\t" + toString(startPosition.size())); m->mothurOutEndLine(); }
345 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
347 m->mothurOutEndLine();
348 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
349 m->mothurOut(summaryFile); m->mothurOutEndLine(); outputNames.push_back(summaryFile); outputTypes["summary"].push_back(summaryFile);
350 m->mothurOutEndLine();
358 catch(exception& e) {
359 m->errorOut(e, "SeqSummaryCommand", "execute");
363 /**************************************************************************************/
364 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* filePos) {
368 m->openOutputFile(sumFile, outSummary);
370 //print header if you are process 0
371 if (filePos->start == 0) {
372 outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs" << endl;
376 m->openInputFile(filename, in);
378 in.seekg(filePos->start);
385 if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
387 Sequence current(in); m->gobble(in);
389 if (current.getName() != "") {
392 if (namefile != "") {
393 //make sure this sequence is in the namefile, else error
394 map<string, int>::iterator it = nameMap.find(current.getName());
396 if (it == nameMap.end()) { m->mothurOut("[ERROR]: " + current.getName() + " is not in your namefile, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
397 else { num = it->second; }
400 //for each sequence this sequence represents
401 for (int i = 0; i < num; i++) {
402 startPosition.push_back(current.getStartPos());
403 endPosition.push_back(current.getEndPos());
404 seqLength.push_back(current.getNumBases());
405 ambigBases.push_back(current.getAmbigBases());
406 longHomoPolymer.push_back(current.getLongHomoPolymer());
410 outSummary << current.getName() << '\t';
411 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
412 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
413 outSummary << current.getLongHomoPolymer() << '\t' << num << endl;
416 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
417 unsigned long int pos = in.tellg();
418 if ((pos == -1) || (pos >= filePos->end)) { break; }
420 if (in.eof()) { break; }
424 //if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
427 //if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
433 catch(exception& e) {
434 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
439 /**************************************************************************************/
440 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) {
445 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
447 for(int i=0;i<num;i++){
449 if (m->control_pressed) { return 0; }
452 int length = MPIPos[start+i+1] - MPIPos[start+i];
454 char* buf4 = new char[length];
455 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
457 string tempBuf = buf4;
458 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
459 istringstream iss (tempBuf,istringstream::in);
462 Sequence current(iss);
464 if (current.getName() != "") {
467 if (namefile != "") {
468 //make sure this sequence is in the namefile, else error
469 map<string, int>::iterator it = nameMap.find(current.getName());
471 if (it == nameMap.end()) { cout << "[ERROR]: " << current.getName() << " is not in your namefile, please correct." << endl; m->control_pressed = true; }
472 else { num = it->second; }
475 //for each sequence this sequence represents
476 for (int i = 0; i < num; i++) {
477 startPosition.push_back(current.getStartPos());
478 endPosition.push_back(current.getEndPos());
479 seqLength.push_back(current.getNumBases());
480 ambigBases.push_back(current.getAmbigBases());
481 longHomoPolymer.push_back(current.getLongHomoPolymer());
484 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
485 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\t" + toString(num) + "\n";
488 length = outputString.length();
489 char* buf3 = new char[length];
490 memcpy(buf3, outputString.c_str(), length);
492 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
499 catch(exception& e) {
500 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
505 /**************************************************************************************************/
506 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
508 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
513 //loop through and create all the processes you want
514 while (process != processors) {
518 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
521 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
523 //pass numSeqs to parent
525 string tempFile = fastafile + toString(getpid()) + ".num.temp";
526 m->openOutputFile(tempFile, out);
529 out << startPosition.size() << endl;
530 for (int k = 0; k < startPosition.size(); k++) { out << startPosition[k] << '\t'; } out << endl;
531 for (int k = 0; k < endPosition.size(); k++) { out << endPosition[k] << '\t'; } out << endl;
532 for (int k = 0; k < seqLength.size(); k++) { out << seqLength[k] << '\t'; } out << endl;
533 for (int k = 0; k < ambigBases.size(); k++) { out << ambigBases[k] << '\t'; } out << endl;
534 for (int k = 0; k < longHomoPolymer.size(); k++) { out << longHomoPolymer[k] << '\t'; } out << endl;
540 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
541 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
546 //force parent to wait until all the processes are done
547 for (int i=0;i<processors;i++) {
548 int temp = processIDS[i];
552 //parent reads in and combine Filter info
553 for (int i = 0; i < processIDS.size(); i++) {
554 string tempFilename = fastafile + toString(processIDS[i]) + ".num.temp";
556 m->openInputFile(tempFilename, in);
559 in >> tempNum; m->gobble(in); num += tempNum;
560 in >> tempNum; m->gobble(in);
561 for (int k = 0; k < tempNum; k++) { in >> temp; startPosition.push_back(temp); } m->gobble(in);
562 for (int k = 0; k < tempNum; k++) { in >> temp; endPosition.push_back(temp); } m->gobble(in);
563 for (int k = 0; k < tempNum; k++) { in >> temp; seqLength.push_back(temp); } m->gobble(in);
564 for (int k = 0; k < tempNum; k++) { in >> temp; ambigBases.push_back(temp); } m->gobble(in);
565 for (int k = 0; k < tempNum; k++) { in >> temp; longHomoPolymer.push_back(temp); } m->gobble(in);
568 remove(tempFilename.c_str());
574 catch(exception& e) {
575 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
579 /**********************************************************************************************************************/