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