]> git.donarmstrong.com Git - mothur.git/blob - chimeracheckcommand.cpp
modified chimera commands to process multiple fasta files and added checks to pintail...
[mothur.git] / chimeracheckcommand.cpp
1 /*
2  *  chimeracheckcommand.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 "chimeracheckcommand.h"
11 #include "chimeracheckrdp.h"
12
13 //***************************************************************************************************************
14
15 ChimeraCheckCommand::ChimeraCheckCommand(string option)  {
16         try {
17                 abort = false;
18                 
19                 //allow user to run help
20                 if(option == "help") { help(); abort = true; }
21                 
22                 else {
23                         //valid paramters for this command
24                         string Array[] =  {"fasta","processors","increment","template","ksize","svg", "name","outputdir","inputdir" };
25                         vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
26                         
27                         OptionParser parser(option);
28                         map<string,string> parameters = parser.getParameters();
29                         
30                         ValidParameters validParameter;
31                         map<string,string>::iterator it;
32                         
33                         //check to make sure all parameters are valid for command
34                         for (it = parameters.begin(); it != parameters.end(); it++) { 
35                                 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) {  abort = true;  }
36                         }
37                         
38                         //if the user changes the input directory command factory will send this info to us in the output parameter 
39                         string inputDir = validParameter.validFile(parameters, "inputdir", false);              
40                         if (inputDir == "not found"){   inputDir = "";          }
41                         else {
42                                 it = parameters.find("template");
43                                 //user has given a template file
44                                 if(it != parameters.end()){ 
45                                         string path = hasPath(it->second);
46                                         //if the user has not given a path then, add inputdir. else leave path alone.
47                                         if (path == "") {       parameters["template"] = inputDir + it->second;         }
48                                 }
49                         }
50                         
51                         //check for required parameters
52                         fastafile = validParameter.validFile(parameters, "fasta", false);
53                         if (fastafile == "not found") { fastafile = ""; m->mothurOut("fasta is a required parameter for the chimera.check command."); m->mothurOutEndLine(); abort = true;  }
54                         else { 
55                                 splitAtDash(fastafile, fastaFileNames);
56                                 
57                                 //go through files and make sure they are good, if not, then disregard them
58                                 for (int i = 0; i < fastaFileNames.size(); i++) {
59                                         if (inputDir != "") {
60                                                 string path = hasPath(fastaFileNames[i]);
61                                                 //if the user has not given a path then, add inputdir. else leave path alone.
62                                                 if (path == "") {       fastaFileNames[i] = inputDir + fastaFileNames[i];               }
63                                         }
64         
65                                         int ableToOpen;
66                                         ifstream in;
67                                         
68                                         #ifdef USE_MPI  
69                                                 int pid;
70                                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
71                                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
72                                 
73                                                 if (pid == 0) {
74                                         #endif
75
76                                         ableToOpen = openInputFile(fastaFileNames[i], in);
77                                         in.close();
78                                         
79                                         #ifdef USE_MPI  
80                                                         for (int j = 1; j < processors; j++) {
81                                                                 MPI_Send(&ableToOpen, 1, MPI_INT, j, 2001, MPI_COMM_WORLD); 
82                                                         }
83                                                 }else{
84                                                         MPI_Status status;
85                                                         MPI_Recv(&ableToOpen, 1, MPI_INT, 0, 2001, MPI_COMM_WORLD, &status);
86                                                 }
87                                                 
88                                         #endif
89
90                                         if (ableToOpen == 1) { 
91                                                 m->mothurOut(fastaFileNames[i] + " will be disregarded."); m->mothurOutEndLine(); 
92                                                 //erase from file list
93                                                 fastaFileNames.erase(fastaFileNames.begin()+i);
94                                                 i--;
95                                         }
96                                 }
97                                 
98                                 //make sure there is at least one valid file left
99                                 if (fastaFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
100                         }
101                         
102                         //if the user changes the output directory command factory will send this info to us in the output parameter 
103                         outputDir = validParameter.validFile(parameters, "outputdir", false);           if (outputDir == "not found"){  outputDir = ""; }
104
105                         templatefile = validParameter.validFile(parameters, "template", true);
106                         if (templatefile == "not open") { abort = true; }
107                         else if (templatefile == "not found") { templatefile = "";  m->mothurOut("template is a required parameter for the chimera.check command."); m->mothurOutEndLine(); abort = true;  }    
108                         
109                         namefile = validParameter.validFile(parameters, "name", false);
110                         if (namefile == "not found") { namefile = ""; }
111                         else { 
112                                 splitAtDash(namefile, nameFileNames);
113                                 
114                                 //go through files and make sure they are good, if not, then disregard them
115                                 for (int i = 0; i < nameFileNames.size(); i++) {
116                                         if (inputDir != "") {
117                                                 string path = hasPath(nameFileNames[i]);
118                                                 //if the user has not given a path then, add inputdir. else leave path alone.
119                                                 if (path == "") {       nameFileNames[i] = inputDir + nameFileNames[i];         }
120                                         }
121         
122                                         int ableToOpen;
123                                         ifstream in;
124                                         
125                                         #ifdef USE_MPI  
126                                                 int pid;
127                                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
128                                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
129                                 
130                                                 if (pid == 0) {
131                                         #endif
132
133                                         ableToOpen = openInputFile(nameFileNames[i], in);
134                                         in.close();
135                                         
136                                         #ifdef USE_MPI  
137                                                         for (int j = 1; j < processors; j++) {
138                                                                 MPI_Send(&ableToOpen, 1, MPI_INT, j, 2001, MPI_COMM_WORLD); 
139                                                         }
140                                                 }else{
141                                                         MPI_Status status;
142                                                         MPI_Recv(&ableToOpen, 1, MPI_INT, 0, 2001, MPI_COMM_WORLD, &status);
143                                                 }
144                                                 
145                                         #endif
146
147                                         if (ableToOpen == 1) { 
148                                                 m->mothurOut(nameFileNames[i] + " will be disregarded."); m->mothurOutEndLine(); 
149                                                 //erase from file list
150                                                 nameFileNames.erase(nameFileNames.begin()+i);
151                                                 i--;
152                                         }
153                                         
154                                 }
155                                 
156                                 //make sure there is at least one valid file left
157                                 if (nameFileNames.size() != 0) {
158                                         if (nameFileNames.size() != fastaFileNames.size()) { 
159                                                  m->mothurOut("Different number of valid name files and fasta files, aborting command."); m->mothurOutEndLine(); 
160                                                  abort = true;
161                                         }
162                                 }
163                         }
164
165                         string temp = validParameter.validFile(parameters, "processors", false);                if (temp == "not found") { temp = "1"; }
166                         convert(temp, processors);
167                         
168                         temp = validParameter.validFile(parameters, "ksize", false);                    if (temp == "not found") { temp = "7"; }
169                         convert(temp, ksize);
170                         
171                         temp = validParameter.validFile(parameters, "svg", false);                              if (temp == "not found") { temp = "F"; }
172                         svg = isTrue(temp);
173                         if (nameFileNames.size() != 0) { svg = true; }
174                         
175                         temp = validParameter.validFile(parameters, "increment", false);                if (temp == "not found") { temp = "10"; }
176                         convert(temp, increment);                       
177                 }
178         }
179         catch(exception& e) {
180                 m->errorOut(e, "ChimeraCheckCommand", "ChimeraCheckCommand");
181                 exit(1);
182         }
183 }
184 //**********************************************************************************************************************
185
186 void ChimeraCheckCommand::help(){
187         try {
188         
189                 m->mothurOut("The chimera.check command reads a fastafile and templatefile and outputs potentially chimeric sequences.\n");
190                 m->mothurOut("This command was created using the algorythms described in CHIMERA_CHECK version 2.7 written by Niels Larsen. \n");
191                 m->mothurOut("The chimera.check command parameters are fasta, template, processors, ksize, increment, svg and name.\n");
192                 m->mothurOut("The fasta parameter allows you to enter the fasta file containing your potentially chimeric sequences, and is required. \n");
193                 m->mothurOut("You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amzon.fasta \n");
194                 m->mothurOut("The template parameter allows you to enter a template file containing known non-chimeric sequences, and is required. \n");
195                 m->mothurOut("The processors parameter allows you to specify how many processors you would like to use.  The default is 1. \n");
196                 #ifdef USE_MPI
197                 m->mothurOut("When using MPI, the processors parameter is set to the number of MPI processes running. \n");
198                 #endif
199                 m->mothurOut("The increment parameter allows you to specify how far you move each window while finding chimeric sequences, default is 10.\n");
200                 m->mothurOut("The ksize parameter allows you to input kmersize, default is 7. \n");
201                 m->mothurOut("The svg parameter allows you to specify whether or not you would like a svg file outputted for each query sequence, default is False.\n");
202                 m->mothurOut("The name parameter allows you to enter a file containing names of sequences you would like .svg files for.\n");
203                 m->mothurOut("You may enter multiple name files by separating their names with dashes. ie. fasta=abrecovery.svg.names-amzon.svg.names \n");
204                 m->mothurOut("The chimera.check command should be in the following format: \n");
205                 m->mothurOut("chimera.check(fasta=yourFastaFile, template=yourTemplateFile, processors=yourProcessors, ksize=yourKmerSize) \n");
206                 m->mothurOut("Example: chimera.check(fasta=AD.fasta, template=core_set_aligned,imputed.fasta, processors=4, ksize=8) \n");
207                 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");     
208         }
209         catch(exception& e) {
210                 m->errorOut(e, "ChimeraCheckCommand", "help");
211                 exit(1);
212         }
213 }
214
215 //***************************************************************************************************************
216
217 ChimeraCheckCommand::~ChimeraCheckCommand(){    /*      do nothing      */      }
218
219 //***************************************************************************************************************
220
221 int ChimeraCheckCommand::execute(){
222         try{
223                 
224                 if (abort == true) { return 0; }
225                 
226                 for (int i = 0; i < fastaFileNames.size(); i++) {
227                                 
228                         m->mothurOut("Checking sequences from " + fastaFileNames[i] + " ..." ); m->mothurOutEndLine();
229                         
230                         int start = time(NULL); 
231                         
232                         string thisNameFile = "";
233                         if (nameFileNames.size() != 0) { thisNameFile = nameFileNames[i]; }
234                         
235                         chimera = new ChimeraCheckRDP(fastaFileNames[i], templatefile, thisNameFile, svg, increment, ksize, outputDir);                 
236
237                         if (m->control_pressed) { delete chimera;       return 0;       }
238                         
239                         string outputFileName = outputDir + getRootName(getSimpleName(fastaFileNames[i]))  + "chimeracheck.chimeras";
240                         outputNames.push_back(outputFileName);
241                         
242                 #ifdef USE_MPI
243                 
244                                 int pid, end, numSeqsPerProcessor; 
245                                 int tag = 2001;
246                                 vector<long> MPIPos;
247                                 
248                                 MPI_Status status; 
249                                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
250                                 MPI_Comm_size(MPI_COMM_WORLD, &processors); 
251
252                                 MPI_File inMPI;
253                                 MPI_File outMPI;
254                                                         
255                                 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY; 
256                                 int inMode=MPI_MODE_RDONLY; 
257                                                         
258                                 char outFilename[1024];
259                                 strcpy(outFilename, outputFileName.c_str());
260                         
261                                 char inFileName[1024];
262                                 strcpy(inFileName, fastaFileNames[i].c_str());
263
264                                 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI);  //comm, filename, mode, info, filepointer
265                                 MPI_File_open(MPI_COMM_WORLD, outFilename, outMode, MPI_INFO_NULL, &outMPI);
266                                 
267                                 if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  for (int j = 0; j < outputNames.size(); j++) {    remove(outputNames[j].c_str()); } delete chimera; return 0;  }
268                                 
269                                 if (pid == 0) { //you are the root process 
270                                         MPIPos = setFilePosFasta(fastaFileNames[i], numSeqs); //fills MPIPos, returns numSeqs
271                                         
272                                         //send file positions to all processes
273                                         for(int j = 1; j < processors; j++) { 
274                                                 MPI_Send(&numSeqs, 1, MPI_INT, j, tag, MPI_COMM_WORLD);
275                                                 MPI_Send(&MPIPos[0], (numSeqs+1), MPI_LONG, j, tag, MPI_COMM_WORLD);
276                                         }       
277                                         
278                                         //figure out how many sequences you have to align
279                                         numSeqsPerProcessor = numSeqs / processors;
280                                         int startIndex =  pid * numSeqsPerProcessor;
281                                         if(pid == (processors - 1)){    numSeqsPerProcessor = numSeqs - pid * numSeqsPerProcessor;      }
282                                         
283                                 
284                                         //align your part
285                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
286                                         
287                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);  for (int j = 0; j < outputNames.size(); j++) {    remove(outputNames[j].c_str()); }   delete chimera; return 0;  }
288                                         
289                                         //wait on chidren
290                                         for(int j = 1; j < processors; j++) { 
291                                                 char buf[4];
292                                                 MPI_Recv(buf, 4, MPI_CHAR, j, tag, MPI_COMM_WORLD, &status); 
293                                         }
294                                 }else{ //you are a child process
295                                         MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
296                                         MPIPos.resize(numSeqs+1);
297                                         MPI_Recv(&MPIPos[0], (numSeqs+1), MPI_LONG, 0, tag, MPI_COMM_WORLD, &status);
298                                         
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;      }
303                                         
304                                         //align your part
305                                         driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPI, MPIPos);
306                                         
307                                         if (m->control_pressed) {  MPI_File_close(&inMPI);  MPI_File_close(&outMPI);   for (int j = 0; j < outputNames.size(); j++) {   remove(outputNames[j].c_str()); }  delete chimera; return 0;  }
308                                         
309                                         //tell parent you are done.
310                                         char buf[4];
311                                         strcpy(buf, "done"); 
312                                         MPI_Send(buf, 4, MPI_CHAR, 0, tag, MPI_COMM_WORLD);
313                                 }
314                                 
315                                 //close files 
316                                 MPI_File_close(&inMPI);
317                                 MPI_File_close(&outMPI);
318                                 MPI_Barrier(MPI_COMM_WORLD); //make everyone wait - just in case
319                 #else
320                         
321                         //break up file
322                         #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
323                                 if(processors == 1){
324                                         ifstream inFASTA;
325                                         openInputFile(fastaFileNames[i], inFASTA);
326                                         getNumSeqs(inFASTA, numSeqs);
327                                         inFASTA.close();
328                                         
329                                         lines.push_back(new linePair(0, numSeqs));
330                                         
331                                         driver(lines[0], outputFileName, fastaFileNames[i]);
332                                         
333                                         if (m->control_pressed) { 
334                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
335                                                 for (int j = 0; j < lines.size(); j++) {  delete lines[j];  }  lines.clear();
336                                                 delete chimera;
337                                                 return 0;
338                                         }
339                                                                         
340                                 }else{
341                                         vector<int> positions;
342                                         processIDS.resize(0);
343                                         
344                                         ifstream inFASTA;
345                                         openInputFile(fastaFileNames[i], inFASTA);
346                                         
347                                         string input;
348                                         while(!inFASTA.eof()){
349                                                 input = getline(inFASTA);
350                                                 if (input.length() != 0) {
351                                                         if(input[0] == '>'){    long int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1);  }
352                                                 }
353                                         }
354                                         inFASTA.close();
355                                         
356                                         numSeqs = positions.size();
357                                         
358                                         int numSeqsPerProcessor = numSeqs / processors;
359                                         
360                                         for (int j = 0; j < processors; j++) {
361                                                 long int startPos = positions[ j * numSeqsPerProcessor ];
362                                                 if(j == processors - 1){
363                                                         numSeqsPerProcessor = numSeqs - j * numSeqsPerProcessor;
364                                                 }
365                                                 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
366                                         }
367                                         
368                                         
369                                         createProcesses(outputFileName, fastaFileNames[i]); 
370                                 
371                                         rename((outputFileName + toString(processIDS[0]) + ".temp").c_str(), outputFileName.c_str());
372                                                 
373                                         //append output files
374                                         for(int j=1;j<processors;j++){
375                                                 appendFiles((outputFileName + toString(processIDS[j]) + ".temp"), outputFileName);
376                                                 remove((outputFileName + toString(processIDS[j]) + ".temp").c_str());
377                                         }
378                                         
379                                         if (m->control_pressed) { 
380                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
381                                                 for (int j = 0; j < lines.size(); j++) {  delete lines[j];  }  lines.clear();
382                                                 delete chimera;
383                                                 return 0;
384                                         }
385                                 }
386
387                         #else
388                                 ifstream inFASTA;
389                                 openInputFile(fastaFileNames[i], inFASTA);
390                                 getNumSeqs(inFASTA, numSeqs);
391                                 inFASTA.close();
392                                 lines.push_back(new linePair(0, numSeqs));
393                                 
394                                 driver(lines[0], outputFileName, fastaFileNames[i]);
395                                 
396                                 if (m->control_pressed) { 
397                                                 for (int j = 0; j < lines.size(); j++) {  delete lines[j];  }  lines.clear();
398                                                 for (int j = 0; j < outputNames.size(); j++) {  remove(outputNames[j].c_str()); } 
399                                                 delete chimera;
400                                                 return 0;
401                                 }
402                         #endif
403                 #endif          
404                         delete chimera;
405                         for (int j = 0; j < lines.size(); j++) {  delete lines[j];  }  lines.clear();
406                         
407                         m->mothurOutEndLine(); m->mothurOut("This method does not determine if a sequence is chimeric, but allows you to make that determination based on the IS values."); m->mothurOutEndLine(); 
408                         m->mothurOutEndLine(); m->mothurOut("It took " + toString(time(NULL) - start) + " secs to check " + toString(numSeqs) + " sequences."); m->mothurOutEndLine(); m->mothurOutEndLine();
409
410                 }
411                 
412                 m->mothurOutEndLine();
413                 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
414                 for (int i = 0; i < outputNames.size(); i++) {  m->mothurOut(outputNames[i]); m->mothurOutEndLine();    }       
415                 m->mothurOutEndLine();
416         
417                 return 0;
418                 
419         }
420         catch(exception& e) {
421                 m->errorOut(e, "ChimeraCheckCommand", "execute");
422                 exit(1);
423         }
424 }
425 //**********************************************************************************************************************
426
427 int ChimeraCheckCommand::driver(linePair* line, string outputFName, string filename){
428         try {
429                 ofstream out;
430                 openOutputFile(outputFName, out);
431                 
432                 ofstream out2;
433                 
434                 ifstream inFASTA;
435                 openInputFile(filename, inFASTA);
436
437                 inFASTA.seekg(line->start);
438                 
439                 for(int i=0;i<line->numSeqs;i++){
440                 
441                         if (m->control_pressed) {       return 1;       }
442                 
443                         Sequence* candidateSeq = new Sequence(inFASTA);  gobble(inFASTA);
444                                 
445                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
446                                 //find chimeras
447                                 chimera->getChimeras(candidateSeq);
448                                 
449                                 if (m->control_pressed) {       delete candidateSeq; return 1;  }
450         
451                                 //print results
452                                 chimera->print(out, out2);
453                         }
454                         delete candidateSeq;
455                         
456                         //report progress
457                         if((i+1) % 100 == 0){   m->mothurOut("Processing sequence: " + toString(i+1)); m->mothurOutEndLine();           }
458                 }
459                 //report progress
460                 if((line->numSeqs) % 100 != 0){ m->mothurOut("Processing sequence: " + toString(line->numSeqs)); m->mothurOutEndLine();         }
461                 
462                 out.close();
463                 inFASTA.close();
464                                 
465                 return 0;
466         }
467         catch(exception& e) {
468                 m->errorOut(e, "ChimeraCheckCommand", "driver");
469                 exit(1);
470         }
471 }
472 //**********************************************************************************************************************
473 #ifdef USE_MPI
474 int ChimeraCheckCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& outMPI, vector<long>& MPIPos){
475         try {
476                 MPI_File outAccMPI;
477                 MPI_Status status; 
478                 int pid;
479                 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
480                 
481                 for(int i=0;i<num;i++){
482                         
483                         if (m->control_pressed) { return 0; }
484                         
485                         //read next sequence
486                         int length = MPIPos[start+i+1] - MPIPos[start+i];
487         
488                         char* buf4 = new char[length];
489                         MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
490                         
491                         string tempBuf = buf4;
492                         if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length);  }
493                         istringstream iss (tempBuf,istringstream::in);
494                         delete buf4;
495
496                         Sequence* candidateSeq = new Sequence(iss);  gobble(iss);
497                                 
498                         if (candidateSeq->getName() != "") { //incase there is a commented sequence at the end of a file
499                                 //find chimeras
500                                 chimera->getChimeras(candidateSeq);
501                                         
502                                 //print results
503                                 chimera->print(outMPI, outAccMPI);
504                         }
505                         delete candidateSeq;
506                         
507                         //report progress
508                         if((i+1) % 100 == 0){  cout << "Processing sequence: " << (i+1) << endl;        m->mothurOutJustToLog("Processing sequence: " + toString(i+1) + "\n");          }
509                 }
510                 //report progress
511                 if(num % 100 != 0){             cout << "Processing sequence: " << num << endl; m->mothurOutJustToLog("Processing sequence: " + toString(num) + "\n");  }
512                 
513                 return 0;
514         }
515         catch(exception& e) {
516                 m->errorOut(e, "ChimeraCheckCommand", "driverMPI");
517                 exit(1);
518         }
519 }
520 #endif
521
522 /**************************************************************************************************/
523
524 int ChimeraCheckCommand::createProcesses(string outputFileName, string filename) {
525         try {
526 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
527                 int process = 0;
528                 //              processIDS.resize(0);
529                 
530                 //loop through and create all the processes you want
531                 while (process != processors) {
532                         int pid = fork();
533                         
534                         if (pid > 0) {
535                                 processIDS.push_back(pid);  //create map from line number to pid so you can append files in correct order later
536                                 process++;
537                         }else if (pid == 0){
538                                 driver(lines[process], outputFileName + toString(getpid()) + ".temp", filename);
539                                 exit(0);
540                         }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
541                 }
542                 
543                 //force parent to wait until all the processes are done
544                 for (int i=0;i<processors;i++) { 
545                         int temp = processIDS[i];
546                         wait(&temp);
547                 }
548                 
549                 return 0;
550 #endif          
551         }
552         catch(exception& e) {
553                 m->errorOut(e, "ChimeraCheckCommand", "createProcesses");
554                 exit(1);
555         }
556 }
557 /**************************************************************************************************/
558
559