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