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