]> git.donarmstrong.com Git - mothur.git/blob - chimerauchimecommand.cpp
adds group parameter to chimera.uchime so you can check for chimeras with template...
[mothur.git] / chimerauchimecommand.cpp
1 /*
2  *  chimerauchimecommand.cpp
3  *  Mothur
4  *
5  *  Created by westcott on 5/13/11.
6  *  Copyright 2011 Schloss Lab. All rights reserved.
7  *
8  */
9
10 #include "chimerauchimecommand.h"
11 #include "deconvolutecommand.h"
12 #include "uc.h"
13 #include "sequence.hpp"
14 #include "referencedb.h"
15
16
17 //**********************************************************************************************************************
18 vector<string> ChimeraUchimeCommand::setParameters(){   
19         try {
20                 CommandParameter ptemplate("reference", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(ptemplate);
21                 CommandParameter pfasta("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pfasta);
22                 CommandParameter pname("name", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(pname);
23                 CommandParameter pgroup("group", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(pgroup);
24                 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
25                 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
26                 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
27                 CommandParameter pabskew("abskew", "Number", "", "1.9", "", "", "",false,false); parameters.push_back(pabskew);
28                 CommandParameter pchimealns("chimealns", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(pchimealns);
29                 CommandParameter pminh("minh", "Number", "", "0.3", "", "", "",false,false); parameters.push_back(pminh);
30                 CommandParameter pmindiv("mindiv", "Number", "", "0.5", "", "", "",false,false); parameters.push_back(pmindiv);
31                 CommandParameter pxn("xn", "Number", "", "8.0", "", "", "",false,false); parameters.push_back(pxn);
32                 CommandParameter pdn("dn", "Number", "", "1.4", "", "", "",false,false); parameters.push_back(pdn);
33                 CommandParameter pxa("xa", "Number", "", "1", "", "", "",false,false); parameters.push_back(pxa);
34                 CommandParameter pchunks("chunks", "Number", "", "4", "", "", "",false,false); parameters.push_back(pchunks);
35                 CommandParameter pminchunk("minchunk", "Number", "", "64", "", "", "",false,false); parameters.push_back(pminchunk);
36                 CommandParameter pidsmoothwindow("idsmoothwindow", "Number", "", "32", "", "", "",false,false); parameters.push_back(pidsmoothwindow);
37                 //CommandParameter pminsmoothid("minsmoothid", "Number", "", "0.95", "", "", "",false,false); parameters.push_back(pminsmoothid);
38                 CommandParameter pmaxp("maxp", "Number", "", "2", "", "", "",false,false); parameters.push_back(pmaxp);
39                 CommandParameter pskipgaps("skipgaps", "Boolean", "", "T", "", "", "",false,false); parameters.push_back(pskipgaps);
40                 CommandParameter pskipgaps2("skipgaps2", "Boolean", "", "T", "", "", "",false,false); parameters.push_back(pskipgaps2);
41                 CommandParameter pminlen("minlen", "Number", "", "10", "", "", "",false,false); parameters.push_back(pminlen);
42                 CommandParameter pmaxlen("maxlen", "Number", "", "10000", "", "", "",false,false); parameters.push_back(pmaxlen);
43                 CommandParameter pucl("ucl", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(pucl);
44                 CommandParameter pqueryfract("queryfract", "Number", "", "0.5", "", "", "",false,false); parameters.push_back(pqueryfract);
45
46                 vector<string> myArray;
47                 for (int i = 0; i < parameters.size(); i++) {   myArray.push_back(parameters[i].name);          }
48                 return myArray;
49         }
50         catch(exception& e) {
51                 m->errorOut(e, "ChimeraUchimeCommand", "setParameters");
52                 exit(1);
53         }
54 }
55 //**********************************************************************************************************************
56 string ChimeraUchimeCommand::getHelpString(){   
57         try {
58                 string helpString = "";
59                 helpString += "The chimera.uchime command reads a fastafile and referencefile and outputs potentially chimeric sequences.\n";
60                 helpString += "This command is a wrapper for uchime written by Robert C. Edgar.\n";
61                 helpString += "The chimera.uchime command parameters are fasta, name, reference, processors, abskew, chimealns, minh, mindiv, xn, dn, xa, chunks, minchunk, idsmoothwindow, minsmoothid, maxp, skipgaps, skipgaps2, minlen, maxlen, ucl and queryfact.\n";
62                 helpString += "The fasta parameter allows you to enter the fasta file containing your potentially chimeric sequences, and is required, unless you have a valid current fasta file. \n";
63                 helpString += "The name parameter allows you to provide a name file, if you are using template=self. \n";
64                 helpString += "You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amazon.fasta \n";
65                 helpString += "The group parameter allows you to provide a group file. The group file can be used with a namesfile and reference=self. When checking sequences, only sequences from the same group as the query sequence will be used as the reference. \n";
66                 helpString += "The reference parameter allows you to enter a reference file containing known non-chimeric sequences, and is required. You may also set template=self, in this case the abundant sequences will be used as potential parents. \n";
67                 helpString += "The processors parameter allows you to specify how many processors you would like to use.  The default is 1. \n";
68                 helpString += "The abskew parameter can only be used with template=self. Minimum abundance skew. Default 1.9. Abundance skew is: min [ abund(parent1), abund(parent2) ] / abund(query).\n";
69                 helpString += "The chimealns parameter allows you to indicate you would like a file containing multiple alignments of query sequences to parents in human readable format. Alignments show columns with differences that support or contradict a chimeric model.\n";
70                 helpString += "The minh parameter - mininum score to report chimera. Default 0.3. Values from 0.1 to 5 might be reasonable. Lower values increase sensitivity but may report more false positives. If you decrease xn you may need to increase minh, and vice versa.\n";
71                 helpString += "The mindiv parameter - minimum divergence ratio, default 0.5. Div ratio is 100%% - %%identity between query sequence and the closest candidate for being a parent. If you don't care about very close chimeras, then you could increase mindiv to, say, 1.0 or 2.0, and also decrease minh, say to 0.1, to increase sensitivity. How well this works will depend on your data. Best is to tune parameters on a good benchmark.\n";
72                 helpString += "The xn parameter - weight of a no vote. Default 8.0. Decreasing this weight to around 3 or 4 may give better performance on denoised data.\n";
73                 helpString += "The dn parameter - pseudo-count prior on number of no votes. Default 1.4. Probably no good reason to change this unless you can retune to a good benchmark for your data. Reasonable values are probably in the range from 0.2 to 2.\n";
74                 helpString += "The xa parameter - weight of an abstain vote. Default 1. So far, results do not seem to be very sensitive to this parameter, but if you have a good training set might be worth trying. Reasonable values might range from 0.1 to 2.\n";
75                 helpString += "The chunks parameter is the number of chunks to extract from the query sequence when searching for parents. Default 4.\n";
76                 helpString += "The minchunk parameter is the minimum length of a chunk. Default 64.\n";
77                 helpString += "The idsmoothwindow parameter is the length of id smoothing window. Default 32.\n";
78                 //helpString += "The minsmoothid parameter - minimum factional identity over smoothed window of candidate parent. Default 0.95.\n";
79                 helpString += "The maxp parameter - maximum number of candidate parents to consider. Default 2. In tests so far, increasing maxp gives only a very small improvement in sensivity but tends to increase the error rate quite a bit.\n";
80                 helpString += "The skipgaps parameter controls how gapped columns affect counting of diffs. If skipgaps is set to T, columns containing gaps do not found as diffs. Default = T.\n";
81                 helpString += "The skipgaps2 parameter controls how gapped columns affect counting of diffs. If skipgaps2 is set to T, if column is immediately adjacent to a column containing a gap, it is not counted as a diff. Default = T.\n";
82                 helpString += "The minlen parameter is the minimum unaligned sequence length. Defaults 10. Applies to both query and reference sequences.\n";
83                 helpString += "The maxlen parameter is the maximum unaligned sequence length. Defaults 10000. Applies to both query and reference sequences.\n";
84                 helpString += "The ucl parameter - use local-X alignments. Default is global-X or false. On tests so far, global-X is always better; this option is retained because it just might work well on some future type of data.\n";
85                 helpString += "The queryfract parameter - minimum fraction of the query sequence that must be covered by a local-X alignment. Default 0.5. Applies only when ucl is true.\n";
86 #ifdef USE_MPI
87                 helpString += "When using MPI, the processors parameter is set to the number of MPI processes running. \n";
88 #endif
89                 helpString += "The chimera.uchime command should be in the following format: \n";
90                 helpString += "chimera.uchime(fasta=yourFastaFile, reference=yourTemplate) \n";
91                 helpString += "Example: chimera.uchime(fasta=AD.align, reference=silva.gold.align) \n";
92                 helpString += "Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n";       
93                 return helpString;
94         }
95         catch(exception& e) {
96                 m->errorOut(e, "ChimeraUchimeCommand", "getHelpString");
97                 exit(1);
98         }
99 }
100 //**********************************************************************************************************************
101 ChimeraUchimeCommand::ChimeraUchimeCommand(){   
102         try {
103                 abort = true; calledHelp = true;
104                 setParameters();
105                 vector<string> tempOutNames;
106                 outputTypes["chimera"] = tempOutNames;
107                 outputTypes["accnos"] = tempOutNames;
108                 outputTypes["alns"] = tempOutNames;
109         }
110         catch(exception& e) {
111                 m->errorOut(e, "ChimeraUchimeCommand", "ChimeraUchimeCommand");
112                 exit(1);
113         }
114 }
115 //***************************************************************************************************************
116 ChimeraUchimeCommand::ChimeraUchimeCommand(string option)  {
117         try {
118                 abort = false; calledHelp = false; 
119                 ReferenceDB* rdb = ReferenceDB::getInstance();
120                 
121                 //allow user to run help
122                 if(option == "help") { help(); abort = true; calledHelp = true; }
123                 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
124                 
125                 else {
126                         vector<string> myArray = setParameters();
127                         
128                         OptionParser parser(option);
129                         map<string,string> parameters = parser.getParameters();
130                         
131                         ValidParameters validParameter("chimera.uchime");
132                         map<string,string>::iterator it;
133                         
134                         //check to make sure all parameters are valid for command
135                         for (it = parameters.begin(); it != parameters.end(); it++) { 
136                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
137                         }
138                         
139                         vector<string> tempOutNames;
140                         outputTypes["chimera"] = tempOutNames;
141                         outputTypes["accnos"] = tempOutNames;
142                         outputTypes["alns"] = tempOutNames;
143                         
144                         //if the user changes the input directory command factory will send this info to us in the output parameter 
145                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
146                         if (inputDir == "not found"){   inputDir = "";          }
147                         
148                         //check for required parameters
149                         fastafile = validParameter.validFile(parameters, "fasta", false);
150                         if (fastafile == "not found") {                                 
151                                 //if there is a current fasta file, use it
152                                 string filename = m->getFastaFile(); 
153                                 if (filename != "") { fastaFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
154                                 else {  m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
155                         }else { 
156                                 m->splitAtDash(fastafile, fastaFileNames);
157                                 
158                                 //go through files and make sure they are good, if not, then disregard them
159                                 for (int i = 0; i < fastaFileNames.size(); i++) {
160                                         
161                                         bool ignore = false;
162                                         if (fastaFileNames[i] == "current") { 
163                                                 fastaFileNames[i] = m->getFastaFile(); 
164                                                 if (fastaFileNames[i] != "") {  m->mothurOut("Using " + fastaFileNames[i] + " as input file for the fasta parameter where you had given current."); m->mothurOutEndLine(); }
165                                                 else {  
166                                                         m->mothurOut("You have no current fastafile, ignoring current."); m->mothurOutEndLine(); ignore=true; 
167                                                         //erase from file list
168                                                         fastaFileNames.erase(fastaFileNames.begin()+i);
169                                                         i--;
170                                                 }
171                                         }
172                                         
173                                         if (!ignore) {
174                                                 
175                                                 if (inputDir != "") {
176                                                         string path = m->hasPath(fastaFileNames[i]);
177                                                         //if the user has not given a path then, add inputdir. else leave path alone.
178                                                         if (path == "") {       fastaFileNames[i] = inputDir + fastaFileNames[i];               }
179                                                 }
180                                                 
181                                                 int ableToOpen;
182                                                 ifstream in;
183                                                 
184                                                 ableToOpen = m->openInputFile(fastaFileNames[i], in, "noerror");
185                                                 
186                                                 //if you can't open it, try default location
187                                                 if (ableToOpen == 1) {
188                                                         if (m->getDefaultPath() != "") { //default path is set
189                                                                 string tryPath = m->getDefaultPath() + m->getSimpleName(fastaFileNames[i]);
190                                                                 m->mothurOut("Unable to open " + fastaFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
191                                                                 ifstream in2;
192                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
193                                                                 in2.close();
194                                                                 fastaFileNames[i] = tryPath;
195                                                         }
196                                                 }
197                                                 
198                                                 if (ableToOpen == 1) {
199                                                         if (m->getOutputDir() != "") { //default path is set
200                                                                 string tryPath = m->getOutputDir() + m->getSimpleName(fastaFileNames[i]);
201                                                                 m->mothurOut("Unable to open " + fastaFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
202                                                                 ifstream in2;
203                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
204                                                                 in2.close();
205                                                                 fastaFileNames[i] = tryPath;
206                                                         }
207                                                 }
208                                                 
209                                                 in.close();
210                                                 
211                                                 if (ableToOpen == 1) { 
212                                                         m->mothurOut("Unable to open " + fastaFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
213                                                         //erase from file list
214                                                         fastaFileNames.erase(fastaFileNames.begin()+i);
215                                                         i--;
216                                                 }else {
217                                                         m->setFastaFile(fastaFileNames[i]);
218                                                 }
219                                         }
220                                 }
221                                 
222                                 //make sure there is at least one valid file left
223                                 if (fastaFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid files."); m->mothurOutEndLine(); abort = true; }
224                         }
225                         
226                         
227                         //check for required parameters
228                         bool hasName = true;
229                         namefile = validParameter.validFile(parameters, "name", false);
230                         if (namefile == "not found") { namefile = "";  hasName = false; }
231                         else { 
232                                 m->splitAtDash(namefile, nameFileNames);
233                                 
234                                 //go through files and make sure they are good, if not, then disregard them
235                                 for (int i = 0; i < nameFileNames.size(); i++) {
236                                         
237                                         bool ignore = false;
238                                         if (nameFileNames[i] == "current") { 
239                                                 nameFileNames[i] = m->getNameFile(); 
240                                                 if (nameFileNames[i] != "") {  m->mothurOut("Using " + nameFileNames[i] + " as input file for the name parameter where you had given current."); m->mothurOutEndLine(); }
241                                                 else {  
242                                                         m->mothurOut("You have no current namefile, ignoring current."); m->mothurOutEndLine(); ignore=true; 
243                                                         //erase from file list
244                                                         nameFileNames.erase(nameFileNames.begin()+i);
245                                                         i--;
246                                                 }
247                                         }
248                                         
249                                         if (!ignore) {
250                                                 
251                                                 if (inputDir != "") {
252                                                         string path = m->hasPath(nameFileNames[i]);
253                                                         //if the user has not given a path then, add inputdir. else leave path alone.
254                                                         if (path == "") {       nameFileNames[i] = inputDir + nameFileNames[i];         }
255                                                 }
256                                                 
257                                                 int ableToOpen;
258                                                 ifstream in;
259                                                 
260                                                 ableToOpen = m->openInputFile(nameFileNames[i], in, "noerror");
261                                                 
262                                                 //if you can't open it, try default location
263                                                 if (ableToOpen == 1) {
264                                                         if (m->getDefaultPath() != "") { //default path is set
265                                                                 string tryPath = m->getDefaultPath() + m->getSimpleName(nameFileNames[i]);
266                                                                 m->mothurOut("Unable to open " + nameFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
267                                                                 ifstream in2;
268                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
269                                                                 in2.close();
270                                                                 nameFileNames[i] = tryPath;
271                                                         }
272                                                 }
273                                                 
274                                                 if (ableToOpen == 1) {
275                                                         if (m->getOutputDir() != "") { //default path is set
276                                                                 string tryPath = m->getOutputDir() + m->getSimpleName(nameFileNames[i]);
277                                                                 m->mothurOut("Unable to open " + nameFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
278                                                                 ifstream in2;
279                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
280                                                                 in2.close();
281                                                                 nameFileNames[i] = tryPath;
282                                                         }
283                                                 }
284                                                 
285                                                 in.close();
286                                                 
287                                                 if (ableToOpen == 1) { 
288                                                         m->mothurOut("Unable to open " + nameFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
289                                                         //erase from file list
290                                                         nameFileNames.erase(nameFileNames.begin()+i);
291                                                         i--;
292                                                 }else {
293                                                         m->setNameFile(nameFileNames[i]);
294                                                 }
295                                         }
296                                 }
297                                 
298                                 //make sure there is at least one valid file left
299                                 if (nameFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid name files."); m->mothurOutEndLine(); abort = true; }
300                         }
301                         
302                         if (hasName && (nameFileNames.size() != fastaFileNames.size())) { m->mothurOut("[ERROR]: The number of namefiles does not match the number of fastafiles, please correct."); m->mothurOutEndLine(); abort=true; }
303                         
304                         bool hasGroup = true;
305                         groupfile = validParameter.validFile(parameters, "group", false);
306                         if (groupfile == "not found") { groupfile = "";  hasGroup = false; }
307                         else { 
308                                 m->splitAtDash(groupfile, groupFileNames);
309                                 
310                                 //go through files and make sure they are good, if not, then disregard them
311                                 for (int i = 0; i < groupFileNames.size(); i++) {
312                                         
313                                         bool ignore = false;
314                                         if (groupFileNames[i] == "current") { 
315                                                 groupFileNames[i] = m->getGroupFile(); 
316                                                 if (groupFileNames[i] != "") {  m->mothurOut("Using " + groupFileNames[i] + " as input file for the group parameter where you had given current."); m->mothurOutEndLine(); }
317                                                 else {  
318                                                         m->mothurOut("You have no current namefile, ignoring current."); m->mothurOutEndLine(); ignore=true; 
319                                                         //erase from file list
320                                                         groupFileNames.erase(groupFileNames.begin()+i);
321                                                         i--;
322                                                 }
323                                         }
324                                         
325                                         if (!ignore) {
326                                                 
327                                                 if (inputDir != "") {
328                                                         string path = m->hasPath(groupFileNames[i]);
329                                                         //if the user has not given a path then, add inputdir. else leave path alone.
330                                                         if (path == "") {       groupFileNames[i] = inputDir + groupFileNames[i];               }
331                                                 }
332                                                 
333                                                 int ableToOpen;
334                                                 ifstream in;
335                                                 
336                                                 ableToOpen = m->openInputFile(groupFileNames[i], in, "noerror");
337                                                 
338                                                 //if you can't open it, try default location
339                                                 if (ableToOpen == 1) {
340                                                         if (m->getDefaultPath() != "") { //default path is set
341                                                                 string tryPath = m->getDefaultPath() + m->getSimpleName(groupFileNames[i]);
342                                                                 m->mothurOut("Unable to open " + groupFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
343                                                                 ifstream in2;
344                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
345                                                                 in2.close();
346                                                                 groupFileNames[i] = tryPath;
347                                                         }
348                                                 }
349                                                 
350                                                 if (ableToOpen == 1) {
351                                                         if (m->getOutputDir() != "") { //default path is set
352                                                                 string tryPath = m->getOutputDir() + m->getSimpleName(groupFileNames[i]);
353                                                                 m->mothurOut("Unable to open " + groupFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
354                                                                 ifstream in2;
355                                                                 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
356                                                                 in2.close();
357                                                                 groupFileNames[i] = tryPath;
358                                                         }
359                                                 }
360                                                 
361                                                 in.close();
362                                                 
363                                                 if (ableToOpen == 1) { 
364                                                         m->mothurOut("Unable to open " + groupFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
365                                                         //erase from file list
366                                                         groupFileNames.erase(groupFileNames.begin()+i);
367                                                         i--;
368                                                 }else {
369                                                         m->setGroupFile(groupFileNames[i]);
370                                                 }
371                                         }
372                                 }
373                                 
374                                 //make sure there is at least one valid file left
375                                 if (groupFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid group files."); m->mothurOutEndLine(); abort = true; }
376                         }
377                         
378                         if (hasGroup && (groupFileNames.size() != fastaFileNames.size())) { m->mothurOut("[ERROR]: The number of groupfiles does not match the number of fastafiles, please correct."); m->mothurOutEndLine(); abort=true; }
379                         
380                         
381                         //if the user changes the output directory command factory will send this info to us in the output parameter 
382                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  outputDir = ""; }
383                         
384                         string path;
385                         it = parameters.find("reference");
386                         //user has given a template file
387                         if(it != parameters.end()){ 
388                                 if (it->second == "self") { templatefile = "self"; }
389                                 else {
390                                         path = m->hasPath(it->second);
391                                         //if the user has not given a path then, add inputdir. else leave path alone.
392                                         if (path == "") {       parameters["reference"] = inputDir + it->second;                }
393                                         
394                                         templatefile = validParameter.validFile(parameters, "reference", true);
395                                         if (templatefile == "not open") { abort = true; }
396                                         else if (templatefile == "not found") { //check for saved reference sequences
397                                                 if (rdb->getSavedReference() != "") {
398                                                         templatefile = rdb->getSavedReference();
399                                                         m->mothurOutEndLine();  m->mothurOut("Using sequences from " + rdb->getSavedReference() + "."); m->mothurOutEndLine();
400                                                 }else {
401                                                         m->mothurOut("[ERROR]: You don't have any saved reference sequences and the reference parameter is a required."); 
402                                                         m->mothurOutEndLine();
403                                                         abort = true; 
404                                                 }
405                                         }
406                                 }
407                         }else if (hasName) {  templatefile = "self"; }
408                         else { 
409                                 if (rdb->getSavedReference() != "") {
410                                         templatefile = rdb->getSavedReference();
411                                         m->mothurOutEndLine();  m->mothurOut("Using sequences from " + rdb->getSavedReference() + "."); m->mothurOutEndLine();
412                                 }else {
413                                         m->mothurOut("[ERROR]: You don't have any saved reference sequences and the reference parameter is a required."); 
414                                         m->mothurOutEndLine();
415                                         templatefile = ""; abort = true; 
416                                 } 
417                         }
418                                 
419                         string temp = validParameter.validFile(parameters, "processors", false);        if (temp == "not found"){       temp = m->getProcessors();      }
420                         m->setProcessors(temp);
421                         convert(temp, processors);
422                         
423                         abskew = validParameter.validFile(parameters, "abskew", false); if (abskew == "not found"){     useAbskew = false;  abskew = "1.9";     }else{  useAbskew = true;  }
424                         if (useAbskew && templatefile != "self") { m->mothurOut("The abskew parameter is only valid with template=self, ignoring."); m->mothurOutEndLine(); useAbskew = false; }
425                         
426                         temp = validParameter.validFile(parameters, "chimealns", false);                        if (temp == "not found") { temp = "f"; }
427                         chimealns = m->isTrue(temp); 
428                         
429                         minh = validParameter.validFile(parameters, "minh", false);                                             if (minh == "not found")                        { useMinH = false; minh = "0.3";                                        }       else{ useMinH = true;                   }
430                         mindiv = validParameter.validFile(parameters, "mindiv", false);                                 if (mindiv == "not found")                      { useMindiv = false; mindiv = "0.5";                            }       else{ useMindiv = true;                 }
431                         xn = validParameter.validFile(parameters, "xn", false);                                                 if (xn == "not found")                          { useXn = false; xn = "8.0";                                            }       else{ useXn = true;                             }
432                         dn = validParameter.validFile(parameters, "dn", false);                                                 if (dn == "not found")                          { useDn = false; dn = "1.4";                                            }       else{ useDn = true;                             }
433                         xa = validParameter.validFile(parameters, "xa", false);                                                 if (xa == "not found")                          { useXa = false; xa = "1";                                                      }       else{ useXa = true;                             }
434                         chunks = validParameter.validFile(parameters, "chunks", false);                                 if (chunks == "not found")                      { useChunks = false; chunks = "4";                                      }       else{ useChunks = true;                 }
435                         minchunk = validParameter.validFile(parameters, "minchunk", false);                             if (minchunk == "not found")            { useMinchunk = false; minchunk = "64";                         }       else{ useMinchunk = true;               }
436                         idsmoothwindow = validParameter.validFile(parameters, "idsmoothwindow", false); if (idsmoothwindow == "not found")      { useIdsmoothwindow = false; idsmoothwindow = "32";     }       else{ useIdsmoothwindow = true; }
437                         //minsmoothid = validParameter.validFile(parameters, "minsmoothid", false);             if (minsmoothid == "not found")         { useMinsmoothid = false; minsmoothid = "0.95";         }       else{ useMinsmoothid = true;    }
438                         maxp = validParameter.validFile(parameters, "maxp", false);                                             if (maxp == "not found")                        { useMaxp = false; maxp = "2";                                          }       else{ useMaxp = true;                   }
439                         minlen = validParameter.validFile(parameters, "minlen", false);                                 if (minlen == "not found")                      { useMinlen = false; minlen = "10";                                     }       else{ useMinlen = true;                 }
440                         maxlen = validParameter.validFile(parameters, "maxlen", false);                                 if (maxlen == "not found")                      { useMaxlen = false; maxlen = "10000";                          }       else{ useMaxlen = true;                 }
441                         
442                         temp = validParameter.validFile(parameters, "ucl", false);                                              if (temp == "not found") { temp = "f"; }
443                         ucl = m->isTrue(temp);
444                         
445                         queryfract = validParameter.validFile(parameters, "queryfract", false);                 if (queryfract == "not found")          { useQueryfract = false; queryfract = "0.5";            }       else{ useQueryfract = true;             }
446                         if (!ucl && useQueryfract) { m->mothurOut("queryfact may only be used when ucl=t, ignoring."); m->mothurOutEndLine(); useQueryfract = false; }
447                         
448                         temp = validParameter.validFile(parameters, "skipgaps", false);                                 if (temp == "not found") { temp = "t"; }
449                         skipgaps = m->isTrue(temp); 
450
451                         temp = validParameter.validFile(parameters, "skipgaps2", false);                                if (temp == "not found") { temp = "t"; }
452                         skipgaps2 = m->isTrue(temp); 
453                         
454                         if (hasName && (templatefile != "self")) { m->mothurOut("You have provided a namefile and the reference parameter is not set to self. I am not sure what reference you are trying to use, aborting."); m->mothurOutEndLine(); abort=true; }
455                         if (hasGroup && (templatefile != "self")) { m->mothurOut("You have provided a group file and the reference parameter is not set to self. I am not sure what reference you are trying to use, aborting."); m->mothurOutEndLine(); abort=true; }
456
457                 }
458         }
459         catch(exception& e) {
460                 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
461                 exit(1);
462         }
463 }
464 //***************************************************************************************************************
465
466 int ChimeraUchimeCommand::execute(){
467         try{
468                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
469                 
470                 m->mothurOut("\nuchime by Robert C. Edgar\nhttp://drive5.com/uchime\nThis code is donated to the public domain.\n\n");
471                 
472                 for (int s = 0; s < fastaFileNames.size(); s++) {
473                         
474                         m->mothurOut("Checking sequences from " + fastaFileNames[s] + " ..." ); m->mothurOutEndLine();
475                         
476                         int start = time(NULL); 
477                         string nameFile = "";
478                         if (outputDir == "") { outputDir = m->hasPath(fastaFileNames[s]);  }//if user entered a file with a path then preserve it                               
479                         string outputFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s])) + "uchime.chimera";
480                         string accnosFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s]))  + "uchime.accnos";
481                         string alnsFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s]))  + "uchime.alns";
482                         string newFasta = m->getRootName(fastaFileNames[s]) + "temp";
483                                 
484                         //you provided a groupfile
485                         string groupFile = "";
486                         if (groupFileNames.size() != 0) { groupFile = groupFileNames[s]; }
487                         
488                         if ((templatefile == "self") && (groupFile == "")) { //you want to run uchime with a reference template
489
490                                 if (processors != 1) { m->mothurOut("When using template=self, mothur can only use 1 processor, continuing."); m->mothurOutEndLine(); processors = 1; }
491                                 if (nameFileNames.size() != 0) { //you provided a namefile and we don't need to create one
492                                         nameFile = nameFileNames[s];
493                                 }else { nameFile = getNamesFile(fastaFileNames[s]); }
494                                                                         
495                                 map<string, string> seqs;  
496                                 readFasta(fastaFileNames[s], seqs);  if (m->control_pressed) { for (int j = 0; j < outputNames.size(); j++) {   m->mothurRemove(outputNames[j]);        }  return 0; }
497
498                                 //read namefile
499                                 vector<seqPriorityNode> nameMapCount;
500                                 int error = m->readNames(nameFile, nameMapCount, seqs); if (m->control_pressed) { for (int j = 0; j < outputNames.size(); j++) {        m->mothurRemove(outputNames[j]);        }  return 0; }
501                                 if (error == 1) { for (int j = 0; j < outputNames.size(); j++) {        m->mothurRemove(outputNames[j]);        }  return 0; }
502                                 if (seqs.size() != nameMapCount.size()) { m->mothurOut( "The number of sequences in your fastafile does not match the number of sequences in your namefile, aborting."); m->mothurOutEndLine(); for (int j = 0; j < outputNames.size(); j++) {  m->mothurRemove(outputNames[j]);        }  return 0; }
503                                 
504                                 printFile(nameMapCount, newFasta);
505                                 fastaFileNames[s] = newFasta;
506                         }
507                         
508                         if (m->control_pressed) {  for (int j = 0; j < outputNames.size(); j++) {       m->mothurRemove(outputNames[j]);        }  return 0;    }                               
509                         
510                         if (groupFile != "") {
511                                 if (nameFileNames.size() != 0) { //you provided a namefile and we don't need to create one
512                                         nameFile = nameFileNames[s];
513                                 }else { nameFile = getNamesFile(fastaFileNames[s]); }
514                                 
515                                 //Parse sequences by group
516                                 SequenceParser parser(groupFile, fastaFileNames[s], nameFile);
517                                 vector<string> groups = parser.getNamesOfGroups();
518                                         
519                                 if (m->control_pressed) { for (int j = 0; j < outputNames.size(); j++) {        m->mothurRemove(outputNames[j]);        }  return 0; }
520                                                                 
521                                 //clears files
522                                 ofstream out, out1, out2;
523                                 m->openOutputFile(outputFileName, out); out.close(); 
524                                 m->openOutputFile(accnosFileName, out1); out1.close();
525                                 if (chimealns) { m->openOutputFile(alnsFileName, out2); out2.close(); }
526                                 int totalSeqs = 0;
527                                 
528         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
529                                 if(processors == 1)     {       totalSeqs = driverGroups(parser, outputFileName, newFasta, accnosFileName, alnsFileName, 0, groups.size(), groups);     }
530                                 else                            {       totalSeqs = createProcessesGroups(parser, outputFileName, newFasta, accnosFileName, alnsFileName, groups);                      }
531         #else
532                                 totalSeqs = driverGroups(parser, outputFileName, newFasta, accnosFileName, alnsFileName, 0, groups.size(), groups);
533         #endif
534                                 if (m->control_pressed) {  for (int j = 0; j < outputNames.size(); j++) {       m->mothurRemove(outputNames[j]);        }  return 0;    }                               
535
536                                 int totalChimeras = deconvoluteResults(parser, outputFileName, accnosFileName, alnsFileName);
537                                 
538                                 m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(totalSeqs) + " sequences. " + toString(totalChimeras) + " chimeras were found.");  m->mothurOutEndLine();
539                                 m->mothurOut("The number of sequences checked may be larger than the number of unique sequences because some sequences are found in several samples."); m->mothurOutEndLine(); 
540                                 
541                                 if (m->control_pressed) {  for (int j = 0; j < outputNames.size(); j++) {       m->mothurRemove(outputNames[j]);        }  return 0;    }                               
542                                         
543                         }else{
544                                 if (m->control_pressed) {  for (int j = 0; j < outputNames.size(); j++) {       m->mothurRemove(outputNames[j]);        }  return 0;    }
545                         
546                                 int numSeqs = 0;
547                                 int numChimeras = 0;
548         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
549                                 if(processors == 1){ numSeqs = driver(outputFileName, fastaFileNames[s], accnosFileName, alnsFileName, numChimeras); }
550                                 else{   numSeqs = createProcesses(outputFileName, fastaFileNames[s], accnosFileName, alnsFileName, numChimeras); }
551         #else
552                                 numSeqs = driver(outputFileName, fastaFileNames[s], accnosFileName, alnsFileName, numChimeras);
553         #endif
554                                 if (m->control_pressed) { for (int j = 0; j < outputNames.size(); j++) {        m->mothurRemove(outputNames[j]);        } return 0; }
555                         
556                                 //remove file made for uchime
557                                 if (templatefile == "self") {  m->mothurRemove(fastaFileNames[s]); }
558                         
559                                 m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences. " + toString(numChimeras) + " chimeras were found.");      m->mothurOutEndLine();
560                         }
561                         
562                         outputNames.push_back(outputFileName); outputTypes["chimera"].push_back(outputFileName);
563                         outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName);
564                         if (chimealns) { outputNames.push_back(alnsFileName); outputTypes["alns"].push_back(alnsFileName); }
565                 }
566         
567                 //set accnos file as new current accnosfile
568                 string current = "";
569                 itTypes = outputTypes.find("accnos");
570                 if (itTypes != outputTypes.end()) {
571                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setAccnosFile(current); }
572                 }
573                 
574                 m->mothurOutEndLine();
575                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
576                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }       
577                 m->mothurOutEndLine();
578                 
579                 return 0;
580                 
581         }
582         catch(exception& e) {
583                 m->errorOut(e, "ChimeraUchimeCommand", "execute");
584                 exit(1);
585         }
586 }
587 //**********************************************************************************************************************
588 int ChimeraUchimeCommand::deconvoluteResults(SequenceParser& parser, string outputFileName, string accnosFileName, string alnsFileName){
589         try {
590                 map<string, string> uniqueNames = parser.getAllSeqsMap();
591                 map<string, string>::iterator itUnique;
592                 int total = 0;
593                 
594                 //edit chimera file
595                 ifstream in; 
596                 m->openInputFile(outputFileName, in);
597                 
598                 ofstream out;
599                 m->openOutputFile(outputFileName+".temp", out); out.setf(ios::fixed, ios::floatfield); out.setf(ios::showpoint);
600                 
601                 float temp1;
602                 string name, rest, parent1, parent2;
603                 set<string> namesInFile; //this is so if a sequence is found to be chimera in several samples we dont write it to the results file more than once
604                 set<string>::iterator itNames;
605                 
606                 //assumptions - in file each read will always look like - if uchime source is updated, revisit this code.
607                 /*
608                  0.000000       F11Fcsw_33372/ab=18/            *       *       *       *       *       *       *       *       *       *       *       *       *       *       N
609                  0.018300       F11Fcsw_14980/ab=16/            F11Fcsw_1915/ab=35/     F11Fcsw_6032/ab=42/     79.9    78.7    78.2    78.7    79.2    3       0       5       11      10      20      1.46    N
610                 */
611                 
612                 while (!in.eof()) {
613                         
614                         if (m->control_pressed) { in.close(); out.close(); m->mothurRemove((outputFileName+".temp")); return 0; }
615                         
616                         in >> temp1;    m->gobble(in);
617                         in >> name;             m->gobble(in);
618                         in >> parent1;  m->gobble(in);
619                         in >> parent2;  m->gobble(in);
620                         rest = m->getline(in); m->gobble(in);
621                         
622                         //parse name - name will look like U68590/ab=1/
623                         string restOfName = "";
624                         int pos = name.find_first_of('/');
625                         if (pos != string::npos) {
626                                 restOfName = name.substr(pos);
627                                 name = name.substr(0, pos);
628                         }
629                         
630                         //find unique name
631                         itUnique = uniqueNames.find(name);
632                         
633                         if (itUnique == uniqueNames.end()) { m->mothurOut("[ERROR]: trouble parsing chimera results. Cannot find "+ name + "."); m->mothurOutEndLine(); m->control_pressed = true; }
634                         else {
635                                 itNames = namesInFile.find((itUnique->second));
636                                 
637                                 if (itNames == namesInFile.end()) {
638                                         out << temp1 << '\t' << itUnique->second << restOfName << '\t';
639                                         namesInFile.insert((itUnique->second));
640                                         
641                                         //parse parent1 names
642                                         if (parent1 != "*") {
643                                                 restOfName = "";
644                                                 pos = parent1.find_first_of('/');
645                                                 if (pos != string::npos) {
646                                                         restOfName = parent1.substr(pos);
647                                                         parent1 = parent1.substr(0, pos);
648                                                 }
649                                                 
650                                                 itUnique = uniqueNames.find(parent1);
651                                                 if (itUnique == uniqueNames.end()) { m->mothurOut("[ERROR]: trouble parsing chimera results. Cannot find parentA "+ parent1 + "."); m->mothurOutEndLine(); m->control_pressed = true; }
652                                                 else {
653                                                         out << itUnique->second << restOfName << '\t';
654                                                 }
655                                         }else { out << parent1 << '\t'; }
656                                         
657                                         //parse parent2 names
658                                         if (parent2 != "*") {
659                                                 restOfName = "";
660                                                 pos = parent2.find_first_of('/');
661                                                 if (pos != string::npos) {
662                                                         restOfName = parent2.substr(pos);
663                                                         parent2 = parent2.substr(0, pos);
664                                                 }
665                                                 
666                                                 itUnique = uniqueNames.find(parent2);
667                                                 if (itUnique == uniqueNames.end()) { m->mothurOut("[ERROR]: trouble parsing chimera results. Cannot find parentB "+ parent2 + "."); m->mothurOutEndLine(); m->control_pressed = true; }
668                                                 else {
669                                                         out << itUnique->second << restOfName << '\t';
670                                                 }
671                                         }else { out << parent2 << '\t'; }
672                                         
673                                         out  << rest << endl;
674                                 }
675                         }
676                 }
677                 in.close();
678                 out.close();
679                 
680                 m->mothurRemove(outputFileName);
681                 rename((outputFileName+".temp").c_str(), outputFileName.c_str());
682                 
683                 //edit accnos file
684                 ifstream in2; 
685                 m->openInputFile(accnosFileName, in2);
686                 
687                 ofstream out2;
688                 m->openOutputFile(accnosFileName+".temp", out2);
689                 
690                 name = "";
691                 namesInFile.clear();
692                 
693                 while (!in2.eof()) {
694                         if (m->control_pressed) { in2.close(); out2.close(); m->mothurRemove(outputFileName); m->mothurRemove((accnosFileName+".temp")); return 0; }
695                         
696                         in2 >> name; m->gobble(in2);
697                         
698                         //find unique name
699                         itUnique = uniqueNames.find(name);
700                         
701                         if (itUnique == uniqueNames.end()) { m->mothurOut("[ERROR]: trouble parsing accnos results. Cannot find "+ name + "."); m->mothurOutEndLine(); m->control_pressed = true; }
702                         else {
703                                 itNames = namesInFile.find((itUnique->second));
704         
705                                 if (itNames == namesInFile.end()) {
706                                         out2 << itUnique->second << endl;
707                                         namesInFile.insert((itUnique->second));
708                                         total++;
709                                 }
710                         }
711                 }
712                 in2.close();
713                 out2.close();
714                 
715                 m->mothurRemove(accnosFileName);
716                 rename((accnosFileName+".temp").c_str(), accnosFileName.c_str());
717                 
718                 //edit anls file
719                 //assumptions - in file each read will always look like - if uchime source is updated, revisit this code.
720                 /*
721                  ------------------------------------------------------------------------
722                  Query   (  179 nt) F21Fcsw_11639/ab=591/
723                  ParentA (  179 nt) F11Fcsw_6529/ab=1625/
724                  ParentB (  181 nt) F21Fcsw_12128/ab=1827/
725                  
726                  A     1 AAGgAAGAtTAATACaagATGgCaTCatgAGtccgCATgTtcAcatGATTAAAG--gTaTtcCGGTagacGATGGGGATG 78
727                  Q     1 AAGTAAGACTAATACCCAATGACGTCTCTAGAAGACATCTGAAAGAGATTAAAG--ATTTATCGGTGATGGATGGGGATG 78
728                  B     1 AAGgAAGAtTAATcCaggATGggaTCatgAGttcACATgTccgcatGATTAAAGgtATTTtcCGGTagacGATGGGGATG 80
729                  Diffs      N    N    A N?N   N N  NNN  N?NB   N ?NaNNN          B B NN    NNNN          
730                  Votes      0    0    + 000   0 0  000  000+   0 00!000            + 00    0000          
731                  Model   AAAAAAAAAAAAAAAAAAAAAAxBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBB
732                  
733                  A    79 CGTtccATTAGaTaGTaGGCGGGGTAACGGCCCACCtAGtCttCGATggaTAGGGGTTCTGAGAGGAAGGTCCCCCACAT 158
734                  Q    79 CGTCTGATTAGCTTGTTGGCGGGGTAACGGCCCACCAAGGCAACGATCAGTAGGGGTTCTGAGAGGAAGGTCCCCCACAT 158
735                  B    81 CGTtccATTAGaTaGTaGGCGGGGTAACGGCCCACCtAGtCAACGATggaTAGGGGTTCTGAGAGGAAGGTCCCCCACAT 160
736                  Diffs      NNN     N N  N                   N  N BB    NNN                              
737                  Votes      000     0 0  0                   0  0 ++    000                              
738                  Model   BBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBBB
739                  
740                  A   159 TGGAACTGAGACACGGTCCAA 179
741                  Q   159 TGGAACTGAGACACGGTCCAA 179
742                  B   161 TGGAACTGAGACACGGTCCAA 181
743                  Diffs                        
744                  Votes                        
745                  Model   BBBBBBBBBBBBBBBBBBBBB
746                  
747                  Ids.  QA 76.6%, QB 77.7%, AB 93.7%, QModel 78.9%, Div. +1.5%
748                  Diffs Left 7: N 0, A 6, Y 1 (14.3%); Right 35: N 1, A 30, Y 4 (11.4%), Score 0.0047
749                 */
750                 if (chimealns) {
751                         ifstream in3; 
752                         m->openInputFile(alnsFileName, in3);
753                 
754                         ofstream out3;
755                         m->openOutputFile(alnsFileName+".temp", out3); out3.setf(ios::fixed, ios::floatfield); out3.setf(ios::showpoint);
756                 
757                         name = "";
758                         namesInFile.clear();
759                         string line = "";
760                         
761                         while (!in3.eof()) {
762                                 if (m->control_pressed) { in3.close(); out3.close(); m->mothurRemove(outputFileName); m->mothurRemove((accnosFileName)); m->mothurRemove((alnsFileName+".temp")); return 0; }
763                                 
764                                 line = "";
765                                 line = m->getline(in3); 
766                                 string temp = "";
767                                 
768                                 if (line != "") {
769                                         istringstream iss(line);
770                                         iss >> temp;
771                                         
772                                         //are you a name line
773                                         if ((temp == "Query") || (temp == "ParentA") || (temp == "ParentB")) {
774                                                 int spot = 0;
775                                                 for (int i = 0; i < line.length(); i++) {
776                                                         spot = i;
777                                                         if (line[i] == ')') { break; }
778                                                         else { out3 << line[i]; }
779                                                 }
780                                                 
781                                                 if (spot == (line.length() - 1)) { m->mothurOut("[ERROR]: could not line sequence name in line " + line + "."); m->mothurOutEndLine(); m->control_pressed = true; }
782                                                 else if ((spot+2) > (line.length() - 1)) { m->mothurOut("[ERROR]: could not line sequence name in line " + line + "."); m->mothurOutEndLine(); m->control_pressed = true; }
783                                                 else {
784                                                         out << line[spot] << line[spot+1];
785                                                         
786                                                         name = line.substr(spot+2);
787                                                         
788                                                         //parse name - name will either look like U68590/ab=1/ or U68590
789                                                         string restOfName = "";
790                                                         int pos = name.find_first_of('/');
791                                                         if (pos != string::npos) {
792                                                                 restOfName = name.substr(pos);
793                                                                 name = name.substr(0, pos);
794                                                         }
795                                                         
796                                                         //find unique name
797                                                         itUnique = uniqueNames.find(name);
798                                                         
799                                                         if (itUnique == uniqueNames.end()) { m->mothurOut("[ERROR]: trouble parsing alns results. Cannot find "+ name + "."); m->mothurOutEndLine();m->control_pressed = true;  }
800                                                         else {
801                                                                 //only limit repeats on query names
802                                                                 if (temp == "Query") {
803                                                                         itNames = namesInFile.find((itUnique->second));
804                                                                         
805                                                                         if (itNames == namesInFile.end()) {
806                                                                                 out << itUnique->second << restOfName << endl;
807                                                                                 namesInFile.insert((itUnique->second));
808                                                                         }
809                                                                 }else { out << itUnique->second << restOfName << endl;  }
810                                                         }
811                                                         
812                                                 }
813                                                 
814                                         }else { //not need to alter line
815                                                 out3 << line << endl;
816                                         }
817                                 }else { out3 << endl; }
818                         }
819                         in3.close();
820                         out3.close();
821                         
822                         m->mothurRemove(alnsFileName);
823                         rename((alnsFileName+".temp").c_str(), alnsFileName.c_str());
824                 }
825                 
826                 return total;
827         }
828         catch(exception& e) {
829                 m->errorOut(e, "ChimeraUchimeCommand", "deconvoluteResults");
830                 exit(1);
831         }
832 }       
833 //**********************************************************************************************************************
834 int ChimeraUchimeCommand::printFile(vector<seqPriorityNode>& nameMapCount, string filename){
835         try {
836                 
837                 sort(nameMapCount.begin(), nameMapCount.end(), compareSeqPriorityNodes);
838                 
839                 ofstream out;
840                 m->openOutputFile(filename, out);
841                 
842                 //print new file in order of
843                 for (int i = 0; i < nameMapCount.size(); i++) {
844                         out << ">" << nameMapCount[i].name  << "/ab=" << nameMapCount[i].numIdentical << "/" << endl << nameMapCount[i].seq << endl;
845                 }
846                 out.close();
847                 
848                 return 0;
849         }
850         catch(exception& e) {
851                 m->errorOut(e, "ChimeraUchimeCommand", "printFile");
852                 exit(1);
853         }
854 }       
855 //**********************************************************************************************************************
856 int ChimeraUchimeCommand::readFasta(string filename, map<string, string>& seqs){
857         try {
858                 //create input file for uchime
859                 //read through fastafile and store info
860                 ifstream in;
861                 m->openInputFile(filename, in);
862                 
863                 while (!in.eof()) {
864                         
865                         if (m->control_pressed) { in.close(); return 0; }
866                         
867                         Sequence seq(in); m->gobble(in);
868                         seqs[seq.getName()] = seq.getAligned();
869                 }
870                 in.close();
871                 
872                 return 0;
873         }
874         catch(exception& e) {
875                 m->errorOut(e, "ChimeraUchimeCommand", "readFasta");
876                 exit(1);
877         }
878 }       
879 //**********************************************************************************************************************
880
881 string ChimeraUchimeCommand::getNamesFile(string& inputFile){
882         try {
883                 string nameFile = "";
884                 
885                 m->mothurOutEndLine(); m->mothurOut("No namesfile given, running unique.seqs command to generate one."); m->mothurOutEndLine(); m->mothurOutEndLine();
886                 
887                 //use unique.seqs to create new name and fastafile
888                 string inputString = "fasta=" + inputFile;
889                 m->mothurOut("/******************************************/"); m->mothurOutEndLine(); 
890                 m->mothurOut("Running command: unique.seqs(" + inputString + ")"); m->mothurOutEndLine(); 
891                 
892                 Command* uniqueCommand = new DeconvoluteCommand(inputString);
893                 uniqueCommand->execute();
894                 
895                 map<string, vector<string> > filenames = uniqueCommand->getOutputFiles();
896                 
897                 delete uniqueCommand;
898                 
899                 m->mothurOut("/******************************************/"); m->mothurOutEndLine(); 
900                 
901                 nameFile = filenames["name"][0];
902                 inputFile = filenames["fasta"][0];
903                 
904                 return nameFile;
905         }
906         catch(exception& e) {
907                 m->errorOut(e, "ChimeraUchimeCommand", "getNamesFile");
908                 exit(1);
909         }
910 }
911 //**********************************************************************************************************************
912 int ChimeraUchimeCommand::driverGroups(SequenceParser& parser, string outputFName, string filename, string accnos, string alns, int start, int end, vector<string> groups){
913         try {
914                 
915                 int totalSeqs = 0;
916                 int numChimeras = 0;
917                 
918                 for (int i = start; i < end; i++) {
919                         int start = time(NULL);  if (m->control_pressed) {  return 0; }
920                         
921                         int error = parser.getSeqs(groups[i], filename, true); if ((error == 1) || m->control_pressed) {  return 0; }
922                         
923                         int numSeqs = driver((outputFName + groups[i]), filename, (accnos+ groups[i]), (alns+ groups[i]), numChimeras);
924                         totalSeqs += numSeqs;
925                         
926                         if (m->control_pressed) { return 0; }
927                         
928                         //remove file made for uchime
929                         m->mothurRemove(filename);
930                         
931                         //append files
932                         m->appendFiles((outputFName+groups[i]), outputFName); m->mothurRemove((outputFName+groups[i]));
933                         m->appendFiles((accnos+groups[i]), accnos); m->mothurRemove((accnos+groups[i]));
934                         if (chimealns) { m->appendFiles((alns+groups[i]), alns); m->mothurRemove((alns+groups[i])); }
935                         
936                         m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences from group " + groups[i] + ".");    m->mothurOutEndLine();                                  
937                 }       
938                 
939                 return totalSeqs;
940                 
941         }
942         catch(exception& e) {
943                 m->errorOut(e, "ChimeraUchimeCommand", "driverGroups");
944                 exit(1);
945         }
946 }       
947 //**********************************************************************************************************************
948
949 int ChimeraUchimeCommand::driver(string outputFName, string filename, string accnos, string alns, int& numChimeras){
950         try {
951                 
952                 vector<char*> cPara;
953                 
954                 char* tempUchime = new char[8];  
955                 strcpy(tempUchime, "./uchime "); 
956                 cPara.push_back(tempUchime);
957                 
958                 char* tempIn = new char[8]; 
959                 *tempIn = '\0'; strncat(tempIn, "--input", 7);
960                 //strcpy(tempIn, "--input"); 
961                 cPara.push_back(tempIn);
962                 char* temp = new char[filename.length()+1];
963                 *temp = '\0'; strncat(temp, filename.c_str(), filename.length());
964                 //strcpy(temp, filename.c_str());
965                 cPara.push_back(temp);
966                 
967                 //are you using a reference file
968                 if (templatefile != "self") {
969                         //add reference file
970                         char* tempRef = new char[5]; 
971                         //strcpy(tempRef, "--db"); 
972                         *tempRef = '\0'; strncat(tempRef, "--db", 4);
973                         cPara.push_back(tempRef);  
974                         char* tempR = new char[templatefile.length()+1];
975                         //strcpy(tempR, templatefile.c_str());
976                         *tempR = '\0'; strncat(tempR, templatefile.c_str(), templatefile.length());
977                         cPara.push_back(tempR);
978                 }
979                 
980                 char* tempO = new char[12]; 
981                 *tempO = '\0'; strncat(tempO, "--uchimeout", 11);
982                 //strcpy(tempO, "--uchimeout"); 
983                 cPara.push_back(tempO);
984                 char* tempout = new char[outputFName.length()+1];
985                 //strcpy(tempout, outputFName.c_str());
986                 *tempout = '\0'; strncat(tempout, outputFName.c_str(), outputFName.length());
987                 cPara.push_back(tempout);
988                 
989                 if (chimealns) {
990                         char* tempA = new char[13]; 
991                         *tempA = '\0'; strncat(tempA, "--uchimealns", 12);
992                         //strcpy(tempA, "--uchimealns"); 
993                         cPara.push_back(tempA);
994                         char* tempa = new char[alns.length()+1];
995                         //strcpy(tempa, alns.c_str());
996                         *tempa = '\0'; strncat(tempa, alns.c_str(), alns.length());
997                         cPara.push_back(tempa);
998                 }
999                 
1000                 if (useAbskew) {
1001                         char* tempskew = new char[9];
1002                         *tempskew = '\0'; strncat(tempskew, "--abskew", 8);
1003                         //strcpy(tempskew, "--abskew"); 
1004                         cPara.push_back(tempskew);
1005                         char* tempSkew = new char[abskew.length()+1];
1006                         //strcpy(tempSkew, abskew.c_str());
1007                         *tempSkew = '\0'; strncat(tempSkew, abskew.c_str(), abskew.length());
1008                         cPara.push_back(tempSkew);
1009                 }
1010                 
1011                 if (useMinH) {
1012                         char* tempminh = new char[7]; 
1013                         *tempminh = '\0'; strncat(tempminh, "--minh", 6);
1014                         //strcpy(tempminh, "--minh"); 
1015                         cPara.push_back(tempminh);
1016                         char* tempMinH = new char[minh.length()+1];
1017                         *tempMinH = '\0'; strncat(tempMinH, minh.c_str(), minh.length());
1018                         //strcpy(tempMinH, minh.c_str());
1019                         cPara.push_back(tempMinH);
1020                 }
1021                 
1022                 if (useMindiv) {
1023                         char* tempmindiv = new char[9]; 
1024                         *tempmindiv = '\0'; strncat(tempmindiv, "--mindiv", 8);
1025                         //strcpy(tempmindiv, "--mindiv"); 
1026                         cPara.push_back(tempmindiv);
1027                         char* tempMindiv = new char[mindiv.length()+1];
1028                         *tempMindiv = '\0'; strncat(tempMindiv, mindiv.c_str(), mindiv.length());
1029                         //strcpy(tempMindiv, mindiv.c_str());
1030                         cPara.push_back(tempMindiv);
1031                 }
1032                 
1033                 if (useXn) {
1034                         char* tempxn = new char[5]; 
1035                         //strcpy(tempxn, "--xn"); 
1036                         *tempxn = '\0'; strncat(tempxn, "--xn", 4);
1037                         cPara.push_back(tempxn);
1038                         char* tempXn = new char[xn.length()+1];
1039                         //strcpy(tempXn, xn.c_str());
1040                         *tempXn = '\0'; strncat(tempXn, xn.c_str(), xn.length());
1041                         cPara.push_back(tempXn);
1042                 }
1043                 
1044                 if (useDn) {
1045                         char* tempdn = new char[5]; 
1046                         //strcpy(tempdn, "--dn"); 
1047                         *tempdn = '\0'; strncat(tempdn, "--dn", 4);
1048                         cPara.push_back(tempdn);
1049                         char* tempDn = new char[dn.length()+1];
1050                         *tempDn = '\0'; strncat(tempDn, dn.c_str(), dn.length());
1051                         //strcpy(tempDn, dn.c_str());
1052                         cPara.push_back(tempDn);
1053                 }
1054                 
1055                 if (useXa) {
1056                         char* tempxa = new char[5]; 
1057                         //strcpy(tempxa, "--xa"); 
1058                         *tempxa = '\0'; strncat(tempxa, "--xa", 4);
1059                         cPara.push_back(tempxa);
1060                         char* tempXa = new char[xa.length()+1];
1061                         *tempXa = '\0'; strncat(tempXa, xa.c_str(), xa.length());
1062                         //strcpy(tempXa, xa.c_str());
1063                         cPara.push_back(tempXa);
1064                 }
1065                 
1066                 if (useChunks) {
1067                         char* tempchunks = new char[9]; 
1068                         //strcpy(tempchunks, "--chunks"); 
1069                         *tempchunks = '\0'; strncat(tempchunks, "--chunks", 8);
1070                         cPara.push_back(tempchunks);
1071                         char* tempChunks = new char[chunks.length()+1];
1072                         *tempChunks = '\0'; strncat(tempChunks, chunks.c_str(), chunks.length());
1073                         //strcpy(tempChunks, chunks.c_str());
1074                         cPara.push_back(tempChunks);
1075                 }
1076                 
1077                 if (useMinchunk) {
1078                         char* tempminchunk = new char[11]; 
1079                         //strcpy(tempminchunk, "--minchunk"); 
1080                         *tempminchunk = '\0'; strncat(tempminchunk, "--minchunk", 10);
1081                         cPara.push_back(tempminchunk);
1082                         char* tempMinchunk = new char[minchunk.length()+1];
1083                         *tempMinchunk = '\0'; strncat(tempMinchunk, minchunk.c_str(), minchunk.length());
1084                         //strcpy(tempMinchunk, minchunk.c_str());
1085                         cPara.push_back(tempMinchunk);
1086                 }
1087                 
1088                 if (useIdsmoothwindow) {
1089                         char* tempidsmoothwindow = new char[17]; 
1090                         *tempidsmoothwindow = '\0'; strncat(tempidsmoothwindow, "--idsmoothwindow", 16);
1091                         //strcpy(tempidsmoothwindow, "--idsmoothwindow"); 
1092                         cPara.push_back(tempidsmoothwindow);
1093                         char* tempIdsmoothwindow = new char[idsmoothwindow.length()+1];
1094                         *tempIdsmoothwindow = '\0'; strncat(tempIdsmoothwindow, idsmoothwindow.c_str(), idsmoothwindow.length());
1095                         //strcpy(tempIdsmoothwindow, idsmoothwindow.c_str());
1096                         cPara.push_back(tempIdsmoothwindow);
1097                 }
1098                 
1099                 /*if (useMinsmoothid) {
1100                         char* tempminsmoothid = new char[14]; 
1101                         //strcpy(tempminsmoothid, "--minsmoothid"); 
1102                         *tempminsmoothid = '\0'; strncat(tempminsmoothid, "--minsmoothid", 13);
1103                         cPara.push_back(tempminsmoothid);
1104                         char* tempMinsmoothid = new char[minsmoothid.length()+1];
1105                         *tempMinsmoothid = '\0'; strncat(tempMinsmoothid, minsmoothid.c_str(), minsmoothid.length());
1106                         //strcpy(tempMinsmoothid, minsmoothid.c_str());
1107                         cPara.push_back(tempMinsmoothid);
1108                 }*/
1109                 
1110                 if (useMaxp) {
1111                         char* tempmaxp = new char[7]; 
1112                         //strcpy(tempmaxp, "--maxp"); 
1113                         *tempmaxp = '\0'; strncat(tempmaxp, "--maxp", 6);
1114                         cPara.push_back(tempmaxp);
1115                         char* tempMaxp = new char[maxp.length()+1];
1116                         *tempMaxp = '\0'; strncat(tempMaxp, maxp.c_str(), maxp.length());
1117                         //strcpy(tempMaxp, maxp.c_str());
1118                         cPara.push_back(tempMaxp);
1119                 }
1120                 
1121                 if (!skipgaps) {
1122                         char* tempskipgaps = new char[13]; 
1123                         //strcpy(tempskipgaps, "--[no]skipgaps");
1124                         *tempskipgaps = '\0'; strncat(tempskipgaps, "--noskipgaps", 12);
1125                         cPara.push_back(tempskipgaps);
1126                 }
1127                 
1128                 if (!skipgaps2) {
1129                         char* tempskipgaps2 = new char[14]; 
1130                         //strcpy(tempskipgaps2, "--[no]skipgaps2"); 
1131                         *tempskipgaps2 = '\0'; strncat(tempskipgaps2, "--noskipgaps2", 13);
1132                         cPara.push_back(tempskipgaps2);
1133                 }
1134                 
1135                 if (useMinlen) {
1136                         char* tempminlen = new char[9]; 
1137                         *tempminlen = '\0'; strncat(tempminlen, "--minlen", 8);
1138                         //strcpy(tempminlen, "--minlen"); 
1139                         cPara.push_back(tempminlen);
1140                         char* tempMinlen = new char[minlen.length()+1];
1141                         //strcpy(tempMinlen, minlen.c_str());
1142                         *tempMinlen = '\0'; strncat(tempMinlen, minlen.c_str(), minlen.length());
1143                         cPara.push_back(tempMinlen);
1144                 }
1145                 
1146                 if (useMaxlen) {
1147                         char* tempmaxlen = new char[9]; 
1148                         //strcpy(tempmaxlen, "--maxlen"); 
1149                         *tempmaxlen = '\0'; strncat(tempmaxlen, "--maxlen", 8);
1150                         cPara.push_back(tempmaxlen);
1151                         char* tempMaxlen = new char[maxlen.length()+1];
1152                         *tempMaxlen = '\0'; strncat(tempMaxlen, maxlen.c_str(), maxlen.length());
1153                         //strcpy(tempMaxlen, maxlen.c_str());
1154                         cPara.push_back(tempMaxlen);
1155                 }
1156                 
1157                 if (ucl) {
1158                         char* tempucl = new char[5]; 
1159                         strcpy(tempucl, "--ucl"); 
1160                         cPara.push_back(tempucl);
1161                 }
1162                 
1163                 if (useQueryfract) {
1164                         char* tempqueryfract = new char[13]; 
1165                         *tempqueryfract = '\0'; strncat(tempqueryfract, "--queryfract", 12);
1166                         //strcpy(tempqueryfract, "--queryfract"); 
1167                         cPara.push_back(tempqueryfract);
1168                         char* tempQueryfract = new char[queryfract.length()+1];
1169                         *tempQueryfract = '\0'; strncat(tempQueryfract, queryfract.c_str(), queryfract.length());
1170                         //strcpy(tempQueryfract, queryfract.c_str());
1171                         cPara.push_back(tempQueryfract);
1172                 }
1173                 
1174                 
1175                 char** uchimeParameters;
1176                 uchimeParameters = new char*[cPara.size()];
1177                 for (int i = 0; i < cPara.size(); i++) {  uchimeParameters[i] = cPara[i];  } 
1178                 int numArgs = cPara.size();
1179                 
1180                 uchime_main(numArgs, uchimeParameters); 
1181                 
1182                 //free memory
1183                 for(int i = 0; i < cPara.size(); i++)  {  delete[] cPara[i];  }
1184                 delete[] uchimeParameters; 
1185                 
1186                 if (m->control_pressed) { return 0; }
1187                 
1188                 //create accnos file from uchime results
1189                 ifstream in; 
1190                 m->openInputFile(outputFName, in);
1191                 
1192                 ofstream out;
1193                 m->openOutputFile(accnos, out);
1194                 
1195                 int num = 0;
1196                 numChimeras = 0;
1197                 while(!in.eof()) {
1198                         
1199                         if (m->control_pressed) { break; }
1200                         
1201                         string name = "";
1202                         string chimeraFlag = "";
1203                         in >> chimeraFlag >> name;
1204                         
1205                         //fix name if needed
1206                         if (templatefile == "self") { 
1207                                 name = name.substr(0, name.length()-1); //rip off last /
1208                                 name = name.substr(0, name.find_last_of('/'));
1209                         }
1210                         
1211                         for (int i = 0; i < 15; i++) {  in >> chimeraFlag; }
1212                         m->gobble(in);
1213                         
1214                         if (chimeraFlag == "Y") {  out << name << endl; numChimeras++; }
1215                         num++;
1216                 }
1217                 in.close();
1218                 out.close();
1219                 
1220                 return num;
1221         }
1222         catch(exception& e) {
1223                 m->errorOut(e, "ChimeraUchimeCommand", "driver");
1224                 exit(1);
1225         }
1226 }
1227 /**************************************************************************************************/
1228
1229 int ChimeraUchimeCommand::createProcesses(string outputFileName, string filename, string accnos, string alns, int& numChimeras) {
1230         try {
1231                 
1232                 processIDS.clear();
1233                 int process = 1;
1234                 int num = 0;
1235 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)           
1236                 //break up file into multiple files
1237                 vector<string> files;
1238                 m->divideFile(filename, processors, files);
1239                 
1240                 if (m->control_pressed) {  return 0;  }
1241                                 
1242                 //loop through and create all the processes you want
1243                 while (process != processors) {
1244                         int pid = fork();
1245                         
1246                         if (pid > 0) {
1247                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
1248                                 process++;
1249                         }else if (pid == 0){
1250                                 num = driver(outputFileName + toString(getpid()) + ".temp", files[process], accnos + toString(getpid()) + ".temp", alns + toString(getpid()) + ".temp", numChimeras);
1251                                 
1252                                 //pass numSeqs to parent
1253                                 ofstream out;
1254                                 string tempFile = outputFileName + toString(getpid()) + ".num.temp";
1255                                 m->openOutputFile(tempFile, out);
1256                                 out << num << endl;
1257                                 out << numChimeras << endl;
1258                                 out.close();
1259                                 
1260                                 exit(0);
1261                         }else { 
1262                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
1263                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
1264                                 exit(0);
1265                         }
1266                 }
1267                 
1268                 //do my part
1269                 num = driver(outputFileName, files[0], accnos, alns, numChimeras);
1270                 
1271                 //force parent to wait until all the processes are done
1272                 for (int i=0;i<processIDS.size();i++) { 
1273                         int temp = processIDS[i];
1274                         wait(&temp);
1275                 }
1276                 
1277                 for (int i = 0; i < processIDS.size(); i++) {
1278                         ifstream in;
1279                         string tempFile =  outputFileName + toString(processIDS[i]) + ".num.temp";
1280                         m->openInputFile(tempFile, in);
1281                         if (!in.eof()) { 
1282                                 int tempNum = 0; 
1283                                 in >> tempNum; m->gobble(in);
1284                                 num += tempNum; 
1285                                 in >> tempNum;
1286                                 numChimeras += tempNum;
1287                         }
1288                         in.close(); m->mothurRemove(tempFile);
1289                 }
1290                 
1291                 
1292                 //append output files
1293                 for(int i=0;i<processIDS[i];i++){
1294                         m->appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
1295                         m->mothurRemove((outputFileName + toString(processIDS[i]) + ".temp"));
1296                         
1297                         m->appendFiles((accnos + toString(processIDS[i]) + ".temp"), accnos);
1298                         m->mothurRemove((accnos + toString(processIDS[i]) + ".temp"));
1299                         
1300                         if (chimealns) {
1301                                 m->appendFiles((alns + toString(processIDS[i]) + ".temp"), alns);
1302                                 m->mothurRemove((alns + toString(processIDS[i]) + ".temp"));
1303                         }
1304                 }
1305                 
1306                 //get rid of the file pieces.
1307                 for (int i = 0; i < files.size(); i++) { m->mothurRemove(files[i]); }
1308 #endif          
1309                 return num;     
1310         }
1311         catch(exception& e) {
1312                 m->errorOut(e, "ChimeraUchimeCommand", "createProcesses");
1313                 exit(1);
1314         }
1315 }
1316 /**************************************************************************************************/
1317
1318 int ChimeraUchimeCommand::createProcessesGroups(SequenceParser& parser, string outputFName, string filename, string accnos, string alns, vector<string> groups) {
1319         try {
1320                 
1321                 processIDS.clear();
1322                 int process = 1;
1323                 int num = 0;
1324                 
1325                 //sanity check
1326                 if (groups.size() < processors) { processors = groups.size(); }
1327                 
1328                 //divide the groups between the processors
1329                 vector<linePair> lines;
1330                 int numGroupsPerProcessor = groups.size() / processors;
1331                 for (int i = 0; i < processors; i++) {
1332                         int startIndex =  i * numGroupsPerProcessor;
1333                         int endIndex = (i+1) * numGroupsPerProcessor;
1334                         if(i == (processors - 1)){      endIndex = groups.size();       }
1335                         lines.push_back(linePair(startIndex, endIndex));
1336                 }
1337                 
1338 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)           
1339                                 
1340                 //loop through and create all the processes you want
1341                 while (process != processors) {
1342                         int pid = fork();
1343                         
1344                         if (pid > 0) {
1345                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
1346                                 process++;
1347                         }else if (pid == 0){
1348                                 num = driverGroups(parser, outputFName + toString(getpid()) + ".temp", filename + toString(getpid()) + ".temp", accnos + toString(getpid()) + ".temp", alns + toString(getpid()) + ".temp", lines[process].start, lines[process].end, groups);
1349                                 
1350                                 //pass numSeqs to parent
1351                                 ofstream out;
1352                                 string tempFile = outputFName + toString(getpid()) + ".num.temp";
1353                                 m->openOutputFile(tempFile, out);
1354                                 out << num << endl;
1355                                 out.close();
1356                                 
1357                                 exit(0);
1358                         }else { 
1359                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
1360                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
1361                                 exit(0);
1362                         }
1363                 }
1364                 
1365                 //do my part
1366                 num = driverGroups(parser, outputFName, filename, accnos, alns, lines[0].start, lines[0].end, groups);
1367                 
1368                 //force parent to wait until all the processes are done
1369                 for (int i=0;i<processIDS.size();i++) { 
1370                         int temp = processIDS[i];
1371                         wait(&temp);
1372                 }
1373 #endif          
1374                 
1375                 for (int i = 0; i < processIDS.size(); i++) {
1376                         ifstream in;
1377                         string tempFile =  outputFName + toString(processIDS[i]) + ".num.temp";
1378                         m->openInputFile(tempFile, in);
1379                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
1380                         in.close(); m->mothurRemove(tempFile);
1381                 }
1382                 
1383                 
1384                 //append output files
1385                 for(int i=0;i<processIDS[i];i++){
1386                         m->appendFiles((outputFName + toString(processIDS[i]) + ".temp"), outputFName);
1387                         m->mothurRemove((outputFName + toString(processIDS[i]) + ".temp"));
1388                         
1389                         m->appendFiles((accnos + toString(processIDS[i]) + ".temp"), accnos);
1390                         m->mothurRemove((accnos + toString(processIDS[i]) + ".temp"));
1391                         
1392                         if (chimealns) {
1393                                 m->appendFiles((alns + toString(processIDS[i]) + ".temp"), alns);
1394                                 m->mothurRemove((alns + toString(processIDS[i]) + ".temp"));
1395                         }
1396                 }
1397                 
1398                 return num;     
1399                 
1400         }
1401         catch(exception& e) {
1402                 m->errorOut(e, "ChimeraUchimeCommand", "createProcessesGroups");
1403                 exit(1);
1404         }
1405 }
1406 /**************************************************************************************************/
1407