5 * Created by Pat Schloss on 12/27/10.
6 * Copyright 2010 Schloss Lab. All rights reserved.
10 #include "shhhercommand.h"
12 #include "readcolumn.h"
13 #include "readmatrix.hpp"
14 #include "rabundvector.hpp"
15 #include "sabundvector.hpp"
16 #include "listvector.hpp"
17 #include "cluster.hpp"
18 #include "sparsematrix.hpp"
21 //**********************************************************************************************************************
22 vector<string> ShhherCommand::setParameters(){
24 CommandParameter pflow("flow", "InputTypes", "", "", "none", "fileflow", "none",false,false); parameters.push_back(pflow);
25 CommandParameter pfile("file", "InputTypes", "", "", "none", "fileflow", "none",false,false); parameters.push_back(pfile);
26 CommandParameter plookup("lookup", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(plookup);
27 CommandParameter pcutoff("cutoff", "Number", "", "0.01", "", "", "",false,false); parameters.push_back(pcutoff);
28 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
29 CommandParameter pmaxiter("maxiter", "Number", "", "1000", "", "", "",false,false); parameters.push_back(pmaxiter);
30 CommandParameter psigma("sigma", "Number", "", "60", "", "", "",false,false); parameters.push_back(psigma);
31 CommandParameter pmindelta("mindelta", "Number", "", "0.000001", "", "", "",false,false); parameters.push_back(pmindelta);
32 CommandParameter porder("order", "String", "", "", "", "", "",false,false); parameters.push_back(porder);
33 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
34 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
36 vector<string> myArray;
37 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
41 m->errorOut(e, "ShhherCommand", "setParameters");
45 //**********************************************************************************************************************
46 string ShhherCommand::getHelpString(){
48 string helpString = "";
49 helpString += "The shhh.seqs command reads a file containing flowgrams and creates a file of corrected sequences.\n";
53 m->errorOut(e, "ShhherCommand", "getHelpString");
57 //**********************************************************************************************************************
59 ShhherCommand::ShhherCommand(){
61 abort = true; calledHelp = true;
64 //initialize outputTypes
65 // vector<string> tempOutNames;
66 // outputTypes["pn.dist"] = tempOutNames;
70 m->errorOut(e, "ShhherCommand", "ShhherCommand");
75 //**********************************************************************************************************************
77 ShhherCommand::ShhherCommand(string option) {
81 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
82 MPI_Comm_size(MPI_COMM_WORLD, &ncpus);
88 abort = false; calledHelp = false;
91 //allow user to run help
92 if(option == "help") { help(); abort = true; calledHelp = true; }
93 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
96 vector<string> myArray = setParameters();
98 OptionParser parser(option);
99 map<string,string> parameters = parser.getParameters();
101 ValidParameters validParameter;
102 map<string,string>::iterator it;
104 //check to make sure all parameters are valid for command
105 for (it = parameters.begin(); it != parameters.end(); it++) {
106 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
109 //initialize outputTypes
110 vector<string> tempOutNames;
111 // outputTypes["pn.dist"] = tempOutNames;
112 // outputTypes["fasta"] = tempOutNames;
114 //if the user changes the input directory command factory will send this info to us in the output parameter
115 string inputDir = validParameter.validFile(parameters, "inputdir", false);
116 if (inputDir == "not found"){ inputDir = ""; }
119 it = parameters.find("flow");
120 //user has given a template file
121 if(it != parameters.end()){
122 path = m->hasPath(it->second);
123 //if the user has not given a path then, add inputdir. else leave path alone.
124 if (path == "") { parameters["flow"] = inputDir + it->second; }
127 it = parameters.find("lookup");
128 //user has given a template file
129 if(it != parameters.end()){
130 path = m->hasPath(it->second);
131 //if the user has not given a path then, add inputdir. else leave path alone.
132 if (path == "") { parameters["lookup"] = inputDir + it->second; }
135 it = parameters.find("file");
136 //user has given a template file
137 if(it != parameters.end()){
138 path = m->hasPath(it->second);
139 //if the user has not given a path then, add inputdir. else leave path alone.
140 if (path == "") { parameters["file"] = inputDir + it->second; }
145 //check for required parameters
146 flowFileName = validParameter.validFile(parameters, "flow", true);
147 flowFilesFileName = validParameter.validFile(parameters, "file", true);
148 if (flowFileName == "not found" && flowFilesFileName == "not found") {
149 m->mothurOut("values for either flow or file must be provided for the shhh.seqs command.");
150 m->mothurOutEndLine();
153 else if (flowFileName == "not open" || flowFilesFileName == "not open") { abort = true; }
155 if(flowFileName != "not found"){
156 compositeFASTAFileName = "";
157 compositeNamesFileName = "";
162 //flow.files = 9 character offset
163 compositeFASTAFileName = flowFilesFileName.substr(0, flowFilesFileName.length()-10) + "shhh.fasta";
164 m->openOutputFile(compositeFASTAFileName, temp);
167 compositeNamesFileName = flowFilesFileName.substr(0, flowFilesFileName.length()-10) + "shhh.names";
168 m->openOutputFile(compositeNamesFileName, temp);
172 //if the user changes the output directory command factory will send this info to us in the output parameter
173 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
175 outputDir += m->hasPath(flowFileName); //if user entered a file with a path then preserve it
179 //check for optional parameter and set defaults
180 // ...at some point should added some additional type checking...
182 temp = validParameter.validFile(parameters, "lookup", true);
183 if (temp == "not found") {
184 lookupFileName = "LookUp_Titanium.pat";
188 ableToOpen = m->openInputFile(lookupFileName, in, "noerror");
191 //if you can't open it, try input location
192 if (ableToOpen == 1) {
193 if (inputDir != "") { //default path is set
194 string tryPath = inputDir + lookupFileName;
195 m->mothurOut("Unable to open " + lookupFileName + ". Trying input directory " + tryPath); m->mothurOutEndLine();
197 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
199 lookupFileName = tryPath;
203 //if you can't open it, try default location
204 if (ableToOpen == 1) {
205 if (m->getDefaultPath() != "") { //default path is set
206 string tryPath = m->getDefaultPath() + m->getSimpleName(lookupFileName);
207 m->mothurOut("Unable to open " + lookupFileName + ". Trying default " + tryPath); m->mothurOutEndLine();
209 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
211 lookupFileName = tryPath;
215 //if you can't open it its not in current working directory or inputDir, try mothur excutable location
216 if (ableToOpen == 1) {
217 string exepath = m->argv;
218 string tempPath = exepath;
219 for (int i = 0; i < exepath.length(); i++) { tempPath[i] = tolower(exepath[i]); }
220 exepath = exepath.substr(0, (tempPath.find_last_of('m')));
222 string tryPath = m->getFullPathName(exepath) + m->getSimpleName(lookupFileName);
223 m->mothurOut("Unable to open " + lookupFileName + ". Trying mothur's executable location " + tryPath); m->mothurOutEndLine();
225 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
227 lookupFileName = tryPath;
230 if (ableToOpen == 1) { m->mothurOut("Unable to open " + lookupFileName + "."); m->mothurOutEndLine(); abort=true; }
232 else if(temp == "not open") {
234 lookupFileName = validParameter.validFile(parameters, "lookup", false);
236 //if you can't open it its not inputDir, try mothur excutable location
237 string exepath = m->argv;
238 string tempPath = exepath;
239 for (int i = 0; i < exepath.length(); i++) { tempPath[i] = tolower(exepath[i]); }
240 exepath = exepath.substr(0, (tempPath.find_last_of('m')));
242 string tryPath = m->getFullPathName(exepath) + lookupFileName;
243 m->mothurOut("Unable to open " + lookupFileName + ". Trying mothur's executable location " + tryPath); m->mothurOutEndLine();
245 int ableToOpen = m->openInputFile(tryPath, in2, "noerror");
247 lookupFileName = tryPath;
249 if (ableToOpen == 1) { m->mothurOut("Unable to open " + lookupFileName + "."); m->mothurOutEndLine(); abort=true; }
250 }else { lookupFileName = temp; }
252 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
253 m->setProcessors(temp);
254 convert(temp, processors);
256 temp = validParameter.validFile(parameters, "cutoff", false); if (temp == "not found"){ temp = "0.01"; }
257 convert(temp, cutoff);
259 temp = validParameter.validFile(parameters, "mindelta", false); if (temp == "not found"){ temp = "0.000001"; }
260 convert(temp, minDelta);
262 temp = validParameter.validFile(parameters, "maxiter", false); if (temp == "not found"){ temp = "1000"; }
263 convert(temp, maxIters);
265 temp = validParameter.validFile(parameters, "sigma", false);if (temp == "not found") { temp = "60"; }
266 convert(temp, sigma);
268 flowOrder = validParameter.validFile(parameters, "order", false);
269 if (flowOrder == "not found"){ flowOrder = "TACG"; }
270 else if(flowOrder.length() != 4){
271 m->mothurOut("The value of the order option must be four bases long\n");
281 catch(exception& e) {
282 m->errorOut(e, "ShhherCommand", "ShhherCommand");
286 //**********************************************************************************************************************
288 int ShhherCommand::execute(){
290 if (abort == true) { if (calledHelp) { return 0; } return 2; }
297 for(int i=1;i<ncpus;i++){
298 MPI_Send(&abort, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
300 if(abort == 1){ return 0; }
304 m->mothurOut("\nGetting preliminary data...\n");
305 getSingleLookUp(); if (m->control_pressed) { return 0; }
306 getJointLookUp(); if (m->control_pressed) { return 0; }
308 vector<string> flowFileVector;
309 if(flowFilesFileName != "not found"){
312 ifstream flowFilesFile;
313 m->openInputFile(flowFilesFileName, flowFilesFile);
314 while(flowFilesFile){
315 fName = m->getline(flowFilesFile);
316 flowFileVector.push_back(fName);
317 m->gobble(flowFilesFile);
321 flowFileVector.push_back(flowFileName);
323 int numFiles = flowFileVector.size();
325 for(int i=1;i<ncpus;i++){
326 MPI_Send(&numFiles, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
329 for(int i=0;i<numFiles;i++){
331 if (m->control_pressed) { break; }
333 double begClock = clock();
334 unsigned long long begTime = time(NULL);
336 flowFileName = flowFileVector[i];
338 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
339 m->mothurOut("Reading flowgrams...\n");
342 if (m->control_pressed) { break; }
344 m->mothurOut("Identifying unique flowgrams...\n");
347 if (m->control_pressed) { break; }
349 m->mothurOut("Calculating distances between flowgrams...\n");
351 strcpy(fileName, flowFileName.c_str());
353 for(int i=1;i<ncpus;i++){
354 MPI_Send(&fileName[0], 1024, MPI_CHAR, i, tag, MPI_COMM_WORLD);
356 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
357 MPI_Send(&numUniques, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
358 MPI_Send(&numFlowCells, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
359 MPI_Send(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, i, tag, MPI_COMM_WORLD);
360 MPI_Send(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
361 MPI_Send(&mapUniqueToSeq[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
362 MPI_Send(&mapSeqToUnique[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
363 MPI_Send(&lengths[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
364 MPI_Send(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
365 MPI_Send(&cutoff, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
368 string distFileName = flowDistMPI(0, int(sqrt(1.0/float(ncpus)) * numUniques));
370 if (m->control_pressed) { break; }
373 for(int i=1;i<ncpus;i++){
374 MPI_Recv(&done, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
376 m->appendFiles((distFileName + ".temp." + toString(i)), distFileName);
377 m->mothurRemove((distFileName + ".temp." + toString(i)));
380 string namesFileName = createNamesFile();
382 if (m->control_pressed) { break; }
384 m->mothurOut("\nClustering flowgrams...\n");
385 string listFileName = cluster(distFileName, namesFileName);
387 if (m->control_pressed) { break; }
389 getOTUData(listFileName);
391 m->mothurRemove(distFileName);
392 m->mothurRemove(namesFileName);
393 m->mothurRemove(listFileName);
395 if (m->control_pressed) { break; }
399 if (m->control_pressed) { break; }
401 for(int i=1;i<ncpus;i++){
402 MPI_Send(&numOTUs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
403 MPI_Send(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
404 MPI_Send(&uniqueFlowgrams[0], numFlowCells * numUniques, MPI_SHORT, i, tag, MPI_COMM_WORLD);
405 MPI_Send(&sigma, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
408 if (m->control_pressed) { break; }
413 int numOTUsOnCPU = numOTUs / ncpus;
414 int numSeqsOnCPU = numSeqs / ncpus;
415 m->mothurOut("\nDenoising flowgrams...\n");
416 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
418 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
420 double cycClock = clock();
421 unsigned long long cycTime = time(NULL);
424 if (m->control_pressed) { break; }
426 int total = singleTau.size();
427 for(int i=1;i<ncpus;i++){
428 MPI_Send(&total, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
429 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
430 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
432 MPI_Send(&singleTau[0], total, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
433 MPI_Send(&seqNumber[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
434 MPI_Send(&seqIndex[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
435 MPI_Send(&nSeqsPerOTU[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
436 MPI_Send(&cumNumSeqs[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
439 calcCentroidsDriver(0, numOTUsOnCPU);
441 for(int i=1;i<ncpus;i++){
442 int otuStart = i * numOTUs / ncpus;
443 int otuStop = (i + 1) * numOTUs / ncpus;
445 vector<int> tempCentroids(numOTUs, 0);
446 vector<short> tempChange(numOTUs, 0);
448 MPI_Recv(&tempCentroids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
449 MPI_Recv(&tempChange[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD, &status);
451 for(int j=otuStart;j<otuStop;j++){
452 centroids[j] = tempCentroids[j];
453 change[j] = tempChange[j];
457 maxDelta = getNewWeights(); if (m->control_pressed) { break; }
458 double nLL = getLikelihood(); if (m->control_pressed) { break; }
459 checkCentroids(); if (m->control_pressed) { break; }
461 for(int i=1;i<ncpus;i++){
462 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
463 MPI_Send(&weight[0], numOTUs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
464 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
467 calcNewDistancesParent(0, numSeqsOnCPU);
469 total = singleTau.size();
471 for(int i=1;i<ncpus;i++){
473 int seqStart = i * numSeqs / ncpus;
474 int seqStop = (i + 1) * numSeqs / ncpus;
476 MPI_Recv(&childTotal, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
478 vector<int> childSeqIndex(childTotal, 0);
479 vector<double> childSingleTau(childTotal, 0);
480 vector<double> childDist(numSeqs * numOTUs, 0);
481 vector<int> otuIndex(childTotal, 0);
483 MPI_Recv(&childSeqIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
484 MPI_Recv(&childSingleTau[0], childTotal, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
485 MPI_Recv(&childDist[0], numOTUs * numSeqs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
486 MPI_Recv(&otuIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
488 int oldTotal = total;
490 singleTau.resize(total, 0);
491 seqIndex.resize(total, 0);
492 seqNumber.resize(total, 0);
496 for(int j=oldTotal;j<total;j++){
497 int otuI = otuIndex[childIndex];
498 int seqI = childSeqIndex[childIndex];
500 singleTau[j] = childSingleTau[childIndex];
502 aaP[otuI][nSeqsPerOTU[otuI]] = j;
503 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
508 int index = seqStart * numOTUs;
509 for(int j=seqStart;j<seqStop;j++){
510 for(int k=0;k<numOTUs;k++){
511 dist[index] = childDist[index];
519 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
521 if((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
523 for(int i=1;i<ncpus;i++){
524 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
529 for(int i=1;i<ncpus;i++){
530 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD); //send kill command
536 if (m->control_pressed) { break; }
538 m->mothurOut("\nFinalizing...\n");
541 if (m->control_pressed) { break; }
545 vector<int> otuCounts(numOTUs, 0);
546 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
547 calcCentroidsDriver(0, numOTUs);
549 if (m->control_pressed) { break; }
551 writeQualities(otuCounts); if (m->control_pressed) { break; }
552 writeSequences(otuCounts); if (m->control_pressed) { break; }
553 writeNames(otuCounts); if (m->control_pressed) { break; }
554 writeClusters(otuCounts); if (m->control_pressed) { break; }
555 writeGroups(); if (m->control_pressed) { break; }
558 m->mothurOut("Total time to process " + toString(flowFileName) + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
564 MPI_Recv(&abort, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
565 if(abort){ return 0; }
568 MPI_Recv(&numFiles, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
570 for(int i=0;i<numFiles;i++){
572 if (m->control_pressed) { break; }
574 //Now into the pyrodist part
578 MPI_Recv(&fileName, 1024, MPI_CHAR, 0, tag, MPI_COMM_WORLD, &status);
579 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
580 MPI_Recv(&numUniques, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
581 MPI_Recv(&numFlowCells, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
583 flowDataIntI.resize(numSeqs * numFlowCells);
584 flowDataPrI.resize(numSeqs * numFlowCells);
585 mapUniqueToSeq.resize(numSeqs);
586 mapSeqToUnique.resize(numSeqs);
587 lengths.resize(numSeqs);
588 jointLookUp.resize(NUMBINS * NUMBINS);
590 MPI_Recv(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
591 MPI_Recv(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
592 MPI_Recv(&mapUniqueToSeq[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
593 MPI_Recv(&mapSeqToUnique[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
594 MPI_Recv(&lengths[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
595 MPI_Recv(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
596 MPI_Recv(&cutoff, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
598 flowFileName = string(fileName);
599 int flowDistStart = int(sqrt(float(pid)/float(ncpus)) * numUniques);
600 int flowDistEnd = int(sqrt(float(pid+1)/float(ncpus)) * numUniques);
602 string distanceStringChild = flowDistMPI(flowDistStart, flowDistEnd);
604 if (m->control_pressed) { break; }
607 MPI_Send(&done, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
609 //Now into the pyronoise part
610 MPI_Recv(&numOTUs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
612 singleLookUp.resize(HOMOPS * NUMBINS);
613 uniqueFlowgrams.resize(numUniques * numFlowCells);
614 weight.resize(numOTUs);
615 centroids.resize(numOTUs);
616 change.resize(numOTUs);
617 dist.assign(numOTUs * numSeqs, 0);
618 nSeqsPerOTU.resize(numOTUs);
619 cumNumSeqs.resize(numOTUs);
621 MPI_Recv(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
622 MPI_Recv(&uniqueFlowgrams[0], uniqueFlowgrams.size(), MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
623 MPI_Recv(&sigma, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
625 int startOTU = pid * numOTUs / ncpus;
626 int endOTU = (pid + 1) * numOTUs / ncpus;
628 int startSeq = pid * numSeqs / ncpus;
629 int endSeq = (pid + 1) * numSeqs /ncpus;
635 if (m->control_pressed) { break; }
637 MPI_Recv(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
638 singleTau.assign(total, 0.0000);
639 seqNumber.assign(total, 0);
640 seqIndex.assign(total, 0);
642 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
643 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
644 MPI_Recv(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
645 MPI_Recv(&seqNumber[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
646 MPI_Recv(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
647 MPI_Recv(&nSeqsPerOTU[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
648 MPI_Recv(&cumNumSeqs[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
650 calcCentroidsDriver(startOTU, endOTU);
652 MPI_Send(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD);
653 MPI_Send(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD);
655 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
656 MPI_Recv(&weight[0], numOTUs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
657 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
659 vector<int> otuIndex(total, 0);
660 calcNewDistancesChildMPI(startSeq, endSeq, otuIndex);
661 total = otuIndex.size();
663 MPI_Send(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
664 MPI_Send(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
665 MPI_Send(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
666 MPI_Send(&dist[0], numOTUs * numSeqs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
667 MPI_Send(&otuIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
669 MPI_Recv(&live, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
674 if (m->control_pressed) { for (int i = 0; i < outputNames.size(); i++) { m->mothurRemove(outputNames[i]); } return 0; }
676 MPI_Barrier(MPI_COMM_WORLD);
679 if(compositeFASTAFileName != ""){
680 outputNames.push_back(compositeFASTAFileName);
681 outputNames.push_back(compositeNamesFileName);
684 m->mothurOutEndLine();
685 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
686 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
687 m->mothurOutEndLine();
692 catch(exception& e) {
693 m->errorOut(e, "ShhherCommand", "execute");
698 /**************************************************************************************************/
700 string ShhherCommand::flowDistMPI(int startSeq, int stopSeq){
702 ostringstream outStream;
703 outStream.setf(ios::fixed, ios::floatfield);
704 outStream.setf(ios::dec, ios::basefield);
705 outStream.setf(ios::showpoint);
706 outStream.precision(6);
708 int begTime = time(NULL);
709 double begClock = clock();
711 for(int i=startSeq;i<stopSeq;i++){
713 if (m->control_pressed) { break; }
715 for(int j=0;j<i;j++){
716 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
718 if(flowDistance < 1e-6){
719 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
721 else if(flowDistance <= cutoff){
722 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
726 m->mothurOut(toString(i) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
730 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.dist";
731 if(pid != 0){ fDistFileName += ".temp." + toString(pid); }
733 if (m->control_pressed) { return fDistFileName; }
735 m->mothurOut(toString(stopSeq) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
737 ofstream distFile(fDistFileName.c_str());
738 distFile << outStream.str();
741 return fDistFileName;
743 catch(exception& e) {
744 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
750 //**********************************************************************************************************************
752 int ShhherCommand::execute(){
754 if (abort == true) { return 0; }
756 getSingleLookUp(); if (m->control_pressed) { return 0; }
757 getJointLookUp(); if (m->control_pressed) { return 0; }
760 vector<string> flowFileVector;
761 if(flowFilesFileName != "not found"){
764 ifstream flowFilesFile;
765 m->openInputFile(flowFilesFileName, flowFilesFile);
766 while(flowFilesFile){
767 fName = m->getline(flowFilesFile);
768 flowFileVector.push_back(fName);
769 m->gobble(flowFilesFile);
773 flowFileVector.push_back(flowFileName);
775 int numFiles = flowFileVector.size();
778 for(int i=0;i<numFiles;i++){
780 if (m->control_pressed) { break; }
782 flowFileName = flowFileVector[i];
784 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
785 m->mothurOut("Reading flowgrams...\n");
788 if (m->control_pressed) { break; }
790 m->mothurOut("Identifying unique flowgrams...\n");
793 if (m->control_pressed) { break; }
795 m->mothurOut("Calculating distances between flowgrams...\n");
796 string distFileName = createDistFile(processors);
797 string namesFileName = createNamesFile();
799 if (m->control_pressed) { break; }
801 m->mothurOut("\nClustering flowgrams...\n");
802 string listFileName = cluster(distFileName, namesFileName);
804 if (m->control_pressed) { break; }
806 getOTUData(listFileName);
808 if (m->control_pressed) { break; }
810 m->mothurRemove(distFileName);
811 m->mothurRemove(namesFileName);
812 m->mothurRemove(listFileName);
816 if (m->control_pressed) { break; }
821 double begClock = clock();
822 unsigned long long begTime = time(NULL);
824 m->mothurOut("\nDenoising flowgrams...\n");
825 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
827 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
829 if (m->control_pressed) { break; }
831 double cycClock = clock();
832 unsigned long long cycTime = time(NULL);
835 if (m->control_pressed) { break; }
839 if (m->control_pressed) { break; }
841 maxDelta = getNewWeights(); if (m->control_pressed) { break; }
842 double nLL = getLikelihood(); if (m->control_pressed) { break; }
845 if (m->control_pressed) { break; }
849 if (m->control_pressed) { break; }
853 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
857 if (m->control_pressed) { break; }
859 m->mothurOut("\nFinalizing...\n");
862 if (m->control_pressed) { break; }
866 if (m->control_pressed) { break; }
868 vector<int> otuCounts(numOTUs, 0);
869 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
871 calcCentroidsDriver(0, numOTUs); if (m->control_pressed) { break; }
872 writeQualities(otuCounts); if (m->control_pressed) { break; }
873 writeSequences(otuCounts); if (m->control_pressed) { break; }
874 writeNames(otuCounts); if (m->control_pressed) { break; }
875 writeClusters(otuCounts); if (m->control_pressed) { break; }
876 writeGroups(); if (m->control_pressed) { break; }
878 m->mothurOut("Total time to process " + flowFileName + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
881 if (m->control_pressed) { for (int i = 0; i < outputNames.size(); i++) { m->mothurRemove(outputNames[i]); } return 0; }
884 if(compositeFASTAFileName != ""){
885 outputNames.push_back(compositeFASTAFileName);
886 outputNames.push_back(compositeNamesFileName);
889 m->mothurOutEndLine();
890 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
891 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
892 m->mothurOutEndLine();
896 catch(exception& e) {
897 m->errorOut(e, "ShhherCommand", "execute");
902 /**************************************************************************************************/
904 void ShhherCommand::getFlowData(){
907 m->openInputFile(flowFileName, flowFile);
910 seqNameVector.clear();
912 flowDataIntI.clear();
916 int currentNumFlowCells;
920 flowFile >> numFlowCells;
921 int index = 0;//pcluster
922 while(!flowFile.eof()){
924 if (m->control_pressed) { break; }
926 flowFile >> seqName >> currentNumFlowCells;
927 lengths.push_back(currentNumFlowCells);
929 seqNameVector.push_back(seqName);
930 nameMap[seqName] = index++;//pcluster
932 for(int i=0;i<numFlowCells;i++){
933 flowFile >> intensity;
934 if(intensity > 9.99) { intensity = 9.99; }
935 int intI = int(100 * intensity + 0.0001);
936 flowDataIntI.push_back(intI);
942 numSeqs = seqNameVector.size();
944 for(int i=0;i<numSeqs;i++){
946 if (m->control_pressed) { break; }
948 int iNumFlowCells = i * numFlowCells;
949 for(int j=lengths[i];j<numFlowCells;j++){
950 flowDataIntI[iNumFlowCells + j] = 0;
955 catch(exception& e) {
956 m->errorOut(e, "ShhherCommand", "getFlowData");
961 /**************************************************************************************************/
963 void ShhherCommand::getSingleLookUp(){
965 // these are the -log probabilities that a signal corresponds to a particular homopolymer length
966 singleLookUp.assign(HOMOPS * NUMBINS, 0);
970 m->openInputFile(lookupFileName, lookUpFile);
972 for(int i=0;i<HOMOPS;i++){
974 if (m->control_pressed) { break; }
977 lookUpFile >> logFracFreq;
979 for(int j=0;j<NUMBINS;j++) {
980 lookUpFile >> singleLookUp[index];
986 catch(exception& e) {
987 m->errorOut(e, "ShhherCommand", "getSingleLookUp");
992 /**************************************************************************************************/
994 void ShhherCommand::getJointLookUp(){
997 // the most likely joint probability (-log) that two intenities have the same polymer length
998 jointLookUp.resize(NUMBINS * NUMBINS, 0);
1000 for(int i=0;i<NUMBINS;i++){
1002 if (m->control_pressed) { break; }
1004 for(int j=0;j<NUMBINS;j++){
1006 double minSum = 100000000;
1008 for(int k=0;k<HOMOPS;k++){
1009 double sum = singleLookUp[k * NUMBINS + i] + singleLookUp[k * NUMBINS + j];
1011 if(sum < minSum) { minSum = sum; }
1013 jointLookUp[i * NUMBINS + j] = minSum;
1017 catch(exception& e) {
1018 m->errorOut(e, "ShhherCommand", "getJointLookUp");
1023 /**************************************************************************************************/
1025 double ShhherCommand::getProbIntensity(int intIntensity){
1028 double minNegLogProb = 100000000;
1031 for(int i=0;i<HOMOPS;i++){//loop signal strength
1033 if (m->control_pressed) { break; }
1035 float negLogProb = singleLookUp[i * NUMBINS + intIntensity];
1036 if(negLogProb < minNegLogProb) { minNegLogProb = negLogProb; }
1039 return minNegLogProb;
1041 catch(exception& e) {
1042 m->errorOut(e, "ShhherCommand", "getProbIntensity");
1047 /**************************************************************************************************/
1049 void ShhherCommand::getUniques(){
1054 uniqueFlowgrams.assign(numFlowCells * numSeqs, -1);
1055 uniqueCount.assign(numSeqs, 0); // anWeights
1056 uniqueLengths.assign(numSeqs, 0);
1057 mapSeqToUnique.assign(numSeqs, -1);
1058 mapUniqueToSeq.assign(numSeqs, -1);
1060 vector<short> uniqueFlowDataIntI(numFlowCells * numSeqs, -1);
1062 for(int i=0;i<numSeqs;i++){
1064 if (m->control_pressed) { break; }
1068 vector<short> current(numFlowCells);
1069 for(int j=0;j<numFlowCells;j++){
1070 current[j] = short(((flowDataIntI[i * numFlowCells + j] + 50.0)/100.0));
1073 for(int j=0;j<numUniques;j++){
1074 int offset = j * numFlowCells;
1078 if(lengths[i] < uniqueLengths[j]) { shorterLength = lengths[i]; }
1079 else { shorterLength = uniqueLengths[j]; }
1081 for(int k=0;k<shorterLength;k++){
1082 if(current[k] != uniqueFlowgrams[offset + k]){
1089 mapSeqToUnique[i] = j;
1092 if(lengths[i] > uniqueLengths[j]) { uniqueLengths[j] = lengths[i]; }
1098 if(index == numUniques){
1099 uniqueLengths[numUniques] = lengths[i];
1100 uniqueCount[numUniques] = 1;
1101 mapSeqToUnique[i] = numUniques;//anMap
1102 mapUniqueToSeq[numUniques] = i;//anF
1104 for(int k=0;k<numFlowCells;k++){
1105 uniqueFlowgrams[numUniques * numFlowCells + k] = current[k];
1106 uniqueFlowDataIntI[numUniques * numFlowCells + k] = flowDataIntI[i * numFlowCells + k];
1112 uniqueFlowDataIntI.resize(numFlowCells * numUniques);
1113 uniqueLengths.resize(numUniques);
1115 flowDataPrI.resize(numSeqs * numFlowCells, 0);
1116 for(int i=0;i<flowDataPrI.size();i++) { if (m->control_pressed) { break; } flowDataPrI[i] = getProbIntensity(flowDataIntI[i]); }
1118 catch(exception& e) {
1119 m->errorOut(e, "ShhherCommand", "getUniques");
1124 /**************************************************************************************************/
1126 float ShhherCommand::calcPairwiseDist(int seqA, int seqB){
1128 int minLength = lengths[mapSeqToUnique[seqA]];
1129 if(lengths[seqB] < minLength){ minLength = lengths[mapSeqToUnique[seqB]]; }
1131 int ANumFlowCells = seqA * numFlowCells;
1132 int BNumFlowCells = seqB * numFlowCells;
1136 for(int i=0;i<minLength;i++){
1138 if (m->control_pressed) { break; }
1140 int flowAIntI = flowDataIntI[ANumFlowCells + i];
1141 float flowAPrI = flowDataPrI[ANumFlowCells + i];
1143 int flowBIntI = flowDataIntI[BNumFlowCells + i];
1144 float flowBPrI = flowDataPrI[BNumFlowCells + i];
1145 dist += jointLookUp[flowAIntI * NUMBINS + flowBIntI] - flowAPrI - flowBPrI;
1148 dist /= (float) minLength;
1151 catch(exception& e) {
1152 m->errorOut(e, "ShhherCommand", "calcPairwiseDist");
1157 /**************************************************************************************************/
1159 void ShhherCommand::flowDistParentFork(string distFileName, int startSeq, int stopSeq){
1162 ostringstream outStream;
1163 outStream.setf(ios::fixed, ios::floatfield);
1164 outStream.setf(ios::dec, ios::basefield);
1165 outStream.setf(ios::showpoint);
1166 outStream.precision(6);
1168 int begTime = time(NULL);
1169 double begClock = clock();
1171 for(int i=startSeq;i<stopSeq;i++){
1173 if (m->control_pressed) { break; }
1175 for(int j=0;j<i;j++){
1176 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
1178 if(flowDistance < 1e-6){
1179 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
1181 else if(flowDistance <= cutoff){
1182 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
1186 m->mothurOut(toString(i) + "\t" + toString(time(NULL) - begTime));
1187 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1188 m->mothurOutEndLine();
1192 ofstream distFile(distFileName.c_str());
1193 distFile << outStream.str();
1196 if (m->control_pressed) {}
1198 m->mothurOut(toString(stopSeq-1) + "\t" + toString(time(NULL) - begTime));
1199 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1200 m->mothurOutEndLine();
1203 catch(exception& e) {
1204 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
1209 /**************************************************************************************************/
1211 string ShhherCommand::createDistFile(int processors){
1213 //////////////////////// until I figure out the shared memory issue //////////////////////
1214 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1218 //////////////////////// until I figure out the shared memory issue //////////////////////
1220 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.dist";
1222 unsigned long long begTime = time(NULL);
1223 double begClock = clock();
1225 if (numSeqs < processors){ processors = 1; }
1227 if(processors == 1) { flowDistParentFork(fDistFileName, 0, numUniques); }
1229 else{ //you have multiple processors
1231 vector<int> start(processors, 0);
1232 vector<int> end(processors, 0);
1235 vector<int> processIDs;
1237 for (int i = 0; i < processors; i++) {
1238 start[i] = int(sqrt(float(i)/float(processors)) * numUniques);
1239 end[i] = int(sqrt(float(i+1)/float(processors)) * numUniques);
1242 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1244 //loop through and create all the processes you want
1245 while (process != processors) {
1249 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1251 }else if (pid == 0){
1252 flowDistParentFork(fDistFileName + toString(getpid()) + ".temp", start[process], end[process]);
1255 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1257 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1262 //parent does its part
1263 flowDistParentFork(fDistFileName, start[0], end[0]);
1265 //force parent to wait until all the processes are done
1266 for (int i=0;i<processIDs.size();i++) {
1267 int temp = processIDs[i];
1271 //////////////////////////////////////////////////////////////////////////////////////////////////////
1272 //Windows version shared memory, so be careful when passing variables through the flowDistParentForkData struct.
1273 //Above fork() will clone, so memory is separate, but that's not the case with windows,
1274 //////////////////////////////////////////////////////////////////////////////////////////////////////
1276 vector<flowDistParentForkData*> pDataArray;
1277 DWORD dwThreadIdArray[processors-1];
1278 HANDLE hThreadArray[processors-1];
1280 //Create processor worker threads.
1281 for(int i = 0; i < processors-1; i++){
1282 // Allocate memory for thread data.
1283 string extension = extension = toString(i) + ".temp";
1285 flowDistParentForkData* tempdist = new flowDistParentForkData((fDistFileName + extension), mapUniqueToSeq, mapSeqToUnique, lengths, flowDataIntI, flowDataPrI, jointLookUp, m, start[i+1], end[i+1], numFlowCells, cutoff, i);
1286 pDataArray.push_back(tempdist);
1287 processIDs.push_back(i);
1289 //MySeqSumThreadFunction is in header. It must be global or static to work with the threads.
1290 //default security attributes, thread function name, argument to thread function, use default creation flags, returns the thread identifier
1291 hThreadArray[i] = CreateThread(NULL, 0, MyflowDistParentForkThreadFunction, pDataArray[i], 0, &dwThreadIdArray[i]);
1294 //parent does its part
1295 flowDistParentFork(fDistFileName, start[0], end[0]);
1297 //Wait until all threads have terminated.
1298 WaitForMultipleObjects(processors-1, hThreadArray, TRUE, INFINITE);
1300 //Close all thread handles and free memory allocations.
1301 for(int i=0; i < pDataArray.size(); i++){
1302 CloseHandle(hThreadArray[i]);
1303 delete pDataArray[i];
1308 //append and remove temp files
1309 for (int i=0;i<processIDs.size();i++) {
1310 m->appendFiles((fDistFileName + toString(processIDs[i]) + ".temp"), fDistFileName);
1311 m->mothurRemove((fDistFileName + toString(processIDs[i]) + ".temp"));
1316 m->mothurOutEndLine();
1317 m->mothurOut("Total time: " + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/CLOCKS_PER_SEC) + '\n');
1319 return fDistFileName;
1321 catch(exception& e) {
1322 m->errorOut(e, "ShhherCommand", "createDistFile");
1328 /**************************************************************************************************/
1330 string ShhherCommand::createNamesFile(){
1333 vector<string> duplicateNames(numUniques, "");
1334 for(int i=0;i<numSeqs;i++){
1335 duplicateNames[mapSeqToUnique[i]] += seqNameVector[i] + ',';
1338 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.names";
1341 m->openOutputFile(nameFileName, nameFile);
1343 for(int i=0;i<numUniques;i++){
1345 if (m->control_pressed) { break; }
1347 // nameFile << seqNameVector[mapUniqueToSeq[i]] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1348 nameFile << mapUniqueToSeq[i] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1352 return nameFileName;
1354 catch(exception& e) {
1355 m->errorOut(e, "ShhherCommand", "createNamesFile");
1360 //**********************************************************************************************************************
1362 string ShhherCommand::cluster(string distFileName, string namesFileName){
1365 ReadMatrix* read = new ReadColumnMatrix(distFileName);
1366 read->setCutoff(cutoff);
1368 NameAssignment* clusterNameMap = new NameAssignment(namesFileName);
1369 clusterNameMap->readMap();
1370 read->read(clusterNameMap);
1372 ListVector* list = read->getListVector();
1373 SparseMatrix* matrix = read->getMatrix();
1376 delete clusterNameMap;
1378 RAbundVector* rabund = new RAbundVector(list->getRAbundVector());
1380 Cluster* cluster = new CompleteLinkage(rabund, list, matrix, cutoff, "furthest");
1381 string tag = cluster->getTag();
1383 double clusterCutoff = cutoff;
1384 while (matrix->getSmallDist() <= clusterCutoff && matrix->getNNodes() > 0){
1386 if (m->control_pressed) { break; }
1388 cluster->update(clusterCutoff);
1391 list->setLabel(toString(cutoff));
1393 string listFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.list";
1395 m->openOutputFile(listFileName, listFile);
1396 list->print(listFile);
1399 delete matrix; delete cluster; delete rabund; delete list;
1401 return listFileName;
1403 catch(exception& e) {
1404 m->errorOut(e, "ShhherCommand", "cluster");
1409 /**************************************************************************************************/
1411 void ShhherCommand::getOTUData(string listFileName){
1415 m->openInputFile(listFileName, listFile);
1418 listFile >> label >> numOTUs;
1420 otuData.assign(numSeqs, 0);
1421 cumNumSeqs.assign(numOTUs, 0);
1422 nSeqsPerOTU.assign(numOTUs, 0);
1423 aaP.clear();aaP.resize(numOTUs);
1429 string singleOTU = "";
1431 for(int i=0;i<numOTUs;i++){
1433 if (m->control_pressed) { break; }
1435 listFile >> singleOTU;
1437 istringstream otuString(singleOTU);
1441 string seqName = "";
1443 for(int j=0;j<singleOTU.length();j++){
1444 char letter = otuString.get();
1450 map<string,int>::iterator nmIt = nameMap.find(seqName);
1451 int index = nmIt->second;
1453 nameMap.erase(nmIt);
1457 aaP[i].push_back(index);
1462 map<string,int>::iterator nmIt = nameMap.find(seqName);
1464 int index = nmIt->second;
1465 nameMap.erase(nmIt);
1469 aaP[i].push_back(index);
1474 sort(aaP[i].begin(), aaP[i].end());
1475 for(int j=0;j<nSeqsPerOTU[i];j++){
1476 seqNumber.push_back(aaP[i][j]);
1478 for(int j=nSeqsPerOTU[i];j<numSeqs;j++){
1479 aaP[i].push_back(0);
1485 for(int i=1;i<numOTUs;i++){
1486 cumNumSeqs[i] = cumNumSeqs[i-1] + nSeqsPerOTU[i-1];
1489 seqIndex = seqNumber;
1494 catch(exception& e) {
1495 m->errorOut(e, "ShhherCommand", "getOTUData");
1500 /**************************************************************************************************/
1502 void ShhherCommand::initPyroCluster(){
1505 if (numOTUs < processors) { processors = 1; }
1507 dist.assign(numSeqs * numOTUs, 0);
1508 change.assign(numOTUs, 1);
1509 centroids.assign(numOTUs, -1);
1510 weight.assign(numOTUs, 0);
1511 singleTau.assign(numSeqs, 1.0);
1513 nSeqsBreaks.assign(processors+1, 0);
1514 nOTUsBreaks.assign(processors+1, 0);
1517 for(int i=0;i<processors;i++){
1518 nSeqsBreaks[i+1] = nSeqsBreaks[i] + (int)((double) numSeqs / (double) processors);
1519 nOTUsBreaks[i+1] = nOTUsBreaks[i] + (int)((double) numOTUs / (double) processors);
1521 nSeqsBreaks[processors] = numSeqs;
1522 nOTUsBreaks[processors] = numOTUs;
1524 catch(exception& e) {
1525 m->errorOut(e, "ShhherCommand", "initPyroCluster");
1530 /**************************************************************************************************/
1532 void ShhherCommand::fill(){
1535 for(int i=0;i<numOTUs;i++){
1537 if (m->control_pressed) { break; }
1539 cumNumSeqs[i] = index;
1540 for(int j=0;j<nSeqsPerOTU[i];j++){
1541 seqNumber[index] = aaP[i][j];
1542 seqIndex[index] = aaI[i][j];
1548 catch(exception& e) {
1549 m->errorOut(e, "ShhherCommand", "fill");
1554 /**************************************************************************************************/
1556 void ShhherCommand::calcCentroids(){
1559 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1561 if(processors == 1) {
1562 calcCentroidsDriver(0, numOTUs);
1564 else{ //you have multiple processors
1565 if (numOTUs < processors){ processors = 1; }
1568 vector<int> processIDs;
1570 //loop through and create all the processes you want
1571 while (process != processors) {
1575 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1577 }else if (pid == 0){
1578 calcCentroidsDriver(nOTUsBreaks[process], nOTUsBreaks[process+1]);
1581 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1583 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1588 //parent does its part
1589 calcCentroidsDriver(nOTUsBreaks[0], nOTUsBreaks[1]);
1591 //force parent to wait until all the processes are done
1592 for (int i=0;i<processIDs.size();i++) {
1593 int temp = processIDs[i];
1599 calcCentroidsDriver(0, numOTUs);
1602 catch(exception& e) {
1603 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1608 /**************************************************************************************************/
1610 void ShhherCommand::calcCentroidsDriver(int start, int finish){
1612 //this function gets the most likely homopolymer length at a flow position for a group of sequences
1617 for(int i=start;i<finish;i++){
1619 if (m->control_pressed) { break; }
1623 int minFlowGram = 100000000;
1624 double minFlowValue = 1e8;
1625 change[i] = 0; //FALSE
1627 for(int j=0;j<nSeqsPerOTU[i];j++){
1628 count += singleTau[seqNumber[cumNumSeqs[i] + j]];
1631 if(nSeqsPerOTU[i] > 0 && count > MIN_COUNT){
1632 vector<double> adF(nSeqsPerOTU[i]);
1633 vector<int> anL(nSeqsPerOTU[i]);
1635 for(int j=0;j<nSeqsPerOTU[i];j++){
1636 int index = cumNumSeqs[i] + j;
1637 int nI = seqIndex[index];
1638 int nIU = mapSeqToUnique[nI];
1641 for(k=0;k<position;k++){
1647 anL[position] = nIU;
1648 adF[position] = 0.0000;
1653 for(int j=0;j<nSeqsPerOTU[i];j++){
1654 int index = cumNumSeqs[i] + j;
1655 int nI = seqIndex[index];
1657 double tauValue = singleTau[seqNumber[index]];
1659 for(int k=0;k<position;k++){
1660 double dist = getDistToCentroid(anL[k], nI, lengths[nI]);
1661 adF[k] += dist * tauValue;
1665 for(int j=0;j<position;j++){
1666 if(adF[j] < minFlowValue){
1668 minFlowValue = adF[j];
1672 if(centroids[i] != anL[minFlowGram]){
1674 centroids[i] = anL[minFlowGram];
1677 else if(centroids[i] != -1){
1683 catch(exception& e) {
1684 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1689 /**************************************************************************************************/
1691 double ShhherCommand::getDistToCentroid(int cent, int flow, int length){
1694 int flowAValue = cent * numFlowCells;
1695 int flowBValue = flow * numFlowCells;
1699 for(int i=0;i<length;i++){
1700 dist += singleLookUp[uniqueFlowgrams[flowAValue] * NUMBINS + flowDataIntI[flowBValue]];
1705 return dist / (double)length;
1707 catch(exception& e) {
1708 m->errorOut(e, "ShhherCommand", "getDistToCentroid");
1713 /**************************************************************************************************/
1715 double ShhherCommand::getNewWeights(){
1718 double maxChange = 0;
1720 for(int i=0;i<numOTUs;i++){
1722 if (m->control_pressed) { break; }
1724 double difference = weight[i];
1727 for(int j=0;j<nSeqsPerOTU[i];j++){
1728 int index = cumNumSeqs[i] + j;
1729 double tauValue = singleTau[seqNumber[index]];
1730 weight[i] += tauValue;
1733 difference = fabs(weight[i] - difference);
1734 if(difference > maxChange){ maxChange = difference; }
1738 catch(exception& e) {
1739 m->errorOut(e, "ShhherCommand", "getNewWeights");
1744 /**************************************************************************************************/
1746 double ShhherCommand::getLikelihood(){
1750 vector<long double> P(numSeqs, 0);
1753 for(int i=0;i<numOTUs;i++){
1754 if(weight[i] > MIN_WEIGHT){
1760 for(int i=0;i<numOTUs;i++){
1762 if (m->control_pressed) { break; }
1764 for(int j=0;j<nSeqsPerOTU[i];j++){
1765 int index = cumNumSeqs[i] + j;
1766 int nI = seqIndex[index];
1767 double singleDist = dist[seqNumber[index]];
1769 P[nI] += weight[i] * exp(-singleDist * sigma);
1773 for(int i=0;i<numSeqs;i++){
1774 if(P[i] == 0){ P[i] = DBL_EPSILON; }
1779 nLL = nLL -(double)numSeqs * log(sigma);
1783 catch(exception& e) {
1784 m->errorOut(e, "ShhherCommand", "getNewWeights");
1789 /**************************************************************************************************/
1791 void ShhherCommand::checkCentroids(){
1793 vector<int> unique(numOTUs, 1);
1795 for(int i=0;i<numOTUs;i++){
1796 if(centroids[i] == -1 || weight[i] < MIN_WEIGHT){
1801 for(int i=0;i<numOTUs;i++){
1803 if (m->control_pressed) { break; }
1806 for(int j=i+1;j<numOTUs;j++){
1809 if(centroids[j] == centroids[i]){
1813 weight[i] += weight[j];
1821 catch(exception& e) {
1822 m->errorOut(e, "ShhherCommand", "checkCentroids");
1827 /**************************************************************************************************/
1829 void ShhherCommand::calcNewDistances(){
1832 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1834 if(processors == 1) {
1835 calcNewDistancesParent(0, numSeqs);
1837 else{ //you have multiple processors
1838 if (numSeqs < processors){ processors = 1; }
1840 vector<vector<int> > child_otuIndex(processors);
1841 vector<vector<int> > child_seqIndex(processors);
1842 vector<vector<double> > child_singleTau(processors);
1843 vector<int> totals(processors);
1846 vector<int> processIDs;
1848 //loop through and create all the processes you want
1849 while (process != processors) {
1853 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1855 }else if (pid == 0){
1856 calcNewDistancesChild(nSeqsBreaks[process], nSeqsBreaks[process+1], child_otuIndex[process], child_seqIndex[process], child_singleTau[process]);
1857 totals[process] = child_otuIndex[process].size();
1861 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1863 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1868 //parent does its part
1869 calcNewDistancesParent(nSeqsBreaks[0], nSeqsBreaks[1]);
1870 int total = seqIndex.size();
1872 //force parent to wait until all the processes are done
1873 for (int i=0;i<processIDs.size();i++) {
1874 int temp = processIDs[i];
1878 for(int i=1;i<processors;i++){
1879 int oldTotal = total;
1882 singleTau.resize(total, 0);
1883 seqIndex.resize(total, 0);
1884 seqNumber.resize(total, 0);
1888 for(int j=oldTotal;j<total;j++){
1889 int otuI = child_otuIndex[i][childIndex];
1890 int seqI = child_seqIndex[i][childIndex];
1892 singleTau[j] = child_singleTau[i][childIndex];
1893 aaP[otuI][nSeqsPerOTU[otuI]] = j;
1894 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
1895 nSeqsPerOTU[otuI]++;
1902 calcNewDistancesParent(0, numSeqs);
1905 catch(exception& e) {
1906 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1911 /**************************************************************************************************/
1913 void ShhherCommand::calcNewDistancesChildMPI(int startSeq, int stopSeq, vector<int>& otuIndex){
1916 vector<double> newTau(numOTUs,0);
1917 vector<double> norms(numSeqs, 0);
1922 for(int i=startSeq;i<stopSeq;i++){
1924 if (m->control_pressed) { break; }
1926 double offset = 1e8;
1927 int indexOffset = i * numOTUs;
1929 for(int j=0;j<numOTUs;j++){
1931 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1932 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1934 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1935 offset = dist[indexOffset + j];
1939 for(int j=0;j<numOTUs;j++){
1940 if(weight[j] > MIN_WEIGHT){
1941 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1942 norms[i] += newTau[j];
1949 for(int j=0;j<numOTUs;j++){
1951 newTau[j] /= norms[i];
1953 if(newTau[j] > MIN_TAU){
1954 otuIndex.push_back(j);
1955 seqIndex.push_back(i);
1956 singleTau.push_back(newTau[j]);
1962 catch(exception& e) {
1963 m->errorOut(e, "ShhherCommand", "calcNewDistancesChildMPI");
1968 /**************************************************************************************************/
1970 void ShhherCommand::calcNewDistancesChild(int startSeq, int stopSeq, vector<int>& child_otuIndex, vector<int>& child_seqIndex, vector<double>& child_singleTau){
1973 vector<double> newTau(numOTUs,0);
1974 vector<double> norms(numSeqs, 0);
1975 child_otuIndex.resize(0);
1976 child_seqIndex.resize(0);
1977 child_singleTau.resize(0);
1979 for(int i=startSeq;i<stopSeq;i++){
1981 if (m->control_pressed) { break; }
1983 double offset = 1e8;
1984 int indexOffset = i * numOTUs;
1987 for(int j=0;j<numOTUs;j++){
1988 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1989 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1992 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1993 offset = dist[indexOffset + j];
1997 for(int j=0;j<numOTUs;j++){
1998 if(weight[j] > MIN_WEIGHT){
1999 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
2000 norms[i] += newTau[j];
2007 for(int j=0;j<numOTUs;j++){
2008 newTau[j] /= norms[i];
2010 if(newTau[j] > MIN_TAU){
2011 child_otuIndex.push_back(j);
2012 child_seqIndex.push_back(i);
2013 child_singleTau.push_back(newTau[j]);
2018 catch(exception& e) {
2019 m->errorOut(e, "ShhherCommand", "calcNewDistancesChild");
2024 /**************************************************************************************************/
2026 void ShhherCommand::calcNewDistancesParent(int startSeq, int stopSeq){
2031 vector<double> newTau(numOTUs,0);
2032 vector<double> norms(numSeqs, 0);
2033 nSeqsPerOTU.assign(numOTUs, 0);
2035 for(int i=startSeq;i<stopSeq;i++){
2037 if (m->control_pressed) { break; }
2039 int indexOffset = i * numOTUs;
2041 double offset = 1e8;
2043 for(int j=0;j<numOTUs;j++){
2045 if(weight[j] > MIN_WEIGHT && change[j] == 1){
2046 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
2049 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
2050 offset = dist[indexOffset + j];
2054 for(int j=0;j<numOTUs;j++){
2055 if(weight[j] > MIN_WEIGHT){
2056 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
2057 norms[i] += newTau[j];
2064 for(int j=0;j<numOTUs;j++){
2065 newTau[j] /= norms[i];
2068 for(int j=0;j<numOTUs;j++){
2069 if(newTau[j] > MIN_TAU){
2071 int oldTotal = total;
2075 singleTau.resize(total, 0);
2076 seqNumber.resize(total, 0);
2077 seqIndex.resize(total, 0);
2079 singleTau[oldTotal] = newTau[j];
2081 aaP[j][nSeqsPerOTU[j]] = oldTotal;
2082 aaI[j][nSeqsPerOTU[j]] = i;
2090 catch(exception& e) {
2091 m->errorOut(e, "ShhherCommand", "calcNewDistancesParent");
2096 /**************************************************************************************************/
2098 void ShhherCommand::setOTUs(){
2101 vector<double> bigTauMatrix(numOTUs * numSeqs, 0.0000);
2103 for(int i=0;i<numOTUs;i++){
2105 if (m->control_pressed) { break; }
2107 for(int j=0;j<nSeqsPerOTU[i];j++){
2108 int index = cumNumSeqs[i] + j;
2109 double tauValue = singleTau[seqNumber[index]];
2110 int sIndex = seqIndex[index];
2111 bigTauMatrix[sIndex * numOTUs + i] = tauValue;
2115 for(int i=0;i<numSeqs;i++){
2116 double maxTau = -1.0000;
2118 for(int j=0;j<numOTUs;j++){
2119 if(bigTauMatrix[i * numOTUs + j] > maxTau){
2120 maxTau = bigTauMatrix[i * numOTUs + j];
2125 otuData[i] = maxOTU;
2128 nSeqsPerOTU.assign(numOTUs, 0);
2130 for(int i=0;i<numSeqs;i++){
2131 int index = otuData[i];
2133 singleTau[i] = 1.0000;
2136 aaP[index][nSeqsPerOTU[index]] = i;
2137 aaI[index][nSeqsPerOTU[index]] = i;
2139 nSeqsPerOTU[index]++;
2143 catch(exception& e) {
2144 m->errorOut(e, "ShhherCommand", "calcNewDistances");
2149 /**************************************************************************************************/
2151 void ShhherCommand::writeQualities(vector<int> otuCounts){
2154 string qualityFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.qual";
2156 ofstream qualityFile;
2157 m->openOutputFile(qualityFileName, qualityFile);
2159 qualityFile.setf(ios::fixed, ios::floatfield);
2160 qualityFile.setf(ios::showpoint);
2161 qualityFile << setprecision(6);
2163 vector<vector<int> > qualities(numOTUs);
2164 vector<double> pr(HOMOPS, 0);
2167 for(int i=0;i<numOTUs;i++){
2169 if (m->control_pressed) { break; }
2174 if(nSeqsPerOTU[i] > 0){
2175 qualities[i].assign(1024, -1);
2177 while(index < numFlowCells){
2178 double maxPrValue = 1e8;
2179 short maxPrIndex = -1;
2180 double count = 0.0000;
2182 pr.assign(HOMOPS, 0);
2184 for(int j=0;j<nSeqsPerOTU[i];j++){
2185 int lIndex = cumNumSeqs[i] + j;
2186 double tauValue = singleTau[seqNumber[lIndex]];
2187 int sequenceIndex = aaI[i][j];
2188 short intensity = flowDataIntI[sequenceIndex * numFlowCells + index];
2192 for(int s=0;s<HOMOPS;s++){
2193 pr[s] += tauValue * singleLookUp[s * NUMBINS + intensity];
2197 maxPrIndex = uniqueFlowgrams[centroids[i] * numFlowCells + index];
2198 maxPrValue = pr[maxPrIndex];
2200 if(count > MIN_COUNT){
2202 double norm = 0.0000;
2204 for(int s=0;s<HOMOPS;s++){
2205 norm += exp(-(pr[s] - maxPrValue));
2208 for(int s=1;s<=maxPrIndex;s++){
2210 double temp = 0.0000;
2212 U += exp(-(pr[s-1]-maxPrValue))/norm;
2220 temp = floor(-10 * temp);
2221 value = (int)floor(temp);
2222 if(value > 100){ value = 100; }
2224 qualities[i][base] = (int)value;
2234 if(otuCounts[i] > 0){
2235 qualityFile << '>' << seqNameVector[mapUniqueToSeq[i]] << endl;
2237 int j=4; //need to get past the first four bases
2238 while(qualities[i][j] != -1){
2239 qualityFile << qualities[i][j] << ' ';
2242 qualityFile << endl;
2245 qualityFile.close();
2246 outputNames.push_back(qualityFileName);
2249 catch(exception& e) {
2250 m->errorOut(e, "ShhherCommand", "writeQualities");
2255 /**************************************************************************************************/
2257 void ShhherCommand::writeSequences(vector<int> otuCounts){
2260 string fastaFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.fasta";
2262 m->openOutputFile(fastaFileName, fastaFile);
2264 vector<string> names(numOTUs, "");
2266 for(int i=0;i<numOTUs;i++){
2268 if (m->control_pressed) { break; }
2270 int index = centroids[i];
2272 if(otuCounts[i] > 0){
2273 fastaFile << '>' << seqNameVector[aaI[i][0]] << endl;
2277 for(int j=0;j<numFlowCells;j++){
2279 char base = flowOrder[j % 4];
2280 for(int k=0;k<uniqueFlowgrams[index * numFlowCells + j];k++){
2285 fastaFile << newSeq.substr(4) << endl;
2290 outputNames.push_back(fastaFileName);
2292 if(compositeFASTAFileName != ""){
2293 m->appendFiles(fastaFileName, compositeFASTAFileName);
2296 catch(exception& e) {
2297 m->errorOut(e, "ShhherCommand", "writeSequences");
2302 /**************************************************************************************************/
2304 void ShhherCommand::writeNames(vector<int> otuCounts){
2306 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.names";
2308 m->openOutputFile(nameFileName, nameFile);
2310 for(int i=0;i<numOTUs;i++){
2312 if (m->control_pressed) { break; }
2314 if(otuCounts[i] > 0){
2315 nameFile << seqNameVector[aaI[i][0]] << '\t' << seqNameVector[aaI[i][0]];
2317 for(int j=1;j<nSeqsPerOTU[i];j++){
2318 nameFile << ',' << seqNameVector[aaI[i][j]];
2325 outputNames.push_back(nameFileName);
2328 if(compositeNamesFileName != ""){
2329 m->appendFiles(nameFileName, compositeNamesFileName);
2332 catch(exception& e) {
2333 m->errorOut(e, "ShhherCommand", "writeNames");
2338 /**************************************************************************************************/
2340 void ShhherCommand::writeGroups(){
2342 string fileRoot = flowFileName.substr(0,flowFileName.find_last_of('.'));
2343 string groupFileName = fileRoot + ".shhh.groups";
2345 m->openOutputFile(groupFileName, groupFile);
2347 for(int i=0;i<numSeqs;i++){
2348 if (m->control_pressed) { break; }
2349 groupFile << seqNameVector[i] << '\t' << fileRoot << endl;
2352 outputNames.push_back(groupFileName);
2355 catch(exception& e) {
2356 m->errorOut(e, "ShhherCommand", "writeGroups");
2361 /**************************************************************************************************/
2363 void ShhherCommand::writeClusters(vector<int> otuCounts){
2365 string otuCountsFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.counts";
2366 ofstream otuCountsFile;
2367 m->openOutputFile(otuCountsFileName, otuCountsFile);
2369 string bases = flowOrder;
2371 for(int i=0;i<numOTUs;i++){
2373 if (m->control_pressed) {
2376 //output the translated version of the centroid sequence for the otu
2377 if(otuCounts[i] > 0){
2378 int index = centroids[i];
2380 otuCountsFile << "ideal\t";
2381 for(int j=8;j<numFlowCells;j++){
2382 char base = bases[j % 4];
2383 for(int s=0;s<uniqueFlowgrams[index * numFlowCells + j];s++){
2384 otuCountsFile << base;
2387 otuCountsFile << endl;
2389 for(int j=0;j<nSeqsPerOTU[i];j++){
2390 int sequence = aaI[i][j];
2391 otuCountsFile << seqNameVector[sequence] << '\t';
2395 for(int k=0;k<lengths[sequence];k++){
2396 char base = bases[k % 4];
2397 int freq = int(0.01 * (double)flowDataIntI[sequence * numFlowCells + k] + 0.5);
2399 for(int s=0;s<freq;s++){
2401 //otuCountsFile << base;
2404 otuCountsFile << newSeq.substr(4) << endl;
2406 otuCountsFile << endl;
2409 otuCountsFile.close();
2410 outputNames.push_back(otuCountsFileName);
2413 catch(exception& e) {
2414 m->errorOut(e, "ShhherCommand", "writeClusters");
2419 //**********************************************************************************************************************