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