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;
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;
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 //send file positions to all processes
149 MPI_Bcast(&numSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD); //send numSeqs
150 MPI_Bcast(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //send file pos
152 //figure out how many sequences you have to do
153 numSeqsPerProcessor = numSeqs / processors;
154 int startIndex = pid * numSeqsPerProcessor;
155 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
158 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
160 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
162 //get the info from the child processes
163 for(int i = 1; i < processors; i++) {
166 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
168 vector<int> temp; temp.resize(size+1);
170 for(int j = 0; j < 5; j++) {
172 MPI_Recv(&temp[0], (size+1), MPI_INT, i, 2001, MPI_COMM_WORLD, &status);
173 int receiveTag = temp[temp.size()-1]; //child process added a int to the end to indicate what count this is for
175 if (receiveTag == startTag) {
176 for (int k = 0; k < size; k++) { startPosition.push_back(temp[k]); }
177 }else if (receiveTag == endTag) {
178 for (int k = 0; k < size; k++) { endPosition.push_back(temp[k]); }
179 }else if (receiveTag == lengthTag) {
180 for (int k = 0; k < size; k++) { seqLength.push_back(temp[k]); }
181 }else if (receiveTag == baseTag) {
182 for (int k = 0; k < size; k++) { ambigBases.push_back(temp[k]); }
183 }else if (receiveTag == lhomoTag) {
184 for (int k = 0; k < size; k++) { longHomoPolymer.push_back(temp[k]); }
191 }else { //i am the child process
193 MPI_Bcast(&numSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD); //get numSeqs
194 MPIPos.resize(numSeqs+1);
195 MPI_Bcast(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //get file positions
197 //figure out how many sequences you have to align
198 numSeqsPerProcessor = numSeqs / processors;
199 int startIndex = pid * numSeqsPerProcessor;
200 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
203 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
205 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); return 0; }
208 int size = startPosition.size();
209 MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
211 startPosition.push_back(startTag);
212 int ierr = MPI_Send(&(startPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
213 endPosition.push_back(endTag);
214 ierr = MPI_Send (&(endPosition[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
215 seqLength.push_back(lengthTag);
216 ierr = MPI_Send(&(seqLength[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
217 ambigBases.push_back(baseTag);
218 ierr = MPI_Send(&(ambigBases[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
219 longHomoPolymer.push_back(lhomoTag);
220 ierr = MPI_Send(&(longHomoPolymer[0]), (size+1), MPI_INT, 0, 2001, MPI_COMM_WORLD);
224 MPI_File_close(&inMPI);
225 MPI_File_close(&outMPI);
228 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
231 openInputFile(fastafile, inFASTA);
232 numSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
235 lines.push_back(new linePair(0, numSeqs));
237 driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
239 numSeqs = setLines(fastafile);
240 createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile);
242 rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
244 for(int i=1;i<processors;i++){
245 appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
246 remove((summaryFile + toString(processIDS[i]) + ".temp").c_str());
250 if (m->control_pressed) { return 0; }
253 openInputFile(fastafile, inFASTA);
254 numSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
257 lines.push_back(new linePair(0, numSeqs));
259 driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
260 if (m->control_pressed) { return 0; }
268 sort(startPosition.begin(), startPosition.end());
269 sort(endPosition.begin(), endPosition.end());
270 sort(seqLength.begin(), seqLength.end());
271 sort(ambigBases.begin(), ambigBases.end());
272 sort(longHomoPolymer.begin(), longHomoPolymer.end());
274 int ptile0_25 = int(numSeqs * 0.025);
275 int ptile25 = int(numSeqs * 0.250);
276 int ptile50 = int(numSeqs * 0.500);
277 int ptile75 = int(numSeqs * 0.750);
278 int ptile97_5 = int(numSeqs * 0.975);
279 int ptile100 = numSeqs - 1;
281 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
282 if (startPosition[0] == -1) { startPosition[0] = 0; }
283 if (endPosition[0] == -1) { endPosition[0] = 0; }
285 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
287 m->mothurOutEndLine();
288 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
289 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();
290 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();
291 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();
292 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();
293 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();
294 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();
295 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();
296 m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine();
298 if (m->control_pressed) { remove(summaryFile.c_str()); return 0; }
300 m->mothurOutEndLine();
301 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
302 m->mothurOut(summaryFile); m->mothurOutEndLine();
303 m->mothurOutEndLine();
311 catch(exception& e) {
312 m->errorOut(e, "SeqSummaryCommand", "execute");
316 /**************************************************************************************/
317 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* line) {
321 openOutputFile(sumFile, outSummary);
323 //print header if you are process 0
324 if (line->start == 0) {
325 outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer" << endl;
329 openInputFile(filename, in);
331 in.seekg(line->start);
333 for(int i=0;i<line->num;i++){
335 if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
337 Sequence current(in);
338 if (current.getName() != "") {
339 startPosition.push_back(current.getStartPos());
340 endPosition.push_back(current.getEndPos());
341 seqLength.push_back(current.getNumBases());
342 ambigBases.push_back(current.getAmbigBases());
343 longHomoPolymer.push_back(current.getLongHomoPolymer());
345 outSummary << current.getName() << '\t';
346 outSummary << current.getStartPos() << '\t' << current.getEndPos() << '\t';
347 outSummary << current.getNumBases() << '\t' << current.getAmbigBases() << '\t';
348 outSummary << current.getLongHomoPolymer() << endl;
356 catch(exception& e) {
357 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
362 /**************************************************************************************/
363 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<long>& MPIPos) {
368 for(int i=0;i<num;i++){
370 if (m->control_pressed) { return 0; }
373 int length = MPIPos[start+i+1] - MPIPos[start+i];
375 char* buf4 = new char[length];
376 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
378 string tempBuf = buf4;
379 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
380 istringstream iss (tempBuf,istringstream::in);
383 Sequence current(iss);
385 if (current.getName() != "") {
386 startPosition.push_back(current.getStartPos());
387 endPosition.push_back(current.getEndPos());
388 seqLength.push_back(current.getNumBases());
389 ambigBases.push_back(current.getAmbigBases());
390 longHomoPolymer.push_back(current.getLongHomoPolymer());
392 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
393 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\n";
396 length = outputString.length();
397 char* buf3 = new char[length];
398 memcpy(buf3, outputString.c_str(), length);
400 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
407 catch(exception& e) {
408 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
413 /**************************************************************************************************/
414 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
416 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
421 //loop through and create all the processes you want
422 while (process != processors) {
426 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
429 driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
431 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
434 //force parent to wait until all the processes are done
435 for (int i=0;i<processors;i++) {
436 int temp = processIDS[i];
443 catch(exception& e) {
444 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
448 /**************************************************************************************************/
450 int SeqSummaryCommand::setLines(string filename) {
453 vector<long int> positions;
456 openInputFile(filename, inFASTA);
459 while(!inFASTA.eof()){
460 input = getline(inFASTA);
462 if (input.length() != 0) {
463 if(input[0] == '>'){ long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
468 int numFastaSeqs = positions.size();
473 //get num bytes in file
474 pFile = fopen (filename.c_str(),"rb");
475 if (pFile==NULL) perror ("Error opening file");
477 fseek (pFile, 0, SEEK_END);
482 int numSeqsPerProcessor = numFastaSeqs / processors;
484 for (int i = 0; i < processors; i++) {
486 long int startPos = positions[ i * numSeqsPerProcessor ];
487 if(i == processors - 1){
488 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
490 long int myEnd = positions[ (i+1) * numSeqsPerProcessor ];
492 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
497 catch(exception& e) {
498 m->errorOut(e, "SeqSummaryCommand", "setLines");
502 //***************************************************************************************************************