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;
82 currentFiles = CurrentFile::getInstance();
84 //allow user to run help
85 if(option == "help") { help(); abort = true; calledHelp = true;}
89 //valid paramters for this command
90 string AlignArray[] = {"template","candidate","search","ksize","align","match","mismatch","gapopen","gapextend", "processors","flip","threshold","outputdir","inputdir"};
91 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
93 OptionParser parser(option);
94 map<string, string> parameters = parser.getParameters();
96 ValidParameters validParameter("align.seqs");
97 map<string, string>::iterator it;
99 //check to make sure all parameters are valid for command
100 for (it = parameters.begin(); it != parameters.end(); it++) {
101 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
104 //initialize outputTypes
105 vector<string> tempOutNames;
106 outputTypes["fasta"] = tempOutNames;
107 outputTypes["alignreport"] = tempOutNames;
108 outputTypes["accnos"] = tempOutNames;
110 //if the user changes the output directory command factory will send this info to us in the output parameter
111 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
114 //if the user changes the input directory command factory will send this info to us in the output parameter
115 string inputDir = validParameter.validFile(parameters, "inputdir", false);
117 if (inputDir == "not found"){ inputDir = ""; }
121 it = parameters.find("template");
123 //user has given a template file
124 if(it != parameters.end()){
125 path = m->hasPath(it->second);
126 //if the user has not given a path then, add inputdir. else leave path alone.
127 if (path == "") { parameters["template"] = inputDir + it->second; }
131 //check for required parameters
132 templateFileName = validParameter.validFile(parameters, "template", true);
134 if (templateFileName == "not found") {
135 m->mothurOut("template is a required parameter for the align.seqs command.");
136 m->mothurOutEndLine();
138 }else if (templateFileName == "not open") { abort = true; }
140 candidateFileName = validParameter.validFile(parameters, "candidate", false);
141 if (candidateFileName == "not found") {
142 //check currentFiles for a fasta file
143 if (currentFiles->getFastaFile() != "") { candidateFileName = currentFiles->getFastaFile(); m->mothurOut("Using " + candidateFileName + " as candidate file."); m->mothurOutEndLine();
144 }else { m->mothurOut("candidate is a required parameter for the align.seqs command."); m->mothurOutEndLine(); abort = true; }
146 m->splitAtDash(candidateFileName, candidateFileNames);
148 //go through files and make sure they are good, if not, then disregard them
149 for (int i = 0; i < candidateFileNames.size(); i++) {
150 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
152 if (inputDir != "") {
153 string path = m->hasPath(candidateFileNames[i]);
154 //if the user has not given a path then, add inputdir. else leave path alone.
155 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
160 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
163 //if you can't open it, try default location
164 if (ableToOpen == 1) {
165 if (m->getDefaultPath() != "") { //default path is set
166 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
167 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
169 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
171 candidateFileNames[i] = tryPath;
175 //if you can't open it, try output location
176 if (ableToOpen == 1) {
177 if (m->getOutputDir() != "") { //default path is set
178 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
179 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
181 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
183 candidateFileNames[i] = tryPath;
189 if (ableToOpen == 1) {
190 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
191 //erase from file list
192 candidateFileNames.erase(candidateFileNames.begin()+i);
198 //make sure there is at least one valid file left
199 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
202 //check for optional parameter and set defaults
203 // ...at some point should added some additional type checking...
205 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
206 convert(temp, kmerSize);
208 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
209 convert(temp, match);
211 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
212 convert(temp, misMatch);
214 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
215 convert(temp, gapOpen);
217 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
218 convert(temp, gapExtend);
220 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
221 convert(temp, processors);
223 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
224 flip = m->isTrue(temp);
226 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
227 convert(temp, threshold);
229 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
231 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
235 catch(exception& e) {
236 m->errorOut(e, "AlignCommand", "AlignCommand");
240 //**********************************************************************************************************************
242 AlignCommand::~AlignCommand(){
244 if (abort == false) {
245 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
251 //**********************************************************************************************************************
253 void AlignCommand::help(){
255 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
256 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen, gapextend and processors.\n");
257 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");
258 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");
259 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");
260 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
261 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
262 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
263 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
264 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
265 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");
266 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");
267 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");
268 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");
269 m->mothurOut("The align.seqs command should be in the following format: \n");
270 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
271 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");
272 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
274 catch(exception& e) {
275 m->errorOut(e, "AlignCommand", "help");
281 //**********************************************************************************************************************
283 int AlignCommand::execute(){
285 if (abort == true) { if (calledHelp) { return 0; } return 2; }
287 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
288 int longestBase = templateDB->getLongestBase();
290 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
291 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
292 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
293 else if(align == "noalign") { alignment = new NoAlign(); }
295 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
296 m->mothurOutEndLine();
297 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
300 for (int s = 0; s < candidateFileNames.size(); s++) {
301 if (m->control_pressed) { outputTypes.clear(); return 0; }
303 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
305 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
306 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
307 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
308 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
309 bool hasAccnos = true;
311 int numFastaSeqs = 0;
312 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
313 int start = time(NULL);
316 int pid, numSeqsPerProcessor;
318 vector<unsigned long int> MPIPos;
319 MPIWroteAccnos = false;
322 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
323 MPI_Comm_size(MPI_COMM_WORLD, &processors);
326 MPI_File outMPIAlign;
327 MPI_File outMPIReport;
328 MPI_File outMPIAccnos;
330 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
331 int inMode=MPI_MODE_RDONLY;
333 char outAlignFilename[1024];
334 strcpy(outAlignFilename, alignFileName.c_str());
336 char outReportFilename[1024];
337 strcpy(outReportFilename, reportFileName.c_str());
339 char outAccnosFilename[1024];
340 strcpy(outAccnosFilename, accnosFileName.c_str());
342 char inFileName[1024];
343 strcpy(inFileName, candidateFileNames[s].c_str());
345 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
346 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
347 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
348 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
350 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
352 if (pid == 0) { //you are the root process
354 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
356 //send file positions to all processes
357 for(int i = 1; i < processors; i++) {
358 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
359 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
362 //figure out how many sequences you have to align
363 numSeqsPerProcessor = numFastaSeqs / processors;
364 int startIndex = pid * numSeqsPerProcessor;
365 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
368 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
370 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
372 for (int i = 1; i < processors; i++) {
374 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
375 if (tempResult != 0) { MPIWroteAccnos = true; }
377 }else{ //you are a child process
378 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
379 MPIPos.resize(numFastaSeqs+1);
380 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
383 //figure out how many sequences you have to align
384 numSeqsPerProcessor = numFastaSeqs / processors;
385 int startIndex = pid * numSeqsPerProcessor;
386 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
390 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
392 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
394 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
398 MPI_File_close(&inMPI);
399 MPI_File_close(&outMPIAlign);
400 MPI_File_close(&outMPIReport);
401 MPI_File_close(&outMPIAccnos);
403 //delete accnos file if blank
405 //delete accnos file if its blank else report to user
406 if (MPIWroteAccnos) {
407 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
409 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
410 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
411 m->mothurOutEndLine();
414 //MPI_File_delete(outAccnosFilename, info);
416 remove(accnosFileName.c_str());
422 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
423 for (int i = 0; i < (positions.size()-1); i++) {
424 lines.push_back(new linePair(positions[i], positions[(i+1)]));
426 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
428 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
430 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
433 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
436 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear(); return 0; }
438 //delete accnos file if its blank else report to user
439 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
441 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
443 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
444 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
445 m->mothurOutEndLine();
452 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
454 if (pid == 0) { //only one process should output to screen
457 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
458 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
459 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
465 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
466 m->mothurOutEndLine();
467 m->mothurOutEndLine();
470 //set align file as new current fastafile
471 string currentFasta = "";
472 itTypes = outputTypes.find("fasta");
473 if (itTypes != outputTypes.end()) {
474 if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; }
476 currentFiles->setFastaFile(currentFasta);
477 cout << "current fasta = " << currentFiles->getFastaFile() << endl;
478 m->mothurOutEndLine();
479 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
480 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
481 m->mothurOutEndLine();
485 catch(exception& e) {
486 m->errorOut(e, "AlignCommand", "execute");
491 //**********************************************************************************************************************
493 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
495 ofstream alignmentFile;
496 m->openOutputFile(alignFName, alignmentFile);
499 m->openOutputFile(accnosFName, accnosFile);
501 NastReport report(reportFName);
504 m->openInputFile(filename, inFASTA);
506 inFASTA.seekg(filePos->start);
513 if (m->control_pressed) { return 0; }
515 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
516 report.setCandidate(candidateSeq);
518 int origNumBases = candidateSeq->getNumBases();
519 string originalUnaligned = candidateSeq->getUnaligned();
520 int numBasesNeeded = origNumBases * threshold;
522 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
523 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
524 alignment->resize(candidateSeq->getUnaligned().length()+1);
527 Sequence temp = templateDB->findClosestSequence(candidateSeq);
528 Sequence* templateSeq = &temp;
530 float searchScore = templateDB->getSearchScore();
532 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
537 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
538 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
539 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
540 //so this bool tells you if you need to delete it
542 //if there is a possibility that this sequence should be reversed
543 if (candidateSeq->getNumBases() < numBasesNeeded) {
545 string wasBetter = "";
546 //if the user wants you to try the reverse
549 //get reverse compliment
550 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
551 copy->reverseComplement();
554 Sequence temp2 = templateDB->findClosestSequence(copy);
555 Sequence* templateSeq2 = &temp2;
557 searchScore = templateDB->getSearchScore();
559 nast2 = new Nast(alignment, copy, templateSeq2);
561 //check if any better
562 if (copy->getNumBases() > candidateSeq->getNumBases()) {
563 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
564 templateSeq = templateSeq2;
567 needToDeleteCopy = true;
568 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
570 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
576 //create accnos file with names
577 accnosFile << candidateSeq->getName() << wasBetter << endl;
580 report.setTemplate(templateSeq);
581 report.setSearchParameters(search, searchScore);
582 report.setAlignmentParameters(align, alignment);
583 report.setNastParameters(*nast);
585 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
589 if (needToDeleteCopy) { delete copy; }
595 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
596 unsigned long int pos = inFASTA.tellg();
597 if ((pos == -1) || (pos >= filePos->end)) { break; }
599 if (inFASTA.eof()) { break; }
603 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
607 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
609 alignmentFile.close();
615 catch(exception& e) {
616 m->errorOut(e, "AlignCommand", "driver");
620 //**********************************************************************************************************************
622 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
624 string outputString = "";
625 MPI_Status statusReport;
626 MPI_Status statusAlign;
627 MPI_Status statusAccnos;
630 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
635 outputString = report.getHeaders();
636 int length = outputString.length();
638 char* buf = new char[length];
639 memcpy(buf, outputString.c_str(), length);
641 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
646 for(int i=0;i<num;i++){
648 if (m->control_pressed) { return 0; }
651 int length = MPIPos[start+i+1] - MPIPos[start+i];
653 char* buf4 = new char[length];
654 //memcpy(buf4, outputString.c_str(), length);
656 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
658 string tempBuf = buf4;
662 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
664 istringstream iss (tempBuf,istringstream::in);
666 Sequence* candidateSeq = new Sequence(iss);
667 report.setCandidate(candidateSeq);
669 int origNumBases = candidateSeq->getNumBases();
670 string originalUnaligned = candidateSeq->getUnaligned();
671 int numBasesNeeded = origNumBases * threshold;
673 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
674 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
675 alignment->resize(candidateSeq->getUnaligned().length()+1);
678 Sequence temp = templateDB->findClosestSequence(candidateSeq);
679 Sequence* templateSeq = &temp;
681 float searchScore = templateDB->getSearchScore();
683 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
687 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
688 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
689 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
690 //so this bool tells you if you need to delete it
692 //if there is a possibility that this sequence should be reversed
693 if (candidateSeq->getNumBases() < numBasesNeeded) {
695 string wasBetter = "";
696 //if the user wants you to try the reverse
698 //get reverse compliment
699 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
700 copy->reverseComplement();
703 Sequence temp2 = templateDB->findClosestSequence(copy);
704 Sequence* templateSeq2 = &temp2;
706 searchScore = templateDB->getSearchScore();
708 nast2 = new Nast(alignment, copy, templateSeq2);
710 //check if any better
711 if (copy->getNumBases() > candidateSeq->getNumBases()) {
712 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
713 templateSeq = templateSeq2;
716 needToDeleteCopy = true;
717 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
719 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
725 //create accnos file with names
726 outputString = candidateSeq->getName() + wasBetter + "\n";
728 //send results to parent
729 int length = outputString.length();
731 char* buf = new char[length];
732 memcpy(buf, outputString.c_str(), length);
734 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
736 MPIWroteAccnos = true;
739 report.setTemplate(templateSeq);
740 report.setSearchParameters(search, searchScore);
741 report.setAlignmentParameters(align, alignment);
742 report.setNastParameters(*nast);
744 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
746 //send results to parent
747 int length = outputString.length();
748 char* buf2 = new char[length];
749 memcpy(buf2, outputString.c_str(), length);
751 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
755 outputString = report.getReport();
757 //send results to parent
758 length = outputString.length();
759 char* buf3 = new char[length];
760 memcpy(buf3, outputString.c_str(), length);
762 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
766 if (needToDeleteCopy) { delete copy; }
771 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
774 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
778 catch(exception& e) {
779 m->errorOut(e, "AlignCommand", "driverMPI");
784 /**************************************************************************************************/
786 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
788 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
789 processIDS.resize(0);
792 // processIDS.resize(0);
794 //loop through and create all the processes you want
795 while (process != processors) {
799 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
802 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
804 //pass numSeqs to parent
806 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
807 m->openOutputFile(tempFile, out);
813 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
814 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
820 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
822 //force parent to wait until all the processes are done
823 for (int i=0;i<processors;i++) {
824 int temp = processIDS[i];
828 vector<string> nonBlankAccnosFiles;
829 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
830 else { remove(accnosFName.c_str()); } //remove so other files can be renamed to it
832 for (int i = 0; i < processIDS.size(); i++) {
834 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
835 m->openInputFile(tempFile, in);
836 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
837 in.close(); remove(tempFile.c_str());
839 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
840 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
842 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
843 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
845 if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
846 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
847 }else { remove((accnosFName + toString(processIDS[i]) + ".temp").c_str()); }
851 //append accnos files
852 if (nonBlankAccnosFiles.size() != 0) {
853 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
855 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
856 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
857 remove(nonBlankAccnosFiles[h].c_str());
859 }else { //recreate the accnosfile if needed
861 m->openOutputFile(accnosFName, out);
868 catch(exception& e) {
869 m->errorOut(e, "AlignCommand", "createProcesses");
873 /**************************************************************************************************/
875 void AlignCommand::appendAlignFiles(string temp, string filename) {
880 m->openOutputFileAppend(filename, output);
881 m->openInputFile(temp, input);
883 while(char c = input.get()){
884 if(input.eof()) { break; }
885 else { output << c; }
891 catch(exception& e) {
892 m->errorOut(e, "AlignCommand", "appendAlignFiles");
896 //**********************************************************************************************************************
898 void AlignCommand::appendReportFiles(string temp, string filename) {
903 m->openOutputFileAppend(filename, output);
904 m->openInputFile(temp, input);
906 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
908 while(char c = input.get()){
909 if(input.eof()) { break; }
910 else { output << c; }
916 catch(exception& e) {
917 m->errorOut(e, "AlignCommand", "appendReportFiles");
921 //**********************************************************************************************************************