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