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","processors","outputdir","inputdir"};
17 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
21 m->errorOut(e, "SeqSummaryCommand", "getValidParameters");
25 //**********************************************************************************************************************
26 SeqSummaryCommand::SeqSummaryCommand(){
29 //initialize outputTypes
30 vector<string> tempOutNames;
31 outputTypes["summary"] = tempOutNames;
34 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
38 //**********************************************************************************************************************
39 vector<string> SeqSummaryCommand::getRequiredParameters(){
41 string Array[] = {"fasta"};
42 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
46 m->errorOut(e, "SeqSummaryCommand", "getRequiredParameters");
50 //**********************************************************************************************************************
51 vector<string> SeqSummaryCommand::getRequiredFiles(){
53 vector<string> myArray;
57 m->errorOut(e, "SeqSummaryCommand", "getRequiredFiles");
61 //***************************************************************************************************************
63 SeqSummaryCommand::SeqSummaryCommand(string option) {
67 //allow user to run help
68 if(option == "help") { help(); abort = true; }
71 //valid paramters for this command
72 string Array[] = {"fasta","processors","outputdir","inputdir"};
73 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
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; }
100 //initialize outputTypes
101 vector<string> tempOutNames;
102 outputTypes["summary"] = tempOutNames;
104 //check for required parameters
105 fastafile = validParameter.validFile(parameters, "fasta", true);
106 if (fastafile == "not open") { abort = true; }
107 else if (fastafile == "not found") { fastafile = ""; m->mothurOut("fasta is a required parameter for the summary.seqs command."); m->mothurOutEndLine(); abort = true; }
109 //if the user changes the output directory command factory will send this info to us in the output parameter
110 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
112 outputDir += m->hasPath(fastafile); //if user entered a file with a path then preserve it
115 string temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
116 convert(temp, processors);
121 catch(exception& e) {
122 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
126 //**********************************************************************************************************************
128 void SeqSummaryCommand::help(){
130 m->mothurOut("The summary.seqs command reads a fastafile and summarizes the sequences.\n");
131 m->mothurOut("The summary.seqs command parameters are fasta and processors, fasta is required.\n");
132 m->mothurOut("The summary.seqs command should be in the following format: \n");
133 m->mothurOut("summary.seqs(fasta=yourFastaFile, processors=2) \n");
134 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");
136 catch(exception& e) {
137 m->errorOut(e, "SeqSummaryCommand", "help");
142 //***************************************************************************************************************
144 SeqSummaryCommand::~SeqSummaryCommand(){ /* do nothing */ }
146 //***************************************************************************************************************
148 int SeqSummaryCommand::execute(){
151 if (abort == true) { return 0; }
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;
164 int pid, numSeqsPerProcessor;
166 int startTag = 1; int endTag = 2; int lengthTag = 3; int baseTag = 4; int lhomoTag = 5;
167 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
168 vector<unsigned long int> MPIPos;
171 MPI_Status statusOut;
174 MPI_Comm_size(MPI_COMM_WORLD, &processors);
175 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
177 char tempFileName[1024];
178 strcpy(tempFileName, fastafile.c_str());
180 char sumFileName[1024];
181 strcpy(sumFileName, summaryFile.c_str());
183 MPI_File_open(MPI_COMM_WORLD, tempFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
184 MPI_File_open(MPI_COMM_WORLD, sumFileName, outMode, MPI_INFO_NULL, &outMPI);
186 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
188 if (pid == 0) { //you are the root process
190 string outputString = "seqname\tstart\tend\tnbases\tambigs\tpolymer\n";
191 int length = outputString.length();
192 char* buf2 = new char[length];
193 memcpy(buf2, outputString.c_str(), length);
195 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
198 MPIPos = m->setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
200 for(int i = 1; i < processors; i++) {
201 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
202 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
205 //figure out how many sequences you have to do
206 numSeqsPerProcessor = numSeqs / processors;
207 int startIndex = pid * numSeqsPerProcessor;
208 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
211 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
213 }else { //i am the child process
215 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
216 MPIPos.resize(numSeqs+1);
217 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
219 //figure out how many sequences you have to align
220 numSeqsPerProcessor = numSeqs / processors;
221 int startIndex = pid * numSeqsPerProcessor;
222 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
225 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
228 MPI_File_close(&inMPI);
229 MPI_File_close(&outMPI);
230 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
233 //get the info from the child processes
234 for(int i = 1; i < processors; i++) {
236 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
238 vector<int> temp; temp.resize(size+1);
240 for(int j = 0; j < 5; j++) {
242 MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
243 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what count this is for
245 if (receiveTag == startTag) {
246 for (int k = 0; k < size; k++) { startPosition.push_back(temp[k]); }
247 }else if (receiveTag == endTag) {
248 for (int k = 0; k < size; k++) { endPosition.push_back(temp[k]); }
249 }else if (receiveTag == lengthTag) {
250 for (int k = 0; k < size; k++) { seqLength.push_back(temp[k]); }
251 }else if (receiveTag == baseTag) {
252 for (int k = 0; k < size; k++) { ambigBases.push_back(temp[k]); }
253 }else if (receiveTag == lhomoTag) {
254 for (int k = 0; k < size; k++) { longHomoPolymer.push_back(temp[k]); }
262 int size = startPosition.size();
263 MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
265 startPosition.push_back(startTag);
266 int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
267 endPosition.push_back(endTag);
268 ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
269 seqLength.push_back(lengthTag);
270 ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
271 ambigBases.push_back(baseTag);
272 ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
273 longHomoPolymer.push_back(lhomoTag);
274 ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
277 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
279 vector<unsigned long int> positions = m->divideFile(fastafile, processors);
281 for (int i = 0; i < (positions.size()-1); i++) {
282 lines.push_back(new linePair(positions[i], positions[(i+1)]));
285 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
287 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
289 numSeqs = createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile);
291 rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
293 for(int i=1;i<processors;i++){
294 m->appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
295 remove((summaryFile + toString(processIDS[i]) + ".temp").c_str());
299 if (m->control_pressed) { return 0; }
301 numSeqs = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
302 if (m->control_pressed) { return 0; }
310 sort(startPosition.begin(), startPosition.end());
311 sort(endPosition.begin(), endPosition.end());
312 sort(seqLength.begin(), seqLength.end());
313 sort(ambigBases.begin(), ambigBases.end());
314 sort(longHomoPolymer.begin(), longHomoPolymer.end());
316 int ptile0_25 = int(numSeqs * 0.025);
317 int ptile25 = int(numSeqs * 0.250);
318 int ptile50 = int(numSeqs * 0.500);
319 int ptile75 = int(numSeqs * 0.750);
320 int ptile97_5 = int(numSeqs * 0.975);
321 int ptile100 = numSeqs - 1;
323 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
324 if (startPosition[0] == -1) { startPosition[0] = 0; }
325 if (endPosition[0] == -1) { endPosition[0] = 0; }
327 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
329 m->mothurOutEndLine();
330 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
331 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();
332 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();
333 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();
334 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();
335 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();
336 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();
337 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();
338 m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine();
340 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
342 m->mothurOutEndLine();
343 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
344 m->mothurOut(summaryFile); m->mothurOutEndLine(); outputNames.push_back(summaryFile); outputTypes["summary"].push_back(summaryFile);
345 m->mothurOutEndLine();
353 catch(exception& e) {
354 m->errorOut(e, "SeqSummaryCommand", "execute");
358 /**************************************************************************************/
359 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* filePos) {
363 m->openOutputFile(sumFile, outSummary);
365 //print header if you are process 0
366 if (filePos->start == 0) {
367 outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer" << endl;
371 m->openInputFile(filename, in);
373 in.seekg(filePos->start);
380 if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
382 Sequence current(in); m->gobble(in);
384 if (current.getName() != "") {
385 startPosition.push_back(current.getStartPos());
386 endPosition.push_back(current.getEndPos());
387 seqLength.push_back(current.getNumBases());
388 ambigBases.push_back(current.getAmbigBases());
389 longHomoPolymer.push_back(current.getLongHomoPolymer());
391 outSummary << current.getName() << '\t';
392 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
393 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
394 outSummary << current.getLongHomoPolymer() << endl;
398 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
399 unsigned long int pos = in.tellg();
400 if ((pos == -1) || (pos >= filePos->end)) { break; }
402 if (in.eof()) { break; }
406 //if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
409 //if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
415 catch(exception& e) {
416 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
421 /**************************************************************************************/
422 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) {
427 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
430 for(int i=0;i<num;i++){
432 if (m->control_pressed) { return 0; }
435 int length = MPIPos[start+i+1] - MPIPos[start+i];
437 char* buf4 = new char[length];
438 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
440 string tempBuf = buf4;
441 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
442 istringstream iss (tempBuf,istringstream::in);
445 Sequence current(iss);
447 if (current.getName() != "") {
448 startPosition.push_back(current.getStartPos());
449 endPosition.push_back(current.getEndPos());
450 seqLength.push_back(current.getNumBases());
451 ambigBases.push_back(current.getAmbigBases());
452 longHomoPolymer.push_back(current.getLongHomoPolymer());
454 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
455 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\n";
458 length = outputString.length();
459 char* buf3 = new char[length];
460 memcpy(buf3, outputString.c_str(), length);
462 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
469 catch(exception& e) {
470 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
475 /**************************************************************************************************/
476 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
478 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
483 //loop through and create all the processes you want
484 while (process != processors) {
488 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
491 num = driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
493 //pass numSeqs to parent
495 string tempFile = fastafile + toString(getpid()) + ".num.temp";
496 m->openOutputFile(tempFile, out);
499 for (int k = 0; k < startPosition.size(); k++) { out << startPosition[k] << '\t'; } out << endl;
500 for (int k = 0; k < endPosition.size(); k++) { out << endPosition[k] << '\t'; } out << endl;
501 for (int k = 0; k < seqLength.size(); k++) { out << seqLength[k] << '\t'; } out << endl;
502 for (int k = 0; k < ambigBases.size(); k++) { out << ambigBases[k] << '\t'; } out << endl;
503 for (int k = 0; k < longHomoPolymer.size(); k++) { out << longHomoPolymer[k] << '\t'; } out << endl;
509 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
510 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
515 //force parent to wait until all the processes are done
516 for (int i=0;i<processors;i++) {
517 int temp = processIDS[i];
521 //parent reads in and combine Filter info
522 for (int i = 0; i < processIDS.size(); i++) {
523 string tempFilename = fastafile + toString(processIDS[i]) + ".num.temp";
525 m->openInputFile(tempFilename, in);
528 in >> tempNum; m->gobble(in); num += tempNum;
529 for (int k = 0; k < tempNum; k++) { in >> temp; startPosition.push_back(temp); } m->gobble(in);
530 for (int k = 0; k < tempNum; k++) { in >> temp; endPosition.push_back(temp); } m->gobble(in);
531 for (int k = 0; k < tempNum; k++) { in >> temp; seqLength.push_back(temp); } m->gobble(in);
532 for (int k = 0; k < tempNum; k++) { in >> temp; ambigBases.push_back(temp); } m->gobble(in);
533 for (int k = 0; k < tempNum; k++) { in >> temp; longHomoPolymer.push_back(temp); } m->gobble(in);
536 remove(tempFilename.c_str());
542 catch(exception& e) {
543 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
547 //***************************************************************************************************************