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++) {
93 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
96 string path = m->hasPath(candidateFileNames[i]);
97 //if the user has not given a path then, add inputdir. else leave path alone.
98 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
104 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
106 //if you can't open it, try default location
107 if (ableToOpen == 1) {
108 if (m->getDefaultPath() != "") { //default path is set
109 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
110 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
111 ableToOpen = m->openInputFile(tryPath, in, "noerror");
112 candidateFileNames[i] = tryPath;
117 if (ableToOpen == 1) {
118 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
119 //erase from file list
120 candidateFileNames.erase(candidateFileNames.begin()+i);
126 //make sure there is at least one valid file left
127 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
130 //check for optional parameter and set defaults
131 // ...at some point should added some additional type checking...
133 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
134 convert(temp, kmerSize);
136 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
137 convert(temp, match);
139 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
140 convert(temp, misMatch);
142 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
143 convert(temp, gapOpen);
145 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
146 convert(temp, gapExtend);
148 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
149 convert(temp, processors);
151 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
152 flip = m->isTrue(temp);
154 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
155 convert(temp, threshold);
157 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
159 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
163 catch(exception& e) {
164 m->errorOut(e, "AlignCommand", "AlignCommand");
169 //**********************************************************************************************************************
171 AlignCommand::~AlignCommand(){
173 if (abort == false) {
174 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
180 //**********************************************************************************************************************
182 void AlignCommand::help(){
184 m->mothurOut("The align.seqs command reads a file containing sequences and creates an alignment file and a report file.\n");
185 m->mothurOut("The align.seqs command parameters are template, candidate, search, ksize, align, match, mismatch, gapopen and gapextend.\n");
186 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");
187 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");
188 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");
189 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
190 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
191 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
192 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
193 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
194 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");
195 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");
196 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");
197 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");
198 m->mothurOut("The align.seqs command should be in the following format: \n");
199 m->mothurOut("align.seqs(template=yourTemplateFile, candidate=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty) \n");
200 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");
201 m->mothurOut("Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).\n\n");
203 catch(exception& e) {
204 m->errorOut(e, "AlignCommand", "help");
210 //**********************************************************************************************************************
212 int AlignCommand::execute(){
214 if (abort == true) { return 0; }
216 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
217 int longestBase = templateDB->getLongestBase();
219 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
220 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
221 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
222 else if(align == "noalign") { alignment = new NoAlign(); }
224 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
225 m->mothurOutEndLine();
226 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
228 vector<string> outputNames;
230 for (int s = 0; s < candidateFileNames.size(); s++) {
231 if (m->control_pressed) { return 0; }
233 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
235 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
236 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
237 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
238 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
239 bool hasAccnos = true;
241 int numFastaSeqs = 0;
242 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
243 int start = time(NULL);
246 int pid, end, numSeqsPerProcessor;
248 vector<unsigned long int> MPIPos;
249 MPIWroteAccnos = false;
252 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
253 MPI_Comm_size(MPI_COMM_WORLD, &processors);
256 MPI_File outMPIAlign;
257 MPI_File outMPIReport;
258 MPI_File outMPIAccnos;
260 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
261 int inMode=MPI_MODE_RDONLY;
263 char outAlignFilename[1024];
264 strcpy(outAlignFilename, alignFileName.c_str());
266 char outReportFilename[1024];
267 strcpy(outReportFilename, reportFileName.c_str());
269 char outAccnosFilename[1024];
270 strcpy(outAccnosFilename, accnosFileName.c_str());
272 char inFileName[1024];
273 strcpy(inFileName, candidateFileNames[s].c_str());
275 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
276 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
277 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
278 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
280 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
282 if (pid == 0) { //you are the root process
284 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
286 //send file positions to all processes
287 for(int i = 1; i < processors; i++) {
288 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
289 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
292 //figure out how many sequences you have to align
293 numSeqsPerProcessor = numFastaSeqs / processors;
294 int startIndex = pid * numSeqsPerProcessor;
295 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
298 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
300 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
302 for (int i = 1; i < processors; i++) {
304 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
305 if (tempResult != 0) { MPIWroteAccnos = true; }
307 }else{ //you are a child process
308 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
309 MPIPos.resize(numFastaSeqs+1);
310 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
313 //figure out how many sequences you have to align
314 numSeqsPerProcessor = numFastaSeqs / processors;
315 int startIndex = pid * numSeqsPerProcessor;
316 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
320 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
322 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); return 0; }
324 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
328 MPI_File_close(&inMPI);
329 MPI_File_close(&outMPIAlign);
330 MPI_File_close(&outMPIReport);
331 MPI_File_close(&outMPIAccnos);
333 //delete accnos file if blank
335 //delete accnos file if its blank else report to user
336 if (MPIWroteAccnos) {
337 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
339 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
340 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
341 m->mothurOutEndLine();
344 //MPI_File_delete(outAccnosFilename, info);
346 remove(accnosFileName.c_str());
352 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
353 for (int i = 0; i < (positions.size()-1); i++) {
354 lines.push_back(new linePair(positions[i], positions[(i+1)]));
356 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
358 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
360 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
362 //delete accnos file if its blank else report to user
363 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
365 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
367 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
368 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
369 m->mothurOutEndLine();
372 processIDS.resize(0);
374 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
376 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
377 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
379 //append alignment and report files
380 for(int i=1;i<processors;i++){
381 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
382 remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
384 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
385 remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
388 vector<string> nonBlankAccnosFiles;
389 //delete blank accnos files generated with multiple processes
390 for(int i=0;i<processors;i++){
391 if (!(m->isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
392 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
393 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
396 //append accnos files
397 if (nonBlankAccnosFiles.size() != 0) {
398 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
400 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
401 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
402 remove(nonBlankAccnosFiles[h].c_str());
404 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
406 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
407 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
408 m->mothurOutEndLine();
409 }else{ hasAccnos = false; }
411 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
414 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
416 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
418 //delete accnos file if its blank else report to user
419 if (m->isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
421 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
423 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
424 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
425 m->mothurOutEndLine();
434 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
436 if (pid == 0) { //only one process should output to screen
439 outputNames.push_back(alignFileName);
440 outputNames.push_back(reportFileName);
441 if (hasAccnos) { outputNames.push_back(accnosFileName); }
447 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
448 m->mothurOutEndLine();
449 m->mothurOutEndLine();
453 m->mothurOutEndLine();
454 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
455 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
456 m->mothurOutEndLine();
460 catch(exception& e) {
461 m->errorOut(e, "AlignCommand", "execute");
466 //**********************************************************************************************************************
468 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
470 ofstream alignmentFile;
471 m->openOutputFile(alignFName, alignmentFile);
474 m->openOutputFile(accnosFName, accnosFile);
476 NastReport report(reportFName);
479 m->openInputFile(filename, inFASTA);
481 inFASTA.seekg(filePos->start);
488 if (m->control_pressed) { return 0; }
490 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
492 int origNumBases = candidateSeq->getNumBases();
493 string originalUnaligned = candidateSeq->getUnaligned();
494 int numBasesNeeded = origNumBases * threshold;
496 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
497 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
498 alignment->resize(candidateSeq->getUnaligned().length()+1);
501 Sequence temp = templateDB->findClosestSequence(candidateSeq);
502 Sequence* templateSeq = &temp;
504 float searchScore = templateDB->getSearchScore();
506 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
510 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
511 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
512 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
513 //so this bool tells you if you need to delete it
515 //if there is a possibility that this sequence should be reversed
516 if (candidateSeq->getNumBases() < numBasesNeeded) {
518 string wasBetter = "";
519 //if the user wants you to try the reverse
521 //get reverse compliment
522 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
523 copy->reverseComplement();
526 Sequence temp2 = templateDB->findClosestSequence(copy);
527 Sequence* templateSeq2 = &temp2;
529 searchScore = templateDB->getSearchScore();
531 nast2 = new Nast(alignment, copy, templateSeq2);
533 //check if any better
534 if (copy->getNumBases() > candidateSeq->getNumBases()) {
535 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
536 templateSeq = templateSeq2;
539 needToDeleteCopy = true;
540 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
542 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
548 //create accnos file with names
549 accnosFile << candidateSeq->getName() << wasBetter << endl;
552 report.setCandidate(candidateSeq);
553 report.setTemplate(templateSeq);
554 report.setSearchParameters(search, searchScore);
555 report.setAlignmentParameters(align, alignment);
556 report.setNastParameters(*nast);
558 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
562 if (needToDeleteCopy) { delete copy; }
568 unsigned long int pos = inFASTA.tellg();
569 if ((pos == -1) || (pos >= filePos->end)) { break; }
572 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
576 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
578 alignmentFile.close();
584 catch(exception& e) {
585 m->errorOut(e, "AlignCommand", "driver");
589 //**********************************************************************************************************************
591 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
593 string outputString = "";
594 MPI_Status statusReport;
595 MPI_Status statusAlign;
596 MPI_Status statusAccnos;
599 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
604 outputString = report.getHeaders();
605 int length = outputString.length();
607 char* buf = new char[length];
608 memcpy(buf, outputString.c_str(), length);
610 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
615 for(int i=0;i<num;i++){
617 if (m->control_pressed) { return 0; }
620 int length = MPIPos[start+i+1] - MPIPos[start+i];
622 char* buf4 = new char[length];
623 //memcpy(buf4, outputString.c_str(), length);
625 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
627 string tempBuf = buf4;
631 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
633 istringstream iss (tempBuf,istringstream::in);
635 Sequence* candidateSeq = new Sequence(iss);
637 int origNumBases = candidateSeq->getNumBases();
638 string originalUnaligned = candidateSeq->getUnaligned();
639 int numBasesNeeded = origNumBases * threshold;
641 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
642 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
643 alignment->resize(candidateSeq->getUnaligned().length()+1);
646 Sequence temp = templateDB->findClosestSequence(candidateSeq);
647 Sequence* templateSeq = &temp;
649 float searchScore = templateDB->getSearchScore();
651 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
655 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
656 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
657 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
658 //so this bool tells you if you need to delete it
660 //if there is a possibility that this sequence should be reversed
661 if (candidateSeq->getNumBases() < numBasesNeeded) {
663 string wasBetter = "";
664 //if the user wants you to try the reverse
666 //get reverse compliment
667 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
668 copy->reverseComplement();
671 Sequence temp2 = templateDB->findClosestSequence(copy);
672 Sequence* templateSeq2 = &temp2;
674 searchScore = templateDB->getSearchScore();
676 nast2 = new Nast(alignment, copy, templateSeq2);
678 //check if any better
679 if (copy->getNumBases() > candidateSeq->getNumBases()) {
680 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
681 templateSeq = templateSeq2;
684 needToDeleteCopy = true;
686 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
692 //create accnos file with names
693 outputString = candidateSeq->getName() + wasBetter + "\n";
695 //send results to parent
696 int length = outputString.length();
698 char* buf = new char[length];
699 memcpy(buf, outputString.c_str(), length);
701 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
703 MPIWroteAccnos = true;
706 report.setCandidate(candidateSeq);
707 report.setTemplate(templateSeq);
708 report.setSearchParameters(search, searchScore);
709 report.setAlignmentParameters(align, alignment);
710 report.setNastParameters(*nast);
712 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
714 //send results to parent
715 int length = outputString.length();
716 char* buf2 = new char[length];
717 memcpy(buf2, outputString.c_str(), length);
719 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
723 outputString = report.getReport();
725 //send results to parent
726 length = outputString.length();
727 char* buf3 = new char[length];
728 memcpy(buf3, outputString.c_str(), length);
730 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
734 if (needToDeleteCopy) { delete copy; }
739 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
742 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
746 catch(exception& e) {
747 m->errorOut(e, "AlignCommand", "driverMPI");
752 /**************************************************************************************************/
754 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
756 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
759 // processIDS.resize(0);
761 //loop through and create all the processes you want
762 while (process != processors) {
766 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
769 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
771 //pass numSeqs to parent
773 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
774 m->openOutputFile(tempFile, out);
779 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
782 //force parent to wait until all the processes are done
783 for (int i=0;i<processors;i++) {
784 int temp = processIDS[i];
788 for (int i = 0; i < processIDS.size(); i++) {
790 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
791 m->openInputFile(tempFile, in);
792 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
793 in.close(); remove(tempFile.c_str());
799 catch(exception& e) {
800 m->errorOut(e, "AlignCommand", "createProcesses");
804 /**************************************************************************************************/
806 void AlignCommand::appendAlignFiles(string temp, string filename) {
811 m->openOutputFileAppend(filename, output);
812 m->openInputFile(temp, input);
814 while(char c = input.get()){
815 if(input.eof()) { break; }
816 else { output << c; }
822 catch(exception& e) {
823 m->errorOut(e, "AlignCommand", "appendAlignFiles");
827 //**********************************************************************************************************************
829 void AlignCommand::appendReportFiles(string temp, string filename) {
834 m->openOutputFileAppend(filename, output);
835 m->openInputFile(temp, input);
837 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
839 while(char c = input.get()){
840 if(input.eof()) { break; }
841 else { output << c; }
847 catch(exception& e) {
848 m->errorOut(e, "AlignCommand", "appendReportFiles");
852 //**********************************************************************************************************************