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 m->mothurConvert(temp, kmerSize);
229 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
230 m->mothurConvert(temp, match);
232 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
233 m->mothurConvert(temp, misMatch);
235 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
236 m->mothurConvert(temp, gapOpen);
238 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
239 m->mothurConvert(temp, gapExtend);
241 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
242 m->setProcessors(temp);
243 m->mothurConvert(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);
251 if (save) { //clear out old references
255 //this has to go after save so that if the user sets save=t and provides no reference we abort
256 templateFileName = validParameter.validFile(parameters, "reference", true);
257 if (templateFileName == "not found") {
258 //check for saved reference sequences
259 if (rdb->referenceSeqs.size() != 0) {
260 templateFileName = "saved";
262 m->mothurOut("[ERROR]: You don't have any saved reference sequences and the reference parameter is a required for the align.seqs command.");
263 m->mothurOutEndLine();
266 }else if (templateFileName == "not open") { abort = true; }
267 else { if (save) { rdb->setSavedReference(templateFileName); } }
269 temp = validParameter.validFile(parameters, "threshold", false); if (temp == "not found"){ temp = "0.50"; }
270 m->mothurConvert(temp, threshold);
272 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
273 if ((search != "suffix") && (search != "kmer") && (search != "blast")) { m->mothurOut("invalid search option: choices are kmer, suffix or blast."); m->mothurOutEndLine(); abort=true; }
275 align = validParameter.validFile(parameters, "align", false); if (align == "not found"){ align = "needleman"; }
276 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; }
281 catch(exception& e) {
282 m->errorOut(e, "AlignCommand", "AlignCommand");
286 //**********************************************************************************************************************
287 AlignCommand::~AlignCommand(){
289 if (abort == false) {
290 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
294 //**********************************************************************************************************************
296 int AlignCommand::execute(){
298 if (abort == true) { if (calledHelp) { return 0; } return 2; }
300 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch, rand());
302 for (int s = 0; s < candidateFileNames.size(); s++) {
303 if (m->control_pressed) { outputTypes.clear(); return 0; }
305 m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
307 if (outputDir == "") { outputDir += m->hasPath(candidateFileNames[s]); }
308 string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
309 string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
310 string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
311 bool hasAccnos = true;
313 int numFastaSeqs = 0;
314 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
315 int start = time(NULL);
318 int pid, numSeqsPerProcessor;
320 vector<unsigned long long> MPIPos;
321 MPIWroteAccnos = false;
324 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
325 MPI_Comm_size(MPI_COMM_WORLD, &processors);
328 MPI_File outMPIAlign;
329 MPI_File outMPIReport;
330 MPI_File outMPIAccnos;
332 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
333 int inMode=MPI_MODE_RDONLY;
335 char outAlignFilename[1024];
336 strcpy(outAlignFilename, alignFileName.c_str());
338 char outReportFilename[1024];
339 strcpy(outReportFilename, reportFileName.c_str());
341 char outAccnosFilename[1024];
342 strcpy(outAccnosFilename, accnosFileName.c_str());
344 char inFileName[1024];
345 strcpy(inFileName, candidateFileNames[s].c_str());
347 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
348 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
349 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
350 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
352 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
354 if (pid == 0) { //you are the root process
356 MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
358 //send file positions to all processes
359 for(int i = 1; i < processors; i++) {
360 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
361 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
364 //figure out how many sequences you have to align
365 numSeqsPerProcessor = numFastaSeqs / processors;
366 int startIndex = pid * numSeqsPerProcessor;
367 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
370 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
372 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
374 for (int i = 1; i < processors; i++) {
376 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
377 if (tempResult != 0) { MPIWroteAccnos = true; }
379 }else{ //you are a child process
380 MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
381 MPIPos.resize(numFastaSeqs+1);
382 MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
385 //figure out how many sequences you have to align
386 numSeqsPerProcessor = numFastaSeqs / processors;
387 int startIndex = pid * numSeqsPerProcessor;
388 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
392 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
394 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPIAlign); MPI_File_close(&outMPIReport); MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
396 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
400 MPI_File_close(&inMPI);
401 MPI_File_close(&outMPIAlign);
402 MPI_File_close(&outMPIReport);
403 MPI_File_close(&outMPIAccnos);
405 //delete accnos file if blank
407 //delete accnos file if its blank else report to user
408 if (MPIWroteAccnos) {
409 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
411 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
412 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
413 m->mothurOutEndLine();
416 //MPI_File_delete(outAccnosFilename, info);
418 m->mothurRemove(accnosFileName);
424 vector<unsigned long long> positions;
425 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
426 positions = m->divideFile(candidateFileNames[s], processors);
427 for (int i = 0; i < (positions.size()-1); i++) { lines.push_back(new linePair(positions[i], positions[(i+1)])); }
429 if (processors == 1) {
430 lines.push_back(new linePair(0, 1000));
432 positions = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs);
433 if (positions.size() < processors) { processors = positions.size(); }
435 //figure out how many sequences you have to process
436 int numSeqsPerProcessor = numFastaSeqs / processors;
437 for (int i = 0; i < processors; i++) {
438 int startIndex = i * numSeqsPerProcessor;
439 if(i == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor; }
440 lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
446 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
448 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
451 if (m->control_pressed) { m->mothurRemove(accnosFileName); m->mothurRemove(alignFileName); m->mothurRemove(reportFileName); outputTypes.clear(); return 0; }
453 //delete accnos file if its blank else report to user
454 if (m->isBlank(accnosFileName)) { m->mothurRemove(accnosFileName); hasAccnos = false; }
456 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
458 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well.");
459 }else{ m->mothurOut(" If the reverse compliment proved to be better it was reported."); }
460 m->mothurOutEndLine();
467 MPI_Comm_rank(MPI_COMM_WORLD, &pid);
469 if (pid == 0) { //only one process should output to screen
472 outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
473 outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
474 if (hasAccnos) { outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName); }
480 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
481 m->mothurOutEndLine();
482 m->mothurOutEndLine();
485 //set align file as new current fastafile
486 string currentFasta = "";
487 itTypes = outputTypes.find("fasta");
488 if (itTypes != outputTypes.end()) {
489 if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; m->setFastaFile(currentFasta); }
492 m->mothurOutEndLine();
493 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
494 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
495 m->mothurOutEndLine();
499 catch(exception& e) {
500 m->errorOut(e, "AlignCommand", "execute");
505 //**********************************************************************************************************************
506 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
508 ofstream alignmentFile;
509 m->openOutputFile(alignFName, alignmentFile);
512 m->openOutputFile(accnosFName, accnosFile);
514 NastReport report(reportFName);
517 m->openInputFile(filename, inFASTA);
519 inFASTA.seekg(filePos->start);
524 //moved this into driver to avoid deep copies in windows paralellized version
525 Alignment* alignment;
526 int longestBase = templateDB->getLongestBase();
527 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
528 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
529 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
530 else if(align == "noalign") { alignment = new NoAlign(); }
532 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
533 m->mothurOutEndLine();
534 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
539 if (m->control_pressed) { break; }
541 Sequence* candidateSeq = new Sequence(inFASTA); m->gobble(inFASTA);
542 report.setCandidate(candidateSeq);
544 int origNumBases = candidateSeq->getNumBases();
545 string originalUnaligned = candidateSeq->getUnaligned();
546 int numBasesNeeded = origNumBases * threshold;
548 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
549 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
550 alignment->resize(candidateSeq->getUnaligned().length()+1);
553 Sequence temp = templateDB->findClosestSequence(candidateSeq);
554 Sequence* templateSeq = &temp;
556 float searchScore = templateDB->getSearchScore();
558 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
563 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
564 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
565 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
566 //so this bool tells you if you need to delete it
568 //if there is a possibility that this sequence should be reversed
569 if (candidateSeq->getNumBases() < numBasesNeeded) {
571 string wasBetter = "";
572 //if the user wants you to try the reverse
575 //get reverse compliment
576 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
577 copy->reverseComplement();
580 Sequence temp2 = templateDB->findClosestSequence(copy);
581 Sequence* templateSeq2 = &temp2;
583 searchScore = templateDB->getSearchScore();
585 nast2 = new Nast(alignment, copy, templateSeq2);
587 //check if any better
588 if (copy->getNumBases() > candidateSeq->getNumBases()) {
589 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
590 templateSeq = templateSeq2;
593 needToDeleteCopy = true;
594 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
596 wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
602 //create accnos file with names
603 accnosFile << candidateSeq->getName() << wasBetter << endl;
606 report.setTemplate(templateSeq);
607 report.setSearchParameters(search, searchScore);
608 report.setAlignmentParameters(align, alignment);
609 report.setNastParameters(*nast);
611 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
615 if (needToDeleteCopy) { delete copy; }
621 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
622 unsigned long long pos = inFASTA.tellg();
623 if ((pos == -1) || (pos >= filePos->end)) { break; }
625 if (inFASTA.eof()) { break; }
629 if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
633 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine(); }
636 alignmentFile.close();
642 catch(exception& e) {
643 m->errorOut(e, "AlignCommand", "driver");
647 //**********************************************************************************************************************
649 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long long>& MPIPos){
651 string outputString = "";
652 MPI_Status statusReport;
653 MPI_Status statusAlign;
654 MPI_Status statusAccnos;
657 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
662 outputString = report.getHeaders();
663 int length = outputString.length();
665 char* buf = new char[length];
666 memcpy(buf, outputString.c_str(), length);
668 MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
673 Alignment* alignment;
674 int longestBase = templateDB->getLongestBase();
675 if(align == "gotoh") { alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase); }
676 else if(align == "needleman") { alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase); }
677 else if(align == "blast") { alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch); }
678 else if(align == "noalign") { alignment = new NoAlign(); }
680 m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
681 m->mothurOutEndLine();
682 alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
686 for(int i=0;i<num;i++){
688 if (m->control_pressed) { delete alignment; return 0; }
691 int length = MPIPos[start+i+1] - MPIPos[start+i];
693 char* buf4 = new char[length];
694 //memcpy(buf4, outputString.c_str(), length);
696 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
698 string tempBuf = buf4;
702 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
704 istringstream iss (tempBuf,istringstream::in);
706 Sequence* candidateSeq = new Sequence(iss);
707 report.setCandidate(candidateSeq);
709 int origNumBases = candidateSeq->getNumBases();
710 string originalUnaligned = candidateSeq->getUnaligned();
711 int numBasesNeeded = origNumBases * threshold;
713 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
714 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
715 alignment->resize(candidateSeq->getUnaligned().length()+1);
718 Sequence temp = templateDB->findClosestSequence(candidateSeq);
719 Sequence* templateSeq = &temp;
721 float searchScore = templateDB->getSearchScore();
723 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
727 bool needToDeleteCopy = false; //this is needed in case you have you enter the ifs below
728 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
729 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
730 //so this bool tells you if you need to delete it
732 //if there is a possibility that this sequence should be reversed
733 if (candidateSeq->getNumBases() < numBasesNeeded) {
735 string wasBetter = "";
736 //if the user wants you to try the reverse
738 //get reverse compliment
739 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
740 copy->reverseComplement();
743 Sequence temp2 = templateDB->findClosestSequence(copy);
744 Sequence* templateSeq2 = &temp2;
746 searchScore = templateDB->getSearchScore();
748 nast2 = new Nast(alignment, copy, templateSeq2);
750 //check if any better
751 if (copy->getNumBases() > candidateSeq->getNumBases()) {
752 candidateSeq->setAligned(copy->getAligned()); //use reverse compliments alignment since its better
753 templateSeq = templateSeq2;
756 needToDeleteCopy = true;
757 wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
759 wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
765 //create accnos file with names
766 outputString = candidateSeq->getName() + wasBetter + "\n";
768 //send results to parent
769 int length = outputString.length();
771 char* buf = new char[length];
772 memcpy(buf, outputString.c_str(), length);
774 MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
776 MPIWroteAccnos = true;
779 report.setTemplate(templateSeq);
780 report.setSearchParameters(search, searchScore);
781 report.setAlignmentParameters(align, alignment);
782 report.setNastParameters(*nast);
784 outputString = ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
786 //send results to parent
787 int length = outputString.length();
788 char* buf2 = new char[length];
789 memcpy(buf2, outputString.c_str(), length);
791 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
795 outputString = report.getReport();
797 //send results to parent
798 length = outputString.length();
799 char* buf3 = new char[length];
800 memcpy(buf3, outputString.c_str(), length);
802 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
806 if (needToDeleteCopy) { delete copy; }
811 if((i+1) % 100 == 0){ cout << (toString(i+1)) << endl; }
814 if((num) % 100 != 0){ cout << (toString(num)) << endl; }
818 catch(exception& e) {
819 m->errorOut(e, "AlignCommand", "driverMPI");
824 /**************************************************************************************************/
826 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
829 processIDS.resize(0);
830 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
833 //loop through and create all the processes you want
834 while (process != processors) {
838 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
841 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
843 //pass numSeqs to parent
845 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
846 m->openOutputFile(tempFile, out);
852 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine();
853 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
859 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
861 //force parent to wait until all the processes are done
862 for (int i=0;i<processIDS.size();i++) {
863 int temp = processIDS[i];
867 vector<string> nonBlankAccnosFiles;
868 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
869 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
871 for (int i = 0; i < processIDS.size(); i++) {
873 string tempFile = alignFileName + toString(processIDS[i]) + ".num.temp";
874 m->openInputFile(tempFile, in);
875 if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
876 in.close(); m->mothurRemove(tempFile);
878 appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
879 m->mothurRemove((alignFileName + toString(processIDS[i]) + ".temp"));
881 appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
882 m->mothurRemove((reportFileName + toString(processIDS[i]) + ".temp"));
884 if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
885 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
886 }else { m->mothurRemove((accnosFName + toString(processIDS[i]) + ".temp")); }
890 //append accnos files
891 if (nonBlankAccnosFiles.size() != 0) {
892 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
894 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
895 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
896 m->mothurRemove(nonBlankAccnosFiles[h]);
898 }else { //recreate the accnosfile if needed
900 m->openOutputFile(accnosFName, out);
904 //////////////////////////////////////////////////////////////////////////////////////////////////////
905 //Windows version shared memory, so be careful when passing variables through the alignData struct.
906 //Above fork() will clone, so memory is separate, but that's not the case with windows,
907 //////////////////////////////////////////////////////////////////////////////////////////////////////
909 vector<alignData*> pDataArray;
910 DWORD dwThreadIdArray[processors-1];
911 HANDLE hThreadArray[processors-1];
913 //Create processor worker threads.
914 for( int i=0; i<processors-1; i++ ){
916 //AlignmentDB* tempDB = new AlignmentDB(*templateDB);
918 // Allocate memory for thread data.
919 string extension = "";
920 if (i != 0) { extension = toString(i) + ".temp"; }
922 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);
923 pDataArray.push_back(tempalign);
924 processIDS.push_back(i);
926 //MySeqSumThreadFunction is in header. It must be global or static to work with the threads.
927 //default security attributes, thread function name, argument to thread function, use default creation flags, returns the thread identifier
928 hThreadArray[i] = CreateThread(NULL, 0, MyAlignThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
931 //need to check for line ending error
933 m->openInputFile(filename, inFASTA);
934 inFASTA.seekg(lines[processors-1]->start-1);
935 char c = inFASTA.peek();
937 if (c != '>') { //we need to move back
938 lines[processors-1]->start--;
941 //using the main process as a worker saves time and memory
942 //do my part - do last piece because windows is looking for eof
943 num = driver(lines[processors-1], (alignFileName + toString(processors-1) + ".temp"), (reportFileName + toString(processors-1) + ".temp"), (accnosFName + toString(processors-1) + ".temp"), filename);
945 //Wait until all threads have terminated.
946 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
948 //Close all thread handles and free memory allocations.
949 for(int i=0; i < pDataArray.size(); i++){
950 num += pDataArray[i]->count;
951 CloseHandle(hThreadArray[i]);
952 delete pDataArray[i];
955 vector<string> nonBlankAccnosFiles;
956 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
957 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
959 for (int i = 1; i < processors; i++) {
960 appendAlignFiles((alignFileName + toString(i) + ".temp"), alignFileName);
961 m->mothurRemove((alignFileName + toString(i) + ".temp"));
963 appendReportFiles((reportFileName + toString(i) + ".temp"), reportFileName);
964 m->mothurRemove((reportFileName + toString(i) + ".temp"));
966 if (!(m->isBlank(accnosFName + toString(i) + ".temp"))) {
967 nonBlankAccnosFiles.push_back(accnosFName + toString(i) + ".temp");
968 }else { m->mothurRemove((accnosFName + toString(i) + ".temp")); }
971 //append accnos files
972 if (nonBlankAccnosFiles.size() != 0) {
973 rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
975 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
976 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
977 m->mothurRemove(nonBlankAccnosFiles[h]);
979 }else { //recreate the accnosfile if needed
981 m->openOutputFile(accnosFName, out);
988 catch(exception& e) {
989 m->errorOut(e, "AlignCommand", "createProcesses");
993 /**************************************************************************************************/
995 void AlignCommand::appendAlignFiles(string temp, string filename) {
1000 m->openOutputFileAppend(filename, output);
1001 m->openInputFile(temp, input);
1003 while(char c = input.get()){
1004 if(input.eof()) { break; }
1005 else { output << c; }
1011 catch(exception& e) {
1012 m->errorOut(e, "AlignCommand", "appendAlignFiles");
1016 //**********************************************************************************************************************
1018 void AlignCommand::appendReportFiles(string temp, string filename) {
1023 m->openOutputFileAppend(filename, output);
1024 m->openInputFile(temp, input);
1026 while (!input.eof()) { char c = input.get(); if (c == 10 || c == 13){ break; } } // get header line
1028 while(char c = input.get()){
1029 if(input.eof()) { break; }
1030 else { output << c; }
1036 catch(exception& e) {
1037 m->errorOut(e, "AlignCommand", "appendReportFiles");
1041 //**********************************************************************************************************************