5 * Created by Sarah Westcott on 5/15/09.
6 * Copyright 2009 Schloss Lab UMASS Amherst. All rights reserved.
8 * This version of nast does everything I think that the greengenes nast server does and then some. I have added the
9 * feature of allowing users to define their database, kmer size for searching, alignment penalty values and alignment
10 * method. This latter feature is perhaps most significant. nastPlus enables a user to use either a Needleman-Wunsch
11 * (non-affine gap penalty) or Gotoh (affine gap penalty) pairwise alignment algorithm. This is significant because it
12 * allows for a global alignment and not the local alignment provided by bLAst. Furthermore, it has the potential to
13 * provide a better alignment because of the banding method employed by blast (I'm not sure about this).
17 #include "aligncommand.h"
18 #include "sequence.hpp"
20 #include "gotohoverlap.hpp"
21 #include "needlemanoverlap.hpp"
22 #include "blastalign.hpp"
23 #include "noalign.hpp"
26 #include "nastreport.hpp"
29 //**********************************************************************************************************************
31 AlignCommand::AlignCommand(string option) {
36 //allow user to run help
37 if(option == "help") { help(); abort = true; }
41 //valid paramters for this command
42 string AlignArray[] = {"template","candidate","search","ksize","align","match","mismatch","gapopen","gapextend", "processors","flip","threshold","outputdir","inputdir"};
43 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
45 OptionParser parser(option);
46 map<string, string> parameters = parser.getParameters();
48 ValidParameters validParameter("align.seqs");
49 map<string, string>::iterator it;
51 //check to make sure all parameters are valid for command
52 for (it = parameters.begin(); it != parameters.end(); it++) {
53 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
56 //if the user changes the output directory command factory will send this info to us in the output parameter
57 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
60 //if the user changes the input directory command factory will send this info to us in the output parameter
61 string inputDir = validParameter.validFile(parameters, "inputdir", false);
63 if (inputDir == "not found"){ inputDir = ""; }
67 it = parameters.find("template");
69 //user has given a template file
70 if(it != parameters.end()){
71 path = hasPath(it->second);
72 //if the user has not given a path then, add inputdir. else leave path alone.
73 if (path == "") { parameters["template"] = inputDir + it->second; }
77 //check for required parameters
78 templateFileName = validParameter.validFile(parameters, "template", true);
80 if (templateFileName == "not found") {
81 m->mothurOut("template is a required parameter for the align.seqs command.");
82 m->mothurOutEndLine();
84 }else if (templateFileName == "not open") { abort = true; }
86 candidateFileName = validParameter.validFile(parameters, "candidate", false);
87 if (candidateFileName == "not found") { m->mothurOut("candidate is a required parameter for the align.seqs command."); m->mothurOutEndLine(); abort = true; }
89 splitAtDash(candidateFileName, candidateFileNames);
91 //go through files and make sure they are good, if not, then disregard them
92 for (int i = 0; i < candidateFileNames.size(); i++) {
94 string path = hasPath(candidateFileNames[i]);
95 //if the user has not given a path then, add inputdir. else leave path alone.
96 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
104 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
105 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
110 ableToOpen = openInputFile(candidateFileNames[i], in, "noerror");
112 //if you can't open it, try default location
113 if (ableToOpen == 1) {
114 if (m->getDefaultPath() != "") { //default path is set
115 string tryPath = m->getDefaultPath() + getSimpleName(candidateFileNames[i]);
116 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
117 ableToOpen = openInputFile(tryPath, in, "noerror");
118 candidateFileNames[i] = tryPath;
123 for (int j = 1; j < processors; j++) {
124 MPI_Send(&ableToOpen, 1, MPI_INT, j, 2001, MPI_COMM_WORLD);
128 MPI_Recv(&ableToOpen, 1, MPI_INT, 0, 2001, MPI_COMM_WORLD, &status);
133 if (ableToOpen == 1) {
134 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
135 //erase from file list
136 candidateFileNames.erase(candidateFileNames.begin()+i);
142 //make sure there is at least one valid file left
143 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
146 //check for optional parameter and set defaults
147 // ...at some point should added some additional type checking...
149 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
150 convert(temp, kmerSize);
152 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
153 convert(temp, match);
155 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
156 convert(temp, misMatch);
158 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
159 convert(temp, gapOpen);
161 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
162 convert(temp, gapExtend);
164 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
165 convert(temp, processors);
167 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
170 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
171 convert(temp, threshold);
173 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
175 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
179 catch(exception& e) {
180 m->errorOut(e, "AlignCommand", "AlignCommand");
185 //**********************************************************************************************************************
187 AlignCommand::~AlignCommand(){
189 if (abort == false) {
190 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
196 //**********************************************************************************************************************
198 void AlignCommand::help(){
200 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
201 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen and gapextend.\n");
202 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");
203 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");
204 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");
205 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
206 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
207 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
208 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
209 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
210 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");
211 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");
212 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");
213 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");
214 m->mothurOut("The align.seqs command should be in the following format: \n");
215 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
216 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");
217 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
219 catch(exception& e) {
220 m->errorOut(e, "AlignCommand", "help");
226 //**********************************************************************************************************************
228 int AlignCommand::execute(){
230 if (abort == true) { return 0; }
232 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
233 int longestBase = templateDB->getLongestBase();
235 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
236 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
237 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
238 else if(align == "noalign") { alignment = new NoAlign(); }
240 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
241 m->mothurOutEndLine();
242 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
244 vector<string> outputNames;
246 for (int s = 0; s < candidateFileNames.size(); s++) {
247 if (m->control_pressed) { return 0; }
249 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
251 if (outputDir == "") { outputDir += hasPath(candidateFileNames[s]); }
252 string alignFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align";
253 string reportFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align.report";
254 string accnosFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "flip.accnos";
255 bool hasAccnos = true;
257 int numFastaSeqs = 0;
258 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
259 int start = time(NULL);
262 int pid, end, numSeqsPerProcessor;
264 vector<unsigned long int> MPIPos;
265 MPIWroteAccnos = false;
268 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
269 MPI_Comm_size(MPI_COMM_WORLD, &processors);
272 MPI_File outMPIAlign;
273 MPI_File outMPIReport;
274 MPI_File outMPIAccnos;
276 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
277 int inMode=MPI_MODE_RDONLY;
279 char outAlignFilename[1024];
280 strcpy(outAlignFilename, alignFileName.c_str());
282 char outReportFilename[1024];
283 strcpy(outReportFilename, reportFileName.c_str());
285 char outAccnosFilename[1024];
286 strcpy(outAccnosFilename, accnosFileName.c_str());
288 char inFileName[1024];
289 strcpy(inFileName, candidateFileNames[s].c_str());
291 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
292 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
293 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
294 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
296 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
298 if (pid == 0) { //you are the root process
300 MPIPos = setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
302 //send file positions to all processes
303 for(int i = 1; i < processors; i++) {
304 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
305 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
308 //figure out how many sequences you have to align
309 numSeqsPerProcessor = numFastaSeqs / processors;
310 int startIndex = pid * numSeqsPerProcessor;
311 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
315 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
317 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
319 for (int i = 1; i < processors; i++) {
321 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
322 if (tempResult != 0) { MPIWroteAccnos = true; }
324 }else{ //you are a child process
325 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
326 MPIPos.resize(numFastaSeqs+1);
327 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
330 //figure out how many sequences you have to align
331 numSeqsPerProcessor = numFastaSeqs / processors;
332 int startIndex = pid * numSeqsPerProcessor;
333 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
337 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
339 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
341 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
345 MPI_File_close(&inMPI);
346 MPI_File_close(&outMPIAlign);
347 MPI_File_close(&outMPIReport);
348 MPI_File_close(&outMPIAccnos);
350 //delete accnos file if blank
352 //delete accnos file if its blank else report to user
353 if (MPIWroteAccnos) {
354 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
356 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
357 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
358 m->mothurOutEndLine();
361 //MPI_File_delete(outAccnosFilename, info);
363 remove(accnosFileName.c_str());
369 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
372 openInputFile(candidateFileNames[s], inFASTA);
373 getNumSeqs(inFASTA, numFastaSeqs);
376 lines.push_back(new linePair(0, numFastaSeqs));
378 driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
380 if (m->control_pressed) {
381 remove(accnosFileName.c_str());
382 remove(alignFileName.c_str());
383 remove(reportFileName.c_str());
387 //delete accnos file if its blank else report to user
388 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
390 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
392 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
393 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
394 m->mothurOutEndLine();
398 vector<unsigned long int> positions;
399 processIDS.resize(0);
402 openInputFile(candidateFileNames[s], inFASTA);
405 while(!inFASTA.eof()){
406 input = getline(inFASTA);
407 if (input.length() != 0) {
408 if(input[0] == '>'){ unsigned long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
413 numFastaSeqs = positions.size();
415 int numSeqsPerProcessor = numFastaSeqs / processors;
417 for (int i = 0; i < processors; i++) {
418 unsigned long int startPos = positions[ i * numSeqsPerProcessor ];
419 if(i == processors - 1){
420 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
422 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
425 createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
427 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
428 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
430 //append alignment and report files
431 for(int i=1;i<processors;i++){
432 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
433 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
435 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
436 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
439 vector<string> nonBlankAccnosFiles;
440 //delete blank accnos files generated with multiple processes
441 for(int i=0;i<processors;i++){
442 if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
443 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
444 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
447 //append accnos files
448 if (nonBlankAccnosFiles.size() != 0) {
449 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
451 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
452 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
453 remove(nonBlankAccnosFiles[h].c_str());
455 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
457 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
458 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
459 m->mothurOutEndLine();
460 }else{ hasAccnos = false; }
462 if (m->control_pressed) {
463 remove(accnosFileName.c_str());
464 remove(alignFileName.c_str());
465 remove(reportFileName.c_str());
471 openInputFile(candidateFileNames[s], inFASTA);
472 getNumSeqs(inFASTA, numFastaSeqs);
475 lines.push_back(new linePair(0, numFastaSeqs));
477 driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
479 if (m->control_pressed) {
480 remove(accnosFileName.c_str());
481 remove(alignFileName.c_str());
482 remove(reportFileName.c_str());
486 //delete accnos file if its blank else report to user
487 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
489 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
491 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
492 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
493 m->mothurOutEndLine();
502 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
504 if (pid == 0) { //only one process should output to screen
507 outputNames.push_back(alignFileName);
508 outputNames.push_back(reportFileName);
509 if (hasAccnos) { outputNames.push_back(accnosFileName); }
515 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
516 m->mothurOutEndLine();
517 m->mothurOutEndLine();
521 m->mothurOutEndLine();
522 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
523 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
524 m->mothurOutEndLine();
528 catch(exception& e) {
529 m->errorOut(e, "AlignCommand", "execute");
534 //**********************************************************************************************************************
536 int AlignCommand::driver(linePair* line, string alignFName, string reportFName, string accnosFName, string filename){
538 ofstream alignmentFile;
539 openOutputFile(alignFName, alignmentFile);
542 openOutputFile(accnosFName, accnosFile);
544 NastReport report(reportFName);
547 openInputFile(filename, inFASTA);
549 inFASTA.seekg(line->start);
551 for(int i=0;i<line->numSeqs;i++){
553 if (m->control_pressed) { return 0; }
555 Sequence* candidateSeq = new Sequence(inFASTA); gobble(inFASTA);
557 int origNumBases = candidateSeq->getNumBases();
558 string originalUnaligned = candidateSeq->getUnaligned();
559 int numBasesNeeded = origNumBases * threshold;
561 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
562 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
563 alignment->resize(candidateSeq->getUnaligned().length()+1);
566 Sequence temp = templateDB->findClosestSequence(candidateSeq);
567 Sequence* templateSeq = &temp;
569 float searchScore = templateDB->getSearchScore();
571 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
575 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
576 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
577 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
578 //so this bool tells you if you need to delete it
580 //if there is a possibility that this sequence should be reversed
581 if (candidateSeq->getNumBases() < numBasesNeeded) {
583 string wasBetter = "";
584 //if the user wants you to try the reverse
586 //get reverse compliment
587 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
588 copy->reverseComplement();
591 Sequence temp2 = templateDB->findClosestSequence(copy);
592 Sequence* templateSeq2 = &temp2;
594 searchScore = templateDB->getSearchScore();
596 nast2 = new Nast(alignment, copy, templateSeq2);
598 //check if any better
599 if (copy->getNumBases() > candidateSeq->getNumBases()) {
600 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
601 templateSeq = templateSeq2;
604 needToDeleteCopy = true;
606 wasBetter = "\treverse complement did NOT produce a better alignment, 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; }
631 if((i+1) % 100 == 0){ m->mothurOut(toString(i+1)); m->mothurOutEndLine(); }
634 if((line->numSeqs) % 100 != 0){ m->mothurOut(toString(line->numSeqs)); m->mothurOutEndLine(); }
636 alignmentFile.close();
642 catch(exception& e) {
643 m->errorOut(e, "AlignCommand", "driver");
647 //**********************************************************************************************************************
649 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
651 string outputString = "";
652 MPI_Status statusReport;
653 MPI_Status statusAlign;
654 MPI_Status statusAccnos;
657 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
662 outputString = report.getHeaders();
663 int length = outputString.length();
665 char* buf = new char[length];
666 memcpy(buf, outputString.c_str(), length);
668 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
673 for(int i=0;i<num;i++){
675 if (m->control_pressed) { return 0; }
678 int length = MPIPos[start+i+1] - MPIPos[start+i];
680 char* buf4 = new char[length];
681 memcpy(buf4, outputString.c_str(), length);
683 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
685 string tempBuf = buf4;
689 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
690 istringstream iss (tempBuf,istringstream::in);
692 Sequence* candidateSeq = new Sequence(iss);
693 int origNumBases = candidateSeq->getNumBases();
694 string originalUnaligned = candidateSeq->getUnaligned();
695 int numBasesNeeded = origNumBases * threshold;
697 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
698 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
699 alignment->resize(candidateSeq->getUnaligned().length()+1);
702 Sequence temp = templateDB->findClosestSequence(candidateSeq);
703 Sequence* templateSeq = &temp;
705 float searchScore = templateDB->getSearchScore();
707 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
711 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
712 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
713 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
714 //so this bool tells you if you need to delete it
716 //if there is a possibility that this sequence should be reversed
717 if (candidateSeq->getNumBases() < numBasesNeeded) {
719 string wasBetter = "";
720 //if the user wants you to try the reverse
722 //get reverse compliment
723 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
724 copy->reverseComplement();
727 Sequence temp2 = templateDB->findClosestSequence(copy);
728 Sequence* templateSeq2 = &temp2;
730 searchScore = templateDB->getSearchScore();
732 nast2 = new Nast(alignment, copy, templateSeq2);
734 //check if any better
735 if (copy->getNumBases() > candidateSeq->getNumBases()) {
736 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
737 templateSeq = templateSeq2;
740 needToDeleteCopy = true;
742 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
748 //create accnos file with names
749 outputString = candidateSeq->getName() + wasBetter + "\n";
751 //send results to parent
752 int length = outputString.length();
754 char* buf = new char[length];
755 memcpy(buf, outputString.c_str(), length);
757 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
759 MPIWroteAccnos = true;
762 report.setCandidate(candidateSeq);
763 report.setTemplate(templateSeq);
764 report.setSearchParameters(search, searchScore);
765 report.setAlignmentParameters(align, alignment);
766 report.setNastParameters(*nast);
768 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
770 //send results to parent
771 int length = outputString.length();
772 char* buf2 = new char[length];
773 memcpy(buf2, outputString.c_str(), length);
775 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
779 outputString = report.getReport();
781 //send results to parent
782 length = outputString.length();
783 char* buf3 = new char[length];
784 memcpy(buf3, outputString.c_str(), length);
786 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
790 if (needToDeleteCopy) { delete copy; }
795 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
798 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
802 catch(exception& e) {
803 m->errorOut(e, "AlignCommand", "driverMPI");
808 /**************************************************************************************************/
810 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
812 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
815 // processIDS.resize(0);
817 //loop through and create all the processes you want
818 while (process != processors) {
822 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
825 exitCommand = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
827 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
830 //force parent to wait until all the processes are done
831 for (int i=0;i<processors;i++) {
832 int temp = processIDS[i];
839 catch(exception& e) {
840 m->errorOut(e, "AlignCommand", "createProcesses");
845 /**************************************************************************************************/
847 void AlignCommand::appendAlignFiles(string temp, string filename) {
852 openOutputFileAppend(filename, output);
853 openInputFile(temp, input);
855 while(char c = input.get()){
856 if(input.eof()) { break; }
857 else { output << c; }
863 catch(exception& e) {
864 m->errorOut(e, "AlignCommand", "appendAlignFiles");
868 //**********************************************************************************************************************
870 void AlignCommand::appendReportFiles(string temp, string filename) {
875 openOutputFileAppend(filename, output);
876 openInputFile(temp, input);
878 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
880 while(char c = input.get()){
881 if(input.eof()) { break; }
882 else { output << c; }
888 catch(exception& e) {
889 m->errorOut(e, "AlignCommand", "appendReportFiles");
893 //**********************************************************************************************************************