]> git.donarmstrong.com Git - mothur.git/blob - aligncommand.cpp
Merge remote-tracking branch 'mothur/master'
[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 "referencedb.h"
19
20 //**********************************************************************************************************************
21 vector<string> AlignCommand::setParameters(){   
22         try {
23                 CommandParameter ptemplate("reference", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(ptemplate);
24                 CommandParameter pcandidate("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pcandidate);
25                 CommandParameter psearch("search", "Multiple", "kmer-blast-suffix", "kmer", "", "", "",false,false); parameters.push_back(psearch);
26                 CommandParameter pksize("ksize", "Number", "", "8", "", "", "",false,false); parameters.push_back(pksize);
27                 CommandParameter pmatch("match", "Number", "", "1.0", "", "", "",false,false); parameters.push_back(pmatch);
28                 CommandParameter palign("align", "Multiple", "needleman-gotoh-blast-noalign", "needleman", "", "", "",false,false); parameters.push_back(palign);
29                 CommandParameter pmismatch("mismatch", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pmismatch);
30                 CommandParameter pgapopen("gapopen", "Number", "", "-2.0", "", "", "",false,false); parameters.push_back(pgapopen);
31                 CommandParameter pgapextend("gapextend", "Number", "", "-1.0", "", "", "",false,false); parameters.push_back(pgapextend);
32                 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
33                 CommandParameter pflip("flip", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(pflip);
34                 CommandParameter psave("save", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(psave);
35                 CommandParameter pthreshold("threshold", "Number", "", "0.50", "", "", "",false,false); parameters.push_back(pthreshold);
36                 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
37                 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
38
39                 vector<string> myArray;
40                 for (int i = 0; i < parameters.size(); i++) {   myArray.push_back(parameters[i].name);          }
41                 return myArray;
42         }
43         catch(exception& e) {
44                 m->errorOut(e, "AlignCommand", "setParameters");
45                 exit(1);
46         }
47 }
48 //**********************************************************************************************************************
49 string AlignCommand::getHelpString(){   
50         try {
51                 string helpString = "";
52                 helpString += "The align.seqs command reads a file containing sequences and creates an alignment file and a report file.";
53                 helpString += "The align.seqs command parameters are reference, fasta, search, ksize, align, match, mismatch, gapopen, gapextend and processors.";
54                 helpString += "The reference and fasta parameters are required. You may leave fasta blank if you have a valid fasta file. You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amzon.fasta.";
55                 helpString += "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.";
56                 helpString += "The align parameter allows you to specify the alignment method to use.  Your options are: gotoh, needleman, blast and noalign. The default is needleman.";
57                 helpString += "The ksize parameter allows you to specify the kmer size for finding most similar template to candidate.  The default is 8.";
58                 helpString += "The match parameter allows you to specify the bonus for having the same base. The default is 1.0.";
59                 helpString += "The mistmatch parameter allows you to specify the penalty for having different bases.  The default is -1.0.";
60                 helpString += "The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.";
61                 helpString += "The gapextend parameter allows you to specify the penalty for extending a gap in an alignment.  The default is -1.0.";
62                 helpString += "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.";
63                 helpString += "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.";
64                 helpString += "If the flip parameter is set to true the reverse complement of the sequence is aligned and the better alignment is reported.";
65                 helpString += "If the save parameter is set to true the reference sequences will be saved in memory, to clear them later you can use the clear.memory command. Default=f.";
66                 helpString += "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.";
67                 helpString += "The align.seqs command should be in the following format:";
68                 helpString += "align.seqs(reference=yourTemplateFile, fasta=yourCandidateFile, align=yourAlignmentMethod, search=yourSearchmethod, ksize=yourKmerSize, match=yourMatchBonus, mismatch=yourMismatchpenalty, gapopen=yourGapopenPenalty, gapextend=yourGapExtendPenalty)";
69                 helpString += "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)";
70                 helpString += "Note: No spaces between parameter labels (i.e. candidate), '=' and parameters (i.e.yourFastaFile).";
71                 return helpString;
72         }
73         catch(exception& e) {
74                 m->errorOut(e, "AlignCommand", "getHelpString");
75                 exit(1);
76         }
77 }
78 //**********************************************************************************************************************
79 string AlignCommand::getOutputFileNameTag(string type, string inputName=""){    
80         try {
81         string tag = "";
82                 map<string, vector<string> >::iterator it;
83         
84         //is this a type this command creates
85         it = outputTypes.find(type);
86         if (it == outputTypes.end()) {  m->mothurOut("[ERROR]: this command doesn't create a " + type + " output file.\n"); }
87         else {
88             if (type == "fasta") {  tag = "align"; }
89             else if (type == "alignreport") {  tag = "align.report"; }
90             else if (type == "accnos") {  tag = "flip.accnos"; }
91             else { m->mothurOut("[ERROR]: No definition for type " + type + " output file tag.\n"); m->control_pressed = true;  }
92         }
93         return tag;
94         }
95         catch(exception& e) {
96                 m->errorOut(e, "AlignCommand", "getOutputFileName");
97                 exit(1);
98         }
99 }
100
101 //**********************************************************************************************************************
102 AlignCommand::AlignCommand(){   
103         try {
104                 abort = true; calledHelp = true; 
105                 setParameters();
106                 vector<string> tempOutNames;
107                 outputTypes["fasta"] = tempOutNames;
108                 outputTypes["alignreport"] = tempOutNames;
109                 outputTypes["accnos"] = tempOutNames;
110         }
111         catch(exception& e) {
112                 m->errorOut(e, "AlignCommand", "AlignCommand");
113                 exit(1);
114         }
115 }
116 //**********************************************************************************************************************
117 AlignCommand::AlignCommand(string option)  {
118         try {
119                 abort = false; calledHelp = false;  
120                 ReferenceDB* rdb = ReferenceDB::getInstance();
121                 
122                 //allow user to run help
123                 if(option == "help") { help(); abort = true; calledHelp = true;}
124                 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
125                 
126                 else {
127                         vector<string> myArray = setParameters();
128                         
129                         OptionParser parser(option);
130                         map<string, string> parameters = parser.getParameters(); 
131                         
132                         ValidParameters validParameter("align.seqs");
133                         map<string, string>::iterator it;
134                         
135                         //check to make sure all parameters are valid for command
136                         for (it = parameters.begin(); it != parameters.end(); it++) { 
137                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
138                         }
139                         
140                         //initialize outputTypes
141                         vector<string> tempOutNames;
142                         outputTypes["fasta"] = tempOutNames;
143                         outputTypes["alignreport"] = tempOutNames;
144                         outputTypes["accnos"] = tempOutNames;
145                         
146                         //if the user changes the output directory command factory will send this info to us in the output parameter 
147                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  outputDir = "";         }
148                         
149
150                         //if the user changes the input directory command factory will send this info to us in the output parameter 
151                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
152                         
153                         if (inputDir == "not found"){   inputDir = "";          }
154                         else {
155                                 string path;
156
157                                 it = parameters.find("reference");
158
159                                 //user has given a template file
160                                 if(it != parameters.end()){ 
161                                         path = m->hasPath(it->second);
162                                         //if the user has not given a path then, add inputdir. else leave path alone.
163                                         if (path == "") {       parameters["reference"] = inputDir + it->second;                }
164                                 }
165                         }
166
167                         candidateFileName = validParameter.validFile(parameters, "fasta", false);
168                         if (candidateFileName == "not found") { 
169                                 //if there is a current fasta file, use it
170                                 string filename = m->getFastaFile(); 
171                                 if (filename != "") { candidateFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
172                                 else {  m->mothurOut("You have no current fastafile and the candidate parameter is required."); m->mothurOutEndLine(); abort = true; }
173                         }else { 
174                                 m->splitAtDash(candidateFileName, candidateFileNames);
175                                 
176                                 //go through files and make sure they are good, if not, then disregard them
177                                 for (int i = 0; i < candidateFileNames.size(); i++) {
178                                         //candidateFileNames[i] = m->getFullPathName(candidateFileNames[i]);
179                                         
180                                         bool ignore = false;
181                                         if (candidateFileNames[i] == "current") { 
182                                                 candidateFileNames[i] = m->getFastaFile(); 
183                                                 if (candidateFileNames[i] != "") {  m->mothurOut("Using " + candidateFileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
184                                                 else {  
185                                                         m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true; 
186                                                         //erase from file list
187                                                         candidateFileNames.erase(candidateFileNames.begin()+i);
188                                                         i--;
189                                                 }
190                                         }
191                                         
192                                         if (!ignore) {
193                                         
194                                                 if (inputDir != "") {
195                                                         string path = m->hasPath(candidateFileNames[i]);
196                                                         //if the user has not given a path then, add inputdir. else leave path alone.
197                                                         if (path == "") {       candidateFileNames[i] = inputDir + candidateFileNames[i];               }
198                                                 }
199                 
200                                                 int ableToOpen;
201                                                 ifstream in;
202                                                 ableToOpen = m->openInputFile(candidateFileNames[i], in, "noerror");
203                                                 in.close();     
204                                                 
205                                                 //if you can't open it, try default location
206                                                 if (ableToOpen == 1) {
207                                                         if (m->getDefaultPath() != "") { //default path is set
208                                                                 string tryPath = m->getDefaultPath() + m->getSimpleName(candidateFileNames[i]);
209                                                                 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
210                                                                 ifstream in2;
211                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
212                                                                 in2.close();
213                                                                 candidateFileNames[i] = tryPath;
214                                                         }
215                                                 }
216                                                 
217                                                 //if you can't open it, try output location
218                                                 if (ableToOpen == 1) {
219                                                         if (m->getOutputDir() != "") { //default path is set
220                                                                 string tryPath = m->getOutputDir() + m->getSimpleName(candidateFileNames[i]);
221                                                                 m->mothurOut("Unable to open " + candidateFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
222                                                                 ifstream in2;
223                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
224                                                                 in2.close();
225                                                                 candidateFileNames[i] = tryPath;
226                                                         }
227                                                 }
228                                                 
229                                                                                 
230
231                                                 if (ableToOpen == 1) { 
232                                                         m->mothurOut("Unable to open " + candidateFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
233                                                         //erase from file list
234                                                         candidateFileNames.erase(candidateFileNames.begin()+i);
235                                                         i--;
236                                                 }else {
237                                                         m->setFastaFile(candidateFileNames[i]);
238                                                 }
239                                         }
240                                 }
241                                 
242                                 //make sure there is at least one valid file left
243                                 if (candidateFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
244                         }
245                 
246                         //check for optional parameter and set defaults
247                         // ...at some point should added some additional type checking...
248                         string temp;
249                         temp = validParameter.validFile(parameters, "ksize", false);            if (temp == "not found"){       temp = "8";                             }
250                         m->mothurConvert(temp, kmerSize); 
251                         
252                         temp = validParameter.validFile(parameters, "match", false);            if (temp == "not found"){       temp = "1.0";                   }
253                         m->mothurConvert(temp, match);  
254                         
255                         temp = validParameter.validFile(parameters, "mismatch", false);         if (temp == "not found"){       temp = "-1.0";                  }
256                         m->mothurConvert(temp, misMatch);  
257                         
258                         temp = validParameter.validFile(parameters, "gapopen", false);          if (temp == "not found"){       temp = "-2.0";                  }
259                         m->mothurConvert(temp, gapOpen);  
260                         
261                         temp = validParameter.validFile(parameters, "gapextend", false);        if (temp == "not found"){       temp = "-1.0";                  }
262                         m->mothurConvert(temp, gapExtend); 
263                         
264                         temp = validParameter.validFile(parameters, "processors", false);       if (temp == "not found"){       temp = m->getProcessors();      }
265                         m->setProcessors(temp);
266                         m->mothurConvert(temp, processors); 
267                         
268                         temp = validParameter.validFile(parameters, "flip", false);                     if (temp == "not found"){       temp = "f";                             }
269                         flip = m->isTrue(temp); 
270                         
271                         temp = validParameter.validFile(parameters, "save", false);                     if (temp == "not found"){       temp = "f";                             }
272                         save = m->isTrue(temp); 
273                         rdb->save = save; 
274                         if (save) { //clear out old references
275                                 rdb->clearMemory();
276                         }
277                         
278                         //this has to go after save so that if the user sets save=t and provides no reference we abort
279                         templateFileName = validParameter.validFile(parameters, "reference", true);
280                         if (templateFileName == "not found") { 
281                                 //check for saved reference sequences
282                                 if (rdb->referenceSeqs.size() != 0) {
283                                         templateFileName = "saved";
284                                 }else {
285                                         m->mothurOut("[ERROR]: You don't have any saved reference sequences and the reference parameter is a required for the align.seqs command."); 
286                                         m->mothurOutEndLine();
287                                         abort = true; 
288                                 }
289                         }else if (templateFileName == "not open") { abort = true; }     
290                         else {  if (save) {     rdb->setSavedReference(templateFileName);       }       }
291                         
292                         temp = validParameter.validFile(parameters, "threshold", false);        if (temp == "not found"){       temp = "0.50";                  }
293                         m->mothurConvert(temp, threshold); 
294                         
295                         search = validParameter.validFile(parameters, "search", false);         if (search == "not found"){     search = "kmer";                }
296                         if ((search != "suffix") && (search != "kmer") && (search != "blast")) { m->mothurOut("invalid search option: choices are kmer, suffix or blast."); m->mothurOutEndLine(); abort=true; }
297                         
298                         align = validParameter.validFile(parameters, "align", false);           if (align == "not found"){      align = "needleman";    }
299                         if ((align != "needleman") && (align != "gotoh") && (align != "blast") && (align != "noalign")) { m->mothurOut("invalid align option: choices are needleman, gotoh, blast or noalign."); m->mothurOutEndLine(); abort=true; }
300
301                 }
302                 
303         }
304         catch(exception& e) {
305                 m->errorOut(e, "AlignCommand", "AlignCommand");
306                 exit(1);
307         }
308 }
309 //**********************************************************************************************************************
310 AlignCommand::~AlignCommand(){  
311
312         if (abort == false) {
313                 for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
314                 delete templateDB;
315         }
316 }
317 //**********************************************************************************************************************
318
319 int AlignCommand::execute(){
320         try {
321                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
322
323                 templateDB = new AlignmentDB(templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch, rand());
324                 
325                 for (int s = 0; s < candidateFileNames.size(); s++) {
326                         if (m->control_pressed) { outputTypes.clear(); return 0; }
327                         
328                         m->mothurOut("Aligning sequences from " + candidateFileNames[s] + " ..." ); m->mothurOutEndLine();
329                         
330                         if (outputDir == "") {  outputDir += m->hasPath(candidateFileNames[s]); }
331                         string alignFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + getOutputFileNameTag("fasta");  
332                         string reportFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + getOutputFileNameTag("alignreport");
333                         string accnosFileName = outputDir + m->getRootName(m->getSimpleName(candidateFileNames[s])) + getOutputFileNameTag("accnos");
334                         bool hasAccnos = true;
335                         
336                         int numFastaSeqs = 0;
337                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
338                         int start = time(NULL);
339                 
340 #ifdef USE_MPI  
341                                 int pid, numSeqsPerProcessor; 
342                                 int tag = 2001;
343                                 vector<unsigned long long> MPIPos;
344                                 MPIWroteAccnos = false;
345                         
346                                 MPI_Status status; 
347                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
348                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); 
349
350                                 MPI_File inMPI;
351                                 MPI_File outMPIAlign;
352                                 MPI_File outMPIReport;
353                                 MPI_File outMPIAccnos;
354                                 
355                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
356                                 int inMode=MPI_MODE_RDONLY; 
357                                 
358                                 char outAlignFilename[1024];
359                                 strcpy(outAlignFilename, alignFileName.c_str());
360                                 
361                                 char outReportFilename[1024];
362                                 strcpy(outReportFilename, reportFileName.c_str());
363                                 
364                                 char outAccnosFilename[1024];
365                                 strcpy(outAccnosFilename, accnosFileName.c_str());
366                                 
367                                 char inFileName[1024];
368                                 strcpy(inFileName, candidateFileNames[s].c_str());
369                                 
370                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
371                                 MPI_File_open(MPI_COMM_WORLD, outAlignFilename, outMode, MPI_INFO_NULL, &outMPIAlign);
372                                 MPI_File_open(MPI_COMM_WORLD, outReportFilename, outMode, MPI_INFO_NULL, &outMPIReport);
373                                 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
374                                 
375                                 if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIAlign);  MPI_File_close(&outMPIReport);  MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
376                                 
377                                 if (pid == 0) { //you are the root process 
378                                         
379                                         MPIPos = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
380                                         
381                                         //send file positions to all processes
382                                         for(int i = 1; i < processors; i++) { 
383                                                 MPI_Send(&numFastaSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
384                                                 MPI_Send(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
385                                         }
386                                         
387                                         //figure out how many sequences you have to align
388                                         numSeqsPerProcessor = numFastaSeqs / processors;
389                                         int startIndex =  pid * numSeqsPerProcessor;
390                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
391                                         
392                                         //align your part
393                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
394                                         
395                                         if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIAlign);  MPI_File_close(&outMPIReport);  MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
396
397                                         for (int i = 1; i < processors; i++) {
398                                                 bool tempResult;
399                                                 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
400                                                 if (tempResult != 0) { MPIWroteAccnos = true; }
401                                         }
402                                 }else{ //you are a child process
403                                         MPI_Recv(&numFastaSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
404                                         MPIPos.resize(numFastaSeqs+1);
405                                         MPI_Recv(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
406
407                                         
408                                         //figure out how many sequences you have to align
409                                         numSeqsPerProcessor = numFastaSeqs / processors;
410                                         int startIndex =  pid * numSeqsPerProcessor;
411                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor;         }
412                                         
413                                         
414                                         //align your part
415                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPIAlign, outMPIReport, outMPIAccnos, MPIPos);
416                                         
417                                         if (m->control_pressed) { MPI_File_close(&inMPI);  MPI_File_close(&outMPIAlign);  MPI_File_close(&outMPIReport);  MPI_File_close(&outMPIAccnos); outputTypes.clear(); return 0; }
418
419                                         MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD); 
420                                 }
421                                 
422                                 //close files 
423                                 MPI_File_close(&inMPI);
424                                 MPI_File_close(&outMPIAlign);
425                                 MPI_File_close(&outMPIReport);
426                                 MPI_File_close(&outMPIAccnos);
427                                 
428                                 //delete accnos file if blank
429                                 if (pid == 0) {
430                                         //delete accnos file if its blank else report to user
431                                         if (MPIWroteAccnos) { 
432                                                 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
433                                                 if (!flip) {
434                                                         m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
435                                                 }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
436                                                 m->mothurOutEndLine();
437                                         }else { 
438                                                 //MPI_Info info;
439                                                 //MPI_File_delete(outAccnosFilename, info);
440                                                 hasAccnos = false;      
441                                                 m->mothurRemove(accnosFileName); 
442                                         }
443                                 }
444                                 
445 #else
446
447                         vector<unsigned long long> positions; 
448                 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
449                         positions = m->divideFile(candidateFileNames[s], processors);
450                         for (int i = 0; i < (positions.size()-1); i++) {        lines.push_back(new linePair(positions[i], positions[(i+1)]));  }
451                 #else
452                         if (processors == 1) {
453                                 lines.push_back(new linePair(0, 1000));
454                         }else {
455                                 positions = m->setFilePosFasta(candidateFileNames[s], numFastaSeqs); 
456                                 if (positions.size() < processors) { processors = positions.size(); }
457                 
458                                 //figure out how many sequences you have to process
459                                 int numSeqsPerProcessor = numFastaSeqs / processors;
460                                 for (int i = 0; i < processors; i++) {
461                                         int startIndex =  i * numSeqsPerProcessor;
462                                         if(i == (processors - 1)){      numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;   }
463                                         lines.push_back(new linePair(positions[startIndex], numSeqsPerProcessor));
464                                 }
465                         }
466                 #endif
467                         
468                         if(processors == 1){
469                                 numFastaSeqs = driver(lines[0], alignFileName, reportFileName, accnosFileName, candidateFileNames[s]);
470                         }else{
471                                 numFastaSeqs = createProcesses(alignFileName, reportFileName, accnosFileName, candidateFileNames[s]); 
472                         }
473                                 
474                         if (m->control_pressed) { m->mothurRemove(accnosFileName); m->mothurRemove(alignFileName); m->mothurRemove(reportFileName); outputTypes.clear();  return 0; }
475                         
476                         //delete accnos file if its blank else report to user
477                         if (m->isBlank(accnosFileName)) {  m->mothurRemove(accnosFileName);  hasAccnos = false; }
478                         else { 
479                                 m->mothurOut("Some of you sequences generated alignments that eliminated too many bases, a list is provided in " + accnosFileName + ".");
480                                 if (!flip) {
481                                         m->mothurOut(" If you set the flip parameter to true mothur will try aligning the reverse compliment as well."); 
482                                 }else{  m->mothurOut(" If the reverse compliment proved to be better it was reported.");  }
483                                 m->mothurOutEndLine();
484                         }
485
486 #endif          
487
488
489                 #ifdef USE_MPI
490                         MPI_Comm_rank(MPI_COMM_WORLD, &pid); 
491                                         
492                         if (pid == 0) { //only one process should output to screen
493                 #endif
494
495                         outputNames.push_back(alignFileName); outputTypes["fasta"].push_back(alignFileName);
496                         outputNames.push_back(reportFileName); outputTypes["alignreport"].push_back(reportFileName);
497                         if (hasAccnos)  {       outputNames.push_back(accnosFileName);  outputTypes["accnos"].push_back(accnosFileName);  }
498                         
499                 #ifdef USE_MPI
500                         }
501                 #endif
502
503                         m->mothurOut("It took " + toString(time(NULL) - start) + " secs to align " + toString(numFastaSeqs) + " sequences.");
504                         m->mothurOutEndLine();
505                         m->mothurOutEndLine();
506                 }
507                 
508                 //set align file as new current fastafile
509                 string currentFasta = "";
510                 itTypes = outputTypes.find("fasta");
511                 if (itTypes != outputTypes.end()) {
512                         if ((itTypes->second).size() != 0) { currentFasta = (itTypes->second)[0]; m->setFastaFile(currentFasta); }
513                 }
514                 
515                 m->mothurOutEndLine();
516                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
517                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }
518                 m->mothurOutEndLine();
519
520                 return 0;
521         }
522         catch(exception& e) {
523                 m->errorOut(e, "AlignCommand", "execute");
524                 exit(1);
525         }
526 }
527
528 //**********************************************************************************************************************
529 int AlignCommand::driver(linePair* filePos, string alignFName, string reportFName, string accnosFName, string filename){
530         try {
531                 ofstream alignmentFile;
532                 m->openOutputFile(alignFName, alignmentFile);
533                 
534                 ofstream accnosFile;
535                 m->openOutputFile(accnosFName, accnosFile);
536                 
537                 NastReport report(reportFName);
538                 
539                 ifstream inFASTA;
540                 m->openInputFile(filename, inFASTA);
541
542                 inFASTA.seekg(filePos->start);
543
544                 bool done = false;
545                 int count = 0;
546                 
547                 //moved this into driver to avoid deep copies in windows paralellized version
548                 Alignment* alignment;
549                 int longestBase = templateDB->getLongestBase();
550                 if(align == "gotoh")                    {       alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase);                 }
551                 else if(align == "needleman")   {       alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);                                }
552                 else if(align == "blast")               {       alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch);            }
553                 else if(align == "noalign")             {       alignment = new NoAlign();                                                                                                      }
554                 else {
555                         m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
556                         m->mothurOutEndLine();
557                         alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
558                 }
559         
560                 while (!done) {
561                         
562                         if (m->control_pressed) {  break; }
563                         
564                         Sequence* candidateSeq = new Sequence(inFASTA);  m->gobble(inFASTA);
565                         report.setCandidate(candidateSeq);
566
567                         int origNumBases = candidateSeq->getNumBases();
568                         string originalUnaligned = candidateSeq->getUnaligned();
569                         int numBasesNeeded = origNumBases * threshold;
570         
571                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
572                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
573                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
574                                 }
575                                                                 
576                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
577                                 Sequence* templateSeq = &temp;
578                                 
579                                 float searchScore = templateDB->getSearchScore();
580                                                                 
581                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
582                 
583                                 Sequence* copy;
584                                 
585                                 Nast* nast2;
586                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
587                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
588                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
589                                                                                                 //so this bool tells you if you need to delete it
590                                                                                                 
591                                 //if there is a possibility that this sequence should be reversed
592                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
593                                         
594                                         string wasBetter =  "";
595                                         //if the user wants you to try the reverse
596                                         if (flip) {
597                                 
598                                                 //get reverse compliment
599                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
600                                                 copy->reverseComplement();
601                                                 
602                                                 //rerun alignment
603                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
604                                                 Sequence* templateSeq2 = &temp2;
605                                                 
606                                                 searchScore = templateDB->getSearchScore();
607                                                 
608                                                 nast2 = new Nast(alignment, copy, templateSeq2);
609                         
610                                                 //check if any better
611                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
612                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
613                                                         templateSeq = templateSeq2; 
614                                                         delete nast;
615                                                         nast = nast2;
616                                                         needToDeleteCopy = true;
617                                                         wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
618                                                 }else{  
619                                                         wasBetter = "\treverse complement did NOT produce a better alignment so it was not used, please check sequence.";
620                                                         delete nast2;
621                                                         delete copy;    
622                                                 }
623                                         }
624                                         
625                                         //create accnos file with names
626                                         accnosFile << candidateSeq->getName() << wasBetter << endl;
627                                 }
628                                 
629                                 report.setTemplate(templateSeq);
630                                 report.setSearchParameters(search, searchScore);
631                                 report.setAlignmentParameters(align, alignment);
632                                 report.setNastParameters(*nast);
633         
634                                 alignmentFile << '>' << candidateSeq->getName() << '\n' << candidateSeq->getAligned() << endl;
635                                 
636                                 report.print();
637                                 delete nast;
638                                 if (needToDeleteCopy) {   delete copy;   }
639                                 
640                                 count++;
641                         }
642                         delete candidateSeq;
643                         
644                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
645                                 unsigned long long pos = inFASTA.tellg();
646                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
647                         #else
648                                 if (inFASTA.eof()) { break; }
649                         #endif
650                         
651                         //report progress
652                         if((count) % 100 == 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
653                         
654                 }
655                 //report progress
656                 if((count) % 100 != 0){ m->mothurOut(toString(count)); m->mothurOutEndLine();           }
657                 
658                 delete alignment;
659                 alignmentFile.close();
660                 inFASTA.close();
661                 accnosFile.close();
662                 
663                 return count;
664         }
665         catch(exception& e) {
666                 m->errorOut(e, "AlignCommand", "driver");
667                 exit(1);
668         }
669 }
670 //**********************************************************************************************************************
671 #ifdef USE_MPI
672 int AlignCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& alignFile, MPI_File& reportFile, MPI_File& accnosFile, vector<unsigned long long>& MPIPos){
673         try {
674                 string outputString = "";
675                 MPI_Status statusReport; 
676                 MPI_Status statusAlign; 
677                 MPI_Status statusAccnos; 
678                 MPI_Status status; 
679                 int pid;
680                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
681         
682                 NastReport report;
683                 
684                 if (pid == 0) {
685                         outputString = report.getHeaders();
686                         int length = outputString.length();
687             
688                         char* buf = new char[length];
689                         memcpy(buf, outputString.c_str(), length);
690                 
691                         MPI_File_write_shared(reportFile, buf, length, MPI_CHAR, &statusReport);
692
693             delete buf;
694                 }
695                 
696                 Alignment* alignment;
697                 int longestBase = templateDB->getLongestBase();
698                 if(align == "gotoh")                    {       alignment = new GotohOverlap(gapOpen, gapExtend, match, misMatch, longestBase);                 }
699                 else if(align == "needleman")   {       alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);                                }
700                 else if(align == "blast")               {       alignment = new BlastAlignment(gapOpen, gapExtend, match, misMatch);            }
701                 else if(align == "noalign")             {       alignment = new NoAlign();                                                                                                      }
702                 else {
703                         m->mothurOut(align + " is not a valid alignment option. I will run the command using needleman.");
704                         m->mothurOutEndLine();
705                         alignment = new NeedlemanOverlap(gapOpen, match, misMatch, longestBase);
706                 }
707                 
708                 
709                 for(int i=0;i<num;i++){
710                 
711                         if (m->control_pressed) { delete alignment; return 0; }
712
713                         //read next sequence
714                         int length = MPIPos[start+i+1] - MPIPos[start+i];
715
716                         char* buf4 = new char[length];
717                         //memcpy(buf4, outputString.c_str(), length);
718
719                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
720                         
721                         string tempBuf = buf4;
722
723                         delete buf4;
724
725                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
726         
727                         istringstream iss (tempBuf,istringstream::in);
728
729                         Sequence* candidateSeq = new Sequence(iss);  
730                         report.setCandidate(candidateSeq);
731
732                         int origNumBases = candidateSeq->getNumBases();
733                         string originalUnaligned = candidateSeq->getUnaligned();
734                         int numBasesNeeded = origNumBases * threshold;
735         
736                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
737                                 if (candidateSeq->getUnaligned().length() > alignment->getnRows()) {
738                                         alignment->resize(candidateSeq->getUnaligned().length()+1);
739                                 }
740                                                                 
741                                 Sequence temp = templateDB->findClosestSequence(candidateSeq);
742                                 Sequence* templateSeq = &temp;
743                                 
744                                 float searchScore = templateDB->getSearchScore();
745                                                                 
746                                 Nast* nast = new Nast(alignment, candidateSeq, templateSeq);
747                                 Sequence* copy;
748                                 
749                                 Nast* nast2;
750                                 bool needToDeleteCopy = false;  //this is needed in case you have you enter the ifs below
751                                                                                                 //since nast does not make a copy of hte sequence passed, and it is used by the reporter below
752                                                                                                 //you can't delete the copy sequence til after you report, but you may choose not to create it in the first place
753                                                                                                 //so this bool tells you if you need to delete it
754                                                                                                 
755                                 //if there is a possibility that this sequence should be reversed
756                                 if (candidateSeq->getNumBases() < numBasesNeeded) {
757                                         
758                                         string wasBetter = "";
759                                         //if the user wants you to try the reverse
760                                         if (flip) {
761                                                 //get reverse compliment
762                                                 copy = new Sequence(candidateSeq->getName(), originalUnaligned);
763                                                 copy->reverseComplement();
764                                                 
765                                                 //rerun alignment
766                                                 Sequence temp2 = templateDB->findClosestSequence(copy);
767                                                 Sequence* templateSeq2 = &temp2;
768                                                 
769                                                 searchScore = templateDB->getSearchScore();
770                                                 
771                                                 nast2 = new Nast(alignment, copy, templateSeq2);
772                         
773                                                 //check if any better
774                                                 if (copy->getNumBases() > candidateSeq->getNumBases()) {
775                                                         candidateSeq->setAligned(copy->getAligned());  //use reverse compliments alignment since its better
776                                                         templateSeq = templateSeq2; 
777                                                         delete nast;
778                                                         nast = nast2;
779                                                         needToDeleteCopy = true;
780                                                         wasBetter = "\treverse complement produced a better alignment, so mothur used the reverse complement.";
781                                                 }else{  
782                                                         wasBetter = "\treverse complement did NOT produce a better alignment, please check sequence.";
783                                                         delete nast2;
784                                                         delete copy;    
785                                                 }
786                                         }
787                                         
788                                         //create accnos file with names
789                                         outputString = candidateSeq->getName() + wasBetter + "\n";
790                                         
791                                         //send results to parent
792                                         int length = outputString.length();
793
794                                         char* buf = new char[length];
795                                         memcpy(buf, outputString.c_str(), length);
796                                 
797                                         MPI_File_write_shared(accnosFile, buf, length, MPI_CHAR, &statusAccnos);
798                                         delete buf;
799                                         MPIWroteAccnos = true;
800                                 }
801                                 
802                                 report.setTemplate(templateSeq);
803                                 report.setSearchParameters(search, searchScore);
804                                 report.setAlignmentParameters(align, alignment);
805                                 report.setNastParameters(*nast);
806         
807                                 outputString =  ">" + candidateSeq->getName() + "\n" + candidateSeq->getAligned() + "\n";
808                                 
809                                 //send results to parent
810                                 int length = outputString.length();
811                                 char* buf2 = new char[length];
812                                 memcpy(buf2, outputString.c_str(), length);
813                                 
814                                 MPI_File_write_shared(alignFile, buf2, length, MPI_CHAR, &statusAlign);
815                                 
816                                 delete buf2;
817
818                                 outputString = report.getReport();
819                                 
820                                 //send results to parent
821                                 length = outputString.length();
822                                 char* buf3 = new char[length];
823                                 memcpy(buf3, outputString.c_str(), length);
824                                 
825                                 MPI_File_write_shared(reportFile, buf3, length, MPI_CHAR, &statusReport);
826                                 
827                                 delete buf3;
828                                 delete nast;
829                                 if (needToDeleteCopy) {   delete copy;   }
830                         }
831                         delete candidateSeq;
832                         
833                         //report progress
834                         if((i+1) % 100 == 0){   cout << (toString(i+1)) << endl;                }
835                 }
836                 //report progress
837                 if((num) % 100 != 0){   cout << (toString(num)) << endl;                }
838                 
839                 return 1;
840         }
841         catch(exception& e) {
842                 m->errorOut(e, "AlignCommand", "driverMPI");
843                 exit(1);
844         }
845 }
846 #endif
847 /**************************************************************************************************/
848
849 int AlignCommand::createProcesses(string alignFileName, string reportFileName, string accnosFName, string filename) {
850         try {
851                 int num = 0;
852                 processIDS.resize(0);
853 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux) || (__linux__) || (__unix__) || (__unix)
854                 int process = 1;
855                 
856                 //loop through and create all the processes you want
857                 while (process != processors) {
858                         int pid = fork();
859                         
860                         if (pid > 0) {
861                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
862                                 process++;
863                         }else if (pid == 0){
864                                 num = driver(lines[process], alignFileName + toString(getpid()) + ".temp", reportFileName + toString(getpid()) + ".temp", accnosFName + toString(getpid()) + ".temp", filename);
865                                 
866                                 //pass numSeqs to parent
867                                 ofstream out;
868                                 string tempFile = alignFileName + toString(getpid()) + ".num.temp";
869                                 m->openOutputFile(tempFile, out);
870                                 out << num << endl;
871                                 out.close();
872                                 
873                                 exit(0);
874                         }else { 
875                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
876                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
877                                 exit(0);
878                         }
879                 }
880                 
881                 //do my part
882                 num = driver(lines[0], alignFileName, reportFileName, accnosFName, filename);
883                 
884                 //force parent to wait until all the processes are done
885                 for (int i=0;i<processIDS.size();i++) { 
886                         int temp = processIDS[i];
887                         wait(&temp);
888                 }
889                 
890                 vector<string> nonBlankAccnosFiles;
891                 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
892                 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
893                         
894                 for (int i = 0; i < processIDS.size(); i++) {
895                         ifstream in;
896                         string tempFile =  alignFileName + toString(processIDS[i]) + ".num.temp";
897                         m->openInputFile(tempFile, in);
898                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
899                         in.close(); m->mothurRemove(tempFile);
900                         
901                         m->appendFiles((alignFileName + toString(processIDS[i]) + ".temp"), alignFileName);
902                         m->mothurRemove((alignFileName + toString(processIDS[i]) + ".temp"));
903                         
904                         appendReportFiles((reportFileName + toString(processIDS[i]) + ".temp"), reportFileName);
905                         m->mothurRemove((reportFileName + toString(processIDS[i]) + ".temp"));
906                         
907                         if (!(m->isBlank(accnosFName + toString(processIDS[i]) + ".temp"))) {
908                                 nonBlankAccnosFiles.push_back(accnosFName + toString(processIDS[i]) + ".temp");
909                         }else { m->mothurRemove((accnosFName + toString(processIDS[i]) + ".temp"));  }
910                         
911                 }
912                 
913                 //append accnos files
914                 if (nonBlankAccnosFiles.size() != 0) { 
915                         rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
916                         
917                         for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
918                                 m->appendFiles(nonBlankAccnosFiles[h], accnosFName);
919                                 m->mothurRemove(nonBlankAccnosFiles[h]);
920                         }
921                 }else { //recreate the accnosfile if needed
922                         ofstream out;
923                         m->openOutputFile(accnosFName, out);
924                         out.close();
925                 }
926 #else
927                 //////////////////////////////////////////////////////////////////////////////////////////////////////
928                 //Windows version shared memory, so be careful when passing variables through the alignData struct. 
929                 //Above fork() will clone, so memory is separate, but that's not the case with windows, 
930                 //////////////////////////////////////////////////////////////////////////////////////////////////////
931                 
932                 vector<alignData*> pDataArray; 
933                 DWORD   dwThreadIdArray[processors-1];
934                 HANDLE  hThreadArray[processors-1]; 
935                 
936                 //Create processor worker threads.
937                 for( int i=0; i<processors-1; i++ ){
938                         //copy templateDb
939                         //AlignmentDB* tempDB = new AlignmentDB(*templateDB);
940                         
941                         // Allocate memory for thread data.
942                         string extension = "";
943                         if (i != 0) { extension = toString(i) + ".temp"; }
944                         
945                         alignData* tempalign = new alignData(templateFileName, (alignFileName + extension), (reportFileName + extension), (accnosFName + extension), filename, align, search, kmerSize, m, lines[i]->start, lines[i]->end, flip, match, misMatch, gapOpen, gapExtend, threshold, i);
946                         pDataArray.push_back(tempalign);
947                         processIDS.push_back(i);
948                                 
949                         //MySeqSumThreadFunction is in header. It must be global or static to work with the threads.
950                         //default security attributes, thread function name, argument to thread function, use default creation flags, returns the thread identifier
951                         hThreadArray[i] = CreateThread(NULL, 0, MyAlignThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);   
952                 }
953                 
954                 //need to check for line ending error
955                 ifstream inFASTA;
956                 m->openInputFile(filename, inFASTA);
957                 inFASTA.seekg(lines[processors-1]->start-1);
958                 char c = inFASTA.peek();
959                 
960                 if (c != '>') { //we need to move back
961                         lines[processors-1]->start--; 
962                 }
963                 
964                 //using the main process as a worker saves time and memory
965                 //do my part - do last piece because windows is looking for eof
966                 num = driver(lines[processors-1], (alignFileName + toString(processors-1) + ".temp"), (reportFileName + toString(processors-1) + ".temp"), (accnosFName + toString(processors-1) + ".temp"), filename);
967                 
968                 //Wait until all threads have terminated.
969                 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
970                 
971                 //Close all thread handles and free memory allocations.
972                 for(int i=0; i < pDataArray.size(); i++){
973                         num += pDataArray[i]->count;
974                         CloseHandle(hThreadArray[i]);
975                         delete pDataArray[i];
976                 }
977                 
978                 vector<string> nonBlankAccnosFiles;
979                 if (!(m->isBlank(accnosFName))) { nonBlankAccnosFiles.push_back(accnosFName); }
980                 else { m->mothurRemove(accnosFName); } //remove so other files can be renamed to it
981                 
982                 for (int i = 1; i < processors; i++) {
983                         m->appendFiles((alignFileName + toString(i) + ".temp"), alignFileName);
984                         m->mothurRemove((alignFileName + toString(i) + ".temp"));
985                         
986                         appendReportFiles((reportFileName + toString(i) + ".temp"), reportFileName);
987                         m->mothurRemove((reportFileName + toString(i) + ".temp"));
988                         
989                         if (!(m->isBlank(accnosFName + toString(i) + ".temp"))) {
990                                 nonBlankAccnosFiles.push_back(accnosFName + toString(i) + ".temp");
991                         }else { m->mothurRemove((accnosFName + toString(i) + ".temp"));  }
992                 }
993                 
994                 //append accnos files
995                 if (nonBlankAccnosFiles.size() != 0) { 
996                         rename(nonBlankAccnosFiles[0].c_str(), accnosFName.c_str());
997                         
998                         for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
999                                 m->appendFiles(nonBlankAccnosFiles[h], accnosFName);
1000                                 m->mothurRemove(nonBlankAccnosFiles[h]);
1001                         }
1002                 }else { //recreate the accnosfile if needed
1003                         ofstream out;
1004                         m->openOutputFile(accnosFName, out);
1005                         out.close();
1006                 }       
1007 #endif  
1008                 
1009                 return num;
1010         }
1011         catch(exception& e) {
1012                 m->errorOut(e, "AlignCommand", "createProcesses");
1013                 exit(1);
1014         }
1015 }
1016 //**********************************************************************************************************************
1017
1018 void AlignCommand::appendReportFiles(string temp, string filename) {
1019         try{
1020                 
1021                 ofstream output;
1022                 ifstream input;
1023                 m->openOutputFileAppend(filename, output);
1024                 m->openInputFile(temp, input);
1025
1026                 while (!input.eof())    {       char c = input.get(); if (c == 10 || c == 13){  break;  }       } // get header line
1027                                 
1028         char buffer[4096];        
1029         while (!input.eof()) {
1030             input.read(buffer, 4096);
1031             output.write(buffer, input.gcount());
1032         }
1033                 
1034                 input.close();
1035                 output.close();
1036         }
1037         catch(exception& e) {
1038                 m->errorOut(e, "AlignCommand", "appendReportFiles");
1039                 exit(1);
1040         }
1041 }
1042 //**********************************************************************************************************************