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