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