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,true); parameters.push_back(ptemplate);
24 CommandParameter pcandidate("fasta", "InputTypes", "", "", "none", "none", "none","fasta-alignreport-accnos",false,true,true); parameters.push_back(pcandidate);
25 CommandParameter psearch("search", "Multiple", "kmer-blast-suffix", "kmer", "", "", "","",false,false,true); 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,true); parameters.push_back(palign);
29 CommandParameter pmismatch("mismatch", "Number", "", "-1.0", "", "", "","",false,false); parameters.push_back(pmismatch);
30 CommandParameter pgapopen("gapopen", "Number", "", "-5.0", "", "", "","",false,false); parameters.push_back(pgapopen);
31 CommandParameter pgapextend("gapextend", "Number", "", "-2.0", "", "", "","",false,false); parameters.push_back(pgapextend);
32 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "","",false,false,true); 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 -5.0.";
61 helpString += "The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -2.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::getOutputPattern(string type) {
83 if (type == "fasta") { pattern = "[filename],align"; } //makes file like: amazon.align
84 else if (type == "alignreport") { pattern = "[filename],align.report"; }
85 else if (type == "accnos") { pattern = "[filename],flip.accnos"; }
86 else { m->mothurOut("[ERROR]: No definition for type " + type + " output pattern.\n"); m->control_pressed = true; }
91 m->errorOut(e, "AlignCommand", "getOutputPattern");
95 //**********************************************************************************************************************
96 AlignCommand::AlignCommand(){
98 abort = true; calledHelp = true;
100 vector<string> tempOutNames;
101 outputTypes["fasta"] = tempOutNames;
102 outputTypes["alignreport"] = tempOutNames;
103 outputTypes["accnos"] = tempOutNames;
105 catch(exception& e) {
106 m->errorOut(e, "AlignCommand", "AlignCommand");
110 //**********************************************************************************************************************
111 AlignCommand::AlignCommand(string option) {
113 abort = false; calledHelp = false;
114 ReferenceDB* rdb = ReferenceDB::getInstance();
116 //allow user to run help
117 if(option == "help") { help(); abort = true; calledHelp = true;}
118 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
121 vector<string> myArray = setParameters();
123 OptionParser parser(option);
124 map<string, string> parameters = parser.getParameters();
126 ValidParameters validParameter("align.seqs");
127 map<string, string>::iterator it;
129 //check to make sure all parameters are valid for command
130 for (it = parameters.begin(); it != parameters.end(); it++) {
131 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
134 //initialize outputTypes
135 vector<string> tempOutNames;
136 outputTypes["fasta"] = tempOutNames;
137 outputTypes["alignreport"] = tempOutNames;
138 outputTypes["accnos"] = tempOutNames;
140 //if the user changes the output directory command factory will send this info to us in the output parameter
141 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
144 //if the user changes the input directory command factory will send this info to us in the output parameter
145 string inputDir = validParameter.validFile(parameters, "inputdir", false);
147 if (inputDir == "not found"){ inputDir = ""; }
151 it = parameters.find("reference");
153 //user has given a template file
154 if(it != parameters.end()){
155 path = m->hasPath(it->second);
156 //if the user has not given a path then, add inputdir. else leave path alone.
157 if (path == "") { parameters["reference"] = inputDir + it->second; }
161 candidateFileName = validParameter.validFile(parameters, "fasta", false);
162 if (candidateFileName == "not found") {
163 //if there is a current fasta file, use it
164 string filename = m->getFastaFile();
165 if (filename != "") { candidateFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
166 else { m->mothurOut("You have no current fastafile and the candidate parameter is required."); m->mothurOutEndLine(); abort = true; }
168 m->splitAtDash(candidateFileName, candidateFileNames);
170 //go through files and make sure they are good, if not, then disregard them
171 for (int i = 0; i < candidateFileNames.size(); i++) {
172 //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
175 if (candidateFileNames[i] == "current") {
176 candidateFileNames[i] = m->getFastaFile();
177 if (candidateFileNames[i] != "") { m->mothurOut("Using " + candidateFileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
179 m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true;
180 //erase from file list
181 candidateFileNames.erase(candidateFileNames.begin()+i);
188 if (inputDir != "") {
189 string path = m->hasPath(candidateFileNames[i]);
190 //if the user has not given a path then, add inputdir. else leave path alone.
191 if (path == "") { candidateFileNames[i] = inputDir + candidateFileNames[i]; }
196 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
199 //if you can't open it, try default location
200 if (ableToOpen == 1) {
201 if (m->getDefaultPath() != "") { //default path is set
202 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
203 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
205 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
207 candidateFileNames[i] = tryPath;
211 //if you can't open it, try output location
212 if (ableToOpen == 1) {
213 if (m->getOutputDir() != "") { //default path is set
214 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
215 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
217 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
219 candidateFileNames[i] = tryPath;
225 if (ableToOpen == 1) {
226 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine();
227 //erase from file list
228 candidateFileNames.erase(candidateFileNames.begin()+i);
231 m->setFastaFile(candidateFileNames[i]);
236 //make sure there is at least one valid file left
237 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
240 //check for optional parameter and set defaults
241 // ...at some point should added some additional type checking...
243 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
244 m->mothurConvert(temp, kmerSize);
246 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
247 m->mothurConvert(temp, match);
249 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
250 m->mothurConvert(temp, misMatch);
252 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-5.0"; }
253 m->mothurConvert(temp, gapOpen);
255 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-2.0"; }
256 m->mothurConvert(temp, gapExtend);
258 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
259 m->setProcessors(temp);
260 m->mothurConvert(temp, processors);
262 temp = validParameter.validFile(parameters, "flip", false); if (temp == "not found"){ temp = "f"; }
263 flip = m->isTrue(temp);
265 temp = validParameter.validFile(parameters, "save", false); if (temp == "not found"){ temp = "f"; }
266 save = m->isTrue(temp);
268 if (save) { //clear out old references
272 //this has to go after save so that if the user sets save=t and provides no reference we abort
273 templateFileName = validParameter.validFile(parameters, "reference", true);
274 if (templateFileName == "not found") {
275 //check for saved reference sequences
276 if (rdb->referenceSeqs.size() != 0) {
277 templateFileName = "saved";
279 m->mothurOut("[ERROR]: You don't have any saved reference sequences and the reference parameter is a required for the align.seqs command.");
280 m->mothurOutEndLine();
283 }else if (templateFileName == "not open") { abort = true; }
284 else { if (save) { rdb->setSavedReference(templateFileName); } }
286 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
287 m->mothurConvert(temp, threshold);
289 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
290 if ((search != "suffix") && (search != "kmer") && (search != "blast")) { m->mothurOut("invalid search option: choices are kmer, suffix or blast."); m->mothurOutEndLine(); abort=true; }
292 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
293 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; }
298 catch(exception& e) {
299 m->errorOut(e, "AlignCommand", "AlignCommand");
303 //**********************************************************************************************************************
304 AlignCommand::~AlignCommand(){
306 if (abort == false) {
307 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
311 //**********************************************************************************************************************
313 int AlignCommand::execute(){
315 if (abort == true) { if (calledHelp) { return 0; } return 2; }
317 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch, rand());
319 for (int s = 0; s < candidateFileNames.size(); s++) {
320 if (m->control_pressed) { outputTypes.clear(); return 0; }
322 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
324 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
325 map<string, string> variables; variables["[filename]"] = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s]));
326 string alignFileName = getOutputFileName("fasta", variables);
327 string reportFileName = getOutputFileName("alignreport", variables);
328 string accnosFileName = getOutputFileName("accnos", variables);
330 bool hasAccnos = true;
332 int numFastaSeqs = 0;
333 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
334 int start = time(NULL);
337 int pid, numSeqsPerProcessor;
339 vector<unsigned long long> MPIPos;
340 MPIWroteAccnos = false;
343 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
344 MPI_Comm_size(MPI_COMM_WORLD, &processors);
347 MPI_File outMPIAlign;
348 MPI_File outMPIReport;
349 MPI_File outMPIAccnos;
351 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
352 int inMode=MPI_MODE_RDONLY;
354 char outAlignFilename[1024];
355 strcpy(outAlignFilename, alignFileName.c_str());
357 char outReportFilename[1024];
358 strcpy(outReportFilename, reportFileName.c_str());
360 char outAccnosFilename[1024];
361 strcpy(outAccnosFilename, accnosFileName.c_str());
363 char inFileName[1024];
364 strcpy(inFileName, candidateFileNames[s].c_str());
366 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
367 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
368 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
369 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
371 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
373 if (pid == 0) { //you are the root process
375 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
377 //send file positions to all processes
378 for(int i = 1; i < processors; i++) {
379 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
380 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
383 //figure out how many sequences you have to align
384 numSeqsPerProcessor = numFastaSeqs / processors;
385 int startIndex = pid * numSeqsPerProcessor;
386 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
389 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
391 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
393 for (int i = 1; i < processors; i++) {
395 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
396 if (tempResult != 0) { MPIWroteAccnos = true; }
398 }else{ //you are a child process
399 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
400 MPIPos.resize(numFastaSeqs+1);
401 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
404 //figure out how many sequences you have to align
405 numSeqsPerProcessor = numFastaSeqs / processors;
406 int startIndex = pid * numSeqsPerProcessor;
407 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
411 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
413 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
415 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
419 MPI_File_close(&inMPI);
420 MPI_File_close(&outMPIAlign);
421 MPI_File_close(&outMPIReport);
422 MPI_File_close(&outMPIAccnos);
424 //delete accnos file if blank
426 //delete accnos file if its blank else report to user
427 if (MPIWroteAccnos) {
428 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
430 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
431 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
432 m->mothurOutEndLine();
435 //MPI_File_delete(outAccnosFilename, info);
437 m->mothurRemove(accnosFileName);
443 vector<unsigned long long> positions;
444 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
445 positions = m->divideFile(candidateFileNames[s], processors);
446 for (int i = 0; i < (positions.size()-1); i++) { lines.push_back(new linePair(positions[i], positions[(i+1)])); }
448 if (processors == 1) {
449 lines.push_back(new linePair(0, 1000));
451 positions = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs);
452 if (positions.size() < processors) { processors = positions.size(); }
454 //figure out how many sequences you have to process
455 int numSeqsPerProcessor = numFastaSeqs / processors;
456 for (int i = 0; i < processors; i++) {
457 int startIndex = i * numSeqsPerProcessor;
458 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
459 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
465 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
467 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
470 if (m->control_pressed) { m->mothurRemove(accnosFileName); m->mothurRemove(alignFileName); m->mothurRemove(reportFileName); outputTypes.clear(); return 0; }
472 //delete accnos file if its blank else report to user
473 if (m->isBlank(accnosFileName)) { m->mothurRemove(accnosFileName); hasAccnos = false; }
475 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
477 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
478 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
479 m->mothurOutEndLine();
486 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
488 if (pid == 0) { //only one process should output to screen
491 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
492 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
493 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
499 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
500 m->mothurOutEndLine();
501 m->mothurOutEndLine();
504 //set align file as new current fastafile
505 string currentFasta = "";
506 itTypes = outputTypes.find("fasta");
507 if (itTypes != outputTypes.end()) {
508 if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; m->setFastaFile(currentFasta); }
511 m->mothurOutEndLine();
512 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
513 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
514 m->mothurOutEndLine();
518 catch(exception& e) {
519 m->errorOut(e, "AlignCommand", "execute");
524 //**********************************************************************************************************************
525 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
527 ofstream alignmentFile;
528 m->openOutputFile(alignFName, alignmentFile);
531 m->openOutputFile(accnosFName, accnosFile);
533 NastReport report(reportFName);
536 m->openInputFile(filename, inFASTA);
538 inFASTA.seekg(filePos->start);
543 //moved this into driver to avoid deep copies in windows paralellized version
544 Alignment* alignment;
545 int longestBase = templateDB->getLongestBase();
546 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
547 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
548 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
549 else if(align == "noalign") { alignment = new NoAlign(); }
551 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
552 m->mothurOutEndLine();
553 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
558 if (m->control_pressed) { break; }
560 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
561 cout << candidateSeq->getAligned() << endl;
562 report.setCandidate(candidateSeq);
564 int origNumBases = candidateSeq->getNumBases();
565 string originalUnaligned = candidateSeq->getUnaligned();
566 int numBasesNeeded = origNumBases * threshold;
568 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
569 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
570 alignment->resize(candidateSeq->getUnaligned().length()+1);
572 Sequence temp = templateDB->findClosestSequence(candidateSeq);
573 Sequence* templateSeq = &temp;
575 float searchScore = templateDB->getSearchScore();
577 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
582 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
583 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
584 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
585 //so this bool tells you if you need to delete it
587 //if there is a possibility that this sequence should be reversed
588 if (candidateSeq->getNumBases() < numBasesNeeded) {
590 string wasBetter = "";
591 //if the user wants you to try the reverse
594 //get reverse compliment
595 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
596 copy->reverseComplement();
599 Sequence temp2 = templateDB->findClosestSequence(copy);
600 Sequence* templateSeq2 = &temp2;
602 searchScore = templateDB->getSearchScore();
604 nast2 = new Nast(alignment, copy, templateSeq2);
606 //check if any better
607 if (copy->getNumBases() > candidateSeq->getNumBases()) {
608 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
609 templateSeq = templateSeq2;
612 needToDeleteCopy = true;
613 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
615 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
621 //create accnos file with names
622 accnosFile << candidateSeq->getName() << wasBetter << endl;
625 report.setTemplate(templateSeq);
626 report.setSearchParameters(search, searchScore);
627 report.setAlignmentParameters(align, alignment);
628 report.setNastParameters(*nast);
630 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
634 if (needToDeleteCopy) { delete copy; }
640 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
641 unsigned long long pos = inFASTA.tellg();
642 if ((pos == -1) || (pos >= filePos->end)) { break; }
644 if (inFASTA.eof()) { break; }
648 if((count) % 100 == 0){ m->mothurOutJustToScreen(toString(count) + "\n"); }
652 if((count) % 100 != 0){ m->mothurOutJustToScreen(toString(count) + "\n"); }
655 alignmentFile.close();
661 catch(exception& e) {
662 m->errorOut(e, "AlignCommand", "driver");
666 //**********************************************************************************************************************
668 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long long>& MPIPos){
670 string outputString = "";
671 MPI_Status statusReport;
672 MPI_Status statusAlign;
673 MPI_Status statusAccnos;
676 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
681 outputString = report.getHeaders();
682 int length = outputString.length();
684 char* buf = new char[length];
685 memcpy(buf, outputString.c_str(), length);
687 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
692 Alignment* alignment;
693 int longestBase = templateDB->getLongestBase();
694 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
695 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
696 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
697 else if(align == "noalign") { alignment = new NoAlign(); }
699 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
700 m->mothurOutEndLine();
701 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
705 for(int i=0;i<num;i++){
707 if (m->control_pressed) { delete alignment; return 0; }
710 int length = MPIPos[start+i+1] - MPIPos[start+i];
712 char* buf4 = new char[length];
713 //memcpy(buf4, outputString.c_str(), length);
715 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
717 string tempBuf = buf4;
721 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
723 istringstream iss (tempBuf,istringstream::in);
725 Sequence* candidateSeq = new Sequence(iss);
726 report.setCandidate(candidateSeq);
728 int origNumBases = candidateSeq->getNumBases();
729 string originalUnaligned = candidateSeq->getUnaligned();
730 int numBasesNeeded = origNumBases * threshold;
732 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
733 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
734 alignment->resize(candidateSeq->getUnaligned().length()+1);
737 Sequence temp = templateDB->findClosestSequence(candidateSeq);
738 Sequence* templateSeq = &temp;
740 float searchScore = templateDB->getSearchScore();
742 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
746 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
747 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
748 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
749 //so this bool tells you if you need to delete it
751 //if there is a possibility that this sequence should be reversed
752 if (candidateSeq->getNumBases() < numBasesNeeded) {
754 string wasBetter = "";
755 //if the user wants you to try the reverse
757 //get reverse compliment
758 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
759 copy->reverseComplement();
762 Sequence temp2 = templateDB->findClosestSequence(copy);
763 Sequence* templateSeq2 = &temp2;
765 searchScore = templateDB->getSearchScore();
767 nast2 = new Nast(alignment, copy, templateSeq2);
769 //check if any better
770 if (copy->getNumBases() > candidateSeq->getNumBases()) {
771 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
772 templateSeq = templateSeq2;
775 needToDeleteCopy = true;
776 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
778 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
784 //create accnos file with names
785 outputString = candidateSeq->getName() + wasBetter + "\n";
787 //send results to parent
788 int length = outputString.length();
790 char* buf = new char[length];
791 memcpy(buf, outputString.c_str(), length);
793 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
795 MPIWroteAccnos = true;
798 report.setTemplate(templateSeq);
799 report.setSearchParameters(search, searchScore);
800 report.setAlignmentParameters(align, alignment);
801 report.setNastParameters(*nast);
803 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
805 //send results to parent
806 int length = outputString.length();
807 char* buf2 = new char[length];
808 memcpy(buf2, outputString.c_str(), length);
810 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
814 outputString = report.getReport();
816 //send results to parent
817 length = outputString.length();
818 char* buf3 = new char[length];
819 memcpy(buf3, outputString.c_str(), length);
821 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
825 if (needToDeleteCopy) { delete copy; }
830 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
833 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
837 catch(exception& e) {
838 m->errorOut(e, "AlignCommand", "driverMPI");
843 /**************************************************************************************************/
845 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
848 processIDS.resize(0);
849 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
852 //loop through and create all the processes you want
853 while (process != processors) {
857 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
860 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
862 //pass numSeqs to parent
864 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
865 m->openOutputFile(tempFile, out);
871 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
872 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
878 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
880 //force parent to wait until all the processes are done
881 for (int i=0;i<processIDS.size();i++) {
882 int temp = processIDS[i];
886 vector<string> nonBlankAccnosFiles;
887 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
888 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
890 for (int i = 0; i < processIDS.size(); i++) {
892 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
893 m->openInputFile(tempFile, in);
894 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
895 in.close(); m->mothurRemove(tempFile);
897 m->appendFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
898 m->mothurRemove((alignFileName + toString(processIDS[i]) + ".temp"));
900 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
901 m->mothurRemove((reportFileName + toString(processIDS[i]) + ".temp"));
903 if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
904 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
905 }else { m->mothurRemove((accnosFName + toString(processIDS[i]) + ".temp")); }
909 //append accnos files
910 if (nonBlankAccnosFiles.size() != 0) {
911 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
913 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
914 m->appendFiles(nonBlankAccnosFiles[h], accnosFName);
915 m->mothurRemove(nonBlankAccnosFiles[h]);
917 }else { //recreate the accnosfile if needed
919 m->openOutputFile(accnosFName, out);
923 //////////////////////////////////////////////////////////////////////////////////////////////////////
924 //Windows version shared memory, so be careful when passing variables through the alignData struct.
925 //Above fork() will clone, so memory is separate, but that's not the case with windows,
926 //////////////////////////////////////////////////////////////////////////////////////////////////////
928 vector<alignData*> pDataArray;
929 DWORD dwThreadIdArray[processors-1];
930 HANDLE hThreadArray[processors-1];
932 //Create processor worker threads.
933 for( int i=0; i<processors-1; i++ ){
935 //AlignmentDB* tempDB = new AlignmentDB(*templateDB);
937 // Allocate memory for thread data.
938 string extension = "";
939 if (i != 0) { extension = toString(i) + ".temp"; }
941 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);
942 pDataArray.push_back(tempalign);
943 processIDS.push_back(i);
945 //MySeqSumThreadFunction is in header. It must be global or static to work with the threads.
946 //default security attributes, thread function name, argument to thread function, use default creation flags, returns the thread identifier
947 hThreadArray[i] = CreateThread(NULL, 0, MyAlignThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
950 //need to check for line ending error
952 m->openInputFile(filename, inFASTA);
953 inFASTA.seekg(lines[processors-1]->start-1);
954 char c = inFASTA.peek();
956 if (c != '>') { //we need to move back
957 lines[processors-1]->start--;
960 //using the main process as a worker saves time and memory
961 //do my part - do last piece because windows is looking for eof
962 num = driver(lines[processors-1], (alignFileName + toString(processors-1) + ".temp"), (reportFileName + toString(processors-1) + ".temp"), (accnosFName + toString(processors-1) + ".temp"), filename);
964 //Wait until all threads have terminated.
965 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
967 //Close all thread handles and free memory allocations.
968 for(int i=0; i < pDataArray.size(); i++){
969 if (pDataArray[i]->count != pDataArray[i]->end) {
970 m->mothurOut("[ERROR]: process " + toString(i) + " only processed " + toString(pDataArray[i]->count) + " of " + toString(pDataArray[i]->end) + " sequences assigned to it, quitting. \n"); m->control_pressed = true;
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 //**********************************************************************************************************************