]> git.donarmstrong.com Git - mothur.git/blob - seqsummarycommand.cpp
some bug fixes
[mothur.git] / seqsummarycommand.cpp
1 /*
2  *  seqcoordcommand.cpp
3  *  Mothur
4  *
5  *  Created by Pat Schloss on 5/30/09.
6  *  Copyright 2009 Patrick D. Schloss. All rights reserved.
7  *
8  */
9
10 #include "seqsummarycommand.h"
11 #include "sequence.hpp"
12
13 //***************************************************************************************************************
14
15 SeqSummaryCommand::SeqSummaryCommand(string option)  {
16         try {
17                 abort = false;
18                 
19                 //allow user to run help
20                 if(option == "help") { help(); abort = true; }
21                 
22                 else {
23                         //valid paramters for this command
24                         string Array[] =  {"fasta","processors","outputdir","inputdir"};
25                         vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
26                         
27                         OptionParser parser(option);
28                         map<string,string> parameters = parser.getParameters();
29                         
30                         ValidParameters validParameter("summary.seqs");
31                         map<string,string>::iterator it;
32                         
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;  }
36                         }
37                         
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 = "";          }
41                         else {
42                                 string path;
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;            }
49                                 }
50                         }
51                         
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;  }       
56                         
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"){  
59                                 outputDir = ""; 
60                                 outputDir += hasPath(fastafile); //if user entered a file with a path then preserve it  
61                         }
62                         
63                         string temp = validParameter.validFile(parameters, "processors", false);        if (temp == "not found"){       temp = "1";                             }
64                         convert(temp, processors); 
65
66
67                 }
68         }
69         catch(exception& e) {
70                 m->errorOut(e, "SeqSummaryCommand", "SeqSummaryCommand");
71                 exit(1);
72         }
73 }
74 //**********************************************************************************************************************
75
76 void SeqSummaryCommand::help(){
77         try {
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");     
83         }
84         catch(exception& e) {
85                 m->errorOut(e, "SeqSummaryCommand", "help");
86                 exit(1);
87         }
88 }
89
90 //***************************************************************************************************************
91
92 SeqSummaryCommand::~SeqSummaryCommand(){        /*      do nothing      */      }
93
94 //***************************************************************************************************************
95
96 int SeqSummaryCommand::execute(){
97         try{
98                 
99                 if (abort == true) { return 0; }
100                 
101                 string summaryFile = outputDir + getSimpleName(fastafile) + ".summary";
102                                 
103                 int numSeqs = 0;
104                 
105                 vector<int> startPosition;
106                 vector<int> endPosition;
107                 vector<int> seqLength;
108                 vector<int> ambigBases;
109                 vector<int> longHomoPolymer;
110                         
111 #ifdef USE_MPI  
112                                 int pid, numSeqsPerProcessor; 
113                                 int tag = 2001;
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;
117                                 
118                                 MPI_Status status; 
119                                 MPI_Status statusOut;
120                                 MPI_File inMPI; 
121                                 MPI_File outMPI; 
122                                 MPI_Comm_size(MPI_COMM_WORLD, &processors);
123                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
124                                                         
125                                 char tempFileName[1024];
126                                 strcpy(tempFileName, fastafile.c_str());
127                                 
128                                 char sumFileName[1024];
129                                 strcpy(sumFileName, summaryFile.c_str());
130                 
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);
133                                 
134                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI); return 0;  }
135                                 
136                                 if (pid == 0) { //you are the root process
137                                                 //print header
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);
142                                         
143                                                 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &statusOut);
144                                                 delete buf2;
145                                                 
146                                                 MPIPos = setFilePosFasta(fastafile, numSeqs); //fills MPIPos, returns numSeqs
147                                         
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);
151                                                 }
152                                                 
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;      }
157                                                 
158                                                 //do your part
159                                                 MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
160                                                 
161                                 }else { //i am the child process
162                         
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);
166                                 
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;      }
171                                 
172                                         //do your part
173                                         MPICreateSummary(startIndex, numSeqsPerProcessor, startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, inMPI, outMPI, MPIPos);
174                                 }
175                                 
176                                 MPI_File_close(&inMPI);
177                                 MPI_File_close(&outMPI);
178                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
179                                 
180                                 if (pid == 0) {
181                                         //get the info from the child processes
182                                         for(int i = 1; i < processors; i++) { 
183                                                 int size;
184                                                 MPI_Recv(&size, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
185
186                                                 vector<int> temp; temp.resize(size+1);
187                                                 
188                                                 for(int j = 0; j < 5; j++) { 
189                                                 
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
192                                                         
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]);     }
203                                                         }
204                                                 } 
205                                         }
206
207                                 }else{
208                                 
209                                         //send my counts
210                                         int size = startPosition.size();
211                                         MPI_Send(&size, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
212                                         
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);
223                                 }
224                                 
225                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
226 #else
227                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
228                                 if(processors == 1){
229                                         ifstream inFASTA;
230                                         openInputFile(fastafile, inFASTA);
231                                         getNumSeqs(inFASTA, numSeqs);
232                                         inFASTA.close();        
233                                                 
234                                         lines.push_back(new linePair(0, numSeqs));
235                                         
236                                         driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
237                                 }else{
238                                         numSeqs = setLines(fastafile);                                  
239                                         createProcessesCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile); 
240                                         
241                                         rename((summaryFile + toString(processIDS[0]) + ".temp").c_str(), summaryFile.c_str());
242                                         //append files
243                                         for(int i=1;i<processors;i++){
244                                                 appendFiles((summaryFile + toString(processIDS[i]) + ".temp"), summaryFile);
245                                                 remove((summaryFile + toString(processIDS[i]) + ".temp").c_str());
246                                         }
247                                 }
248                                 
249                                 if (m->control_pressed) {  return 0; }
250                 #else
251                                 ifstream inFASTA;
252                                 openInputFile(fastafile, inFASTA);
253                                 getNumSeqs(inFASTA, numSeqs);
254                                 inFASTA.close();
255                                 
256                                 lines.push_back(new linePair(0, numSeqs));
257                                 
258                                 driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, summaryFile, lines[0]);
259                                 if (m->control_pressed) {  return 0; }
260                 #endif
261 #endif
262                         
263                 #ifdef USE_MPI
264                         if (pid == 0) { 
265                 #endif
266                 
267                 sort(startPosition.begin(), startPosition.end());
268                 sort(endPosition.begin(), endPosition.end());
269                 sort(seqLength.begin(), seqLength.end());
270                 sort(ambigBases.begin(), ambigBases.end());
271                 sort(longHomoPolymer.begin(), longHomoPolymer.end());
272                 
273                 int ptile0_25   = int(numSeqs * 0.025);
274                 int ptile25             = int(numSeqs * 0.250);
275                 int ptile50             = int(numSeqs * 0.500);
276                 int ptile75             = int(numSeqs * 0.750);
277                 int ptile97_5   = int(numSeqs * 0.975);
278                 int ptile100    = numSeqs - 1;
279                 
280                 //to compensate for blank sequences that would result in startPosition and endPostion equalling -1
281                 if (startPosition[0] == -1) {  startPosition[0] = 0;    }
282                 if (endPosition[0] == -1)       {  endPosition[0] = 0;          }
283                 
284                 if (m->control_pressed) {  remove(summaryFile.c_str()); return 0; }
285                 
286                 m->mothurOutEndLine();
287                 m->mothurOut("\t\tStart\tEnd\tNBases\tAmbigs\tPolymer"); m->mothurOutEndLine();
288                 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();
289                 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();
290                 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();
291                 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();
292                 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();
293                 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();
294                 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();
295                 m->mothurOut("# of Seqs:\t" + toString(numSeqs)); m->mothurOutEndLine();
296                 
297                 if (m->control_pressed) {  remove(summaryFile.c_str()); return 0; }
298                 
299                 m->mothurOutEndLine();
300                 m->mothurOut("Output File Name: "); m->mothurOutEndLine();
301                 m->mothurOut(summaryFile); m->mothurOutEndLine();       
302                 m->mothurOutEndLine();
303                 
304                 #ifdef USE_MPI
305                         }
306                 #endif
307
308                 return 0;
309         }
310         catch(exception& e) {
311                 m->errorOut(e, "SeqSummaryCommand", "execute");
312                 exit(1);
313         }
314 }
315 /**************************************************************************************/
316 int SeqSummaryCommand::driverCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile, linePair* line) {      
317         try {
318                 
319                 ofstream outSummary;
320                 openOutputFile(sumFile, outSummary);
321                 
322                 //print header if you are process 0
323                 if (line->start == 0) {
324                         outSummary << "seqname\tstart\tend\tnbases\tambigs\tpolymer" << endl;   
325                 }
326                                 
327                 ifstream in;
328                 openInputFile(filename, in);
329                                 
330                 in.seekg(line->start);
331                 
332                 for(int i=0;i<line->num;i++){
333                                 
334                         if (m->control_pressed) { in.close(); outSummary.close(); return 1; }
335                                         
336                         Sequence current(in);
337         
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());
344                                 
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;
349                         }
350                         gobble(in);
351                 }
352                 in.close();
353                 
354                 return 0;
355         }
356         catch(exception& e) {
357                 m->errorOut(e, "SeqSummaryCommand", "driverCreateSummary");
358                 exit(1);
359         }
360 }
361 #ifdef USE_MPI
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<unsigned long int>& MPIPos) {        
364         try {
365                 
366                 int pid;
367                 MPI_Status status; 
368                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
369
370                 
371                 for(int i=0;i<num;i++){
372                         
373                         if (m->control_pressed) { return 0; }
374                         
375                         //read next sequence
376                         int length = MPIPos[start+i+1] - MPIPos[start+i];
377         
378                         char* buf4 = new char[length];
379                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
380                         
381                         string tempBuf = buf4;
382                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
383                         istringstream iss (tempBuf,istringstream::in);
384                         delete buf4;
385
386                         Sequence current(iss);  
387
388                         if (current.getName() != "") {
389                                 startPosition.push_back(current.getStartPos());
390                                 endPosition.push_back(current.getEndPos());
391                                 seqLength.push_back(current.getNumBases());
392                                 ambigBases.push_back(current.getAmbigBases());
393                                 longHomoPolymer.push_back(current.getLongHomoPolymer());
394                                 
395                                 string outputString = current.getName() + "\t" + toString(current.getStartPos()) + "\t" + toString(current.getEndPos()) + "\t";
396                                 outputString += toString(current.getNumBases()) + "\t" + toString(current.getAmbigBases()) + "\t" + toString(current.getLongHomoPolymer()) + "\n";
397                                 
398                                 //output to file
399                                 length = outputString.length();
400                                 char* buf3 = new char[length];
401                                 memcpy(buf3, outputString.c_str(), length);
402                                         
403                                 MPI_File_write_shared(outMPI, buf3, length, MPI_CHAR, &status);
404                                 delete buf3;
405                         }       
406                 }
407                 
408                 return 0;
409         }
410         catch(exception& e) {
411                 m->errorOut(e, "SeqSummaryCommand", "MPICreateSummary");
412                 exit(1);
413         }
414 }
415 #endif
416 /**************************************************************************************************/
417 int SeqSummaryCommand::createProcessesCreateSummary(vector<int>& startPosition, vector<int>& endPosition, vector<int>& seqLength, vector<int>& ambigBases, vector<int>& longHomoPolymer, string filename, string sumFile) {
418         try {
419 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
420                 int process = 0;
421                 int exitCommand = 1;
422                 processIDS.clear();
423                 
424                 //loop through and create all the processes you want
425                 while (process != processors) {
426                         int pid = vfork();
427                         
428                         if (pid > 0) {
429                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
430                                 process++;
431                         }else if (pid == 0){
432                                 driverCreateSummary(startPosition, endPosition, seqLength, ambigBases, longHomoPolymer, fastafile, sumFile + toString(getpid()) + ".temp", lines[process]);
433                                 exit(0);
434                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
435                 }
436                 
437                 //force parent to wait until all the processes are done
438                 for (int i=0;i<processors;i++) { 
439                         int temp = processIDS[i];
440                         wait(&temp);
441                 }
442                 
443                 return exitCommand;
444 #endif          
445         }
446         catch(exception& e) {
447                 m->errorOut(e, "SeqSummaryCommand", "createProcessesCreateSummary");
448                 exit(1);
449         }
450 }
451 /**************************************************************************************************/
452
453 int SeqSummaryCommand::setLines(string filename) {
454         try {
455                 
456                 vector<unsigned long int> positions;
457                 
458                 ifstream inFASTA;
459                 openInputFile(filename, inFASTA);
460                         
461                 string input;
462                 while(!inFASTA.eof()){  
463                         input = getline(inFASTA);
464
465                         if (input.length() != 0) {
466                                 if(input[0] == '>'){ unsigned long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);    }
467                         }
468                 }
469                 inFASTA.close();
470                 
471                 int numFastaSeqs = positions.size();
472         
473                 FILE * pFile;
474                 unsigned long int size;
475                 
476                 //get num bytes in file
477                 pFile = fopen (filename.c_str(),"rb");
478                 if (pFile==NULL) perror ("Error opening file");
479                 else{
480                         fseek (pFile, 0, SEEK_END);
481                         size=ftell (pFile);
482                         fclose (pFile);
483                 }
484                 
485                 int numSeqsPerProcessor = numFastaSeqs / processors;
486                 
487                 for (int i = 0; i < processors; i++) {
488
489                         unsigned long int startPos = positions[ i * numSeqsPerProcessor ];
490                         if(i == processors - 1){
491                                 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
492                         }else{  
493                                 unsigned long int myEnd = positions[ (i+1) * numSeqsPerProcessor ];
494                         }
495                         lines.push_back(new linePair(startPos, numSeqsPerProcessor));
496                 }
497                 
498                 return numFastaSeqs;
499         }
500         catch(exception& e) {
501                 m->errorOut(e, "SeqSummaryCommand", "setLines");
502                 exit(1);
503         }
504 }
505 //***************************************************************************************************************
506
507