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 //**********************************************************************************************************************
26 #define MIN_WEIGHT 0.1
27 #define MIN_TAU 0.0001
30 //**********************************************************************************************************************
32 vector<string> ShhherCommand::getValidParameters(){
35 "file", "flow", "lookup", "cutoff", "sigma", "outputdir","inputdir", "processors", "maxiter", "mindelta"
38 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
42 m->errorOut(e, "ShhherCommand", "getValidParameters");
47 //**********************************************************************************************************************
49 ShhherCommand::ShhherCommand(){
51 abort = true; calledHelp = true;
53 //initialize outputTypes
54 vector<string> tempOutNames;
55 outputTypes["pn.dist"] = tempOutNames;
59 m->errorOut(e, "ShhherCommand", "ShhherCommand");
64 //**********************************************************************************************************************
66 vector<string> ShhherCommand::getRequiredParameters(){
68 string Array[] = {"flow"};
69 vector<string> myArray (Array, Array+(sizeof(Array)/sizeof(string)));
73 m->errorOut(e, "ShhherCommand", "getRequiredParameters");
78 //**********************************************************************************************************************
80 vector<string> ShhherCommand::getRequiredFiles(){
82 vector<string> myArray;
86 m->errorOut(e, "ShhherCommand", "getRequiredFiles");
91 //**********************************************************************************************************************
93 ShhherCommand::ShhherCommand(string option) {
97 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
98 MPI_Comm_size(MPI_COMM_WORLD, &ncpus);
104 abort = false; calledHelp = false;
107 //allow user to run help
108 if(option == "help") { help(); abort = true; calledHelp = true; }
112 //valid paramters for this command
113 string AlignArray[] = {
114 "file", "flow", "lookup", "cutoff", "sigma", "outputdir","inputdir", "processors", "maxiter", "mindelta"
117 vector<string> myArray (AlignArray, AlignArray+(sizeof(AlignArray)/sizeof(string)));
119 OptionParser parser(option);
120 map<string,string> parameters = parser.getParameters();
122 ValidParameters validParameter;
123 map<string,string>::iterator it;
125 //check to make sure all parameters are valid for command
126 for (it = parameters.begin(); it != parameters.end(); it++) {
127 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
130 //initialize outputTypes
131 vector<string> tempOutNames;
132 outputTypes["pn.dist"] = tempOutNames;
133 // outputTypes["fasta"] = tempOutNames;
135 //if the user changes the input directory command factory will send this info to us in the output parameter
136 string inputDir = validParameter.validFile(parameters, "inputdir", false);
137 if (inputDir == "not found"){ inputDir = ""; }
140 it = parameters.find("flow");
141 //user has given a template file
142 if(it != parameters.end()){
143 path = m->hasPath(it->second);
144 //if the user has not given a path then, add inputdir. else leave path alone.
145 if (path == "") { parameters["flow"] = inputDir + it->second; }
148 it = parameters.find("lookup");
149 //user has given a template file
150 if(it != parameters.end()){
151 path = m->hasPath(it->second);
152 //if the user has not given a path then, add inputdir. else leave path alone.
153 if (path == "") { parameters["lookup"] = inputDir + it->second; }
156 it = parameters.find("file");
157 //user has given a template file
158 if(it != parameters.end()){
159 path = m->hasPath(it->second);
160 //if the user has not given a path then, add inputdir. else leave path alone.
161 if (path == "") { parameters["file"] = inputDir + it->second; }
166 //check for required parameters
167 flowFileName = validParameter.validFile(parameters, "flow", true);
168 flowFilesFileName = validParameter.validFile(parameters, "file", true);
169 if (flowFileName == "not found" && flowFilesFileName == "not found") {
170 m->mothurOut("values for either flow or file must be provided for the shhh.seqs command.");
171 m->mothurOutEndLine();
174 else if (flowFileName == "not open" || flowFilesFileName == "not open") { abort = true; }
176 if(flowFileName != "not found"){ compositeFASTAFileName = ""; }
178 compositeFASTAFileName = flowFilesFileName.substr(0, flowFilesFileName.length()-10) + "pn.fasta";
180 m->openOutputFile(compositeFASTAFileName, temp);
184 //if the user changes the output directory command factory will send this info to us in the output parameter
185 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
187 outputDir += m->hasPath(flowFileName); //if user entered a file with a path then preserve it
191 //check for optional parameter and set defaults
192 // ...at some point should added some additional type checking...
194 temp = validParameter.validFile(parameters, "lookup", true);
195 if (temp == "not found") { lookupFileName = "LookUp_Titanium.pat"; }
196 else if(temp == "not open") { abort = true; }
197 else { lookupFileName = temp; }
199 temp = validParameter.validFile(parameters, "processors", false);if (temp == "not found"){ temp = "1"; }
200 convert(temp, processors);
202 temp = validParameter.validFile(parameters, "cutoff", false); if (temp == "not found"){ temp = "0.01"; }
203 convert(temp, cutoff);
205 temp = validParameter.validFile(parameters, "mindelta", false); if (temp == "not found"){ temp = "0.000001"; }
206 convert(temp, minDelta);
208 temp = validParameter.validFile(parameters, "maxiter", false); if (temp == "not found"){ temp = "1000"; }
209 convert(temp, maxIters);
211 temp = validParameter.validFile(parameters, "sigma", false);if (temp == "not found") { temp = "60"; }
212 convert(temp, sigma);
214 globaldata = GlobalData::getInstance();
222 catch(exception& e) {
223 m->errorOut(e, "ShhherCommand", "ShhherCommand");
228 //**********************************************************************************************************************
230 ShhherCommand::~ShhherCommand(){}
232 //**********************************************************************************************************************
234 void ShhherCommand::help(){
236 m->mothurOut("The shhher command reads a file containing flowgrams and creates a file of corrected sequences.\n");
238 catch(exception& e) {
239 m->errorOut(e, "ShhherCommand", "help");
244 //**********************************************************************************************************************
246 int ShhherCommand::execute(){
248 if (abort == true) { if (calledHelp) { return 0; } return 2; }
253 double begClock = clock();
254 unsigned long int begTime = time(NULL);
258 for(int i=1;i<ncpus;i++){
259 MPI_Send(&abort, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
261 if(abort == 1){ return 0; }
265 m->mothurOut("\nGetting preliminary data...\n");
269 vector<string> flowFileVector;
270 if(flowFilesFileName != "not found"){
273 ifstream flowFilesFile;
274 m->openInputFile(flowFilesFileName, flowFilesFile);
275 while(flowFilesFile){
276 flowFilesFile >> fName;
277 flowFileVector.push_back(fName);
278 m->gobble(flowFilesFile);
282 flowFileVector.push_back(flowFileName);
284 int numFiles = flowFileVector.size();
286 for(int i=1;i<ncpus;i++){
287 MPI_Send(&numFiles, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
290 for(int i=0;i<numFiles;i++){
291 flowFileName = flowFileVector[i];
295 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
296 m->mothurOut("Reading flowgrams...\n");
299 m->mothurOut("Identifying unique flowgrams...\n");
302 m->mothurOut("Calculating distances between flowgrams...\n");
304 strcpy(fileName, flowFileName.c_str());
306 for(int i=1;i<ncpus;i++){
307 MPI_Send(&fileName[0], 1024, MPI_CHAR, i, tag, MPI_COMM_WORLD);
309 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
310 MPI_Send(&numUniques, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
311 MPI_Send(&numFlowCells, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
312 MPI_Send(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, i, tag, MPI_COMM_WORLD);
313 MPI_Send(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
314 MPI_Send(&mapUniqueToSeq[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
315 MPI_Send(&mapSeqToUnique[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
316 MPI_Send(&lengths[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
317 MPI_Send(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
318 MPI_Send(&cutoff, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
321 string distFileName = flowDistMPI(0, int(sqrt(1.0/float(ncpus)) * numUniques));
324 for(int i=1;i<ncpus;i++){
325 MPI_Recv(&done, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
327 m->appendFiles((distFileName + ".temp." + toString(i)), distFileName);
328 remove((distFileName + ".temp." + toString(i)).c_str());
331 string namesFileName = createNamesFile();
333 m->mothurOut("\nClustering flowgrams...\n");
334 string listFileName = cluster(distFileName, namesFileName);
336 getOTUData(listFileName);
339 for(int i=1;i<ncpus;i++){
340 MPI_Send(&numOTUs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
341 MPI_Send(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
342 MPI_Send(&uniqueFlowgrams[0], numFlowCells * numUniques, MPI_SHORT, i, tag, MPI_COMM_WORLD);
343 MPI_Send(&sigma, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
350 int numOTUsOnCPU = numOTUs / ncpus;
351 int numSeqsOnCPU = numSeqs / ncpus;
352 m->mothurOut("\nDenoising flowgrams...\n");
353 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
355 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
357 double cycClock = clock();
358 unsigned long int cycTime = time(NULL);
361 int total = singleTau.size();
362 for(int i=1;i<ncpus;i++){
363 MPI_Send(&total, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
364 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
365 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
367 MPI_Send(&singleTau[0], total, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
368 MPI_Send(&seqNumber[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
369 MPI_Send(&seqIndex[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
370 MPI_Send(&nSeqsPerOTU[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
371 MPI_Send(&cumNumSeqs[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
374 calcCentroidsDriver(0, numOTUsOnCPU);
376 for(int i=1;i<ncpus;i++){
377 int otuStart = i * numOTUs / ncpus;
378 int otuStop = (i + 1) * numOTUs / ncpus;
380 vector<int> tempCentroids(numOTUs, 0);
381 vector<short> tempChange(numOTUs, 0);
383 MPI_Recv(&tempCentroids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
384 MPI_Recv(&tempChange[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD, &status);
386 for(int j=otuStart;j<otuStop;j++){
387 centroids[j] = tempCentroids[j];
388 change[j] = tempChange[j];
392 maxDelta = getNewWeights();
393 double nLL = getLikelihood();
396 for(int i=1;i<ncpus;i++){
397 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
398 MPI_Send(&weight[0], numOTUs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
399 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
402 calcNewDistancesParent(0, numSeqsOnCPU);
404 total = singleTau.size();
406 for(int i=1;i<ncpus;i++){
408 int seqStart = i * numSeqs / ncpus;
409 int seqStop = (i + 1) * numSeqs / ncpus;
411 MPI_Recv(&childTotal, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
413 vector<int> childSeqIndex(childTotal, 0);
414 vector<double> childSingleTau(childTotal, 0);
415 vector<double> childDist(numSeqs * numOTUs, 0);
416 vector<int> otuIndex(childTotal, 0);
418 MPI_Recv(&childSeqIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
419 MPI_Recv(&childSingleTau[0], childTotal, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
420 MPI_Recv(&childDist[0], numOTUs * numSeqs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
421 MPI_Recv(&otuIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
423 int oldTotal = total;
425 singleTau.resize(total, 0);
426 seqIndex.resize(total, 0);
427 seqNumber.resize(total, 0);
431 for(int j=oldTotal;j<total;j++){
432 int otuI = otuIndex[childIndex];
433 int seqI = childSeqIndex[childIndex];
435 singleTau[j] = childSingleTau[childIndex];
437 aaP[otuI][nSeqsPerOTU[otuI]] = j;
438 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
443 int index = seqStart * numOTUs;
444 for(int j=seqStart;j<seqStop;j++){
445 for(int k=0;k<numOTUs;k++){
446 dist[index] = childDist[index];
454 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
456 if((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
458 for(int i=1;i<ncpus;i++){
459 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
464 for(int i=1;i<ncpus;i++){
465 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD); //send kill command
471 m->mothurOut("\nFinalizing...\n");
474 vector<int> otuCounts(numOTUs, 0);
475 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
476 calcCentroidsDriver(0, numOTUs);
477 writeQualities(otuCounts);
478 writeSequences(otuCounts);
479 writeNames(otuCounts);
480 writeClusters(otuCounts);
483 remove(distFileName.c_str());
484 remove(namesFileName.c_str());
485 remove(listFileName.c_str());
487 m->mothurOut("Total time to process " + toString(flowFileName) + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
494 MPI_Recv(&abort, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
495 if(abort){ return 0; }
498 MPI_Recv(&numFiles, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
500 for(int i=0;i<numFiles;i++){
501 //Now into the pyrodist part
503 MPI_Recv(&fileName, 1024, MPI_CHAR, 0, tag, MPI_COMM_WORLD, &status);
504 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
505 MPI_Recv(&numUniques, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
506 MPI_Recv(&numFlowCells, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
508 flowDataIntI.resize(numSeqs * numFlowCells);
509 flowDataPrI.resize(numSeqs * numFlowCells);
510 mapUniqueToSeq.resize(numSeqs);
511 mapSeqToUnique.resize(numSeqs);
512 lengths.resize(numSeqs);
513 jointLookUp.resize(NUMBINS * NUMBINS);
515 MPI_Recv(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
516 MPI_Recv(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
517 MPI_Recv(&mapUniqueToSeq[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
518 MPI_Recv(&mapSeqToUnique[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
519 MPI_Recv(&lengths[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
520 MPI_Recv(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
521 MPI_Recv(&cutoff, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
523 flowFileName = string(fileName);
524 int flowDistStart = int(sqrt(float(pid)/float(ncpus)) * numUniques);
525 int flowDistEnd = int(sqrt(float(pid+1)/float(ncpus)) * numUniques);
527 string distanceStringChild = flowDistMPI(flowDistStart, flowDistEnd);
530 MPI_Send(&done, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
532 //Now into the pyronoise part
533 MPI_Recv(&numOTUs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
535 singleLookUp.resize(HOMOPS * NUMBINS);
536 uniqueFlowgrams.resize(numUniques * numFlowCells);
537 weight.resize(numOTUs);
538 centroids.resize(numOTUs);
539 change.resize(numOTUs);
540 dist.assign(numOTUs * numSeqs, 0);
541 nSeqsPerOTU.resize(numOTUs);
542 cumNumSeqs.resize(numOTUs);
544 MPI_Recv(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
545 MPI_Recv(&uniqueFlowgrams[0], uniqueFlowgrams.size(), MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
546 MPI_Recv(&sigma, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
548 int startOTU = pid * numOTUs / ncpus;
549 int endOTU = (pid + 1) * numOTUs / ncpus;
551 int startSeq = pid * numSeqs / ncpus;
552 int endSeq = (pid + 1) * numSeqs /ncpus;
558 MPI_Recv(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
559 singleTau.assign(total, 0.0000);
560 seqNumber.assign(total, 0);
561 seqIndex.assign(total, 0);
563 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
564 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
565 MPI_Recv(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
566 MPI_Recv(&seqNumber[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
567 MPI_Recv(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
568 MPI_Recv(&nSeqsPerOTU[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
569 MPI_Recv(&cumNumSeqs[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
571 calcCentroidsDriver(startOTU, endOTU);
573 MPI_Send(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD);
574 MPI_Send(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD);
577 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
578 MPI_Recv(&weight[0], numOTUs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
579 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
581 vector<int> otuIndex(total, 0);
582 calcNewDistancesChildMPI(startSeq, endSeq, otuIndex);
583 total = otuIndex.size();
585 MPI_Send(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
586 MPI_Send(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
587 MPI_Send(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
588 MPI_Send(&dist[0], numOTUs * numSeqs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
589 MPI_Send(&otuIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
591 MPI_Recv(&live, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
596 MPI_Barrier(MPI_COMM_WORLD);
600 catch(exception& e) {
601 m->errorOut(e, "ShhherCommand", "execute");
606 /**************************************************************************************************/
608 string ShhherCommand::flowDistMPI(int startSeq, int stopSeq){
610 ostringstream outStream;
611 outStream.setf(ios::fixed, ios::floatfield);
612 outStream.setf(ios::dec, ios::basefield);
613 outStream.setf(ios::showpoint);
614 outStream.precision(6);
616 int begTime = time(NULL);
617 double begClock = clock();
619 for(int i=startSeq;i<stopSeq;i++){
620 for(int j=0;j<i;j++){
621 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
623 if(flowDistance < 1e-6){
624 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
626 else if(flowDistance <= cutoff){
627 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
631 m->mothurOut(toString(i) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
635 m->mothurOut(toString(stopSeq) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
637 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
638 if(pid != 0){ fDistFileName += ".temp." + toString(pid); }
640 ofstream distFile(fDistFileName.c_str());
641 distFile << outStream.str();
644 return fDistFileName;
646 catch(exception& e) {
647 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
653 //**********************************************************************************************************************
655 int ShhherCommand::execute(){
657 if (abort == true) { return 0; }
662 vector<string> flowFileVector;
663 if(flowFilesFileName != "not found"){
666 ifstream flowFilesFile;
667 m->openInputFile(flowFilesFileName, flowFilesFile);
668 while(flowFilesFile){
669 flowFilesFile >> fName;
670 flowFileVector.push_back(fName);
671 m->gobble(flowFilesFile);
675 flowFileVector.push_back(flowFileName);
677 int numFiles = flowFileVector.size();
680 for(int i=0;i<numFiles;i++){
681 flowFileName = flowFileVector[i];
683 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
684 m->mothurOut("Reading flowgrams...\n");
687 m->mothurOut("Identifying unique flowgrams...\n");
691 m->mothurOut("Calculating distances between flowgrams...\n");
692 string distFileName = createDistFile(processors);
693 string namesFileName = createNamesFile();
695 m->mothurOut("\nClustering flowgrams...\n");
696 string listFileName = cluster(distFileName, namesFileName);
697 getOTUData(listFileName);
704 double begClock = clock();
705 unsigned long int begTime = time(NULL);
708 cout << numOTUs << endl;
711 m->mothurOut("\nDenoising flowgrams...\n");
712 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
714 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
716 double cycClock = clock();
717 unsigned long int cycTime = time(NULL);
722 maxDelta = getNewWeights();
723 double nLL = getLikelihood();
730 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
734 m->mothurOut("\nFinalizing...\n");
738 vector<int> otuCounts(numOTUs, 0);
739 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
741 calcCentroidsDriver(0, numOTUs);
742 writeQualities(otuCounts);
743 writeSequences(otuCounts);
744 writeNames(otuCounts);
745 writeClusters(otuCounts);
748 remove(distFileName.c_str());
749 remove(namesFileName.c_str());
750 remove(listFileName.c_str());
752 m->mothurOut("Total time to process " + flowFileName + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
756 catch(exception& e) {
757 m->errorOut(e, "ShhherCommand", "execute");
762 /**************************************************************************************************/
764 void ShhherCommand::getFlowData(){
767 m->openInputFile(flowFileName, flowFile);
770 seqNameVector.clear();
772 flowDataIntI.clear();
776 int currentNumFlowCells;
780 flowFile >> numFlowCells;
781 int index = 0;//pcluster
782 while(!flowFile.eof()){
783 flowFile >> seqName >> currentNumFlowCells;
784 lengths.push_back(currentNumFlowCells);
786 seqNameVector.push_back(seqName);
787 nameMap[seqName] = index++;//pcluster
789 for(int i=0;i<numFlowCells;i++){
790 flowFile >> intensity;
791 if(intensity > 9.99) { intensity = 9.99; }
792 int intI = int(100 * intensity + 0.0001);
793 flowDataIntI.push_back(intI);
799 numSeqs = seqNameVector.size();
801 for(int i=0;i<numSeqs;i++){
802 int iNumFlowCells = i * numFlowCells;
803 for(int j=lengths[i];j<numFlowCells;j++){
804 flowDataIntI[iNumFlowCells + j] = 0;
809 catch(exception& e) {
810 m->errorOut(e, "ShhherCommand", "getFlowData");
815 /**************************************************************************************************/
817 void ShhherCommand::getSingleLookUp(){
819 // these are the -log probabilities that a signal corresponds to a particular homopolymer length
820 singleLookUp.assign(HOMOPS * NUMBINS, 0);
824 m->openInputFile(lookupFileName, lookUpFile);
826 for(int i=0;i<HOMOPS;i++){
828 lookUpFile >> logFracFreq;
830 for(int j=0;j<NUMBINS;j++) {
831 lookUpFile >> singleLookUp[index];
837 catch(exception& e) {
838 m->errorOut(e, "ShhherCommand", "getSingleLookUp");
843 /**************************************************************************************************/
845 void ShhherCommand::getJointLookUp(){
848 // the most likely joint probability (-log) that two intenities have the same polymer length
849 jointLookUp.resize(NUMBINS * NUMBINS, 0);
851 for(int i=0;i<NUMBINS;i++){
852 for(int j=0;j<NUMBINS;j++){
854 double minSum = 100000000;
856 for(int k=0;k<HOMOPS;k++){
857 double sum = singleLookUp[k * NUMBINS + i] + singleLookUp[k * NUMBINS + j];
859 if(sum < minSum) { minSum = sum; }
861 jointLookUp[i * NUMBINS + j] = minSum;
865 catch(exception& e) {
866 m->errorOut(e, "ShhherCommand", "getJointLookUp");
871 /**************************************************************************************************/
873 double ShhherCommand::getProbIntensity(int intIntensity){
876 double minNegLogProb = 100000000;
879 for(int i=0;i<HOMOPS;i++){//loop signal strength
880 float negLogProb = singleLookUp[i * NUMBINS + intIntensity];
881 if(negLogProb < minNegLogProb) { minNegLogProb = negLogProb; }
884 return minNegLogProb;
886 catch(exception& e) {
887 m->errorOut(e, "ShhherCommand", "getProbIntensity");
892 /**************************************************************************************************/
894 void ShhherCommand::getUniques(){
899 uniqueFlowgrams.assign(numFlowCells * numSeqs, -1);
900 uniqueCount.assign(numSeqs, 0); // anWeights
901 uniqueLengths.assign(numSeqs, 0);
902 mapSeqToUnique.assign(numSeqs, -1);
903 mapUniqueToSeq.assign(numSeqs, -1);
905 vector<short> uniqueFlowDataIntI(numFlowCells * numSeqs, -1);
907 for(int i=0;i<numSeqs;i++){
910 vector<short> current(numFlowCells);
911 for(int j=0;j<numFlowCells;j++){
912 current[j] = short(((flowDataIntI[i * numFlowCells + j] + 50.0)/100.0));
915 for(int j=0;j<numUniques;j++){
916 int offset = j * numFlowCells;
920 if(lengths[i] < uniqueLengths[j]) { shorterLength = lengths[i]; }
921 else { shorterLength = uniqueLengths[j]; }
923 for(int k=0;k<shorterLength;k++){
924 if(current[k] != uniqueFlowgrams[offset + k]){
931 mapSeqToUnique[i] = j;
934 if(lengths[i] > uniqueLengths[j]) { uniqueLengths[j] = lengths[i]; }
940 if(index == numUniques){
941 uniqueLengths[numUniques] = lengths[i];
942 uniqueCount[numUniques] = 1;
943 mapSeqToUnique[i] = numUniques;//anMap
944 mapUniqueToSeq[numUniques] = i;//anF
946 for(int k=0;k<numFlowCells;k++){
947 uniqueFlowgrams[numUniques * numFlowCells + k] = current[k];
948 uniqueFlowDataIntI[numUniques * numFlowCells + k] = flowDataIntI[i * numFlowCells + k];
954 uniqueFlowDataIntI.resize(numFlowCells * numUniques);
955 uniqueLengths.resize(numUniques);
957 flowDataPrI.resize(numSeqs * numFlowCells, 0);
958 for(int i=0;i<flowDataPrI.size();i++) { flowDataPrI[i] = getProbIntensity(flowDataIntI[i]); }
960 catch(exception& e) {
961 m->errorOut(e, "ShhherCommand", "getUniques");
966 /**************************************************************************************************/
968 float ShhherCommand::calcPairwiseDist(int seqA, int seqB){
970 int minLength = lengths[mapSeqToUnique[seqA]];
971 if(lengths[seqB] < minLength){ minLength = lengths[mapSeqToUnique[seqB]]; }
973 int ANumFlowCells = seqA * numFlowCells;
974 int BNumFlowCells = seqB * numFlowCells;
978 for(int i=0;i<minLength;i++){
979 int flowAIntI = flowDataIntI[ANumFlowCells + i];
980 float flowAPrI = flowDataPrI[ANumFlowCells + i];
982 int flowBIntI = flowDataIntI[BNumFlowCells + i];
983 float flowBPrI = flowDataPrI[BNumFlowCells + i];
984 dist += jointLookUp[flowAIntI * NUMBINS + flowBIntI] - flowAPrI - flowBPrI;
987 dist /= (float) minLength;
990 catch(exception& e) {
991 m->errorOut(e, "ShhherCommand", "calcPairwiseDist");
996 /**************************************************************************************************/
998 void ShhherCommand::flowDistParentFork(string distFileName, int startSeq, int stopSeq){
1001 ostringstream outStream;
1002 outStream.setf(ios::fixed, ios::floatfield);
1003 outStream.setf(ios::dec, ios::basefield);
1004 outStream.setf(ios::showpoint);
1005 outStream.precision(6);
1007 int begTime = time(NULL);
1008 double begClock = clock();
1010 for(int i=startSeq;i<stopSeq;i++){
1011 for(int j=0;j<i;j++){
1012 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
1014 if(flowDistance < 1e-6){
1015 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
1017 else if(flowDistance <= cutoff){
1018 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
1022 m->mothurOut(toString(i) + "\t" + toString(time(NULL) - begTime));
1023 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1024 m->mothurOutEndLine();
1027 m->mothurOut(toString(stopSeq-1) + "\t" + toString(time(NULL) - begTime));
1028 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1029 m->mothurOutEndLine();
1031 ofstream distFile(distFileName.c_str());
1032 distFile << outStream.str();
1035 catch(exception& e) {
1036 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
1041 /**************************************************************************************************/
1043 string ShhherCommand::createDistFile(int processors){
1045 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
1047 unsigned long int begTime = time(NULL);
1048 double begClock = clock();
1053 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1054 if(processors == 1) { flowDistParentFork(fDistFileName, 0, numUniques); }
1055 else{ //you have multiple processors
1057 if (numSeqs < processors){ processors = 1; }
1059 vector<int> start(processors, 0);
1060 vector<int> end(processors, 0);
1062 for (int i = 0; i < processors; i++) {
1063 start[i] = int(sqrt(float(i)/float(processors)) * numUniques);
1064 end[i] = int(sqrt(float(i+1)/float(processors)) * numUniques);
1068 vector<int> processIDs;
1070 //loop through and create all the processes you want
1071 while (process != processors) {
1075 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1077 }else if (pid == 0){
1078 flowDistParentFork(fDistFileName + toString(getpid()) + ".temp", start[process], end[process]);
1081 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1083 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1088 //parent does its part
1089 flowDistParentFork(fDistFileName, start[0], end[0]);
1091 //force parent to wait until all the processes are done
1092 for (int i=0;i<processIDs.size();i++) {
1093 int temp = processIDs[i];
1097 //append and remove temp files
1098 for (int i=0;i<processIDs.size();i++) {
1099 m->appendFiles((fDistFileName + toString(processIDs[i]) + ".temp"), fDistFileName);
1100 remove((fDistFileName + toString(processIDs[i]) + ".temp").c_str());
1106 flowDistParentFork(fDistFileName, 0, numUniques);
1109 m->mothurOutEndLine();
1111 m->mothurOut("Total time: " + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/CLOCKS_PER_SEC) + '\n');
1114 return fDistFileName;
1116 catch(exception& e) {
1117 m->errorOut(e, "ShhherCommand", "createDistFile");
1123 /**************************************************************************************************/
1125 string ShhherCommand::createNamesFile(){
1128 vector<string> duplicateNames(numUniques, "");
1129 for(int i=0;i<numSeqs;i++){
1130 duplicateNames[mapSeqToUnique[i]] += seqNameVector[i] + ',';
1133 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.names";
1136 m->openOutputFile(nameFileName, nameFile);
1138 for(int i=0;i<numUniques;i++){
1139 // nameFile << seqNameVector[mapUniqueToSeq[i]] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1140 nameFile << mapUniqueToSeq[i] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1144 return nameFileName;
1146 catch(exception& e) {
1147 m->errorOut(e, "ShhherCommand", "createNamesFile");
1152 //**********************************************************************************************************************
1154 string ShhherCommand::cluster(string distFileName, string namesFileName){
1158 globaldata->setNameFile(namesFileName);
1159 globaldata->setColumnFile(distFileName);
1160 globaldata->setFormat("column");
1162 ReadMatrix* read = new ReadColumnMatrix(distFileName);
1163 read->setCutoff(cutoff);
1165 NameAssignment* clusterNameMap = new NameAssignment(namesFileName);
1166 clusterNameMap->readMap();
1167 read->read(clusterNameMap);
1169 ListVector* list = read->getListVector();
1170 SparseMatrix* matrix = read->getMatrix();
1173 delete clusterNameMap;
1175 RAbundVector* rabund = new RAbundVector(list->getRAbundVector());
1177 Cluster* cluster = new CompleteLinkage(rabund, list, matrix, cutoff, "furthest");
1178 string tag = cluster->getTag();
1180 double clusterCutoff = cutoff;
1181 while (matrix->getSmallDist() <= clusterCutoff && matrix->getNNodes() > 0){
1182 cluster->update(clusterCutoff);
1185 list->setLabel(toString(cutoff));
1187 string listFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.list";
1189 m->openOutputFile(listFileName, listFile);
1190 list->print(listFile);
1193 delete matrix; delete cluster; delete rabund; delete list;
1195 return listFileName;
1197 catch(exception& e) {
1198 m->errorOut(e, "ShhherCommand", "cluster");
1203 /**************************************************************************************************/
1205 void ShhherCommand::getOTUData(string listFileName){
1209 m->openInputFile(listFileName, listFile);
1212 listFile >> label >> numOTUs;
1214 otuData.assign(numSeqs, 0);
1215 cumNumSeqs.assign(numOTUs, 0);
1216 nSeqsPerOTU.assign(numOTUs, 0);
1217 aaP.clear();aaP.resize(numOTUs);
1223 string singleOTU = "";
1225 for(int i=0;i<numOTUs;i++){
1227 listFile >> singleOTU;
1229 istringstream otuString(singleOTU);
1233 string seqName = "";
1235 for(int j=0;j<singleOTU.length();j++){
1236 char letter = otuString.get();
1242 map<string,int>::iterator nmIt = nameMap.find(seqName);
1243 int index = nmIt->second;
1245 nameMap.erase(nmIt);
1249 aaP[i].push_back(index);
1254 map<string,int>::iterator nmIt = nameMap.find(seqName);
1256 int index = nmIt->second;
1257 nameMap.erase(nmIt);
1261 aaP[i].push_back(index);
1266 sort(aaP[i].begin(), aaP[i].end());
1267 for(int j=0;j<nSeqsPerOTU[i];j++){
1268 seqNumber.push_back(aaP[i][j]);
1270 for(int j=nSeqsPerOTU[i];j<numSeqs;j++){
1271 aaP[i].push_back(0);
1277 for(int i=1;i<numOTUs;i++){
1278 cumNumSeqs[i] = cumNumSeqs[i-1] + nSeqsPerOTU[i-1];
1281 seqIndex = seqNumber;
1285 for(int i=0;i<seqNumber.size();i++){
1286 cout << seqNumber[i] << ' ';
1292 catch(exception& e) {
1293 m->errorOut(e, "ShhherCommand", "getOTUData");
1298 /**************************************************************************************************/
1300 void ShhherCommand::initPyroCluster(){
1302 dist.assign(numSeqs * numOTUs, 0);
1303 change.assign(numOTUs, 1);
1304 centroids.assign(numOTUs, -1);
1305 weight.assign(numOTUs, 0);
1306 singleTau.assign(numSeqs, 1.0);
1308 nSeqsBreaks.assign(processors+1, 0);
1309 nOTUsBreaks.assign(processors+1, 0);
1311 cout << numSeqs << '\t' << numOTUs << '\t' << processors << endl;
1314 for(int i=0;i<processors;i++){
1315 nSeqsBreaks[i+1] = nSeqsBreaks[i] + (int)((double) numSeqs / (double) processors);
1316 nOTUsBreaks[i+1] = nOTUsBreaks[i] + (int)((double) numOTUs / (double) processors);
1318 nSeqsBreaks[processors] = numSeqs;
1319 nOTUsBreaks[processors] = numOTUs;
1321 catch(exception& e) {
1322 m->errorOut(e, "ShhherCommand", "initPyroCluster");
1327 /**************************************************************************************************/
1329 void ShhherCommand::fill(){
1332 for(int i=0;i<numOTUs;i++){
1333 cumNumSeqs[i] = index;
1334 for(int j=0;j<nSeqsPerOTU[i];j++){
1335 seqNumber[index] = aaP[i][j];
1336 seqIndex[index] = aaI[i][j];
1342 catch(exception& e) {
1343 m->errorOut(e, "ShhherCommand", "fill");
1348 /**************************************************************************************************/
1350 void ShhherCommand::calcCentroids(){
1353 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1355 if(processors == 1) {
1356 calcCentroidsDriver(0, numOTUs);
1358 else{ //you have multiple processors
1359 if (numOTUs < processors){ processors = 1; }
1362 vector<int> processIDs;
1364 //loop through and create all the processes you want
1365 while (process != processors) {
1369 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1371 }else if (pid == 0){
1372 calcCentroidsDriver(nOTUsBreaks[process], nOTUsBreaks[process+1]);
1375 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1377 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1382 //parent does its part
1383 calcCentroidsDriver(nOTUsBreaks[0], nOTUsBreaks[1]);
1385 //force parent to wait until all the processes are done
1386 for (int i=0;i<processIDs.size();i++) {
1387 int temp = processIDs[i];
1393 calcCentroidsDriver(0, numOTUs);
1396 catch(exception& e) {
1397 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1402 /**************************************************************************************************/
1404 void ShhherCommand::calcCentroidsDriver(int start, int finish){
1406 //this function gets the most likely homopolymer length at a flow position for a group of sequences
1411 // for(int i=0;i<seqNumber.size();i++){
1412 // cout << seqNumber[i] << ' ';
1415 for(int i=start;i<finish;i++){
1419 int minFlowGram = 100000000;
1420 double minFlowValue = 1e8;
1421 change[i] = 0; //FALSE
1423 for(int j=0;j<nSeqsPerOTU[i];j++){
1424 count += singleTau[seqNumber[cumNumSeqs[i] + j]];
1427 if(nSeqsPerOTU[i] > 0 && count > MIN_COUNT){
1428 vector<double> adF(nSeqsPerOTU[i]);
1429 vector<int> anL(nSeqsPerOTU[i]);
1431 for(int j=0;j<nSeqsPerOTU[i];j++){
1432 int index = cumNumSeqs[i] + j;
1433 int nI = seqIndex[index];
1434 int nIU = mapSeqToUnique[nI];
1437 for(k=0;k<position;k++){
1443 anL[position] = nIU;
1444 adF[position] = 0.0000;
1449 for(int j=0;j<nSeqsPerOTU[i];j++){
1450 int index = cumNumSeqs[i] + j;
1451 int nI = seqIndex[index];
1453 double tauValue = singleTau[seqNumber[index]];
1455 for(int k=0;k<position;k++){
1456 double dist = getDistToCentroid(anL[k], nI, lengths[nI]);
1457 adF[k] += dist * tauValue;
1461 for(int j=0;j<position;j++){
1462 if(adF[j] < minFlowValue){
1464 minFlowValue = adF[j];
1468 if(centroids[i] != anL[minFlowGram]){
1470 centroids[i] = anL[minFlowGram];
1473 else if(centroids[i] != -1){
1479 catch(exception& e) {
1480 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1485 /**************************************************************************************************/
1487 double ShhherCommand::getDistToCentroid(int cent, int flow, int length){
1490 int flowAValue = cent * numFlowCells;
1491 int flowBValue = flow * numFlowCells;
1495 for(int i=0;i<length;i++){
1496 dist += singleLookUp[uniqueFlowgrams[flowAValue] * NUMBINS + flowDataIntI[flowBValue]];
1501 return dist / (double)length;
1503 catch(exception& e) {
1504 m->errorOut(e, "ShhherCommand", "getDistToCentroid");
1509 /**************************************************************************************************/
1511 double ShhherCommand::getNewWeights(){
1514 double maxChange = 0;
1516 for(int i=0;i<numOTUs;i++){
1518 double difference = weight[i];
1521 for(int j=0;j<nSeqsPerOTU[i];j++){
1522 int index = cumNumSeqs[i] + j;
1523 double tauValue = singleTau[seqNumber[index]];
1524 weight[i] += tauValue;
1527 difference = fabs(weight[i] - difference);
1528 if(difference > maxChange){ maxChange = difference; }
1532 catch(exception& e) {
1533 m->errorOut(e, "ShhherCommand", "getNewWeights");
1538 /**************************************************************************************************/
1540 double ShhherCommand::getLikelihood(){
1544 vector<long double> P(numSeqs, 0);
1547 for(int i=0;i<numOTUs;i++){
1548 if(weight[i] > MIN_WEIGHT){
1554 for(int i=0;i<numOTUs;i++){
1555 for(int j=0;j<nSeqsPerOTU[i];j++){
1556 int index = cumNumSeqs[i] + j;
1557 int nI = seqIndex[index];
1558 double singleDist = dist[seqNumber[index]];
1560 P[nI] += weight[i] * exp(-singleDist * sigma);
1564 for(int i=0;i<numSeqs;i++){
1565 if(P[i] == 0){ P[i] = DBL_EPSILON; }
1570 nLL = nLL -(double)numSeqs * log(sigma);
1574 catch(exception& e) {
1575 m->errorOut(e, "ShhherCommand", "getNewWeights");
1580 /**************************************************************************************************/
1582 void ShhherCommand::checkCentroids(){
1584 vector<int> unique(numOTUs, 1);
1586 for(int i=0;i<numOTUs;i++){
1587 if(centroids[i] == -1 || weight[i] < MIN_WEIGHT){
1592 for(int i=0;i<numOTUs;i++){
1594 for(int j=i+1;j<numOTUs;j++){
1597 if(centroids[j] == centroids[i]){
1601 weight[i] += weight[j];
1609 catch(exception& e) {
1610 m->errorOut(e, "ShhherCommand", "checkCentroids");
1615 /**************************************************************************************************/
1617 void ShhherCommand::calcNewDistances(){
1620 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1622 if(processors == 1) {
1623 calcNewDistancesParent(0, numSeqs);
1625 else{ //you have multiple processors
1626 if (numSeqs < processors){ processors = 1; }
1628 vector<vector<int> > child_otuIndex(processors);
1629 vector<vector<int> > child_seqIndex(processors);
1630 vector<vector<double> > child_singleTau(processors);
1631 vector<int> totals(processors);
1634 vector<int> processIDs;
1636 //loop through and create all the processes you want
1637 while (process != processors) {
1641 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1643 }else if (pid == 0){
1644 calcNewDistancesChild(nSeqsBreaks[process], nSeqsBreaks[process+1], child_otuIndex[process], child_seqIndex[process], child_singleTau[process]);
1645 totals[process] = child_otuIndex[process].size();
1649 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1651 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1656 //parent does its part
1657 calcNewDistancesParent(nSeqsBreaks[0], nSeqsBreaks[1]);
1658 int total = seqIndex.size();
1660 //force parent to wait until all the processes are done
1661 for (int i=0;i<processIDs.size();i++) {
1662 int temp = processIDs[i];
1666 for(int i=1;i<processors;i++){
1667 int oldTotal = total;
1670 singleTau.resize(total, 0);
1671 seqIndex.resize(total, 0);
1672 seqNumber.resize(total, 0);
1676 for(int j=oldTotal;j<total;j++){
1677 int otuI = child_otuIndex[i][childIndex];
1678 int seqI = child_seqIndex[i][childIndex];
1680 singleTau[j] = child_singleTau[i][childIndex];
1681 aaP[otuI][nSeqsPerOTU[otuI]] = j;
1682 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
1683 nSeqsPerOTU[otuI]++;
1690 calcNewDistancesParent(0, numSeqs);
1693 catch(exception& e) {
1694 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1699 /**************************************************************************************************/
1701 void ShhherCommand::calcNewDistancesChildMPI(int startSeq, int stopSeq, vector<int>& otuIndex){
1704 vector<double> newTau(numOTUs,0);
1705 vector<double> norms(numSeqs, 0);
1712 for(int i=startSeq;i<stopSeq;i++){
1713 double offset = 1e8;
1714 int indexOffset = i * numOTUs;
1716 for(int j=0;j<numOTUs;j++){
1718 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1719 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1721 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1722 offset = dist[indexOffset + j];
1726 for(int j=0;j<numOTUs;j++){
1727 if(weight[j] > MIN_WEIGHT){
1728 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1729 norms[i] += newTau[j];
1736 for(int j=0;j<numOTUs;j++){
1738 newTau[j] /= norms[i];
1740 if(newTau[j] > MIN_TAU){
1741 otuIndex.push_back(j);
1742 seqIndex.push_back(i);
1743 singleTau.push_back(newTau[j]);
1749 catch(exception& e) {
1750 m->errorOut(e, "ShhherCommand", "calcNewDistancesChildMPI");
1755 /**************************************************************************************************/
1757 void ShhherCommand::calcNewDistancesChild(int startSeq, int stopSeq, vector<int>& child_otuIndex, vector<int>& child_seqIndex, vector<double>& child_singleTau){
1760 vector<double> newTau(numOTUs,0);
1761 vector<double> norms(numSeqs, 0);
1762 child_otuIndex.resize(0);
1763 child_seqIndex.resize(0);
1764 child_singleTau.resize(0);
1766 for(int i=startSeq;i<stopSeq;i++){
1767 double offset = 1e8;
1768 int indexOffset = i * numOTUs;
1771 for(int j=0;j<numOTUs;j++){
1772 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1773 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1776 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1777 offset = dist[indexOffset + j];
1781 for(int j=0;j<numOTUs;j++){
1782 if(weight[j] > MIN_WEIGHT){
1783 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1784 norms[i] += newTau[j];
1791 for(int j=0;j<numOTUs;j++){
1792 newTau[j] /= norms[i];
1794 if(newTau[j] > MIN_TAU){
1795 child_otuIndex.push_back(j);
1796 child_seqIndex.push_back(i);
1797 child_singleTau.push_back(newTau[j]);
1802 catch(exception& e) {
1803 m->errorOut(e, "ShhherCommand", "calcNewDistancesChild");
1808 /**************************************************************************************************/
1810 void ShhherCommand::calcNewDistancesParent(int startSeq, int stopSeq){
1815 vector<double> newTau(numOTUs,0);
1816 vector<double> norms(numSeqs, 0);
1817 nSeqsPerOTU.assign(numOTUs, 0);
1819 for(int i=startSeq;i<stopSeq;i++){
1820 int indexOffset = i * numOTUs;
1822 double offset = 1e8;
1824 for(int j=0;j<numOTUs;j++){
1825 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1826 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1829 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1830 offset = dist[indexOffset + j];
1834 for(int j=0;j<numOTUs;j++){
1835 if(weight[j] > MIN_WEIGHT){
1836 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1837 norms[i] += newTau[j];
1844 for(int j=0;j<numOTUs;j++){
1845 newTau[j] /= norms[i];
1848 for(int j=0;j<numOTUs;j++){
1849 if(newTau[j] > MIN_TAU){
1851 int oldTotal = total;
1855 singleTau.resize(total, 0);
1856 seqNumber.resize(total, 0);
1857 seqIndex.resize(total, 0);
1859 singleTau[oldTotal] = newTau[j];
1861 aaP[j][nSeqsPerOTU[j]] = oldTotal;
1862 aaI[j][nSeqsPerOTU[j]] = i;
1868 catch(exception& e) {
1869 m->errorOut(e, "ShhherCommand", "calcNewDistancesParent");
1874 /**************************************************************************************************/
1876 void ShhherCommand::setOTUs(){
1879 vector<double> bigTauMatrix(numOTUs * numSeqs, 0.0000);
1881 for(int i=0;i<numOTUs;i++){
1882 for(int j=0;j<nSeqsPerOTU[i];j++){
1883 int index = cumNumSeqs[i] + j;
1884 double tauValue = singleTau[seqNumber[index]];
1885 int sIndex = seqIndex[index];
1886 bigTauMatrix[sIndex * numOTUs + i] = tauValue;
1890 for(int i=0;i<numSeqs;i++){
1891 double maxTau = -1.0000;
1893 for(int j=0;j<numOTUs;j++){
1894 if(bigTauMatrix[i * numOTUs + j] > maxTau){
1895 maxTau = bigTauMatrix[i * numOTUs + j];
1900 otuData[i] = maxOTU;
1903 nSeqsPerOTU.assign(numOTUs, 0);
1905 for(int i=0;i<numSeqs;i++){
1906 int index = otuData[i];
1908 singleTau[i] = 1.0000;
1911 aaP[index][nSeqsPerOTU[index]] = i;
1912 aaI[index][nSeqsPerOTU[index]] = i;
1914 nSeqsPerOTU[index]++;
1918 catch(exception& e) {
1919 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1924 /**************************************************************************************************/
1926 void ShhherCommand::writeQualities(vector<int> otuCounts){
1929 string qualityFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.qual";
1931 ofstream qualityFile;
1932 m->openOutputFile(qualityFileName, qualityFile);
1934 qualityFile.setf(ios::fixed, ios::floatfield);
1935 qualityFile.setf(ios::showpoint);
1936 qualityFile << setprecision(6);
1938 vector<vector<int> > qualities(numOTUs);
1939 vector<double> pr(HOMOPS, 0);
1942 for(int i=0;i<numOTUs;i++){
1946 if(nSeqsPerOTU[i] > 0){
1947 qualities[i].assign(1024, -1);
1949 while(index < numFlowCells){
1950 double maxPrValue = 1e8;
1951 short maxPrIndex = -1;
1952 double count = 0.0000;
1954 pr.assign(HOMOPS, 0);
1956 for(int j=0;j<nSeqsPerOTU[i];j++){
1957 int lIndex = cumNumSeqs[i] + j;
1958 double tauValue = singleTau[seqNumber[lIndex]];
1959 int sequenceIndex = aaI[i][j];
1960 short intensity = flowDataIntI[sequenceIndex * numFlowCells + index];
1964 for(int s=0;s<HOMOPS;s++){
1965 pr[s] += tauValue * singleLookUp[s * NUMBINS + intensity];
1969 maxPrIndex = uniqueFlowgrams[centroids[i] * numFlowCells + index];
1970 maxPrValue = pr[maxPrIndex];
1972 if(count > MIN_COUNT){
1974 double norm = 0.0000;
1976 for(int s=0;s<HOMOPS;s++){
1977 norm += exp(-(pr[s] - maxPrValue));
1980 for(int s=1;s<=maxPrIndex;s++){
1982 double temp = 0.0000;
1984 U += exp(-(pr[s-1]-maxPrValue))/norm;
1992 temp = floor(-10 * temp);
1993 value = (int)floor(temp);
1994 if(value > 100){ value = 100; }
1996 qualities[i][base] = (int)value;
2006 if(otuCounts[i] > 0){
2007 qualityFile << '>' << seqNameVector[mapUniqueToSeq[i]] << endl;
2009 int j=4; //need to get past the first four bases
2010 while(qualities[i][j] != -1){
2011 qualityFile << qualities[i][j] << ' ';
2014 qualityFile << endl;
2017 qualityFile.close();
2020 catch(exception& e) {
2021 m->errorOut(e, "ShhherCommand", "writeQualities");
2026 /**************************************************************************************************/
2028 void ShhherCommand::writeSequences(vector<int> otuCounts){
2030 string bases = "TACG";
2032 string fastaFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.fasta";
2034 m->openOutputFile(fastaFileName, fastaFile);
2036 vector<string> names(numOTUs, "");
2038 for(int i=0;i<numOTUs;i++){
2039 int index = centroids[i];
2041 if(otuCounts[i] > 0){
2042 fastaFile << '>' << seqNameVector[aaI[i][0]] << endl;
2044 for(int j=8;j<numFlowCells;j++){
2046 char base = bases[j % 4];
2047 for(int k=0;k<uniqueFlowgrams[index * numFlowCells + j];k++){
2056 if(compositeFASTAFileName != ""){
2057 m->appendFiles(fastaFileName, compositeFASTAFileName);
2060 catch(exception& e) {
2061 m->errorOut(e, "ShhherCommand", "writeSequences");
2066 /**************************************************************************************************/
2068 void ShhherCommand::writeNames(vector<int> otuCounts){
2070 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.final.names";
2072 m->openOutputFile(nameFileName, nameFile);
2074 for(int i=0;i<numOTUs;i++){
2075 if(otuCounts[i] > 0){
2076 nameFile << seqNameVector[aaI[i][0]] << '\t' << seqNameVector[aaI[i][0]];
2078 for(int j=1;j<nSeqsPerOTU[i];j++){
2079 nameFile << ',' << seqNameVector[aaI[i][j]];
2087 catch(exception& e) {
2088 m->errorOut(e, "ShhherCommand", "writeNames");
2093 /**************************************************************************************************/
2095 void ShhherCommand::writeGroups(){
2097 string fileRoot = flowFileName.substr(0,flowFileName.find_last_of('.'));
2098 string groupFileName = fileRoot + ".pn.groups";
2100 m->openOutputFile(groupFileName, groupFile);
2102 for(int i=0;i<numSeqs;i++){
2103 groupFile << seqNameVector[i] << '\t' << fileRoot << endl;
2107 catch(exception& e) {
2108 m->errorOut(e, "ShhherCommand", "writeGroups");
2113 /**************************************************************************************************/
2115 void ShhherCommand::writeClusters(vector<int> otuCounts){
2117 string otuCountsFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.counts";
2118 ofstream otuCountsFile;
2119 m->openOutputFile(otuCountsFileName, otuCountsFile);
2121 string bases = "TACG";
2123 for(int i=0;i<numOTUs;i++){
2124 //output the translated version of the centroid sequence for the otu
2125 if(otuCounts[i] > 0){
2126 int index = centroids[i];
2128 otuCountsFile << "ideal\t";
2129 for(int j=8;j<numFlowCells;j++){
2130 char base = bases[j % 4];
2131 for(int s=0;s<uniqueFlowgrams[index * numFlowCells + j];s++){
2132 otuCountsFile << base;
2135 otuCountsFile << endl;
2137 for(int j=0;j<nSeqsPerOTU[i];j++){
2138 int sequence = aaI[i][j];
2139 otuCountsFile << seqNameVector[sequence] << '\t';
2141 for(int k=8;k<lengths[sequence];k++){
2142 char base = bases[k % 4];
2143 int freq = int(0.01 * (double)flowDataIntI[sequence * numFlowCells + k] + 0.5);
2145 for(int s=0;s<freq;s++){
2146 otuCountsFile << base;
2149 otuCountsFile << endl;
2151 otuCountsFile << endl;
2154 otuCountsFile.close();
2156 catch(exception& e) {
2157 m->errorOut(e, "ShhherCommand", "writeClusters");
2162 //**********************************************************************************************************************