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