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"
27 #include "referencedb.h"
29 //**********************************************************************************************************************
30 vector<string> AlignCommand::setParameters(){
32 CommandParameter ptemplate("reference", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(ptemplate);
33 CommandParameter pcandidate("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pcandidate);
34 CommandParameter psearch("search", "Multiple", "kmer-blast-suffix", "kmer", "", "", "",false,false); parameters.push_back(psearch);
35 CommandParameter pksize("ksize", "Number", "", "8", "", "", "",false,false); parameters.push_back(pksize);
36 CommandParameter pmatch("match", "Number", "", "1.0", "", "", "",false,false); parameters.push_back(pmatch);
37 CommandParameter palign("align", "Multiple", "needleman-gotoh-blast-noalign", "needleman", "", "", "",false,false); parameters.push_back(palign);
38 CommandParameter pmismatch("mismatch", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pmismatch);
39 CommandParameter pgapopen("gapopen", "Number", "", "-2.0", "", "", "",false,false); parameters.push_back(pgapopen);
40 CommandParameter pgapextend("gapextend", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pgapextend);
41 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
42 CommandParameter pflip("flip", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(pflip);
43 CommandParameter psave("save", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(psave);
44 CommandParameter pthreshold("threshold", "Number", "", "0.50", "", "", "",false,false); parameters.push_back(pthreshold);
45 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
46 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
48 vector<string> myArray;
49 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
53 m->errorOut(e, "AlignCommand", "setParameters");
57 //**********************************************************************************************************************
58 string AlignCommand::getHelpString(){
60 string helpString = "";
61 helpString += "The align.seqs command reads a file containing sequences and creates an alignment file and a report file.";
62 helpString += "The align.seqs command parameters are reference, fasta, search, ksize, align, match, mismatch, gapopen, gapextend and processors.";
63 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.";
64 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.";
65 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.";
66 helpString += "The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.";
67 helpString += "The match parameter allows you to specify the bonus for having the same base. The default is 1.0.";
68 helpString += "The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.";
69 helpString += "The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.";
70 helpString += "The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.";
71 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.";
72 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.";
73 helpString += "If the flip parameter is set to true the reverse complement of the sequence is aligned and the better alignment is reported.";
74 helpString += "If the save parameter is set to true the reference sequences will be saved in memory, to clear them later you can use the clear.memory command. Default=f.";
75 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.";
76 helpString += "The align.seqs command should be in the following format:";
77 helpString += "align.seqs(reference=yourTemplateFile, fasta=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty)";
78 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)";
79 helpString += "Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).";
83 m->errorOut(e, "AlignCommand", "getHelpString");
87 //**********************************************************************************************************************
88 AlignCommand::AlignCommand(){
90 abort = true; calledHelp = true;
92 vector<string> tempOutNames;
93 outputTypes["fasta"] = tempOutNames;
94 outputTypes["alignreport"] = tempOutNames;
95 outputTypes["accnos"] = tempOutNames;
98 m->errorOut(e, "AlignCommand", "AlignCommand");
102 //**********************************************************************************************************************
103 AlignCommand::AlignCommand(string option) {
105 abort = false; calledHelp = false;
106 ReferenceDB* rdb = ReferenceDB::getInstance();
108 //allow user to run help
109 if(option == "help") { help(); abort = true; calledHelp = true;}
110 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
113 vector<string> myArray = setParameters();
115 OptionParser parser(option);
116 map<string, string> parameters = parser.getParameters();
118 ValidParameters validParameter("align.seqs");
119 map<string, string>::iterator it;
121 //check to make sure all parameters are valid for command
122 for (it = parameters.begin(); it != parameters.end(); it++) {
123 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
126 //initialize outputTypes
127 vector<string> tempOutNames;
128 outputTypes["fasta"] = tempOutNames;
129 outputTypes["alignreport"] = tempOutNames;
130 outputTypes["accnos"] = tempOutNames;
132 //if the user changes the output directory command factory will send this info to us in the output parameter
133 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
136 //if the user changes the input directory command factory will send this info to us in the output parameter
137 string inputDir = validParameter.validFile(parameters, "inputdir", false);
139 if (inputDir == "not found"){ inputDir = ""; }
143 it = parameters.find("reference");
145 //user has given a template file
146 if(it != parameters.end()){
147 path = m->hasPath(it->second);
148 //if the user has not given a path then, add inputdir. else leave path alone.
149 if (path == "") { parameters["reference"] = inputDir + it->second; }
153 candidateFileName = validParameter.validFile(parameters, "fasta", false);
154 if (candidateFileName == "not found") {
155 //if there is a current fasta file, use it
156 string filename = m->getFastaFile();
157 if (filename != "") { candidateFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
158 else { m->mothurOut("You have no current fastafile and the candidate parameter is required."); m->mothurOutEndLine(); abort = true; }
160 m->splitAtDash(candidateFileName, candidateFileNames);
162 //go through files and make sure they are good, if not, then disregard them
163 for (int i = 0; i < candidateFileNames.size(); i++) {
164 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
167 if (candidateFileNames[i] == "current") {
168 candidateFileNames[i] = m->getFastaFile();
169 if (candidateFileNames[i] != "") { m->mothurOut("Using " + candidateFileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
171 m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true;
172 //erase from file list
173 candidateFileNames.erase(candidateFileNames.begin()+i);
180 if (inputDir != "") {
181 string path = m->hasPath(candidateFileNames[i]);
182 //if the user has not given a path then, add inputdir. else leave path alone.
183 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
188 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
191 //if you can't open it, try default location
192 if (ableToOpen == 1) {
193 if (m->getDefaultPath() != "") { //default path is set
194 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
195 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
197 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
199 candidateFileNames[i] = tryPath;
203 //if you can't open it, try output location
204 if (ableToOpen == 1) {
205 if (m->getOutputDir() != "") { //default path is set
206 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
207 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
209 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
211 candidateFileNames[i] = tryPath;
217 if (ableToOpen == 1) {
218 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
219 //erase from file list
220 candidateFileNames.erase(candidateFileNames.begin()+i);
223 m->setFastaFile(candidateFileNames[i]);
228 //make sure there is at least one valid file left
229 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
232 //check for optional parameter and set defaults
233 // ...at some point should added some additional type checking...
235 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
236 convert(temp, kmerSize);
238 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
239 convert(temp, match);
241 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
242 convert(temp, misMatch);
244 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
245 convert(temp, gapOpen);
247 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
248 convert(temp, gapExtend);
250 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
251 m->setProcessors(temp);
252 convert(temp, processors);
254 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
255 flip = m->isTrue(temp);
257 temp = validParameter.validFile(parameters, "save", false); if (temp == "not found"){ temp = "f"; }
258 save = m->isTrue(temp);
260 if (save) { //clear out old references
264 //this has to go after save so that if the user sets save=t and provides no reference we abort
265 templateFileName = validParameter.validFile(parameters, "reference", true);
266 if (templateFileName == "not found") {
267 //check for saved reference sequences
268 if (rdb->referenceSeqs.size() != 0) {
269 templateFileName = "saved";
271 m->mothurOut("[ERROR]: You don't have any saved reference sequences and the reference parameter is a required for the align.seqs command.");
272 m->mothurOutEndLine();
275 }else if (templateFileName == "not open") { abort = true; }
276 else { if (save) { rdb->setSavedReference(templateFileName); } }
278 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
279 convert(temp, threshold);
281 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
282 if ((search != "suffix") && (search != "kmer") && (search != "blast")) { m->mothurOut("invalid search option: choices are kmer, suffix or blast."); m->mothurOutEndLine(); abort=true; }
284 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
285 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; }
290 catch(exception& e) {
291 m->errorOut(e, "AlignCommand", "AlignCommand");
295 //**********************************************************************************************************************
296 AlignCommand::~AlignCommand(){
298 if (abort == false) {
299 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
304 //**********************************************************************************************************************
306 int AlignCommand::execute(){
308 if (abort == true) { if (calledHelp) { return 0; } return 2; }
310 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
311 int longestBase = templateDB->getLongestBase();
313 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
314 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
315 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
316 else if(align == "noalign") { alignment = new NoAlign(); }
318 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
319 m->mothurOutEndLine();
320 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
323 for (int s = 0; s < candidateFileNames.size(); s++) {
324 if (m->control_pressed) { outputTypes.clear(); return 0; }
326 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
328 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
329 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
330 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
331 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
332 bool hasAccnos = true;
334 int numFastaSeqs = 0;
335 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
336 int start = time(NULL);
339 int pid, numSeqsPerProcessor;
341 vector<unsigned long int> MPIPos;
342 MPIWroteAccnos = false;
345 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
346 MPI_Comm_size(MPI_COMM_WORLD, &processors);
349 MPI_File outMPIAlign;
350 MPI_File outMPIReport;
351 MPI_File outMPIAccnos;
353 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
354 int inMode=MPI_MODE_RDONLY;
356 char outAlignFilename[1024];
357 strcpy(outAlignFilename, alignFileName.c_str());
359 char outReportFilename[1024];
360 strcpy(outReportFilename, reportFileName.c_str());
362 char outAccnosFilename[1024];
363 strcpy(outAccnosFilename, accnosFileName.c_str());
365 char inFileName[1024];
366 strcpy(inFileName, candidateFileNames[s].c_str());
368 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
369 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
370 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
371 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
373 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
375 if (pid == 0) { //you are the root process
377 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
379 //send file positions to all processes
380 for(int i = 1; i < processors; i++) {
381 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
382 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
385 //figure out how many sequences you have to align
386 numSeqsPerProcessor = numFastaSeqs / processors;
387 int startIndex = pid * numSeqsPerProcessor;
388 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
391 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
393 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
395 for (int i = 1; i < processors; i++) {
397 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
398 if (tempResult != 0) { MPIWroteAccnos = true; }
400 }else{ //you are a child process
401 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
402 MPIPos.resize(numFastaSeqs+1);
403 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
406 //figure out how many sequences you have to align
407 numSeqsPerProcessor = numFastaSeqs / processors;
408 int startIndex = pid * numSeqsPerProcessor;
409 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
413 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
415 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
417 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
421 MPI_File_close(&inMPI);
422 MPI_File_close(&outMPIAlign);
423 MPI_File_close(&outMPIReport);
424 MPI_File_close(&outMPIAccnos);
426 //delete accnos file if blank
428 //delete accnos file if its blank else report to user
429 if (MPIWroteAccnos) {
430 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
432 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
433 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
434 m->mothurOutEndLine();
437 //MPI_File_delete(outAccnosFilename, info);
439 remove(accnosFileName.c_str());
445 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
446 for (int i = 0; i < (positions.size()-1); i++) {
447 lines.push_back(new linePair(positions[i], positions[(i+1)]));
449 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
451 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
453 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
456 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
459 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
461 //delete accnos file if its blank else report to user
462 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
464 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
466 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
467 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
468 m->mothurOutEndLine();
475 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
477 if (pid == 0) { //only one process should output to screen
480 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
481 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
482 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
488 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
489 m->mothurOutEndLine();
490 m->mothurOutEndLine();
493 //set align file as new current fastafile
494 string currentFasta = "";
495 itTypes = outputTypes.find("fasta");
496 if (itTypes != outputTypes.end()) {
497 if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; m->setFastaFile(currentFasta); }
500 m->mothurOutEndLine();
501 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
502 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
503 m->mothurOutEndLine();
507 catch(exception& e) {
508 m->errorOut(e, "AlignCommand", "execute");
513 //**********************************************************************************************************************
515 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
517 ofstream alignmentFile;
518 m->openOutputFile(alignFName, alignmentFile);
521 m->openOutputFile(accnosFName, accnosFile);
523 NastReport report(reportFName);
526 m->openInputFile(filename, inFASTA);
528 inFASTA.seekg(filePos->start);
535 if (m->control_pressed) { return 0; }
537 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
538 report.setCandidate(candidateSeq);
540 int origNumBases = candidateSeq->getNumBases();
541 string originalUnaligned = candidateSeq->getUnaligned();
542 int numBasesNeeded = origNumBases * threshold;
544 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
545 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
546 alignment->resize(candidateSeq->getUnaligned().length()+1);
549 Sequence temp = templateDB->findClosestSequence(candidateSeq);
550 Sequence* templateSeq = &temp;
552 float searchScore = templateDB->getSearchScore();
554 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
559 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
560 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
561 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
562 //so this bool tells you if you need to delete it
564 //if there is a possibility that this sequence should be reversed
565 if (candidateSeq->getNumBases() < numBasesNeeded) {
567 string wasBetter = "";
568 //if the user wants you to try the reverse
571 //get reverse compliment
572 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
573 copy->reverseComplement();
576 Sequence temp2 = templateDB->findClosestSequence(copy);
577 Sequence* templateSeq2 = &temp2;
579 searchScore = templateDB->getSearchScore();
581 nast2 = new Nast(alignment, copy, templateSeq2);
583 //check if any better
584 if (copy->getNumBases() > candidateSeq->getNumBases()) {
585 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
586 templateSeq = templateSeq2;
589 needToDeleteCopy = true;
590 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
592 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
598 //create accnos file with names
599 accnosFile << candidateSeq->getName() << wasBetter << endl;
602 report.setTemplate(templateSeq);
603 report.setSearchParameters(search, searchScore);
604 report.setAlignmentParameters(align, alignment);
605 report.setNastParameters(*nast);
607 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
611 if (needToDeleteCopy) { delete copy; }
617 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
618 unsigned long int pos = inFASTA.tellg();
619 if ((pos == -1) || (pos >= filePos->end)) { break; }
621 if (inFASTA.eof()) { break; }
625 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
629 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
631 alignmentFile.close();
637 catch(exception& e) {
638 m->errorOut(e, "AlignCommand", "driver");
642 //**********************************************************************************************************************
644 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
646 string outputString = "";
647 MPI_Status statusReport;
648 MPI_Status statusAlign;
649 MPI_Status statusAccnos;
652 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
657 outputString = report.getHeaders();
658 int length = outputString.length();
660 char* buf = new char[length];
661 memcpy(buf, outputString.c_str(), length);
663 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
668 for(int i=0;i<num;i++){
670 if (m->control_pressed) { return 0; }
673 int length = MPIPos[start+i+1] - MPIPos[start+i];
675 char* buf4 = new char[length];
676 //memcpy(buf4, outputString.c_str(), length);
678 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
680 string tempBuf = buf4;
684 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
686 istringstream iss (tempBuf,istringstream::in);
688 Sequence* candidateSeq = new Sequence(iss);
689 report.setCandidate(candidateSeq);
691 int origNumBases = candidateSeq->getNumBases();
692 string originalUnaligned = candidateSeq->getUnaligned();
693 int numBasesNeeded = origNumBases * threshold;
695 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
696 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
697 alignment->resize(candidateSeq->getUnaligned().length()+1);
700 Sequence temp = templateDB->findClosestSequence(candidateSeq);
701 Sequence* templateSeq = &temp;
703 float searchScore = templateDB->getSearchScore();
705 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
709 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
710 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
711 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
712 //so this bool tells you if you need to delete it
714 //if there is a possibility that this sequence should be reversed
715 if (candidateSeq->getNumBases() < numBasesNeeded) {
717 string wasBetter = "";
718 //if the user wants you to try the reverse
720 //get reverse compliment
721 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
722 copy->reverseComplement();
725 Sequence temp2 = templateDB->findClosestSequence(copy);
726 Sequence* templateSeq2 = &temp2;
728 searchScore = templateDB->getSearchScore();
730 nast2 = new Nast(alignment, copy, templateSeq2);
732 //check if any better
733 if (copy->getNumBases() > candidateSeq->getNumBases()) {
734 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
735 templateSeq = templateSeq2;
738 needToDeleteCopy = true;
739 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
741 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
747 //create accnos file with names
748 outputString = candidateSeq->getName() + wasBetter + "\n";
750 //send results to parent
751 int length = outputString.length();
753 char* buf = new char[length];
754 memcpy(buf, outputString.c_str(), length);
756 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
758 MPIWroteAccnos = true;
761 report.setTemplate(templateSeq);
762 report.setSearchParameters(search, searchScore);
763 report.setAlignmentParameters(align, alignment);
764 report.setNastParameters(*nast);
766 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
768 //send results to parent
769 int length = outputString.length();
770 char* buf2 = new char[length];
771 memcpy(buf2, outputString.c_str(), length);
773 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
777 outputString = report.getReport();
779 //send results to parent
780 length = outputString.length();
781 char* buf3 = new char[length];
782 memcpy(buf3, outputString.c_str(), length);
784 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
788 if (needToDeleteCopy) { delete copy; }
793 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
796 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
800 catch(exception& e) {
801 m->errorOut(e, "AlignCommand", "driverMPI");
806 /**************************************************************************************************/
808 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
810 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
811 processIDS.resize(0);
814 // processIDS.resize(0);
816 //loop through and create all the processes you want
817 while (process != processors) {
821 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
824 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
826 //pass numSeqs to parent
828 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
829 m->openOutputFile(tempFile, out);
835 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
836 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
842 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
844 //force parent to wait until all the processes are done
845 for (int i=0;i<processIDS.size();i++) {
846 int temp = processIDS[i];
850 vector<string> nonBlankAccnosFiles;
851 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
852 else { remove(accnosFName.c_str()); } //remove so other files can be renamed to it
854 for (int i = 0; i < processIDS.size(); i++) {
856 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
857 m->openInputFile(tempFile, in);
858 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
859 in.close(); remove(tempFile.c_str());
861 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
862 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
864 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
865 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
867 if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
868 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
869 }else { remove((accnosFName + toString(processIDS[i]) + ".temp").c_str()); }
873 //append accnos files
874 if (nonBlankAccnosFiles.size() != 0) {
875 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
877 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
878 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
879 remove(nonBlankAccnosFiles[h].c_str());
881 }else { //recreate the accnosfile if needed
883 m->openOutputFile(accnosFName, out);
890 catch(exception& e) {
891 m->errorOut(e, "AlignCommand", "createProcesses");
895 /**************************************************************************************************/
897 void AlignCommand::appendAlignFiles(string temp, string filename) {
902 m->openOutputFileAppend(filename, output);
903 m->openInputFile(temp, input);
905 while(char c = input.get()){
906 if(input.eof()) { break; }
907 else { output << c; }
913 catch(exception& e) {
914 m->errorOut(e, "AlignCommand", "appendAlignFiles");
918 //**********************************************************************************************************************
920 void AlignCommand::appendReportFiles(string temp, string filename) {
925 m->openOutputFileAppend(filename, output);
926 m->openInputFile(temp, input);
928 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
930 while(char c = input.get()){
931 if(input.eof()) { break; }
932 else { output << c; }
938 catch(exception& e) {
939 m->errorOut(e, "AlignCommand", "appendReportFiles");
943 //**********************************************************************************************************************