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