2 * classifyseqscommand.cpp
5 * Created by westcott on 11/2/09.
6 * Copyright 2009 Schloss Lab. All rights reserved.
10 #include "classifyseqscommand.h"
11 #include "sequence.hpp"
13 #include "phylotree.h"
16 //**********************************************************************************************************************
18 ClassifySeqsCommand::ClassifySeqsCommand(string option) {
22 //allow user to run help
23 if(option == "help") { help(); abort = true; }
27 //valid paramters for this command
28 string AlignArray[] = {"template","fasta","name","search","ksize","method","processors","taxonomy","match","mismatch","gapopen","gapextend","numwanted","cutoff","probs","iters", "outputdir","inputdir"};
29 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
31 OptionParser parser(option);
32 map<string, string> parameters = parser.getParameters();
34 ValidParameters validParameter;
35 map<string, string>::iterator it;
37 //check to make sure all parameters are valid for command
38 for (it = parameters.begin(); it != parameters.end(); it++) {
39 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
42 //if the user changes the output directory command factory will send this info to us in the output parameter
43 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){ outputDir = ""; }
45 //if the user changes the input directory command factory will send this info to us in the output parameter
46 string inputDir = validParameter.validFile(parameters, "inputdir", false);
47 if (inputDir == "not found"){ inputDir = ""; }
50 it = parameters.find("template");
51 //user has given a template file
52 if(it != parameters.end()){
53 path = hasPath(it->second);
54 //if the user has not given a path then, add inputdir. else leave path alone.
55 if (path == "") { parameters["template"] = inputDir + it->second; }
58 it = parameters.find("taxonomy");
59 //user has given a template file
60 if(it != parameters.end()){
61 path = hasPath(it->second);
62 //if the user has not given a path then, add inputdir. else leave path alone.
63 if (path == "") { parameters["taxonomy"] = inputDir + it->second; }
67 //check for required parameters
68 templateFileName = validParameter.validFile(parameters, "template", true);
69 if (templateFileName == "not found") {
70 m->mothurOut("template is a required parameter for the classify.seqs command.");
71 m->mothurOutEndLine();
74 else if (templateFileName == "not open") { abort = true; }
76 fastaFileName = validParameter.validFile(parameters, "fasta", false);
77 if (fastaFileName == "not found") { m->mothurOut("fasta is a required parameter for the classify.seqs command."); m->mothurOutEndLine(); abort = true; }
79 splitAtDash(fastaFileName, fastaFileNames);
81 //go through files and make sure they are good, if not, then disregard them
82 for (int i = 0; i < fastaFileNames.size(); i++) {
84 string path = hasPath(fastaFileNames[i]);
85 //if the user has not given a path then, add inputdir. else leave path alone.
86 if (path == "") { fastaFileNames[i] = inputDir + fastaFileNames[i]; }
93 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
94 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
100 ableToOpen = openInputFile(fastaFileNames[i], in);
104 for (int j = 1; j < processors; j++) {
105 MPI_Send(&ableToOpen, 1, MPI_INT, j, 2001, MPI_COMM_WORLD);
109 MPI_Recv(&ableToOpen, 1, MPI_INT, 0, 2001, MPI_COMM_WORLD, &status);
114 if (ableToOpen == 1) {
115 m->mothurOut(fastaFileNames[i] + " will be disregarded."); m->mothurOutEndLine();
116 //erase from file list
117 fastaFileNames.erase(fastaFileNames.begin()+i);
123 //make sure there is at least one valid file left
124 if (fastaFileNames.size() == 0) { m->mothurOut("no valid files."); m->mothurOutEndLine(); abort = true; }
128 taxonomyFileName = validParameter.validFile(parameters, "taxonomy", true);
129 if (taxonomyFileName == "not found") {
130 m->mothurOut("taxonomy is a required parameter for the classify.seqs command.");
131 m->mothurOutEndLine();
134 else if (taxonomyFileName == "not open") { abort = true; }
137 namefile = validParameter.validFile(parameters, "name", false);
138 if (namefile == "not found") { namefile = ""; }
141 splitAtDash(namefile, namefileNames);
143 //go through files and make sure they are good, if not, then disregard them
144 for (int i = 0; i < namefileNames.size(); i++) {
145 if (inputDir != "") {
146 string path = hasPath(namefileNames[i]);
147 //if the user has not given a path then, add inputdir. else leave path alone.
148 if (path == "") { namefileNames[i] = inputDir + namefileNames[i]; }
154 MPI_Comm_size(MPI_COMM_WORLD, &processors); //set processors to the number of mpi processes running
155 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
161 ableToOpen = openInputFile(namefileNames[i], in);
165 for (int j = 1; j < processors; j++) {
166 MPI_Send(&ableToOpen, 1, MPI_INT, j, 2001, MPI_COMM_WORLD);
170 MPI_Recv(&ableToOpen, 1, MPI_INT, 0, 2001, MPI_COMM_WORLD, &status);
174 if (ableToOpen == 1) { m->mothurOut("Unable to match name file with fasta file."); m->mothurOutEndLine(); abort = true; }
179 if (namefile != "") {
180 if (namefileNames.size() != fastaFileNames.size()) { abort = true; m->mothurOut("If you provide a name file, you must have one for each fasta file."); m->mothurOutEndLine(); }
183 //check for optional parameter and set defaults
184 // ...at some point should added some additional type checking...
186 temp = validParameter.validFile(parameters, "ksize", false); if (temp == "not found"){ temp = "8"; }
187 convert(temp, kmerSize);
189 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = "1"; }
190 convert(temp, processors);
192 search = validParameter.validFile(parameters, "search", false); if (search == "not found"){ search = "kmer"; }
194 method = validParameter.validFile(parameters, "method", false); if (method == "not found"){ method = "bayesian"; }
196 temp = validParameter.validFile(parameters, "match", false); if (temp == "not found"){ temp = "1.0"; }
197 convert(temp, match);
199 temp = validParameter.validFile(parameters, "mismatch", false); if (temp == "not found"){ temp = "-1.0"; }
200 convert(temp, misMatch);
202 temp = validParameter.validFile(parameters, "gapopen", false); if (temp == "not found"){ temp = "-2.0"; }
203 convert(temp, gapOpen);
205 temp = validParameter.validFile(parameters, "gapextend", false); if (temp == "not found"){ temp = "-1.0"; }
206 convert(temp, gapExtend);
208 temp = validParameter.validFile(parameters, "numwanted", false); if (temp == "not found"){ temp = "10"; }
209 convert(temp, numWanted);
211 temp = validParameter.validFile(parameters, "cutoff", false); if (temp == "not found"){ temp = "0"; }
212 convert(temp, cutoff);
214 temp = validParameter.validFile(parameters, "probs", false); if (temp == "not found"){ temp = "true"; }
215 probs = isTrue(temp);
217 temp = validParameter.validFile(parameters, "iters", false); if (temp == "not found") { temp = "100"; }
218 convert(temp, iters);
222 if ((method == "bayesian") && (search != "kmer")) {
223 m->mothurOut("The bayesian method requires the kmer search." + search + "will be disregarded." ); m->mothurOutEndLine();
229 catch(exception& e) {
230 m->errorOut(e, "ClassifySeqsCommand", "ClassifySeqsCommand");
235 //**********************************************************************************************************************
237 ClassifySeqsCommand::~ClassifySeqsCommand(){
239 if (abort == false) {
240 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
244 //**********************************************************************************************************************
246 void ClassifySeqsCommand::help(){
248 m->mothurOut("The classify.seqs command reads a fasta file containing sequences and creates a .taxonomy file and a .tax.summary file.\n");
249 m->mothurOut("The classify.seqs command parameters are template, fasta, name, search, ksize, method, taxonomy, processors, match, mismatch, gapopen, gapextend, numwanted and probs.\n");
250 m->mothurOut("The template, fasta and taxonomy parameters are required. You may enter multiple fasta files by separating their names with dashes. ie. fasta=abrecovery.fasta-amzon.fasta \n");
251 m->mothurOut("The search parameter allows you to specify the method to find most similar template. Your options are: suffix, kmer, blast and distance. The default is kmer.\n");
252 m->mothurOut("The name parameter allows you add a names file with your fasta file, if you enter multiple fasta files, you must enter matching names files for them.\n");
253 m->mothurOut("The method parameter allows you to specify classification method to use. Your options are: bayesian and knn. The default is bayesian.\n");
254 m->mothurOut("The ksize parameter allows you to specify the kmer size for finding most similar template to candidate. The default is 8.\n");
255 m->mothurOut("The processors parameter allows you to specify the number of processors to use. The default is 1.\n");
257 m->mothurOut("When using MPI, the processors parameter is set to the number of MPI processes running. \n");
259 m->mothurOut("The match parameter allows you to specify the bonus for having the same base. The default is 1.0.\n");
260 m->mothurOut("The mistmatch parameter allows you to specify the penalty for having different bases. The default is -1.0.\n");
261 m->mothurOut("The gapopen parameter allows you to specify the penalty for opening a gap in an alignment. The default is -2.0.\n");
262 m->mothurOut("The gapextend parameter allows you to specify the penalty for extending a gap in an alignment. The default is -1.0.\n");
263 m->mothurOut("The numwanted parameter allows you to specify the number of sequence matches you want with the knn method. The default is 10.\n");
264 m->mothurOut("The cutoff parameter allows you to specify a bootstrap confidence threshold for your taxonomy. The default is 0.\n");
265 m->mothurOut("The probs parameter shut off the bootstrapping results for the bayesian method. The default is true, meaning you want the bootstrapping to be run.\n");
266 m->mothurOut("The iters parameter allows you to specify how many iterations to do when calculating the bootstrap confidence score for your taxonomy with the bayesian method. The default is 100.\n");
267 m->mothurOut("The classify.seqs command should be in the following format: \n");
268 m->mothurOut("classify.seqs(template=yourTemplateFile, fasta=yourFastaFile, method=yourClassificationMethod, search=yourSearchmethod, ksize=yourKmerSize, taxonomy=yourTaxonomyFile, processors=yourProcessors) \n");
269 m->mothurOut("Example classify.seqs(fasta=amazon.fasta, template=core.filtered, method=knn, search=gotoh, ksize=8, processors=2)\n");
270 m->mothurOut("The .taxonomy file consists of 2 columns: 1 = your sequence name, 2 = the taxonomy for your sequence. \n");
271 m->mothurOut("The .tax.summary is a summary of the different taxonomies represented in your fasta file. \n");
272 m->mothurOut("Note: No spaces between parameter labels (i.e. fasta), '=' and parameters (i.e.yourFastaFile).\n\n");
274 catch(exception& e) {
275 m->errorOut(e, "ClassifySeqsCommand", "help");
281 //**********************************************************************************************************************
283 int ClassifySeqsCommand::execute(){
285 if (abort == true) { return 0; }
287 if(method == "bayesian"){ classify = new Bayesian(taxonomyFileName, templateFileName, search, kmerSize, cutoff, iters); }
288 else if(method == "knn"){ classify = new Knn(taxonomyFileName, templateFileName, search, kmerSize, gapOpen, gapExtend, match, misMatch, numWanted); }
290 m->mothurOut(search + " is not a valid method option. I will run the command using bayesian.");
291 m->mothurOutEndLine();
292 classify = new Bayesian(taxonomyFileName, templateFileName, search, kmerSize, cutoff, iters);
295 if (m->control_pressed) { delete classify; return 0; }
297 vector<string> outputNames;
299 for (int s = 0; s < fastaFileNames.size(); s++) {
301 m->mothurOut("Classifying sequences from " + fastaFileNames[s] + " ..." ); m->mothurOutEndLine();
303 if (outputDir == "") { outputDir += hasPath(fastaFileNames[s]); }
304 string newTaxonomyFile = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + getRootName(getSimpleName(taxonomyFileName)) + "taxonomy";
305 string tempTaxonomyFile = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + "taxonomy.temp";
306 string taxSummary = outputDir + getRootName(getSimpleName(fastaFileNames[s])) + getRootName(getSimpleName(taxonomyFileName)) + "tax.summary";
308 outputNames.push_back(newTaxonomyFile);
309 outputNames.push_back(taxSummary);
311 int start = time(NULL);
312 int numFastaSeqs = 0;
313 for (int i = 0; i < lines.size(); i++) { delete lines[i]; } lines.clear();
316 int pid, end, numSeqsPerProcessor;
321 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
322 MPI_Comm_size(MPI_COMM_WORLD, &processors);
325 MPI_File outMPINewTax;
326 MPI_File outMPITempTax;
328 int outMode=MPI_MODE_CREATE|MPI_MODE_WRONLY;
329 int inMode=MPI_MODE_RDONLY;
331 char* outNewTax = new char[newTaxonomyFile.length()];
\r
332 memcpy(outNewTax, newTaxonomyFile.c_str(), newTaxonomyFile.length());
334 char* outTempTax = new char[tempTaxonomyFile.length()];
\r
335 memcpy(outTempTax, tempTaxonomyFile.c_str(), tempTaxonomyFile.length());
337 char* inFileName = new char[fastaFileNames[s].length()];
\r
338 memcpy(inFileName, fastaFileNames[s].c_str(), fastaFileNames[s].length());
340 MPI_File_open(MPI_COMM_WORLD, inFileName, inMode, MPI_INFO_NULL, &inMPI); //comm, filename, mode, info, filepointer
341 MPI_File_open(MPI_COMM_WORLD, outNewTax, outMode, MPI_INFO_NULL, &outMPINewTax);
342 MPI_File_open(MPI_COMM_WORLD, outTempTax, outMode, MPI_INFO_NULL, &outMPITempTax);
348 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPINewTax); MPI_File_close(&outMPITempTax); delete classify; return 0; }
350 if(namefile != "") { MPIReadNamesFile(namefileNames[s]); }
352 if (pid == 0) { //you are the root process
354 MPIPos = setFilePosFasta(fastaFileNames[s], numFastaSeqs); //fills MPIPos, returns numSeqs
356 //send file positions to all processes
357 MPI_Bcast(&numFastaSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD); //send numSeqs
358 MPI_Bcast(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //send file pos
360 //figure out how many sequences you have to align
361 numSeqsPerProcessor = numFastaSeqs / processors;
362 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
363 int startIndex = pid * numSeqsPerProcessor;
366 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPINewTax, outMPITempTax, MPIPos);
368 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPINewTax); MPI_File_close(&outMPITempTax); for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } delete classify; return 0; }
370 for (int i = 1; i < processors; i++) {
372 MPI_Recv(&done, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
374 }else{ //you are a child process
375 MPI_Bcast(&numFastaSeqs, 1, MPI_INT, 0, MPI_COMM_WORLD); //get numSeqs
376 MPIPos.resize(numFastaSeqs+1);
377 MPI_Bcast(&MPIPos[0], (numFastaSeqs+1), MPI_LONG, 0, MPI_COMM_WORLD); //get file positions
379 //figure out how many sequences you have to align
380 numSeqsPerProcessor = numFastaSeqs / processors;
381 if(pid == (processors - 1)){ numSeqsPerProcessor = numFastaSeqs - pid * numSeqsPerProcessor; }
382 int startIndex = pid * numSeqsPerProcessor;
385 driverMPI(startIndex, numSeqsPerProcessor, inMPI, outMPINewTax, outMPITempTax, MPIPos);
387 if (m->control_pressed) { MPI_File_close(&inMPI); MPI_File_close(&outMPINewTax); MPI_File_close(&outMPITempTax); delete classify; return 0; }
390 MPI_Send(&done, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
394 MPI_File_close(&inMPI);
395 MPI_File_close(&outMPINewTax);
396 MPI_File_close(&outMPITempTax);
401 nameMap.clear(); //remove old names
404 openInputFile(namefileNames[s], inNames);
406 string firstCol, secondCol;
407 while(!inNames.eof()) {
408 inNames >> firstCol >> secondCol; gobble(inNames);
409 nameMap[firstCol] = getNumNames(secondCol); //ex. seq1 seq1,seq3,seq5 -> seq1 = 3.
414 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
417 openInputFile(fastaFileNames[s], inFASTA);
418 numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
421 lines.push_back(new linePair(0, numFastaSeqs));
423 driver(lines[0], newTaxonomyFile, tempTaxonomyFile, fastaFileNames[s]);
426 vector<int> positions;
427 processIDS.resize(0);
430 openInputFile(fastaFileNames[s], inFASTA);
433 while(!inFASTA.eof()){
434 input = getline(inFASTA);
435 if (input.length() != 0) {
436 if(input[0] == '>'){ int pos = inFASTA.tellg(); positions.push_back(pos - input.length() - 1); }
441 numFastaSeqs = positions.size();
443 int numSeqsPerProcessor = numFastaSeqs / processors;
445 for (int i = 0; i < processors; i++) {
446 int startPos = positions[ i * numSeqsPerProcessor ];
447 if(i == processors - 1){
448 numSeqsPerProcessor = numFastaSeqs - i * numSeqsPerProcessor;
450 lines.push_back(new linePair(startPos, numSeqsPerProcessor));
452 createProcesses(newTaxonomyFile, tempTaxonomyFile, fastaFileNames[s]);
454 rename((newTaxonomyFile + toString(processIDS[0]) + ".temp").c_str(), newTaxonomyFile.c_str());
455 rename((tempTaxonomyFile + toString(processIDS[0]) + ".temp").c_str(), tempTaxonomyFile.c_str());
457 for(int i=1;i<processors;i++){
458 appendTaxFiles((newTaxonomyFile + toString(processIDS[i]) + ".temp"), newTaxonomyFile);
459 appendTaxFiles((tempTaxonomyFile + toString(processIDS[i]) + ".temp"), tempTaxonomyFile);
460 remove((newTaxonomyFile + toString(processIDS[i]) + ".temp").c_str());
461 remove((tempTaxonomyFile + toString(processIDS[i]) + ".temp").c_str());
467 openInputFile(fastaFileNames[s], inFASTA);
468 numFastaSeqs=count(istreambuf_iterator<char>(inFASTA),istreambuf_iterator<char>(), '>');
471 lines.push_back(new linePair(0, numFastaSeqs));
473 driver(lines[0], newTaxonomyFile, tempTaxonomyFile, fastaFileNames[s]);
478 if (pid == 0) { //this part does not need to be paralellized
481 //make taxonomy tree from new taxonomy file
482 PhyloTree taxaBrowser;
484 if (m->control_pressed) { for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } delete classify; return 0; }
487 openInputFile(tempTaxonomyFile, in);
489 //read in users taxonomy file and add sequences to tree
492 in >> name >> taxon; gobble(in);
494 if (m->control_pressed) { for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } remove(tempTaxonomyFile.c_str()); delete classify; return 0; }
496 if (namefile != "") {
497 itNames = nameMap.find(name);
499 if (itNames == nameMap.end()) {
500 m->mothurOut(name + " is not in your name file please correct."); m->mothurOutEndLine(); exit(1);
502 for (int i = 0; i < itNames->second; i++) {
503 taxaBrowser.addSeqToTree(name+toString(i), taxon); //add it as many times as there are identical seqs
506 }else { taxaBrowser.addSeqToTree(name, taxon); } //add it once
510 taxaBrowser.assignHeirarchyIDs(0);
512 if (m->control_pressed) { for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } remove(tempTaxonomyFile.c_str()); delete classify; return 0; }
514 taxaBrowser.binUnclassified();
516 remove(tempTaxonomyFile.c_str());
518 if (m->control_pressed) { for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } delete classify; return 0; }
523 openOutputFile(taxSummary, outTaxTree);
524 taxaBrowser.print(outTaxTree);
527 //output taxonomy with the unclassified bins added
529 openInputFile(newTaxonomyFile, inTax);
532 string unclass = newTaxonomyFile + ".unclass.temp";
533 openOutputFile(unclass, outTax);
535 //get maxLevel from phylotree so you know how many 'unclassified's to add
536 int maxLevel = taxaBrowser.getMaxLevel();
538 //read taxfile - this reading and rewriting is done to preserve the confidence scores.
539 while (!inTax.eof()) {
540 if (m->control_pressed) { for (int i = 0; i < outputNames.size(); i++) { remove(outputNames[i].c_str()); } remove(unclass.c_str()); delete classify; return 0; }
542 inTax >> name >> taxon; gobble(inTax);
544 string newTax = addUnclassifieds(taxon, maxLevel);
546 outTax << name << '\t' << newTax << endl;
551 remove(newTaxonomyFile.c_str());
552 rename(unclass.c_str(), newTaxonomyFile.c_str());
558 m->mothurOutEndLine();
559 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
560 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
561 m->mothurOutEndLine();
564 m->mothurOutEndLine();
565 m->mothurOut("It took " + toString(time(NULL) - start) + " secs to classify " + toString(numFastaSeqs) + " sequences."); m->mothurOutEndLine(); m->mothurOutEndLine();
571 catch(exception& e) {
572 m->errorOut(e, "ClassifySeqsCommand", "execute");
577 /**************************************************************************************************/
578 string ClassifySeqsCommand::addUnclassifieds(string tax, int maxlevel) {
580 string newTax, taxon;
583 //keep what you have counting the levels
584 while (tax.find_first_of(';') != -1) {
586 taxon = tax.substr(0,tax.find_first_of(';'))+';';
587 tax = tax.substr(tax.find_first_of(';')+1, tax.length());
592 //add "unclassified" until you reach maxLevel
593 while (level < maxlevel) {
594 newTax += "unclassified;";
600 catch(exception& e) {
601 m->errorOut(e, "ClassifySeqsCommand", "addUnclassifieds");
606 /**************************************************************************************************/
608 void ClassifySeqsCommand::createProcesses(string taxFileName, string tempTaxFile, string filename) {
610 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
612 // processIDS.resize(0);
614 //loop through and create all the processes you want
615 while (process != processors) {
619 processIDS.push_back(pid); //create map from line number to pid so you can append files in correct order later
622 driver(lines[process], taxFileName + toString(getpid()) + ".temp", tempTaxFile + toString(getpid()) + ".temp", filename);
624 }else { m->mothurOut("unable to spawn the necessary processes."); m->mothurOutEndLine(); exit(0); }
627 //force parent to wait until all the processes are done
628 for (int i=0;i<processors;i++) {
629 int temp = processIDS[i];
634 catch(exception& e) {
635 m->errorOut(e, "ClassifySeqsCommand", "createProcesses");
639 /**************************************************************************************************/
641 void ClassifySeqsCommand::appendTaxFiles(string temp, string filename) {
646 openOutputFileAppend(filename, output);
647 openInputFile(temp, input);
649 while(char c = input.get()){
650 if(input.eof()) { break; }
651 else { output << c; }
657 catch(exception& e) {
658 m->errorOut(e, "ClassifySeqsCommand", "appendTaxFiles");
663 //**********************************************************************************************************************
665 int ClassifySeqsCommand::driver(linePair* line, string taxFName, string tempTFName, string filename){
668 openOutputFile(taxFName, outTax);
670 ofstream outTaxSimple;
671 openOutputFile(tempTFName, outTaxSimple);
674 openInputFile(filename, inFASTA);
676 inFASTA.seekg(line->start);
680 for(int i=0;i<line->numSeqs;i++){
681 if (m->control_pressed) { return 0; }
683 Sequence* candidateSeq = new Sequence(inFASTA);
685 if (candidateSeq->getName() != "") {
686 taxonomy = classify->getTaxonomy(candidateSeq);
688 if (m->control_pressed) { delete candidateSeq; return 0; }
690 if (taxonomy != "bad seq") {
691 //output confidence scores or not
693 outTax << candidateSeq->getName() << '\t' << taxonomy << endl;
695 outTax << candidateSeq->getName() << '\t' << classify->getSimpleTax() << endl;
698 outTaxSimple << candidateSeq->getName() << '\t' << classify->getSimpleTax() << endl;
703 if((i+1) % 100 == 0){
704 m->mothurOut("Classifying sequence " + toString(i+1)); m->mothurOutEndLine();
710 outTaxSimple.close();
714 catch(exception& e) {
715 m->errorOut(e, "ClassifySeqsCommand", "driver");
719 //**********************************************************************************************************************
721 int ClassifySeqsCommand::driverMPI(int start, int num, MPI_File& inMPI, MPI_File& newFile, MPI_File& tempFile, vector<long>& MPIPos){
723 MPI_Status statusNew;
724 MPI_Status statusTemp;
728 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
733 for(int i=0;i<num;i++){
735 if (m->control_pressed) { return 0; }
738 int length = MPIPos[start+i+1] - MPIPos[start+i];
739 char* buf4 = new char[length];
740 MPI_File_read_at(inMPI, MPIPos[start+i], buf4, length, MPI_CHAR, &status);
742 string tempBuf = buf4;
743 if (tempBuf.length() > length) { tempBuf = tempBuf.substr(0, length); }
744 istringstream iss (tempBuf,istringstream::in);
747 Sequence* candidateSeq = new Sequence(iss);
749 if (candidateSeq->getName() != "") {
750 taxonomy = classify->getTaxonomy(candidateSeq);
752 if (taxonomy != "bad seq") {
753 //output confidence scores or not
755 outputString = candidateSeq->getName() + "\t" + taxonomy + "\n";
757 outputString = candidateSeq->getName() + "\t" + classify->getSimpleTax() + "\n";
760 int length = outputString.length();
761 char* buf2 = new char[length];
\r
762 memcpy(buf2, outputString.c_str(), length);
764 MPI_File_write_shared(newFile, buf2, length, MPI_CHAR, &statusNew);
767 outputString = candidateSeq->getName() + "\t" + classify->getSimpleTax() + "\n";
768 length = outputString.length();
769 char* buf = new char[length];
\r
770 memcpy(buf, outputString.c_str(), length);
772 MPI_File_write_shared(tempFile, buf, length, MPI_CHAR, &statusTemp);
778 if((i+1) % 100 == 0){ cout << "Classifying sequence " << (i+1) << endl; }
781 if(num % 100 != 0){ cout << "Classifying sequence " << (num) << endl; }
786 catch(exception& e) {
787 m->errorOut(e, "ClassifySeqsCommand", "driverMPI");
792 //**********************************************************************************************************************
793 int ClassifySeqsCommand::MPIReadNamesFile(string nameFilename){
796 nameMap.clear(); //remove old names
802 char* inFileName = new char[nameFilename.length()];
\r
803 memcpy(inFileName, nameFilename.c_str(), nameFilename.length());
805 MPI_File_open(MPI_COMM_WORLD, inFileName, MPI_MODE_RDONLY, MPI_INFO_NULL, &inMPI);
806 MPI_File_get_size(inMPI, &size);
809 char* buffer = new char[size];
810 MPI_File_read(inMPI, buffer, size, MPI_CHAR, &status);
812 string tempBuf = buffer;
813 if (tempBuf.length() > size) { tempBuf = tempBuf.substr(0, size); }
814 istringstream iss (tempBuf,istringstream::in);
817 string firstCol, secondCol;
819 iss >> firstCol >> secondCol; gobble(iss);
820 nameMap[firstCol] = getNumNames(secondCol); //ex. seq1 seq1,seq3,seq5 -> seq1 = 3.
823 MPI_File_close(&inMPI);
827 catch(exception& e) {
828 m->errorOut(e, "ClassifySeqsCommand", "MPIReadNamesFile");
833 /**************************************************************************************************/