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 = m->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 m->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 = m->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]; }
102 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
104 //if you can't open it, try default location
105 if (ableToOpen == 1) {
106 if (m->getDefaultPath() != "") { //default path is set
107 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
108 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
109 ableToOpen = m->openInputFile(tryPath, in, "noerror");
110 candidateFileNames[i] = tryPath;
115 if (ableToOpen == 1) {
116 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
117 //erase from file list
118 candidateFileNames.erase(candidateFileNames.begin()+i);
124 //make sure there is at least one valid file left
125 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
128 //check for optional parameter and set defaults
129 // ...at some point should added some additional type checking...
131 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
132 convert(temp, kmerSize);
134 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
135 convert(temp, match);
137 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
138 convert(temp, misMatch);
140 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
141 convert(temp, gapOpen);
143 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
144 convert(temp, gapExtend);
146 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
147 convert(temp, processors);
149 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
150 flip = m->isTrue(temp);
152 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
153 convert(temp, threshold);
155 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
157 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
161 catch(exception& e) {
162 m->errorOut(e, "AlignCommand", "AlignCommand");
167 //**********************************************************************************************************************
169 AlignCommand::~AlignCommand(){
171 if (abort == false) {
172 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
178 //**********************************************************************************************************************
180 void AlignCommand::help(){
182 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
183 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen and gapextend.\n");
184 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");
185 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");
186 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");
187 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
188 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
189 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
190 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
191 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
192 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");
193 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");
194 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");
195 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");
196 m->mothurOut("The align.seqs command should be in the following format: \n");
197 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
198 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");
199 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
201 catch(exception& e) {
202 m->errorOut(e, "AlignCommand", "help");
208 //**********************************************************************************************************************
210 int AlignCommand::execute(){
212 if (abort == true) { return 0; }
214 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
215 int longestBase = templateDB->getLongestBase();
217 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
218 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
219 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
220 else if(align == "noalign") { alignment = new NoAlign(); }
222 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
223 m->mothurOutEndLine();
224 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
226 vector<string> outputNames;
228 for (int s = 0; s < candidateFileNames.size(); s++) {
229 if (m->control_pressed) { return 0; }
231 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
233 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
234 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
235 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
236 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
237 bool hasAccnos = true;
239 int numFastaSeqs = 0;
240 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
241 int start = time(NULL);
244 int pid, end, numSeqsPerProcessor;
246 vector<unsigned long int> MPIPos;
247 MPIWroteAccnos = false;
250 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
251 MPI_Comm_size(MPI_COMM_WORLD, &processors);
254 MPI_File outMPIAlign;
255 MPI_File outMPIReport;
256 MPI_File outMPIAccnos;
258 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
259 int inMode=MPI_MODE_RDONLY;
261 char outAlignFilename[1024];
262 strcpy(outAlignFilename, alignFileName.c_str());
264 char outReportFilename[1024];
265 strcpy(outReportFilename, reportFileName.c_str());
267 char outAccnosFilename[1024];
268 strcpy(outAccnosFilename, accnosFileName.c_str());
270 char inFileName[1024];
271 strcpy(inFileName, candidateFileNames[s].c_str());
273 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
274 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
275 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
276 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
278 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
280 if (pid == 0) { //you are the root process
282 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
284 //send file positions to all processes
285 for(int i = 1; i < processors; i++) {
286 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
287 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
290 //figure out how many sequences you have to align
291 numSeqsPerProcessor = numFastaSeqs / processors;
292 int startIndex = pid * numSeqsPerProcessor;
293 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
296 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
298 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
300 for (int i = 1; i < processors; i++) {
302 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
303 if (tempResult != 0) { MPIWroteAccnos = true; }
305 }else{ //you are a child process
306 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
307 MPIPos.resize(numFastaSeqs+1);
308 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
311 //figure out how many sequences you have to align
312 numSeqsPerProcessor = numFastaSeqs / processors;
313 int startIndex = pid * numSeqsPerProcessor;
314 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
318 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
320 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
322 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
326 MPI_File_close(&inMPI);
327 MPI_File_close(&outMPIAlign);
328 MPI_File_close(&outMPIReport);
329 MPI_File_close(&outMPIAccnos);
331 //delete accnos file if blank
333 //delete accnos file if its blank else report to user
334 if (MPIWroteAccnos) {
335 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
337 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
338 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
339 m->mothurOutEndLine();
342 //MPI_File_delete(outAccnosFilename, info);
344 remove(accnosFileName.c_str());
350 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
352 for (int i = 0; i < (positions.size()-1); i++) {
353 lines.push_back(new linePair(positions[i], positions[(i+1)]));
355 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
357 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
359 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
361 //delete accnos file if its blank else report to user
362 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
364 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
366 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
367 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
368 m->mothurOutEndLine();
371 processIDS.resize(0);
373 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
375 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
376 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
378 //append alignment and report files
379 for(int i=1;i<processors;i++){
380 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
381 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
383 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
384 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
387 vector<string> nonBlankAccnosFiles;
388 //delete blank accnos files generated with multiple processes
389 for(int i=0;i<processors;i++){
390 if (!(m->isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
391 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
392 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
395 //append accnos files
396 if (nonBlankAccnosFiles.size() != 0) {
397 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
399 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
400 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
401 remove(nonBlankAccnosFiles[h].c_str());
403 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
405 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
406 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
407 m->mothurOutEndLine();
408 }else{ hasAccnos = false; }
410 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
413 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
415 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
417 //delete accnos file if its blank else report to user
418 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
420 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
422 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
423 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
424 m->mothurOutEndLine();
433 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
435 if (pid == 0) { //only one process should output to screen
438 outputNames.push_back(alignFileName);
439 outputNames.push_back(reportFileName);
440 if (hasAccnos) { outputNames.push_back(accnosFileName); }
446 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
447 m->mothurOutEndLine();
448 m->mothurOutEndLine();
452 m->mothurOutEndLine();
453 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
454 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
455 m->mothurOutEndLine();
459 catch(exception& e) {
460 m->errorOut(e, "AlignCommand", "execute");
465 //**********************************************************************************************************************
467 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
469 ofstream alignmentFile;
470 m->openOutputFile(alignFName, alignmentFile);
473 m->openOutputFile(accnosFName, accnosFile);
475 NastReport report(reportFName);
478 m->openInputFile(filename, inFASTA);
480 inFASTA.seekg(filePos->start);
487 if (m->control_pressed) { return 0; }
489 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
491 int origNumBases = candidateSeq->getNumBases();
492 string originalUnaligned = candidateSeq->getUnaligned();
493 int numBasesNeeded = origNumBases * threshold;
495 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
496 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
497 alignment->resize(candidateSeq->getUnaligned().length()+1);
500 Sequence temp = templateDB->findClosestSequence(candidateSeq);
501 Sequence* templateSeq = &temp;
503 float searchScore = templateDB->getSearchScore();
505 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
509 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
510 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
511 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
512 //so this bool tells you if you need to delete it
514 //if there is a possibility that this sequence should be reversed
515 if (candidateSeq->getNumBases() < numBasesNeeded) {
517 string wasBetter = "";
518 //if the user wants you to try the reverse
520 //get reverse compliment
521 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
522 copy->reverseComplement();
525 Sequence temp2 = templateDB->findClosestSequence(copy);
526 Sequence* templateSeq2 = &temp2;
528 searchScore = templateDB->getSearchScore();
530 nast2 = new Nast(alignment, copy, templateSeq2);
532 //check if any better
533 if (copy->getNumBases() > candidateSeq->getNumBases()) {
534 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
535 templateSeq = templateSeq2;
538 needToDeleteCopy = true;
539 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
541 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
547 //create accnos file with names
548 accnosFile << candidateSeq->getName() << wasBetter << endl;
551 report.setCandidate(candidateSeq);
552 report.setTemplate(templateSeq);
553 report.setSearchParameters(search, searchScore);
554 report.setAlignmentParameters(align, alignment);
555 report.setNastParameters(*nast);
557 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
561 if (needToDeleteCopy) { delete copy; }
567 unsigned long int pos = inFASTA.tellg();
568 if ((pos == -1) || (pos >= filePos->end)) { break; }
571 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
575 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
577 alignmentFile.close();
583 catch(exception& e) {
584 m->errorOut(e, "AlignCommand", "driver");
588 //**********************************************************************************************************************
590 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
592 string outputString = "";
593 MPI_Status statusReport;
594 MPI_Status statusAlign;
595 MPI_Status statusAccnos;
598 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
603 outputString = report.getHeaders();
604 int length = outputString.length();
606 char* buf = new char[length];
607 memcpy(buf, outputString.c_str(), length);
609 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
614 for(int i=0;i<num;i++){
616 if (m->control_pressed) { return 0; }
619 int length = MPIPos[start+i+1] - MPIPos[start+i];
621 char* buf4 = new char[length];
622 //memcpy(buf4, outputString.c_str(), length);
624 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
626 string tempBuf = buf4;
630 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
632 istringstream iss (tempBuf,istringstream::in);
634 Sequence* candidateSeq = new Sequence(iss);
636 int origNumBases = candidateSeq->getNumBases();
637 string originalUnaligned = candidateSeq->getUnaligned();
638 int numBasesNeeded = origNumBases * threshold;
640 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
641 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
642 alignment->resize(candidateSeq->getUnaligned().length()+1);
645 Sequence temp = templateDB->findClosestSequence(candidateSeq);
646 Sequence* templateSeq = &temp;
648 float searchScore = templateDB->getSearchScore();
650 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
654 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
655 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
656 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
657 //so this bool tells you if you need to delete it
659 //if there is a possibility that this sequence should be reversed
660 if (candidateSeq->getNumBases() < numBasesNeeded) {
662 string wasBetter = "";
663 //if the user wants you to try the reverse
665 //get reverse compliment
666 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
667 copy->reverseComplement();
670 Sequence temp2 = templateDB->findClosestSequence(copy);
671 Sequence* templateSeq2 = &temp2;
673 searchScore = templateDB->getSearchScore();
675 nast2 = new Nast(alignment, copy, templateSeq2);
677 //check if any better
678 if (copy->getNumBases() > candidateSeq->getNumBases()) {
679 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
680 templateSeq = templateSeq2;
683 needToDeleteCopy = true;
685 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
691 //create accnos file with names
692 outputString = candidateSeq->getName() + wasBetter + "\n";
694 //send results to parent
695 int length = outputString.length();
697 char* buf = new char[length];
698 memcpy(buf, outputString.c_str(), length);
700 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
702 MPIWroteAccnos = true;
705 report.setCandidate(candidateSeq);
706 report.setTemplate(templateSeq);
707 report.setSearchParameters(search, searchScore);
708 report.setAlignmentParameters(align, alignment);
709 report.setNastParameters(*nast);
711 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
713 //send results to parent
714 int length = outputString.length();
715 char* buf2 = new char[length];
716 memcpy(buf2, outputString.c_str(), length);
718 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
722 outputString = report.getReport();
724 //send results to parent
725 length = outputString.length();
726 char* buf3 = new char[length];
727 memcpy(buf3, outputString.c_str(), length);
729 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
733 if (needToDeleteCopy) { delete copy; }
738 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
741 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
745 catch(exception& e) {
746 m->errorOut(e, "AlignCommand", "driverMPI");
751 /**************************************************************************************************/
753 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
755 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
758 // processIDS.resize(0);
760 //loop through and create all the processes you want
761 while (process != processors) {
765 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
768 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
770 //pass numSeqs to parent
772 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
773 m->openOutputFile(tempFile, out);
778 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
781 //force parent to wait until all the processes are done
782 for (int i=0;i<processors;i++) {
783 int temp = processIDS[i];
787 for (int i = 0; i < processIDS.size(); i++) {
789 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
790 m->openInputFile(tempFile, in);
791 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
792 in.close(); remove(tempFile.c_str());
798 catch(exception& e) {
799 m->errorOut(e, "AlignCommand", "createProcesses");
803 /**************************************************************************************************/
805 void AlignCommand::appendAlignFiles(string temp, string filename) {
810 m->openOutputFileAppend(filename, output);
811 m->openInputFile(temp, input);
813 while(char c = input.get()){
814 if(input.eof()) { break; }
815 else { output << c; }
821 catch(exception& e) {
822 m->errorOut(e, "AlignCommand", "appendAlignFiles");
826 //**********************************************************************************************************************
828 void AlignCommand::appendReportFiles(string temp, string filename) {
833 m->openOutputFileAppend(filename, output);
834 m->openInputFile(temp, input);
836 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
838 while(char c = input.get()){
839 if(input.eof()) { break; }
840 else { output << c; }
846 catch(exception& e) {
847 m->errorOut(e, "AlignCommand", "appendReportFiles");
851 //**********************************************************************************************************************