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 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
218 if (outputDir == "") { outputDir += hasPath(candidateFileNames[s]); }
219 string alignFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align";
220 string reportFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align.report";
221 string accnosFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "flip.accnos";
222 bool hasAccnos = true;
224 int numFastaSeqs = 0;
225 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
226 int start = time(NULL);
228 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
231 openInputFile(candidateFileNames[s], inFASTA);
232 numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
235 lines.push_back(new linePair(0, numFastaSeqs));
237 int exitCommand = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
238 if (exitCommand == 0) {
239 remove(accnosFileName.c_str());
240 remove(alignFileName.c_str());
241 remove(reportFileName.c_str());
245 //delete accnos file if its blank else report to user
246 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
248 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
250 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
251 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
252 m->mothurOutEndLine();
256 vector<int> positions;
257 processIDS.resize(0);
260 openInputFile(candidateFileNames[s], inFASTA);
263 while(!inFASTA.eof()){
264 input = getline(inFASTA);
265 if (input.length() != 0) {
266 if(input[0] == '>'){ long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
271 numFastaSeqs = positions.size();
273 int numSeqsPerProcessor = numFastaSeqs / processors;
275 for (int i = 0; i < processors; i++) {
276 long int startPos = positions[ i * numSeqsPerProcessor ];
277 if(i == processors - 1){
278 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
280 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
283 int exitCommand = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
285 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
286 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
288 //append alignment and report files
289 for(int i=1;i<processors;i++){
290 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
291 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
293 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
294 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
297 vector<string> nonBlankAccnosFiles;
298 //delete blank accnos files generated with multiple processes
299 for(int i=0;i<processors;i++){
300 if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
301 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
302 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
305 //append accnos files
306 if (nonBlankAccnosFiles.size() != 0) {
307 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
309 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
310 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
311 remove(nonBlankAccnosFiles[h].c_str());
313 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
315 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
316 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
317 m->mothurOutEndLine();
318 }else{ hasAccnos = false; }
320 if (exitCommand == 0) {
321 remove(accnosFileName.c_str());
322 remove(alignFileName.c_str());
323 remove(reportFileName.c_str());
329 openInputFile(candidateFileNames[s], inFASTA);
330 numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
333 lines.push_back(new linePair(0, numFastaSeqs));
335 int exitCommand = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
336 if (exitCommand == 0) {
337 remove(accnosFileName.c_str());
338 remove(alignFileName.c_str());
339 remove(reportFileName.c_str());
343 //delete accnos file if its blank else report to user
344 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
346 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
348 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
349 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
350 m->mothurOutEndLine();
355 outputNames.push_back(alignFileName);
356 outputNames.push_back(reportFileName);
357 if (hasAccnos) { outputNames.push_back(accnosFileName); }
359 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
360 m->mothurOutEndLine();
361 m->mothurOutEndLine();
365 m->mothurOutEndLine();
366 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
367 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
368 m->mothurOutEndLine();
372 catch(exception& e) {
373 m->errorOut(e, "AlignCommand", "execute");
378 //**********************************************************************************************************************
380 int AlignCommand::driver(linePair* line, string alignFName, string reportFName, string accnosFName, string filename){
382 ofstream alignmentFile;
383 openOutputFile(alignFName, alignmentFile);
386 openOutputFile(accnosFName, accnosFile);
388 NastReport report(reportFName);
391 openInputFile(filename, inFASTA);
393 inFASTA.seekg(line->start);
395 for(int i=0;i<line->numSeqs;i++){
397 if (m->control_pressed) { return 0; }
399 Sequence* candidateSeq = new Sequence(inFASTA); gobble(inFASTA);
400 int origNumBases = candidateSeq->getNumBases();
401 string originalUnaligned = candidateSeq->getUnaligned();
402 int numBasesNeeded = origNumBases * threshold;
404 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
405 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
406 alignment->resize(candidateSeq->getUnaligned().length()+1);
409 Sequence temp = templateDB->findClosestSequence(candidateSeq);
410 Sequence* templateSeq = &temp;
412 float searchScore = templateDB->getSearchScore();
414 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
418 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
419 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
420 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
421 //so this bool tells you if you need to delete it
423 //if there is a possibility that this sequence should be reversed
424 if (candidateSeq->getNumBases() < numBasesNeeded) {
426 string wasBetter = "";
427 //if the user wants you to try the reverse
429 //get reverse compliment
430 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
431 copy->reverseComplement();
434 Sequence temp2 = templateDB->findClosestSequence(copy);
435 Sequence* templateSeq2 = &temp2;
437 searchScore = templateDB->getSearchScore();
439 nast2 = new Nast(alignment, copy, templateSeq2);
441 //check if any better
442 if (copy->getNumBases() > candidateSeq->getNumBases()) {
443 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
444 templateSeq = templateSeq2;
447 needToDeleteCopy = true;
449 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
455 //create accnos file with names
456 accnosFile << candidateSeq->getName() << wasBetter << endl;
459 report.setCandidate(candidateSeq);
460 report.setTemplate(templateSeq);
461 report.setSearchParameters(search, searchScore);
462 report.setAlignmentParameters(align, alignment);
463 report.setNastParameters(*nast);
465 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
469 if (needToDeleteCopy) { delete copy; }
474 if((i+1) % 100 == 0){ m->mothurOut(toString(i+1)); m->mothurOutEndLine(); }
477 if((line->numSeqs) % 100 != 0){ m->mothurOut(toString(line->numSeqs)); m->mothurOutEndLine(); }
479 alignmentFile.close();
485 catch(exception& e) {
486 m->errorOut(e, "AlignCommand", "driver");
491 /**************************************************************************************************/
493 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
495 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
498 // processIDS.resize(0);
500 //loop through and create all the processes you want
501 while (process != processors) {
505 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
508 exitCommand = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
510 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
513 //force parent to wait until all the processes are done
514 for (int i=0;i<processors;i++) {
515 int temp = processIDS[i];
522 catch(exception& e) {
523 m->errorOut(e, "AlignCommand", "createProcesses");
528 /**************************************************************************************************/
530 void AlignCommand::appendAlignFiles(string temp, string filename) {
535 openOutputFileAppend(filename, output);
536 openInputFile(temp, input);
538 while(char c = input.get()){
539 if(input.eof()) { break; }
540 else { output << c; }
546 catch(exception& e) {
547 m->errorOut(e, "AlignCommand", "appendAlignFiles");
551 //**********************************************************************************************************************
553 void AlignCommand::appendReportFiles(string temp, string filename) {
558 openOutputFileAppend(filename, output);
559 openInputFile(temp, input);
561 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
563 while(char c = input.get()){
564 if(input.eof()) { break; }
565 else { output << c; }
571 catch(exception& e) {
572 m->errorOut(e, "AlignCommand", "appendReportFiles");
577 //**********************************************************************************************************************