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"
28 //**********************************************************************************************************************
29 vector<string> AlignCommand::setParameters(){
31 CommandParameter ptemplate("reference", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(ptemplate);
32 CommandParameter pcandidate("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pcandidate);
33 CommandParameter psearch("search", "Multiple", "kmer-blast-suffix", "kmer", "", "", "",false,false); parameters.push_back(psearch);
34 CommandParameter pksize("ksize", "Number", "", "8", "", "", "",false,false); parameters.push_back(pksize);
35 CommandParameter pmatch("match", "Number", "", "1.0", "", "", "",false,false); parameters.push_back(pmatch);
36 CommandParameter palign("align", "Multiple", "needleman-gotoh-blast-noalign", "needleman", "", "", "",false,false); parameters.push_back(palign);
37 CommandParameter pmismatch("mismatch", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pmismatch);
38 CommandParameter pgapopen("gapopen", "Number", "", "-2.0", "", "", "",false,false); parameters.push_back(pgapopen);
39 CommandParameter pgapextend("gapextend", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pgapextend);
40 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
41 CommandParameter pflip("flip", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(pflip);
42 CommandParameter pthreshold("threshold", "Number", "", "0.50", "", "", "",false,false); parameters.push_back(pthreshold);
43 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
44 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
46 vector<string> myArray;
47 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
51 m->errorOut(e, "AlignCommand", "setParameters");
55 //**********************************************************************************************************************
56 string AlignCommand::getHelpString(){
58 string helpString = "";
59 helpString += "The align.seqs command reads a file containing sequences and creates an alignment file and a report file.";
60 helpString += "The align.seqs command parameters are reference, fasta, search, ksize, align, match, mismatch, gapopen, gapextend and processors.";
61 helpString += "The reference and fasta parameters are required. You may leave fasta blank if you have a valid fasta file. You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amzon.fasta.";
62 helpString += "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.";
63 helpString += "The align parameter allows you to specify the alignment method to use. Your options are: gotoh, needleman, blast and noalign. The default is needleman.";
64 helpString += "The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.";
65 helpString += "The match parameter allows you to specify the bonus for having the same base. The default is 1.0.";
66 helpString += "The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.";
67 helpString += "The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.";
68 helpString += "The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.";
69 helpString += "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.";
70 helpString += "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.";
71 helpString += "If the flip parameter is set to true the reverse complement of the sequence is aligned and the better alignment is reported.";
72 helpString += "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.";
73 helpString += "The align.seqs command should be in the following format:";
74 helpString += "align.seqs(reference=yourTemplateFile, fasta=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty)";
75 helpString += "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)";
76 helpString += "Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).";
80 m->errorOut(e, "AlignCommand", "getHelpString");
84 //**********************************************************************************************************************
85 AlignCommand::AlignCommand(){
87 abort = true; calledHelp = true;
89 vector<string> tempOutNames;
90 outputTypes["fasta"] = tempOutNames;
91 outputTypes["alignreport"] = tempOutNames;
92 outputTypes["accnos"] = tempOutNames;
95 m->errorOut(e, "AlignCommand", "AlignCommand");
99 //**********************************************************************************************************************
100 AlignCommand::AlignCommand(string option) {
102 abort = false; calledHelp = false;
104 //allow user to run help
105 if(option == "help") { help(); abort = true; calledHelp = true;}
108 vector<string> myArray = setParameters();
110 OptionParser parser(option);
111 map<string, string> parameters = parser.getParameters();
113 ValidParameters validParameter("align.seqs");
114 map<string, string>::iterator it;
116 //check to make sure all parameters are valid for command
117 for (it = parameters.begin(); it != parameters.end(); it++) {
118 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
121 //initialize outputTypes
122 vector<string> tempOutNames;
123 outputTypes["fasta"] = tempOutNames;
124 outputTypes["alignreport"] = tempOutNames;
125 outputTypes["accnos"] = tempOutNames;
127 //if the user changes the output directory command factory will send this info to us in the output parameter
128 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
131 //if the user changes the input directory command factory will send this info to us in the output parameter
132 string inputDir = validParameter.validFile(parameters, "inputdir", false);
134 if (inputDir == "not found"){ inputDir = ""; }
138 it = parameters.find("reference");
140 //user has given a template file
141 if(it != parameters.end()){
142 path = m->hasPath(it->second);
143 //if the user has not given a path then, add inputdir. else leave path alone.
144 if (path == "") { parameters["reference"] = inputDir + it->second; }
148 //check for required parameters
149 templateFileName = validParameter.validFile(parameters, "reference", true);
151 if (templateFileName == "not found") {
152 m->mothurOut("reference is a required parameter for the align.seqs command.");
153 m->mothurOutEndLine();
155 }else if (templateFileName == "not open") { abort = true; }
157 candidateFileName = validParameter.validFile(parameters, "fasta", false);
158 if (candidateFileName == "not found") {
159 //if there is a current fasta file, use it
160 string filename = m->getFastaFile();
161 if (filename != "") { candidateFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
162 else { m->mothurOut("You have no current fastafile and the candidate parameter is required."); m->mothurOutEndLine(); abort = true; }
164 m->splitAtDash(candidateFileName, candidateFileNames);
166 //go through files and make sure they are good, if not, then disregard them
167 for (int i = 0; i < candidateFileNames.size(); i++) {
168 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
171 if (candidateFileNames[i] == "current") {
172 candidateFileNames[i] = m->getFastaFile();
173 if (candidateFileNames[i] != "") { m->mothurOut("Using " + candidateFileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
175 m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true;
176 //erase from file list
177 candidateFileNames.erase(candidateFileNames.begin()+i);
184 if (inputDir != "") {
185 string path = m->hasPath(candidateFileNames[i]);
186 //if the user has not given a path then, add inputdir. else leave path alone.
187 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
192 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
195 //if you can't open it, try default location
196 if (ableToOpen == 1) {
197 if (m->getDefaultPath() != "") { //default path is set
198 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
199 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
201 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
203 candidateFileNames[i] = tryPath;
207 //if you can't open it, try output location
208 if (ableToOpen == 1) {
209 if (m->getOutputDir() != "") { //default path is set
210 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
211 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
213 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
215 candidateFileNames[i] = tryPath;
221 if (ableToOpen == 1) {
222 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
223 //erase from file list
224 candidateFileNames.erase(candidateFileNames.begin()+i);
230 //make sure there is at least one valid file left
231 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
234 //check for optional parameter and set defaults
235 // ...at some point should added some additional type checking...
237 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
238 convert(temp, kmerSize);
240 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
241 convert(temp, match);
243 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
244 convert(temp, misMatch);
246 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
247 convert(temp, gapOpen);
249 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
250 convert(temp, gapExtend);
252 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
253 m->setProcessors(temp);
254 convert(temp, processors);
256 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
257 flip = m->isTrue(temp);
259 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
260 convert(temp, threshold);
262 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
263 if ((search != "suffix") && (search != "kmer") && (search != "blast")) { m->mothurOut("invalid search option: choices are kmer, suffix or blast."); m->mothurOutEndLine(); abort=true; }
265 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
266 if ((align != "needleman") && (align != "gotoh") && (align != "blast") && (align != "noalign")) { m->mothurOut("invalid align option: choices are needleman, gotoh, blast or noalign."); m->mothurOutEndLine(); abort=true; }
271 catch(exception& e) {
272 m->errorOut(e, "AlignCommand", "AlignCommand");
276 //**********************************************************************************************************************
277 AlignCommand::~AlignCommand(){
279 if (abort == false) {
280 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
285 //**********************************************************************************************************************
287 int AlignCommand::execute(){
289 if (abort == true) { if (calledHelp) { return 0; } return 2; }
291 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
292 int longestBase = templateDB->getLongestBase();
294 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
295 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
296 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
297 else if(align == "noalign") { alignment = new NoAlign(); }
299 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
300 m->mothurOutEndLine();
301 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
304 for (int s = 0; s < candidateFileNames.size(); s++) {
305 if (m->control_pressed) { outputTypes.clear(); return 0; }
307 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
309 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
310 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
311 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
312 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
313 bool hasAccnos = true;
315 int numFastaSeqs = 0;
316 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
317 int start = time(NULL);
320 int pid, numSeqsPerProcessor;
322 vector<unsigned long int> MPIPos;
323 MPIWroteAccnos = false;
326 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
327 MPI_Comm_size(MPI_COMM_WORLD, &processors);
330 MPI_File outMPIAlign;
331 MPI_File outMPIReport;
332 MPI_File outMPIAccnos;
334 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
335 int inMode=MPI_MODE_RDONLY;
337 char outAlignFilename[1024];
338 strcpy(outAlignFilename, alignFileName.c_str());
340 char outReportFilename[1024];
341 strcpy(outReportFilename, reportFileName.c_str());
343 char outAccnosFilename[1024];
344 strcpy(outAccnosFilename, accnosFileName.c_str());
346 char inFileName[1024];
347 strcpy(inFileName, candidateFileNames[s].c_str());
349 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
350 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
351 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
352 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
354 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
356 if (pid == 0) { //you are the root process
358 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
360 //send file positions to all processes
361 for(int i = 1; i < processors; i++) {
362 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
363 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
366 //figure out how many sequences you have to align
367 numSeqsPerProcessor = numFastaSeqs / processors;
368 int startIndex = pid * numSeqsPerProcessor;
369 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
372 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
374 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
376 for (int i = 1; i < processors; i++) {
378 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
379 if (tempResult != 0) { MPIWroteAccnos = true; }
381 }else{ //you are a child process
382 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
383 MPIPos.resize(numFastaSeqs+1);
384 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
387 //figure out how many sequences you have to align
388 numSeqsPerProcessor = numFastaSeqs / processors;
389 int startIndex = pid * numSeqsPerProcessor;
390 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
394 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
396 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
398 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
402 MPI_File_close(&inMPI);
403 MPI_File_close(&outMPIAlign);
404 MPI_File_close(&outMPIReport);
405 MPI_File_close(&outMPIAccnos);
407 //delete accnos file if blank
409 //delete accnos file if its blank else report to user
410 if (MPIWroteAccnos) {
411 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
413 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
414 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
415 m->mothurOutEndLine();
418 //MPI_File_delete(outAccnosFilename, info);
420 remove(accnosFileName.c_str());
426 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
427 for (int i = 0; i < (positions.size()-1); i++) {
428 lines.push_back(new linePair(positions[i], positions[(i+1)]));
430 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
432 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
434 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
437 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
440 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
442 //delete accnos file if its blank else report to user
443 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
445 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
447 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
448 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
449 m->mothurOutEndLine();
456 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
458 if (pid == 0) { //only one process should output to screen
461 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
462 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
463 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
469 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
470 m->mothurOutEndLine();
471 m->mothurOutEndLine();
474 //set align file as new current fastafile
475 string currentFasta = "";
476 itTypes = outputTypes.find("fasta");
477 if (itTypes != outputTypes.end()) {
478 if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; m->setFastaFile(currentFasta); }
481 m->mothurOutEndLine();
482 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
483 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
484 m->mothurOutEndLine();
488 catch(exception& e) {
489 m->errorOut(e, "AlignCommand", "execute");
494 //**********************************************************************************************************************
496 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
498 ofstream alignmentFile;
499 m->openOutputFile(alignFName, alignmentFile);
502 m->openOutputFile(accnosFName, accnosFile);
504 NastReport report(reportFName);
507 m->openInputFile(filename, inFASTA);
509 inFASTA.seekg(filePos->start);
516 if (m->control_pressed) { return 0; }
518 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
519 report.setCandidate(candidateSeq);
521 int origNumBases = candidateSeq->getNumBases();
522 string originalUnaligned = candidateSeq->getUnaligned();
523 int numBasesNeeded = origNumBases * threshold;
525 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
526 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
527 alignment->resize(candidateSeq->getUnaligned().length()+1);
530 Sequence temp = templateDB->findClosestSequence(candidateSeq);
531 Sequence* templateSeq = &temp;
533 float searchScore = templateDB->getSearchScore();
535 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
540 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
541 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
542 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
543 //so this bool tells you if you need to delete it
545 //if there is a possibility that this sequence should be reversed
546 if (candidateSeq->getNumBases() < numBasesNeeded) {
548 string wasBetter = "";
549 //if the user wants you to try the reverse
552 //get reverse compliment
553 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
554 copy->reverseComplement();
557 Sequence temp2 = templateDB->findClosestSequence(copy);
558 Sequence* templateSeq2 = &temp2;
560 searchScore = templateDB->getSearchScore();
562 nast2 = new Nast(alignment, copy, templateSeq2);
564 //check if any better
565 if (copy->getNumBases() > candidateSeq->getNumBases()) {
566 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
567 templateSeq = templateSeq2;
570 needToDeleteCopy = true;
571 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
573 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
579 //create accnos file with names
580 accnosFile << candidateSeq->getName() << wasBetter << endl;
583 report.setTemplate(templateSeq);
584 report.setSearchParameters(search, searchScore);
585 report.setAlignmentParameters(align, alignment);
586 report.setNastParameters(*nast);
588 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
592 if (needToDeleteCopy) { delete copy; }
598 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
599 unsigned long int pos = inFASTA.tellg();
600 if ((pos == -1) || (pos >= filePos->end)) { break; }
602 if (inFASTA.eof()) { break; }
606 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
610 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
612 alignmentFile.close();
618 catch(exception& e) {
619 m->errorOut(e, "AlignCommand", "driver");
623 //**********************************************************************************************************************
625 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
627 string outputString = "";
628 MPI_Status statusReport;
629 MPI_Status statusAlign;
630 MPI_Status statusAccnos;
633 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
638 outputString = report.getHeaders();
639 int length = outputString.length();
641 char* buf = new char[length];
642 memcpy(buf, outputString.c_str(), length);
644 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
649 for(int i=0;i<num;i++){
651 if (m->control_pressed) { return 0; }
654 int length = MPIPos[start+i+1] - MPIPos[start+i];
656 char* buf4 = new char[length];
657 //memcpy(buf4, outputString.c_str(), length);
659 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
661 string tempBuf = buf4;
665 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
667 istringstream iss (tempBuf,istringstream::in);
669 Sequence* candidateSeq = new Sequence(iss);
670 report.setCandidate(candidateSeq);
672 int origNumBases = candidateSeq->getNumBases();
673 string originalUnaligned = candidateSeq->getUnaligned();
674 int numBasesNeeded = origNumBases * threshold;
676 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
677 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
678 alignment->resize(candidateSeq->getUnaligned().length()+1);
681 Sequence temp = templateDB->findClosestSequence(candidateSeq);
682 Sequence* templateSeq = &temp;
684 float searchScore = templateDB->getSearchScore();
686 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
690 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
691 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
692 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
693 //so this bool tells you if you need to delete it
695 //if there is a possibility that this sequence should be reversed
696 if (candidateSeq->getNumBases() < numBasesNeeded) {
698 string wasBetter = "";
699 //if the user wants you to try the reverse
701 //get reverse compliment
702 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
703 copy->reverseComplement();
706 Sequence temp2 = templateDB->findClosestSequence(copy);
707 Sequence* templateSeq2 = &temp2;
709 searchScore = templateDB->getSearchScore();
711 nast2 = new Nast(alignment, copy, templateSeq2);
713 //check if any better
714 if (copy->getNumBases() > candidateSeq->getNumBases()) {
715 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
716 templateSeq = templateSeq2;
719 needToDeleteCopy = true;
720 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
722 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
728 //create accnos file with names
729 outputString = candidateSeq->getName() + wasBetter + "\n";
731 //send results to parent
732 int length = outputString.length();
734 char* buf = new char[length];
735 memcpy(buf, outputString.c_str(), length);
737 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
739 MPIWroteAccnos = true;
742 report.setTemplate(templateSeq);
743 report.setSearchParameters(search, searchScore);
744 report.setAlignmentParameters(align, alignment);
745 report.setNastParameters(*nast);
747 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
749 //send results to parent
750 int length = outputString.length();
751 char* buf2 = new char[length];
752 memcpy(buf2, outputString.c_str(), length);
754 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
758 outputString = report.getReport();
760 //send results to parent
761 length = outputString.length();
762 char* buf3 = new char[length];
763 memcpy(buf3, outputString.c_str(), length);
765 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
769 if (needToDeleteCopy) { delete copy; }
774 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
777 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
781 catch(exception& e) {
782 m->errorOut(e, "AlignCommand", "driverMPI");
787 /**************************************************************************************************/
789 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
791 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
792 processIDS.resize(0);
795 // processIDS.resize(0);
797 //loop through and create all the processes you want
798 while (process != processors) {
802 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
805 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
807 //pass numSeqs to parent
809 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
810 m->openOutputFile(tempFile, out);
816 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
817 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
823 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
825 //force parent to wait until all the processes are done
826 for (int i=0;i<processIDS.size();i++) {
827 int temp = processIDS[i];
831 vector<string> nonBlankAccnosFiles;
832 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
833 else { remove(accnosFName.c_str()); } //remove so other files can be renamed to it
835 for (int i = 0; i < processIDS.size(); i++) {
837 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
838 m->openInputFile(tempFile, in);
839 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
840 in.close(); remove(tempFile.c_str());
842 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
843 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
845 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
846 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
848 if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
849 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
850 }else { remove((accnosFName + toString(processIDS[i]) + ".temp").c_str()); }
854 //append accnos files
855 if (nonBlankAccnosFiles.size() != 0) {
856 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
858 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
859 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
860 remove(nonBlankAccnosFiles[h].c_str());
862 }else { //recreate the accnosfile if needed
864 m->openOutputFile(accnosFName, out);
871 catch(exception& e) {
872 m->errorOut(e, "AlignCommand", "createProcesses");
876 /**************************************************************************************************/
878 void AlignCommand::appendAlignFiles(string temp, string filename) {
883 m->openOutputFileAppend(filename, output);
884 m->openInputFile(temp, input);
886 while(char c = input.get()){
887 if(input.eof()) { break; }
888 else { output << c; }
894 catch(exception& e) {
895 m->errorOut(e, "AlignCommand", "appendAlignFiles");
899 //**********************************************************************************************************************
901 void AlignCommand::appendReportFiles(string temp, string filename) {
906 m->openOutputFileAppend(filename, output);
907 m->openInputFile(temp, input);
909 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
911 while(char c = input.get()){
912 if(input.eof()) { break; }
913 else { output << c; }
919 catch(exception& e) {
920 m->errorOut(e, "AlignCommand", "appendReportFiles");
924 //**********************************************************************************************************************