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"
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"
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);
707 m->mothurOut("\nDenoising flowgrams...\n");
708 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
710 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
712 double cycClock = clock();
713 unsigned long int cycTime = time(NULL);
718 maxDelta = getNewWeights();
719 double nLL = getLikelihood();
726 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
730 m->mothurOut("\nFinalizing...\n");
734 vector<int> otuCounts(numOTUs, 0);
735 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
737 calcCentroidsDriver(0, numOTUs);
738 writeQualities(otuCounts);
739 writeSequences(otuCounts);
740 writeNames(otuCounts);
741 writeClusters(otuCounts);
744 remove(distFileName.c_str());
745 remove(namesFileName.c_str());
746 remove(listFileName.c_str());
748 m->mothurOut("Total time to process " + flowFileName + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
752 catch(exception& e) {
753 m->errorOut(e, "ShhherCommand", "execute");
758 /**************************************************************************************************/
760 void ShhherCommand::getFlowData(){
763 m->openInputFile(flowFileName, flowFile);
766 seqNameVector.clear();
768 flowDataIntI.clear();
772 int currentNumFlowCells;
776 flowFile >> numFlowCells;
777 int index = 0;//pcluster
778 while(!flowFile.eof()){
779 flowFile >> seqName >> currentNumFlowCells;
780 lengths.push_back(currentNumFlowCells);
782 seqNameVector.push_back(seqName);
783 nameMap[seqName] = index++;//pcluster
785 for(int i=0;i<numFlowCells;i++){
786 flowFile >> intensity;
787 if(intensity > 9.99) { intensity = 9.99; }
788 int intI = int(100 * intensity + 0.0001);
789 flowDataIntI.push_back(intI);
795 numSeqs = seqNameVector.size();
797 for(int i=0;i<numSeqs;i++){
798 int iNumFlowCells = i * numFlowCells;
799 for(int j=lengths[i];j<numFlowCells;j++){
800 flowDataIntI[iNumFlowCells + j] = 0;
805 catch(exception& e) {
806 m->errorOut(e, "ShhherCommand", "getFlowData");
811 /**************************************************************************************************/
813 void ShhherCommand::getSingleLookUp(){
815 // these are the -log probabilities that a signal corresponds to a particular homopolymer length
816 singleLookUp.assign(HOMOPS * NUMBINS, 0);
820 m->openInputFile(lookupFileName, lookUpFile);
822 for(int i=0;i<HOMOPS;i++){
824 lookUpFile >> logFracFreq;
826 for(int j=0;j<NUMBINS;j++) {
827 lookUpFile >> singleLookUp[index];
833 catch(exception& e) {
834 m->errorOut(e, "ShhherCommand", "getSingleLookUp");
839 /**************************************************************************************************/
841 void ShhherCommand::getJointLookUp(){
844 // the most likely joint probability (-log) that two intenities have the same polymer length
845 jointLookUp.resize(NUMBINS * NUMBINS, 0);
847 for(int i=0;i<NUMBINS;i++){
848 for(int j=0;j<NUMBINS;j++){
850 double minSum = 100000000;
852 for(int k=0;k<HOMOPS;k++){
853 double sum = singleLookUp[k * NUMBINS + i] + singleLookUp[k * NUMBINS + j];
855 if(sum < minSum) { minSum = sum; }
857 jointLookUp[i * NUMBINS + j] = minSum;
861 catch(exception& e) {
862 m->errorOut(e, "ShhherCommand", "getJointLookUp");
867 /**************************************************************************************************/
869 double ShhherCommand::getProbIntensity(int intIntensity){
872 double minNegLogProb = 100000000;
875 for(int i=0;i<HOMOPS;i++){//loop signal strength
876 float negLogProb = singleLookUp[i * NUMBINS + intIntensity];
877 if(negLogProb < minNegLogProb) { minNegLogProb = negLogProb; }
880 return minNegLogProb;
882 catch(exception& e) {
883 m->errorOut(e, "ShhherCommand", "getProbIntensity");
888 /**************************************************************************************************/
890 void ShhherCommand::getUniques(){
895 uniqueFlowgrams.assign(numFlowCells * numSeqs, -1);
896 uniqueCount.assign(numSeqs, 0); // anWeights
897 uniqueLengths.assign(numSeqs, 0);
898 mapSeqToUnique.assign(numSeqs, -1);
899 mapUniqueToSeq.assign(numSeqs, -1);
901 vector<short> uniqueFlowDataIntI(numFlowCells * numSeqs, -1);
903 for(int i=0;i<numSeqs;i++){
906 vector<short> current(numFlowCells);
907 for(int j=0;j<numFlowCells;j++){
908 current[j] = short(((flowDataIntI[i * numFlowCells + j] + 50.0)/100.0));
911 for(int j=0;j<numUniques;j++){
912 int offset = j * numFlowCells;
916 if(lengths[i] < uniqueLengths[j]) { shorterLength = lengths[i]; }
917 else { shorterLength = uniqueLengths[j]; }
919 for(int k=0;k<shorterLength;k++){
920 if(current[k] != uniqueFlowgrams[offset + k]){
927 mapSeqToUnique[i] = j;
930 if(lengths[i] > uniqueLengths[j]) { uniqueLengths[j] = lengths[i]; }
936 if(index == numUniques){
937 uniqueLengths[numUniques] = lengths[i];
938 uniqueCount[numUniques] = 1;
939 mapSeqToUnique[i] = numUniques;//anMap
940 mapUniqueToSeq[numUniques] = i;//anF
942 for(int k=0;k<numFlowCells;k++){
943 uniqueFlowgrams[numUniques * numFlowCells + k] = current[k];
944 uniqueFlowDataIntI[numUniques * numFlowCells + k] = flowDataIntI[i * numFlowCells + k];
950 uniqueFlowDataIntI.resize(numFlowCells * numUniques);
951 uniqueLengths.resize(numUniques);
953 flowDataPrI.resize(numSeqs * numFlowCells, 0);
954 for(int i=0;i<flowDataPrI.size();i++) { flowDataPrI[i] = getProbIntensity(flowDataIntI[i]); }
956 catch(exception& e) {
957 m->errorOut(e, "ShhherCommand", "getUniques");
962 /**************************************************************************************************/
964 float ShhherCommand::calcPairwiseDist(int seqA, int seqB){
966 int minLength = lengths[mapSeqToUnique[seqA]];
967 if(lengths[seqB] < minLength){ minLength = lengths[mapSeqToUnique[seqB]]; }
969 int ANumFlowCells = seqA * numFlowCells;
970 int BNumFlowCells = seqB * numFlowCells;
974 for(int i=0;i<minLength;i++){
975 int flowAIntI = flowDataIntI[ANumFlowCells + i];
976 float flowAPrI = flowDataPrI[ANumFlowCells + i];
978 int flowBIntI = flowDataIntI[BNumFlowCells + i];
979 float flowBPrI = flowDataPrI[BNumFlowCells + i];
980 dist += jointLookUp[flowAIntI * NUMBINS + flowBIntI] - flowAPrI - flowBPrI;
983 dist /= (float) minLength;
986 catch(exception& e) {
987 m->errorOut(e, "ShhherCommand", "calcPairwiseDist");
992 /**************************************************************************************************/
994 void ShhherCommand::flowDistParentFork(string distFileName, int startSeq, int stopSeq){
997 ostringstream outStream;
998 outStream.setf(ios::fixed, ios::floatfield);
999 outStream.setf(ios::dec, ios::basefield);
1000 outStream.setf(ios::showpoint);
1001 outStream.precision(6);
1003 int begTime = time(NULL);
1004 double begClock = clock();
1006 for(int i=startSeq;i<stopSeq;i++){
1007 for(int j=0;j<i;j++){
1008 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
1010 if(flowDistance < 1e-6){
1011 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
1013 else if(flowDistance <= cutoff){
1014 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
1018 m->mothurOut(toString(i) + "\t" + toString(time(NULL) - begTime));
1019 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1020 m->mothurOutEndLine();
1023 m->mothurOut(toString(stopSeq-1) + "\t" + toString(time(NULL) - begTime));
1024 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1025 m->mothurOutEndLine();
1027 ofstream distFile(distFileName.c_str());
1028 distFile << outStream.str();
1031 catch(exception& e) {
1032 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
1037 /**************************************************************************************************/
1039 string ShhherCommand::createDistFile(int processors){
1041 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
1043 unsigned long int begTime = time(NULL);
1044 double begClock = clock();
1049 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1050 if(processors == 1) { flowDistParentFork(fDistFileName, 0, numUniques); }
1051 else{ //you have multiple processors
1053 if (numSeqs < processors){ processors = 1; }
1055 vector<int> start(processors, 0);
1056 vector<int> end(processors, 0);
1058 for (int i = 0; i < processors; i++) {
1059 start[i] = int(sqrt(float(i)/float(processors)) * numUniques);
1060 end[i] = int(sqrt(float(i+1)/float(processors)) * numUniques);
1064 vector<int> processIDs;
1066 //loop through and create all the processes you want
1067 while (process != processors) {
1071 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1073 }else if (pid == 0){
1074 flowDistParentFork(fDistFileName + toString(getpid()) + ".temp", start[process], end[process]);
1077 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1079 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1084 //parent does its part
1085 flowDistParentFork(fDistFileName, start[0], end[0]);
1087 //force parent to wait until all the processes are done
1088 for (int i=0;i<processIDs.size();i++) {
1089 int temp = processIDs[i];
1093 //append and remove temp files
1094 for (int i=0;i<processIDs.size();i++) {
1095 m->appendFiles((fDistFileName + toString(processIDs[i]) + ".temp"), fDistFileName);
1096 remove((fDistFileName + toString(processIDs[i]) + ".temp").c_str());
1102 flowDistParentFork(fDistFileName, 0, numUniques);
1105 m->mothurOutEndLine();
1107 m->mothurOut("Total time: " + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/CLOCKS_PER_SEC) + '\n');
1110 return fDistFileName;
1112 catch(exception& e) {
1113 m->errorOut(e, "ShhherCommand", "createDistFile");
1119 /**************************************************************************************************/
1121 string ShhherCommand::createNamesFile(){
1124 vector<string> duplicateNames(numUniques, "");
1125 for(int i=0;i<numSeqs;i++){
1126 duplicateNames[mapSeqToUnique[i]] += seqNameVector[i] + ',';
1129 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.names";
1132 m->openOutputFile(nameFileName, nameFile);
1134 for(int i=0;i<numUniques;i++){
1135 // nameFile << seqNameVector[mapUniqueToSeq[i]] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1136 nameFile << mapUniqueToSeq[i] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1140 return nameFileName;
1142 catch(exception& e) {
1143 m->errorOut(e, "ShhherCommand", "createNamesFile");
1148 //**********************************************************************************************************************
1150 string ShhherCommand::cluster(string distFileName, string namesFileName){
1154 globaldata->setNameFile(namesFileName);
1155 globaldata->setColumnFile(distFileName);
1156 globaldata->setFormat("column");
1158 ReadMatrix* read = new ReadColumnMatrix(distFileName);
1159 read->setCutoff(cutoff);
1161 NameAssignment* clusterNameMap = new NameAssignment(namesFileName);
1162 clusterNameMap->readMap();
1163 read->read(clusterNameMap);
1165 ListVector* list = read->getListVector();
1166 SparseMatrix* matrix = read->getMatrix();
1169 delete clusterNameMap;
1171 RAbundVector* rabund = new RAbundVector(list->getRAbundVector());
1173 Cluster* cluster = new CompleteLinkage(rabund, list, matrix, cutoff, "furthest");
1174 string tag = cluster->getTag();
1176 double clusterCutoff = cutoff;
1177 while (matrix->getSmallDist() <= clusterCutoff && matrix->getNNodes() > 0){
1178 cluster->update(clusterCutoff);
1181 list->setLabel(toString(cutoff));
1183 string listFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.list";
1185 m->openOutputFile(listFileName, listFile);
1186 list->print(listFile);
1189 delete matrix; delete cluster; delete rabund; delete list;
1191 return listFileName;
1193 catch(exception& e) {
1194 m->errorOut(e, "ShhherCommand", "cluster");
1199 /**************************************************************************************************/
1201 void ShhherCommand::getOTUData(string listFileName){
1205 m->openInputFile(listFileName, listFile);
1208 listFile >> label >> numOTUs;
1210 otuData.assign(numSeqs, 0);
1211 cumNumSeqs.assign(numOTUs, 0);
1212 nSeqsPerOTU.assign(numOTUs, 0);
1213 aaP.resize(numOTUs);
1219 string singleOTU = "";
1221 for(int i=0;i<numOTUs;i++){
1223 listFile >> singleOTU;
1225 istringstream otuString(singleOTU);
1229 string seqName = "";
1231 for(int j=0;j<singleOTU.length();j++){
1232 char letter = otuString.get();
1238 map<string,int>::iterator nmIt = nameMap.find(seqName);
1239 int index = nmIt->second;
1241 nameMap.erase(nmIt);
1245 aaP[i].push_back(index);
1250 map<string,int>::iterator nmIt = nameMap.find(seqName);
1252 int index = nmIt->second;
1253 nameMap.erase(nmIt);
1257 aaP[i].push_back(index);
1262 sort(aaP[i].begin(), aaP[i].end());
1263 for(int j=0;j<nSeqsPerOTU[i];j++){
1264 seqNumber.push_back(aaP[i][j]);
1266 for(int j=nSeqsPerOTU[i];j<numSeqs;j++){
1267 aaP[i].push_back(0);
1271 for(int i=1;i<numOTUs;i++){
1272 cumNumSeqs[i] = cumNumSeqs[i-1] + nSeqsPerOTU[i-1];
1275 seqIndex = seqNumber;
1279 catch(exception& e) {
1280 m->errorOut(e, "ShhherCommand", "getOTUData");
1285 /**************************************************************************************************/
1287 void ShhherCommand::initPyroCluster(){
1289 dist.assign(numSeqs * numOTUs, 0);
1290 change.assign(numOTUs, 1);
1291 centroids.assign(numOTUs, -1);
1292 weight.assign(numOTUs, 0);
1293 singleTau.assign(numSeqs, 1.0);
1295 nSeqsBreaks.assign(processors+1, 0);
1296 nOTUsBreaks.assign(processors+1, 0);
1299 for(int i=0;i<processors;i++){
1300 nSeqsBreaks[i+1] = nSeqsBreaks[i] + (int)((double) numSeqs / (double) processors);
1301 nOTUsBreaks[i+1] = nOTUsBreaks[i] + (int)((double) numOTUs / (double) processors);
1303 nSeqsBreaks[processors] = numSeqs;
1304 nOTUsBreaks[processors] = numOTUs;
1306 catch(exception& e) {
1307 m->errorOut(e, "ShhherCommand", "initPyroCluster");
1312 /**************************************************************************************************/
1314 void ShhherCommand::fill(){
1317 for(int i=0;i<numOTUs;i++){
1318 cumNumSeqs[i] = index;
1319 for(int j=0;j<nSeqsPerOTU[i];j++){
1320 seqNumber[index] = aaP[i][j];
1321 seqIndex[index] = aaI[i][j];
1327 catch(exception& e) {
1328 m->errorOut(e, "ShhherCommand", "fill");
1333 /**************************************************************************************************/
1335 void ShhherCommand::calcCentroids(){
1338 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1340 if(processors == 1) {
1341 calcCentroidsDriver(0, numOTUs);
1343 else{ //you have multiple processors
1344 if (numOTUs < processors){ processors = 1; }
1347 vector<int> processIDs;
1349 //loop through and create all the processes you want
1350 while (process != processors) {
1354 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1356 }else if (pid == 0){
1357 calcCentroidsDriver(nOTUsBreaks[process], nOTUsBreaks[process+1]);
1360 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1362 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1367 //parent does its part
1368 calcCentroidsDriver(nOTUsBreaks[0], nOTUsBreaks[1]);
1370 //force parent to wait until all the processes are done
1371 for (int i=0;i<processIDs.size();i++) {
1372 int temp = processIDs[i];
1378 calcCentroidsDriver(0, numOTUs);
1381 catch(exception& e) {
1382 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1387 /**************************************************************************************************/
1389 void ShhherCommand::calcCentroidsDriver(int start, int finish){
1391 //this function gets the most likely homopolymer length at a flow position for a group of sequences
1396 for(int i=start;i<finish;i++){
1400 int minFlowGram = 100000000;
1401 double minFlowValue = 1e8;
1402 change[i] = 0; //FALSE
1404 for(int j=0;j<nSeqsPerOTU[i];j++){
1405 count += singleTau[seqNumber[cumNumSeqs[i] + j]];
1408 if(nSeqsPerOTU[i] > 0 && count > MIN_COUNT){
1409 vector<double> adF(nSeqsPerOTU[i]);
1410 vector<int> anL(nSeqsPerOTU[i]);
1412 for(int j=0;j<nSeqsPerOTU[i];j++){
1413 int index = cumNumSeqs[i] + j;
1414 int nI = seqIndex[index];
1415 int nIU = mapSeqToUnique[nI];
1418 for(k=0;k<position;k++){
1424 anL[position] = nIU;
1425 adF[position] = 0.0000;
1430 for(int j=0;j<nSeqsPerOTU[i];j++){
1431 int index = cumNumSeqs[i] + j;
1432 int nI = seqIndex[index];
1434 double tauValue = singleTau[seqNumber[index]];
1436 for(int k=0;k<position;k++){
1437 double dist = getDistToCentroid(anL[k], nI, lengths[nI]);
1438 adF[k] += dist * tauValue;
1442 for(int j=0;j<position;j++){
1443 if(adF[j] < minFlowValue){
1445 minFlowValue = adF[j];
1449 if(centroids[i] != anL[minFlowGram]){
1451 centroids[i] = anL[minFlowGram];
1454 else if(centroids[i] != -1){
1460 catch(exception& e) {
1461 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1466 /**************************************************************************************************/
1468 double ShhherCommand::getDistToCentroid(int cent, int flow, int length){
1471 int flowAValue = cent * numFlowCells;
1472 int flowBValue = flow * numFlowCells;
1476 for(int i=0;i<length;i++){
1477 dist += singleLookUp[uniqueFlowgrams[flowAValue] * NUMBINS + flowDataIntI[flowBValue]];
1482 return dist / (double)length;
1484 catch(exception& e) {
1485 m->errorOut(e, "ShhherCommand", "getDistToCentroid");
1490 /**************************************************************************************************/
1492 double ShhherCommand::getNewWeights(){
1495 double maxChange = 0;
1497 for(int i=0;i<numOTUs;i++){
1499 double difference = weight[i];
1502 for(int j=0;j<nSeqsPerOTU[i];j++){
1503 int index = cumNumSeqs[i] + j;
1504 double tauValue = singleTau[seqNumber[index]];
1505 weight[i] += tauValue;
1508 difference = fabs(weight[i] - difference);
1509 if(difference > maxChange){ maxChange = difference; }
1513 catch(exception& e) {
1514 m->errorOut(e, "ShhherCommand", "getNewWeights");
1519 /**************************************************************************************************/
1521 double ShhherCommand::getLikelihood(){
1525 vector<long double> P(numSeqs, 0);
1528 for(int i=0;i<numOTUs;i++){
1529 if(weight[i] > MIN_WEIGHT){
1535 for(int i=0;i<numOTUs;i++){
1536 for(int j=0;j<nSeqsPerOTU[i];j++){
1537 int index = cumNumSeqs[i] + j;
1538 int nI = seqIndex[index];
1539 double singleDist = dist[seqNumber[index]];
1541 P[nI] += weight[i] * exp(-singleDist * sigma);
1545 for(int i=0;i<numSeqs;i++){
1546 if(P[i] == 0){ P[i] = DBL_EPSILON; }
1551 nLL = nLL -(double)numSeqs * log(sigma);
1555 catch(exception& e) {
1556 m->errorOut(e, "ShhherCommand", "getNewWeights");
1561 /**************************************************************************************************/
1563 void ShhherCommand::checkCentroids(){
1565 vector<int> unique(numOTUs, 1);
1567 for(int i=0;i<numOTUs;i++){
1568 if(centroids[i] == -1 || weight[i] < MIN_WEIGHT){
1573 for(int i=0;i<numOTUs;i++){
1575 for(int j=i+1;j<numOTUs;j++){
1578 if(centroids[j] == centroids[i]){
1582 weight[i] += weight[j];
1590 catch(exception& e) {
1591 m->errorOut(e, "ShhherCommand", "checkCentroids");
1596 /**************************************************************************************************/
1598 void ShhherCommand::calcNewDistances(){
1601 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1603 if(processors == 1) {
1604 calcNewDistancesParent(0, numSeqs);
1606 else{ //you have multiple processors
1607 if (numSeqs < processors){ processors = 1; }
1609 vector<vector<int> > child_otuIndex(processors);
1610 vector<vector<int> > child_seqIndex(processors);
1611 vector<vector<double> > child_singleTau(processors);
1612 vector<int> totals(processors);
1615 vector<int> processIDs;
1617 //loop through and create all the processes you want
1618 while (process != processors) {
1622 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1624 }else if (pid == 0){
1625 calcNewDistancesChild(nSeqsBreaks[process], nSeqsBreaks[process+1], child_otuIndex[process], child_seqIndex[process], child_singleTau[process]);
1626 totals[process] = child_otuIndex[process].size();
1630 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1632 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1637 //parent does its part
1638 calcNewDistancesParent(nSeqsBreaks[0], nSeqsBreaks[1]);
1639 int total = seqIndex.size();
1641 //force parent to wait until all the processes are done
1642 for (int i=0;i<processIDs.size();i++) {
1643 int temp = processIDs[i];
1647 for(int i=1;i<processors;i++){
1648 int oldTotal = total;
1651 singleTau.resize(total, 0);
1652 seqIndex.resize(total, 0);
1653 seqNumber.resize(total, 0);
1657 for(int j=oldTotal;j<total;j++){
1658 int otuI = child_otuIndex[i][childIndex];
1659 int seqI = child_seqIndex[i][childIndex];
1661 singleTau[j] = child_singleTau[i][childIndex];
1662 aaP[otuI][nSeqsPerOTU[otuI]] = j;
1663 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
1664 nSeqsPerOTU[otuI]++;
1671 calcNewDistancesParent(0, numSeqs);
1674 catch(exception& e) {
1675 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1680 /**************************************************************************************************/
1682 void ShhherCommand::calcNewDistancesChildMPI(int startSeq, int stopSeq, vector<int>& otuIndex){
1685 vector<double> newTau(numOTUs,0);
1686 vector<double> norms(numSeqs, 0);
1693 for(int i=startSeq;i<stopSeq;i++){
1694 double offset = 1e8;
1695 int indexOffset = i * numOTUs;
1697 for(int j=0;j<numOTUs;j++){
1699 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1700 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1702 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1703 offset = dist[indexOffset + j];
1707 for(int j=0;j<numOTUs;j++){
1708 if(weight[j] > MIN_WEIGHT){
1709 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1710 norms[i] += newTau[j];
1717 for(int j=0;j<numOTUs;j++){
1719 newTau[j] /= norms[i];
1721 if(newTau[j] > MIN_TAU){
1722 otuIndex.push_back(j);
1723 seqIndex.push_back(i);
1724 singleTau.push_back(newTau[j]);
1730 catch(exception& e) {
1731 m->errorOut(e, "ShhherCommand", "calcNewDistancesChildMPI");
1736 /**************************************************************************************************/
1738 void ShhherCommand::calcNewDistancesChild(int startSeq, int stopSeq, vector<int>& child_otuIndex, vector<int>& child_seqIndex, vector<double>& child_singleTau){
1741 vector<double> newTau(numOTUs,0);
1742 vector<double> norms(numSeqs, 0);
1743 child_otuIndex.resize(0);
1744 child_seqIndex.resize(0);
1745 child_singleTau.resize(0);
1747 for(int i=startSeq;i<stopSeq;i++){
1748 double offset = 1e8;
1749 int indexOffset = i * numOTUs;
1752 for(int j=0;j<numOTUs;j++){
1753 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1754 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1757 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1758 offset = dist[indexOffset + j];
1762 for(int j=0;j<numOTUs;j++){
1763 if(weight[j] > MIN_WEIGHT){
1764 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1765 norms[i] += newTau[j];
1772 for(int j=0;j<numOTUs;j++){
1773 newTau[j] /= norms[i];
1775 if(newTau[j] > MIN_TAU){
1776 child_otuIndex.push_back(j);
1777 child_seqIndex.push_back(i);
1778 child_singleTau.push_back(newTau[j]);
1783 catch(exception& e) {
1784 m->errorOut(e, "ShhherCommand", "calcNewDistancesChild");
1789 /**************************************************************************************************/
1791 void ShhherCommand::calcNewDistancesParent(int startSeq, int stopSeq){
1796 vector<double> newTau(numOTUs,0);
1797 vector<double> norms(numSeqs, 0);
1798 nSeqsPerOTU.assign(numOTUs, 0);
1800 for(int i=startSeq;i<stopSeq;i++){
1801 int indexOffset = i * numOTUs;
1803 double offset = 1e8;
1805 for(int j=0;j<numOTUs;j++){
1806 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1807 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1810 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1811 offset = dist[indexOffset + j];
1815 for(int j=0;j<numOTUs;j++){
1816 if(weight[j] > MIN_WEIGHT){
1817 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1818 norms[i] += newTau[j];
1825 for(int j=0;j<numOTUs;j++){
1826 newTau[j] /= norms[i];
1829 for(int j=0;j<numOTUs;j++){
1830 if(newTau[j] > MIN_TAU){
1832 int oldTotal = total;
1836 singleTau.resize(total, 0);
1837 seqNumber.resize(total, 0);
1838 seqIndex.resize(total, 0);
1840 singleTau[oldTotal] = newTau[j];
1842 aaP[j][nSeqsPerOTU[j]] = oldTotal;
1843 aaI[j][nSeqsPerOTU[j]] = i;
1849 catch(exception& e) {
1850 m->errorOut(e, "ShhherCommand", "calcNewDistancesParent");
1855 /**************************************************************************************************/
1857 void ShhherCommand::setOTUs(){
1860 vector<double> bigTauMatrix(numOTUs * numSeqs, 0.0000);
1862 for(int i=0;i<numOTUs;i++){
1863 for(int j=0;j<nSeqsPerOTU[i];j++){
1864 int index = cumNumSeqs[i] + j;
1865 double tauValue = singleTau[seqNumber[index]];
1866 int sIndex = seqIndex[index];
1867 bigTauMatrix[sIndex * numOTUs + i] = tauValue;
1871 for(int i=0;i<numSeqs;i++){
1872 double maxTau = -1.0000;
1874 for(int j=0;j<numOTUs;j++){
1875 if(bigTauMatrix[i * numOTUs + j] > maxTau){
1876 maxTau = bigTauMatrix[i * numOTUs + j];
1881 otuData[i] = maxOTU;
1884 nSeqsPerOTU.assign(numOTUs, 0);
1886 for(int i=0;i<numSeqs;i++){
1887 int index = otuData[i];
1889 singleTau[i] = 1.0000;
1892 aaP[index][nSeqsPerOTU[index]] = i;
1893 aaI[index][nSeqsPerOTU[index]] = i;
1895 nSeqsPerOTU[index]++;
1899 catch(exception& e) {
1900 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1905 /**************************************************************************************************/
1907 void ShhherCommand::writeQualities(vector<int> otuCounts){
1910 string qualityFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.qual";
1912 ofstream qualityFile;
1913 m->openOutputFile(qualityFileName, qualityFile);
1915 qualityFile.setf(ios::fixed, ios::floatfield);
1916 qualityFile.setf(ios::showpoint);
1917 qualityFile << setprecision(6);
1919 vector<vector<int> > qualities(numOTUs);
1920 vector<double> pr(HOMOPS, 0);
1923 for(int i=0;i<numOTUs;i++){
1927 if(nSeqsPerOTU[i] > 0){
1928 qualities[i].assign(1024, -1);
1930 while(index < numFlowCells){
1931 double maxPrValue = 1e8;
1932 short maxPrIndex = -1;
1933 double count = 0.0000;
1935 pr.assign(HOMOPS, 0);
1937 for(int j=0;j<nSeqsPerOTU[i];j++){
1938 int lIndex = cumNumSeqs[i] + j;
1939 double tauValue = singleTau[seqNumber[lIndex]];
1940 int sequenceIndex = aaI[i][j];
1941 short intensity = flowDataIntI[sequenceIndex * numFlowCells + index];
1945 for(int s=0;s<HOMOPS;s++){
1946 pr[s] += tauValue * singleLookUp[s * NUMBINS + intensity];
1950 maxPrIndex = uniqueFlowgrams[centroids[i] * numFlowCells + index];
1951 maxPrValue = pr[maxPrIndex];
1953 if(count > MIN_COUNT){
1955 double norm = 0.0000;
1957 for(int s=0;s<HOMOPS;s++){
1958 norm += exp(-(pr[s] - maxPrValue));
1961 for(int s=1;s<=maxPrIndex;s++){
1963 double temp = 0.0000;
1965 U += exp(-(pr[s-1]-maxPrValue))/norm;
1973 temp = floor(-10 * temp);
1974 value = (int)floor(temp);
1975 if(value > 100){ value = 100; }
1977 qualities[i][base] = (int)value;
1987 if(otuCounts[i] > 0){
1988 qualityFile << '>' << seqNameVector[mapUniqueToSeq[i]] << endl;
1990 int j=4; //need to get past the first four bases
1991 while(qualities[i][j] != -1){
1992 qualityFile << qualities[i][j] << ' ';
1995 qualityFile << endl;
1998 qualityFile.close();
2001 catch(exception& e) {
2002 m->errorOut(e, "ShhherCommand", "writeQualities");
2007 /**************************************************************************************************/
2009 void ShhherCommand::writeSequences(vector<int> otuCounts){
2011 string bases = "TACG";
2013 string fastaFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.fasta";
2015 m->openOutputFile(fastaFileName, fastaFile);
2017 vector<string> names(numOTUs, "");
2019 for(int i=0;i<numOTUs;i++){
2020 int index = centroids[i];
2022 if(otuCounts[i] > 0){
2023 fastaFile << '>' << seqNameVector[aaI[i][0]] << endl;
2025 for(int j=8;j<numFlowCells;j++){
2027 char base = bases[j % 4];
2028 for(int k=0;k<uniqueFlowgrams[index * numFlowCells + j];k++){
2037 if(compositeFASTAFileName != ""){
2038 m->appendFiles(fastaFileName, compositeFASTAFileName);
2041 catch(exception& e) {
2042 m->errorOut(e, "ShhherCommand", "writeSequences");
2047 /**************************************************************************************************/
2049 void ShhherCommand::writeNames(vector<int> otuCounts){
2051 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.final.names";
2053 m->openOutputFile(nameFileName, nameFile);
2055 for(int i=0;i<numOTUs;i++){
2056 if(otuCounts[i] > 0){
2057 nameFile << seqNameVector[aaI[i][0]] << '\t' << seqNameVector[aaI[i][0]];
2059 for(int j=1;j<nSeqsPerOTU[i];j++){
2060 nameFile << ',' << seqNameVector[aaI[i][j]];
2068 catch(exception& e) {
2069 m->errorOut(e, "ShhherCommand", "writeNames");
2074 /**************************************************************************************************/
2076 void ShhherCommand::writeGroups(){
2078 string fileRoot = flowFileName.substr(0,flowFileName.find_last_of('.'));
2079 string groupFileName = fileRoot + ".pn.groups";
2081 m->openOutputFile(groupFileName, groupFile);
2083 for(int i=0;i<numSeqs;i++){
2084 groupFile << seqNameVector[i] << '\t' << fileRoot << endl;
2088 catch(exception& e) {
2089 m->errorOut(e, "ShhherCommand", "writeGroups");
2094 /**************************************************************************************************/
2096 void ShhherCommand::writeClusters(vector<int> otuCounts){
2098 string otuCountsFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.counts";
2099 ofstream otuCountsFile;
2100 m->openOutputFile(otuCountsFileName, otuCountsFile);
2102 string bases = "TACG";
2104 for(int i=0;i<numOTUs;i++){
2105 //output the translated version of the centroid sequence for the otu
2106 if(otuCounts[i] > 0){
2107 int index = centroids[i];
2109 otuCountsFile << "ideal\t";
2110 for(int j=8;j<numFlowCells;j++){
2111 char base = bases[j % 4];
2112 for(int s=0;s<uniqueFlowgrams[index * numFlowCells + j];s++){
2113 otuCountsFile << base;
2116 otuCountsFile << endl;
2118 for(int j=0;j<nSeqsPerOTU[i];j++){
2119 int sequence = aaI[i][j];
2120 otuCountsFile << seqNameVector[sequence] << '\t';
2122 for(int k=8;k<lengths[sequence];k++){
2123 char base = bases[k % 4];
2124 int freq = int(0.01 * (double)flowDataIntI[sequence * numFlowCells + k] + 0.5);
2126 for(int s=0;s<freq;s++){
2127 otuCountsFile << base;
2130 otuCountsFile << endl;
2132 otuCountsFile << endl;
2135 otuCountsFile.close();
2137 catch(exception& e) {
2138 m->errorOut(e, "ShhherCommand", "writeClusters");
2143 //**********************************************************************************************************************