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; }
255 for(int i=1;i<ncpus;i++){
256 MPI_Send(&abort, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
258 if(abort == 1){ return 0; }
262 m->mothurOut("\nGetting preliminary data...\n");
266 vector<string> flowFileVector;
267 if(flowFilesFileName != "not found"){
270 ifstream flowFilesFile;
271 m->openInputFile(flowFilesFileName, flowFilesFile);
272 while(flowFilesFile){
273 flowFilesFile >> fName;
274 flowFileVector.push_back(fName);
275 m->gobble(flowFilesFile);
279 flowFileVector.push_back(flowFileName);
281 int numFiles = flowFileVector.size();
283 for(int i=1;i<ncpus;i++){
284 MPI_Send(&numFiles, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
287 for(int i=0;i<numFiles;i++){
288 double begClock = clock();
289 unsigned long int begTime = time(NULL);
291 flowFileName = flowFileVector[i];
293 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
294 m->mothurOut("Reading flowgrams...\n");
297 m->mothurOut("Identifying unique flowgrams...\n");
300 m->mothurOut("Calculating distances between flowgrams...\n");
302 strcpy(fileName, flowFileName.c_str());
304 for(int i=1;i<ncpus;i++){
305 MPI_Send(&fileName[0], 1024, MPI_CHAR, i, tag, MPI_COMM_WORLD);
307 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
308 MPI_Send(&numUniques, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
309 MPI_Send(&numFlowCells, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
310 MPI_Send(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, i, tag, MPI_COMM_WORLD);
311 MPI_Send(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
312 MPI_Send(&mapUniqueToSeq[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
313 MPI_Send(&mapSeqToUnique[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
314 MPI_Send(&lengths[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
315 MPI_Send(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
316 MPI_Send(&cutoff, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
319 string distFileName = flowDistMPI(0, int(sqrt(1.0/float(ncpus)) * numUniques));
322 for(int i=1;i<ncpus;i++){
323 MPI_Recv(&done, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
325 m->appendFiles((distFileName + ".temp." + toString(i)), distFileName);
326 remove((distFileName + ".temp." + toString(i)).c_str());
329 string namesFileName = createNamesFile();
331 m->mothurOut("\nClustering flowgrams...\n");
332 string listFileName = cluster(distFileName, namesFileName);
334 getOTUData(listFileName);
337 for(int i=1;i<ncpus;i++){
338 MPI_Send(&numOTUs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
339 MPI_Send(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
340 MPI_Send(&uniqueFlowgrams[0], numFlowCells * numUniques, MPI_SHORT, i, tag, MPI_COMM_WORLD);
341 MPI_Send(&sigma, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
348 int numOTUsOnCPU = numOTUs / ncpus;
349 int numSeqsOnCPU = numSeqs / ncpus;
350 m->mothurOut("\nDenoising flowgrams...\n");
351 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
353 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
355 double cycClock = clock();
356 unsigned long int cycTime = time(NULL);
359 int total = singleTau.size();
360 for(int i=1;i<ncpus;i++){
361 MPI_Send(&total, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
362 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
363 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
365 MPI_Send(&singleTau[0], total, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
366 MPI_Send(&seqNumber[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
367 MPI_Send(&seqIndex[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
368 MPI_Send(&nSeqsPerOTU[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
369 MPI_Send(&cumNumSeqs[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
372 calcCentroidsDriver(0, numOTUsOnCPU);
374 for(int i=1;i<ncpus;i++){
375 int otuStart = i * numOTUs / ncpus;
376 int otuStop = (i + 1) * numOTUs / ncpus;
378 vector<int> tempCentroids(numOTUs, 0);
379 vector<short> tempChange(numOTUs, 0);
381 MPI_Recv(&tempCentroids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
382 MPI_Recv(&tempChange[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD, &status);
384 for(int j=otuStart;j<otuStop;j++){
385 centroids[j] = tempCentroids[j];
386 change[j] = tempChange[j];
390 maxDelta = getNewWeights();
391 double nLL = getLikelihood();
394 for(int i=1;i<ncpus;i++){
395 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
396 MPI_Send(&weight[0], numOTUs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
397 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
400 calcNewDistancesParent(0, numSeqsOnCPU);
402 total = singleTau.size();
404 for(int i=1;i<ncpus;i++){
406 int seqStart = i * numSeqs / ncpus;
407 int seqStop = (i + 1) * numSeqs / ncpus;
409 MPI_Recv(&childTotal, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
411 vector<int> childSeqIndex(childTotal, 0);
412 vector<double> childSingleTau(childTotal, 0);
413 vector<double> childDist(numSeqs * numOTUs, 0);
414 vector<int> otuIndex(childTotal, 0);
416 MPI_Recv(&childSeqIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
417 MPI_Recv(&childSingleTau[0], childTotal, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
418 MPI_Recv(&childDist[0], numOTUs * numSeqs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
419 MPI_Recv(&otuIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
421 int oldTotal = total;
423 singleTau.resize(total, 0);
424 seqIndex.resize(total, 0);
425 seqNumber.resize(total, 0);
429 for(int j=oldTotal;j<total;j++){
430 int otuI = otuIndex[childIndex];
431 int seqI = childSeqIndex[childIndex];
433 singleTau[j] = childSingleTau[childIndex];
435 aaP[otuI][nSeqsPerOTU[otuI]] = j;
436 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
441 int index = seqStart * numOTUs;
442 for(int j=seqStart;j<seqStop;j++){
443 for(int k=0;k<numOTUs;k++){
444 dist[index] = childDist[index];
452 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
454 if((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
456 for(int i=1;i<ncpus;i++){
457 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
462 for(int i=1;i<ncpus;i++){
463 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD); //send kill command
469 m->mothurOut("\nFinalizing...\n");
472 vector<int> otuCounts(numOTUs, 0);
473 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
474 calcCentroidsDriver(0, numOTUs);
475 writeQualities(otuCounts);
476 writeSequences(otuCounts);
477 writeNames(otuCounts);
478 writeClusters(otuCounts);
481 remove(distFileName.c_str());
482 remove(namesFileName.c_str());
483 remove(listFileName.c_str());
485 m->mothurOut("Total time to process " + toString(flowFileName) + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
491 MPI_Recv(&abort, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
492 if(abort){ return 0; }
495 MPI_Recv(&numFiles, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
497 for(int i=0;i<numFiles;i++){
498 //Now into the pyrodist part
502 MPI_Recv(&fileName, 1024, MPI_CHAR, 0, tag, MPI_COMM_WORLD, &status);
503 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
504 MPI_Recv(&numUniques, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
505 MPI_Recv(&numFlowCells, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
507 flowDataIntI.resize(numSeqs * numFlowCells);
508 flowDataPrI.resize(numSeqs * numFlowCells);
509 mapUniqueToSeq.resize(numSeqs);
510 mapSeqToUnique.resize(numSeqs);
511 lengths.resize(numSeqs);
512 jointLookUp.resize(NUMBINS * NUMBINS);
514 MPI_Recv(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
515 MPI_Recv(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
516 MPI_Recv(&mapUniqueToSeq[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
517 MPI_Recv(&mapSeqToUnique[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
518 MPI_Recv(&lengths[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
519 MPI_Recv(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
520 MPI_Recv(&cutoff, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
522 flowFileName = string(fileName);
523 int flowDistStart = int(sqrt(float(pid)/float(ncpus)) * numUniques);
524 int flowDistEnd = int(sqrt(float(pid+1)/float(ncpus)) * numUniques);
526 string distanceStringChild = flowDistMPI(flowDistStart, flowDistEnd);
529 MPI_Send(&done, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
531 //Now into the pyronoise part
532 MPI_Recv(&numOTUs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
534 singleLookUp.resize(HOMOPS * NUMBINS);
535 uniqueFlowgrams.resize(numUniques * numFlowCells);
536 weight.resize(numOTUs);
537 centroids.resize(numOTUs);
538 change.resize(numOTUs);
539 dist.assign(numOTUs * numSeqs, 0);
540 nSeqsPerOTU.resize(numOTUs);
541 cumNumSeqs.resize(numOTUs);
543 MPI_Recv(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
544 MPI_Recv(&uniqueFlowgrams[0], uniqueFlowgrams.size(), MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
545 MPI_Recv(&sigma, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
547 int startOTU = pid * numOTUs / ncpus;
548 int endOTU = (pid + 1) * numOTUs / ncpus;
550 int startSeq = pid * numSeqs / ncpus;
551 int endSeq = (pid + 1) * numSeqs /ncpus;
557 MPI_Recv(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
558 singleTau.assign(total, 0.0000);
559 seqNumber.assign(total, 0);
560 seqIndex.assign(total, 0);
562 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
563 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
564 MPI_Recv(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
565 MPI_Recv(&seqNumber[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
566 MPI_Recv(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
567 MPI_Recv(&nSeqsPerOTU[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
568 MPI_Recv(&cumNumSeqs[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
570 calcCentroidsDriver(startOTU, endOTU);
572 MPI_Send(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD);
573 MPI_Send(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD);
575 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
576 MPI_Recv(&weight[0], numOTUs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
577 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
579 vector<int> otuIndex(total, 0);
580 calcNewDistancesChildMPI(startSeq, endSeq, otuIndex);
581 total = otuIndex.size();
583 MPI_Send(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
584 MPI_Send(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
585 MPI_Send(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
586 MPI_Send(&dist[0], numOTUs * numSeqs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
587 MPI_Send(&otuIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
589 MPI_Recv(&live, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
593 MPI_Barrier(MPI_COMM_WORLD);
598 catch(exception& e) {
599 m->errorOut(e, "ShhherCommand", "execute");
604 /**************************************************************************************************/
606 string ShhherCommand::flowDistMPI(int startSeq, int stopSeq){
608 ostringstream outStream;
609 outStream.setf(ios::fixed, ios::floatfield);
610 outStream.setf(ios::dec, ios::basefield);
611 outStream.setf(ios::showpoint);
612 outStream.precision(6);
614 int begTime = time(NULL);
615 double begClock = clock();
617 for(int i=startSeq;i<stopSeq;i++){
618 for(int j=0;j<i;j++){
619 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
621 if(flowDistance < 1e-6){
622 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
624 else if(flowDistance <= cutoff){
625 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
629 m->mothurOut(toString(i) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
633 m->mothurOut(toString(stopSeq) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
635 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
636 if(pid != 0){ fDistFileName += ".temp." + toString(pid); }
638 ofstream distFile(fDistFileName.c_str());
639 distFile << outStream.str();
642 return fDistFileName;
644 catch(exception& e) {
645 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
651 //**********************************************************************************************************************
653 int ShhherCommand::execute(){
655 if (abort == true) { return 0; }
660 vector<string> flowFileVector;
661 if(flowFilesFileName != "not found"){
664 ifstream flowFilesFile;
665 m->openInputFile(flowFilesFileName, flowFilesFile);
666 while(flowFilesFile){
667 flowFilesFile >> fName;
668 flowFileVector.push_back(fName);
669 m->gobble(flowFilesFile);
673 flowFileVector.push_back(flowFileName);
675 int numFiles = flowFileVector.size();
678 for(int i=0;i<numFiles;i++){
679 flowFileName = flowFileVector[i];
681 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
682 m->mothurOut("Reading flowgrams...\n");
685 m->mothurOut("Identifying unique flowgrams...\n");
689 m->mothurOut("Calculating distances between flowgrams...\n");
690 string distFileName = createDistFile(processors);
691 string namesFileName = createNamesFile();
693 m->mothurOut("\nClustering flowgrams...\n");
694 string listFileName = cluster(distFileName, namesFileName);
695 getOTUData(listFileName);
702 double begClock = clock();
703 unsigned long int begTime = time(NULL);
705 m->mothurOut("\nDenoising flowgrams...\n");
706 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
708 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
710 double cycClock = clock();
711 unsigned long int cycTime = time(NULL);
716 maxDelta = getNewWeights();
717 double nLL = getLikelihood();
724 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
728 m->mothurOut("\nFinalizing...\n");
732 vector<int> otuCounts(numOTUs, 0);
733 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
735 calcCentroidsDriver(0, numOTUs);
736 writeQualities(otuCounts);
737 writeSequences(otuCounts);
738 writeNames(otuCounts);
739 writeClusters(otuCounts);
742 remove(distFileName.c_str());
743 remove(namesFileName.c_str());
744 remove(listFileName.c_str());
746 m->mothurOut("Total time to process " + flowFileName + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
750 catch(exception& e) {
751 m->errorOut(e, "ShhherCommand", "execute");
756 /**************************************************************************************************/
758 void ShhherCommand::getFlowData(){
761 m->openInputFile(flowFileName, flowFile);
764 seqNameVector.clear();
766 flowDataIntI.clear();
770 int currentNumFlowCells;
774 flowFile >> numFlowCells;
775 int index = 0;//pcluster
776 while(!flowFile.eof()){
777 flowFile >> seqName >> currentNumFlowCells;
778 lengths.push_back(currentNumFlowCells);
780 seqNameVector.push_back(seqName);
781 nameMap[seqName] = index++;//pcluster
783 for(int i=0;i<numFlowCells;i++){
784 flowFile >> intensity;
785 if(intensity > 9.99) { intensity = 9.99; }
786 int intI = int(100 * intensity + 0.0001);
787 flowDataIntI.push_back(intI);
793 numSeqs = seqNameVector.size();
795 for(int i=0;i<numSeqs;i++){
796 int iNumFlowCells = i * numFlowCells;
797 for(int j=lengths[i];j<numFlowCells;j++){
798 flowDataIntI[iNumFlowCells + j] = 0;
803 catch(exception& e) {
804 m->errorOut(e, "ShhherCommand", "getFlowData");
809 /**************************************************************************************************/
811 void ShhherCommand::getSingleLookUp(){
813 // these are the -log probabilities that a signal corresponds to a particular homopolymer length
814 singleLookUp.assign(HOMOPS * NUMBINS, 0);
818 m->openInputFile(lookupFileName, lookUpFile);
820 for(int i=0;i<HOMOPS;i++){
822 lookUpFile >> logFracFreq;
824 for(int j=0;j<NUMBINS;j++) {
825 lookUpFile >> singleLookUp[index];
831 catch(exception& e) {
832 m->errorOut(e, "ShhherCommand", "getSingleLookUp");
837 /**************************************************************************************************/
839 void ShhherCommand::getJointLookUp(){
842 // the most likely joint probability (-log) that two intenities have the same polymer length
843 jointLookUp.resize(NUMBINS * NUMBINS, 0);
845 for(int i=0;i<NUMBINS;i++){
846 for(int j=0;j<NUMBINS;j++){
848 double minSum = 100000000;
850 for(int k=0;k<HOMOPS;k++){
851 double sum = singleLookUp[k * NUMBINS + i] + singleLookUp[k * NUMBINS + j];
853 if(sum < minSum) { minSum = sum; }
855 jointLookUp[i * NUMBINS + j] = minSum;
859 catch(exception& e) {
860 m->errorOut(e, "ShhherCommand", "getJointLookUp");
865 /**************************************************************************************************/
867 double ShhherCommand::getProbIntensity(int intIntensity){
870 double minNegLogProb = 100000000;
873 for(int i=0;i<HOMOPS;i++){//loop signal strength
874 float negLogProb = singleLookUp[i * NUMBINS + intIntensity];
875 if(negLogProb < minNegLogProb) { minNegLogProb = negLogProb; }
878 return minNegLogProb;
880 catch(exception& e) {
881 m->errorOut(e, "ShhherCommand", "getProbIntensity");
886 /**************************************************************************************************/
888 void ShhherCommand::getUniques(){
893 uniqueFlowgrams.assign(numFlowCells * numSeqs, -1);
894 uniqueCount.assign(numSeqs, 0); // anWeights
895 uniqueLengths.assign(numSeqs, 0);
896 mapSeqToUnique.assign(numSeqs, -1);
897 mapUniqueToSeq.assign(numSeqs, -1);
899 vector<short> uniqueFlowDataIntI(numFlowCells * numSeqs, -1);
901 for(int i=0;i<numSeqs;i++){
904 vector<short> current(numFlowCells);
905 for(int j=0;j<numFlowCells;j++){
906 current[j] = short(((flowDataIntI[i * numFlowCells + j] + 50.0)/100.0));
909 for(int j=0;j<numUniques;j++){
910 int offset = j * numFlowCells;
914 if(lengths[i] < uniqueLengths[j]) { shorterLength = lengths[i]; }
915 else { shorterLength = uniqueLengths[j]; }
917 for(int k=0;k<shorterLength;k++){
918 if(current[k] != uniqueFlowgrams[offset + k]){
925 mapSeqToUnique[i] = j;
928 if(lengths[i] > uniqueLengths[j]) { uniqueLengths[j] = lengths[i]; }
934 if(index == numUniques){
935 uniqueLengths[numUniques] = lengths[i];
936 uniqueCount[numUniques] = 1;
937 mapSeqToUnique[i] = numUniques;//anMap
938 mapUniqueToSeq[numUniques] = i;//anF
940 for(int k=0;k<numFlowCells;k++){
941 uniqueFlowgrams[numUniques * numFlowCells + k] = current[k];
942 uniqueFlowDataIntI[numUniques * numFlowCells + k] = flowDataIntI[i * numFlowCells + k];
948 uniqueFlowDataIntI.resize(numFlowCells * numUniques);
949 uniqueLengths.resize(numUniques);
951 flowDataPrI.resize(numSeqs * numFlowCells, 0);
952 for(int i=0;i<flowDataPrI.size();i++) { flowDataPrI[i] = getProbIntensity(flowDataIntI[i]); }
954 catch(exception& e) {
955 m->errorOut(e, "ShhherCommand", "getUniques");
960 /**************************************************************************************************/
962 float ShhherCommand::calcPairwiseDist(int seqA, int seqB){
964 int minLength = lengths[mapSeqToUnique[seqA]];
965 if(lengths[seqB] < minLength){ minLength = lengths[mapSeqToUnique[seqB]]; }
967 int ANumFlowCells = seqA * numFlowCells;
968 int BNumFlowCells = seqB * numFlowCells;
972 for(int i=0;i<minLength;i++){
973 int flowAIntI = flowDataIntI[ANumFlowCells + i];
974 float flowAPrI = flowDataPrI[ANumFlowCells + i];
976 int flowBIntI = flowDataIntI[BNumFlowCells + i];
977 float flowBPrI = flowDataPrI[BNumFlowCells + i];
978 dist += jointLookUp[flowAIntI * NUMBINS + flowBIntI] - flowAPrI - flowBPrI;
981 dist /= (float) minLength;
984 catch(exception& e) {
985 m->errorOut(e, "ShhherCommand", "calcPairwiseDist");
990 /**************************************************************************************************/
992 void ShhherCommand::flowDistParentFork(string distFileName, int startSeq, int stopSeq){
995 ostringstream outStream;
996 outStream.setf(ios::fixed, ios::floatfield);
997 outStream.setf(ios::dec, ios::basefield);
998 outStream.setf(ios::showpoint);
999 outStream.precision(6);
1001 int begTime = time(NULL);
1002 double begClock = clock();
1004 for(int i=startSeq;i<stopSeq;i++){
1005 for(int j=0;j<i;j++){
1006 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
1008 if(flowDistance < 1e-6){
1009 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
1011 else if(flowDistance <= cutoff){
1012 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
1016 m->mothurOut(toString(i) + "\t" + toString(time(NULL) - begTime));
1017 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1018 m->mothurOutEndLine();
1021 m->mothurOut(toString(stopSeq-1) + "\t" + toString(time(NULL) - begTime));
1022 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1023 m->mothurOutEndLine();
1025 ofstream distFile(distFileName.c_str());
1026 distFile << outStream.str();
1029 catch(exception& e) {
1030 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
1035 /**************************************************************************************************/
1037 string ShhherCommand::createDistFile(int processors){
1039 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
1041 unsigned long int begTime = time(NULL);
1042 double begClock = clock();
1047 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1048 if(processors == 1) { flowDistParentFork(fDistFileName, 0, numUniques); }
1049 else{ //you have multiple processors
1051 if (numSeqs < processors){ processors = 1; }
1053 vector<int> start(processors, 0);
1054 vector<int> end(processors, 0);
1056 for (int i = 0; i < processors; i++) {
1057 start[i] = int(sqrt(float(i)/float(processors)) * numUniques);
1058 end[i] = int(sqrt(float(i+1)/float(processors)) * numUniques);
1062 vector<int> processIDs;
1064 //loop through and create all the processes you want
1065 while (process != processors) {
1069 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1071 }else if (pid == 0){
1072 flowDistParentFork(fDistFileName + toString(getpid()) + ".temp", start[process], end[process]);
1075 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1077 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1082 //parent does its part
1083 flowDistParentFork(fDistFileName, start[0], end[0]);
1085 //force parent to wait until all the processes are done
1086 for (int i=0;i<processIDs.size();i++) {
1087 int temp = processIDs[i];
1091 //append and remove temp files
1092 for (int i=0;i<processIDs.size();i++) {
1093 m->appendFiles((fDistFileName + toString(processIDs[i]) + ".temp"), fDistFileName);
1094 remove((fDistFileName + toString(processIDs[i]) + ".temp").c_str());
1100 flowDistParentFork(fDistFileName, 0, numUniques);
1103 m->mothurOutEndLine();
1105 m->mothurOut("Total time: " + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/CLOCKS_PER_SEC) + '\n');
1108 return fDistFileName;
1110 catch(exception& e) {
1111 m->errorOut(e, "ShhherCommand", "createDistFile");
1117 /**************************************************************************************************/
1119 string ShhherCommand::createNamesFile(){
1122 vector<string> duplicateNames(numUniques, "");
1123 for(int i=0;i<numSeqs;i++){
1124 duplicateNames[mapSeqToUnique[i]] += seqNameVector[i] + ',';
1127 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.names";
1130 m->openOutputFile(nameFileName, nameFile);
1132 for(int i=0;i<numUniques;i++){
1133 // nameFile << seqNameVector[mapUniqueToSeq[i]] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1134 nameFile << mapUniqueToSeq[i] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1138 return nameFileName;
1140 catch(exception& e) {
1141 m->errorOut(e, "ShhherCommand", "createNamesFile");
1146 //**********************************************************************************************************************
1148 string ShhherCommand::cluster(string distFileName, string namesFileName){
1152 globaldata->setNameFile(namesFileName);
1153 globaldata->setColumnFile(distFileName);
1154 globaldata->setFormat("column");
1156 ReadMatrix* read = new ReadColumnMatrix(distFileName);
1157 read->setCutoff(cutoff);
1159 NameAssignment* clusterNameMap = new NameAssignment(namesFileName);
1160 clusterNameMap->readMap();
1161 read->read(clusterNameMap);
1163 ListVector* list = read->getListVector();
1164 SparseMatrix* matrix = read->getMatrix();
1167 delete clusterNameMap;
1169 RAbundVector* rabund = new RAbundVector(list->getRAbundVector());
1171 Cluster* cluster = new CompleteLinkage(rabund, list, matrix, cutoff, "furthest");
1172 string tag = cluster->getTag();
1174 double clusterCutoff = cutoff;
1175 while (matrix->getSmallDist() <= clusterCutoff && matrix->getNNodes() > 0){
1176 cluster->update(clusterCutoff);
1179 list->setLabel(toString(cutoff));
1181 string listFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.list";
1183 m->openOutputFile(listFileName, listFile);
1184 list->print(listFile);
1187 delete matrix; delete cluster; delete rabund; delete list;
1189 return listFileName;
1191 catch(exception& e) {
1192 m->errorOut(e, "ShhherCommand", "cluster");
1197 /**************************************************************************************************/
1199 void ShhherCommand::getOTUData(string listFileName){
1203 m->openInputFile(listFileName, listFile);
1206 listFile >> label >> numOTUs;
1208 otuData.assign(numSeqs, 0);
1209 cumNumSeqs.assign(numOTUs, 0);
1210 nSeqsPerOTU.assign(numOTUs, 0);
1211 aaP.clear();aaP.resize(numOTUs);
1217 string singleOTU = "";
1219 for(int i=0;i<numOTUs;i++){
1221 listFile >> singleOTU;
1223 istringstream otuString(singleOTU);
1227 string seqName = "";
1229 for(int j=0;j<singleOTU.length();j++){
1230 char letter = otuString.get();
1236 map<string,int>::iterator nmIt = nameMap.find(seqName);
1237 int index = nmIt->second;
1239 nameMap.erase(nmIt);
1243 aaP[i].push_back(index);
1248 map<string,int>::iterator nmIt = nameMap.find(seqName);
1250 int index = nmIt->second;
1251 nameMap.erase(nmIt);
1255 aaP[i].push_back(index);
1260 sort(aaP[i].begin(), aaP[i].end());
1261 for(int j=0;j<nSeqsPerOTU[i];j++){
1262 seqNumber.push_back(aaP[i][j]);
1264 for(int j=nSeqsPerOTU[i];j<numSeqs;j++){
1265 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;
1280 catch(exception& e) {
1281 m->errorOut(e, "ShhherCommand", "getOTUData");
1286 /**************************************************************************************************/
1288 void ShhherCommand::initPyroCluster(){
1290 dist.assign(numSeqs * numOTUs, 0);
1291 change.assign(numOTUs, 1);
1292 centroids.assign(numOTUs, -1);
1293 weight.assign(numOTUs, 0);
1294 singleTau.assign(numSeqs, 1.0);
1296 nSeqsBreaks.assign(processors+1, 0);
1297 nOTUsBreaks.assign(processors+1, 0);
1300 for(int i=0;i<processors;i++){
1301 nSeqsBreaks[i+1] = nSeqsBreaks[i] + (int)((double) numSeqs / (double) processors);
1302 nOTUsBreaks[i+1] = nOTUsBreaks[i] + (int)((double) numOTUs / (double) processors);
1304 nSeqsBreaks[processors] = numSeqs;
1305 nOTUsBreaks[processors] = numOTUs;
1307 catch(exception& e) {
1308 m->errorOut(e, "ShhherCommand", "initPyroCluster");
1313 /**************************************************************************************************/
1315 void ShhherCommand::fill(){
1318 for(int i=0;i<numOTUs;i++){
1319 cumNumSeqs[i] = index;
1320 for(int j=0;j<nSeqsPerOTU[i];j++){
1321 seqNumber[index] = aaP[i][j];
1322 seqIndex[index] = aaI[i][j];
1328 catch(exception& e) {
1329 m->errorOut(e, "ShhherCommand", "fill");
1334 /**************************************************************************************************/
1336 void ShhherCommand::calcCentroids(){
1339 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1341 if(processors == 1) {
1342 calcCentroidsDriver(0, numOTUs);
1344 else{ //you have multiple processors
1345 if (numOTUs < processors){ processors = 1; }
1348 vector<int> processIDs;
1350 //loop through and create all the processes you want
1351 while (process != processors) {
1355 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1357 }else if (pid == 0){
1358 calcCentroidsDriver(nOTUsBreaks[process], nOTUsBreaks[process+1]);
1361 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1363 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1368 //parent does its part
1369 calcCentroidsDriver(nOTUsBreaks[0], nOTUsBreaks[1]);
1371 //force parent to wait until all the processes are done
1372 for (int i=0;i<processIDs.size();i++) {
1373 int temp = processIDs[i];
1379 calcCentroidsDriver(0, numOTUs);
1382 catch(exception& e) {
1383 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1388 /**************************************************************************************************/
1390 void ShhherCommand::calcCentroidsDriver(int start, int finish){
1392 //this function gets the most likely homopolymer length at a flow position for a group of sequences
1398 for(int i=start;i<finish;i++){
1402 int minFlowGram = 100000000;
1403 double minFlowValue = 1e8;
1404 change[i] = 0; //FALSE
1406 for(int j=0;j<nSeqsPerOTU[i];j++){
1407 count += singleTau[seqNumber[cumNumSeqs[i] + j]];
1410 if(nSeqsPerOTU[i] > 0 && count > MIN_COUNT){
1411 vector<double> adF(nSeqsPerOTU[i]);
1412 vector<int> anL(nSeqsPerOTU[i]);
1414 for(int j=0;j<nSeqsPerOTU[i];j++){
1415 int index = cumNumSeqs[i] + j;
1416 int nI = seqIndex[index];
1417 int nIU = mapSeqToUnique[nI];
1420 for(k=0;k<position;k++){
1426 anL[position] = nIU;
1427 adF[position] = 0.0000;
1432 for(int j=0;j<nSeqsPerOTU[i];j++){
1433 int index = cumNumSeqs[i] + j;
1434 int nI = seqIndex[index];
1436 double tauValue = singleTau[seqNumber[index]];
1438 for(int k=0;k<position;k++){
1439 double dist = getDistToCentroid(anL[k], nI, lengths[nI]);
1440 adF[k] += dist * tauValue;
1444 for(int j=0;j<position;j++){
1445 if(adF[j] < minFlowValue){
1447 minFlowValue = adF[j];
1451 if(centroids[i] != anL[minFlowGram]){
1453 centroids[i] = anL[minFlowGram];
1456 else if(centroids[i] != -1){
1462 catch(exception& e) {
1463 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1468 /**************************************************************************************************/
1470 double ShhherCommand::getDistToCentroid(int cent, int flow, int length){
1473 int flowAValue = cent * numFlowCells;
1474 int flowBValue = flow * numFlowCells;
1478 for(int i=0;i<length;i++){
1479 dist += singleLookUp[uniqueFlowgrams[flowAValue] * NUMBINS + flowDataIntI[flowBValue]];
1484 return dist / (double)length;
1486 catch(exception& e) {
1487 m->errorOut(e, "ShhherCommand", "getDistToCentroid");
1492 /**************************************************************************************************/
1494 double ShhherCommand::getNewWeights(){
1497 double maxChange = 0;
1499 for(int i=0;i<numOTUs;i++){
1501 double difference = weight[i];
1504 for(int j=0;j<nSeqsPerOTU[i];j++){
1505 int index = cumNumSeqs[i] + j;
1506 double tauValue = singleTau[seqNumber[index]];
1507 weight[i] += tauValue;
1510 difference = fabs(weight[i] - difference);
1511 if(difference > maxChange){ maxChange = difference; }
1515 catch(exception& e) {
1516 m->errorOut(e, "ShhherCommand", "getNewWeights");
1521 /**************************************************************************************************/
1523 double ShhherCommand::getLikelihood(){
1527 vector<long double> P(numSeqs, 0);
1530 for(int i=0;i<numOTUs;i++){
1531 if(weight[i] > MIN_WEIGHT){
1537 for(int i=0;i<numOTUs;i++){
1538 for(int j=0;j<nSeqsPerOTU[i];j++){
1539 int index = cumNumSeqs[i] + j;
1540 int nI = seqIndex[index];
1541 double singleDist = dist[seqNumber[index]];
1543 P[nI] += weight[i] * exp(-singleDist * sigma);
1547 for(int i=0;i<numSeqs;i++){
1548 if(P[i] == 0){ P[i] = DBL_EPSILON; }
1553 nLL = nLL -(double)numSeqs * log(sigma);
1557 catch(exception& e) {
1558 m->errorOut(e, "ShhherCommand", "getNewWeights");
1563 /**************************************************************************************************/
1565 void ShhherCommand::checkCentroids(){
1567 vector<int> unique(numOTUs, 1);
1569 for(int i=0;i<numOTUs;i++){
1570 if(centroids[i] == -1 || weight[i] < MIN_WEIGHT){
1575 for(int i=0;i<numOTUs;i++){
1577 for(int j=i+1;j<numOTUs;j++){
1580 if(centroids[j] == centroids[i]){
1584 weight[i] += weight[j];
1592 catch(exception& e) {
1593 m->errorOut(e, "ShhherCommand", "checkCentroids");
1598 /**************************************************************************************************/
1600 void ShhherCommand::calcNewDistances(){
1603 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1605 if(processors == 1) {
1606 calcNewDistancesParent(0, numSeqs);
1608 else{ //you have multiple processors
1609 if (numSeqs < processors){ processors = 1; }
1611 vector<vector<int> > child_otuIndex(processors);
1612 vector<vector<int> > child_seqIndex(processors);
1613 vector<vector<double> > child_singleTau(processors);
1614 vector<int> totals(processors);
1617 vector<int> processIDs;
1619 //loop through and create all the processes you want
1620 while (process != processors) {
1624 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1626 }else if (pid == 0){
1627 calcNewDistancesChild(nSeqsBreaks[process], nSeqsBreaks[process+1], child_otuIndex[process], child_seqIndex[process], child_singleTau[process]);
1628 totals[process] = child_otuIndex[process].size();
1632 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1634 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1639 //parent does its part
1640 calcNewDistancesParent(nSeqsBreaks[0], nSeqsBreaks[1]);
1641 int total = seqIndex.size();
1643 //force parent to wait until all the processes are done
1644 for (int i=0;i<processIDs.size();i++) {
1645 int temp = processIDs[i];
1649 for(int i=1;i<processors;i++){
1650 int oldTotal = total;
1653 singleTau.resize(total, 0);
1654 seqIndex.resize(total, 0);
1655 seqNumber.resize(total, 0);
1659 for(int j=oldTotal;j<total;j++){
1660 int otuI = child_otuIndex[i][childIndex];
1661 int seqI = child_seqIndex[i][childIndex];
1663 singleTau[j] = child_singleTau[i][childIndex];
1664 aaP[otuI][nSeqsPerOTU[otuI]] = j;
1665 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
1666 nSeqsPerOTU[otuI]++;
1673 calcNewDistancesParent(0, numSeqs);
1676 catch(exception& e) {
1677 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1682 /**************************************************************************************************/
1684 void ShhherCommand::calcNewDistancesChildMPI(int startSeq, int stopSeq, vector<int>& otuIndex){
1687 vector<double> newTau(numOTUs,0);
1688 vector<double> norms(numSeqs, 0);
1695 for(int i=startSeq;i<stopSeq;i++){
1696 double offset = 1e8;
1697 int indexOffset = i * numOTUs;
1699 for(int j=0;j<numOTUs;j++){
1701 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1702 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1704 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1705 offset = dist[indexOffset + j];
1709 for(int j=0;j<numOTUs;j++){
1710 if(weight[j] > MIN_WEIGHT){
1711 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1712 norms[i] += newTau[j];
1719 for(int j=0;j<numOTUs;j++){
1721 newTau[j] /= norms[i];
1723 if(newTau[j] > MIN_TAU){
1724 otuIndex.push_back(j);
1725 seqIndex.push_back(i);
1726 singleTau.push_back(newTau[j]);
1732 catch(exception& e) {
1733 m->errorOut(e, "ShhherCommand", "calcNewDistancesChildMPI");
1738 /**************************************************************************************************/
1740 void ShhherCommand::calcNewDistancesChild(int startSeq, int stopSeq, vector<int>& child_otuIndex, vector<int>& child_seqIndex, vector<double>& child_singleTau){
1743 vector<double> newTau(numOTUs,0);
1744 vector<double> norms(numSeqs, 0);
1745 child_otuIndex.resize(0);
1746 child_seqIndex.resize(0);
1747 child_singleTau.resize(0);
1749 for(int i=startSeq;i<stopSeq;i++){
1750 double offset = 1e8;
1751 int indexOffset = i * numOTUs;
1754 for(int j=0;j<numOTUs;j++){
1755 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1756 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1759 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1760 offset = dist[indexOffset + j];
1764 for(int j=0;j<numOTUs;j++){
1765 if(weight[j] > MIN_WEIGHT){
1766 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1767 norms[i] += newTau[j];
1774 for(int j=0;j<numOTUs;j++){
1775 newTau[j] /= norms[i];
1777 if(newTau[j] > MIN_TAU){
1778 child_otuIndex.push_back(j);
1779 child_seqIndex.push_back(i);
1780 child_singleTau.push_back(newTau[j]);
1785 catch(exception& e) {
1786 m->errorOut(e, "ShhherCommand", "calcNewDistancesChild");
1791 /**************************************************************************************************/
1793 void ShhherCommand::calcNewDistancesParent(int startSeq, int stopSeq){
1798 vector<double> newTau(numOTUs,0);
1799 vector<double> norms(numSeqs, 0);
1800 nSeqsPerOTU.assign(numOTUs, 0);
1802 for(int i=startSeq;i<stopSeq;i++){
1803 int indexOffset = i * numOTUs;
1805 double offset = 1e8;
1807 for(int j=0;j<numOTUs;j++){
1808 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1809 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1812 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1813 offset = dist[indexOffset + j];
1817 for(int j=0;j<numOTUs;j++){
1818 if(weight[j] > MIN_WEIGHT){
1819 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1820 norms[i] += newTau[j];
1827 for(int j=0;j<numOTUs;j++){
1828 newTau[j] /= norms[i];
1831 for(int j=0;j<numOTUs;j++){
1832 if(newTau[j] > MIN_TAU){
1834 int oldTotal = total;
1838 singleTau.resize(total, 0);
1839 seqNumber.resize(total, 0);
1840 seqIndex.resize(total, 0);
1842 singleTau[oldTotal] = newTau[j];
1844 aaP[j][nSeqsPerOTU[j]] = oldTotal;
1845 aaI[j][nSeqsPerOTU[j]] = i;
1851 catch(exception& e) {
1852 m->errorOut(e, "ShhherCommand", "calcNewDistancesParent");
1857 /**************************************************************************************************/
1859 void ShhherCommand::setOTUs(){
1862 vector<double> bigTauMatrix(numOTUs * numSeqs, 0.0000);
1864 for(int i=0;i<numOTUs;i++){
1865 for(int j=0;j<nSeqsPerOTU[i];j++){
1866 int index = cumNumSeqs[i] + j;
1867 double tauValue = singleTau[seqNumber[index]];
1868 int sIndex = seqIndex[index];
1869 bigTauMatrix[sIndex * numOTUs + i] = tauValue;
1873 for(int i=0;i<numSeqs;i++){
1874 double maxTau = -1.0000;
1876 for(int j=0;j<numOTUs;j++){
1877 if(bigTauMatrix[i * numOTUs + j] > maxTau){
1878 maxTau = bigTauMatrix[i * numOTUs + j];
1883 otuData[i] = maxOTU;
1886 nSeqsPerOTU.assign(numOTUs, 0);
1888 for(int i=0;i<numSeqs;i++){
1889 int index = otuData[i];
1891 singleTau[i] = 1.0000;
1894 aaP[index][nSeqsPerOTU[index]] = i;
1895 aaI[index][nSeqsPerOTU[index]] = i;
1897 nSeqsPerOTU[index]++;
1901 catch(exception& e) {
1902 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1907 /**************************************************************************************************/
1909 void ShhherCommand::writeQualities(vector<int> otuCounts){
1912 string qualityFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.qual";
1914 ofstream qualityFile;
1915 m->openOutputFile(qualityFileName, qualityFile);
1917 qualityFile.setf(ios::fixed, ios::floatfield);
1918 qualityFile.setf(ios::showpoint);
1919 qualityFile << setprecision(6);
1921 vector<vector<int> > qualities(numOTUs);
1922 vector<double> pr(HOMOPS, 0);
1925 for(int i=0;i<numOTUs;i++){
1929 if(nSeqsPerOTU[i] > 0){
1930 qualities[i].assign(1024, -1);
1932 while(index < numFlowCells){
1933 double maxPrValue = 1e8;
1934 short maxPrIndex = -1;
1935 double count = 0.0000;
1937 pr.assign(HOMOPS, 0);
1939 for(int j=0;j<nSeqsPerOTU[i];j++){
1940 int lIndex = cumNumSeqs[i] + j;
1941 double tauValue = singleTau[seqNumber[lIndex]];
1942 int sequenceIndex = aaI[i][j];
1943 short intensity = flowDataIntI[sequenceIndex * numFlowCells + index];
1947 for(int s=0;s<HOMOPS;s++){
1948 pr[s] += tauValue * singleLookUp[s * NUMBINS + intensity];
1952 maxPrIndex = uniqueFlowgrams[centroids[i] * numFlowCells + index];
1953 maxPrValue = pr[maxPrIndex];
1955 if(count > MIN_COUNT){
1957 double norm = 0.0000;
1959 for(int s=0;s<HOMOPS;s++){
1960 norm += exp(-(pr[s] - maxPrValue));
1963 for(int s=1;s<=maxPrIndex;s++){
1965 double temp = 0.0000;
1967 U += exp(-(pr[s-1]-maxPrValue))/norm;
1975 temp = floor(-10 * temp);
1976 value = (int)floor(temp);
1977 if(value > 100){ value = 100; }
1979 qualities[i][base] = (int)value;
1989 if(otuCounts[i] > 0){
1990 qualityFile << '>' << seqNameVector[mapUniqueToSeq[i]] << endl;
1992 int j=4; //need to get past the first four bases
1993 while(qualities[i][j] != -1){
1994 qualityFile << qualities[i][j] << ' ';
1997 qualityFile << endl;
2000 qualityFile.close();
2003 catch(exception& e) {
2004 m->errorOut(e, "ShhherCommand", "writeQualities");
2009 /**************************************************************************************************/
2011 void ShhherCommand::writeSequences(vector<int> otuCounts){
2013 string bases = "TACG";
2015 string fastaFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.fasta";
2017 m->openOutputFile(fastaFileName, fastaFile);
2019 vector<string> names(numOTUs, "");
2021 for(int i=0;i<numOTUs;i++){
2022 int index = centroids[i];
2024 if(otuCounts[i] > 0){
2025 fastaFile << '>' << seqNameVector[aaI[i][0]] << endl;
2027 for(int j=8;j<numFlowCells;j++){
2029 char base = bases[j % 4];
2030 for(int k=0;k<uniqueFlowgrams[index * numFlowCells + j];k++){
2039 if(compositeFASTAFileName != ""){
2040 m->appendFiles(fastaFileName, compositeFASTAFileName);
2043 catch(exception& e) {
2044 m->errorOut(e, "ShhherCommand", "writeSequences");
2049 /**************************************************************************************************/
2051 void ShhherCommand::writeNames(vector<int> otuCounts){
2053 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.final.names";
2055 m->openOutputFile(nameFileName, nameFile);
2057 for(int i=0;i<numOTUs;i++){
2058 if(otuCounts[i] > 0){
2059 nameFile << seqNameVector[aaI[i][0]] << '\t' << seqNameVector[aaI[i][0]];
2061 for(int j=1;j<nSeqsPerOTU[i];j++){
2062 nameFile << ',' << seqNameVector[aaI[i][j]];
2070 catch(exception& e) {
2071 m->errorOut(e, "ShhherCommand", "writeNames");
2076 /**************************************************************************************************/
2078 void ShhherCommand::writeGroups(){
2080 string fileRoot = flowFileName.substr(0,flowFileName.find_last_of('.'));
2081 string groupFileName = fileRoot + ".pn.groups";
2083 m->openOutputFile(groupFileName, groupFile);
2085 for(int i=0;i<numSeqs;i++){
2086 groupFile << seqNameVector[i] << '\t' << fileRoot << endl;
2090 catch(exception& e) {
2091 m->errorOut(e, "ShhherCommand", "writeGroups");
2096 /**************************************************************************************************/
2098 void ShhherCommand::writeClusters(vector<int> otuCounts){
2100 string otuCountsFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.counts";
2101 ofstream otuCountsFile;
2102 m->openOutputFile(otuCountsFileName, otuCountsFile);
2104 string bases = "TACG";
2106 for(int i=0;i<numOTUs;i++){
2107 //output the translated version of the centroid sequence for the otu
2108 if(otuCounts[i] > 0){
2109 int index = centroids[i];
2111 otuCountsFile << "ideal\t";
2112 for(int j=8;j<numFlowCells;j++){
2113 char base = bases[j % 4];
2114 for(int s=0;s<uniqueFlowgrams[index * numFlowCells + j];s++){
2115 otuCountsFile << base;
2118 otuCountsFile << endl;
2120 for(int j=0;j<nSeqsPerOTU[i];j++){
2121 int sequence = aaI[i][j];
2122 otuCountsFile << seqNameVector[sequence] << '\t';
2124 for(int k=8;k<lengths[sequence];k++){
2125 char base = bases[k % 4];
2126 int freq = int(0.01 * (double)flowDataIntI[sequence * numFlowCells + k] + 0.5);
2128 for(int s=0;s<freq;s++){
2129 otuCountsFile << base;
2132 otuCountsFile << endl;
2134 otuCountsFile << endl;
2137 otuCountsFile.close();
2139 catch(exception& e) {
2140 m->errorOut(e, "ShhherCommand", "writeClusters");
2145 //**********************************************************************************************************************