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(){
67 //initialize outputTypes
68 vector<string> tempOutNames;
69 outputTypes["fasta"] = tempOutNames;
70 outputTypes["alignreport"] = tempOutNames;
71 outputTypes["accnos"] = tempOutNames;
74 m->errorOut(e, "AlignCommand", "AlignCommand");
78 //**********************************************************************************************************************
79 AlignCommand::AlignCommand(string option) {
83 //allow user to run help
84 if(option == "help") { help(); abort = true; }
88 //valid paramters for this command
89 string AlignArray[] = {"template","candidate","search","ksize","align","match","mismatch","gapopen","gapextend", "processors","flip","threshold","outputdir","inputdir"};
90 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
92 OptionParser parser(option);
93 map<string, string> parameters = parser.getParameters();
95 ValidParameters validParameter("align.seqs");
96 map<string, string>::iterator it;
98 //check to make sure all parameters are valid for command
99 for (it = parameters.begin(); it != parameters.end(); it++) {
100 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
103 //initialize outputTypes
104 vector<string> tempOutNames;
105 outputTypes["fasta"] = tempOutNames;
106 outputTypes["alignreport"] = tempOutNames;
107 outputTypes["accnos"] = tempOutNames;
109 //if the user changes the output directory command factory will send this info to us in the output parameter
110 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
113 //if the user changes the input directory command factory will send this info to us in the output parameter
114 string inputDir = validParameter.validFile(parameters, "inputdir", false);
116 if (inputDir == "not found"){ inputDir = ""; }
120 it = parameters.find("template");
122 //user has given a template file
123 if(it != parameters.end()){
124 path = m->hasPath(it->second);
125 //if the user has not given a path then, add inputdir. else leave path alone.
126 if (path == "") { parameters["template"] = inputDir + it->second; }
130 //check for required parameters
131 templateFileName = validParameter.validFile(parameters, "template", true);
133 if (templateFileName == "not found") {
134 m->mothurOut("template is a required parameter for the align.seqs command.");
135 m->mothurOutEndLine();
137 }else if (templateFileName == "not open") { abort = true; }
139 candidateFileName = validParameter.validFile(parameters, "candidate", false);
140 if (candidateFileName == "not found") { m->mothurOut("candidate is a required parameter for the align.seqs command."); m->mothurOutEndLine(); abort = true; }
142 m->splitAtDash(candidateFileName, candidateFileNames);
144 //go through files and make sure they are good, if not, then disregard them
145 for (int i = 0; i < candidateFileNames.size(); i++) {
146 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
148 if (inputDir != "") {
149 string path = m->hasPath(candidateFileNames[i]);
150 //if the user has not given a path then, add inputdir. else leave path alone.
151 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
157 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
159 //if you can't open it, try default location
160 if (ableToOpen == 1) {
161 if (m->getDefaultPath() != "") { //default path is set
162 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
163 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
164 ableToOpen = m->openInputFile(tryPath, in, "noerror");
165 candidateFileNames[i] = tryPath;
169 //if you can't open it, try default location
170 if (ableToOpen == 1) {
171 if (m->getOutputDir() != "") { //default path is set
172 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
173 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
174 ableToOpen = m->openInputFile(tryPath, in, "noerror");
175 candidateFileNames[i] = tryPath;
181 if (ableToOpen == 1) {
182 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
183 //erase from file list
184 candidateFileNames.erase(candidateFileNames.begin()+i);
190 //make sure there is at least one valid file left
191 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
194 //check for optional parameter and set defaults
195 // ...at some point should added some additional type checking...
197 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
198 convert(temp, kmerSize);
200 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
201 convert(temp, match);
203 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
204 convert(temp, misMatch);
206 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
207 convert(temp, gapOpen);
209 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
210 convert(temp, gapExtend);
212 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
213 convert(temp, processors);
215 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
216 flip = m->isTrue(temp);
218 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
219 convert(temp, threshold);
221 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
223 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
227 catch(exception& e) {
228 m->errorOut(e, "AlignCommand", "AlignCommand");
232 //**********************************************************************************************************************
234 AlignCommand::~AlignCommand(){
236 if (abort == false) {
237 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
243 //**********************************************************************************************************************
245 void AlignCommand::help(){
247 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
248 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen, gapextend and processors.\n");
249 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");
250 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");
251 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");
252 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
253 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
254 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
255 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
256 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
257 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");
258 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");
259 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");
260 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");
261 m->mothurOut("The align.seqs command should be in the following format: \n");
262 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
263 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");
264 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
266 catch(exception& e) {
267 m->errorOut(e, "AlignCommand", "help");
273 //**********************************************************************************************************************
275 int AlignCommand::execute(){
277 if (abort == true) { return 0; }
279 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
280 int longestBase = templateDB->getLongestBase();
282 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
283 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
284 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
285 else if(align == "noalign") { alignment = new NoAlign(); }
287 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
288 m->mothurOutEndLine();
289 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
292 for (int s = 0; s < candidateFileNames.size(); s++) {
293 if (m->control_pressed) { outputTypes.clear(); return 0; }
295 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
297 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
298 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
299 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
300 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
301 bool hasAccnos = true;
303 int numFastaSeqs = 0;
304 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
305 int start = time(NULL);
308 int pid, end, numSeqsPerProcessor;
310 vector<unsigned long int> MPIPos;
311 MPIWroteAccnos = false;
314 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
315 MPI_Comm_size(MPI_COMM_WORLD, &processors);
318 MPI_File outMPIAlign;
319 MPI_File outMPIReport;
320 MPI_File outMPIAccnos;
322 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
323 int inMode=MPI_MODE_RDONLY;
325 char outAlignFilename[1024];
326 strcpy(outAlignFilename, alignFileName.c_str());
328 char outReportFilename[1024];
329 strcpy(outReportFilename, reportFileName.c_str());
331 char outAccnosFilename[1024];
332 strcpy(outAccnosFilename, accnosFileName.c_str());
334 char inFileName[1024];
335 strcpy(inFileName, candidateFileNames[s].c_str());
337 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
338 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
339 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
340 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
342 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
344 if (pid == 0) { //you are the root process
346 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
348 //send file positions to all processes
349 for(int i = 1; i < processors; i++) {
350 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
351 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
354 //figure out how many sequences you have to align
355 numSeqsPerProcessor = numFastaSeqs / processors;
356 int startIndex = pid * numSeqsPerProcessor;
357 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
360 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
362 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
364 for (int i = 1; i < processors; i++) {
366 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
367 if (tempResult != 0) { MPIWroteAccnos = true; }
369 }else{ //you are a child process
370 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
371 MPIPos.resize(numFastaSeqs+1);
372 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
375 //figure out how many sequences you have to align
376 numSeqsPerProcessor = numFastaSeqs / processors;
377 int startIndex = pid * numSeqsPerProcessor;
378 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
382 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
384 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
386 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
390 MPI_File_close(&inMPI);
391 MPI_File_close(&outMPIAlign);
392 MPI_File_close(&outMPIReport);
393 MPI_File_close(&outMPIAccnos);
395 //delete accnos file if blank
397 //delete accnos file if its blank else report to user
398 if (MPIWroteAccnos) {
399 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
401 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
402 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
403 m->mothurOutEndLine();
406 //MPI_File_delete(outAccnosFilename, info);
408 remove(accnosFileName.c_str());
414 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
415 for (int i = 0; i < (positions.size()-1); i++) {
416 lines.push_back(new linePair(positions[i], positions[(i+1)]));
418 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
420 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
422 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
424 //delete accnos file if its blank else report to user
425 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
427 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
429 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
430 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
431 m->mothurOutEndLine();
434 processIDS.resize(0);
436 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
438 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
439 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
441 //append alignment and report files
442 for(int i=1;i<processors;i++){
443 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
444 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
446 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
447 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
450 vector<string> nonBlankAccnosFiles;
451 //delete blank accnos files generated with multiple processes
452 for(int i=0;i<processors;i++){
453 if (!(m->isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
454 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
455 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
458 //append accnos files
459 if (nonBlankAccnosFiles.size() != 0) {
460 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
462 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
463 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
464 remove(nonBlankAccnosFiles[h].c_str());
466 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
468 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
469 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
470 m->mothurOutEndLine();
471 }else{ hasAccnos = false; }
473 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
476 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
478 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
480 //delete accnos file if its blank else report to user
481 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
483 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
485 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
486 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
487 m->mothurOutEndLine();
496 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
498 if (pid == 0) { //only one process should output to screen
501 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
502 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
503 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
509 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
510 m->mothurOutEndLine();
511 m->mothurOutEndLine();
515 m->mothurOutEndLine();
516 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
517 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
518 m->mothurOutEndLine();
522 catch(exception& e) {
523 m->errorOut(e, "AlignCommand", "execute");
528 //**********************************************************************************************************************
530 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
532 ofstream alignmentFile;
533 m->openOutputFile(alignFName, alignmentFile);
536 m->openOutputFile(accnosFName, accnosFile);
538 NastReport report(reportFName);
541 m->openInputFile(filename, inFASTA);
543 inFASTA.seekg(filePos->start);
550 if (m->control_pressed) { return 0; }
552 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
554 int origNumBases = candidateSeq->getNumBases();
555 string originalUnaligned = candidateSeq->getUnaligned();
556 int numBasesNeeded = origNumBases * threshold;
558 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
559 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
560 alignment->resize(candidateSeq->getUnaligned().length()+1);
563 Sequence temp = templateDB->findClosestSequence(candidateSeq);
564 Sequence* templateSeq = &temp;
566 float searchScore = templateDB->getSearchScore();
568 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
573 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
574 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
575 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
576 //so this bool tells you if you need to delete it
578 //if there is a possibility that this sequence should be reversed
579 if (candidateSeq->getNumBases() < numBasesNeeded) {
581 string wasBetter = "";
582 //if the user wants you to try the reverse
585 //get reverse compliment
586 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
587 copy->reverseComplement();
590 Sequence temp2 = templateDB->findClosestSequence(copy);
591 Sequence* templateSeq2 = &temp2;
593 searchScore = templateDB->getSearchScore();
595 nast2 = new Nast(alignment, copy, templateSeq2);
597 //check if any better
598 if (copy->getNumBases() > candidateSeq->getNumBases()) {
599 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
600 templateSeq = templateSeq2;
603 needToDeleteCopy = true;
604 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
606 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
612 //create accnos file with names
613 accnosFile << candidateSeq->getName() << wasBetter << endl;
616 report.setCandidate(candidateSeq);
617 report.setTemplate(templateSeq);
618 report.setSearchParameters(search, searchScore);
619 report.setAlignmentParameters(align, alignment);
620 report.setNastParameters(*nast);
622 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
626 if (needToDeleteCopy) { delete copy; }
632 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
633 unsigned long int pos = inFASTA.tellg();
634 if ((pos == -1) || (pos >= filePos->end)) { break; }
636 if (inFASTA.eof()) { break; }
640 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
644 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
646 alignmentFile.close();
652 catch(exception& e) {
653 m->errorOut(e, "AlignCommand", "driver");
657 //**********************************************************************************************************************
659 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
661 string outputString = "";
662 MPI_Status statusReport;
663 MPI_Status statusAlign;
664 MPI_Status statusAccnos;
667 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
672 outputString = report.getHeaders();
673 int length = outputString.length();
675 char* buf = new char[length];
676 memcpy(buf, outputString.c_str(), length);
678 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
683 for(int i=0;i<num;i++){
685 if (m->control_pressed) { return 0; }
688 int length = MPIPos[start+i+1] - MPIPos[start+i];
690 char* buf4 = new char[length];
691 //memcpy(buf4, outputString.c_str(), length);
693 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
695 string tempBuf = buf4;
699 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
701 istringstream iss (tempBuf,istringstream::in);
703 Sequence* candidateSeq = new Sequence(iss);
705 int origNumBases = candidateSeq->getNumBases();
706 string originalUnaligned = candidateSeq->getUnaligned();
707 int numBasesNeeded = origNumBases * threshold;
709 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
710 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
711 alignment->resize(candidateSeq->getUnaligned().length()+1);
714 Sequence temp = templateDB->findClosestSequence(candidateSeq);
715 Sequence* templateSeq = &temp;
717 float searchScore = templateDB->getSearchScore();
719 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
723 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
724 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
725 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
726 //so this bool tells you if you need to delete it
728 //if there is a possibility that this sequence should be reversed
729 if (candidateSeq->getNumBases() < numBasesNeeded) {
731 string wasBetter = "";
732 //if the user wants you to try the reverse
734 //get reverse compliment
735 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
736 copy->reverseComplement();
739 Sequence temp2 = templateDB->findClosestSequence(copy);
740 Sequence* templateSeq2 = &temp2;
742 searchScore = templateDB->getSearchScore();
744 nast2 = new Nast(alignment, copy, templateSeq2);
746 //check if any better
747 if (copy->getNumBases() > candidateSeq->getNumBases()) {
748 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
749 templateSeq = templateSeq2;
752 needToDeleteCopy = true;
754 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
760 //create accnos file with names
761 outputString = candidateSeq->getName() + wasBetter + "\n";
763 //send results to parent
764 int length = outputString.length();
766 char* buf = new char[length];
767 memcpy(buf, outputString.c_str(), length);
769 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
771 MPIWroteAccnos = true;
774 report.setCandidate(candidateSeq);
775 report.setTemplate(templateSeq);
776 report.setSearchParameters(search, searchScore);
777 report.setAlignmentParameters(align, alignment);
778 report.setNastParameters(*nast);
780 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
782 //send results to parent
783 int length = outputString.length();
784 char* buf2 = new char[length];
785 memcpy(buf2, outputString.c_str(), length);
787 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
791 outputString = report.getReport();
793 //send results to parent
794 length = outputString.length();
795 char* buf3 = new char[length];
796 memcpy(buf3, outputString.c_str(), length);
798 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
802 if (needToDeleteCopy) { delete copy; }
807 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
810 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
814 catch(exception& e) {
815 m->errorOut(e, "AlignCommand", "driverMPI");
820 /**************************************************************************************************/
822 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
824 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
827 // processIDS.resize(0);
829 //loop through and create all the processes you want
830 while (process != processors) {
834 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
837 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
839 //pass numSeqs to parent
841 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
842 m->openOutputFile(tempFile, out);
847 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
850 //force parent to wait until all the processes are done
851 for (int i=0;i<processors;i++) {
852 int temp = processIDS[i];
856 for (int i = 0; i < processIDS.size(); i++) {
858 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
859 m->openInputFile(tempFile, in);
860 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
861 in.close(); remove(tempFile.c_str());
867 catch(exception& e) {
868 m->errorOut(e, "AlignCommand", "createProcesses");
872 /**************************************************************************************************/
874 void AlignCommand::appendAlignFiles(string temp, string filename) {
879 m->openOutputFileAppend(filename, output);
880 m->openInputFile(temp, input);
882 while(char c = input.get()){
883 if(input.eof()) { break; }
884 else { output << c; }
890 catch(exception& e) {
891 m->errorOut(e, "AlignCommand", "appendAlignFiles");
895 //**********************************************************************************************************************
897 void AlignCommand::appendReportFiles(string temp, string filename) {
902 m->openOutputFileAppend(filename, output);
903 m->openInputFile(temp, input);
905 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
907 while(char c = input.get()){
908 if(input.eof()) { break; }
909 else { output << c; }
915 catch(exception& e) {
916 m->errorOut(e, "AlignCommand", "appendReportFiles");
920 //**********************************************************************************************************************