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;
49 map<string, string>::iterator it;
51 //check to make sure all parameters are valid for command
52 for (it = parameters.begin(); it != parameters.end(); it++) {
53 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
56 //if the user changes the output directory command factory will send this info to us in the output parameter
57 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
60 //if the user changes the input directory command factory will send this info to us in the output parameter
61 string inputDir = validParameter.validFile(parameters, "inputdir", false);
63 if (inputDir == "not found"){ inputDir = ""; }
67 it = parameters.find("template");
69 //user has given a template file
70 if(it != parameters.end()){
71 path = hasPath(it->second);
72 //if the user has not given a path then, add inputdir. else leave path alone.
73 if (path == "") { parameters["template"] = inputDir + it->second; }
77 //check for required parameters
78 templateFileName = validParameter.validFile(parameters, "template", true);
80 if (templateFileName == "not found") {
81 m->mothurOut("template is a required parameter for the align.seqs command.");
82 m->mothurOutEndLine();
84 }else if (templateFileName == "not open") { abort = true; }
86 candidateFileName = validParameter.validFile(parameters, "candidate", false);
87 if (candidateFileName == "not found") { m->mothurOut("candidate is a required parameter for the align.seqs command."); m->mothurOutEndLine(); abort = true; }
89 splitAtDash(candidateFileName, candidateFileNames);
91 //go through files and make sure they are good, if not, then disregard them
92 for (int i = 0; i < candidateFileNames.size(); i++) {
94 string path = hasPath(candidateFileNames[i]);
95 //if the user has not given a path then, add inputdir. else leave path alone.
96 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
104 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
105 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
110 ableToOpen = openInputFile(candidateFileNames[i], in);
114 for (int j = 1; j < processors; j++) {
115 MPI_Send(&ableToOpen, 1, MPI_INT, j, 2001, MPI_COMM_WORLD);
119 MPI_Recv(&ableToOpen, 1, MPI_INT, 0, 2001, MPI_COMM_WORLD, &status);
124 if (ableToOpen == 1) {
125 m->mothurOut(candidateFileNames[i] + " will be disregarded."); m->mothurOutEndLine();
126 //erase from file list
127 candidateFileNames.erase(candidateFileNames.begin()+i);
133 //make sure there is at least one valid file left
134 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
137 //check for optional parameter and set defaults
138 // ...at some point should added some additional type checking...
140 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
141 convert(temp, kmerSize);
143 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
144 convert(temp, match);
146 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
147 convert(temp, misMatch);
149 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
150 convert(temp, gapOpen);
152 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
153 convert(temp, gapExtend);
155 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
156 convert(temp, processors);
158 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
161 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
162 convert(temp, threshold);
164 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
166 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
170 catch(exception& e) {
171 m->errorOut(e, "AlignCommand", "AlignCommand");
176 //**********************************************************************************************************************
178 AlignCommand::~AlignCommand(){
180 if (abort == false) {
181 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
187 //**********************************************************************************************************************
189 void AlignCommand::help(){
191 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
192 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen and gapextend.\n");
193 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");
194 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");
195 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");
196 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
197 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
198 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
199 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
200 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
201 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");
202 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");
203 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");
204 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");
205 m->mothurOut("The align.seqs command should be in the following format: \n");
206 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
207 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");
208 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
210 catch(exception& e) {
211 m->errorOut(e, "AlignCommand", "help");
217 //**********************************************************************************************************************
219 int AlignCommand::execute(){
221 if (abort == true) { return 0; }
223 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
224 int longestBase = templateDB->getLongestBase();
226 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
227 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
228 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
229 else if(align == "noalign") { alignment = new NoAlign(); }
231 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
232 m->mothurOutEndLine();
233 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
235 vector<string> outputNames;
237 for (int s = 0; s < candidateFileNames.size(); s++) {
238 if (m->control_pressed) { return 0; }
240 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
242 if (outputDir == "") { outputDir += hasPath(candidateFileNames[s]); }
243 string alignFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align";
244 string reportFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align.report";
245 string accnosFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "flip.accnos";
246 bool hasAccnos = true;
248 int numFastaSeqs = 0;
249 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
250 int start = time(NULL);
253 int pid, end, numSeqsPerProcessor;
256 MPIWroteAccnos = false;
259 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
260 MPI_Comm_size(MPI_COMM_WORLD, &processors);
263 MPI_File outMPIAlign;
264 MPI_File outMPIReport;
265 MPI_File outMPIAccnos;
267 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
268 int inMode=MPI_MODE_RDONLY;
270 //char* outAlignFilename = new char[alignFileName.length()];
271 //memcpy(outAlignFilename, alignFileName.c_str(), alignFileName.length());
273 char outAlignFilename[1024];
274 strcpy(outAlignFilename, alignFileName.c_str());
276 //char* outReportFilename = new char[reportFileName.length()];
277 //memcpy(outReportFilename, reportFileName.c_str(), reportFileName.length());
279 char outReportFilename[1024];
280 strcpy(outReportFilename, reportFileName.c_str());
282 //char* outAccnosFilename = new char[accnosFileName.length()];
283 //memcpy(outAccnosFilename, accnosFileName.c_str(), accnosFileName.length());
285 char outAccnosFilename[1024];
286 strcpy(outAccnosFilename, accnosFileName.c_str());
288 //char* inFileName = new char[candidateFileNames[s].length()];
289 //memcpy(inFileName, candidateFileNames[s].c_str(), candidateFileNames[s].length());
291 char inFileName[1024];
292 strcpy(inFileName, candidateFileNames[s].c_str());
294 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
295 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
296 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
297 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
299 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
301 if (pid == 0) { //you are the root process
303 MPIPos = setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
305 //send file positions to all processes
306 for(int i = 1; i < processors; i++) {
307 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
308 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
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 for (int i = 1; i < processors; i++) {
324 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
325 if (tempResult != 0) { MPIWroteAccnos = true; }
327 }else{ //you are a child process
328 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
329 MPIPos.resize(numFastaSeqs+1);
330 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
333 //figure out how many sequences you have to align
334 numSeqsPerProcessor = numFastaSeqs / processors;
335 int startIndex = pid * numSeqsPerProcessor;
336 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
340 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
342 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
344 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
348 MPI_File_close(&inMPI);
349 MPI_File_close(&outMPIAlign);
350 MPI_File_close(&outMPIReport);
351 MPI_File_close(&outMPIAccnos);
353 //delete accnos file if blank
355 //delete accnos file if its blank else report to user
356 if (MPIWroteAccnos) {
357 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
359 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
360 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
361 m->mothurOutEndLine();
364 //MPI_File_delete(outAccnosFilename, info);
366 remove(accnosFileName.c_str());
372 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
375 openInputFile(candidateFileNames[s], inFASTA);
376 getNumSeqs(inFASTA, numFastaSeqs);
379 lines.push_back(new linePair(0, numFastaSeqs));
381 driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
383 if (m->control_pressed) {
384 remove(accnosFileName.c_str());
385 remove(alignFileName.c_str());
386 remove(reportFileName.c_str());
390 //delete accnos file if its blank else report to user
391 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
393 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
395 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
396 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
397 m->mothurOutEndLine();
401 vector<int> positions;
402 processIDS.resize(0);
405 openInputFile(candidateFileNames[s], inFASTA);
408 while(!inFASTA.eof()){
409 input = getline(inFASTA);
410 if (input.length() != 0) {
411 if(input[0] == '>'){ long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
416 numFastaSeqs = positions.size();
418 int numSeqsPerProcessor = numFastaSeqs / processors;
420 for (int i = 0; i < processors; i++) {
421 long int startPos = positions[ i * numSeqsPerProcessor ];
422 if(i == processors - 1){
423 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
425 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
428 createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
430 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
431 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
433 //append alignment and report files
434 for(int i=1;i<processors;i++){
435 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
436 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
438 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
439 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
442 vector<string> nonBlankAccnosFiles;
443 //delete blank accnos files generated with multiple processes
444 for(int i=0;i<processors;i++){
445 if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
446 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
447 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
450 //append accnos files
451 if (nonBlankAccnosFiles.size() != 0) {
452 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
454 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
455 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
456 remove(nonBlankAccnosFiles[h].c_str());
458 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
460 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
461 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
462 m->mothurOutEndLine();
463 }else{ hasAccnos = false; }
465 if (m->control_pressed) {
466 remove(accnosFileName.c_str());
467 remove(alignFileName.c_str());
468 remove(reportFileName.c_str());
474 openInputFile(candidateFileNames[s], inFASTA);
475 getNumSeqs(inFASTA, numFastaSeqs);
478 lines.push_back(new linePair(0, numFastaSeqs));
480 driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
482 if (m->control_pressed) {
483 remove(accnosFileName.c_str());
484 remove(alignFileName.c_str());
485 remove(reportFileName.c_str());
489 //delete accnos file if its blank else report to user
490 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
492 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
494 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
495 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
496 m->mothurOutEndLine();
505 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
507 if (pid == 0) { //only one process should output to screen
510 outputNames.push_back(alignFileName);
511 outputNames.push_back(reportFileName);
512 if (hasAccnos) { outputNames.push_back(accnosFileName); }
518 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
519 m->mothurOutEndLine();
520 m->mothurOutEndLine();
524 m->mothurOutEndLine();
525 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
526 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
527 m->mothurOutEndLine();
531 catch(exception& e) {
532 m->errorOut(e, "AlignCommand", "execute");
537 //**********************************************************************************************************************
539 int AlignCommand::driver(linePair* line, string alignFName, string reportFName, string accnosFName, string filename){
541 ofstream alignmentFile;
542 openOutputFile(alignFName, alignmentFile);
545 openOutputFile(accnosFName, accnosFile);
547 NastReport report(reportFName);
550 openInputFile(filename, inFASTA);
552 inFASTA.seekg(line->start);
554 for(int i=0;i<line->numSeqs;i++){
556 if (m->control_pressed) { return 0; }
558 Sequence* candidateSeq = new Sequence(inFASTA); gobble(inFASTA);
560 int origNumBases = candidateSeq->getNumBases();
561 string originalUnaligned = candidateSeq->getUnaligned();
562 int numBasesNeeded = origNumBases * threshold;
564 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
565 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
566 alignment->resize(candidateSeq->getUnaligned().length()+1);
569 Sequence temp = templateDB->findClosestSequence(candidateSeq);
570 Sequence* templateSeq = &temp;
572 float searchScore = templateDB->getSearchScore();
574 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
578 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
579 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
580 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
581 //so this bool tells you if you need to delete it
583 //if there is a possibility that this sequence should be reversed
584 if (candidateSeq->getNumBases() < numBasesNeeded) {
586 string wasBetter = "";
587 //if the user wants you to try the reverse
589 //get reverse compliment
590 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
591 copy->reverseComplement();
594 Sequence temp2 = templateDB->findClosestSequence(copy);
595 Sequence* templateSeq2 = &temp2;
597 searchScore = templateDB->getSearchScore();
599 nast2 = new Nast(alignment, copy, templateSeq2);
601 //check if any better
602 if (copy->getNumBases() > candidateSeq->getNumBases()) {
603 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
604 templateSeq = templateSeq2;
607 needToDeleteCopy = true;
609 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
615 //create accnos file with names
616 accnosFile << candidateSeq->getName() << wasBetter << endl;
619 report.setCandidate(candidateSeq);
620 report.setTemplate(templateSeq);
621 report.setSearchParameters(search, searchScore);
622 report.setAlignmentParameters(align, alignment);
623 report.setNastParameters(*nast);
625 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
629 if (needToDeleteCopy) { delete copy; }
634 if((i+1) % 100 == 0){ m->mothurOut(toString(i+1)); m->mothurOutEndLine(); }
637 if((line->numSeqs) % 100 != 0){ m->mothurOut(toString(line->numSeqs)); m->mothurOutEndLine(); }
639 alignmentFile.close();
645 catch(exception& e) {
646 m->errorOut(e, "AlignCommand", "driver");
650 //**********************************************************************************************************************
652 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<long>& MPIPos){
654 string outputString = "";
655 MPI_Status statusReport;
656 MPI_Status statusAlign;
657 MPI_Status statusAccnos;
660 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
665 outputString = report.getHeaders();
666 int length = outputString.length();
668 char* buf = new char[length];
669 memcpy(buf, outputString.c_str(), length);
671 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
676 for(int i=0;i<num;i++){
678 if (m->control_pressed) { return 0; }
681 int length = MPIPos[start+i+1] - MPIPos[start+i];
683 char* buf4 = new char[length];
684 memcpy(buf4, outputString.c_str(), length);
686 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
688 string tempBuf = buf4;
692 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
693 istringstream iss (tempBuf,istringstream::in);
695 Sequence* candidateSeq = new Sequence(iss);
696 int origNumBases = candidateSeq->getNumBases();
697 string originalUnaligned = candidateSeq->getUnaligned();
698 int numBasesNeeded = origNumBases * threshold;
700 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
701 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
702 alignment->resize(candidateSeq->getUnaligned().length()+1);
705 Sequence temp = templateDB->findClosestSequence(candidateSeq);
706 Sequence* templateSeq = &temp;
708 float searchScore = templateDB->getSearchScore();
710 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
714 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
715 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
716 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
717 //so this bool tells you if you need to delete it
719 //if there is a possibility that this sequence should be reversed
720 if (candidateSeq->getNumBases() < numBasesNeeded) {
722 string wasBetter = "";
723 //if the user wants you to try the reverse
725 //get reverse compliment
726 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
727 copy->reverseComplement();
730 Sequence temp2 = templateDB->findClosestSequence(copy);
731 Sequence* templateSeq2 = &temp2;
733 searchScore = templateDB->getSearchScore();
735 nast2 = new Nast(alignment, copy, templateSeq2);
737 //check if any better
738 if (copy->getNumBases() > candidateSeq->getNumBases()) {
739 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
740 templateSeq = templateSeq2;
743 needToDeleteCopy = true;
745 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
751 //create accnos file with names
752 outputString = candidateSeq->getName() + wasBetter + "\n";
754 //send results to parent
755 int length = outputString.length();
757 char* buf = new char[length];
758 memcpy(buf, outputString.c_str(), length);
760 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
762 MPIWroteAccnos = true;
765 report.setCandidate(candidateSeq);
766 report.setTemplate(templateSeq);
767 report.setSearchParameters(search, searchScore);
768 report.setAlignmentParameters(align, alignment);
769 report.setNastParameters(*nast);
771 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
773 //send results to parent
774 int length = outputString.length();
775 char* buf2 = new char[length];
776 memcpy(buf2, outputString.c_str(), length);
778 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
782 outputString = report.getReport();
784 //send results to parent
785 length = outputString.length();
786 char* buf3 = new char[length];
787 memcpy(buf3, outputString.c_str(), length);
789 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
793 if (needToDeleteCopy) { delete copy; }
798 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
801 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
805 catch(exception& e) {
806 m->errorOut(e, "AlignCommand", "driverMPI");
811 /**************************************************************************************************/
813 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
815 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
818 // processIDS.resize(0);
820 //loop through and create all the processes you want
821 while (process != processors) {
825 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
828 exitCommand = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
830 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
833 //force parent to wait until all the processes are done
834 for (int i=0;i<processors;i++) {
835 int temp = processIDS[i];
842 catch(exception& e) {
843 m->errorOut(e, "AlignCommand", "createProcesses");
848 /**************************************************************************************************/
850 void AlignCommand::appendAlignFiles(string temp, string filename) {
855 openOutputFileAppend(filename, output);
856 openInputFile(temp, input);
858 while(char c = input.get()){
859 if(input.eof()) { break; }
860 else { output << c; }
866 catch(exception& e) {
867 m->errorOut(e, "AlignCommand", "appendAlignFiles");
871 //**********************************************************************************************************************
873 void AlignCommand::appendReportFiles(string temp, string filename) {
878 openOutputFileAppend(filename, output);
879 openInputFile(temp, input);
881 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
883 while(char c = input.get()){
884 if(input.eof()) { break; }
885 else { output << c; }
891 catch(exception& e) {
892 m->errorOut(e, "AlignCommand", "appendReportFiles");
896 //**********************************************************************************************************************