]> git.donarmstrong.com Git - mothur.git/blob - aligncommand.cpp
worked on hcluster. made .single command run using a sharedfile. and various other...
[mothur.git] / aligncommand.cpp
1 /*
2  *  aligncommand.cpp
3  *  Mothur
4  *
5  *  Created by Sarah Westcott on 5/15/09.
6  *  Copyright 2009 Schloss Lab UMASS Amherst. All rights reserved.
7  *
8  *      This version of nast does everything I think that the greengenes nast server does and then some.  I have added the 
9  *      feature of allowing users to define their database, kmer size for searching, alignment penalty values and alignment 
10  *      method.  This latter feature is perhaps most significant.  nastPlus enables a user to use either a Needleman-Wunsch 
11  *      (non-affine gap penalty) or Gotoh (affine gap penalty) pairwise alignment algorithm.  This is significant because it
12  *      allows for a global alignment and not the local alignment provided by bLAst.  Furthermore, it has the potential to
13  *      provide a better alignment because of the banding method employed by blast (I'm not sure about this).
14  *
15  */
16
17 #include "aligncommand.h"
18 #include "sequence.hpp"
19
20 #include "gotohoverlap.hpp"
21 #include "needlemanoverlap.hpp"
22 #include "blastalign.hpp"
23 #include "noalign.hpp"
24
25 #include "nast.hpp"
26 #include "nastreport.hpp"
27
28
29 //**********************************************************************************************************************
30
31 AlignCommand::AlignCommand(string option){
32         try {
33                 //              globaldata = GlobalData::getInstance();
34                 abort = false;
35                 
36                 //allow user to run help
37                 if(option == "help") { help(); abort = true; }
38                 
39                 else {
40                         
41                         //valid paramters for this command
42                         string AlignArray[] =  {"template","candidate","search","ksize","align","match","mismatch","gapopen","gapextend", "processors","flip","threshold"};
43                         vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
44                         
45                         OptionParser parser(option);
46                         map<string, string> parameters = parser.getParameters(); 
47                         
48                         ValidParameters validParameter;
49                         
50                         //check to make sure all parameters are valid for command
51                         for (map<string, string>::iterator it = parameters.begin(); it != parameters.end(); it++) { 
52                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
53                         }
54                         
55                         //check for required parameters
56                         templateFileName = validParameter.validFile(parameters, "template", true);
57                         if (templateFileName == "not found") { 
58                                 mothurOut("template is a required parameter for the align.seqs command."); 
59                                 mothurOutEndLine();
60                                 abort = true; 
61                         }
62                         else if (templateFileName == "not open") { abort = true; }      
63                         
64                         candidateFileName = validParameter.validFile(parameters, "candidate", false);
65                         if (candidateFileName == "not found") { mothurOut("candidate is a required parameter for the align.seqs command."); mothurOutEndLine(); abort = true;  }
66                         else { 
67                                 splitAtDash(candidateFileName, candidateFileNames);
68                                 
69                                 //go through files and make sure they are good, if not, then disregard them
70                                 for (int i = 0; i < candidateFileNames.size(); i++) {
71                                         int ableToOpen;
72                                         ifstream in;
73                                         ableToOpen = openInputFile(candidateFileNames[i], in);
74                                         if (ableToOpen == 1) { 
75                                                 mothurOut(candidateFileNames[i] + " will be disregarded."); mothurOutEndLine(); 
76                                                 //erase from file list
77                                                 candidateFileNames.erase(candidateFileNames.begin()+i);
78                                                 i--;
79                                         }
80                                         in.close();
81                                 }
82                                 
83                                 //make sure there is at least one valid file left
84                                 if (candidateFileNames.size() == 0) { mothurOut("no valid files."); mothurOutEndLine(); abort = true; }
85                         }
86                                 
87                         
88                         //check for optional parameter and set defaults
89                         // ...at some point should added some additional type checking...
90                         string temp;
91                         temp = validParameter.validFile(parameters, "ksize", false);            if (temp == "not found"){       temp = "8";                             }
92                         convert(temp, kmerSize); 
93                         
94                         temp = validParameter.validFile(parameters, "match", false);            if (temp == "not found"){       temp = "1.0";                   }
95                         convert(temp, match);  
96                         
97                         temp = validParameter.validFile(parameters, "mismatch", false);         if (temp == "not found"){       temp = "-1.0";                  }
98                         convert(temp, misMatch);  
99                         
100                         temp = validParameter.validFile(parameters, "gapopen", false);          if (temp == "not found"){       temp = "-2.0";                  }
101                         convert(temp, gapOpen);  
102                         
103                         temp = validParameter.validFile(parameters, "gapextend", false);        if (temp == "not found"){       temp = "-1.0";                  }
104                         convert(temp, gapExtend); 
105                         
106                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = "1";                             }
107                         convert(temp, processors); 
108                         
109                         temp = validParameter.validFile(parameters, "flip", false);                     if (temp == "not found"){       temp = "f";                             }
110                         flip = isTrue(temp); 
111                         
112                         temp = validParameter.validFile(parameters, "threshold", false);        if (temp == "not found"){       temp = "0.50";                  }
113                         convert(temp, threshold); 
114                         
115                         search = validParameter.validFile(parameters, "search", false);         if (search == "not found"){     search = "kmer";                }
116                         
117                         align = validParameter.validFile(parameters, "align", false);           if (align == "not found"){      align = "needleman";    }
118                 }
119                 
120         }
121         catch(exception& e) {
122                 errorOut(e, "AlignCommand", "AlignCommand");
123                 exit(1);
124         }
125 }
126
127 //**********************************************************************************************************************
128
129 AlignCommand::~AlignCommand(){  
130
131         if (abort == false) {
132                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
133                 delete templateDB;
134                 delete alignment;
135         }
136 }
137
138 //**********************************************************************************************************************
139
140 void AlignCommand::help(){
141         try {
142                 mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
143                 mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen and gapextend.\n");
144                 mothurOut("The template and candidate parameters are required. You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amzon.fasta \n");
145                 mothurOut("The search parameter allows you to specify the method to find most similar template.  Your options are: suffix, kmer and blast. The default is kmer.\n");
146                 mothurOut("The align parameter allows you to specify the alignment method to use.  Your options are: gotoh, needleman, blast and noalign. The default is needleman.\n");
147                 mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate.  The default is 8.\n");
148                 mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
149                 mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases.  The default is -1.0.\n");
150                 mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
151                 mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment.  The default is -1.0.\n");
152                 mothurOut("The flip parameter is used to specify whether or not you want mothur to try the reverse compement if a sequence falls below the threshold.  The default is false.\n");
153                 mothurOut("The threshold is used to specify a cutoff at which an alignment is deemed 'bad' and the reverse complement may be tried. \n");
154                 mothurOut("If the flip parameter is set to true the reverse complement of the sequence is aligned and the better alignment is reported.\n");
155                 mothurOut("The default for the threshold parameter is 0.50, meaning at least 50% of the bases must remain or the sequence is reported as potentially reversed.\n");
156                 mothurOut("The align.seqs command should be in the following format: \n");
157                 mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
158                 mothurOut("Example align.seqs(candidate=candidate.fasta, template=core.filtered, align=kmer, search=gotoh, ksize=8, match=2.0, mismatch=3.0, gapopen=-2.0, gapextend=-1.0)\n");
159                 mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
160         }
161         catch(exception& e) {
162                 errorOut(e, "AlignCommand", "help");
163                 exit(1);
164         }
165 }
166
167
168 //**********************************************************************************************************************
169
170 int AlignCommand::execute(){
171         try {
172                 if (abort == true) {    return 0;       }
173                 
174                 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
175                 int longestBase = templateDB->getLongestBase();
176         
177                 if(align == "gotoh")                    {       alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase);                 }
178                 else if(align == "needleman")   {       alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);                                }
179                 else if(align == "blast")               {       alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch);            }
180                 else if(align == "noalign")             {       alignment = new NoAlign();                                                                                                      }
181                 else {
182                         mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
183                         mothurOutEndLine();
184                         alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
185                 }
186                 
187                 for (int s = 0; s < candidateFileNames.size(); s++) {
188                         mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); mothurOutEndLine();
189                         
190                         string alignFileName = candidateFileNames[s].substr(0,candidateFileNames[s].find_last_of(".")+1) + "align";
191                         string reportFileName = candidateFileNames[s].substr(0,candidateFileNames[s].find_last_of(".")+1) + "align.report";
192                         string accnosFileName = candidateFileNames[s].substr(0,candidateFileNames[s].find_last_of(".")+1) + "flip.accnos";
193                         
194                         int numFastaSeqs = 0;
195                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
196                         int start = time(NULL);
197                         
198 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
199                         if(processors == 1){
200                                 ifstream inFASTA;
201                                 openInputFile(candidateFileNames[s], inFASTA);
202                                 numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
203                                 inFASTA.close();
204                                 
205                                 lines.push_back(new linePair(0, numFastaSeqs));
206                                 
207                                 driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
208                                 
209                                 //delete accnos file if its blank else report to user
210                                 if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  }
211                                 else { 
212                                         mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
213                                         if (!flip) {
214                                                 mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
215                                         }else{  mothurOut(" If the reverse compliment proved to be better it was reported.");  }
216                                         mothurOutEndLine();
217                                 }
218                         }
219                         else{
220                                 vector<int> positions;
221                                 processIDS.resize(0);
222                                 
223                                 ifstream inFASTA;
224                                 openInputFile(candidateFileNames[s], inFASTA);
225                                 
226                                 string input;
227                                 while(!inFASTA.eof()){
228                                         input = getline(inFASTA);
229                                         if (input.length() != 0) {
230                                                 if(input[0] == '>'){    long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);  }
231                                         }
232                                 }
233                                 inFASTA.close();
234                                 
235                                 numFastaSeqs = positions.size();
236                                 
237                                 int numSeqsPerProcessor = numFastaSeqs / processors;
238                                 
239                                 for (int i = 0; i < processors; i++) {
240                                         long int startPos = positions[ i * numSeqsPerProcessor ];
241                                         if(i == processors - 1){
242                                                 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
243                                         }
244                                         lines.push_back(new linePair(startPos, numSeqsPerProcessor));
245                                 }
246                                 
247                                 createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]); 
248                                 
249                                 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
250                                 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
251                                 
252                                 //append alignment and report files
253                                 for(int i=1;i<processors;i++){
254                                         appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
255                                         remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
256                                         
257                                         appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
258                                         remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
259                                 }
260                                 
261                                 vector<string> nonBlankAccnosFiles;
262                                 //delete blank accnos files generated with multiple processes
263                                 for(int i=0;i<processors;i++){  
264                                         if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
265                                                 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
266                                         }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());  }
267                                 }
268                                 
269                                 //append accnos files
270                                 if (nonBlankAccnosFiles.size() != 0) { 
271                                         rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
272                                         
273                                         for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
274                                                 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
275                                                 remove(nonBlankAccnosFiles[h].c_str());
276                                         }
277                                         mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
278                                         if (!flip) {
279                                                 mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
280                                         }else{  mothurOut(" If the reverse compliment proved to be better it was reported.");  }
281                                         mothurOutEndLine();
282                                 }
283                         }
284 #else
285                         ifstream inFASTA;
286                         openInputFile(candidateFileName[s], inFASTA);
287                         numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
288                         inFASTA.close();
289                         
290                         lines.push_back(new linePair(0, numFastaSeqs));
291                         
292                         driver(lines[0], alignFileName, reportFileName, accnosFileName);
293                         
294                         //delete accnos file if its blank else report to user
295                         if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  }
296                         else { 
297                                 mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
298                                 if (!flip) {
299                                         mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
300                                 }else{  mothurOut(" If the reverse compliment proved to be better it was reported.");  }
301                                 mothurOutEndLine();
302                         }
303                         
304 #endif
305                         
306                         
307                         
308                         mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
309                         mothurOutEndLine();
310                         mothurOutEndLine();
311                 }
312                 
313                 return 0;
314         }
315         catch(exception& e) {
316                 errorOut(e, "AlignCommand", "execute");
317                 exit(1);
318         }
319 }
320
321 //**********************************************************************************************************************
322
323 int AlignCommand::driver(linePair* line, string alignFName, string reportFName, string accnosFName, string filename){
324         try {
325                 ofstream alignmentFile;
326                 openOutputFile(alignFName, alignmentFile);
327                 
328                 ofstream accnosFile;
329                 openOutputFile(accnosFName, accnosFile);
330                 
331                 NastReport report(reportFName);
332                 
333                 ifstream inFASTA;
334                 openInputFile(filename, inFASTA);
335
336                 inFASTA.seekg(line->start);
337                 
338                 for(int i=0;i<line->numSeqs;i++){
339                 
340                         Sequence* candidateSeq = new Sequence(inFASTA);  gobble(inFASTA);
341                         int origNumBases = candidateSeq->getNumBases();
342                         string originalUnaligned = candidateSeq->getUnaligned();
343                         int numBasesNeeded = origNumBases * threshold;
344         
345                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
346                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
347                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
348                                 }
349                                                                 
350                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
351                                 Sequence* templateSeq = &temp;
352                                 
353                                 float searchScore = templateDB->getSearchScore();
354                                                                 
355                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
356                                 Sequence* copy;
357                                 
358                                 Nast* nast2;
359                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
360                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
361                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
362                                                                                                 //so this bool tells you if you need to delete it
363                                                                                                 
364                                 //if there is a possibility that this sequence should be reversed
365                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
366                                         
367                                         string wasBetter = "";
368                                         //if the user wants you to try the reverse
369                                         if (flip) {
370                                                 //get reverse compliment
371                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
372                                                 copy->reverseComplement();
373                                                 
374                                                 //rerun alignment
375                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
376                                                 Sequence* templateSeq2 = &temp2;
377                                                 
378                                                 searchScore = templateDB->getSearchScore();
379                                                 
380                                                 nast2 = new Nast(alignment, copy, templateSeq2);
381                         
382                                                 //check if any better
383                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
384                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
385                                                         templateSeq = templateSeq2; 
386                                                         delete nast;
387                                                         nast = nast2;
388                                                         needToDeleteCopy = true;
389                                                 }else{  
390                                                         wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
391                                                         delete nast2;
392                                                         delete copy;    
393                                                 }
394                                         }
395                                         
396                                         //create accnos file with names
397                                         accnosFile << candidateSeq->getName() << wasBetter << endl;
398                                 }
399                                 
400                                 report.setCandidate(candidateSeq);
401                                 report.setTemplate(templateSeq);
402                                 report.setSearchParameters(search, searchScore);
403                                 report.setAlignmentParameters(align, alignment);
404                                 report.setNastParameters(*nast);
405         
406                                 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
407                                 
408                                 report.print();
409                                 delete nast;
410                                 if (needToDeleteCopy) {   delete copy;   }
411                         }
412                         delete candidateSeq;
413                         
414                         //report progress
415                         if((i+1) % 100 == 0){   mothurOut(toString(i+1)); mothurOutEndLine();           }
416                 }
417                 //report progress
418                 if((line->numSeqs) % 100 != 0){ mothurOut(toString(line->numSeqs)); mothurOutEndLine();         }
419                 
420                 alignmentFile.close();
421                 inFASTA.close();
422                 accnosFile.close();
423                 
424                 return 1;
425         }
426         catch(exception& e) {
427                 errorOut(e, "AlignCommand", "driver");
428                 exit(1);
429         }
430 }
431
432 /**************************************************************************************************/
433
434 void AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
435         try {
436 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
437                 int process = 0;
438                 //              processIDS.resize(0);
439                 
440                 //loop through and create all the processes you want
441                 while (process != processors) {
442                         int pid = fork();
443                         
444                         if (pid > 0) {
445                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
446                                 process++;
447                         }else if (pid == 0){
448                                 driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
449                                 exit(0);
450                         }else { mothurOut("unable to spawn the necessary processes."); mothurOutEndLine(); exit(0); }
451                 }
452                 
453                 //force parent to wait until all the processes are done
454                 for (int i=0;i<processors;i++) { 
455                         int temp = processIDS[i];
456                         wait(&temp);
457                 }
458 #endif          
459         }
460         catch(exception& e) {
461                 errorOut(e, "AlignCommand", "createProcesses");
462                 exit(1);
463         }
464 }
465
466 /**************************************************************************************************/
467
468 void AlignCommand::appendAlignFiles(string temp, string filename) {
469         try{
470                 
471                 ofstream output;
472                 ifstream input;
473                 openOutputFileAppend(filename, output);
474                 openInputFile(temp, input);
475                 
476                 while(char c = input.get()){
477                         if(input.eof())         {       break;                  }
478                         else                            {       output << c;    }
479                 }
480                 
481                 input.close();
482                 output.close();
483         }
484         catch(exception& e) {
485                 errorOut(e, "AlignCommand", "appendAlignFiles");
486                 exit(1);
487         }
488 }
489 //**********************************************************************************************************************
490
491 void AlignCommand::appendReportFiles(string temp, string filename) {
492         try{
493                 
494                 ofstream output;
495                 ifstream input;
496                 openOutputFileAppend(filename, output);
497                 openInputFile(temp, input);
498
499                 while (!input.eof())    {       char c = input.get(); if (c == 10 || c == 13){  break;  }       } // get header line
500                                 
501                 while(char c = input.get()){
502                         if(input.eof())         {       break;                  }
503                         else                            {       output << c;    }
504                 }
505                 
506                 input.close();
507                 output.close();
508         }
509         catch(exception& e) {
510                 errorOut(e, "AlignCommand", "appendReportFiles");
511                 exit(1);
512         }
513 }
514
515 //**********************************************************************************************************************