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