]> git.donarmstrong.com Git - mothur.git/blob - chimeraslayercommand.cpp
added getCommandInfoCommand for gui
[mothur.git] / chimeraslayercommand.cpp
1 /*
2  *  chimeraslayercommand.cpp
3  *  Mothur
4  *
5  *  Created by westcott on 3/31/10.
6  *  Copyright 2010 Schloss Lab. All rights reserved.
7  *
8  */
9
10 #include "chimeraslayercommand.h"
11 #include "chimeraslayer.h"
12 #include "deconvolutecommand.h"
13
14 //**********************************************************************************************************************
15 vector<string> ChimeraSlayerCommand::setParameters(){   
16         try {
17                 CommandParameter ptemplate("reference", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(ptemplate);
18                 CommandParameter pfasta("fasta", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(pfasta);
19                 CommandParameter pname("name", "InputTypes", "", "", "none", "none", "none",false,false); parameters.push_back(pname);
20                 CommandParameter pwindow("window", "Number", "", "50", "", "", "",false,false); parameters.push_back(pwindow);
21                 CommandParameter pksize("ksize", "Number", "", "7", "", "", "",false,false); parameters.push_back(pksize);
22                 CommandParameter pmatch("match", "Number", "", "5.0", "", "", "",false,false); parameters.push_back(pmatch);
23                 CommandParameter pmismatch("mismatch", "Number", "", "-4.0", "", "", "",false,false); parameters.push_back(pmismatch);
24                 CommandParameter pminsim("minsim", "Number", "", "90", "", "", "",false,false); parameters.push_back(pminsim);
25                 CommandParameter pmincov("mincov", "Number", "", "70", "", "", "",false,false); parameters.push_back(pmincov);
26                 CommandParameter pminsnp("minsnp", "Number", "", "100", "", "", "",false,false); parameters.push_back(pminsnp);
27                 CommandParameter pminbs("minbs", "Number", "", "90", "", "", "",false,false); parameters.push_back(pminbs);
28                 CommandParameter psearch("search", "Multiple", "kmer-blast-distance", "distance", "", "", "",false,false); parameters.push_back(psearch);
29                 CommandParameter pinclude("include", "Multiple", "greater-greaterequal-all", "greater", "", "", "",false,false); parameters.push_back(pinclude);
30                 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
31                 CommandParameter prealign("realign", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(prealign);
32                 CommandParameter ptrim("trim", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(ptrim);
33                 CommandParameter psplit("split", "Boolean", "", "F", "", "", "",false,false); parameters.push_back(psplit);
34                 CommandParameter pnumwanted("numwanted", "Number", "", "15", "", "", "",false,false); parameters.push_back(pnumwanted);
35                 CommandParameter piters("iters", "Number", "", "100", "", "", "",false,false); parameters.push_back(piters);
36                 CommandParameter pdivergence("divergence", "Number", "", "1.007", "", "", "",false,false); parameters.push_back(pdivergence);
37                 CommandParameter pparents("parents", "Number", "", "3", "", "", "",false,false); parameters.push_back(pparents);
38                 CommandParameter pincrement("increment", "Number", "", "5", "", "", "",false,false); parameters.push_back(pincrement);
39                 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
40                 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
41                 
42                 vector<string> myArray;
43                 for (int i = 0; i < parameters.size(); i++) {   myArray.push_back(parameters[i].name);          }
44                 return myArray;
45         }
46         catch(exception& e) {
47                 m->errorOut(e, "ChimeraSlayerCommand", "setParameters");
48                 exit(1);
49         }
50 }
51 //**********************************************************************************************************************
52 string ChimeraSlayerCommand::getHelpString(){   
53         try {
54                 string helpString = "";
55                 helpString += "The chimera.slayer command reads a fastafile and referencefile and outputs potentially chimeric sequences.\n";
56                 helpString += "This command was modeled after the chimeraSlayer written by the Broad Institute.\n";
57                 helpString += "The chimera.slayer command parameters are fasta, name, template, processors, trim, ksize, window, match, mismatch, divergence. minsim, mincov, minbs, minsnp, parents, search, iters, increment and numwanted.\n"; //realign,
58                 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";
59                 helpString += "The name parameter allows you to provide a name file, if you are using template=self. \n";
60                 helpString += "You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amazon.fasta \n";
61                 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";
62                 helpString += "The include parameter is used when template=self and allows you to choose which sequences will make up the \"template\". Options are greater, greaterequal and all, default=greater, meaning sequences with greater abundance than the query sequence. \n";
63                 helpString += "The processors parameter allows you to specify how many processors you would like to use.  The default is 1. \n";
64 #ifdef USE_MPI
65                 helpString += "When using MPI, the processors parameter is set to the number of MPI processes running. \n";
66 #endif
67                 helpString += "The trim parameter allows you to output a new fasta file containing your sequences with the chimeric ones trimmed to include only their longest piece, default=F. \n";
68                 helpString += "The split parameter allows you to check both pieces of non-chimeric sequence for chimeras, thus looking for trimeras and quadmeras. default=F. \n";
69                 helpString += "The window parameter allows you to specify the window size for searching for chimeras, default=50. \n";
70                 helpString += "The increment parameter allows you to specify how far you move each window while finding chimeric sequences, default=5.\n";
71                 helpString += "The numwanted parameter allows you to specify how many sequences you would each query sequence compared with, default=15.\n";
72                 helpString += "The ksize parameter allows you to input kmersize, default is 7, used if search is kmer. \n";
73                 helpString += "The match parameter allows you to reward matched bases in blast search, default is 5. \n";
74                 helpString += "The parents parameter allows you to select the number of potential parents to investigate from the numwanted best matches after rating them, default is 3. \n";
75                 helpString += "The mismatch parameter allows you to penalize mismatched bases in blast search, default is -4. \n";
76                 helpString += "The divergence parameter allows you to set a cutoff for chimera determination, default is 1.007. \n";
77                 helpString += "The iters parameter allows you to specify the number of bootstrap iters to do with the chimeraslayer method, default=100.\n";
78                 helpString += "The minsim parameter allows you to specify a minimum similarity with the parent fragments, default=90. \n";
79                 helpString += "The mincov parameter allows you to specify minimum coverage by closest matches found in template. Default is 70, meaning 70%. \n";
80                 helpString += "The minbs parameter allows you to specify minimum bootstrap support for calling a sequence chimeric. Default is 90, meaning 90%. \n";
81                 helpString += "The minsnp parameter allows you to specify percent of SNPs to sample on each side of breakpoint for computing bootstrap support (default: 100) \n";
82                 helpString += "The search parameter allows you to specify search method for finding the closest parent. Choices are distance, blast, and kmer, default distance. \n";
83                 helpString += "The realign parameter allows you to realign the query to the potential parents. Choices are true or false, default false.  \n";
84                 helpString += "The chimera.slayer command should be in the following format: \n";
85                 helpString += "chimera.slayer(fasta=yourFastaFile, reference=yourTemplate, search=yourSearch) \n";
86                 helpString += "Example: chimera.slayer(fasta=AD.align, reference=core_set_aligned.imputed.fasta, search=kmer) \n";
87                 helpString += "Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n";       
88                 return helpString;
89         }
90         catch(exception& e) {
91                 m->errorOut(e, "ChimeraSlayerCommand", "getHelpString");
92                 exit(1);
93         }
94 }
95 //**********************************************************************************************************************
96 ChimeraSlayerCommand::ChimeraSlayerCommand(){   
97         try {
98                 abort = true; calledHelp = true;
99                 setParameters();
100                 vector<string> tempOutNames;
101                 outputTypes["chimera"] = tempOutNames;
102                 outputTypes["accnos"] = tempOutNames;
103                 outputTypes["fasta"] = tempOutNames;
104         }
105         catch(exception& e) {
106                 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
107                 exit(1);
108         }
109 }
110 //***************************************************************************************************************
111 ChimeraSlayerCommand::ChimeraSlayerCommand(string option)  {
112         try {
113                 abort = false; calledHelp = false;   
114                 
115                 //allow user to run help
116                 if(option == "help") { help(); abort = true; calledHelp = true; }
117                 
118                 else {
119                         vector<string> myArray = setParameters();
120                         
121                         OptionParser parser(option);
122                         map<string,string> parameters = parser.getParameters();
123                         
124                         ValidParameters validParameter("chimera.slayer");
125                         map<string,string>::iterator it;
126                         
127                         //check to make sure all parameters are valid for command
128                         for (it = parameters.begin(); it != parameters.end(); it++) { 
129                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
130                         }
131                         
132                         vector<string> tempOutNames;
133                         outputTypes["chimera"] = tempOutNames;
134                         outputTypes["accnos"] = tempOutNames;
135                         outputTypes["fasta"] = tempOutNames;
136                 
137                         //if the user changes the input directory command factory will send this info to us in the output parameter 
138                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
139                         if (inputDir == "not found"){   inputDir = "";          }
140                                                 
141                         //check for required parameters
142                         fastafile = validParameter.validFile(parameters, "fasta", false);
143                         if (fastafile == "not found") {                                 
144                                 //if there is a current fasta file, use it
145                                 string filename = m->getFastaFile(); 
146                                 if (filename != "") { fastaFileNames.push_back(filename); m->mothurOut("Using " + filename + " as input file for the fasta parameter."); m->mothurOutEndLine(); }
147                                 else {  m->mothurOut("You have no current fastafile and the fasta parameter is required."); m->mothurOutEndLine(); abort = true; }
148                         }else { 
149                                 m->splitAtDash(fastafile, fastaFileNames);
150                                 
151                                 //go through files and make sure they are good, if not, then disregard them
152                                 for (int i = 0; i < fastaFileNames.size(); i++) {
153                                         if (inputDir != "") {
154                                                 string path = m->hasPath(fastaFileNames[i]);
155                                                 //if the user has not given a path then, add inputdir. else leave path alone.
156                                                 if (path == "") {       fastaFileNames[i] = inputDir + fastaFileNames[i];               }
157                                         }
158         
159                                         int ableToOpen;
160                                         ifstream in;
161                                         
162                                         ableToOpen = m->openInputFile(fastaFileNames[i], in, "noerror");
163                                 
164                                         //if you can't open it, try default location
165                                         if (ableToOpen == 1) {
166                                                 if (m->getDefaultPath() != "") { //default path is set
167                                                         string tryPath = m->getDefaultPath() + m->getSimpleName(fastaFileNames[i]);
168                                                         m->mothurOut("Unable to open " + fastaFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
169                                                         ifstream in2;
170                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
171                                                         in2.close();
172                                                         fastaFileNames[i] = tryPath;
173                                                 }
174                                         }
175                                         
176                                         if (ableToOpen == 1) {
177                                                 if (m->getOutputDir() != "") { //default path is set
178                                                         string tryPath = m->getOutputDir() + m->getSimpleName(fastaFileNames[i]);
179                                                         m->mothurOut("Unable to open " + fastaFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
180                                                         ifstream in2;
181                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
182                                                         in2.close();
183                                                         fastaFileNames[i] = tryPath;
184                                                 }
185                                         }
186                                         
187                                         in.close();
188                                         
189                                         if (ableToOpen == 1) { 
190                                                 m->mothurOut("Unable to open " + fastaFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
191                                                 //erase from file list
192                                                 fastaFileNames.erase(fastaFileNames.begin()+i);
193                                                 i--;
194                                         }
195                                 }
196                                 
197                                 //make sure there is at least one valid file left
198                                 if (fastaFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid files."); m->mothurOutEndLine(); abort = true; }
199                         }
200                         
201                         
202                         //check for required parameters
203                         bool hasName = true;
204                         namefile = validParameter.validFile(parameters, "name", false);
205                         if (namefile == "not found") { namefile = "";  hasName = false; }
206                         else { 
207                                 m->splitAtDash(namefile, nameFileNames);
208                                 
209                                 //go through files and make sure they are good, if not, then disregard them
210                                 for (int i = 0; i < nameFileNames.size(); i++) {
211                                         if (inputDir != "") {
212                                                 string path = m->hasPath(nameFileNames[i]);
213                                                 //if the user has not given a path then, add inputdir. else leave path alone.
214                                                 if (path == "") {       nameFileNames[i] = inputDir + nameFileNames[i];         }
215                                         }
216                                         
217                                         int ableToOpen;
218                                         ifstream in;
219                                         
220                                         ableToOpen = m->openInputFile(nameFileNames[i], in, "noerror");
221                                         
222                                         //if you can't open it, try default location
223                                         if (ableToOpen == 1) {
224                                                 if (m->getDefaultPath() != "") { //default path is set
225                                                         string tryPath = m->getDefaultPath() + m->getSimpleName(nameFileNames[i]);
226                                                         m->mothurOut("Unable to open " + nameFileNames[i] + ". Trying default " + tryPath); m->mothurOutEndLine();
227                                                         ifstream in2;
228                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
229                                                         in2.close();
230                                                         nameFileNames[i] = tryPath;
231                                                 }
232                                         }
233                                         
234                                         if (ableToOpen == 1) {
235                                                 if (m->getOutputDir() != "") { //default path is set
236                                                         string tryPath = m->getOutputDir() + m->getSimpleName(nameFileNames[i]);
237                                                         m->mothurOut("Unable to open " + nameFileNames[i] + ". Trying output directory " + tryPath); m->mothurOutEndLine();
238                                                         ifstream in2;
239                                                         ableToOpen = m->openInputFile(tryPath, in2, "noerror");
240                                                         in2.close();
241                                                         nameFileNames[i] = tryPath;
242                                                 }
243                                         }
244                                         
245                                         in.close();
246                                         
247                                         if (ableToOpen == 1) { 
248                                                 m->mothurOut("Unable to open " + nameFileNames[i] + ". It will be disregarded."); m->mothurOutEndLine(); 
249                                                 //erase from file list
250                                                 nameFileNames.erase(nameFileNames.begin()+i);
251                                                 i--;
252                                         }
253                                 }
254                                 
255                                 //make sure there is at least one valid file left
256                                 if (nameFileNames.size() == 0) { m->mothurOut("[ERROR]: no valid name files."); m->mothurOutEndLine(); abort = true; }
257                         }
258                         
259                         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; }
260                         
261                         //if the user changes the output directory command factory will send this info to us in the output parameter 
262                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  outputDir = ""; }
263                         
264                         
265                         string path;
266                         it = parameters.find("reference");
267                         //user has given a template file
268                         if(it != parameters.end()){ 
269                                 if (it->second == "self") { templatefile = "self"; }
270                                 else {
271                                         path = m->hasPath(it->second);
272                                         //if the user has not given a path then, add inputdir. else leave path alone.
273                                         if (path == "") {       parameters["reference"] = inputDir + it->second;                }
274                                         
275                                         templatefile = validParameter.validFile(parameters, "reference", true);
276                                         if (templatefile == "not open") { abort = true; }
277                                         else if (templatefile == "not found") { templatefile = "";  m->mothurOut("reference is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true;  }  
278                                 }
279                         }
280                         
281                         string temp = validParameter.validFile(parameters, "processors", false);        if (temp == "not found"){       temp = m->getProcessors();      }
282                         m->setProcessors(temp);
283                         convert(temp, processors);
284                         
285                         includeAbunds = validParameter.validFile(parameters, "include", false);         if (includeAbunds == "not found") { includeAbunds = "greater"; }
286                         if ((includeAbunds != "greater") && (includeAbunds != "greaterequal") && (includeAbunds != "all")) { includeAbunds = "greater"; m->mothurOut("Invalid include setting. options are greater, greaterequal or all. using greater."); m->mothurOutEndLine(); }
287                         
288                         temp = validParameter.validFile(parameters, "ksize", false);                    if (temp == "not found") { temp = "7"; }
289                         convert(temp, ksize);
290                                                 
291                         temp = validParameter.validFile(parameters, "window", false);                   if (temp == "not found") { temp = "50"; }                       
292                         convert(temp, window);
293                         
294                         temp = validParameter.validFile(parameters, "match", false);                    if (temp == "not found") { temp = "5"; }
295                         convert(temp, match);
296                         
297                         temp = validParameter.validFile(parameters, "mismatch", false);                 if (temp == "not found") { temp = "-4"; }
298                         convert(temp, mismatch);
299                         
300                         temp = validParameter.validFile(parameters, "divergence", false);               if (temp == "not found") { temp = "1.007"; }
301                         convert(temp, divR);
302                         
303                         temp = validParameter.validFile(parameters, "minsim", false);                   if (temp == "not found") { temp = "90"; }
304                         convert(temp, minSimilarity);
305                         
306                         temp = validParameter.validFile(parameters, "mincov", false);                   if (temp == "not found") { temp = "70"; }
307                         convert(temp, minCoverage);
308                         
309                         temp = validParameter.validFile(parameters, "minbs", false);                    if (temp == "not found") { temp = "90"; }
310                         convert(temp, minBS);
311                         
312                         temp = validParameter.validFile(parameters, "minsnp", false);                   if (temp == "not found") { temp = "100"; }
313                         convert(temp, minSNP);
314
315                         temp = validParameter.validFile(parameters, "parents", false);                  if (temp == "not found") { temp = "3"; }
316                         convert(temp, parents); 
317                         
318                         temp = validParameter.validFile(parameters, "realign", false);                  if (temp == "not found") { temp = "f"; }
319                         realign = m->isTrue(temp); 
320                         
321                         temp = validParameter.validFile(parameters, "trim", false);                             if (temp == "not found") { temp = "f"; }
322                         trim = m->isTrue(temp); 
323                         
324                         temp = validParameter.validFile(parameters, "split", false);                            if (temp == "not found") { temp = "f"; }
325                         trimera = m->isTrue(temp); 
326                         
327                         search = validParameter.validFile(parameters, "search", false);                 if (search == "not found") { search = "distance"; }
328                         
329                         temp = validParameter.validFile(parameters, "iters", false);                    if (temp == "not found") { temp = "100"; }              
330                         convert(temp, iters); 
331                          
332                         temp = validParameter.validFile(parameters, "increment", false);                if (temp == "not found") { temp = "5"; }
333                         convert(temp, increment);
334                         
335                         temp = validParameter.validFile(parameters, "numwanted", false);                if (temp == "not found") { temp = "15"; }               
336                         convert(temp, numwanted);
337
338                         if ((search != "distance") && (search != "blast") && (search != "kmer")) { m->mothurOut(search + " is not a valid search."); m->mothurOutEndLine(); abort = true;  }
339                 }
340         }
341         catch(exception& e) {
342                 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
343                 exit(1);
344         }
345 }
346 //***************************************************************************************************************
347
348 int ChimeraSlayerCommand::execute(){
349         try{
350                 
351                 if (abort == true) { if (calledHelp) { return 0; }  return 2;   }
352                 
353                 for (int s = 0; s < fastaFileNames.size(); s++) {
354                                 
355                         m->mothurOut("Checking sequences from " + fastaFileNames[s] + " ..." ); m->mothurOutEndLine();
356                 
357                         int start = time(NULL); 
358                         
359                         if (templatefile != "self") { //you want to run slayer with a refernce template
360                                 chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, trim, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);     
361                         }else {
362                                 if (nameFileNames.size() != 0) { //you provided a namefile and we don't need to create one
363                                         chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, trim, nameFileNames[s], search, includeAbunds, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);    
364                                 }else {
365                                         
366                                         m->mothurOutEndLine(); m->mothurOut("No namesfile given, running unique.seqs command to generate one."); m->mothurOutEndLine(); m->mothurOutEndLine();
367                                         
368                                         //use unique.seqs to create new name and fastafile
369                                         string inputString = "fasta=" + fastaFileNames[s];
370                                         m->mothurOut("/******************************************/"); m->mothurOutEndLine(); 
371                                         m->mothurOut("Running command: unique.seqs(" + inputString + ")"); m->mothurOutEndLine(); 
372                                                                  
373                                         Command* uniqueCommand = new DeconvoluteCommand(inputString);
374                                         uniqueCommand->execute();
375                                         
376                                         map<string, vector<string> > filenames = uniqueCommand->getOutputFiles();
377                                         
378                                         delete uniqueCommand;
379                                         
380                                         m->mothurOut("/******************************************/"); m->mothurOutEndLine(); 
381                                         
382                                         string nameFile = filenames["name"][0];
383                                         fastaFileNames[s] = filenames["fasta"][0];
384                         
385                                         chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, trim, nameFile, search, includeAbunds, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);    
386                                 }
387                         }
388                                 
389                         if (outputDir == "") { outputDir = m->hasPath(fastaFileNames[s]);  }//if user entered a file with a path then preserve it                               
390                         string outputFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s])) + "slayer.chimera";
391                         string accnosFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s]))  + "slayer.accnos";
392                         string trimFastaFileName = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s]))  + "slayer.fasta";
393                         
394                         if (m->control_pressed) { delete chimera; for (int j = 0; j < outputNames.size(); j++) {        remove(outputNames[j].c_str()); }  return 0;    }
395                         
396                         if (chimera->getUnaligned()) { 
397                                 m->mothurOut("Your template sequences are different lengths, please correct."); m->mothurOutEndLine(); 
398                                 delete chimera;
399                                 return 0; 
400                         }
401                         templateSeqsLength = chimera->getLength();
402                         
403                 #ifdef USE_MPI  
404                         int pid, numSeqsPerProcessor; 
405                                 int tag = 2001;
406                                 vector<unsigned long int> MPIPos;
407                                 
408                                 MPI_Status status; 
409                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
410                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); 
411
412                                 MPI_File inMPI;
413                                 MPI_File outMPI;
414                                 MPI_File outMPIAccnos;
415                                 MPI_File outMPIFasta;
416                                 
417                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
418                                 int inMode=MPI_MODE_RDONLY; 
419                                 
420                                 char outFilename[1024];
421                                 strcpy(outFilename, outputFileName.c_str());
422                                 
423                                 char outAccnosFilename[1024];
424                                 strcpy(outAccnosFilename, accnosFileName.c_str());
425                         
426                                 char outFastaFilename[1024];
427                                 strcpy(outFastaFilename, trimFastaFileName.c_str());
428                                 
429                                 char inFileName[1024];
430                                 strcpy(inFileName, fastaFileNames[s].c_str());
431
432                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
433                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
434                                 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
435                                 if (trim) { MPI_File_open(MPI_COMM_WORLD, outFastaFilename, outMode, MPI_INFO_NULL, &outMPIFasta); }
436
437                         if (m->control_pressed) { outputTypes.clear();  MPI_File_close(&inMPI);  MPI_File_close(&outMPI); if (trim) {  MPI_File_close(&outMPIFasta);  } MPI_File_close(&outMPIAccnos); for (int j = 0; j < outputNames.size(); j++) {   remove(outputNames[j].c_str()); }   delete chimera; return 0;  }
438                         
439                                 if (pid == 0) { //you are the root process 
440                                         m->mothurOutEndLine();
441                                         m->mothurOut("Only reporting sequence supported by " + toString(minBS) + "% of bootstrapped results.");
442                                         m->mothurOutEndLine();
443                 
444                                         string outTemp = "Name\tLeftParent\tRightParent\tDivQLAQRB\tPerIDQLAQRB\tBootStrapA\tDivQLBQRA\tPerIDQLBQRA\tBootStrapB\tFlag\tLeftWindow\tRightWindow\n";
445                                         
446                                         //print header
447                                         int length = outTemp.length();
448                                         char* buf2 = new char[length];
449                                         memcpy(buf2, outTemp.c_str(), length);
450
451                                         MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &status);
452                                         delete buf2;
453
454                                         MPIPos = m->setFilePosFasta(fastaFileNames[s], numSeqs); //fills MPIPos, returns numSeqs
455                                         
456                                         //send file positions to all processes
457                                         for(int i = 1; i < processors; i++) { 
458                                                 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
459                                                 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
460                                         }
461                                         
462                                         //figure out how many sequences you have to align
463                                         numSeqsPerProcessor = numSeqs / processors;
464                                         int startIndex =  pid * numSeqsPerProcessor;
465                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
466                                 
467                                         //do your part
468                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, outMPIFasta, MPIPos);
469                                         
470                                         if (m->control_pressed) { outputTypes.clear();  MPI_File_close(&inMPI);  MPI_File_close(&outMPI); if (trim) { MPI_File_close(&outMPIFasta); }  MPI_File_close(&outMPIAccnos);  for (int j = 0; j < outputNames.size(); j++) {   remove(outputNames[j].c_str()); }  remove(outputFileName.c_str());  remove(accnosFileName.c_str());  delete chimera; return 0;  }
471
472                                 }else{ //you are a child process
473                                         MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
474                                         MPIPos.resize(numSeqs+1);
475                                         MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
476                                         
477                                         //figure out how many sequences you have to align
478                                         numSeqsPerProcessor = numSeqs / processors;
479                                         int startIndex =  pid * numSeqsPerProcessor;
480                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
481                                         
482                                         //do your part
483                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, outMPIFasta, MPIPos);
484                                         
485                                         if (m->control_pressed) { outputTypes.clear();  MPI_File_close(&inMPI);  MPI_File_close(&outMPI); if (trim) { MPI_File_close(&outMPIFasta); }  MPI_File_close(&outMPIAccnos);  for (int j = 0; j < outputNames.size(); j++) {   remove(outputNames[j].c_str()); }  delete chimera; return 0;  }
486                                 }
487                                 
488                                 //close files 
489                                 MPI_File_close(&inMPI);
490                                 MPI_File_close(&outMPI);
491                                 MPI_File_close(&outMPIAccnos); 
492                                 if (trim) { MPI_File_close(&outMPIFasta); }
493                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
494                                 
495                 #else
496                         ofstream outHeader;
497                         string tempHeader = outputDir + m->getRootName(m->getSimpleName(fastaFileNames[s])) + "slayer.chimeras.tempHeader";
498                         m->openOutputFile(tempHeader, outHeader);
499                         
500                         chimera->printHeader(outHeader);
501                         outHeader.close();
502                         
503                         vector<unsigned long int> positions = m->divideFile(fastaFileNames[s], processors);
504                                 
505                         for (int i = 0; i < (positions.size()-1); i++) {
506                                 lines.push_back(new linePair(positions[i], positions[(i+1)]));
507                         }       
508
509                         //break up file
510                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
511                                 if(processors == 1){
512                                         numSeqs = driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName);
513                                         
514                                         if (m->control_pressed) { outputTypes.clear(); if (trim) { remove(trimFastaFileName.c_str()); } remove(outputFileName.c_str()); remove(tempHeader.c_str()); remove(accnosFileName.c_str()); for (int j = 0; j < outputNames.size(); j++) {      remove(outputNames[j].c_str()); } for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear(); delete chimera; return 0; }
515                                         
516                                 }else{
517                                         processIDS.resize(0);
518                                         
519                                         numSeqs = createProcesses(outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName); 
520                                 
521                                         rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
522                                         rename((accnosFileName + toString(processIDS[0]) + ".temp").c_str(), accnosFileName.c_str());
523                                         if (trim) {  rename((trimFastaFileName + toString(processIDS[0]) + ".temp").c_str(), trimFastaFileName.c_str()); }
524                                                 
525                                         //append output files
526                                         for(int i=1;i<processors;i++){
527                                                 m->appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
528                                                 remove((outputFileName + toString(processIDS[i]) + ".temp").c_str());
529                                         }
530                                         
531                                         //append output files
532                                         for(int i=1;i<processors;i++){
533                                                 m->appendFiles((accnosFileName + toString(processIDS[i]) + ".temp"), accnosFileName);
534                                                 remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str());
535                                         }
536                                         
537                                         if (trim) {
538                                                 for(int i=1;i<processors;i++){
539                                                         m->appendFiles((trimFastaFileName + toString(processIDS[i]) + ".temp"), trimFastaFileName);
540                                                         remove((trimFastaFileName + toString(processIDS[i]) + ".temp").c_str());
541                                                 }
542                                         }
543                                         
544                                         if (m->control_pressed) { outputTypes.clear(); if (trim) { remove(trimFastaFileName.c_str()); } remove(outputFileName.c_str()); remove(accnosFileName.c_str()); for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear(); delete chimera; return 0; }
545                                 }
546
547                         #else
548                                 numSeqs = driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName, trimFastaFileName);
549                                 
550                                 if (m->control_pressed) { outputTypes.clear(); if (trim) { remove(trimFastaFileName.c_str()); } remove(outputFileName.c_str()); remove(tempHeader.c_str()); remove(accnosFileName.c_str()); for (int j = 0; j < outputNames.size(); j++) {      remove(outputNames[j].c_str()); } for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear(); delete chimera; return 0; }
551                                 
552                         #endif
553                         
554                         m->appendFiles(outputFileName, tempHeader);
555                 
556                         remove(outputFileName.c_str());
557                         rename(tempHeader.c_str(), outputFileName.c_str());
558                         
559                 #endif
560                         delete chimera;
561                         
562                         
563                         for (int i = 0; i < lines.size(); i++) {  delete lines[i];  }  lines.clear();
564                         
565                         outputNames.push_back(outputFileName); outputTypes["chimera"].push_back(outputFileName);
566                         outputNames.push_back(accnosFileName); outputTypes["accnos"].push_back(accnosFileName);
567                         if (trim) {  outputNames.push_back(trimFastaFileName); outputTypes["fasta"].push_back(trimFastaFileName); }
568                         
569                         m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine();
570                 }
571                 
572                 //set accnos file as new current accnosfile
573                 string current = "";
574                 itTypes = outputTypes.find("accnos");
575                 if (itTypes != outputTypes.end()) {
576                         if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setAccnosFile(current); }
577                 }
578                 
579                 if (trim) {
580                         itTypes = outputTypes.find("fasta");
581                         if (itTypes != outputTypes.end()) {
582                                 if ((itTypes->second).size() != 0) { current = (itTypes->second)[0]; m->setFastaFile(current); }
583                         }
584                 }
585                 
586                 m->mothurOutEndLine();
587                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
588                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }       
589                 m->mothurOutEndLine();
590
591                 return 0;
592                 
593         }
594         catch(exception& e) {
595                 m->errorOut(e, "ChimeraSlayerCommand", "execute");
596                 exit(1);
597         }
598 }
599 //**********************************************************************************************************************
600
601 int ChimeraSlayerCommand::driver(linePair* filePos, string outputFName, string filename, string accnos, string fasta){
602         try {
603                 ofstream out;
604                 m->openOutputFile(outputFName, out);
605                 
606                 ofstream out2;
607                 m->openOutputFile(accnos, out2);
608                 
609                 ofstream out3;
610                 if (trim) {  m->openOutputFile(fasta, out3); }
611                 
612                 ifstream inFASTA;
613                 m->openInputFile(filename, inFASTA);
614
615                 inFASTA.seekg(filePos->start);
616
617                 bool done = false;
618                 int count = 0;
619         
620                 while (!done) {
621                 
622                         if (m->control_pressed) {       out.close(); out2.close(); if (trim) { out3.close(); } inFASTA.close(); return 1;       }
623                 
624                         Sequence* candidateSeq = new Sequence(inFASTA);  m->gobble(inFASTA);
625                         string candidateAligned = candidateSeq->getAligned();
626                                 
627                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
628                                 
629                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
630                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
631                                 }else{
632                                         //find chimeras
633                                         chimera->getChimeras(candidateSeq);
634                                         
635                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
636                                                 
637                                         //if you are not chimeric, then check each half
638                                         data_results wholeResults = chimera->getResults();
639                                         
640                                         //determine if we need to split
641                                         bool isChimeric = false;
642                                         
643                                         if (wholeResults.flag == "yes") {
644                                                 string chimeraFlag = "no";
645                                                 if(  (wholeResults.results[0].bsa >= minBS && wholeResults.results[0].divr_qla_qrb >= divR)
646                                                    ||
647                                                    (wholeResults.results[0].bsb >= minBS && wholeResults.results[0].divr_qlb_qra >= divR) ) { chimeraFlag = "yes"; }
648                                                 
649                                                 
650                                                 if (chimeraFlag == "yes") {     
651                                                         if ((wholeResults.results[0].bsa >= minBS) || (wholeResults.results[0].bsb >= minBS)) { isChimeric = true; }
652                                                 }
653                                         }
654                                         
655                                         if ((!isChimeric) && trimera) {
656                                                 
657                                                 //split sequence in half by bases
658                                                 string leftQuery, rightQuery;
659                                                 Sequence tempSeq(candidateSeq->getName(), candidateAligned);
660                                                 divideInHalf(tempSeq, leftQuery, rightQuery);
661                                                 
662                                                 //run chimeraSlayer on each piece
663                                                 Sequence* left = new Sequence(candidateSeq->getName(), leftQuery);
664                                                 Sequence* right = new Sequence(candidateSeq->getName(), rightQuery);
665                                                 
666                                                 //find chimeras
667                                                 chimera->getChimeras(left);
668                                                 data_results leftResults = chimera->getResults();
669                                                 
670                                                 chimera->getChimeras(right);
671                                                 data_results rightResults = chimera->getResults();
672                                                 
673                                                 //if either piece is chimeric then report
674                                                 Sequence* trimmed = chimera->print(out, out2, leftResults, rightResults);
675                                                 if (trim) { trimmed->printSequence(out3); delete trimmed; }
676                                                 
677                                                 delete left; delete right;
678                                                 
679                                         }else { //already chimeric
680                                                 //print results
681                                                 Sequence* trimmed = chimera->print(out, out2);
682                                                 if (trim) { trimmed->printSequence(out3); delete trimmed; }
683                                         }
684                                         
685                                         
686                                 }
687                                 count++;
688                         }
689                         delete candidateSeq;
690                         
691                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
692                                 unsigned long int pos = inFASTA.tellg();
693                                 if ((pos == -1) || (pos >= filePos->end)) { break; }
694                         #else
695                                 if (inFASTA.eof()) { break; }
696                         #endif
697                         
698                         //report progress
699                         if((count) % 100 == 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
700                 }
701                 //report progress
702                 if((count) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(count)); m->mothurOutEndLine();         }
703                 
704                 out.close();
705                 out2.close();
706                 if (trim) { out3.close(); }
707                 inFASTA.close();
708                                 
709                 return count;
710         }
711         catch(exception& e) {
712                 m->errorOut(e, "ChimeraSlayerCommand", "driver");
713                 exit(1);
714         }
715 }
716 //**********************************************************************************************************************
717 #ifdef USE_MPI
718 int ChimeraSlayerCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, MPI_File& outAccMPI, MPI_File& outFastaMPI, vector<unsigned long int>& MPIPos){
719         try {                           
720                 MPI_Status status; 
721                 int pid;
722                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
723                 
724                 for(int i=0;i<num;i++){
725                         
726                         if (m->control_pressed) {       return 1;       }
727                         
728                         //read next sequence
729                         int length = MPIPos[start+i+1] - MPIPos[start+i];
730
731                         char* buf4 = new char[length];
732                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
733         
734                         string tempBuf = buf4;
735                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
736                         istringstream iss (tempBuf,istringstream::in);
737
738                         delete buf4;
739
740                         Sequence* candidateSeq = new Sequence(iss);  m->gobble(iss);
741                         string candidateAligned = candidateSeq->getAligned();
742                 
743                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
744                                 
745                                 if (candidateSeq->getAligned().length() != templateSeqsLength) {  
746                                         m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
747                                 }else{
748                 
749                                         //find chimeras
750                                         chimera->getChimeras(candidateSeq);
751                         
752                                         if (m->control_pressed) {       delete candidateSeq; return 1;  }
753                                         
754                                         //if you are not chimeric, then check each half
755                                         data_results wholeResults = chimera->getResults();
756                                         
757                                         //determine if we need to split
758                                         bool isChimeric = false;
759                                         
760                                         if (wholeResults.flag == "yes") {
761                                                 string chimeraFlag = "no";
762                                                 if(  (wholeResults.results[0].bsa >= minBS && wholeResults.results[0].divr_qla_qrb >= divR)
763                                                    ||
764                                                    (wholeResults.results[0].bsb >= minBS && wholeResults.results[0].divr_qlb_qra >= divR) ) { chimeraFlag = "yes"; }
765                                                 
766                                                 
767                                                 if (chimeraFlag == "yes") {     
768                                                         if ((wholeResults.results[0].bsa >= minBS) || (wholeResults.results[0].bsb >= minBS)) { isChimeric = true; }
769                                                 }
770                                         }
771                                         
772                                         if ((!isChimeric) && trimera) {                                                 
773                                                 //split sequence in half by bases
774                                                 string leftQuery, rightQuery;
775                                                 Sequence tempSeq(candidateSeq->getName(), candidateAligned);
776                                                 divideInHalf(tempSeq, leftQuery, rightQuery);
777                                                 
778                                                 //run chimeraSlayer on each piece
779                                                 Sequence* left = new Sequence(candidateSeq->getName(), leftQuery);
780                                                 Sequence* right = new Sequence(candidateSeq->getName(), rightQuery);
781                                                 
782                                                 //find chimeras
783                                                 chimera->getChimeras(left);
784                                                 data_results leftResults = chimera->getResults();
785                                                 
786                                                 chimera->getChimeras(right);
787                                                 data_results rightResults = chimera->getResults();
788                                                 
789                                                 //if either piece is chimeric then report
790                                                 Sequence* trimmed = chimera->print(outMPI, outAccMPI, leftResults, rightResults);
791                                                 if (trim) {  
792                                                         string outputString = ">" + trimmed->getName() + "\n" + trimmed->getAligned() + "\n";
793                                                         delete trimmed;
794                                                         
795                                                         //write to accnos file
796                                                         int length = outputString.length();
797                                                         char* buf2 = new char[length];
798                                                         memcpy(buf2, outputString.c_str(), length);
799                                                         
800                                                         MPI_File_write_shared(outFastaMPI, buf2, length, MPI_CHAR, &status);
801                                                         delete buf2;
802                                                 }
803                                                 
804                                                 delete left; delete right;
805                                                 
806                                         }else { 
807                                                 //print results
808                                                 Sequence* trimmed = chimera->print(outMPI, outAccMPI);
809                                                 
810                                                 if (trim) {  
811                                                         string outputString = ">" + trimmed->getName() + "\n" + trimmed->getAligned() + "\n";
812                                                         delete trimmed;
813                                                         
814                                                         //write to accnos file
815                                                         int length = outputString.length();
816                                                         char* buf2 = new char[length];
817                                                         memcpy(buf2, outputString.c_str(), length);
818                                                         
819                                                         MPI_File_write_shared(outFastaMPI, buf2, length, MPI_CHAR, &status);
820                                                         delete buf2;
821                                                 }
822                                         }
823                                         
824                                 }
825                         }
826                         delete candidateSeq;
827                         
828                         //report progress
829                         if((i+1) % 100 == 0){  cout << "Processing sequence: " << (i+1) << endl;        m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n");          }
830                 }
831                 //report progress
832                 if(num % 100 != 0){             cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n");  }
833                 
834                                 
835                 return 0;
836         }
837         catch(exception& e) {
838                 m->errorOut(e, "ChimeraSlayerCommand", "driverMPI");
839                 exit(1);
840         }
841 }
842 #endif
843
844 /**************************************************************************************************/
845
846 int ChimeraSlayerCommand::createProcesses(string outputFileName, string filename, string accnos, string fasta) {
847         try {
848 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
849                 int process = 0;
850                 int num = 0;
851                 
852                 //loop through and create all the processes you want
853                 while (process != processors) {
854                         int pid = fork();
855                         
856                         if (pid > 0) {
857                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
858                                 process++;
859                         }else if (pid == 0){
860                                 num = driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename, accnos + toString(getpid()) + ".temp", fasta + toString(getpid()) + ".temp");
861                                 
862                                 //pass numSeqs to parent
863                                 ofstream out;
864                                 string tempFile = outputFileName + toString(getpid()) + ".num.temp";
865                                 m->openOutputFile(tempFile, out);
866                                 out << num << endl;
867                                 out.close();
868                                 
869                                 exit(0);
870                         }else { 
871                                 m->mothurOut("[ERROR]: unable to spawn the necessary processes."); m->mothurOutEndLine(); 
872                                 for (int i = 0; i < processIDS.size(); i++) { kill (processIDS[i], SIGINT); }
873                                 exit(0);
874                         }
875                 }
876                 
877                 //force parent to wait until all the processes are done
878                 for (int i=0;i<processors;i++) { 
879                         int temp = processIDS[i];
880                         wait(&temp);
881                 }
882                 
883                 for (int i = 0; i < processIDS.size(); i++) {
884                         ifstream in;
885                         string tempFile =  outputFileName + toString(processIDS[i]) + ".num.temp";
886                         m->openInputFile(tempFile, in);
887                         if (!in.eof()) { int tempNum = 0; in >> tempNum; num += tempNum; }
888                         in.close(); remove(tempFile.c_str());
889                 }
890                 
891                 return num;
892 #endif          
893         }
894         catch(exception& e) {
895                 m->errorOut(e, "ChimeraSlayerCommand", "createProcesses");
896                 exit(1);
897         }
898 }
899
900 /**************************************************************************************************/
901
902 int ChimeraSlayerCommand::divideInHalf(Sequence querySeq, string& leftQuery, string& rightQuery) {
903         try {
904                 
905                 string queryUnAligned = querySeq.getUnaligned();
906                 int numBases = int(queryUnAligned.length() * 0.5);
907                 
908                 string queryAligned = querySeq.getAligned();
909                 leftQuery = querySeq.getAligned();
910                 rightQuery = querySeq.getAligned();
911                 
912                 int baseCount = 0;
913                 int leftSpot = 0;
914                 for (int i = 0; i < queryAligned.length(); i++) {
915                         //if you are a base
916                         if (isalpha(queryAligned[i])) {         
917                                 baseCount++; 
918                         }
919                         
920                         //if you have half
921                         if (baseCount >= numBases) {  leftSpot = i; break; } //first half
922                 }
923                 
924                 //blank out right side
925                 for (int i = leftSpot; i < leftQuery.length(); i++) { leftQuery[i] = '.'; }
926                 
927                 //blank out left side
928                 for (int i = 0; i < leftSpot; i++) { rightQuery[i] = '.'; }
929                 
930                 return 0;
931                 
932         }
933         catch(exception& e) {
934                 m->errorOut(e, "ChimeraSlayerCommand", "divideInHalf");
935                 exit(1);
936         }
937 }
938
939 /**************************************************************************************************/
940