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 "referencedb.h"
20 //**********************************************************************************************************************
21 vector<string> AlignCommand::setParameters(){
23 CommandParameter ptemplate("reference", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(ptemplate);
24 CommandParameter pcandidate("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pcandidate);
25 CommandParameter psearch("search", "Multiple", "kmer-blast-suffix", "kmer", "", "", "",false,false); parameters.push_back(psearch);
26 CommandParameter pksize("ksize", "Number", "", "8", "", "", "",false,false); parameters.push_back(pksize);
27 CommandParameter pmatch("match", "Number", "", "1.0", "", "", "",false,false); parameters.push_back(pmatch);
28 CommandParameter palign("align", "Multiple", "needleman-gotoh-blast-noalign", "needleman", "", "", "",false,false); parameters.push_back(palign);
29 CommandParameter pmismatch("mismatch", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pmismatch);
30 CommandParameter pgapopen("gapopen", "Number", "", "-2.0", "", "", "",false,false); parameters.push_back(pgapopen);
31 CommandParameter pgapextend("gapextend", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pgapextend);
32 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
33 CommandParameter pflip("flip", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(pflip);
34 CommandParameter psave("save", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(psave);
35 CommandParameter pthreshold("threshold", "Number", "", "0.50", "", "", "",false,false); parameters.push_back(pthreshold);
36 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
37 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
39 vector<string> myArray;
40 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
44 m->errorOut(e, "AlignCommand", "setParameters");
48 //**********************************************************************************************************************
49 string AlignCommand::getHelpString(){
51 string helpString = "";
52 helpString += "The align.seqs command reads a file containing sequences and creates an alignment file and a report file.";
53 helpString += "The align.seqs command parameters are reference, fasta, search, ksize, align, match, mismatch, gapopen, gapextend and processors.";
54 helpString += "The reference and fasta parameters are required. You may leave fasta blank if you have a valid fasta file. You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amzon.fasta.";
55 helpString += "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.";
56 helpString += "The align parameter allows you to specify the alignment method to use. Your options are: gotoh, needleman, blast and noalign. The default is needleman.";
57 helpString += "The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.";
58 helpString += "The match parameter allows you to specify the bonus for having the same base. The default is 1.0.";
59 helpString += "The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.";
60 helpString += "The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.";
61 helpString += "The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.";
62 helpString += "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.";
63 helpString += "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.";
64 helpString += "If the flip parameter is set to true the reverse complement of the sequence is aligned and the better alignment is reported.";
65 helpString += "If the save parameter is set to true the reference sequences will be saved in memory, to clear them later you can use the clear.memory command. Default=f.";
66 helpString += "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.";
67 helpString += "The align.seqs command should be in the following format:";
68 helpString += "align.seqs(reference=yourTemplateFile, fasta=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty)";
69 helpString += "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)";
70 helpString += "Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).";
74 m->errorOut(e, "AlignCommand", "getHelpString");
78 //**********************************************************************************************************************
79 AlignCommand::AlignCommand(){
81 abort = true; calledHelp = true;
83 vector<string> tempOutNames;
84 outputTypes["fasta"] = tempOutNames;
85 outputTypes["alignreport"] = tempOutNames;
86 outputTypes["accnos"] = tempOutNames;
89 m->errorOut(e, "AlignCommand", "AlignCommand");
93 //**********************************************************************************************************************
94 AlignCommand::AlignCommand(string option) {
96 abort = false; calledHelp = false;
97 ReferenceDB* rdb = ReferenceDB::getInstance();
99 //allow user to run help
100 if(option == "help") { help(); abort = true; calledHelp = true;}
101 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
104 vector<string> myArray = setParameters();
106 OptionParser parser(option);
107 map<string, string> parameters = parser.getParameters();
109 ValidParameters validParameter("align.seqs");
110 map<string, string>::iterator it;
112 //check to make sure all parameters are valid for command
113 for (it = parameters.begin(); it != parameters.end(); it++) {
114 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
117 //initialize outputTypes
118 vector<string> tempOutNames;
119 outputTypes["fasta"] = tempOutNames;
120 outputTypes["alignreport"] = tempOutNames;
121 outputTypes["accnos"] = tempOutNames;
123 //if the user changes the output directory command factory will send this info to us in the output parameter
124 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
127 //if the user changes the input directory command factory will send this info to us in the output parameter
128 string inputDir = validParameter.validFile(parameters, "inputdir", false);
130 if (inputDir == "not found"){ inputDir = ""; }
134 it = parameters.find("reference");
136 //user has given a template file
137 if(it != parameters.end()){
138 path = m->hasPath(it->second);
139 //if the user has not given a path then, add inputdir. else leave path alone.
140 if (path == "") { parameters["reference"] = inputDir + it->second; }
144 candidateFileName = validParameter.validFile(parameters, "fasta", false);
145 if (candidateFileName == "not found") {
146 //if there is a current fasta file, use it
147 string filename = m->getFastaFile();
148 if (filename != "") { candidateFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
149 else { m->mothurOut("You have no current fastafile and the candidate parameter is required."); m->mothurOutEndLine(); abort = true; }
151 m->splitAtDash(candidateFileName, candidateFileNames);
153 //go through files and make sure they are good, if not, then disregard them
154 for (int i = 0; i < candidateFileNames.size(); i++) {
155 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
158 if (candidateFileNames[i] == "current") {
159 candidateFileNames[i] = m->getFastaFile();
160 if (candidateFileNames[i] != "") { m->mothurOut("Using " + candidateFileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
162 m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true;
163 //erase from file list
164 candidateFileNames.erase(candidateFileNames.begin()+i);
171 if (inputDir != "") {
172 string path = m->hasPath(candidateFileNames[i]);
173 //if the user has not given a path then, add inputdir. else leave path alone.
174 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
179 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
182 //if you can't open it, try default location
183 if (ableToOpen == 1) {
184 if (m->getDefaultPath() != "") { //default path is set
185 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
186 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
188 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
190 candidateFileNames[i] = tryPath;
194 //if you can't open it, try output location
195 if (ableToOpen == 1) {
196 if (m->getOutputDir() != "") { //default path is set
197 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
198 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
200 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
202 candidateFileNames[i] = tryPath;
208 if (ableToOpen == 1) {
209 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
210 //erase from file list
211 candidateFileNames.erase(candidateFileNames.begin()+i);
214 m->setFastaFile(candidateFileNames[i]);
219 //make sure there is at least one valid file left
220 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
223 //check for optional parameter and set defaults
224 // ...at some point should added some additional type checking...
226 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
227 convert(temp, kmerSize);
229 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
230 convert(temp, match);
232 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
233 convert(temp, misMatch);
235 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
236 convert(temp, gapOpen);
238 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
239 convert(temp, gapExtend);
241 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
242 m->setProcessors(temp);
243 convert(temp, processors);
245 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
246 flip = m->isTrue(temp);
248 temp = validParameter.validFile(parameters, "save", false); if (temp == "not found"){ temp = "f"; }
249 save = m->isTrue(temp);
250 //this is so the threads can quickly load the reference data
251 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
253 if (processors != 1) { save = true; }
256 if (save) { //clear out old references
260 //this has to go after save so that if the user sets save=t and provides no reference we abort
261 templateFileName = validParameter.validFile(parameters, "reference", true);
262 if (templateFileName == "not found") {
263 //check for saved reference sequences
264 if (rdb->referenceSeqs.size() != 0) {
265 templateFileName = "saved";
267 m->mothurOut("[ERROR]: You don't have any saved reference sequences and the reference parameter is a required for the align.seqs command.");
268 m->mothurOutEndLine();
271 }else if (templateFileName == "not open") { abort = true; }
272 else { if (save) { rdb->setSavedReference(templateFileName); } }
274 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
275 convert(temp, threshold);
277 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
278 if ((search != "suffix") && (search != "kmer") && (search != "blast")) { m->mothurOut("invalid search option: choices are kmer, suffix or blast."); m->mothurOutEndLine(); abort=true; }
280 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
281 if ((align != "needleman") && (align != "gotoh") && (align != "blast") && (align != "noalign")) { m->mothurOut("invalid align option: choices are needleman, gotoh, blast or noalign."); m->mothurOutEndLine(); abort=true; }
286 catch(exception& e) {
287 m->errorOut(e, "AlignCommand", "AlignCommand");
291 //**********************************************************************************************************************
292 AlignCommand::~AlignCommand(){
294 if (abort == false) {
295 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
299 //**********************************************************************************************************************
301 int AlignCommand::execute(){
303 if (abort == true) { if (calledHelp) { return 0; } return 2; }
305 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch, rand());
307 for (int s = 0; s < candidateFileNames.size(); s++) {
308 if (m->control_pressed) { outputTypes.clear(); return 0; }
310 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
312 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
313 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
314 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
315 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
316 bool hasAccnos = true;
318 int numFastaSeqs = 0;
319 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
320 int start = time(NULL);
323 int pid, numSeqsPerProcessor;
325 vector<unsigned long int> MPIPos;
326 MPIWroteAccnos = false;
329 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
330 MPI_Comm_size(MPI_COMM_WORLD, &processors);
333 MPI_File outMPIAlign;
334 MPI_File outMPIReport;
335 MPI_File outMPIAccnos;
337 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
338 int inMode=MPI_MODE_RDONLY;
340 char outAlignFilename[1024];
341 strcpy(outAlignFilename, alignFileName.c_str());
343 char outReportFilename[1024];
344 strcpy(outReportFilename, reportFileName.c_str());
346 char outAccnosFilename[1024];
347 strcpy(outAccnosFilename, accnosFileName.c_str());
349 char inFileName[1024];
350 strcpy(inFileName, candidateFileNames[s].c_str());
352 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
353 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
354 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
355 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
357 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
359 if (pid == 0) { //you are the root process
361 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
363 //send file positions to all processes
364 for(int i = 1; i < processors; i++) {
365 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
366 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
369 //figure out how many sequences you have to align
370 numSeqsPerProcessor = numFastaSeqs / processors;
371 int startIndex = pid * numSeqsPerProcessor;
372 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
375 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
377 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
379 for (int i = 1; i < processors; i++) {
381 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
382 if (tempResult != 0) { MPIWroteAccnos = true; }
384 }else{ //you are a child process
385 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
386 MPIPos.resize(numFastaSeqs+1);
387 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
390 //figure out how many sequences you have to align
391 numSeqsPerProcessor = numFastaSeqs / processors;
392 int startIndex = pid * numSeqsPerProcessor;
393 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
397 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
399 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
401 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
405 MPI_File_close(&inMPI);
406 MPI_File_close(&outMPIAlign);
407 MPI_File_close(&outMPIReport);
408 MPI_File_close(&outMPIAccnos);
410 //delete accnos file if blank
412 //delete accnos file if its blank else report to user
413 if (MPIWroteAccnos) {
414 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
416 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
417 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
418 m->mothurOutEndLine();
421 //MPI_File_delete(outAccnosFilename, info);
423 m->mothurRemove(accnosFileName);
429 vector<unsigned long int> positions;
430 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
431 positions = m->divideFile(candidateFileNames[s], processors);
432 for (int i = 0; i < (positions.size()-1); i++) { lines.push_back(new linePair(positions[i], positions[(i+1)])); }
434 if (processors == 1) {
435 lines.push_back(new linePair(0, 1000));
437 positions = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs);
439 //figure out how many sequences you have to process
440 int numSeqsPerProcessor = numFastaSeqs / processors;
441 for (int i = 0; i < processors; i++) {
442 int startIndex = i * numSeqsPerProcessor;
443 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
444 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
450 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
452 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
455 if (m->control_pressed) { m->mothurRemove(accnosFileName); m->mothurRemove(alignFileName); m->mothurRemove(reportFileName); outputTypes.clear(); return 0; }
457 //delete accnos file if its blank else report to user
458 if (m->isBlank(accnosFileName)) { m->mothurRemove(accnosFileName); hasAccnos = false; }
460 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
462 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
463 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
464 m->mothurOutEndLine();
471 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
473 if (pid == 0) { //only one process should output to screen
476 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
477 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
478 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
484 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
485 m->mothurOutEndLine();
486 m->mothurOutEndLine();
489 //set align file as new current fastafile
490 string currentFasta = "";
491 itTypes = outputTypes.find("fasta");
492 if (itTypes != outputTypes.end()) {
493 if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; m->setFastaFile(currentFasta); }
496 m->mothurOutEndLine();
497 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
498 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
499 m->mothurOutEndLine();
503 catch(exception& e) {
504 m->errorOut(e, "AlignCommand", "execute");
509 //**********************************************************************************************************************
510 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
512 ofstream alignmentFile;
513 m->openOutputFile(alignFName, alignmentFile);
516 m->openOutputFile(accnosFName, accnosFile);
518 NastReport report(reportFName);
521 m->openInputFile(filename, inFASTA);
523 inFASTA.seekg(filePos->start);
528 //moved this into driver to avoid deep copies in windows paralellized version
529 Alignment* alignment;
530 int longestBase = templateDB->getLongestBase();
531 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
532 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
533 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
534 else if(align == "noalign") { alignment = new NoAlign(); }
536 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
537 m->mothurOutEndLine();
538 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
543 if (m->control_pressed) { break; }
545 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
546 report.setCandidate(candidateSeq);
548 int origNumBases = candidateSeq->getNumBases();
549 string originalUnaligned = candidateSeq->getUnaligned();
550 int numBasesNeeded = origNumBases * threshold;
552 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
553 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
554 alignment->resize(candidateSeq->getUnaligned().length()+1);
557 Sequence temp = templateDB->findClosestSequence(candidateSeq);
558 Sequence* templateSeq = &temp;
560 float searchScore = templateDB->getSearchScore();
562 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
567 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
568 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
569 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
570 //so this bool tells you if you need to delete it
572 //if there is a possibility that this sequence should be reversed
573 if (candidateSeq->getNumBases() < numBasesNeeded) {
575 string wasBetter = "";
576 //if the user wants you to try the reverse
579 //get reverse compliment
580 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
581 copy->reverseComplement();
584 Sequence temp2 = templateDB->findClosestSequence(copy);
585 Sequence* templateSeq2 = &temp2;
587 searchScore = templateDB->getSearchScore();
589 nast2 = new Nast(alignment, copy, templateSeq2);
591 //check if any better
592 if (copy->getNumBases() > candidateSeq->getNumBases()) {
593 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
594 templateSeq = templateSeq2;
597 needToDeleteCopy = true;
598 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
600 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
606 //create accnos file with names
607 accnosFile << candidateSeq->getName() << wasBetter << endl;
610 report.setTemplate(templateSeq);
611 report.setSearchParameters(search, searchScore);
612 report.setAlignmentParameters(align, alignment);
613 report.setNastParameters(*nast);
615 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
619 if (needToDeleteCopy) { delete copy; }
625 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
626 unsigned long int pos = inFASTA.tellg();
627 if ((pos == -1) || (pos >= filePos->end)) { break; }
629 if (inFASTA.eof()) { break; }
633 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
637 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
640 alignmentFile.close();
646 catch(exception& e) {
647 m->errorOut(e, "AlignCommand", "driver");
651 //**********************************************************************************************************************
653 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
655 string outputString = "";
656 MPI_Status statusReport;
657 MPI_Status statusAlign;
658 MPI_Status statusAccnos;
661 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
666 outputString = report.getHeaders();
667 int length = outputString.length();
669 char* buf = new char[length];
670 memcpy(buf, outputString.c_str(), length);
672 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
677 Alignment* alignment;
678 int longestBase = templateDB->getLongestBase();
679 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
680 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
681 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
682 else if(align == "noalign") { alignment = new NoAlign(); }
684 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
685 m->mothurOutEndLine();
686 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
690 for(int i=0;i<num;i++){
692 if (m->control_pressed) { delete alignment; return 0; }
695 int length = MPIPos[start+i+1] - MPIPos[start+i];
697 char* buf4 = new char[length];
698 //memcpy(buf4, outputString.c_str(), length);
700 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
702 string tempBuf = buf4;
706 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
708 istringstream iss (tempBuf,istringstream::in);
710 Sequence* candidateSeq = new Sequence(iss);
711 report.setCandidate(candidateSeq);
713 int origNumBases = candidateSeq->getNumBases();
714 string originalUnaligned = candidateSeq->getUnaligned();
715 int numBasesNeeded = origNumBases * threshold;
717 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
718 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
719 alignment->resize(candidateSeq->getUnaligned().length()+1);
722 Sequence temp = templateDB->findClosestSequence(candidateSeq);
723 Sequence* templateSeq = &temp;
725 float searchScore = templateDB->getSearchScore();
727 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
731 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
732 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
733 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
734 //so this bool tells you if you need to delete it
736 //if there is a possibility that this sequence should be reversed
737 if (candidateSeq->getNumBases() < numBasesNeeded) {
739 string wasBetter = "";
740 //if the user wants you to try the reverse
742 //get reverse compliment
743 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
744 copy->reverseComplement();
747 Sequence temp2 = templateDB->findClosestSequence(copy);
748 Sequence* templateSeq2 = &temp2;
750 searchScore = templateDB->getSearchScore();
752 nast2 = new Nast(alignment, copy, templateSeq2);
754 //check if any better
755 if (copy->getNumBases() > candidateSeq->getNumBases()) {
756 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
757 templateSeq = templateSeq2;
760 needToDeleteCopy = true;
761 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
763 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
769 //create accnos file with names
770 outputString = candidateSeq->getName() + wasBetter + "\n";
772 //send results to parent
773 int length = outputString.length();
775 char* buf = new char[length];
776 memcpy(buf, outputString.c_str(), length);
778 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
780 MPIWroteAccnos = true;
783 report.setTemplate(templateSeq);
784 report.setSearchParameters(search, searchScore);
785 report.setAlignmentParameters(align, alignment);
786 report.setNastParameters(*nast);
788 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
790 //send results to parent
791 int length = outputString.length();
792 char* buf2 = new char[length];
793 memcpy(buf2, outputString.c_str(), length);
795 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
799 outputString = report.getReport();
801 //send results to parent
802 length = outputString.length();
803 char* buf3 = new char[length];
804 memcpy(buf3, outputString.c_str(), length);
806 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
810 if (needToDeleteCopy) { delete copy; }
815 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
818 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
822 catch(exception& e) {
823 m->errorOut(e, "AlignCommand", "driverMPI");
828 /**************************************************************************************************/
830 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
833 processIDS.resize(0);
834 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
837 //loop through and create all the processes you want
838 while (process != processors) {
842 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
845 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
847 //pass numSeqs to parent
849 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
850 m->openOutputFile(tempFile, out);
856 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
857 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
863 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
865 //force parent to wait until all the processes are done
866 for (int i=0;i<processIDS.size();i++) {
867 int temp = processIDS[i];
871 vector<string> nonBlankAccnosFiles;
872 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
873 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
875 for (int i = 0; i < processIDS.size(); i++) {
877 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
878 m->openInputFile(tempFile, in);
879 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
880 in.close(); m->mothurRemove(tempFile);
882 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
883 m->mothurRemove((alignFileName + toString(processIDS[i]) + ".temp"));
885 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
886 m->mothurRemove((reportFileName + toString(processIDS[i]) + ".temp"));
888 if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
889 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
890 }else { m->mothurRemove((accnosFName + toString(processIDS[i]) + ".temp")); }
894 //append accnos files
895 if (nonBlankAccnosFiles.size() != 0) {
896 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
898 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
899 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
900 m->mothurRemove(nonBlankAccnosFiles[h]);
902 }else { //recreate the accnosfile if needed
904 m->openOutputFile(accnosFName, out);
908 //////////////////////////////////////////////////////////////////////////////////////////////////////
909 //Windows version shared memory, so be careful when passing variables through the alignData struct.
910 //Above fork() will clone, so memory is separate, but that's not the case with windows,
911 //////////////////////////////////////////////////////////////////////////////////////////////////////
913 vector<alignData*> pDataArray;
914 DWORD dwThreadIdArray[processors-1];
915 HANDLE hThreadArray[processors-1];
917 //Create processor worker threads.
918 for( int i=0; i<processors-1; i++ ){
920 //AlignmentDB* tempDB = new AlignmentDB(*templateDB);
922 // Allocate memory for thread data.
923 string extension = "";
924 if (i != 0) { extension = toString(i) + ".temp"; }
926 alignData* tempalign = new alignData((alignFileName + extension), (reportFileName + extension), (accnosFName + extension), filename, align, search, kmerSize, m, lines[i]->start, lines[i]->end, flip, match, misMatch, gapOpen, gapExtend, threshold, i);
927 pDataArray.push_back(tempalign);
928 processIDS.push_back(i);
930 //MySeqSumThreadFunction is in header. It must be global or static to work with the threads.
931 //default security attributes, thread function name, argument to thread function, use default creation flags, returns the thread identifier
932 hThreadArray[i] = CreateThread(NULL, 0, MyAlignThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
935 //need to check for line ending error
937 m->openInputFile(filename, inFASTA);
938 inFASTA.seekg(lines[processors-1]->start-1);
939 char c = inFASTA.peek();
941 if (c != '>') { //we need to move back
942 lines[processors-1]->start--;
945 //using the main process as a worker saves time and memory
946 //do my part - do last piece because windows is looking for eof
947 num = driver(lines[processors-1], (alignFileName + toString(processors-1) + ".temp"), (reportFileName + toString(processors-1) + ".temp"), (accnosFName + toString(processors-1) + ".temp"), filename);
949 //Wait until all threads have terminated.
950 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
952 //Close all thread handles and free memory allocations.
953 for(int i=0; i < pDataArray.size(); i++){
954 num += pDataArray[i]->count;
955 CloseHandle(hThreadArray[i]);
956 delete pDataArray[i];
959 vector<string> nonBlankAccnosFiles;
960 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
961 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
963 for (int i = 1; i < processors; i++) {
964 appendAlignFiles((alignFileName + toString(i) + ".temp"), alignFileName);
965 m->mothurRemove((alignFileName + toString(i) + ".temp"));
967 appendReportFiles((reportFileName + toString(i) + ".temp"), reportFileName);
968 m->mothurRemove((reportFileName + toString(i) + ".temp"));
970 if (!(m->isBlank(accnosFName + toString(i) + ".temp"))) {
971 nonBlankAccnosFiles.push_back(accnosFName + toString(i) + ".temp");
972 }else { m->mothurRemove((accnosFName + toString(i) + ".temp")); }
975 //append accnos files
976 if (nonBlankAccnosFiles.size() != 0) {
977 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
979 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
980 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
981 m->mothurRemove(nonBlankAccnosFiles[h]);
983 }else { //recreate the accnosfile if needed
985 m->openOutputFile(accnosFName, out);
992 catch(exception& e) {
993 m->errorOut(e, "AlignCommand", "createProcesses");
997 /**************************************************************************************************/
999 void AlignCommand::appendAlignFiles(string temp, string filename) {
1004 m->openOutputFileAppend(filename, output);
1005 m->openInputFile(temp, input);
1007 while(char c = input.get()){
1008 if(input.eof()) { break; }
1009 else { output << c; }
1015 catch(exception& e) {
1016 m->errorOut(e, "AlignCommand", "appendAlignFiles");
1020 //**********************************************************************************************************************
1022 void AlignCommand::appendReportFiles(string temp, string filename) {
1027 m->openOutputFileAppend(filename, output);
1028 m->openInputFile(temp, input);
1030 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
1032 while(char c = input.get()){
1033 if(input.eof()) { break; }
1034 else { output << c; }
1040 catch(exception& e) {
1041 m->errorOut(e, "AlignCommand", "appendReportFiles");
1045 //**********************************************************************************************************************