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