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 string AlignCommand::getOutputFileNameTag(string type, string inputName=""){
82 map<string, vector<string> >::iterator it;
84 //is this a type this command creates
85 it = outputTypes.find(type);
86 if (it == outputTypes.end()) { m->mothurOut("[ERROR]: this command doesn't create a " + type + " output file.\n"); }
88 if (type == "fasta") { tag = "align"; }
89 else if (type == "alignreport") { tag = "align.report"; }
90 else if (type == "accnos") { tag = "flip.accnos"; }
91 else { m->mothurOut("[ERROR]: No definition for type " + type + " output file tag.\n"); m->control_pressed = true; }
96 m->errorOut(e, "AlignCommand", "getOutputFileName");
101 //**********************************************************************************************************************
102 AlignCommand::AlignCommand(){
104 abort = true; calledHelp = true;
106 vector<string> tempOutNames;
107 outputTypes["fasta"] = tempOutNames;
108 outputTypes["alignreport"] = tempOutNames;
109 outputTypes["accnos"] = tempOutNames;
111 catch(exception& e) {
112 m->errorOut(e, "AlignCommand", "AlignCommand");
116 //**********************************************************************************************************************
117 AlignCommand::AlignCommand(string option) {
119 abort = false; calledHelp = false;
120 ReferenceDB* rdb = ReferenceDB::getInstance();
122 //allow user to run help
123 if(option == "help") { help(); abort = true; calledHelp = true;}
124 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
127 vector<string> myArray = setParameters();
129 OptionParser parser(option);
130 map<string, string> parameters = parser.getParameters();
132 ValidParameters validParameter("align.seqs");
133 map<string, string>::iterator it;
135 //check to make sure all parameters are valid for command
136 for (it = parameters.begin(); it != parameters.end(); it++) {
137 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
140 //initialize outputTypes
141 vector<string> tempOutNames;
142 outputTypes["fasta"] = tempOutNames;
143 outputTypes["alignreport"] = tempOutNames;
144 outputTypes["accnos"] = tempOutNames;
146 //if the user changes the output directory command factory will send this info to us in the output parameter
147 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
150 //if the user changes the input directory command factory will send this info to us in the output parameter
151 string inputDir = validParameter.validFile(parameters, "inputdir", false);
153 if (inputDir == "not found"){ inputDir = ""; }
157 it = parameters.find("reference");
159 //user has given a template file
160 if(it != parameters.end()){
161 path = m->hasPath(it->second);
162 //if the user has not given a path then, add inputdir. else leave path alone.
163 if (path == "") { parameters["reference"] = inputDir + it->second; }
167 candidateFileName = validParameter.validFile(parameters, "fasta", false);
168 if (candidateFileName == "not found") {
169 //if there is a current fasta file, use it
170 string filename = m->getFastaFile();
171 if (filename != "") { candidateFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
172 else { m->mothurOut("You have no current fastafile and the candidate parameter is required."); m->mothurOutEndLine(); abort = true; }
174 m->splitAtDash(candidateFileName, candidateFileNames);
176 //go through files and make sure they are good, if not, then disregard them
177 for (int i = 0; i < candidateFileNames.size(); i++) {
178 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
181 if (candidateFileNames[i] == "current") {
182 candidateFileNames[i] = m->getFastaFile();
183 if (candidateFileNames[i] != "") { m->mothurOut("Using " + candidateFileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
185 m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true;
186 //erase from file list
187 candidateFileNames.erase(candidateFileNames.begin()+i);
194 if (inputDir != "") {
195 string path = m->hasPath(candidateFileNames[i]);
196 //if the user has not given a path then, add inputdir. else leave path alone.
197 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
202 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
205 //if you can't open it, try default location
206 if (ableToOpen == 1) {
207 if (m->getDefaultPath() != "") { //default path is set
208 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
209 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
211 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
213 candidateFileNames[i] = tryPath;
217 //if you can't open it, try output location
218 if (ableToOpen == 1) {
219 if (m->getOutputDir() != "") { //default path is set
220 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
221 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
223 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
225 candidateFileNames[i] = tryPath;
231 if (ableToOpen == 1) {
232 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
233 //erase from file list
234 candidateFileNames.erase(candidateFileNames.begin()+i);
237 m->setFastaFile(candidateFileNames[i]);
242 //make sure there is at least one valid file left
243 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
246 //check for optional parameter and set defaults
247 // ...at some point should added some additional type checking...
249 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
250 m->mothurConvert(temp, kmerSize);
252 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
253 m->mothurConvert(temp, match);
255 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
256 m->mothurConvert(temp, misMatch);
258 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
259 m->mothurConvert(temp, gapOpen);
261 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
262 m->mothurConvert(temp, gapExtend);
264 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
265 m->setProcessors(temp);
266 m->mothurConvert(temp, processors);
268 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
269 flip = m->isTrue(temp);
271 temp = validParameter.validFile(parameters, "save", false); if (temp == "not found"){ temp = "f"; }
272 save = m->isTrue(temp);
274 if (save) { //clear out old references
278 //this has to go after save so that if the user sets save=t and provides no reference we abort
279 templateFileName = validParameter.validFile(parameters, "reference", true);
280 if (templateFileName == "not found") {
281 //check for saved reference sequences
282 if (rdb->referenceSeqs.size() != 0) {
283 templateFileName = "saved";
285 m->mothurOut("[ERROR]: You don't have any saved reference sequences and the reference parameter is a required for the align.seqs command.");
286 m->mothurOutEndLine();
289 }else if (templateFileName == "not open") { abort = true; }
290 else { if (save) { rdb->setSavedReference(templateFileName); } }
292 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
293 m->mothurConvert(temp, threshold);
295 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
296 if ((search != "suffix") && (search != "kmer") && (search != "blast")) { m->mothurOut("invalid search option: choices are kmer, suffix or blast."); m->mothurOutEndLine(); abort=true; }
298 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
299 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; }
304 catch(exception& e) {
305 m->errorOut(e, "AlignCommand", "AlignCommand");
309 //**********************************************************************************************************************
310 AlignCommand::~AlignCommand(){
312 if (abort == false) {
313 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
317 //**********************************************************************************************************************
319 int AlignCommand::execute(){
321 if (abort == true) { if (calledHelp) { return 0; } return 2; }
323 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch, rand());
325 for (int s = 0; s < candidateFileNames.size(); s++) {
326 if (m->control_pressed) { outputTypes.clear(); return 0; }
328 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
330 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
331 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + getOutputFileNameTag("fasta");
332 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + getOutputFileNameTag("alignreport");
333 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + getOutputFileNameTag("accnos");
334 bool hasAccnos = true;
336 int numFastaSeqs = 0;
337 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
338 int start = time(NULL);
341 int pid, numSeqsPerProcessor;
343 vector<unsigned long long> MPIPos;
344 MPIWroteAccnos = false;
347 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
348 MPI_Comm_size(MPI_COMM_WORLD, &processors);
351 MPI_File outMPIAlign;
352 MPI_File outMPIReport;
353 MPI_File outMPIAccnos;
355 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
356 int inMode=MPI_MODE_RDONLY;
358 char outAlignFilename[1024];
359 strcpy(outAlignFilename, alignFileName.c_str());
361 char outReportFilename[1024];
362 strcpy(outReportFilename, reportFileName.c_str());
364 char outAccnosFilename[1024];
365 strcpy(outAccnosFilename, accnosFileName.c_str());
367 char inFileName[1024];
368 strcpy(inFileName, candidateFileNames[s].c_str());
370 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
371 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
372 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
373 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
375 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
377 if (pid == 0) { //you are the root process
379 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
381 //send file positions to all processes
382 for(int i = 1; i < processors; i++) {
383 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
384 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
387 //figure out how many sequences you have to align
388 numSeqsPerProcessor = numFastaSeqs / processors;
389 int startIndex = pid * numSeqsPerProcessor;
390 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
393 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
395 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
397 for (int i = 1; i < processors; i++) {
399 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
400 if (tempResult != 0) { MPIWroteAccnos = true; }
402 }else{ //you are a child process
403 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
404 MPIPos.resize(numFastaSeqs+1);
405 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
408 //figure out how many sequences you have to align
409 numSeqsPerProcessor = numFastaSeqs / processors;
410 int startIndex = pid * numSeqsPerProcessor;
411 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
415 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
417 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
419 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
423 MPI_File_close(&inMPI);
424 MPI_File_close(&outMPIAlign);
425 MPI_File_close(&outMPIReport);
426 MPI_File_close(&outMPIAccnos);
428 //delete accnos file if blank
430 //delete accnos file if its blank else report to user
431 if (MPIWroteAccnos) {
432 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
434 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
435 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
436 m->mothurOutEndLine();
439 //MPI_File_delete(outAccnosFilename, info);
441 m->mothurRemove(accnosFileName);
447 vector<unsigned long long> positions;
448 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
449 positions = m->divideFile(candidateFileNames[s], processors);
450 for (int i = 0; i < (positions.size()-1); i++) { lines.push_back(new linePair(positions[i], positions[(i+1)])); }
452 if (processors == 1) {
453 lines.push_back(new linePair(0, 1000));
455 positions = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs);
456 if (positions.size() < processors) { processors = positions.size(); }
458 //figure out how many sequences you have to process
459 int numSeqsPerProcessor = numFastaSeqs / processors;
460 for (int i = 0; i < processors; i++) {
461 int startIndex = i * numSeqsPerProcessor;
462 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
463 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
469 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
471 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
474 if (m->control_pressed) { m->mothurRemove(accnosFileName); m->mothurRemove(alignFileName); m->mothurRemove(reportFileName); outputTypes.clear(); return 0; }
476 //delete accnos file if its blank else report to user
477 if (m->isBlank(accnosFileName)) { m->mothurRemove(accnosFileName); hasAccnos = false; }
479 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
481 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
482 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
483 m->mothurOutEndLine();
490 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
492 if (pid == 0) { //only one process should output to screen
495 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
496 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
497 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
503 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
504 m->mothurOutEndLine();
505 m->mothurOutEndLine();
508 //set align file as new current fastafile
509 string currentFasta = "";
510 itTypes = outputTypes.find("fasta");
511 if (itTypes != outputTypes.end()) {
512 if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; m->setFastaFile(currentFasta); }
515 m->mothurOutEndLine();
516 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
517 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
518 m->mothurOutEndLine();
522 catch(exception& e) {
523 m->errorOut(e, "AlignCommand", "execute");
528 //**********************************************************************************************************************
529 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
531 ofstream alignmentFile;
532 m->openOutputFile(alignFName, alignmentFile);
535 m->openOutputFile(accnosFName, accnosFile);
537 NastReport report(reportFName);
540 m->openInputFile(filename, inFASTA);
542 inFASTA.seekg(filePos->start);
547 //moved this into driver to avoid deep copies in windows paralellized version
548 Alignment* alignment;
549 int longestBase = templateDB->getLongestBase();
550 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
551 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
552 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
553 else if(align == "noalign") { alignment = new NoAlign(); }
555 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
556 m->mothurOutEndLine();
557 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
562 if (m->control_pressed) { break; }
564 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
565 report.setCandidate(candidateSeq);
567 int origNumBases = candidateSeq->getNumBases();
568 string originalUnaligned = candidateSeq->getUnaligned();
569 int numBasesNeeded = origNumBases * threshold;
571 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
572 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
573 alignment->resize(candidateSeq->getUnaligned().length()+1);
575 Sequence temp = templateDB->findClosestSequence(candidateSeq);
576 Sequence* templateSeq = &temp;
578 float searchScore = templateDB->getSearchScore();
580 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
585 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
586 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
587 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
588 //so this bool tells you if you need to delete it
590 //if there is a possibility that this sequence should be reversed
591 if (candidateSeq->getNumBases() < numBasesNeeded) {
593 string wasBetter = "";
594 //if the user wants you to try the reverse
597 //get reverse compliment
598 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
599 copy->reverseComplement();
602 Sequence temp2 = templateDB->findClosestSequence(copy);
603 Sequence* templateSeq2 = &temp2;
605 searchScore = templateDB->getSearchScore();
607 nast2 = new Nast(alignment, copy, templateSeq2);
609 //check if any better
610 if (copy->getNumBases() > candidateSeq->getNumBases()) {
611 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
612 templateSeq = templateSeq2;
615 needToDeleteCopy = true;
616 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
618 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
624 //create accnos file with names
625 accnosFile << candidateSeq->getName() << wasBetter << endl;
628 report.setTemplate(templateSeq);
629 report.setSearchParameters(search, searchScore);
630 report.setAlignmentParameters(align, alignment);
631 report.setNastParameters(*nast);
633 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
637 if (needToDeleteCopy) { delete copy; }
643 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
644 unsigned long long pos = inFASTA.tellg();
645 if ((pos == -1) || (pos >= filePos->end)) { break; }
647 if (inFASTA.eof()) { break; }
651 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
655 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
658 alignmentFile.close();
664 catch(exception& e) {
665 m->errorOut(e, "AlignCommand", "driver");
669 //**********************************************************************************************************************
671 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long long>& MPIPos){
673 string outputString = "";
674 MPI_Status statusReport;
675 MPI_Status statusAlign;
676 MPI_Status statusAccnos;
679 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
684 outputString = report.getHeaders();
685 int length = outputString.length();
687 char* buf = new char[length];
688 memcpy(buf, outputString.c_str(), length);
690 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
695 Alignment* alignment;
696 int longestBase = templateDB->getLongestBase();
697 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
698 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
699 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
700 else if(align == "noalign") { alignment = new NoAlign(); }
702 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
703 m->mothurOutEndLine();
704 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
708 for(int i=0;i<num;i++){
710 if (m->control_pressed) { delete alignment; return 0; }
713 int length = MPIPos[start+i+1] - MPIPos[start+i];
715 char* buf4 = new char[length];
716 //memcpy(buf4, outputString.c_str(), length);
718 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
720 string tempBuf = buf4;
724 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
726 istringstream iss (tempBuf,istringstream::in);
728 Sequence* candidateSeq = new Sequence(iss);
729 report.setCandidate(candidateSeq);
731 int origNumBases = candidateSeq->getNumBases();
732 string originalUnaligned = candidateSeq->getUnaligned();
733 int numBasesNeeded = origNumBases * threshold;
735 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
736 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
737 alignment->resize(candidateSeq->getUnaligned().length()+1);
740 Sequence temp = templateDB->findClosestSequence(candidateSeq);
741 Sequence* templateSeq = &temp;
743 float searchScore = templateDB->getSearchScore();
745 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
749 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
750 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
751 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
752 //so this bool tells you if you need to delete it
754 //if there is a possibility that this sequence should be reversed
755 if (candidateSeq->getNumBases() < numBasesNeeded) {
757 string wasBetter = "";
758 //if the user wants you to try the reverse
760 //get reverse compliment
761 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
762 copy->reverseComplement();
765 Sequence temp2 = templateDB->findClosestSequence(copy);
766 Sequence* templateSeq2 = &temp2;
768 searchScore = templateDB->getSearchScore();
770 nast2 = new Nast(alignment, copy, templateSeq2);
772 //check if any better
773 if (copy->getNumBases() > candidateSeq->getNumBases()) {
774 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
775 templateSeq = templateSeq2;
778 needToDeleteCopy = true;
779 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
781 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
787 //create accnos file with names
788 outputString = candidateSeq->getName() + wasBetter + "\n";
790 //send results to parent
791 int length = outputString.length();
793 char* buf = new char[length];
794 memcpy(buf, outputString.c_str(), length);
796 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
798 MPIWroteAccnos = true;
801 report.setTemplate(templateSeq);
802 report.setSearchParameters(search, searchScore);
803 report.setAlignmentParameters(align, alignment);
804 report.setNastParameters(*nast);
806 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
808 //send results to parent
809 int length = outputString.length();
810 char* buf2 = new char[length];
811 memcpy(buf2, outputString.c_str(), length);
813 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
817 outputString = report.getReport();
819 //send results to parent
820 length = outputString.length();
821 char* buf3 = new char[length];
822 memcpy(buf3, outputString.c_str(), length);
824 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
828 if (needToDeleteCopy) { delete copy; }
833 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
836 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
840 catch(exception& e) {
841 m->errorOut(e, "AlignCommand", "driverMPI");
846 /**************************************************************************************************/
848 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
851 processIDS.resize(0);
852 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
855 //loop through and create all the processes you want
856 while (process != processors) {
860 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
863 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
865 //pass numSeqs to parent
867 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
868 m->openOutputFile(tempFile, out);
874 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
875 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
881 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
883 //force parent to wait until all the processes are done
884 for (int i=0;i<processIDS.size();i++) {
885 int temp = processIDS[i];
889 vector<string> nonBlankAccnosFiles;
890 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
891 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
893 for (int i = 0; i < processIDS.size(); i++) {
895 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
896 m->openInputFile(tempFile, in);
897 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
898 in.close(); m->mothurRemove(tempFile);
900 m->appendFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
901 m->mothurRemove((alignFileName + toString(processIDS[i]) + ".temp"));
903 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
904 m->mothurRemove((reportFileName + toString(processIDS[i]) + ".temp"));
906 if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
907 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
908 }else { m->mothurRemove((accnosFName + toString(processIDS[i]) + ".temp")); }
912 //append accnos files
913 if (nonBlankAccnosFiles.size() != 0) {
914 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
916 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
917 m->appendFiles(nonBlankAccnosFiles[h], accnosFName);
918 m->mothurRemove(nonBlankAccnosFiles[h]);
920 }else { //recreate the accnosfile if needed
922 m->openOutputFile(accnosFName, out);
926 //////////////////////////////////////////////////////////////////////////////////////////////////////
927 //Windows version shared memory, so be careful when passing variables through the alignData struct.
928 //Above fork() will clone, so memory is separate, but that's not the case with windows,
929 //////////////////////////////////////////////////////////////////////////////////////////////////////
931 vector<alignData*> pDataArray;
932 DWORD dwThreadIdArray[processors-1];
933 HANDLE hThreadArray[processors-1];
935 //Create processor worker threads.
936 for( int i=0; i<processors-1; i++ ){
938 //AlignmentDB* tempDB = new AlignmentDB(*templateDB);
940 // Allocate memory for thread data.
941 string extension = "";
942 if (i != 0) { extension = toString(i) + ".temp"; }
944 alignData* tempalign = new alignData(templateFileName, (alignFileName + extension), (reportFileName + extension), (accnosFName + extension), filename, align, search, kmerSize, m, lines[i]->start, lines[i]->end, flip, match, misMatch, gapOpen, gapExtend, threshold, i);
945 pDataArray.push_back(tempalign);
946 processIDS.push_back(i);
948 //MySeqSumThreadFunction is in header. It must be global or static to work with the threads.
949 //default security attributes, thread function name, argument to thread function, use default creation flags, returns the thread identifier
950 hThreadArray[i] = CreateThread(NULL, 0, MyAlignThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
953 //need to check for line ending error
955 m->openInputFile(filename, inFASTA);
956 inFASTA.seekg(lines[processors-1]->start-1);
957 char c = inFASTA.peek();
959 if (c != '>') { //we need to move back
960 lines[processors-1]->start--;
963 //using the main process as a worker saves time and memory
964 //do my part - do last piece because windows is looking for eof
965 num = driver(lines[processors-1], (alignFileName + toString(processors-1) + ".temp"), (reportFileName + toString(processors-1) + ".temp"), (accnosFName + toString(processors-1) + ".temp"), filename);
967 //Wait until all threads have terminated.
968 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
970 //Close all thread handles and free memory allocations.
971 for(int i=0; i < pDataArray.size(); i++){
972 num += pDataArray[i]->count;
973 CloseHandle(hThreadArray[i]);
974 delete pDataArray[i];
977 vector<string> nonBlankAccnosFiles;
978 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
979 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
981 for (int i = 1; i < processors; i++) {
982 m->appendFiles((alignFileName + toString(i) + ".temp"), alignFileName);
983 m->mothurRemove((alignFileName + toString(i) + ".temp"));
985 appendReportFiles((reportFileName + toString(i) + ".temp"), reportFileName);
986 m->mothurRemove((reportFileName + toString(i) + ".temp"));
988 if (!(m->isBlank(accnosFName + toString(i) + ".temp"))) {
989 nonBlankAccnosFiles.push_back(accnosFName + toString(i) + ".temp");
990 }else { m->mothurRemove((accnosFName + toString(i) + ".temp")); }
993 //append accnos files
994 if (nonBlankAccnosFiles.size() != 0) {
995 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
997 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
998 m->appendFiles(nonBlankAccnosFiles[h], accnosFName);
999 m->mothurRemove(nonBlankAccnosFiles[h]);
1001 }else { //recreate the accnosfile if needed
1003 m->openOutputFile(accnosFName, out);
1010 catch(exception& e) {
1011 m->errorOut(e, "AlignCommand", "createProcesses");
1015 //**********************************************************************************************************************
1017 void AlignCommand::appendReportFiles(string temp, string filename) {
1022 m->openOutputFileAppend(filename, output);
1023 m->openInputFile(temp, input);
1025 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
1028 while (!input.eof()) {
1029 input.read(buffer, 4096);
1030 output.write(buffer, input.gcount());
1036 catch(exception& e) {
1037 m->errorOut(e, "AlignCommand", "appendReportFiles");
1041 //**********************************************************************************************************************