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