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; }
119 }else { m->setFastaFile(fastafile); }
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); }
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"){
129 outputDir += m->hasPath(fastafile); //if user entered a file with a path then preserve it
132 string temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
133 m->setProcessors(temp);
134 convert(temp, processors);
139 catch(exception& e) {
140 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
144 //***************************************************************************************************************
146 int SeqSummaryCommand::execute(){
149 if (abort == true) { if (calledHelp) { return 0; } return 2; }
151 //set current fasta to fastafile
152 m->setFastaFile(fastafile);
154 string summaryFile = outputDir + m->getSimpleName(fastafile) + ".summary";
158 vector<int> startPosition;
159 vector<int> endPosition;
160 vector<int> seqLength;
161 vector<int> ambigBases;
162 vector<int> longHomoPolymer;
164 if (namefile != "") { nameMap = m->readNames(namefile); }
166 if (m->control_pressed) { return 0; }
169 int pid, numSeqsPerProcessor;
171 int startTag = 1; int endTag = 2; int lengthTag = 3; int baseTag = 4; int lhomoTag = 5;
172 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
173 vector<unsigned long int> MPIPos;
176 MPI_Status statusOut;
179 MPI_Comm_size(MPI_COMM_WORLD, &processors);
180 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
182 char tempFileName[1024];
183 strcpy(tempFileName, fastafile.c_str());
185 char sumFileName[1024];
186 strcpy(sumFileName, summaryFile.c_str());
188 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
189 MPI_File_open(MPI_COMM_WORLD, sumFileName, outMode, MPI_INFO_NULL, &outMPI);
191 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
193 if (pid == 0) { //you are the root process
195 string outputString = "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs\n";
196 int length = outputString.length();
197 char* buf2 = new char[length];
198 memcpy(buf2, outputString.c_str(), length);
200 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
203 MPIPos = m->setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
205 for(int i = 1; i < processors; i++) {
206 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
207 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
210 //figure out how many sequences you have to do
211 numSeqsPerProcessor = numSeqs / processors;
212 int startIndex = pid * numSeqsPerProcessor;
213 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
216 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
218 }else { //i am the child process
220 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
221 MPIPos.resize(numSeqs+1);
222 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
224 //figure out how many sequences you have to align
225 numSeqsPerProcessor = numSeqs / processors;
226 int startIndex = pid * numSeqsPerProcessor;
227 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
230 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
233 MPI_File_close(&inMPI);
234 MPI_File_close(&outMPI);
235 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
238 //get the info from the child processes
239 for(int i = 1; i < processors; i++) {
241 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
243 vector<int> temp; temp.resize(size+1);
245 for(int j = 0; j < 5; j++) {
247 MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
248 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what count this is for
250 if (receiveTag == startTag) {
251 for (int k = 0; k < size; k++) { startPosition.push_back(temp[k]); }
252 }else if (receiveTag == endTag) {
253 for (int k = 0; k < size; k++) { endPosition.push_back(temp[k]); }
254 }else if (receiveTag == lengthTag) {
255 for (int k = 0; k < size; k++) { seqLength.push_back(temp[k]); }
256 }else if (receiveTag == baseTag) {
257 for (int k = 0; k < size; k++) { ambigBases.push_back(temp[k]); }
258 }else if (receiveTag == lhomoTag) {
259 for (int k = 0; k < size; k++) { longHomoPolymer.push_back(temp[k]); }
267 int size = startPosition.size();
268 MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
270 startPosition.push_back(startTag);
271 int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
272 endPosition.push_back(endTag);
273 ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
274 seqLength.push_back(lengthTag);
275 ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
276 ambigBases.push_back(baseTag);
277 ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
278 longHomoPolymer.push_back(lhomoTag);
279 ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
282 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
284 vector<unsigned long int> positions = m->divideFile(fastafile, processors);
286 for (int i = 0; i < (positions.size()-1); i++) {
287 lines.push_back(new linePair(positions[i], positions[(i+1)]));
290 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
292 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
294 numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile);
296 rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
298 for(int i=1;i<processors;i++){
299 m->appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
300 m->mothurRemove((summaryFile + toString(processIDS[i]) + ".temp"));
304 if (m->control_pressed) { return 0; }
306 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
307 if (m->control_pressed) { return 0; }
315 sort(startPosition.begin(), startPosition.end());
316 sort(endPosition.begin(), endPosition.end());
317 sort(seqLength.begin(), seqLength.end());
318 sort(ambigBases.begin(), ambigBases.end());
319 sort(longHomoPolymer.begin(), longHomoPolymer.end());
320 int size = startPosition.size();
322 int ptile0_25 = int(size * 0.025);
323 int ptile25 = int(size * 0.250);
324 int ptile50 = int(size * 0.500);
325 int ptile75 = int(size * 0.750);
326 int ptile97_5 = int(size * 0.975);
327 int ptile100 = size - 1;
329 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
330 if (startPosition[0] == -1) { startPosition[0] = 0; }
331 if (endPosition[0] == -1) { endPosition[0] = 0; }
333 if (m->control_pressed) { m->mothurRemove(summaryFile); return 0; }
335 m->mothurOutEndLine();
336 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
337 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();
338 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();
339 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();
340 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();
341 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();
342 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();
343 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();
344 if (namefile == "") { m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); }
345 else { m->mothurOut("# of unique seqs:\t" + toString(numSeqs)); m->mothurOutEndLine(); m->mothurOut("total # of seqs:\t" + toString(startPosition.size())); m->mothurOutEndLine(); }
347 if (m->control_pressed) { m->mothurRemove(summaryFile); return 0; }
349 m->mothurOutEndLine();
350 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
351 m->mothurOut(summaryFile); m->mothurOutEndLine(); outputNames.push_back(summaryFile); outputTypes["summary"].push_back(summaryFile);
352 m->mothurOutEndLine();
360 catch(exception& e) {
361 m->errorOut(e, "SeqSummaryCommand", "execute");
365 /**************************************************************************************/
366 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* filePos) {
370 m->openOutputFile(sumFile, outSummary);
372 //print header if you are process 0
373 if (filePos->start == 0) {
374 outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer\tnumSeqs" << endl;
378 m->openInputFile(filename, in);
380 in.seekg(filePos->start);
387 if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
389 Sequence current(in); m->gobble(in);
391 if (current.getName() != "") {
394 if (namefile != "") {
395 //make sure this sequence is in the namefile, else error
396 map<string, int>::iterator it = nameMap.find(current.getName());
398 if (it == nameMap.end()) { m->mothurOut("[ERROR]: " + current.getName() + " is not in your namefile, please correct."); m->mothurOutEndLine(); m->control_pressed = true; }
399 else { num = it->second; }
402 //for each sequence this sequence represents
403 for (int i = 0; i < num; i++) {
404 startPosition.push_back(current.getStartPos());
405 endPosition.push_back(current.getEndPos());
406 seqLength.push_back(current.getNumBases());
407 ambigBases.push_back(current.getAmbigBases());
408 longHomoPolymer.push_back(current.getLongHomoPolymer());
412 outSummary << current.getName() << '\t';
413 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
414 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
415 outSummary << current.getLongHomoPolymer() << '\t' << num << endl;
418 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
419 unsigned long int pos = in.tellg();
420 if ((pos == -1) || (pos >= filePos->end)) { break; }
422 if (in.eof()) { break; }
426 //if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
429 //if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
435 catch(exception& e) {
436 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
441 /**************************************************************************************/
442 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) {
447 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
449 for(int i=0;i<num;i++){
451 if (m->control_pressed) { return 0; }
454 int length = MPIPos[start+i+1] - MPIPos[start+i];
456 char* buf4 = new char[length];
457 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
459 string tempBuf = buf4;
460 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
461 istringstream iss (tempBuf,istringstream::in);
464 Sequence current(iss);
466 if (current.getName() != "") {
469 if (namefile != "") {
470 //make sure this sequence is in the namefile, else error
471 map<string, int>::iterator it = nameMap.find(current.getName());
473 if (it == nameMap.end()) { cout << "[ERROR]: " << current.getName() << " is not in your namefile, please correct." << endl; m->control_pressed = true; }
474 else { num = it->second; }
477 //for each sequence this sequence represents
478 for (int i = 0; i < num; i++) {
479 startPosition.push_back(current.getStartPos());
480 endPosition.push_back(current.getEndPos());
481 seqLength.push_back(current.getNumBases());
482 ambigBases.push_back(current.getAmbigBases());
483 longHomoPolymer.push_back(current.getLongHomoPolymer());
486 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
487 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\t" + toString(num) + "\n";
490 length = outputString.length();
491 char* buf3 = new char[length];
492 memcpy(buf3, outputString.c_str(), length);
494 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
501 catch(exception& e) {
502 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
507 /**************************************************************************************************/
508 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
510 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
515 //loop through and create all the processes you want
516 while (process != processors) {
520 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
523 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
525 //pass numSeqs to parent
527 string tempFile = fastafile + toString(getpid()) + ".num.temp";
528 m->openOutputFile(tempFile, out);
531 out << startPosition.size() << endl;
532 for (int k = 0; k < startPosition.size(); k++) { out << startPosition[k] << '\t'; } out << endl;
533 for (int k = 0; k < endPosition.size(); k++) { out << endPosition[k] << '\t'; } out << endl;
534 for (int k = 0; k < seqLength.size(); k++) { out << seqLength[k] << '\t'; } out << endl;
535 for (int k = 0; k < ambigBases.size(); k++) { out << ambigBases[k] << '\t'; } out << endl;
536 for (int k = 0; k < longHomoPolymer.size(); k++) { out << longHomoPolymer[k] << '\t'; } out << endl;
542 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
543 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
548 //force parent to wait until all the processes are done
549 for (int i=0;i<processors;i++) {
550 int temp = processIDS[i];
554 //parent reads in and combine Filter info
555 for (int i = 0; i < processIDS.size(); i++) {
556 string tempFilename = fastafile + toString(processIDS[i]) + ".num.temp";
558 m->openInputFile(tempFilename, in);
561 in >> tempNum; m->gobble(in); num += tempNum;
562 in >> tempNum; m->gobble(in);
563 for (int k = 0; k < tempNum; k++) { in >> temp; startPosition.push_back(temp); } m->gobble(in);
564 for (int k = 0; k < tempNum; k++) { in >> temp; endPosition.push_back(temp); } m->gobble(in);
565 for (int k = 0; k < tempNum; k++) { in >> temp; seqLength.push_back(temp); } m->gobble(in);
566 for (int k = 0; k < tempNum; k++) { in >> temp; ambigBases.push_back(temp); } m->gobble(in);
567 for (int k = 0; k < tempNum; k++) { in >> temp; longHomoPolymer.push_back(temp); } m->gobble(in);
570 m->mothurRemove(tempFilename);
576 catch(exception& e) {
577 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
581 /**********************************************************************************************************************/