2 * chimeraslayercommand.cpp
5 * Created by westcott on 3/31/10.
6 * Copyright 2010 Schloss Lab. All rights reserved.
10 #include "chimeraslayercommand.h"
11 #include "chimeraslayer.h"
14 //***************************************************************************************************************
16 ChimeraSlayerCommand::ChimeraSlayerCommand(string option) {
20 //allow user to run help
21 if(option == "help") { help(); abort = true; }
24 //valid paramters for this command
25 string Array[] = {"fasta", "processors", "window", "template","numwanted", "ksize", "match","mismatch",
26 "divergence", "minsim","mincov","minbs", "minsnp","parents", "iters","outputdir","inputdir", "search","realign" };
27 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
29 OptionParser parser(option);
30 map<string,string> parameters = parser.getParameters();
32 ValidParameters validParameter;
33 map<string,string>::iterator it;
35 //check to make sure all parameters are valid for command
36 for (it = parameters.begin(); it != parameters.end(); it++) {
37 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
40 //if the user changes the input directory command factory will send this info to us in the output parameter
41 string inputDir = validParameter.validFile(parameters, "inputdir", false);
42 if (inputDir == "not found"){ inputDir = ""; }
45 it = parameters.find("template");
46 //user has given a template file
47 if(it != parameters.end()){
48 path = hasPath(it->second);
49 //if the user has not given a path then, add inputdir. else leave path alone.
50 if (path == "") { parameters["template"] = inputDir + it->second; }
55 //check for required parameters
56 fastafile = validParameter.validFile(parameters, "fasta", false);
57 if (fastafile == "not found") { fastafile = ""; m->mothurOut("fasta is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true; }
59 splitAtDash(fastafile, fastaFileNames);
61 //go through files and make sure they are good, if not, then disregard them
62 for (int i = 0; i < fastaFileNames.size(); i++) {
64 string path = hasPath(fastaFileNames[i]);
65 //if the user has not given a path then, add inputdir. else leave path alone.
66 if (path == "") { fastaFileNames[i] = inputDir + fastaFileNames[i]; }
74 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
75 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
80 ableToOpen = openInputFile(fastaFileNames[i], in);
84 for (int j = 1; j < processors; j++) {
85 MPI_Send(&ableToOpen, 1, MPI_INT, j, 2001, MPI_COMM_WORLD);
89 MPI_Recv(&ableToOpen, 1, MPI_INT, 0, 2001, MPI_COMM_WORLD, &status);
94 if (ableToOpen == 1) {
95 m->mothurOut(fastaFileNames[i] + " will be disregarded."); m->mothurOutEndLine();
96 //erase from file list
97 fastaFileNames.erase(fastaFileNames.begin()+i);
102 //make sure there is at least one valid file left
103 if (fastaFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
106 //if the user changes the output directory command factory will send this info to us in the output parameter
107 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
109 outputDir += hasPath(fastafile); //if user entered a file with a path then preserve it
112 templatefile = validParameter.validFile(parameters, "template", true);
113 if (templatefile == "not open") { abort = true; }
114 else if (templatefile == "not found") { templatefile = ""; m->mothurOut("template is a required parameter for the chimera.slayer command."); m->mothurOutEndLine(); abort = true; }
116 string temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found") { temp = "1"; }
117 convert(temp, processors);
119 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found") { temp = "7"; }
120 convert(temp, ksize);
122 temp = validParameter.validFile(parameters, "window", false); if (temp == "not found") { temp = "50"; }
123 convert(temp, window);
125 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found") { temp = "5"; }
126 convert(temp, match);
128 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found") { temp = "-4"; }
129 convert(temp, mismatch);
131 temp = validParameter.validFile(parameters, "divergence", false); if (temp == "not found") { temp = "1.007"; }
134 temp = validParameter.validFile(parameters, "minsim", false); if (temp == "not found") { temp = "90"; }
135 convert(temp, minSimilarity);
137 temp = validParameter.validFile(parameters, "mincov", false); if (temp == "not found") { temp = "70"; }
138 convert(temp, minCoverage);
140 temp = validParameter.validFile(parameters, "minbs", false); if (temp == "not found") { temp = "90"; }
141 convert(temp, minBS);
143 temp = validParameter.validFile(parameters, "minsnp", false); if (temp == "not found") { temp = "10"; }
144 convert(temp, minSNP);
146 temp = validParameter.validFile(parameters, "parents", false); if (temp == "not found") { temp = "3"; }
147 convert(temp, parents);
149 temp = validParameter.validFile(parameters, "realign", false); if (temp == "not found") { temp = "f"; }
150 realign = isTrue(temp);
152 search = validParameter.validFile(parameters, "search", false); if (search == "not found") { search = "distance"; }
154 temp = validParameter.validFile(parameters, "iters", false); if (temp == "not found") { temp = "100"; }
155 convert(temp, iters);
157 temp = validParameter.validFile(parameters, "increment", false); if (temp == "not found") { temp = "5"; }
158 convert(temp, increment);
160 temp = validParameter.validFile(parameters, "numwanted", false); if (temp == "not found") { temp = "15"; }
161 convert(temp, numwanted);
163 if ((search != "distance") && (search != "blast") && (search != "kmer")) { m->mothurOut(search + " is not a valid search."); m->mothurOutEndLine(); abort = true; }
166 catch(exception& e) {
167 m->errorOut(e, "ChimeraSlayerCommand", "ChimeraSlayerCommand");
171 //**********************************************************************************************************************
173 void ChimeraSlayerCommand::help(){
176 m->mothurOut("The chimera.slayer command reads a fastafile and templatefile and outputs potentially chimeric sequences.\n");
177 m->mothurOut("This command was modeled after the chimeraSlayer written by the Broad Institute.\n");
178 m->mothurOut("The chimera.slayer command parameters are fasta, template, processors, ksize, window, match, mismatch, divergence. minsim, mincov, minbs, minsnp, parents, search, iters, increment and numwanted.\n"); //realign,
179 m->mothurOut("The fasta parameter allows you to enter the fasta file containing your potentially chimeric sequences, and is required. \n");
180 m->mothurOut("You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amzon.fasta \n");
181 m->mothurOut("The template parameter allows you to enter a template file containing known non-chimeric sequences, and is required. \n");
182 m->mothurOut("The processors parameter allows you to specify how many processors you would like to use. The default is 1. \n");
184 m->mothurOut("When using MPI, the processors parameter is set to the number of MPI processes running. \n");
186 m->mothurOut("The window parameter allows you to specify the window size for searching for chimeras, default=50. \n");
187 m->mothurOut("The increment parameter allows you to specify how far you move each window while finding chimeric sequences, default=5.\n");
188 m->mothurOut("The numwanted parameter allows you to specify how many sequences you would each query sequence compared with, default=15.\n");
189 m->mothurOut("The ksize parameter allows you to input kmersize, default is 7, used if search is kmer. \n");
190 m->mothurOut("The match parameter allows you to reward matched bases in blast search, default is 5. \n");
191 m->mothurOut("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");
192 m->mothurOut("The mismatch parameter allows you to penalize mismatched bases in blast search, default is -4. \n");
193 m->mothurOut("The divergence parameter allows you to set a cutoff for chimera determination, default is 1.007. \n");
194 m->mothurOut("The iters parameter allows you to specify the number of bootstrap iters to do with the chimeraslayer method, default=100.\n");
195 m->mothurOut("The minsim parameter allows you to specify a minimum similarity with the parent fragments, default=90. \n");
196 m->mothurOut("The mincov parameter allows you to specify minimum coverage by closest matches found in template. Default is 70, meaning 70%. \n");
197 m->mothurOut("The minbs parameter allows you to specify minimum bootstrap support for calling a sequence chimeric. Default is 90, meaning 90%. \n");
198 m->mothurOut("The minsnp parameter allows you to specify percent of SNPs to sample on each side of breakpoint for computing bootstrap support (default: 10) \n");
199 m->mothurOut("The search parameter allows you to specify search method for finding the closest parent. Choices are distance, blast, and kmer, default distance. \n");
200 m->mothurOut("The realign parameter allows you to realign the query to the potential parents. Choices are true or false, default false. \n");
201 m->mothurOut("The chimera.slayer command should be in the following format: \n");
202 m->mothurOut("chimera.slayer(fasta=yourFastaFile, template=yourTemplate, search=yourSearch) \n");
203 m->mothurOut("Example: chimera.slayer(fasta=AD.align, template=core_set_aligned.imputed.fasta, search=kmer) \n");
204 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");
206 catch(exception& e) {
207 m->errorOut(e, "ChimeraSlayerCommand", "help");
212 //***************************************************************************************************************
214 ChimeraSlayerCommand::~ChimeraSlayerCommand(){ /* do nothing */ }
216 //***************************************************************************************************************
218 int ChimeraSlayerCommand::execute(){
221 if (abort == true) { return 0; }
223 for (int s = 0; s < fastaFileNames.size(); s++) {
225 m->mothurOut("Checking sequences from " + fastaFileNames[s] + " ..." ); m->mothurOutEndLine();
227 int start = time(NULL);
229 chimera = new ChimeraSlayer(fastaFileNames[s], templatefile, search, ksize, match, mismatch, window, divR, minSimilarity, minCoverage, minBS, minSNP, parents, iters, increment, numwanted, realign);
231 string outputFileName = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + "slayer.chimeras";
232 string accnosFileName = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + "slayer.accnos";
233 bool hasAccnos = true;
235 if (m->control_pressed) { delete chimera; for (int j = 0; j < outputNames.size(); j++) { remove(outputNames[j].c_str()); } return 0; }
237 if (chimera->getUnaligned()) {
238 m->mothurOut("Your template sequences are different lengths, please correct."); m->mothurOutEndLine();
242 templateSeqsLength = chimera->getLength();
245 int pid, end, numSeqsPerProcessor;
248 MPIWroteAccnos = false;
251 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
252 MPI_Comm_size(MPI_COMM_WORLD, &processors);
256 MPI_File outMPIAccnos;
258 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
259 int inMode=MPI_MODE_RDONLY;
261 char outFilename[1024];
262 strcpy(outFilename, outputFileName.c_str());
264 char outAccnosFilename[1024];
265 strcpy(outAccnosFilename, accnosFileName.c_str());
267 char inFileName[1024];
268 strcpy(inFileName, fastaFileNames[s].c_str());
270 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
271 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
272 MPI_File_open(MPI_COMM_WORLD, outAccnosFilename, outMode, MPI_INFO_NULL, &outMPIAccnos);
274 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); MPI_File_close(&outMPIAccnos); for (int j = 0; j < outputNames.size(); j++) { remove(outputNames[j].c_str()); } delete chimera; return 0; }
276 if (pid == 0) { //you are the root process
277 m->mothurOutEndLine();
278 m->mothurOut("Only reporting sequence supported by " + toString(minBS) + "% of bootstrapped results.");
279 m->mothurOutEndLine();
281 string outTemp = "Name\tLeftParent\tRightParent\tDivQLAQRB\tPerIDQLAQRB\tBootStrapA\tDivQLBQRA\tPerIDQLBQRA\tBootStrapB\tFlag\tLeftWindow\tRightWindow\n";
284 int length = outTemp.length();
285 char* buf2 = new char[length];
286 memcpy(buf2, outTemp.c_str(), length);
288 MPI_File_write_shared(outMPI, buf2, length, MPI_CHAR, &status);
291 MPIPos = setFilePosFasta(fastaFileNames[s], numSeqs); //fills MPIPos, returns numSeqs
293 //send file positions to all processes
294 for(int i = 1; i < processors; i++) {
295 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
296 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, i, tag, MPI_COMM_WORLD);
299 //figure out how many sequences you have to align
300 numSeqsPerProcessor = numSeqs / processors;
301 int startIndex = pid * numSeqsPerProcessor;
302 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
305 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
307 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); 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; }
309 for (int i = 1; i < processors; i++) {
311 MPI_Recv(&tempResult, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
312 if (tempResult != 0) { MPIWroteAccnos = true; }
314 }else{ //you are a child process
315 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
316 MPIPos.resize(numSeqs+1);
317 MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
319 //figure out how many sequences you have to align
320 numSeqsPerProcessor = numSeqs / processors;
321 int startIndex = pid * numSeqsPerProcessor;
322 if(pid == (processors - 1)){ numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor; }
325 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, outMPIAccnos, MPIPos);
327 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPI); MPI_File_close(&outMPIAccnos); for (int j = 0; j < outputNames.size(); j++) { remove(outputNames[j].c_str()); } delete chimera; return 0; }
329 MPI_Send(&MPIWroteAccnos, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
333 MPI_File_close(&inMPI);
334 MPI_File_close(&outMPI);
335 MPI_File_close(&outMPIAccnos);
336 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
338 //delete accnos file if blank
340 if (!MPIWroteAccnos) {
342 remove(accnosFileName.c_str());
348 string tempHeader = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + "slayer.chimeras.tempHeader";
349 openOutputFile(tempHeader, outHeader);
351 chimera->printHeader(outHeader);
355 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
358 openInputFile(fastaFileNames[s], inFASTA);
359 getNumSeqs(inFASTA, numSeqs);
362 lines.push_back(new linePair(0, numSeqs));
364 driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName);
366 if (m->control_pressed) {
367 remove(outputFileName.c_str());
368 remove(tempHeader.c_str());
369 remove(accnosFileName.c_str());
370 for (int j = 0; j < outputNames.size(); j++) { remove(outputNames[j].c_str()); }
371 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
376 //delete accnos file if its blank
377 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
380 vector<int> positions;
381 processIDS.resize(0);
384 openInputFile(fastaFileNames[s], inFASTA);
387 while(!inFASTA.eof()){
388 input = getline(inFASTA);
389 if (input.length() != 0) {
390 if(input[0] == '>'){ long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
395 numSeqs = positions.size();
397 int numSeqsPerProcessor = numSeqs / processors;
399 for (int i = 0; i < processors; i++) {
400 long int startPos = positions[ i * numSeqsPerProcessor ];
401 if(i == processors - 1){
402 numSeqsPerProcessor = numSeqs - i * numSeqsPerProcessor;
404 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
407 createProcesses(outputFileName, fastaFileNames[s], accnosFileName);
409 rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
411 //append output files
412 for(int i=1;i<processors;i++){
413 appendFiles((outputFileName + toString(processIDS[i]) + ".temp"), outputFileName);
414 remove((outputFileName + toString(processIDS[i]) + ".temp").c_str());
417 vector<string> nonBlankAccnosFiles;
418 //delete blank accnos files generated with multiple processes
419 for(int i=0;i<processors;i++){
420 if (!(isBlank(accnosFileName + toString(processIDS[i]) + ".temp"))) {
421 nonBlankAccnosFiles.push_back(accnosFileName + toString(processIDS[i]) + ".temp");
422 }else { remove((accnosFileName + toString(processIDS[i]) + ".temp").c_str()); }
425 //append accnos files
426 if (nonBlankAccnosFiles.size() != 0) {
427 rename(nonBlankAccnosFiles[0].c_str(), accnosFileName.c_str());
429 for (int h=1; h < nonBlankAccnosFiles.size(); h++) {
430 appendFiles(nonBlankAccnosFiles[h], accnosFileName);
431 remove(nonBlankAccnosFiles[h].c_str());
433 }else{ hasAccnos = false; }
435 if (m->control_pressed) {
436 remove(outputFileName.c_str());
437 remove(accnosFileName.c_str());
438 for (int j = 0; j < outputNames.size(); j++) { remove(outputNames[j].c_str()); }
439 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
448 openInputFile(fastaFileNames[s], inFASTA);
449 getNumSeqs(inFASTA, numSeqs);
451 lines.push_back(new linePair(0, numSeqs));
453 driver(lines[0], outputFileName, fastaFileNames[s], accnosFileName);
455 if (m->control_pressed) {
456 remove(outputFileName.c_str());
457 remove(tempHeader.c_str());
458 remove(accnosFileName.c_str());
459 for (int j = 0; j < outputNames.size(); j++) { remove(outputNames[j].c_str()); }
460 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
465 //delete accnos file if its blank
466 if (isBlank(accnosFileName)) { remove(accnosFileName.c_str()); hasAccnos = false; }
469 appendFiles(outputFileName, tempHeader);
471 remove(outputFileName.c_str());
472 rename(tempHeader.c_str(), outputFileName.c_str());
478 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
480 outputNames.push_back(outputFileName);
481 if (hasAccnos) { outputNames.push_back(accnosFileName); }
483 m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine();
486 m->mothurOutEndLine();
487 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
488 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
489 m->mothurOutEndLine();
494 catch(exception& e) {
495 m->errorOut(e, "ChimeraSlayerCommand", "execute");
499 //**********************************************************************************************************************
501 int ChimeraSlayerCommand::driver(linePair* line, string outputFName, string filename, string accnos){
504 openOutputFile(outputFName, out);
507 openOutputFile(accnos, out2);
510 openInputFile(filename, inFASTA);
512 inFASTA.seekg(line->start);
514 for(int i=0;i<line->numSeqs;i++){
516 if (m->control_pressed) { return 1; }
518 Sequence* candidateSeq = new Sequence(inFASTA); gobble(inFASTA);
520 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
522 if (candidateSeq->getAligned().length() != templateSeqsLength) {
523 m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
526 chimera->getChimeras(candidateSeq);
528 if (m->control_pressed) { delete candidateSeq; return 1; }
531 chimera->print(out, out2);
537 if((i+1) % 100 == 0){ m->mothurOut("Processing sequence: " + toString(i+1)); m->mothurOutEndLine(); }
540 if((line->numSeqs) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(line->numSeqs)); m->mothurOutEndLine(); }
548 catch(exception& e) {
549 m->errorOut(e, "ChimeraSlayerCommand", "driver");
553 //**********************************************************************************************************************
555 int ChimeraSlayerCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, MPI_File& outAccMPI, vector<long>& MPIPos){
559 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
561 for(int i=0;i<num;i++){
563 if (m->control_pressed) { return 1; }
566 int length = MPIPos[start+i+1] - MPIPos[start+i];
568 char* buf4 = new char[length];
569 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
571 string tempBuf = buf4;
572 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
573 istringstream iss (tempBuf,istringstream::in);
577 Sequence* candidateSeq = new Sequence(iss); gobble(iss);
579 if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
581 if (candidateSeq->getAligned().length() != templateSeqsLength) {
582 m->mothurOut(candidateSeq->getName() + " is not the same length as the template sequences. Skipping."); m->mothurOutEndLine();
586 chimera->getChimeras(candidateSeq);
588 if (m->control_pressed) { delete candidateSeq; return 1; }
589 //cout << "about to print" << endl;
591 bool isChimeric = chimera->print(outMPI, outAccMPI);
592 if (isChimeric) { MPIWroteAccnos = true; }
599 if((i+1) % 100 == 0){ cout << "Processing sequence: " << (i+1) << endl; m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n"); }
602 if(num % 100 != 0){ cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n"); }
607 catch(exception& e) {
608 m->errorOut(e, "ChimeraSlayerCommand", "driverMPI");
614 /**************************************************************************************************/
616 int ChimeraSlayerCommand::createProcesses(string outputFileName, string filename, string accnos) {
618 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
620 // processIDS.resize(0);
622 //loop through and create all the processes you want
623 while (process != processors) {
627 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
630 driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename, accnos + toString(getpid()) + ".temp");
632 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
635 //force parent to wait until all the processes are done
636 for (int i=0;i<processors;i++) {
637 int temp = processIDS[i];
644 catch(exception& e) {
645 m->errorOut(e, "ChimeraSlayerCommand", "createProcesses");
650 /**************************************************************************************************/