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