]> git.donarmstrong.com Git - mothur.git/blob - aligncommand.cpp
added getCommandInfoCommand for gui
[mothur.git] / aligncommand.cpp
1 /*
2  *  aligncommand.cpp
3  *  Mothur
4  *
5  *  Created by Sarah Westcott on 5/15/09.
6  *  Copyright 2009 Schloss Lab UMASS Amherst. All rights reserved.
7  *
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).
14  *
15  */
16
17 #include "aligncommand.h"
18 #include "sequence.hpp"
19
20 #include "gotohoverlap.hpp"
21 #include "needlemanoverlap.hpp"
22 #include "blastalign.hpp"
23 #include "noalign.hpp"
24
25 #include "nast.hpp"
26 #include "nastreport.hpp"
27
28 //**********************************************************************************************************************
29 vector<string> AlignCommand::setParameters(){   
30         try {
31                 CommandParameter ptemplate("reference", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(ptemplate);
32                 CommandParameter pcandidate("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pcandidate);
33                 CommandParameter psearch("search", "Multiple", "kmer-blast-suffix", "kmer", "", "", "",false,false); parameters.push_back(psearch);
34                 CommandParameter pksize("ksize", "Number", "", "8", "", "", "",false,false); parameters.push_back(pksize);
35                 CommandParameter pmatch("match", "Number", "", "1.0", "", "", "",false,false); parameters.push_back(pmatch);
36                 CommandParameter palign("align", "Multiple", "needleman-gotoh-blast-noalign", "needleman", "", "", "",false,false); parameters.push_back(palign);
37                 CommandParameter pmismatch("mismatch", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pmismatch);
38                 CommandParameter pgapopen("gapopen", "Number", "", "-2.0", "", "", "",false,false); parameters.push_back(pgapopen);
39                 CommandParameter pgapextend("gapextend", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pgapextend);
40                 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
41                 CommandParameter pflip("flip", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(pflip);
42                 CommandParameter pthreshold("threshold", "Number", "", "0.50", "", "", "",false,false); parameters.push_back(pthreshold);
43                 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
44                 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
45
46                 vector<string> myArray;
47                 for (int i = 0; i < parameters.size(); i++) {   myArray.push_back(parameters[i].name);          }
48                 return myArray;
49         }
50         catch(exception& e) {
51                 m->errorOut(e, "AlignCommand", "setParameters");
52                 exit(1);
53         }
54 }
55 //**********************************************************************************************************************
56 string AlignCommand::getHelpString(){   
57         try {
58                 string helpString = "";
59                 helpString += "The align.seqs command reads a file containing sequences and creates an alignment file and a report file.";
60                 helpString += "The align.seqs command parameters are reference, fasta, search, ksize, align, match, mismatch, gapopen, gapextend and processors.";
61                 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.";
62                 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.";
63                 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.";
64                 helpString += "The ksize parameter allows you to specify the kmer size for finding most similar template to candidate.  The default is 8.";
65                 helpString += "The match parameter allows you to specify the bonus for having the same base. The default is 1.0.";
66                 helpString += "The mistmatch parameter allows you to specify the penalty for having different bases.  The default is -1.0.";
67                 helpString += "The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.";
68                 helpString += "The gapextend parameter allows you to specify the penalty for extending a gap in an alignment.  The default is -1.0.";
69                 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.";
70                 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.";
71                 helpString += "If the flip parameter is set to true the reverse complement of the sequence is aligned and the better alignment is reported.";
72                 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.";
73                 helpString += "The align.seqs command should be in the following format:";
74                 helpString += "align.seqs(reference=yourTemplateFile, fasta=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty)";
75                 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)";
76                 helpString += "Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).";
77                 return helpString;
78         }
79         catch(exception& e) {
80                 m->errorOut(e, "AlignCommand", "getHelpString");
81                 exit(1);
82         }
83 }
84 //**********************************************************************************************************************
85 AlignCommand::AlignCommand(){   
86         try {
87                 abort = true; calledHelp = true; 
88                 setParameters();
89                 vector<string> tempOutNames;
90                 outputTypes["fasta"] = tempOutNames;
91                 outputTypes["alignreport"] = tempOutNames;
92                 outputTypes["accnos"] = tempOutNames;
93         }
94         catch(exception& e) {
95                 m->errorOut(e, "AlignCommand", "AlignCommand");
96                 exit(1);
97         }
98 }
99 //**********************************************************************************************************************
100 AlignCommand::AlignCommand(string option)  {
101         try {
102                 abort = false; calledHelp = false;   
103         
104                 //allow user to run help
105                 if(option == "help") { help(); abort = true; calledHelp = true;}
106                 
107                 else {
108                         vector<string> myArray = setParameters();
109                         
110                         OptionParser parser(option);
111                         map<string, string> parameters = parser.getParameters(); 
112                         
113                         ValidParameters validParameter("align.seqs");
114                         map<string, string>::iterator it;
115                         
116                         //check to make sure all parameters are valid for command
117                         for (it = parameters.begin(); it != parameters.end(); it++) { 
118                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
119                         }
120                         
121                         //initialize outputTypes
122                         vector<string> tempOutNames;
123                         outputTypes["fasta"] = tempOutNames;
124                         outputTypes["alignreport"] = tempOutNames;
125                         outputTypes["accnos"] = tempOutNames;
126                         
127                         //if the user changes the output directory command factory will send this info to us in the output parameter 
128                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  outputDir = "";         }
129                         
130
131                         //if the user changes the input directory command factory will send this info to us in the output parameter 
132                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
133                         
134                         if (inputDir == "not found"){   inputDir = "";          }
135                         else {
136                                 string path;
137
138                                 it = parameters.find("reference");
139
140                                 //user has given a template file
141                                 if(it != parameters.end()){ 
142                                         path = m->hasPath(it->second);
143                                         //if the user has not given a path then, add inputdir. else leave path alone.
144                                         if (path == "") {       parameters["reference"] = inputDir + it->second;                }
145                                 }
146                         }
147
148                         //check for required parameters
149                         templateFileName = validParameter.validFile(parameters, "reference", true);
150                         
151                         if (templateFileName == "not found") { 
152                                 m->mothurOut("reference is a required parameter for the align.seqs command."); 
153                                 m->mothurOutEndLine();
154                                 abort = true; 
155                         }else if (templateFileName == "not open") { abort = true; }     
156                         
157                         candidateFileName = validParameter.validFile(parameters, "fasta", false);
158                         if (candidateFileName == "not found") { 
159                                 //if there is a current fasta file, use it
160                                 string filename = m->getFastaFile(); 
161                                 if (filename != "") { candidateFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
162                                 else {  m->mothurOut("You have no current fastafile and the candidate parameter is required."); m->mothurOutEndLine(); abort = true; }
163                         }else { 
164                                 m->splitAtDash(candidateFileName, candidateFileNames);
165                                 
166                                 //go through files and make sure they are good, if not, then disregard them
167                                 for (int i = 0; i < candidateFileNames.size(); i++) {
168                                         //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
169                                         
170                                         if (inputDir != "") {
171                                                 string path = m->hasPath(candidateFileNames[i]);
172                                                 //if the user has not given a path then, add inputdir. else leave path alone.
173                                                 if (path == "") {       candidateFileNames[i] = inputDir + candidateFileNames[i];               }
174                                         }
175         
176                                         int ableToOpen;
177                                         ifstream in;
178                                         ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
179                                         in.close();     
180                                         
181                                         //if you can't open it, try default location
182                                         if (ableToOpen == 1) {
183                                                 if (m->getDefaultPath() != "") { //default path is set
184                                                         string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
185                                                         m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
186                                                         ifstream in2;
187                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
188                                                         in2.close();
189                                                         candidateFileNames[i] = tryPath;
190                                                 }
191                                         }
192                                         
193                                         //if you can't open it, try output location
194                                         if (ableToOpen == 1) {
195                                                 if (m->getOutputDir() != "") { //default path is set
196                                                         string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
197                                                         m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
198                                                         ifstream in2;
199                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
200                                                         in2.close();
201                                                         candidateFileNames[i] = tryPath;
202                                                 }
203                                         }
204                                         
205                                                                         
206
207                                         if (ableToOpen == 1) { 
208                                                 m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
209                                                 //erase from file list
210                                                 candidateFileNames.erase(candidateFileNames.begin()+i);
211                                                 i--;
212                                         }
213                                         
214                                 }
215                                 
216                                 //make sure there is at least one valid file left
217                                 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
218                         }
219                 
220                         //check for optional parameter and set defaults
221                         // ...at some point should added some additional type checking...
222                         string temp;
223                         temp = validParameter.validFile(parameters, "ksize", false);            if (temp == "not found"){       temp = "8";                             }
224                         convert(temp, kmerSize); 
225                         
226                         temp = validParameter.validFile(parameters, "match", false);            if (temp == "not found"){       temp = "1.0";                   }
227                         convert(temp, match);  
228                         
229                         temp = validParameter.validFile(parameters, "mismatch", false);         if (temp == "not found"){       temp = "-1.0";                  }
230                         convert(temp, misMatch);  
231                         
232                         temp = validParameter.validFile(parameters, "gapopen", false);          if (temp == "not found"){       temp = "-2.0";                  }
233                         convert(temp, gapOpen);  
234                         
235                         temp = validParameter.validFile(parameters, "gapextend", false);        if (temp == "not found"){       temp = "-1.0";                  }
236                         convert(temp, gapExtend); 
237                         
238                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = m->getProcessors();      }
239                         m->setProcessors(temp);
240                         convert(temp, processors); 
241                         
242                         temp = validParameter.validFile(parameters, "flip", false);                     if (temp == "not found"){       temp = "f";                             }
243                         flip = m->isTrue(temp); 
244                         
245                         temp = validParameter.validFile(parameters, "threshold", false);        if (temp == "not found"){       temp = "0.50";                  }
246                         convert(temp, threshold); 
247                         
248                         search = validParameter.validFile(parameters, "search", false);         if (search == "not found"){     search = "kmer";                }
249                         if ((search != "suffix") && (search != "kmer") && (search != "blast")) { m->mothurOut("invalid search option: choices are kmer, suffix or blast."); m->mothurOutEndLine(); abort=true; }
250                         
251                         align = validParameter.validFile(parameters, "align", false);           if (align == "not found"){      align = "needleman";    }
252                         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; }
253
254                 }
255                 
256         }
257         catch(exception& e) {
258                 m->errorOut(e, "AlignCommand", "AlignCommand");
259                 exit(1);
260         }
261 }
262 //**********************************************************************************************************************
263 AlignCommand::~AlignCommand(){  
264
265         if (abort == false) {
266                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
267                 delete templateDB;
268                 delete alignment;
269         }
270 }
271 //**********************************************************************************************************************
272
273 int AlignCommand::execute(){
274         try {
275                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
276
277                 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch);
278                 int longestBase = templateDB->getLongestBase();
279                 
280                 if(align == "gotoh")                    {       alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase);                 }
281                 else if(align == "needleman")   {       alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);                                }
282                 else if(align == "blast")               {       alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch);            }
283                 else if(align == "noalign")             {       alignment = new NoAlign();                                                                                                      }
284                 else {
285                         m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
286                         m->mothurOutEndLine();
287                         alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
288                 }
289                 
290                 for (int s = 0; s < candidateFileNames.size(); s++) {
291                         if (m->control_pressed) { outputTypes.clear(); return 0; }
292                         
293                         m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
294                         
295                         if (outputDir == "") {  outputDir += m->hasPath(candidateFileNames[s]); }
296                         string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align";
297                         string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "align.report";
298                         string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + "flip.accnos";
299                         bool hasAccnos = true;
300                         
301                         int numFastaSeqs = 0;
302                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
303                         int start = time(NULL);
304                 
305 #ifdef USE_MPI  
306                                 int pid, numSeqsPerProcessor; 
307                                 int tag = 2001;
308                                 vector<unsigned long int> MPIPos;
309                                 MPIWroteAccnos = false;
310                         
311                                 MPI_Status status; 
312                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
313                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); 
314
315                                 MPI_File inMPI;
316                                 MPI_File outMPIAlign;
317                                 MPI_File outMPIReport;
318                                 MPI_File outMPIAccnos;
319                                 
320                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
321                                 int inMode=MPI_MODE_RDONLY; 
322                                 
323                                 char outAlignFilename[1024];
324                                 strcpy(outAlignFilename, alignFileName.c_str());
325                                 
326                                 char outReportFilename[1024];
327                                 strcpy(outReportFilename, reportFileName.c_str());
328                                 
329                                 char outAccnosFilename[1024];
330                                 strcpy(outAccnosFilename, accnosFileName.c_str());
331                                 
332                                 char inFileName[1024];
333                                 strcpy(inFileName, candidateFileNames[s].c_str());
334                                 
335                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
336                                 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
337                                 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
338                                 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
339                                 
340                                 if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIAlign);  MPI_File_close(&outMPIReport);  MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
341                                 
342                                 if (pid == 0) { //you are the root process 
343                                         
344                                         MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
345                                         
346                                         //send file positions to all processes
347                                         for(int i = 1; i < processors; i++) { 
348                                                 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
349                                                 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
350                                         }
351                                         
352                                         //figure out how many sequences you have to align
353                                         numSeqsPerProcessor = numFastaSeqs / processors;
354                                         int startIndex =  pid * numSeqsPerProcessor;
355                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
356                                         
357                                         //align your part
358                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
359                                         
360                                         if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIAlign);  MPI_File_close(&outMPIReport);  MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
361
362                                         for (int i = 1; i < processors; i++) {
363                                                 bool tempResult;
364                                                 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
365                                                 if (tempResult != 0) { MPIWroteAccnos = true; }
366                                         }
367                                 }else{ //you are a child process
368                                         MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
369                                         MPIPos.resize(numFastaSeqs+1);
370                                         MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
371
372                                         
373                                         //figure out how many sequences you have to align
374                                         numSeqsPerProcessor = numFastaSeqs / processors;
375                                         int startIndex =  pid * numSeqsPerProcessor;
376                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
377                                         
378                                         
379                                         //align your part
380                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
381                                         
382                                         if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIAlign);  MPI_File_close(&outMPIReport);  MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
383
384                                         MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD); 
385                                 }
386                                 
387                                 //close files 
388                                 MPI_File_close(&inMPI);
389                                 MPI_File_close(&outMPIAlign);
390                                 MPI_File_close(&outMPIReport);
391                                 MPI_File_close(&outMPIAccnos);
392                                 
393                                 //delete accnos file if blank
394                                 if (pid == 0) {
395                                         //delete accnos file if its blank else report to user
396                                         if (MPIWroteAccnos) { 
397                                                 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
398                                                 if (!flip) {
399                                                         m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
400                                                 }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
401                                                 m->mothurOutEndLine();
402                                         }else { 
403                                                 //MPI_Info info;
404                                                 //MPI_File_delete(outAccnosFilename, info);
405                                                 hasAccnos = false;      
406                                                 remove(accnosFileName.c_str()); 
407                                         }
408                                 }
409                                 
410 #else
411
412                 vector<unsigned long int> positions = m->divideFile(candidateFileNames[s], processors);
413                 for (int i = 0; i < (positions.size()-1); i++) {
414                         lines.push_back(new linePair(positions[i], positions[(i+1)]));
415                 }       
416         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
417                         if(processors == 1){
418                                 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
419                         }else{
420                                 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]); 
421                         }
422         #else
423                         numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
424         #endif
425                         
426                         if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); outputTypes.clear();  return 0; }
427                         
428                         //delete accnos file if its blank else report to user
429                         if (m->isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
430                         else { 
431                                 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
432                                 if (!flip) {
433                                         m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
434                                 }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
435                                 m->mothurOutEndLine();
436                         }
437
438 #endif          
439
440
441                 #ifdef USE_MPI
442                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
443                                         
444                         if (pid == 0) { //only one process should output to screen
445                 #endif
446
447                         outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
448                         outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
449                         if (hasAccnos)  {       outputNames.push_back(accnosFileName);  outputTypes["accnos"].push_back(accnosFileName);  }
450                         
451                 #ifdef USE_MPI
452                         }
453                 #endif
454
455                         m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
456                         m->mothurOutEndLine();
457                         m->mothurOutEndLine();
458                 }
459                 
460                 //set align file as new current fastafile
461                 string currentFasta = "";
462                 itTypes = outputTypes.find("fasta");
463                 if (itTypes != outputTypes.end()) {
464                         if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; m->setFastaFile(currentFasta); }
465                 }
466                 
467                 m->mothurOutEndLine();
468                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
469                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }
470                 m->mothurOutEndLine();
471
472                 return 0;
473         }
474         catch(exception& e) {
475                 m->errorOut(e, "AlignCommand", "execute");
476                 exit(1);
477         }
478 }
479
480 //**********************************************************************************************************************
481
482 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
483         try {
484                 ofstream alignmentFile;
485                 m->openOutputFile(alignFName, alignmentFile);
486                 
487                 ofstream accnosFile;
488                 m->openOutputFile(accnosFName, accnosFile);
489                 
490                 NastReport report(reportFName);
491                 
492                 ifstream inFASTA;
493                 m->openInputFile(filename, inFASTA);
494
495                 inFASTA.seekg(filePos->start);
496
497                 bool done = false;
498                 int count = 0;
499         
500                 while (!done) {
501                         
502                         if (m->control_pressed) {  return 0; }
503                         
504                         Sequence* candidateSeq = new Sequence(inFASTA);  m->gobble(inFASTA);
505                         report.setCandidate(candidateSeq);
506
507                         int origNumBases = candidateSeq->getNumBases();
508                         string originalUnaligned = candidateSeq->getUnaligned();
509                         int numBasesNeeded = origNumBases * threshold;
510         
511                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
512                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
513                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
514                                 }
515                                                                 
516                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
517                                 Sequence* templateSeq = &temp;
518                         
519                                 float searchScore = templateDB->getSearchScore();
520                                                                 
521                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
522                 
523                                 Sequence* copy;
524                                 
525                                 Nast* nast2;
526                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
527                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
528                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
529                                                                                                 //so this bool tells you if you need to delete it
530                                                                                                 
531                                 //if there is a possibility that this sequence should be reversed
532                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
533                                         
534                                         string wasBetter =  "";
535                                         //if the user wants you to try the reverse
536                                         if (flip) {
537                                 
538                                                 //get reverse compliment
539                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
540                                                 copy->reverseComplement();
541                                                 
542                                                 //rerun alignment
543                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
544                                                 Sequence* templateSeq2 = &temp2;
545                                                 
546                                                 searchScore = templateDB->getSearchScore();
547                                                 
548                                                 nast2 = new Nast(alignment, copy, templateSeq2);
549                         
550                                                 //check if any better
551                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
552                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
553                                                         templateSeq = templateSeq2; 
554                                                         delete nast;
555                                                         nast = nast2;
556                                                         needToDeleteCopy = true;
557                                                         wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
558                                                 }else{  
559                                                         wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
560                                                         delete nast2;
561                                                         delete copy;    
562                                                 }
563                                         }
564                                         
565                                         //create accnos file with names
566                                         accnosFile << candidateSeq->getName() << wasBetter << endl;
567                                 }
568                                 
569                                 report.setTemplate(templateSeq);
570                                 report.setSearchParameters(search, searchScore);
571                                 report.setAlignmentParameters(align, alignment);
572                                 report.setNastParameters(*nast);
573         
574                                 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
575                                 
576                                 report.print();
577                                 delete nast;
578                                 if (needToDeleteCopy) {   delete copy;   }
579                                 
580                                 count++;
581                         }
582                         delete candidateSeq;
583                         
584                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
585                                 unsigned long int pos = inFASTA.tellg();
586                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
587                         #else
588                                 if (inFASTA.eof()) { break; }
589                         #endif
590                         
591                         //report progress
592                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
593                         
594                 }
595                 //report progress
596                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
597                 
598                 alignmentFile.close();
599                 inFASTA.close();
600                 accnosFile.close();
601                 
602                 return count;
603         }
604         catch(exception& e) {
605                 m->errorOut(e, "AlignCommand", "driver");
606                 exit(1);
607         }
608 }
609 //**********************************************************************************************************************
610 #ifdef USE_MPI
611 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
612         try {
613                 string outputString = "";
614                 MPI_Status statusReport; 
615                 MPI_Status statusAlign; 
616                 MPI_Status statusAccnos; 
617                 MPI_Status status; 
618                 int pid;
619                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
620         
621                 NastReport report;
622                 
623                 if (pid == 0) {
624                         outputString = report.getHeaders();
625                         int length = outputString.length();
626             
627                         char* buf = new char[length];
628                         memcpy(buf, outputString.c_str(), length);
629                 
630                         MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
631
632             delete buf;
633                 }
634                 
635                 for(int i=0;i<num;i++){
636                 
637                         if (m->control_pressed) {  return 0; }
638
639                         //read next sequence
640                         int length = MPIPos[start+i+1] - MPIPos[start+i];
641
642                         char* buf4 = new char[length];
643                         //memcpy(buf4, outputString.c_str(), length);
644
645                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
646                         
647                         string tempBuf = buf4;
648
649                         delete buf4;
650
651                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
652         
653                         istringstream iss (tempBuf,istringstream::in);
654
655                         Sequence* candidateSeq = new Sequence(iss);  
656                         report.setCandidate(candidateSeq);
657
658                         int origNumBases = candidateSeq->getNumBases();
659                         string originalUnaligned = candidateSeq->getUnaligned();
660                         int numBasesNeeded = origNumBases * threshold;
661         
662                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
663                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
664                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
665                                 }
666                                                                 
667                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
668                                 Sequence* templateSeq = &temp;
669                                 
670                                 float searchScore = templateDB->getSearchScore();
671                                                                 
672                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
673                                 Sequence* copy;
674                                 
675                                 Nast* nast2;
676                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
677                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
678                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
679                                                                                                 //so this bool tells you if you need to delete it
680                                                                                                 
681                                 //if there is a possibility that this sequence should be reversed
682                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
683                                         
684                                         string wasBetter = "";
685                                         //if the user wants you to try the reverse
686                                         if (flip) {
687                                                 //get reverse compliment
688                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
689                                                 copy->reverseComplement();
690                                                 
691                                                 //rerun alignment
692                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
693                                                 Sequence* templateSeq2 = &temp2;
694                                                 
695                                                 searchScore = templateDB->getSearchScore();
696                                                 
697                                                 nast2 = new Nast(alignment, copy, templateSeq2);
698                         
699                                                 //check if any better
700                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
701                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
702                                                         templateSeq = templateSeq2; 
703                                                         delete nast;
704                                                         nast = nast2;
705                                                         needToDeleteCopy = true;
706                                                         wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
707                                                 }else{  
708                                                         wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
709                                                         delete nast2;
710                                                         delete copy;    
711                                                 }
712                                         }
713                                         
714                                         //create accnos file with names
715                                         outputString = candidateSeq->getName() + wasBetter + "\n";
716                                         
717                                         //send results to parent
718                                         int length = outputString.length();
719
720                                         char* buf = new char[length];
721                                         memcpy(buf, outputString.c_str(), length);
722                                 
723                                         MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
724                                         delete buf;
725                                         MPIWroteAccnos = true;
726                                 }
727                                 
728                                 report.setTemplate(templateSeq);
729                                 report.setSearchParameters(search, searchScore);
730                                 report.setAlignmentParameters(align, alignment);
731                                 report.setNastParameters(*nast);
732         
733                                 outputString =  ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
734                                 
735                                 //send results to parent
736                                 int length = outputString.length();
737                                 char* buf2 = new char[length];
738                                 memcpy(buf2, outputString.c_str(), length);
739                                 
740                                 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
741                                 
742                                 delete buf2;
743
744                                 outputString = report.getReport();
745                                 
746                                 //send results to parent
747                                 length = outputString.length();
748                                 char* buf3 = new char[length];
749                                 memcpy(buf3, outputString.c_str(), length);
750                                 
751                                 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
752                                 
753                                 delete buf3;
754                                 delete nast;
755                                 if (needToDeleteCopy) {   delete copy;   }
756                         }
757                         delete candidateSeq;
758                         
759                         //report progress
760                         if((i+1) % 100 == 0){   cout << (toString(i+1)) << endl;                }
761                 }
762                 //report progress
763                 if((num) % 100 != 0){   cout << (toString(num)) << endl;                }
764                 
765                 return 1;
766         }
767         catch(exception& e) {
768                 m->errorOut(e, "AlignCommand", "driverMPI");
769                 exit(1);
770         }
771 }
772 #endif
773 /**************************************************************************************************/
774
775 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
776         try {
777 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
778                 processIDS.resize(0);
779                 int process = 1;
780                 int num = 0;
781                 //              processIDS.resize(0);
782                 
783                 //loop through and create all the processes you want
784                 while (process != processors) {
785                         int pid = fork();
786                         
787                         if (pid > 0) {
788                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
789                                 process++;
790                         }else if (pid == 0){
791                                 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
792                                 
793                                 //pass numSeqs to parent
794                                 ofstream out;
795                                 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
796                                 m->openOutputFile(tempFile, out);
797                                 out << num << endl;
798                                 out.close();
799                                 
800                                 exit(0);
801                         }else { 
802                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
803                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
804                                 exit(0);
805                         }
806                 }
807                 
808                 //do my part
809                 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
810                 
811                 //force parent to wait until all the processes are done
812                 for (int i=0;i<processIDS.size();i++) { 
813                         int temp = processIDS[i];
814                         wait(&temp);
815                 }
816                 
817                 vector<string> nonBlankAccnosFiles;
818                 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
819                 else { remove(accnosFName.c_str()); } //remove so other files can be renamed to it
820                         
821                 for (int i = 0; i < processIDS.size(); i++) {
822                         ifstream in;
823                         string tempFile =  alignFileName + toString(processIDS[i]) + ".num.temp";
824                         m->openInputFile(tempFile, in);
825                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
826                         in.close(); remove(tempFile.c_str());
827                         
828                         appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
829                         remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
830                         
831                         appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
832                         remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
833                         
834                         if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
835                                 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
836                         }else { remove((accnosFName + toString(processIDS[i]) + ".temp").c_str());  }
837                         
838                 }
839                 
840                 //append accnos files
841                 if (nonBlankAccnosFiles.size() != 0) { 
842                         rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
843                         
844                         for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
845                                 appendAlignFiles(nonBlankAccnosFiles[h], accnosFName);
846                                 remove(nonBlankAccnosFiles[h].c_str());
847                         }
848                 }else { //recreate the accnosfile if needed
849                         ofstream out;
850                         m->openOutputFile(accnosFName, out);
851                         out.close();
852                 }
853                 
854                 return num;
855 #endif          
856         }
857         catch(exception& e) {
858                 m->errorOut(e, "AlignCommand", "createProcesses");
859                 exit(1);
860         }
861 }
862 /**************************************************************************************************/
863
864 void AlignCommand::appendAlignFiles(string temp, string filename) {
865         try{
866                 
867                 ofstream output;
868                 ifstream input;
869                 m->openOutputFileAppend(filename, output);
870                 m->openInputFile(temp, input);
871                 
872                 while(char c = input.get()){
873                         if(input.eof())         {       break;                  }
874                         else                            {       output << c;    }
875                 }
876                 
877                 input.close();
878                 output.close();
879         }
880         catch(exception& e) {
881                 m->errorOut(e, "AlignCommand", "appendAlignFiles");
882                 exit(1);
883         }
884 }
885 //**********************************************************************************************************************
886
887 void AlignCommand::appendReportFiles(string temp, string filename) {
888         try{
889                 
890                 ofstream output;
891                 ifstream input;
892                 m->openOutputFileAppend(filename, output);
893                 m->openInputFile(temp, input);
894
895                 while (!input.eof())    {       char c = input.get(); if (c == 10 || c == 13){  break;  }       } // get header line
896                                 
897                 while(char c = input.get()){
898                         if(input.eof())         {       break;                  }
899                         else                            {       output << c;    }
900                 }
901                 
902                 input.close();
903                 output.close();
904         }
905         catch(exception& e) {
906                 m->errorOut(e, "AlignCommand", "appendReportFiles");
907                 exit(1);
908         }
909 }
910 //**********************************************************************************************************************