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 //***************************************************************************************************************
15 SeqSummaryCommand::SeqSummaryCommand(string option) {
19 //allow user to run help
20 if(option == "help") { help(); abort = true; }
23 //valid paramters for this command
24 string Array[] = {"fasta","processors","outputdir","inputdir"};
25 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
27 OptionParser parser(option);
28 map<string,string> parameters = parser.getParameters();
30 ValidParameters validParameter("summary.seqs");
31 map<string,string>::iterator it;
33 //check to make sure all parameters are valid for command
34 for (it = parameters.begin(); it != parameters.end(); it++) {
35 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
38 //if the user changes the input directory command factory will send this info to us in the output parameter
39 string inputDir = validParameter.validFile(parameters, "inputdir", false);
40 if (inputDir == "not found"){ inputDir = ""; }
43 it = parameters.find("fasta");
44 //user has given a template file
45 if(it != parameters.end()){
46 path = hasPath(it->second);
47 //if the user has not given a path then, add inputdir. else leave path alone.
48 if (path == "") { parameters["fasta"] = inputDir + it->second; }
52 //check for required parameters
53 fastafile = validParameter.validFile(parameters, "fasta", true);
54 if (fastafile == "not open") { abort = true; }
55 else if (fastafile == "not found") { fastafile = ""; m->mothurOut("fasta is a required parameter for the summary.seqs command."); m->mothurOutEndLine(); abort = true; }
57 //if the user changes the output directory command factory will send this info to us in the output parameter
58 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
60 outputDir += hasPath(fastafile); //if user entered a file with a path then preserve it
63 string temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
64 convert(temp, processors);
70 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
74 //**********************************************************************************************************************
76 void SeqSummaryCommand::help(){
78 m->mothurOut("The summary.seqs command reads a fastafile and summarizes the sequences.\n");
79 m->mothurOut("The summary.seqs command parameters are fasta and processors, fasta is required.\n");
80 m->mothurOut("The summary.seqs command should be in the following format: \n");
81 m->mothurOut("summary.seqs(fasta=yourFastaFile, processors=2) \n");
82 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");
85 m->errorOut(e, "SeqSummaryCommand", "help");
90 //***************************************************************************************************************
92 SeqSummaryCommand::~SeqSummaryCommand(){ /* do nothing */ }
94 //***************************************************************************************************************
96 int SeqSummaryCommand::execute(){
99 if (abort == true) { return 0; }
101 string summaryFile = outputDir + getSimpleName(fastafile) + ".summary";
105 vector<int> startPosition;
106 vector<int> endPosition;
107 vector<int> seqLength;
108 vector<int> ambigBases;
109 vector<int> longHomoPolymer;
112 int pid, numSeqsPerProcessor;
114 int startTag = 1; int endTag = 2; int lengthTag = 3; int baseTag = 4; int lhomoTag = 5;
115 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
116 vector<unsigned long int> MPIPos;
119 MPI_Status statusOut;
122 MPI_Comm_size(MPI_COMM_WORLD, &processors);
123 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
125 char tempFileName[1024];
126 strcpy(tempFileName, fastafile.c_str());
128 char sumFileName[1024];
129 strcpy(sumFileName, summaryFile.c_str());
131 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
132 MPI_File_open(MPI_COMM_WORLD, sumFileName, outMode, MPI_INFO_NULL, &outMPI);
134 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
136 if (pid == 0) { //you are the root process
138 string outputString = "seqname\tstart\tend\tnbases\tambigs\tpolymer\n";
139 int length = outputString.length();
140 char* buf2 = new char[length];
141 memcpy(buf2, outputString.c_str(), length);
143 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
146 MPIPos = setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
148 for(int i = 1; i < processors; i++) {
149 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
150 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
153 //figure out how many sequences you have to do
154 numSeqsPerProcessor = numSeqs / processors;
155 int startIndex = pid * numSeqsPerProcessor;
156 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
159 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
161 }else { //i am the child process
163 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
164 MPIPos.resize(numSeqs+1);
165 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
167 //figure out how many sequences you have to align
168 numSeqsPerProcessor = numSeqs / processors;
169 int startIndex = pid * numSeqsPerProcessor;
170 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
173 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
176 MPI_File_close(&inMPI);
177 MPI_File_close(&outMPI);
178 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
181 //get the info from the child processes
182 for(int i = 1; i < processors; i++) {
184 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
186 vector<int> temp; temp.resize(size+1);
188 for(int j = 0; j < 5; j++) {
190 MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
191 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what count this is for
193 if (receiveTag == startTag) {
194 for (int k = 0; k < size; k++) { startPosition.push_back(temp[k]); }
195 }else if (receiveTag == endTag) {
196 for (int k = 0; k < size; k++) { endPosition.push_back(temp[k]); }
197 }else if (receiveTag == lengthTag) {
198 for (int k = 0; k < size; k++) { seqLength.push_back(temp[k]); }
199 }else if (receiveTag == baseTag) {
200 for (int k = 0; k < size; k++) { ambigBases.push_back(temp[k]); }
201 }else if (receiveTag == lhomoTag) {
202 for (int k = 0; k < size; k++) { longHomoPolymer.push_back(temp[k]); }
210 int size = startPosition.size();
211 MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
213 startPosition.push_back(startTag);
214 int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
215 endPosition.push_back(endTag);
216 ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
217 seqLength.push_back(lengthTag);
218 ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
219 ambigBases.push_back(baseTag);
220 ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
221 longHomoPolymer.push_back(lhomoTag);
222 ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
225 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
227 vector<unsigned long int> positions = divideFile(fastafile, processors);
229 for (int i = 0; i < (positions.size()-1); i++) {
230 lines.push_back(new linePair(positions[i], positions[(i+1)]));
233 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
235 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
237 numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile);
239 rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
241 for(int i=1;i<processors;i++){
242 appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
243 remove((summaryFile + toString(processIDS[i]) + ".temp").c_str());
247 if (m->control_pressed) { return 0; }
249 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
250 if (m->control_pressed) { return 0; }
258 sort(startPosition.begin(), startPosition.end());
259 sort(endPosition.begin(), endPosition.end());
260 sort(seqLength.begin(), seqLength.end());
261 sort(ambigBases.begin(), ambigBases.end());
262 sort(longHomoPolymer.begin(), longHomoPolymer.end());
264 int ptile0_25 = int(numSeqs * 0.025);
265 int ptile25 = int(numSeqs * 0.250);
266 int ptile50 = int(numSeqs * 0.500);
267 int ptile75 = int(numSeqs * 0.750);
268 int ptile97_5 = int(numSeqs * 0.975);
269 int ptile100 = numSeqs - 1;
271 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
272 if (startPosition[0] == -1) { startPosition[0] = 0; }
273 if (endPosition[0] == -1) { endPosition[0] = 0; }
275 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
277 m->mothurOutEndLine();
278 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
279 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();
280 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();
281 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();
282 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();
283 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();
284 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();
285 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();
286 m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine();
288 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
290 m->mothurOutEndLine();
291 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
292 m->mothurOut(summaryFile); m->mothurOutEndLine();
293 m->mothurOutEndLine();
301 catch(exception& e) {
302 m->errorOut(e, "SeqSummaryCommand", "execute");
306 /**************************************************************************************/
307 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* filePos) {
311 openOutputFile(sumFile, outSummary);
313 //print header if you are process 0
314 if (filePos->start == 0) {
315 outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer" << endl;
319 openInputFile(filename, in);
321 in.seekg(filePos->start);
328 if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
330 Sequence current(in); gobble(in);
332 if (current.getName() != "") {
333 startPosition.push_back(current.getStartPos());
334 endPosition.push_back(current.getEndPos());
335 seqLength.push_back(current.getNumBases());
336 ambigBases.push_back(current.getAmbigBases());
337 longHomoPolymer.push_back(current.getLongHomoPolymer());
339 outSummary << current.getName() << '\t';
340 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
341 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
342 outSummary << current.getLongHomoPolymer() << endl;
346 unsigned long int pos = in.tellg();
347 if ((pos == -1) || (pos >= filePos->end)) { break; }
350 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
353 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
359 catch(exception& e) {
360 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
365 /**************************************************************************************/
366 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) {
371 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
374 for(int i=0;i<num;i++){
376 if (m->control_pressed) { return 0; }
379 int length = MPIPos[start+i+1] - MPIPos[start+i];
381 char* buf4 = new char[length];
382 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
384 string tempBuf = buf4;
385 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
386 istringstream iss (tempBuf,istringstream::in);
389 Sequence current(iss);
391 if (current.getName() != "") {
392 startPosition.push_back(current.getStartPos());
393 endPosition.push_back(current.getEndPos());
394 seqLength.push_back(current.getNumBases());
395 ambigBases.push_back(current.getAmbigBases());
396 longHomoPolymer.push_back(current.getLongHomoPolymer());
398 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
399 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\n";
402 length = outputString.length();
403 char* buf3 = new char[length];
404 memcpy(buf3, outputString.c_str(), length);
406 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
413 catch(exception& e) {
414 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
419 /**************************************************************************************************/
420 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
422 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
427 //loop through and create all the processes you want
428 while (process != processors) {
432 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
435 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
437 //pass numSeqs to parent
439 string tempFile = toString(getpid()) + ".temp";
440 openOutputFile(tempFile, out);
443 for (int k = 0; k < startPosition.size(); k++) { out << startPosition[k] << '\t'; } out << endl;
444 for (int k = 0; k < endPosition.size(); k++) { out << endPosition[k] << '\t'; } out << endl;
445 for (int k = 0; k < seqLength.size(); k++) { out << seqLength[k] << '\t'; } out << endl;
446 for (int k = 0; k < ambigBases.size(); k++) { out << ambigBases[k] << '\t'; } out << endl;
447 for (int k = 0; k < longHomoPolymer.size(); k++) { out << longHomoPolymer[k] << '\t'; } out << endl;
452 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
455 //force parent to wait until all the processes are done
456 for (int i=0;i<processors;i++) {
457 int temp = processIDS[i];
461 //parent reads in and combine Filter info
462 for (int i = 0; i < processIDS.size(); i++) {
463 string tempFilename = toString(processIDS[i]) + ".temp";
465 openInputFile(tempFilename, in);
468 in >> tempNum; gobble(in); num += tempNum;
469 for (int k = 0; k < tempNum; k++) { in >> temp; startPosition.push_back(temp); } gobble(in);
470 for (int k = 0; k < tempNum; k++) { in >> temp; endPosition.push_back(temp); } gobble(in);
471 for (int k = 0; k < tempNum; k++) { in >> temp; seqLength.push_back(temp); } gobble(in);
472 for (int k = 0; k < tempNum; k++) { in >> temp; ambigBases.push_back(temp); } gobble(in);
473 for (int k = 0; k < tempNum; k++) { in >> temp; longHomoPolymer.push_back(temp); } gobble(in);
476 remove(tempFilename.c_str());
482 catch(exception& e) {
483 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
487 //***************************************************************************************************************