]> git.donarmstrong.com Git - mothur.git/blob - aligncommand.cpp
added MPI code, broke up chimera.seqs into 5 separated commands, added parse.sff...
[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[alignFileName.length()];
271                                 strcpy(outAlignFilename, alignFileName.c_str());
272                                 
273                                 char outReportFilename[reportFileName.length()];
274                                 strcpy(outReportFilename, reportFileName.c_str());
275                                 
276                                 char outAccnosFilename[accnosFileName.length()];
277                                 strcpy(outAccnosFilename, accnosFileName.c_str());
278                                 
279                                 char inFileName[candidateFileNames[s].length()];
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                                         MPI_Bcast(&numFastaSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD);  //send numSeqs
295                                         MPI_Bcast(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //send file pos   
296                                         
297                                         //figure out how many sequences you have to align
298                                         numSeqsPerProcessor = numFastaSeqs / processors;
299                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
300                                         int startIndex =  pid * numSeqsPerProcessor;
301                                 
302                                         //align your part
303                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
304                                         
305                                         if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIAlign);  MPI_File_close(&outMPIReport);  MPI_File_close(&outMPIAccnos); return 0; }
306
307                                         for (int i = 1; i < processors; i++) {
308                                                 bool tempResult;
309                                                 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
310                                                 if (tempResult != 0) { MPIWroteAccnos = true; }
311                                         }
312                                 }else{ //you are a child process
313                                         MPI_Bcast(&numFastaSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD); //get numSeqs
314                                         MPIPos.resize(numFastaSeqs+1);
315                                         MPI_Bcast(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //get file positions
316                                         
317                                         //figure out how many sequences you have to align
318                                         numSeqsPerProcessor = numFastaSeqs / processors;
319                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
320                                         int startIndex =  pid * numSeqsPerProcessor;
321                                         
322                                         //align your part
323                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
324                                         
325                                         if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIAlign);  MPI_File_close(&outMPIReport);  MPI_File_close(&outMPIAccnos); return 0; }
326
327                                         MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD); 
328                                 }
329                                 
330                                 //close files 
331                                 MPI_File_close(&inMPI);
332                                 MPI_File_close(&outMPIAlign);
333                                 MPI_File_close(&outMPIReport);
334                                 MPI_File_close(&outMPIAccnos);
335                                 
336                                 //delete accnos file if blank
337                                 if (pid == 0) {
338                                         //delete accnos file if its blank else report to user
339                                         if (MPIWroteAccnos) { 
340                                                 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
341                                                 if (!flip) {
342                                                         m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
343                                                 }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
344                                                 m->mothurOutEndLine();
345                                         }else { 
346                                                 //MPI_Info info;
347                                                 //MPI_File_delete(outAccnosFilename, info);
348                                                 hasAccnos = false;      
349                                                 remove(accnosFileName.c_str()); 
350                                         }
351                                 }
352                                 
353 #else
354                         
355         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
356                         if(processors == 1){
357                                 ifstream inFASTA;
358                                 openInputFile(candidateFileNames[s], inFASTA);
359                                 numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
360                                 inFASTA.close();
361                                 
362                                 lines.push_back(new linePair(0, numFastaSeqs));
363                                 
364                                 driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
365                                 
366                                 if (m->control_pressed) { 
367                                         remove(accnosFileName.c_str()); 
368                                         remove(alignFileName.c_str()); 
369                                         remove(reportFileName.c_str()); 
370                                         return 0; 
371                                 }
372                                 
373                                 //delete accnos file if its blank else report to user
374                                 if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
375                                 else { 
376                                         m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
377                                         if (!flip) {
378                                                 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
379                                         }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
380                                         m->mothurOutEndLine();
381                                 }
382                         }
383                         else{
384                                 vector<int> positions;
385                                 processIDS.resize(0);
386                                 
387                                 ifstream inFASTA;
388                                 openInputFile(candidateFileNames[s], inFASTA);
389                                 
390                                 string input;
391                                 while(!inFASTA.eof()){
392                                         input = getline(inFASTA);
393                                         if (input.length() != 0) {
394                                                 if(input[0] == '>'){    long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);  }
395                                         }
396                                 }
397                                 inFASTA.close();
398                                 
399                                 numFastaSeqs = positions.size();
400                                 
401                                 int numSeqsPerProcessor = numFastaSeqs / processors;
402                                 
403                                 for (int i = 0; i < processors; i++) {
404                                         long int startPos = positions[ i * numSeqsPerProcessor ];
405                                         if(i == processors - 1){
406                                                 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
407                                         }
408                                         lines.push_back(new linePair(startPos, numSeqsPerProcessor));
409                                 }
410                                 
411                                 createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]); 
412                                 
413                                 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
414                                 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
415                                 
416                                 //append alignment and report files
417                                 for(int i=1;i<processors;i++){
418                                         appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
419                                         remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
420                                         
421                                         appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
422                                         remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
423                                 }
424                                 
425                                 vector<string> nonBlankAccnosFiles;
426                                 //delete blank accnos files generated with multiple processes
427                                 for(int i=0;i<processors;i++){  
428                                         if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
429                                                 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
430                                         }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());  }
431                                 }
432                                 
433                                 //append accnos files
434                                 if (nonBlankAccnosFiles.size() != 0) { 
435                                         rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
436                                         
437                                         for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
438                                                 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
439                                                 remove(nonBlankAccnosFiles[h].c_str());
440                                         }
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                                 }else{ hasAccnos = false;  }
447                                 
448                                 if (m->control_pressed) { 
449                                         remove(accnosFileName.c_str()); 
450                                         remove(alignFileName.c_str()); 
451                                         remove(reportFileName.c_str()); 
452                                         return 0; 
453                                 }
454                         }
455         #else
456                         ifstream inFASTA;
457                         openInputFile(candidateFileNames[s], inFASTA);
458                         numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
459                         inFASTA.close();
460                         
461                         lines.push_back(new linePair(0, numFastaSeqs));
462                         
463                         driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
464                         
465                         if (m->control_pressed) { 
466                                 remove(accnosFileName.c_str()); 
467                                 remove(alignFileName.c_str()); 
468                                 remove(reportFileName.c_str()); 
469                                 return 0; 
470                         }
471                         
472                         //delete accnos file if its blank else report to user
473                         if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
474                         else { 
475                                 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
476                                 if (!flip) {
477                                         m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
478                                 }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
479                                 m->mothurOutEndLine();
480                         }
481                         
482         #endif
483
484 #endif          
485
486
487                 #ifdef USE_MPI
488                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
489                                         
490                         if (pid == 0) { //only one process should output to screen
491                 #endif
492
493                         outputNames.push_back(alignFileName);
494                         outputNames.push_back(reportFileName);
495                         if (hasAccnos)  {       outputNames.push_back(accnosFileName);          }
496                         
497                 #ifdef USE_MPI
498                         }
499                 #endif
500
501                         m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
502                         m->mothurOutEndLine();
503                         m->mothurOutEndLine();
504                 }
505                 
506                 
507                 m->mothurOutEndLine();
508                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
509                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }
510                 m->mothurOutEndLine();
511
512                 return 0;
513         }
514         catch(exception& e) {
515                 m->errorOut(e, "AlignCommand", "execute");
516                 exit(1);
517         }
518 }
519
520 //**********************************************************************************************************************
521
522 int AlignCommand::driver(linePair* line, string alignFName, string reportFName, string accnosFName, string filename){
523         try {
524                 ofstream alignmentFile;
525                 openOutputFile(alignFName, alignmentFile);
526                 
527                 ofstream accnosFile;
528                 openOutputFile(accnosFName, accnosFile);
529                 
530                 NastReport report(reportFName);
531                 
532                 ifstream inFASTA;
533                 openInputFile(filename, inFASTA);
534
535                 inFASTA.seekg(line->start);
536         
537                 for(int i=0;i<line->numSeqs;i++){
538                         
539                         if (m->control_pressed) {  return 0; }
540                         
541                         Sequence* candidateSeq = new Sequence(inFASTA);  gobble(inFASTA);
542         
543                         int origNumBases = candidateSeq->getNumBases();
544                         string originalUnaligned = candidateSeq->getUnaligned();
545                         int numBasesNeeded = origNumBases * threshold;
546         
547                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
548                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
549                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
550                                 }
551                                                                 
552                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
553                                 Sequence* templateSeq = &temp;
554                                 
555                                 float searchScore = templateDB->getSearchScore();
556                                                                 
557                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
558                                 Sequence* copy;
559                                 
560                                 Nast* nast2;
561                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
562                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
563                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
564                                                                                                 //so this bool tells you if you need to delete it
565                                                                                                 
566                                 //if there is a possibility that this sequence should be reversed
567                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
568                                         
569                                         string wasBetter = "";
570                                         //if the user wants you to try the reverse
571                                         if (flip) {
572                                                 //get reverse compliment
573                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
574                                                 copy->reverseComplement();
575                                                 
576                                                 //rerun alignment
577                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
578                                                 Sequence* templateSeq2 = &temp2;
579                                                 
580                                                 searchScore = templateDB->getSearchScore();
581                                                 
582                                                 nast2 = new Nast(alignment, copy, templateSeq2);
583                         
584                                                 //check if any better
585                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
586                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
587                                                         templateSeq = templateSeq2; 
588                                                         delete nast;
589                                                         nast = nast2;
590                                                         needToDeleteCopy = true;
591                                                 }else{  
592                                                         wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
593                                                         delete nast2;
594                                                         delete copy;    
595                                                 }
596                                         }
597                                         
598                                         //create accnos file with names
599                                         accnosFile << candidateSeq->getName() << wasBetter << endl;
600                                 }
601                                 
602                                 report.setCandidate(candidateSeq);
603                                 report.setTemplate(templateSeq);
604                                 report.setSearchParameters(search, searchScore);
605                                 report.setAlignmentParameters(align, alignment);
606                                 report.setNastParameters(*nast);
607         
608                                 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
609                                 
610                                 report.print();
611                                 delete nast;
612                                 if (needToDeleteCopy) {   delete copy;   }
613                         }
614                         delete candidateSeq;
615                         
616                         //report progress
617                         if((i+1) % 100 == 0){   m->mothurOut(toString(i+1)); m->mothurOutEndLine();             }
618                 }
619                 //report progress
620                 if((line->numSeqs) % 100 != 0){ m->mothurOut(toString(line->numSeqs)); m->mothurOutEndLine();           }
621                 
622                 alignmentFile.close();
623                 inFASTA.close();
624                 accnosFile.close();
625                 
626                 return 1;
627         }
628         catch(exception& e) {
629                 m->errorOut(e, "AlignCommand", "driver");
630                 exit(1);
631         }
632 }
633 //**********************************************************************************************************************
634 #ifdef USE_MPI
635 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<long>& MPIPos){
636         try {
637                 string outputString = "";
638                 MPI_Status statusReport; 
639                 MPI_Status statusAlign; 
640                 MPI_Status statusAccnos; 
641                 MPI_Status status; 
642                 int pid;
643                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
644         
645                 NastReport report;
646                 
647                 if (pid == 0) {
648                         outputString = report.getHeaders();
649                         int length = outputString.length();
650                         char buf[length];
651                         strcpy(buf, outputString.c_str()); 
652                 
653                         MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
654                 }
655                 
656                 for(int i=0;i<num;i++){
657                 
658                         if (m->control_pressed) {  return 0; }
659
660                         //read next sequence
661                         int length = MPIPos[start+i+1] - MPIPos[start+i];
662                         char buf4[length];
663                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
664                         
665                         string tempBuf = buf4;
666                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
667                         istringstream iss (tempBuf,istringstream::in);
668         
669                         Sequence* candidateSeq = new Sequence(iss);  
670                         int origNumBases = candidateSeq->getNumBases();
671                         string originalUnaligned = candidateSeq->getUnaligned();
672                         int numBasesNeeded = origNumBases * threshold;
673         
674                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
675                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
676                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
677                                 }
678                                                                 
679                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
680                                 Sequence* templateSeq = &temp;
681                                 
682                                 float searchScore = templateDB->getSearchScore();
683                                                                 
684                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
685                                 Sequence* copy;
686                                 
687                                 Nast* nast2;
688                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
689                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
690                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
691                                                                                                 //so this bool tells you if you need to delete it
692                                                                                                 
693                                 //if there is a possibility that this sequence should be reversed
694                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
695                                         
696                                         string wasBetter = "";
697                                         //if the user wants you to try the reverse
698                                         if (flip) {
699                                                 //get reverse compliment
700                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
701                                                 copy->reverseComplement();
702                                                 
703                                                 //rerun alignment
704                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
705                                                 Sequence* templateSeq2 = &temp2;
706                                                 
707                                                 searchScore = templateDB->getSearchScore();
708                                                 
709                                                 nast2 = new Nast(alignment, copy, templateSeq2);
710                         
711                                                 //check if any better
712                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
713                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
714                                                         templateSeq = templateSeq2; 
715                                                         delete nast;
716                                                         nast = nast2;
717                                                         needToDeleteCopy = true;
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                                         char buf[length];
731                                         strcpy(buf, outputString.c_str()); 
732                                 
733                                         MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
734                                         MPIWroteAccnos = true;
735                                 }
736                                 
737                                 report.setCandidate(candidateSeq);
738                                 report.setTemplate(templateSeq);
739                                 report.setSearchParameters(search, searchScore);
740                                 report.setAlignmentParameters(align, alignment);
741                                 report.setNastParameters(*nast);
742         
743                                 outputString =  ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
744                                 
745                                 //send results to parent
746                                 int length = outputString.length();
747                                 char buf2[length];
748                                 strcpy(buf2, outputString.c_str()); 
749                                 
750                                 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
751                                 
752                                 outputString = report.getReport();
753                                 
754                                 //send results to parent
755                                 length = outputString.length();
756                                 char buf3[length];
757                                 strcpy(buf3, outputString.c_str()); 
758                                 
759                                 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
760
761                                 delete nast;
762                                 if (needToDeleteCopy) {   delete copy;   }
763                         }
764                         delete candidateSeq;
765                         
766                         //report progress
767                         if((i+1) % 100 == 0){   cout << (toString(i+1)) << endl;                }
768                 }
769                 //report progress
770                 if((num) % 100 != 0){   cout << (toString(num)) << endl;                }
771                 
772                 return 1;
773         }
774         catch(exception& e) {
775                 m->errorOut(e, "AlignCommand", "driverMPI");
776                 exit(1);
777         }
778 }
779 #endif
780 /**************************************************************************************************/
781
782 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
783         try {
784 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
785                 int process = 0;
786                 int exitCommand = 1;
787                 //              processIDS.resize(0);
788                 
789                 //loop through and create all the processes you want
790                 while (process != processors) {
791                         int pid = fork();
792                         
793                         if (pid > 0) {
794                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
795                                 process++;
796                         }else if (pid == 0){
797                                 exitCommand = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
798                                 exit(0);
799                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
800                 }
801                 
802                 //force parent to wait until all the processes are done
803                 for (int i=0;i<processors;i++) { 
804                         int temp = processIDS[i];
805                         wait(&temp);
806                 }
807                 
808                 return exitCommand;
809 #endif          
810         }
811         catch(exception& e) {
812                 m->errorOut(e, "AlignCommand", "createProcesses");
813                 exit(1);
814         }
815 }
816
817 /**************************************************************************************************/
818
819 void AlignCommand::appendAlignFiles(string temp, string filename) {
820         try{
821                 
822                 ofstream output;
823                 ifstream input;
824                 openOutputFileAppend(filename, output);
825                 openInputFile(temp, input);
826                 
827                 while(char c = input.get()){
828                         if(input.eof())         {       break;                  }
829                         else                            {       output << c;    }
830                 }
831                 
832                 input.close();
833                 output.close();
834         }
835         catch(exception& e) {
836                 m->errorOut(e, "AlignCommand", "appendAlignFiles");
837                 exit(1);
838         }
839 }
840 //**********************************************************************************************************************
841
842 void AlignCommand::appendReportFiles(string temp, string filename) {
843         try{
844                 
845                 ofstream output;
846                 ifstream input;
847                 openOutputFileAppend(filename, output);
848                 openInputFile(temp, input);
849
850                 while (!input.eof())    {       char c = input.get(); if (c == 10 || c == 13){  break;  }       } // get header line
851                                 
852                 while(char c = input.get()){
853                         if(input.eof())         {       break;                  }
854                         else                            {       output << c;    }
855                 }
856                 
857                 input.close();
858                 output.close();
859         }
860         catch(exception& e) {
861                 m->errorOut(e, "AlignCommand", "appendReportFiles");
862                 exit(1);
863         }
864 }
865 //**********************************************************************************************************************