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 abort = true; calledHelp = true;
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) {
81 abort = false; calledHelp = false;
83 //allow user to run help
84 if(option == "help") { help(); abort = true; calledHelp = 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") {
141 //check currentFiles for a fasta file
142 if (currentFiles->getFastaFile() != "") { candidateFileName = currentFiles->getFastaFile(); m->mothurOut("Using " + candidateFileName + " as candidate file."); m->mothurOutEndLine();
143 }else { m->mothurOut("candidate is a required parameter for the align.seqs command."); m->mothurOutEndLine(); abort = true; }
145 m->splitAtDash(candidateFileName, candidateFileNames);
147 //go through files and make sure they are good, if not, then disregard them
148 for (int i = 0; i < candidateFileNames.size(); i++) {
149 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
151 if (inputDir != "") {
152 string path = m->hasPath(candidateFileNames[i]);
153 //if the user has not given a path then, add inputdir. else leave path alone.
154 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
159 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
162 //if you can't open it, try default location
163 if (ableToOpen == 1) {
164 if (m->getDefaultPath() != "") { //default path is set
165 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
166 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
168 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
170 candidateFileNames[i] = tryPath;
174 //if you can't open it, try output location
175 if (ableToOpen == 1) {
176 if (m->getOutputDir() != "") { //default path is set
177 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
178 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
180 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
182 candidateFileNames[i] = tryPath;
188 if (ableToOpen == 1) {
189 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
190 //erase from file list
191 candidateFileNames.erase(candidateFileNames.begin()+i);
197 //make sure there is at least one valid file left
198 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
201 //check for optional parameter and set defaults
202 // ...at some point should added some additional type checking...
204 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
205 convert(temp, kmerSize);
207 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
208 convert(temp, match);
210 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
211 convert(temp, misMatch);
213 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
214 convert(temp, gapOpen);
216 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
217 convert(temp, gapExtend);
219 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
220 convert(temp, processors);
222 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
223 flip = m->isTrue(temp);
225 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
226 convert(temp, threshold);
228 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
230 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
234 catch(exception& e) {
235 m->errorOut(e, "AlignCommand", "AlignCommand");
239 //**********************************************************************************************************************
241 AlignCommand::~AlignCommand(){
243 if (abort == false) {
244 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
250 //**********************************************************************************************************************
252 void AlignCommand::help(){
254 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
255 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen, gapextend and processors.\n");
256 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");
257 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");
258 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");
259 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
260 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
261 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
262 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
263 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
264 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");
265 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");
266 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");
267 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");
268 m->mothurOut("The align.seqs command should be in the following format: \n");
269 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
270 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");
271 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
273 catch(exception& e) {
274 m->errorOut(e, "AlignCommand", "help");
280 //**********************************************************************************************************************
282 int AlignCommand::execute(){
284 if (abort == true) { if (calledHelp) { return 0; } return 2; }
286 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
287 int longestBase = templateDB->getLongestBase();
289 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
290 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
291 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
292 else if(align == "noalign") { alignment = new NoAlign(); }
294 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
295 m->mothurOutEndLine();
296 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
299 for (int s = 0; s < candidateFileNames.size(); s++) {
300 if (m->control_pressed) { outputTypes.clear(); return 0; }
302 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
304 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
305 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
306 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
307 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
308 bool hasAccnos = true;
310 int numFastaSeqs = 0;
311 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
312 int start = time(NULL);
315 int pid, numSeqsPerProcessor;
317 vector<unsigned long int> MPIPos;
318 MPIWroteAccnos = false;
321 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
322 MPI_Comm_size(MPI_COMM_WORLD, &processors);
325 MPI_File outMPIAlign;
326 MPI_File outMPIReport;
327 MPI_File outMPIAccnos;
329 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
330 int inMode=MPI_MODE_RDONLY;
332 char outAlignFilename[1024];
333 strcpy(outAlignFilename, alignFileName.c_str());
335 char outReportFilename[1024];
336 strcpy(outReportFilename, reportFileName.c_str());
338 char outAccnosFilename[1024];
339 strcpy(outAccnosFilename, accnosFileName.c_str());
341 char inFileName[1024];
342 strcpy(inFileName, candidateFileNames[s].c_str());
344 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
345 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
346 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
347 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
349 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
351 if (pid == 0) { //you are the root process
353 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
355 //send file positions to all processes
356 for(int i = 1; i < processors; i++) {
357 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
358 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
361 //figure out how many sequences you have to align
362 numSeqsPerProcessor = numFastaSeqs / processors;
363 int startIndex = pid * numSeqsPerProcessor;
364 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
367 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
369 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
371 for (int i = 1; i < processors; i++) {
373 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
374 if (tempResult != 0) { MPIWroteAccnos = true; }
376 }else{ //you are a child process
377 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
378 MPIPos.resize(numFastaSeqs+1);
379 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
382 //figure out how many sequences you have to align
383 numSeqsPerProcessor = numFastaSeqs / processors;
384 int startIndex = pid * numSeqsPerProcessor;
385 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
389 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
391 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
393 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
397 MPI_File_close(&inMPI);
398 MPI_File_close(&outMPIAlign);
399 MPI_File_close(&outMPIReport);
400 MPI_File_close(&outMPIAccnos);
402 //delete accnos file if blank
404 //delete accnos file if its blank else report to user
405 if (MPIWroteAccnos) {
406 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
408 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
409 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
410 m->mothurOutEndLine();
413 //MPI_File_delete(outAccnosFilename, info);
415 remove(accnosFileName.c_str());
421 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
422 for (int i = 0; i < (positions.size()-1); i++) {
423 lines.push_back(new linePair(positions[i], positions[(i+1)]));
425 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
427 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
429 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
432 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
435 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
437 //delete accnos file if its blank else report to user
438 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
440 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
442 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
443 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
444 m->mothurOutEndLine();
451 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
453 if (pid == 0) { //only one process should output to screen
456 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
457 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
458 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
464 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
465 m->mothurOutEndLine();
466 m->mothurOutEndLine();
469 //set align file as new current fastafile
470 string currentFasta = "";
471 itTypes = outputTypes.find("fasta");
472 if (itTypes != outputTypes.end()) {
473 if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; }
475 currentFiles->setFastaFile(currentFasta);
476 cout << "current fasta = " << currentFiles->getFastaFile() << endl;
477 m->mothurOutEndLine();
478 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
479 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
480 m->mothurOutEndLine();
484 catch(exception& e) {
485 m->errorOut(e, "AlignCommand", "execute");
490 //**********************************************************************************************************************
492 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
494 ofstream alignmentFile;
495 m->openOutputFile(alignFName, alignmentFile);
498 m->openOutputFile(accnosFName, accnosFile);
500 NastReport report(reportFName);
503 m->openInputFile(filename, inFASTA);
505 inFASTA.seekg(filePos->start);
512 if (m->control_pressed) { return 0; }
514 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
515 report.setCandidate(candidateSeq);
517 int origNumBases = candidateSeq->getNumBases();
518 string originalUnaligned = candidateSeq->getUnaligned();
519 int numBasesNeeded = origNumBases * threshold;
521 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
522 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
523 alignment->resize(candidateSeq->getUnaligned().length()+1);
526 Sequence temp = templateDB->findClosestSequence(candidateSeq);
527 Sequence* templateSeq = &temp;
529 float searchScore = templateDB->getSearchScore();
531 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
536 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
537 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
538 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
539 //so this bool tells you if you need to delete it
541 //if there is a possibility that this sequence should be reversed
542 if (candidateSeq->getNumBases() < numBasesNeeded) {
544 string wasBetter = "";
545 //if the user wants you to try the reverse
548 //get reverse compliment
549 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
550 copy->reverseComplement();
553 Sequence temp2 = templateDB->findClosestSequence(copy);
554 Sequence* templateSeq2 = &temp2;
556 searchScore = templateDB->getSearchScore();
558 nast2 = new Nast(alignment, copy, templateSeq2);
560 //check if any better
561 if (copy->getNumBases() > candidateSeq->getNumBases()) {
562 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
563 templateSeq = templateSeq2;
566 needToDeleteCopy = true;
567 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
569 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
575 //create accnos file with names
576 accnosFile << candidateSeq->getName() << wasBetter << endl;
579 report.setTemplate(templateSeq);
580 report.setSearchParameters(search, searchScore);
581 report.setAlignmentParameters(align, alignment);
582 report.setNastParameters(*nast);
584 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
588 if (needToDeleteCopy) { delete copy; }
594 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
595 unsigned long int pos = inFASTA.tellg();
596 if ((pos == -1) || (pos >= filePos->end)) { break; }
598 if (inFASTA.eof()) { break; }
602 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
606 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
608 alignmentFile.close();
614 catch(exception& e) {
615 m->errorOut(e, "AlignCommand", "driver");
619 //**********************************************************************************************************************
621 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
623 string outputString = "";
624 MPI_Status statusReport;
625 MPI_Status statusAlign;
626 MPI_Status statusAccnos;
629 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
634 outputString = report.getHeaders();
635 int length = outputString.length();
637 char* buf = new char[length];
638 memcpy(buf, outputString.c_str(), length);
640 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
645 for(int i=0;i<num;i++){
647 if (m->control_pressed) { return 0; }
650 int length = MPIPos[start+i+1] - MPIPos[start+i];
652 char* buf4 = new char[length];
653 //memcpy(buf4, outputString.c_str(), length);
655 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
657 string tempBuf = buf4;
661 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
663 istringstream iss (tempBuf,istringstream::in);
665 Sequence* candidateSeq = new Sequence(iss);
666 report.setCandidate(candidateSeq);
668 int origNumBases = candidateSeq->getNumBases();
669 string originalUnaligned = candidateSeq->getUnaligned();
670 int numBasesNeeded = origNumBases * threshold;
672 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
673 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
674 alignment->resize(candidateSeq->getUnaligned().length()+1);
677 Sequence temp = templateDB->findClosestSequence(candidateSeq);
678 Sequence* templateSeq = &temp;
680 float searchScore = templateDB->getSearchScore();
682 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
686 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
687 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
688 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
689 //so this bool tells you if you need to delete it
691 //if there is a possibility that this sequence should be reversed
692 if (candidateSeq->getNumBases() < numBasesNeeded) {
694 string wasBetter = "";
695 //if the user wants you to try the reverse
697 //get reverse compliment
698 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
699 copy->reverseComplement();
702 Sequence temp2 = templateDB->findClosestSequence(copy);
703 Sequence* templateSeq2 = &temp2;
705 searchScore = templateDB->getSearchScore();
707 nast2 = new Nast(alignment, copy, templateSeq2);
709 //check if any better
710 if (copy->getNumBases() > candidateSeq->getNumBases()) {
711 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
712 templateSeq = templateSeq2;
715 needToDeleteCopy = true;
716 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
718 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
724 //create accnos file with names
725 outputString = candidateSeq->getName() + wasBetter + "\n";
727 //send results to parent
728 int length = outputString.length();
730 char* buf = new char[length];
731 memcpy(buf, outputString.c_str(), length);
733 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
735 MPIWroteAccnos = true;
738 report.setTemplate(templateSeq);
739 report.setSearchParameters(search, searchScore);
740 report.setAlignmentParameters(align, alignment);
741 report.setNastParameters(*nast);
743 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
745 //send results to parent
746 int length = outputString.length();
747 char* buf2 = new char[length];
748 memcpy(buf2, outputString.c_str(), length);
750 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
754 outputString = report.getReport();
756 //send results to parent
757 length = outputString.length();
758 char* buf3 = new char[length];
759 memcpy(buf3, outputString.c_str(), length);
761 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
765 if (needToDeleteCopy) { delete copy; }
770 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
773 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
777 catch(exception& e) {
778 m->errorOut(e, "AlignCommand", "driverMPI");
783 /**************************************************************************************************/
785 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
787 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
788 processIDS.resize(0);
791 // processIDS.resize(0);
793 //loop through and create all the processes you want
794 while (process != processors) {
798 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
801 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
803 //pass numSeqs to parent
805 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
806 m->openOutputFile(tempFile, out);
812 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
813 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
819 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
821 //force parent to wait until all the processes are done
822 for (int i=0;i<processors;i++) {
823 int temp = processIDS[i];
827 vector<string> nonBlankAccnosFiles;
828 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
829 else { remove(accnosFName.c_str()); } //remove so other files can be renamed to it
831 for (int i = 0; i < processIDS.size(); i++) {
833 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
834 m->openInputFile(tempFile, in);
835 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
836 in.close(); remove(tempFile.c_str());
838 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
839 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
841 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
842 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
844 if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
845 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
846 }else { remove((accnosFName + toString(processIDS[i]) + ".temp").c_str()); }
850 //append accnos files
851 if (nonBlankAccnosFiles.size() != 0) {
852 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
854 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
855 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
856 remove(nonBlankAccnosFiles[h].c_str());
858 }else { //recreate the accnosfile if needed
860 m->openOutputFile(accnosFName, out);
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 //**********************************************************************************************************************