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 //**********************************************************************************************************************
30 vector<string> AlignCommand::getValidParameters(){
32 string AlignArray[] = {"template","candidate","search","ksize","align","match","mismatch","gapopen","gapextend", "processors","flip","threshold","outputdir","inputdir"};
33 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
37 m->errorOut(e, "AlignCommand", "getValidParameters");
41 //**********************************************************************************************************************
42 vector<string> AlignCommand::getRequiredParameters(){
44 string AlignArray[] = {"template","candidate"};
45 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
49 m->errorOut(e, "AlignCommand", "getRequiredParameters");
53 //**********************************************************************************************************************
54 vector<string> AlignCommand::getRequiredFiles(){
56 vector<string> myArray;
60 m->errorOut(e, "AlignCommand", "getRequiredFiles");
64 //**********************************************************************************************************************
65 AlignCommand::AlignCommand(){
68 //initialize outputTypes
69 vector<string> tempOutNames;
70 outputTypes["fasta"] = tempOutNames;
71 outputTypes["alignreport"] = tempOutNames;
72 outputTypes["accnos"] = tempOutNames;
75 m->errorOut(e, "AlignCommand", "AlignCommand");
79 //**********************************************************************************************************************
80 AlignCommand::AlignCommand(string option) {
84 //allow user to run help
85 if(option == "help") { help(); abort = true; }
89 //valid paramters for this command
90 string AlignArray[] = {"template","candidate","search","ksize","align","match","mismatch","gapopen","gapextend", "processors","flip","threshold","outputdir","inputdir"};
91 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
93 OptionParser parser(option);
94 map<string, string> parameters = parser.getParameters();
96 ValidParameters validParameter("align.seqs");
97 map<string, string>::iterator it;
99 //check to make sure all parameters are valid for command
100 for (it = parameters.begin(); it != parameters.end(); it++) {
101 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
104 //initialize outputTypes
105 vector<string> tempOutNames;
106 outputTypes["fasta"] = tempOutNames;
107 outputTypes["alignreport"] = tempOutNames;
108 outputTypes["accnos"] = tempOutNames;
110 //if the user changes the output directory command factory will send this info to us in the output parameter
111 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
114 //if the user changes the input directory command factory will send this info to us in the output parameter
115 string inputDir = validParameter.validFile(parameters, "inputdir", false);
117 if (inputDir == "not found"){ inputDir = ""; }
121 it = parameters.find("template");
123 //user has given a template file
124 if(it != parameters.end()){
125 path = m->hasPath(it->second);
126 //if the user has not given a path then, add inputdir. else leave path alone.
127 if (path == "") { parameters["template"] = inputDir + it->second; }
131 //check for required parameters
132 templateFileName = validParameter.validFile(parameters, "template", true);
134 if (templateFileName == "not found") {
135 m->mothurOut("template is a required parameter for the align.seqs command.");
136 m->mothurOutEndLine();
138 }else if (templateFileName == "not open") { abort = true; }
140 candidateFileName = validParameter.validFile(parameters, "candidate", false);
141 if (candidateFileName == "not found") { m->mothurOut("candidate is a required parameter for the align.seqs command."); m->mothurOutEndLine(); abort = true; }
143 m->splitAtDash(candidateFileName, candidateFileNames);
145 //go through files and make sure they are good, if not, then disregard them
146 for (int i = 0; i < candidateFileNames.size(); i++) {
147 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
149 if (inputDir != "") {
150 string path = m->hasPath(candidateFileNames[i]);
151 //if the user has not given a path then, add inputdir. else leave path alone.
152 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
157 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
160 //if you can't open it, try default location
161 if (ableToOpen == 1) {
162 if (m->getDefaultPath() != "") { //default path is set
163 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
164 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
166 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
168 candidateFileNames[i] = tryPath;
172 //if you can't open it, try default location
173 if (ableToOpen == 1) {
174 if (m->getOutputDir() != "") { //default path is set
175 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
176 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
178 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
180 candidateFileNames[i] = tryPath;
186 if (ableToOpen == 1) {
187 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
188 //erase from file list
189 candidateFileNames.erase(candidateFileNames.begin()+i);
195 //make sure there is at least one valid file left
196 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
199 //check for optional parameter and set defaults
200 // ...at some point should added some additional type checking...
202 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
203 convert(temp, kmerSize);
205 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
206 convert(temp, match);
208 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
209 convert(temp, misMatch);
211 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
212 convert(temp, gapOpen);
214 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
215 convert(temp, gapExtend);
217 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
218 convert(temp, processors);
220 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
221 flip = m->isTrue(temp);
223 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
224 convert(temp, threshold);
226 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
228 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
232 catch(exception& e) {
233 m->errorOut(e, "AlignCommand", "AlignCommand");
237 //**********************************************************************************************************************
239 AlignCommand::~AlignCommand(){
241 if (abort == false) {
242 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
248 //**********************************************************************************************************************
250 void AlignCommand::help(){
252 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
253 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen, gapextend and processors.\n");
254 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");
255 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");
256 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");
257 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
258 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
259 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
260 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
261 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
262 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");
263 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");
264 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");
265 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");
266 m->mothurOut("The align.seqs command should be in the following format: \n");
267 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
268 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");
269 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
271 catch(exception& e) {
272 m->errorOut(e, "AlignCommand", "help");
278 //**********************************************************************************************************************
280 int AlignCommand::execute(){
282 if (abort == true) { return 0; }
284 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
285 int longestBase = templateDB->getLongestBase();
287 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
288 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
289 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
290 else if(align == "noalign") { alignment = new NoAlign(); }
292 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
293 m->mothurOutEndLine();
294 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
297 for (int s = 0; s < candidateFileNames.size(); s++) {
298 if (m->control_pressed) { outputTypes.clear(); return 0; }
300 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
302 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
303 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
304 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
305 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
306 bool hasAccnos = true;
308 int numFastaSeqs = 0;
309 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
310 int start = time(NULL);
313 int pid, end, numSeqsPerProcessor;
315 vector<unsigned long int> MPIPos;
316 MPIWroteAccnos = false;
319 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
320 MPI_Comm_size(MPI_COMM_WORLD, &processors);
323 MPI_File outMPIAlign;
324 MPI_File outMPIReport;
325 MPI_File outMPIAccnos;
327 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
328 int inMode=MPI_MODE_RDONLY;
330 char outAlignFilename[1024];
331 strcpy(outAlignFilename, alignFileName.c_str());
333 char outReportFilename[1024];
334 strcpy(outReportFilename, reportFileName.c_str());
336 char outAccnosFilename[1024];
337 strcpy(outAccnosFilename, accnosFileName.c_str());
339 char inFileName[1024];
340 strcpy(inFileName, candidateFileNames[s].c_str());
342 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
343 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
344 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
345 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
347 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
349 if (pid == 0) { //you are the root process
351 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
353 //send file positions to all processes
354 for(int i = 1; i < processors; i++) {
355 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
356 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
359 //figure out how many sequences you have to align
360 numSeqsPerProcessor = numFastaSeqs / processors;
361 int startIndex = pid * numSeqsPerProcessor;
362 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
365 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
367 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
369 for (int i = 1; i < processors; i++) {
371 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
372 if (tempResult != 0) { MPIWroteAccnos = true; }
374 }else{ //you are a child process
375 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
376 MPIPos.resize(numFastaSeqs+1);
377 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
380 //figure out how many sequences you have to align
381 numSeqsPerProcessor = numFastaSeqs / processors;
382 int startIndex = pid * numSeqsPerProcessor;
383 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
387 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
389 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
391 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
395 MPI_File_close(&inMPI);
396 MPI_File_close(&outMPIAlign);
397 MPI_File_close(&outMPIReport);
398 MPI_File_close(&outMPIAccnos);
400 //delete accnos file if blank
402 //delete accnos file if its blank else report to user
403 if (MPIWroteAccnos) {
404 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
406 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
407 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
408 m->mothurOutEndLine();
411 //MPI_File_delete(outAccnosFilename, info);
413 remove(accnosFileName.c_str());
419 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
420 for (int i = 0; i < (positions.size()-1); i++) {
421 lines.push_back(new linePair(positions[i], positions[(i+1)]));
423 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
425 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
427 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
429 //delete accnos file if its blank else report to user
430 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
432 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
434 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
435 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
436 m->mothurOutEndLine();
439 processIDS.resize(0);
441 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
443 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
444 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
446 //append alignment and report files
447 for(int i=1;i<processors;i++){
448 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
449 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
451 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
452 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
455 vector<string> nonBlankAccnosFiles;
456 //delete blank accnos files generated with multiple processes
457 for(int i=0;i<processors;i++){
458 if (!(m->isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
459 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
460 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
463 //append accnos files
464 if (nonBlankAccnosFiles.size() != 0) {
465 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
467 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
468 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
469 remove(nonBlankAccnosFiles[h].c_str());
471 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
473 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
474 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
475 m->mothurOutEndLine();
476 }else{ hasAccnos = false; }
478 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
481 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
483 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
485 //delete accnos file if its blank else report to user
486 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
488 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
490 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
491 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
492 m->mothurOutEndLine();
501 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
503 if (pid == 0) { //only one process should output to screen
506 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
507 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
508 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
514 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
515 m->mothurOutEndLine();
516 m->mothurOutEndLine();
520 m->mothurOutEndLine();
521 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
522 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
523 m->mothurOutEndLine();
527 catch(exception& e) {
528 m->errorOut(e, "AlignCommand", "execute");
533 //**********************************************************************************************************************
535 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
537 ofstream alignmentFile;
538 m->openOutputFile(alignFName, alignmentFile);
541 m->openOutputFile(accnosFName, accnosFile);
543 NastReport report(reportFName);
546 m->openInputFile(filename, inFASTA);
548 inFASTA.seekg(filePos->start);
555 if (m->control_pressed) { return 0; }
557 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
559 int origNumBases = candidateSeq->getNumBases();
560 string originalUnaligned = candidateSeq->getUnaligned();
561 int numBasesNeeded = origNumBases * threshold;
563 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
564 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
565 alignment->resize(candidateSeq->getUnaligned().length()+1);
568 Sequence temp = templateDB->findClosestSequence(candidateSeq);
569 Sequence* templateSeq = &temp;
571 float searchScore = templateDB->getSearchScore();
573 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
578 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
579 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
580 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
581 //so this bool tells you if you need to delete it
583 //if there is a possibility that this sequence should be reversed
584 if (candidateSeq->getNumBases() < numBasesNeeded) {
586 string wasBetter = "";
587 //if the user wants you to try the reverse
590 //get reverse compliment
591 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
592 copy->reverseComplement();
595 Sequence temp2 = templateDB->findClosestSequence(copy);
596 Sequence* templateSeq2 = &temp2;
598 searchScore = templateDB->getSearchScore();
600 nast2 = new Nast(alignment, copy, templateSeq2);
602 //check if any better
603 if (copy->getNumBases() > candidateSeq->getNumBases()) {
604 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
605 templateSeq = templateSeq2;
608 needToDeleteCopy = true;
609 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
611 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
617 //create accnos file with names
618 accnosFile << candidateSeq->getName() << wasBetter << endl;
621 report.setCandidate(candidateSeq);
622 report.setTemplate(templateSeq);
623 report.setSearchParameters(search, searchScore);
624 report.setAlignmentParameters(align, alignment);
625 report.setNastParameters(*nast);
627 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
631 if (needToDeleteCopy) { delete copy; }
637 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
638 unsigned long int pos = inFASTA.tellg();
639 if ((pos == -1) || (pos >= filePos->end)) { break; }
641 if (inFASTA.eof()) { break; }
645 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
649 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
651 alignmentFile.close();
657 catch(exception& e) {
658 m->errorOut(e, "AlignCommand", "driver");
662 //**********************************************************************************************************************
664 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
666 string outputString = "";
667 MPI_Status statusReport;
668 MPI_Status statusAlign;
669 MPI_Status statusAccnos;
672 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
677 outputString = report.getHeaders();
678 int length = outputString.length();
680 char* buf = new char[length];
681 memcpy(buf, outputString.c_str(), length);
683 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
688 for(int i=0;i<num;i++){
690 if (m->control_pressed) { return 0; }
693 int length = MPIPos[start+i+1] - MPIPos[start+i];
695 char* buf4 = new char[length];
696 //memcpy(buf4, outputString.c_str(), length);
698 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
700 string tempBuf = buf4;
704 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
706 istringstream iss (tempBuf,istringstream::in);
708 Sequence* candidateSeq = new Sequence(iss);
710 int origNumBases = candidateSeq->getNumBases();
711 string originalUnaligned = candidateSeq->getUnaligned();
712 int numBasesNeeded = origNumBases * threshold;
714 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
715 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
716 alignment->resize(candidateSeq->getUnaligned().length()+1);
719 Sequence temp = templateDB->findClosestSequence(candidateSeq);
720 Sequence* templateSeq = &temp;
722 float searchScore = templateDB->getSearchScore();
724 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
728 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
729 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
730 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
731 //so this bool tells you if you need to delete it
733 //if there is a possibility that this sequence should be reversed
734 if (candidateSeq->getNumBases() < numBasesNeeded) {
736 string wasBetter = "";
737 //if the user wants you to try the reverse
739 //get reverse compliment
740 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
741 copy->reverseComplement();
744 Sequence temp2 = templateDB->findClosestSequence(copy);
745 Sequence* templateSeq2 = &temp2;
747 searchScore = templateDB->getSearchScore();
749 nast2 = new Nast(alignment, copy, templateSeq2);
751 //check if any better
752 if (copy->getNumBases() > candidateSeq->getNumBases()) {
753 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
754 templateSeq = templateSeq2;
757 needToDeleteCopy = true;
758 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
760 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
766 //create accnos file with names
767 outputString = candidateSeq->getName() + wasBetter + "\n";
769 //send results to parent
770 int length = outputString.length();
772 char* buf = new char[length];
773 memcpy(buf, outputString.c_str(), length);
775 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
777 MPIWroteAccnos = true;
780 report.setCandidate(candidateSeq);
781 report.setTemplate(templateSeq);
782 report.setSearchParameters(search, searchScore);
783 report.setAlignmentParameters(align, alignment);
784 report.setNastParameters(*nast);
786 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
788 //send results to parent
789 int length = outputString.length();
790 char* buf2 = new char[length];
791 memcpy(buf2, outputString.c_str(), length);
793 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
797 outputString = report.getReport();
799 //send results to parent
800 length = outputString.length();
801 char* buf3 = new char[length];
802 memcpy(buf3, outputString.c_str(), length);
804 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
808 if (needToDeleteCopy) { delete copy; }
813 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
816 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
820 catch(exception& e) {
821 m->errorOut(e, "AlignCommand", "driverMPI");
826 /**************************************************************************************************/
828 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
830 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
833 // processIDS.resize(0);
835 //loop through and create all the processes you want
836 while (process != processors) {
840 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
843 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
845 //pass numSeqs to parent
847 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
848 m->openOutputFile(tempFile, out);
853 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
856 //force parent to wait until all the processes are done
857 for (int i=0;i<processors;i++) {
858 int temp = processIDS[i];
862 for (int i = 0; i < processIDS.size(); i++) {
864 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
865 m->openInputFile(tempFile, in);
866 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
867 in.close(); remove(tempFile.c_str());
873 catch(exception& e) {
874 m->errorOut(e, "AlignCommand", "createProcesses");
878 /**************************************************************************************************/
880 void AlignCommand::appendAlignFiles(string temp, string filename) {
885 m->openOutputFileAppend(filename, output);
886 m->openInputFile(temp, input);
888 while(char c = input.get()){
889 if(input.eof()) { break; }
890 else { output << c; }
896 catch(exception& e) {
897 m->errorOut(e, "AlignCommand", "appendAlignFiles");
901 //**********************************************************************************************************************
903 void AlignCommand::appendReportFiles(string temp, string filename) {
908 m->openOutputFileAppend(filename, output);
909 m->openInputFile(temp, input);
911 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
913 while(char c = input.get()){
914 if(input.eof()) { break; }
915 else { output << c; }
921 catch(exception& e) {
922 m->errorOut(e, "AlignCommand", "appendReportFiles");
926 //**********************************************************************************************************************