]> git.donarmstrong.com Git - mothur.git/blob - aligncommand.cpp
numSeqs passing fix
[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 = hasPath(it->second);
72                                         //if the user has not given a path then, add inputdir. else leave path alone.
73                                         if (path == "") {       parameters["template"] = inputDir + it->second;         }
74                                 }
75                         }
76
77                         //check for required parameters
78                         templateFileName = validParameter.validFile(parameters, "template", true);
79                         
80                         if (templateFileName == "not found") { 
81                                 m->mothurOut("template is a required parameter for the align.seqs command."); 
82                                 m->mothurOutEndLine();
83                                 abort = true; 
84                         }else if (templateFileName == "not open") { abort = true; }     
85                         
86                         candidateFileName = validParameter.validFile(parameters, "candidate", false);
87                         if (candidateFileName == "not found") { m->mothurOut("candidate is a required parameter for the align.seqs command."); m->mothurOutEndLine(); abort = true;  }
88                         else { 
89                                 splitAtDash(candidateFileName, candidateFileNames);
90                                 
91                                 //go through files and make sure they are good, if not, then disregard them
92                                 for (int i = 0; i < candidateFileNames.size(); i++) {
93                                         if (inputDir != "") {
94                                                 string path = hasPath(candidateFileNames[i]);
95                                                 //if the user has not given a path then, add inputdir. else leave path alone.
96                                                 if (path == "") {       candidateFileNames[i] = inputDir + candidateFileNames[i];               }
97                                         }
98         
99                                         int ableToOpen;
100                                         ifstream in;
101
102                                         ableToOpen = 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() + getSimpleName(candidateFileNames[i]);
108                                                         m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
109                                                         ableToOpen = 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 = 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 += hasPath(candidateFileNames[s]); }
234                         string alignFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align";
235                         string reportFileName = outputDir + getRootName(getSimpleName(candidateFileNames[s])) + "align.report";
236                         string accnosFileName = outputDir + getRootName(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 = 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                 vector<unsigned long int> positions = divideFile(candidateFileNames[s], processors);
350                                 
351                 for (int i = 0; i < (positions.size()-1); i++) {
352                         lines.push_back(new linePair(positions[i], positions[(i+1)]));
353                 }       
354         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
355                         if(processors == 1){
356                                 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
357                                 
358                                 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
359                                 
360                                 //delete accnos file if its blank else report to user
361                                 if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
362                                 else { 
363                                         m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
364                                         if (!flip) {
365                                                 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
366                                         }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
367                                         m->mothurOutEndLine();
368                                 }
369                         }else{
370                                 processIDS.resize(0);
371                                 
372                                 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]); 
373                                 
374                                 rename((alignFileName + toString(processIDS[0]) + ".temp").c_str(), alignFileName.c_str());
375                                 rename((reportFileName + toString(processIDS[0]) + ".temp").c_str(), reportFileName.c_str());
376                                 
377                                 //append alignment and report files
378                                 for(int i=1;i<processors;i++){
379                                         appendAlignFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
380                                         remove((alignFileName + toString(processIDS[i]) + ".temp").c_str());
381                                         
382                                         appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
383                                         remove((reportFileName + toString(processIDS[i]) + ".temp").c_str());
384                                 }
385                                 
386                                 vector<string> nonBlankAccnosFiles;
387                                 //delete blank accnos files generated with multiple processes
388                                 for(int i=0;i<processors;i++){  
389                                         if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
390                                                 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
391                                         }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());  }
392                                 }
393                                 
394                                 //append accnos files
395                                 if (nonBlankAccnosFiles.size() != 0) { 
396                                         rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
397                                         
398                                         for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
399                                                 appendAlignFiles(nonBlankAccnosFiles[h], accnosFileName);
400                                                 remove(nonBlankAccnosFiles[h].c_str());
401                                         }
402                                         m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
403                                         if (!flip) {
404                                                 m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
405                                         }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
406                                         m->mothurOutEndLine();
407                                 }else{ hasAccnos = false;  }
408                                 
409                                 if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
410                         }
411         #else
412                         numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
413                         
414                         if (m->control_pressed) { remove(accnosFileName.c_str()); remove(alignFileName.c_str()); remove(reportFileName.c_str()); return 0; }
415                         
416                         //delete accnos file if its blank else report to user
417                         if (isBlank(accnosFileName)) {  remove(accnosFileName.c_str());  hasAccnos = false; }
418                         else { 
419                                 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
420                                 if (!flip) {
421                                         m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
422                                 }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
423                                 m->mothurOutEndLine();
424                         }
425                         
426         #endif
427
428 #endif          
429
430
431                 #ifdef USE_MPI
432                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
433                                         
434                         if (pid == 0) { //only one process should output to screen
435                 #endif
436
437                         outputNames.push_back(alignFileName);
438                         outputNames.push_back(reportFileName);
439                         if (hasAccnos)  {       outputNames.push_back(accnosFileName);          }
440                         
441                 #ifdef USE_MPI
442                         }
443                 #endif
444
445                         m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
446                         m->mothurOutEndLine();
447                         m->mothurOutEndLine();
448                 }
449                 
450                 
451                 m->mothurOutEndLine();
452                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
453                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }
454                 m->mothurOutEndLine();
455
456                 return 0;
457         }
458         catch(exception& e) {
459                 m->errorOut(e, "AlignCommand", "execute");
460                 exit(1);
461         }
462 }
463
464 //**********************************************************************************************************************
465
466 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
467         try {
468                 ofstream alignmentFile;
469                 openOutputFile(alignFName, alignmentFile);
470                 
471                 ofstream accnosFile;
472                 openOutputFile(accnosFName, accnosFile);
473                 
474                 NastReport report(reportFName);
475                 
476                 ifstream inFASTA;
477                 openInputFile(filename, inFASTA);
478
479                 inFASTA.seekg(filePos->start);
480
481                 bool done = false;
482                 int count = 0;
483         
484                 while (!done) {
485                         
486                         if (m->control_pressed) {  return 0; }
487                         
488                         Sequence* candidateSeq = new Sequence(inFASTA);  gobble(inFASTA);
489         
490                         int origNumBases = candidateSeq->getNumBases();
491                         string originalUnaligned = candidateSeq->getUnaligned();
492                         int numBasesNeeded = origNumBases * threshold;
493         
494                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
495                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
496                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
497                                 }
498                                                                 
499                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
500                                 Sequence* templateSeq = &temp;
501                                 
502                                 float searchScore = templateDB->getSearchScore();
503                                                                 
504                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
505                                 Sequence* copy;
506                                 
507                                 Nast* nast2;
508                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
509                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
510                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
511                                                                                                 //so this bool tells you if you need to delete it
512                                                                                                 
513                                 //if there is a possibility that this sequence should be reversed
514                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
515                                         
516                                         string wasBetter =  "";
517                                         //if the user wants you to try the reverse
518                                         if (flip) {
519                                                 //get reverse compliment
520                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
521                                                 copy->reverseComplement();
522                                                 
523                                                 //rerun alignment
524                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
525                                                 Sequence* templateSeq2 = &temp2;
526                                                 
527                                                 searchScore = templateDB->getSearchScore();
528                                                 
529                                                 nast2 = new Nast(alignment, copy, templateSeq2);
530                         
531                                                 //check if any better
532                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
533                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
534                                                         templateSeq = templateSeq2; 
535                                                         delete nast;
536                                                         nast = nast2;
537                                                         needToDeleteCopy = true;
538                                                         wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
539                                                 }else{  
540                                                         wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
541                                                         delete nast2;
542                                                         delete copy;    
543                                                 }
544                                         }
545                                         
546                                         //create accnos file with names
547                                         accnosFile << candidateSeq->getName() << wasBetter << endl;
548                                 }
549                                 
550                                 report.setCandidate(candidateSeq);
551                                 report.setTemplate(templateSeq);
552                                 report.setSearchParameters(search, searchScore);
553                                 report.setAlignmentParameters(align, alignment);
554                                 report.setNastParameters(*nast);
555         
556                                 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
557                                 
558                                 report.print();
559                                 delete nast;
560                                 if (needToDeleteCopy) {   delete copy;   }
561                                 
562                                 count++;
563                         }
564                         delete candidateSeq;
565                         
566                         unsigned long int pos = inFASTA.tellg();
567                         if ((pos == -1) || (pos >= filePos->end)) { break; }
568                         
569                         //report progress
570                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
571                         
572                 }
573                 //report progress
574                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
575                 
576                 alignmentFile.close();
577                 inFASTA.close();
578                 accnosFile.close();
579                 
580                 return count;
581         }
582         catch(exception& e) {
583                 m->errorOut(e, "AlignCommand", "driver");
584                 exit(1);
585         }
586 }
587 //**********************************************************************************************************************
588 #ifdef USE_MPI
589 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long int>& MPIPos){
590         try {
591                 string outputString = "";
592                 MPI_Status statusReport; 
593                 MPI_Status statusAlign; 
594                 MPI_Status statusAccnos; 
595                 MPI_Status status; 
596                 int pid;
597                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
598         
599                 NastReport report;
600                 
601                 if (pid == 0) {
602                         outputString = report.getHeaders();
603                         int length = outputString.length();
604             
605                         char* buf = new char[length];
606                         memcpy(buf, outputString.c_str(), length);
607                 
608                         MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
609
610             delete buf;
611                 }
612                 
613                 for(int i=0;i<num;i++){
614                 
615                         if (m->control_pressed) {  return 0; }
616
617                         //read next sequence
618                         int length = MPIPos[start+i+1] - MPIPos[start+i];
619
620                         char* buf4 = new char[length];
621                         //memcpy(buf4, outputString.c_str(), length);
622
623                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
624                         
625                         string tempBuf = buf4;
626
627                         delete buf4;
628
629                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
630         
631                         istringstream iss (tempBuf,istringstream::in);
632
633                         Sequence* candidateSeq = new Sequence(iss);  
634         
635                         int origNumBases = candidateSeq->getNumBases();
636                         string originalUnaligned = candidateSeq->getUnaligned();
637                         int numBasesNeeded = origNumBases * threshold;
638         
639                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
640                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
641                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
642                                 }
643                                                                 
644                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
645                                 Sequence* templateSeq = &temp;
646                                 
647                                 float searchScore = templateDB->getSearchScore();
648                                                                 
649                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
650                                 Sequence* copy;
651                                 
652                                 Nast* nast2;
653                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
654                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
655                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
656                                                                                                 //so this bool tells you if you need to delete it
657                                                                                                 
658                                 //if there is a possibility that this sequence should be reversed
659                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
660                                         
661                                         string wasBetter = "";
662                                         //if the user wants you to try the reverse
663                                         if (flip) {
664                                                 //get reverse compliment
665                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
666                                                 copy->reverseComplement();
667                                                 
668                                                 //rerun alignment
669                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
670                                                 Sequence* templateSeq2 = &temp2;
671                                                 
672                                                 searchScore = templateDB->getSearchScore();
673                                                 
674                                                 nast2 = new Nast(alignment, copy, templateSeq2);
675                         
676                                                 //check if any better
677                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
678                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
679                                                         templateSeq = templateSeq2; 
680                                                         delete nast;
681                                                         nast = nast2;
682                                                         needToDeleteCopy = true;
683                                                 }else{  
684                                                         wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
685                                                         delete nast2;
686                                                         delete copy;    
687                                                 }
688                                         }
689                                         
690                                         //create accnos file with names
691                                         outputString = candidateSeq->getName() + wasBetter + "\n";
692                                         
693                                         //send results to parent
694                                         int length = outputString.length();
695
696                                         char* buf = new char[length];
697                                         memcpy(buf, outputString.c_str(), length);
698                                 
699                                         MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
700                                         delete buf;
701                                         MPIWroteAccnos = true;
702                                 }
703                                 
704                                 report.setCandidate(candidateSeq);
705                                 report.setTemplate(templateSeq);
706                                 report.setSearchParameters(search, searchScore);
707                                 report.setAlignmentParameters(align, alignment);
708                                 report.setNastParameters(*nast);
709         
710                                 outputString =  ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
711                                 
712                                 //send results to parent
713                                 int length = outputString.length();
714                                 char* buf2 = new char[length];
715                                 memcpy(buf2, outputString.c_str(), length);
716                                 
717                                 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
718                                 
719                                 delete buf2;
720
721                                 outputString = report.getReport();
722                                 
723                                 //send results to parent
724                                 length = outputString.length();
725                                 char* buf3 = new char[length];
726                                 memcpy(buf3, outputString.c_str(), length);
727                                 
728                                 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
729                                 
730                                 delete buf3;
731                                 delete nast;
732                                 if (needToDeleteCopy) {   delete copy;   }
733                         }
734                         delete candidateSeq;
735                         
736                         //report progress
737                         if((i+1) % 100 == 0){   cout << (toString(i+1)) << endl;                }
738                 }
739                 //report progress
740                 if((num) % 100 != 0){   cout << (toString(num)) << endl;                }
741                 
742                 return 1;
743         }
744         catch(exception& e) {
745                 m->errorOut(e, "AlignCommand", "driverMPI");
746                 exit(1);
747         }
748 }
749 #endif
750 /**************************************************************************************************/
751
752 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
753         try {
754 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
755                 int process = 0;
756                 int num = 0;
757                 //              processIDS.resize(0);
758                 
759                 //loop through and create all the processes you want
760                 while (process != processors) {
761                         int pid = fork();
762                         
763                         if (pid > 0) {
764                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
765                                 process++;
766                         }else if (pid == 0){
767                                 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
768                                 
769                                 //pass numSeqs to parent
770                                 ofstream out;
771                                 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
772                                 openOutputFile(tempFile, out);
773                                 out << num << endl;
774                                 out.close();
775                                 
776                                 exit(0);
777                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
778                 }
779                 
780                 //force parent to wait until all the processes are done
781                 for (int i=0;i<processors;i++) { 
782                         int temp = processIDS[i];
783                         wait(&temp);
784                 }
785                 
786                 for (int i = 0; i < processIDS.size(); i++) {
787                         ifstream in;
788                         string tempFile =  alignFileName + toString(processIDS[i]) + ".num.temp";
789                         openInputFile(tempFile, in);
790                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
791                         in.close(); remove(tempFile.c_str());
792                 }
793                 
794                 return num;
795 #endif          
796         }
797         catch(exception& e) {
798                 m->errorOut(e, "AlignCommand", "createProcesses");
799                 exit(1);
800         }
801 }
802 /**************************************************************************************************/
803
804 void AlignCommand::appendAlignFiles(string temp, string filename) {
805         try{
806                 
807                 ofstream output;
808                 ifstream input;
809                 openOutputFileAppend(filename, output);
810                 openInputFile(temp, input);
811                 
812                 while(char c = input.get()){
813                         if(input.eof())         {       break;                  }
814                         else                            {       output << c;    }
815                 }
816                 
817                 input.close();
818                 output.close();
819         }
820         catch(exception& e) {
821                 m->errorOut(e, "AlignCommand", "appendAlignFiles");
822                 exit(1);
823         }
824 }
825 //**********************************************************************************************************************
826
827 void AlignCommand::appendReportFiles(string temp, string filename) {
828         try{
829                 
830                 ofstream output;
831                 ifstream input;
832                 openOutputFileAppend(filename, output);
833                 openInputFile(temp, input);
834
835                 while (!input.eof())    {       char c = input.get(); if (c == 10 || c == 13){  break;  }       } // get header line
836                                 
837                 while(char c = input.get()){
838                         if(input.eof())         {       break;                  }
839                         else                            {       output << c;    }
840                 }
841                 
842                 input.close();
843                 output.close();
844         }
845         catch(exception& e) {
846                 m->errorOut(e, "AlignCommand", "appendReportFiles");
847                 exit(1);
848         }
849 }
850 //**********************************************************************************************************************