5 * Created by Sarah Westcott on 5/15/09.
6 * Copyright 2009 Schloss Lab UMASS Amherst. All rights reserved.
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).
17 #include "aligncommand.h"
18 #include "sequence.hpp"
20 #include "gotohoverlap.hpp"
21 #include "needlemanoverlap.hpp"
22 #include "blastalign.hpp"
23 #include "noalign.hpp"
26 #include "nastreport.hpp"
29 //**********************************************************************************************************************
31 AlignCommand::AlignCommand(string option) {
36 //allow user to run help
37 if(option == "help") { help(); abort = true; }
41 //valid paramters for this command
42 string AlignArray[] = {"template","candidate","search","ksize","align","match","mismatch","gapopen","gapextend", "processors","flip","threshold","outputdir","inputdir"};
43 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
45 OptionParser parser(option);
46 map<string, string> parameters = parser.getParameters();
48 ValidParameters validParameter;
49 map<string, string>::iterator it;
51 //check to make sure all parameters are valid for command
52 for (it = parameters.begin(); it != parameters.end(); it++) {
53 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
56 //if the user changes the output directory command factory will send this info to us in the output parameter
57 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
60 //if the user changes the input directory command factory will send this info to us in the output parameter
61 string inputDir = validParameter.validFile(parameters, "inputdir", false);
63 if (inputDir == "not found"){ inputDir = ""; }
67 it = parameters.find("template");
69 //user has given a template file
70 if(it != parameters.end()){
71 path = hasPath(it->second);
72 //if the user has not given a path then, add inputdir. else leave path alone.
73 if (path == "") { parameters["template"] = inputDir + it->second; }
77 //check for required parameters
78 templateFileName = validParameter.validFile(parameters, "template", true);
80 if (templateFileName == "not found") {
81 m->mothurOut("template is a required parameter for the align.seqs command.");
82 m->mothurOutEndLine();
84 }else if (templateFileName == "not open") { abort = true; }
86 candidateFileName = validParameter.validFile(parameters, "candidate", false);
87 if (candidateFileName == "not found") { m->mothurOut("candidate is a required parameter for the align.seqs command."); m->mothurOutEndLine(); abort = true; }
89 splitAtDash(candidateFileName, candidateFileNames);
91 //go through files and make sure they are good, if not, then disregard them
92 for (int i = 0; i < candidateFileNames.size(); i++) {
94 string path = hasPath(candidateFileNames[i]);
95 //if the user has not given a path then, add inputdir. else leave path alone.
96 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
101 ableToOpen = openInputFile(candidateFileNames[i], in);
102 if (ableToOpen == 1) {
103 m->mothurOut(candidateFileNames[i] + " will be disregarded."); m->mothurOutEndLine();
104 //erase from file list
105 candidateFileNames.erase(candidateFileNames.begin()+i);
111 //make sure there is at least one valid file left
112 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
115 //check for optional parameter and set defaults
116 // ...at some point should added some additional type checking...
118 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
119 convert(temp, kmerSize);
121 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
122 convert(temp, match);
124 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
125 convert(temp, misMatch);
127 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
128 convert(temp, gapOpen);
130 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
131 convert(temp, gapExtend);
133 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
134 convert(temp, processors);
136 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
139 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
140 convert(temp, threshold);
142 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
144 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
148 catch(exception& e) {
149 m->errorOut(e, "AlignCommand", "AlignCommand");
154 //**********************************************************************************************************************
156 AlignCommand::~AlignCommand(){
158 if (abort == false) {
159 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
165 //**********************************************************************************************************************
167 void AlignCommand::help(){
169 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
170 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen and gapextend.\n");
171 m->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");
172 m->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");
173 m->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");
174 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
175 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
176 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
177 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
178 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
179 m->mothurOut("The flip parameter is used to specify whether or not you want mothur to try the reverse complement if a sequence falls below the threshold. The default is false.\n");
180 m->mothurOut("The threshold is used to specify a cutoff at which an alignment is deemed 'bad' and the reverse complement may be tried. The default threshold is 0.50, meaning 50% of the bases are removed in the alignment.\n");
181 m->mothurOut("If the flip parameter is set to true the reverse complement of the sequence is aligned and the better alignment is reported.\n");
182 m->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");
183 m->mothurOut("The align.seqs command should be in the following format: \n");
184 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
185 m->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");
186 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
188 catch(exception& e) {
189 m->errorOut(e, "AlignCommand", "help");
195 //**********************************************************************************************************************
197 int AlignCommand::execute(){
199 if (abort == true) { return 0; }
201 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
202 int longestBase = templateDB->getLongestBase();
204 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
205 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
206 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
207 else if(align == "noalign") { alignment = new NoAlign(); }
209 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
210 m->mothurOutEndLine();
211 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
213 vector<string> outputNames;
215 for (int s = 0; s < candidateFileNames.size(); s++) {
216 if (m->control_pressed) { return 0; }
218 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
220 if (outputDir == "") { outputDir += hasPath(candidateFileNames[s]); }
221 string alignFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align";
222 string reportFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align.report";
223 string accnosFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "flip.accnos";
224 bool hasAccnos = true;
226 int numFastaSeqs = 0;
227 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
228 int start = time(NULL);
230 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
233 openInputFile(candidateFileNames[s], inFASTA);
234 numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
237 lines.push_back(new linePair(0, numFastaSeqs));
239 driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
241 if (m->control_pressed) {
242 remove(accnosFileName.c_str());
243 remove(alignFileName.c_str());
244 remove(reportFileName.c_str());
248 //delete accnos file if its blank else report to user
249 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
251 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
253 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
254 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
255 m->mothurOutEndLine();
259 vector<int> positions;
260 processIDS.resize(0);
263 openInputFile(candidateFileNames[s], inFASTA);
266 while(!inFASTA.eof()){
267 input = getline(inFASTA);
268 if (input.length() != 0) {
269 if(input[0] == '>'){ long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
274 numFastaSeqs = positions.size();
276 int numSeqsPerProcessor = numFastaSeqs / processors;
278 for (int i = 0; i < processors; i++) {
279 long int startPos = positions[ i * numSeqsPerProcessor ];
280 if(i == processors - 1){
281 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
283 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
286 createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
288 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
289 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
291 //append alignment and report files
292 for(int i=1;i<processors;i++){
293 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
294 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
296 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
297 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
300 vector<string> nonBlankAccnosFiles;
301 //delete blank accnos files generated with multiple processes
302 for(int i=0;i<processors;i++){
303 if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
304 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
305 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
308 //append accnos files
309 if (nonBlankAccnosFiles.size() != 0) {
310 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
312 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
313 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
314 remove(nonBlankAccnosFiles[h].c_str());
316 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
318 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
319 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
320 m->mothurOutEndLine();
321 }else{ hasAccnos = false; }
323 if (m->control_pressed) {
324 remove(accnosFileName.c_str());
325 remove(alignFileName.c_str());
326 remove(reportFileName.c_str());
332 openInputFile(candidateFileNames[s], inFASTA);
333 numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
336 lines.push_back(new linePair(0, numFastaSeqs));
338 driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
340 if (m->control_pressed) {
341 remove(accnosFileName.c_str());
342 remove(alignFileName.c_str());
343 remove(reportFileName.c_str());
347 //delete accnos file if its blank else report to user
348 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
350 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
352 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
353 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
354 m->mothurOutEndLine();
359 outputNames.push_back(alignFileName);
360 outputNames.push_back(reportFileName);
361 if (hasAccnos) { outputNames.push_back(accnosFileName); }
363 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
364 m->mothurOutEndLine();
365 m->mothurOutEndLine();
369 m->mothurOutEndLine();
370 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
371 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
372 m->mothurOutEndLine();
376 catch(exception& e) {
377 m->errorOut(e, "AlignCommand", "execute");
382 //**********************************************************************************************************************
384 int AlignCommand::driver(linePair* line, string alignFName, string reportFName, string accnosFName, string filename){
386 ofstream alignmentFile;
387 openOutputFile(alignFName, alignmentFile);
390 openOutputFile(accnosFName, accnosFile);
392 NastReport report(reportFName);
395 openInputFile(filename, inFASTA);
397 inFASTA.seekg(line->start);
399 for(int i=0;i<line->numSeqs;i++){
401 if (m->control_pressed) { return 0; }
403 Sequence* candidateSeq = new Sequence(inFASTA); gobble(inFASTA);
404 int origNumBases = candidateSeq->getNumBases();
405 string originalUnaligned = candidateSeq->getUnaligned();
406 int numBasesNeeded = origNumBases * threshold;
408 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
409 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
410 alignment->resize(candidateSeq->getUnaligned().length()+1);
413 Sequence temp = templateDB->findClosestSequence(candidateSeq);
414 Sequence* templateSeq = &temp;
416 float searchScore = templateDB->getSearchScore();
418 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
422 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
423 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
424 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
425 //so this bool tells you if you need to delete it
427 //if there is a possibility that this sequence should be reversed
428 if (candidateSeq->getNumBases() < numBasesNeeded) {
430 string wasBetter = "";
431 //if the user wants you to try the reverse
433 //get reverse compliment
434 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
435 copy->reverseComplement();
438 Sequence temp2 = templateDB->findClosestSequence(copy);
439 Sequence* templateSeq2 = &temp2;
441 searchScore = templateDB->getSearchScore();
443 nast2 = new Nast(alignment, copy, templateSeq2);
445 //check if any better
446 if (copy->getNumBases() > candidateSeq->getNumBases()) {
447 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
448 templateSeq = templateSeq2;
451 needToDeleteCopy = true;
453 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
459 //create accnos file with names
460 accnosFile << candidateSeq->getName() << wasBetter << endl;
463 report.setCandidate(candidateSeq);
464 report.setTemplate(templateSeq);
465 report.setSearchParameters(search, searchScore);
466 report.setAlignmentParameters(align, alignment);
467 report.setNastParameters(*nast);
469 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
473 if (needToDeleteCopy) { delete copy; }
478 if((i+1) % 100 == 0){ m->mothurOut(toString(i+1)); m->mothurOutEndLine(); }
481 if((line->numSeqs) % 100 != 0){ m->mothurOut(toString(line->numSeqs)); m->mothurOutEndLine(); }
483 alignmentFile.close();
489 catch(exception& e) {
490 m->errorOut(e, "AlignCommand", "driver");
495 /**************************************************************************************************/
497 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
499 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
502 // processIDS.resize(0);
504 //loop through and create all the processes you want
505 while (process != processors) {
509 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
512 exitCommand = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
514 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
517 //force parent to wait until all the processes are done
518 for (int i=0;i<processors;i++) {
519 int temp = processIDS[i];
526 catch(exception& e) {
527 m->errorOut(e, "AlignCommand", "createProcesses");
532 /**************************************************************************************************/
534 void AlignCommand::appendAlignFiles(string temp, string filename) {
539 openOutputFileAppend(filename, output);
540 openInputFile(temp, input);
542 while(char c = input.get()){
543 if(input.eof()) { break; }
544 else { output << c; }
550 catch(exception& e) {
551 m->errorOut(e, "AlignCommand", "appendAlignFiles");
555 //**********************************************************************************************************************
557 void AlignCommand::appendReportFiles(string temp, string filename) {
562 openOutputFileAppend(filename, output);
563 openInputFile(temp, input);
565 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
567 while(char c = input.get()){
568 if(input.eof()) { break; }
569 else { output << c; }
575 catch(exception& e) {
576 m->errorOut(e, "AlignCommand", "appendReportFiles");
581 //**********************************************************************************************************************