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