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(){
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);
107 //allow user to run help
108 if(option == "help") { help(); abort = 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 the user changes the output directory command factory will send this info to us in the output parameter
177 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
179 outputDir += m->hasPath(flowFileName); //if user entered a file with a path then preserve it
183 //check for optional parameter and set defaults
184 // ...at some point should added some additional type checking...
186 temp = validParameter.validFile(parameters, "lookup", true);
187 if (temp == "not found") { lookupFileName = "LookUp_Titanium.pat"; }
188 else if(temp == "not open") { abort = true; }
189 else { lookupFileName = temp; }
191 temp = validParameter.validFile(parameters, "processors", false);if (temp == "not found"){ temp = "1"; }
192 convert(temp, processors);
194 temp = validParameter.validFile(parameters, "cutoff", false); if (temp == "not found"){ temp = "0.01"; }
195 convert(temp, cutoff);
197 temp = validParameter.validFile(parameters, "mindelta", false); if (temp == "not found"){ temp = "0.000001"; }
198 convert(temp, minDelta);
200 temp = validParameter.validFile(parameters, "maxiter", false); if (temp == "not found"){ temp = "1000"; }
201 convert(temp, maxIters);
203 temp = validParameter.validFile(parameters, "sigma", false);if (temp == "not found") { temp = "60"; }
204 convert(temp, sigma);
206 globaldata = GlobalData::getInstance();
214 catch(exception& e) {
215 m->errorOut(e, "ShhherCommand", "ShhherCommand");
220 //**********************************************************************************************************************
222 ShhherCommand::~ShhherCommand(){}
224 //**********************************************************************************************************************
226 void ShhherCommand::help(){
228 m->mothurOut("The shhher command reads a file containing flowgrams and creates a file of corrected sequences.\n");
230 catch(exception& e) {
231 m->errorOut(e, "ShhherCommand", "help");
236 //**********************************************************************************************************************
238 int ShhherCommand::execute(){
243 double begClock = clock();
244 unsigned long int begTime = time(NULL);
246 cout.setf(ios::fixed, ios::floatfield);
247 cout.setf(ios::showpoint);
248 cout << setprecision(2);
253 for(int i=1;i<ncpus;i++){
254 MPI_Send(&abort, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
256 if(abort == 1){ return 0; }
260 cout << "\nGetting preliminary data..." << endl;
264 vector<string> flowFileVector;
265 if(flowFilesFileName != "not found"){
268 ifstream flowFilesFile;
269 m->openInputFile(flowFilesFileName, flowFilesFile);
270 while(flowFilesFile){
271 flowFilesFile >> fName;
272 flowFileVector.push_back(fName);
273 m->gobble(flowFilesFile);
277 flowFileVector.push_back(flowFileName);
279 int numFiles = flowFileVector.size();
281 for(int i=1;i<ncpus;i++){
282 MPI_Send(&numFiles, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
285 for(int i=0;i<numFiles;i++){
286 flowFileName = flowFileVector[i];
288 cout << "\n>>>>>\tProcessing " << flowFileName << " (file " << i+1 << " of " << numFiles << ")\t<<<<<" << endl;
289 cout << "Reading flowgrams..." << endl;
291 cout << "Identifying unique flowgrams..." << endl;
294 cout << "Calculating distances between flowgrams..." << endl;
296 strcpy(fileName, flowFileName.c_str());
298 for(int i=1;i<ncpus;i++){
299 MPI_Send(&fileName[0], 1024, MPI_CHAR, i, tag, MPI_COMM_WORLD);
301 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
302 MPI_Send(&numUniques, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
303 MPI_Send(&numFlowCells, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
304 MPI_Send(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, i, tag, MPI_COMM_WORLD);
305 MPI_Send(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
306 MPI_Send(&mapUniqueToSeq[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
307 MPI_Send(&mapSeqToUnique[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
308 MPI_Send(&lengths[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
309 MPI_Send(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
310 MPI_Send(&cutoff, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
313 string distFileName = flowDistMPI(0, int(sqrt(1.0/float(ncpus)) * numUniques));
316 for(int i=1;i<ncpus;i++){
317 MPI_Recv(&done, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
319 m->appendFiles((distFileName + ".temp." + toString(i)), distFileName);
320 remove((distFileName + ".temp." + toString(i)).c_str());
323 string namesFileName = createNamesFile();
325 cout << "\nClustering flowgrams..." << endl;
326 string listFileName = cluster(distFileName, namesFileName);
327 // string listFileName = "PriestPot_C7.pn.list";
328 // string listFileName = "test.mock_rep3.v69.pn.list";
330 getOTUData(listFileName);
333 for(int i=1;i<ncpus;i++){
334 MPI_Send(&numOTUs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
335 MPI_Send(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
336 MPI_Send(&uniqueFlowgrams[0], numFlowCells * numUniques, MPI_SHORT, i, tag, MPI_COMM_WORLD);
337 MPI_Send(&sigma, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
344 int numOTUsOnCPU = numOTUs / ncpus;
345 int numSeqsOnCPU = numSeqs / ncpus;
347 cout << "\nDenoising flowgrams..." << endl;
348 cout << "iter\tmaxDelta\tnLL\t\tcycletime" << endl;
350 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
352 double cycClock = clock();
353 unsigned long int cycTime = time(NULL);
356 int total = singleTau.size();
357 for(int i=1;i<ncpus;i++){
358 MPI_Send(&total, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
359 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
360 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
362 MPI_Send(&singleTau[0], total, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
363 MPI_Send(&seqNumber[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
364 MPI_Send(&seqIndex[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
365 MPI_Send(&nSeqsPerOTU[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
366 MPI_Send(&cumNumSeqs[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
369 calcCentroidsDriver(0, numOTUsOnCPU);
371 for(int i=1;i<ncpus;i++){
372 int otuStart = i * numOTUs / ncpus;
373 int otuStop = (i + 1) * numOTUs / ncpus;
375 vector<int> tempCentroids(numOTUs, 0);
376 vector<short> tempChange(numOTUs, 0);
378 MPI_Recv(&tempCentroids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
379 MPI_Recv(&tempChange[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD, &status);
381 for(int j=otuStart;j<otuStop;j++){
382 centroids[j] = tempCentroids[j];
383 change[j] = tempChange[j];
387 maxDelta = getNewWeights();
388 double nLL = getLikelihood();
391 for(int i=1;i<ncpus;i++){
392 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
393 MPI_Send(&weight[0], numOTUs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
394 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
397 calcNewDistancesParent(0, numSeqsOnCPU);
399 total = singleTau.size();
401 for(int i=1;i<ncpus;i++){
403 int seqStart = i * numSeqs / ncpus;
404 int seqStop = (i + 1) * numSeqs / ncpus;
406 MPI_Recv(&childTotal, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
408 vector<int> childSeqIndex(childTotal, 0);
409 vector<double> childSingleTau(childTotal, 0);
410 vector<double> childDist(numSeqs * numOTUs, 0);
411 vector<int> otuIndex(childTotal, 0);
413 MPI_Recv(&childSeqIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
414 MPI_Recv(&childSingleTau[0], childTotal, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
415 MPI_Recv(&childDist[0], numOTUs * numSeqs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
416 MPI_Recv(&otuIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
418 int oldTotal = total;
420 singleTau.resize(total, 0);
421 seqIndex.resize(total, 0);
422 seqNumber.resize(total, 0);
426 for(int j=oldTotal;j<total;j++){
427 int otuI = otuIndex[childIndex];
428 int seqI = childSeqIndex[childIndex];
430 singleTau[j] = childSingleTau[childIndex];
432 aaP[otuI][nSeqsPerOTU[otuI]] = j;
433 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
438 int index = seqStart * numOTUs;
439 for(int j=seqStart;j<seqStop;j++){
440 for(int k=0;k<numOTUs;k++){
441 dist[index] = childDist[index];
449 cout << iter << '\t' << maxDelta << '\t' << setprecision(2) << nLL << '\t' << time(NULL) - cycTime << '\t' << setprecision(6) << (clock() - cycClock)/(double)CLOCKS_PER_SEC << endl;
451 if((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
453 for(int i=1;i<ncpus;i++){
454 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
459 for(int i=1;i<ncpus;i++){
460 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD); //send kill command
466 cout << "\nFinalizing..." << endl;
469 vector<int> otuCounts(numOTUs, 0);
470 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
471 calcCentroidsDriver(0, numOTUs);
472 writeQualities(otuCounts);
473 writeSequences(otuCounts);
474 writeNames(otuCounts);
475 writeClusters(otuCounts);
478 remove(distFileName.c_str());
479 remove(namesFileName.c_str());
480 remove(listFileName.c_str());
482 cout << "Total time to process " << flowFileName << ":\t" << time(NULL) - begTime << '\t' << setprecision(6) << (clock() - begClock)/(double)CLOCKS_PER_SEC << endl;
490 MPI_Recv(&abort, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
491 if(abort){ return 0; }
494 MPI_Recv(&numFiles, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
496 for(int i=0;i<numFiles;i++){
497 //Now into the pyrodist part
499 MPI_Recv(&fileName, 1024, MPI_CHAR, 0, tag, MPI_COMM_WORLD, &status);
500 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
501 MPI_Recv(&numUniques, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
502 MPI_Recv(&numFlowCells, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
504 flowDataIntI.resize(numSeqs * numFlowCells);
505 flowDataPrI.resize(numSeqs * numFlowCells);
506 mapUniqueToSeq.resize(numSeqs);
507 mapSeqToUnique.resize(numSeqs);
508 lengths.resize(numSeqs);
509 jointLookUp.resize(NUMBINS * NUMBINS);
511 MPI_Recv(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
512 MPI_Recv(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
513 MPI_Recv(&mapUniqueToSeq[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
514 MPI_Recv(&mapSeqToUnique[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
515 MPI_Recv(&lengths[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
516 MPI_Recv(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
517 MPI_Recv(&cutoff, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
519 flowFileName = string(fileName);
520 int flowDistStart = int(sqrt(float(pid)/float(ncpus)) * numUniques);
521 int flowDistEnd = int(sqrt(float(pid+1)/float(ncpus)) * numUniques);
523 string distanceStringChild = flowDistMPI(flowDistStart, flowDistEnd);
526 MPI_Send(&done, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
528 //Now into the pyronoise part
529 MPI_Recv(&numOTUs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
531 singleLookUp.resize(HOMOPS * NUMBINS);
532 uniqueFlowgrams.resize(numUniques * numFlowCells);
533 weight.resize(numOTUs);
534 centroids.resize(numOTUs);
535 change.resize(numOTUs);
536 dist.assign(numOTUs * numSeqs, 0);
537 nSeqsPerOTU.resize(numOTUs);
538 cumNumSeqs.resize(numOTUs);
540 MPI_Recv(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
541 MPI_Recv(&uniqueFlowgrams[0], uniqueFlowgrams.size(), MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
542 MPI_Recv(&sigma, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
544 int startOTU = pid * numOTUs / ncpus;
545 int endOTU = (pid + 1) * numOTUs / ncpus;
547 int startSeq = pid * numSeqs / ncpus;
548 int endSeq = (pid + 1) * numSeqs /ncpus;
554 MPI_Recv(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
555 singleTau.assign(total, 0.0000);
556 seqNumber.assign(total, 0);
557 seqIndex.assign(total, 0);
559 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
560 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
561 MPI_Recv(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
562 MPI_Recv(&seqNumber[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
563 MPI_Recv(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
564 MPI_Recv(&nSeqsPerOTU[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
565 MPI_Recv(&cumNumSeqs[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
567 calcCentroidsDriver(startOTU, endOTU);
569 MPI_Send(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD);
570 MPI_Send(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD);
573 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
574 MPI_Recv(&weight[0], numOTUs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
575 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
577 vector<int> otuIndex(total, 0);
578 calcNewDistancesChildMPI(startSeq, endSeq, otuIndex);
579 total = otuIndex.size();
581 MPI_Send(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
582 MPI_Send(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
583 MPI_Send(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
584 MPI_Send(&dist[0], numOTUs * numSeqs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
585 MPI_Send(&otuIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
587 MPI_Recv(&live, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
592 MPI_Barrier(MPI_COMM_WORLD);
596 catch(exception& e) {
597 m->errorOut(e, "ShhherCommand", "execute");
602 /**************************************************************************************************/
604 string ShhherCommand::flowDistMPI(int startSeq, int stopSeq){
606 ostringstream outStream;
607 outStream.setf(ios::fixed, ios::floatfield);
608 outStream.setf(ios::dec, ios::basefield);
609 outStream.setf(ios::showpoint);
610 outStream.precision(6);
612 int begTime = time(NULL);
613 double begClock = clock();
615 for(int i=startSeq;i<stopSeq;i++){
616 for(int j=0;j<i;j++){
617 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
619 if(flowDistance < 1e-6){
620 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
622 else if(flowDistance <= cutoff){
623 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
627 cout << i << "\t" << (time(NULL) - begTime) << "\t" << (clock()-begClock)/CLOCKS_PER_SEC << endl;
630 cout << stopSeq << "\t" << (time(NULL) - begTime) << "\t" << (clock()-begClock)/CLOCKS_PER_SEC << endl;
632 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
633 if(pid != 0){ fDistFileName += ".temp." + toString(pid); }
635 ofstream distFile(fDistFileName.c_str());
636 distFile << outStream.str();
639 return fDistFileName;
641 catch(exception& e) {
642 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
648 //**********************************************************************************************************************
650 int ShhherCommand::execute(){
652 if (abort == true) { return 0; }
654 cout.setf(ios::fixed, ios::floatfield);
655 cout.setf(ios::showpoint);
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 cout << "\n>>>>>\tProcessing " << flowFileName << " (file " << i+1 << " of " << numFiles << ")\t<<<<<" << endl;
682 cout << "Reading flowgrams..." << endl;
684 cout << "Identifying unique flowgrams..." << endl;
688 cout << "Calculating distances between flowgrams..." << endl;
689 string distFileName = createDistFile(processors);
690 string namesFileName = createNamesFile();
692 cout << "\nClustering flowgrams..." << endl;
693 string listFileName = cluster(distFileName, namesFileName);
694 getOTUData(listFileName);
701 double begClock = clock();
702 unsigned long int begTime = time(NULL);
704 cout << "\nDenoising flowgrams..." << endl;
705 cout << "iter\tmaxDelta\tnLL\t\tcycletime" << endl;
707 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
709 double cycClock = clock();
710 unsigned long int cycTime = time(NULL);
715 maxDelta = getNewWeights();
716 double nLL = getLikelihood();
723 cout << iter << '\t' << maxDelta << '\t' << setprecision(2) << nLL << '\t' << time(NULL) - cycTime << '\t' << setprecision(6) << (clock() - cycClock)/(double)CLOCKS_PER_SEC << endl;
726 cout << "\nFinalizing..." << endl;
730 vector<int> otuCounts(numOTUs, 0);
731 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
733 calcCentroidsDriver(0, numOTUs);
734 writeQualities(otuCounts);
735 writeSequences(otuCounts);
736 writeNames(otuCounts);
737 writeClusters(otuCounts);
740 remove(distFileName.c_str());
741 remove(namesFileName.c_str());
742 remove(listFileName.c_str());
744 cout << "Total time to process " << flowFileName << ":\t" << time(NULL) - begTime << '\t' << setprecision(6) << (clock() - begClock)/(double)CLOCKS_PER_SEC << endl;
748 catch(exception& e) {
749 m->errorOut(e, "ShhherCommand", "execute");
754 /**************************************************************************************************/
756 void ShhherCommand::getFlowData(){
759 m->openInputFile(flowFileName, flowFile);
763 int currentNumFlowCells;
767 flowFile >> numFlowCells;
768 int index = 0;//pcluster
769 while(!flowFile.eof()){
770 flowFile >> seqName >> currentNumFlowCells;
771 lengths.push_back(currentNumFlowCells);
773 seqNameVector.push_back(seqName);
774 nameMap[seqName] = index++;//pcluster
776 for(int i=0;i<numFlowCells;i++){
777 flowFile >> intensity;
778 if(intensity > 9.99) { intensity = 9.99; }
779 int intI = int(100 * intensity + 0.0001);
780 flowDataIntI.push_back(intI);
786 numSeqs = seqNameVector.size();
788 for(int i=0;i<numSeqs;i++){
789 int iNumFlowCells = i * numFlowCells;
790 for(int j=lengths[i];j<numFlowCells;j++){
791 flowDataIntI[iNumFlowCells + j] = 0;
796 catch(exception& e) {
797 m->errorOut(e, "ShhherCommand", "getFlowData");
802 /**************************************************************************************************/
804 void ShhherCommand::getSingleLookUp(){
806 // these are the -log probabilities that a signal corresponds to a particular homopolymer length
807 singleLookUp.assign(HOMOPS * NUMBINS, 0);
811 m->openInputFile(lookupFileName, lookUpFile);
813 for(int i=0;i<HOMOPS;i++){
815 lookUpFile >> logFracFreq;
817 for(int j=0;j<NUMBINS;j++) {
818 lookUpFile >> singleLookUp[index];
824 catch(exception& e) {
825 m->errorOut(e, "ShhherCommand", "getSingleLookUp");
830 /**************************************************************************************************/
832 void ShhherCommand::getJointLookUp(){
835 // the most likely joint probability (-log) that two intenities have the same polymer length
836 jointLookUp.resize(NUMBINS * NUMBINS, 0);
838 for(int i=0;i<NUMBINS;i++){
839 for(int j=0;j<NUMBINS;j++){
841 double minSum = 10000000000;
843 for(int k=0;k<HOMOPS;k++){
844 double sum = singleLookUp[k * NUMBINS + i] + singleLookUp[k * NUMBINS + j];
846 if(sum < minSum) { minSum = sum; }
848 jointLookUp[i * NUMBINS + j] = minSum;
852 catch(exception& e) {
853 m->errorOut(e, "ShhherCommand", "getJointLookUp");
858 /**************************************************************************************************/
860 float ShhherCommand::getProbIntensity(int intIntensity){
862 double minNegLogProb = 10000000000;
864 for(int i=0;i<HOMOPS;i++){//loop signal strength
865 float negLogProb = singleLookUp[i * NUMBINS + intIntensity];
866 if(negLogProb < minNegLogProb) { minNegLogProb = negLogProb; }
869 return minNegLogProb;
871 catch(exception& e) {
872 m->errorOut(e, "ShhherCommand", "getProbIntensity");
877 /**************************************************************************************************/
879 void ShhherCommand::getUniques(){
884 uniqueFlowgrams.assign(numFlowCells * numSeqs, -1);
885 uniqueCount.assign(numSeqs, 0); // anWeights
886 uniqueLengths.assign(numSeqs, 0);
887 mapSeqToUnique.assign(numSeqs, -1);
888 mapUniqueToSeq.assign(numSeqs, -1);
890 vector<short> uniqueFlowDataIntI(numFlowCells * numSeqs, -1);
892 for(int i=0;i<numSeqs;i++){
895 vector<short> current(numFlowCells);
896 for(int j=0;j<numFlowCells;j++){ current[j] = short(((flowDataIntI[i * numFlowCells + j] + 50.0)/100.0)); }
898 for(int j=0;j<numUniques;j++){
899 int offset = j * numFlowCells;
902 for(int k=0;k<numFlowCells;k++){
903 if(current[k] != uniqueFlowgrams[offset + k]){
910 mapSeqToUnique[i] = j;
918 if(index == numUniques){
919 uniqueLengths[numUniques] = lengths[i];
920 uniqueCount[numUniques] = 1;
921 mapSeqToUnique[i] = numUniques;//anMap
922 mapUniqueToSeq[numUniques] = i;//anF
924 for(int k=0;k<numFlowCells;k++){
925 uniqueFlowgrams[numUniques * numFlowCells + k] = current[k];
926 uniqueFlowDataIntI[numUniques * numFlowCells + k] = flowDataIntI[i * numFlowCells + k];
932 uniqueFlowDataIntI.resize(numFlowCells * numUniques);
933 uniqueLengths.resize(numUniques);
935 flowDataPrI.assign(numSeqs * numFlowCells, 0);
936 for(int i=0;i<flowDataPrI.size();i++) { flowDataPrI[i] = getProbIntensity(flowDataIntI[i]); }
938 catch(exception& e) {
939 m->errorOut(e, "ShhherCommand", "getUniques");
944 /**************************************************************************************************/
946 float ShhherCommand::calcPairwiseDist(int seqA, int seqB){
948 int minLength = lengths[mapSeqToUnique[seqA]];
949 if(lengths[seqB] < minLength){ minLength = lengths[mapSeqToUnique[seqB]]; }
951 int ANumFlowCells = seqA * numFlowCells;
952 int BNumFlowCells = seqB * numFlowCells;
956 for(int i=0;i<minLength;i++){
957 int flowAIntI = flowDataIntI[ANumFlowCells + i];
958 float flowAPrI = flowDataPrI[ANumFlowCells + i];
960 int flowBIntI = flowDataIntI[BNumFlowCells + i];
961 float flowBPrI = flowDataPrI[BNumFlowCells + i];
962 dist += jointLookUp[flowAIntI * NUMBINS + flowBIntI] - flowAPrI - flowBPrI;
965 dist /= (float) minLength;
968 catch(exception& e) {
969 m->errorOut(e, "ShhherCommand", "calcPairwiseDist");
974 /**************************************************************************************************/
976 void ShhherCommand::flowDistParentFork(string distFileName, int startSeq, int stopSeq){
979 ostringstream outStream;
980 outStream.setf(ios::fixed, ios::floatfield);
981 outStream.setf(ios::dec, ios::basefield);
982 outStream.setf(ios::showpoint);
983 outStream.precision(6);
985 int begTime = time(NULL);
986 double begClock = clock();
988 for(int i=startSeq;i<stopSeq;i++){
989 for(int j=0;j<i;j++){
990 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
992 if(flowDistance < 1e-6){
993 outStream << seqNameVector[mapUniqueToSeq[i]] << '\t' << seqNameVector[mapUniqueToSeq[j]] << '\t' << 0.000000 << endl;
995 else if(flowDistance <= cutoff){
996 outStream << seqNameVector[mapUniqueToSeq[i]] << '\t' << seqNameVector[mapUniqueToSeq[j]] << '\t' << flowDistance << endl;
1000 m->mothurOut(toString(i) + "\t" + toString(time(NULL) - begTime));
1001 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1002 m->mothurOutEndLine();
1005 m->mothurOut(toString(stopSeq-1) + "\t" + toString(time(NULL) - begTime));
1006 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1007 m->mothurOutEndLine();
1009 ofstream distFile(distFileName.c_str());
1010 distFile << outStream.str();
1013 catch(exception& e) {
1014 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
1019 /**************************************************************************************************/
1021 string ShhherCommand::createDistFile(int processors){
1023 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.dist";
1025 unsigned long int begTime = time(NULL);
1026 double begClock = clock();
1031 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1032 if(processors == 1) { flowDistParentFork(fDistFileName, 0, numUniques); }
1033 else{ //you have multiple processors
1035 if (numSeqs < processors){ processors = 1; }
1037 vector<int> start(processors, 0);
1038 vector<int> end(processors, 0);
1040 for (int i = 0; i < processors; i++) {
1041 start[i] = int(sqrt(float(i)/float(processors)) * numUniques);
1042 end[i] = int(sqrt(float(i+1)/float(processors)) * numUniques);
1046 vector<int> processIDs;
1048 //loop through and create all the processes you want
1049 while (process != processors) {
1053 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1055 }else if (pid == 0){
1056 flowDistParentFork(fDistFileName + toString(getpid()) + ".temp", start[process], end[process]);
1059 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1061 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1066 //parent does its part
1067 flowDistParentFork(fDistFileName, start[0], end[0]);
1069 //force parent to wait until all the processes are done
1070 for (int i=0;i<processIDs.size();i++) {
1071 int temp = processIDs[i];
1075 //append and remove temp files
1076 for (int i=0;i<processIDs.size();i++) {
1077 m->appendFiles((fDistFileName + toString(processIDs[i]) + ".temp"), fDistFileName);
1078 remove((fDistFileName + toString(processIDs[i]) + ".temp").c_str());
1084 flowDistParentFork(fDistFileName, 0, numUniques);
1087 m->mothurOutEndLine();
1089 cout << "Total time: " << (time(NULL) - begTime) << "\t" << (clock() - begClock)/CLOCKS_PER_SEC << endl;;
1091 return fDistFileName;
1093 catch(exception& e) {
1094 m->errorOut(e, "ShhherCommand", "createDistFile");
1100 /**************************************************************************************************/
1102 string ShhherCommand::createNamesFile(){
1105 vector<string> duplicateNames(numUniques, "");
1106 for(int i=0;i<numSeqs;i++){
1107 duplicateNames[mapSeqToUnique[i]] += seqNameVector[i] + ',';
1110 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.names";
1113 m->openOutputFile(nameFileName, nameFile);
1115 for(int i=0;i<numUniques;i++){
1116 // nameFile << seqNameVector[mapUniqueToSeq[i]] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1117 nameFile << mapUniqueToSeq[i] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1121 return nameFileName;
1123 catch(exception& e) {
1124 m->errorOut(e, "ShhherCommand", "createNamesFile");
1129 //**********************************************************************************************************************
1131 string ShhherCommand::cluster(string distFileName, string namesFileName){
1134 SparseMatrix* matrix;
1136 RAbundVector* rabund;
1138 globaldata->setNameFile(namesFileName);
1139 globaldata->setColumnFile(distFileName);
1140 globaldata->setFormat("column");
1142 ReadMatrix* read = new ReadColumnMatrix(distFileName);
1143 read->setCutoff(cutoff);
1145 NameAssignment* clusterNameMap = new NameAssignment(namesFileName);
1146 clusterNameMap->readMap();
1147 read->read(clusterNameMap);
1149 list = read->getListVector();
1150 matrix = read->getMatrix();
1153 delete clusterNameMap;
1155 rabund = new RAbundVector(list->getRAbundVector());
1157 Cluster* cluster = new CompleteLinkage(rabund, list, matrix, cutoff, "furthest");
1158 string tag = cluster->getTag();
1160 double clusterCutoff = cutoff;
1161 while (matrix->getSmallDist() <= clusterCutoff && matrix->getNNodes() > 0){
1162 cluster->update(clusterCutoff);
1165 list->setLabel(toString(cutoff));
1167 string listFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.list";
1169 m->openOutputFile(listFileName, listFile);
1170 list->print(listFile);
1173 delete matrix; delete cluster; delete rabund; delete list;
1175 return listFileName;
1177 catch(exception& e) {
1178 m->errorOut(e, "ShhherCommand", "cluster");
1183 /**************************************************************************************************/
1185 void ShhherCommand::getOTUData(string listFileName){
1189 m->openInputFile(listFileName, listFile);
1192 listFile >> label >> numOTUs;
1194 otuData.assign(numSeqs, 0);
1195 cumNumSeqs.assign(numOTUs, 0);
1196 nSeqsPerOTU.assign(numOTUs, 0);
1197 aaP.resize(numOTUs);
1199 string singleOTU = "";
1201 for(int i=0;i<numOTUs;i++){
1203 listFile >> singleOTU;
1205 istringstream otuString(singleOTU);
1209 string seqName = "";
1211 for(int j=0;j<singleOTU.length();j++){
1212 char letter = otuString.get();
1218 map<string,int>::iterator nmIt = nameMap.find(seqName);
1219 int index = nmIt->second;
1221 nameMap.erase(nmIt);
1225 aaP[i].push_back(index);
1230 map<string,int>::iterator nmIt = nameMap.find(seqName);
1232 int index = nmIt->second;
1233 nameMap.erase(nmIt);
1237 aaP[i].push_back(index);
1242 sort(aaP[i].begin(), aaP[i].end());
1243 for(int j=0;j<nSeqsPerOTU[i];j++){
1244 seqNumber.push_back(aaP[i][j]);
1246 for(int j=nSeqsPerOTU[i];j<numSeqs;j++){
1247 aaP[i].push_back(0);
1251 for(int i=1;i<numOTUs;i++){
1252 cumNumSeqs[i] = cumNumSeqs[i-1] + nSeqsPerOTU[i-1];
1255 seqIndex = seqNumber;
1259 catch(exception& e) {
1260 m->errorOut(e, "ShhherCommand", "getOTUData");
1265 /**************************************************************************************************/
1267 void ShhherCommand::initPyroCluster(){
1269 dist.assign(numSeqs * numOTUs, 0);
1270 change.assign(numOTUs, 1);
1271 centroids.assign(numOTUs, -1);
1272 weight.assign(numOTUs, 0);
1273 singleTau.assign(numSeqs, 1.0);
1275 nSeqsBreaks.assign(processors+1, 0);
1276 nOTUsBreaks.assign(processors+1, 0);
1279 for(int i=0;i<processors;i++){
1280 nSeqsBreaks[i+1] = nSeqsBreaks[i] + (int)((double) numSeqs / (double) processors);
1281 nOTUsBreaks[i+1] = nOTUsBreaks[i] + (int)((double) numOTUs / (double) processors);
1283 nSeqsBreaks[processors] = numSeqs;
1284 nOTUsBreaks[processors] = numOTUs;
1286 catch(exception& e) {
1287 m->errorOut(e, "ShhherCommand", "initPyroCluster");
1292 /**************************************************************************************************/
1294 void ShhherCommand::fill(){
1297 for(int i=0;i<numOTUs;i++){
1298 cumNumSeqs[i] = index;
1299 for(int j=0;j<nSeqsPerOTU[i];j++){
1300 seqNumber[index] = aaP[i][j];
1301 seqIndex[index] = aaI[i][j];
1307 catch(exception& e) {
1308 m->errorOut(e, "ShhherCommand", "fill");
1313 /**************************************************************************************************/
1315 void ShhherCommand::calcCentroids(){
1318 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1320 if(processors == 1) {
1321 calcCentroidsDriver(0, numOTUs);
1323 else{ //you have multiple processors
1324 if (numOTUs < processors){ processors = 1; }
1327 vector<int> processIDs;
1329 //loop through and create all the processes you want
1330 while (process != processors) {
1334 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1336 }else if (pid == 0){
1337 calcCentroidsDriver(nOTUsBreaks[process], nOTUsBreaks[process+1]);
1340 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1342 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1347 //parent does its part
1348 calcCentroidsDriver(nOTUsBreaks[0], nOTUsBreaks[1]);
1350 //force parent to wait until all the processes are done
1351 for (int i=0;i<processIDs.size();i++) {
1352 int temp = processIDs[i];
1358 calcCentroidsDriver(0, numOTUs);
1361 catch(exception& e) {
1362 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1367 /**************************************************************************************************/
1369 void ShhherCommand::calcCentroidsDriver(int start, int finish){
1371 //this function gets the most likely homopolymer length at a flow position for a group of sequences
1376 for(int i=start;i<finish;i++){
1380 int minFlowGram = 100000000;
1381 double minFlowValue = 1e8;
1382 change[i] = 0; //FALSE
1384 for(int j=0;j<nSeqsPerOTU[i];j++){
1385 count += singleTau[seqNumber[cumNumSeqs[i] + j]];
1388 if(nSeqsPerOTU[i] > 0 && count > MIN_COUNT){
1389 vector<double> adF(nSeqsPerOTU[i]);
1390 vector<int> anL(nSeqsPerOTU[i]);
1392 for(int j=0;j<nSeqsPerOTU[i];j++){
1393 int index = cumNumSeqs[i] + j;
1394 int nI = seqIndex[index];
1395 int nIU = mapSeqToUnique[nI];
1398 for(k=0;k<position;k++){
1404 anL[position] = nIU;
1405 adF[position] = 0.0000;
1410 for(int j=0;j<nSeqsPerOTU[i];j++){
1411 int index = cumNumSeqs[i] + j;
1412 int nI = seqIndex[index];
1414 double tauValue = singleTau[seqNumber[index]];
1416 for(int k=0;k<position;k++){
1417 double dist = getDistToCentroid(anL[k], nI, lengths[nI]);
1418 adF[k] += dist * tauValue;
1422 for(int j=0;j<position;j++){
1423 if(adF[j] < minFlowValue){
1425 minFlowValue = adF[j];
1429 if(centroids[i] != anL[minFlowGram]){
1431 centroids[i] = anL[minFlowGram];
1434 else if(centroids[i] != -1){
1440 catch(exception& e) {
1441 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1446 /**************************************************************************************************/
1448 double ShhherCommand::getDistToCentroid(int cent, int flow, int length){
1451 int flowAValue = cent * numFlowCells;
1452 int flowBValue = flow * numFlowCells;
1456 for(int i=0;i<length;i++){
1457 dist += singleLookUp[uniqueFlowgrams[flowAValue] * NUMBINS + flowDataIntI[flowBValue]];
1462 return dist / (double)length;
1464 catch(exception& e) {
1465 m->errorOut(e, "ShhherCommand", "getDistToCentroid");
1470 /**************************************************************************************************/
1472 double ShhherCommand::getNewWeights(){
1475 double maxChange = 0;
1477 for(int i=0;i<numOTUs;i++){
1479 double difference = weight[i];
1482 for(int j=0;j<nSeqsPerOTU[i];j++){
1483 int index = cumNumSeqs[i] + j;
1484 double tauValue = singleTau[seqNumber[index]];
1485 weight[i] += tauValue;
1488 difference = fabs(weight[i] - difference);
1489 if(difference > maxChange){ maxChange = difference; }
1493 catch(exception& e) {
1494 m->errorOut(e, "ShhherCommand", "getNewWeights");
1499 /**************************************************************************************************/
1501 double ShhherCommand::getLikelihood(){
1505 vector<long double> P(numSeqs, 0);
1508 for(int i=0;i<numOTUs;i++){
1509 if(weight[i] > MIN_WEIGHT){
1515 for(int i=0;i<numOTUs;i++){
1516 for(int j=0;j<nSeqsPerOTU[i];j++){
1517 int index = cumNumSeqs[i] + j;
1518 int nI = seqIndex[index];
1519 double singleDist = dist[seqNumber[index]];
1521 P[nI] += weight[i] * exp(-singleDist * sigma);
1525 for(int i=0;i<numSeqs;i++){
1526 if(P[i] == 0){ P[i] = DBL_EPSILON; }
1531 nLL = nLL -(double)numSeqs * log(sigma);
1535 catch(exception& e) {
1536 m->errorOut(e, "ShhherCommand", "getNewWeights");
1541 /**************************************************************************************************/
1543 void ShhherCommand::checkCentroids(){
1545 vector<int> unique(numOTUs, 1);
1547 for(int i=0;i<numOTUs;i++){
1548 if(centroids[i] == -1 || weight[i] < MIN_WEIGHT){
1553 for(int i=0;i<numOTUs;i++){
1555 for(int j=i+1;j<numOTUs;j++){
1558 if(centroids[j] == centroids[i]){
1562 weight[i] += weight[j];
1570 catch(exception& e) {
1571 m->errorOut(e, "ShhherCommand", "checkCentroids");
1576 /**************************************************************************************************/
1578 void ShhherCommand::calcNewDistances(){
1581 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1583 if(processors == 1) {
1584 calcNewDistancesParent(0, numSeqs);
1586 else{ //you have multiple processors
1587 if (numSeqs < processors){ processors = 1; }
1589 vector<vector<int> > child_otuIndex(processors);
1590 vector<vector<int> > child_seqIndex(processors);
1591 vector<vector<double> > child_singleTau(processors);
1592 vector<int> totals(processors);
1595 vector<int> processIDs;
1597 //loop through and create all the processes you want
1598 while (process != processors) {
1602 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1604 }else if (pid == 0){
1605 calcNewDistancesChild(nSeqsBreaks[process], nSeqsBreaks[process+1], child_otuIndex[process], child_seqIndex[process], child_singleTau[process]);
1606 totals[process] = child_otuIndex[process].size();
1610 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1612 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1617 //parent does its part
1618 calcNewDistancesParent(nSeqsBreaks[0], nSeqsBreaks[1]);
1619 int total = seqIndex.size();
1621 //force parent to wait until all the processes are done
1622 for (int i=0;i<processIDs.size();i++) {
1623 int temp = processIDs[i];
1627 for(int i=1;i<processors;i++){
1628 int oldTotal = total;
1631 singleTau.resize(total, 0);
1632 seqIndex.resize(total, 0);
1633 seqNumber.resize(total, 0);
1637 for(int j=oldTotal;j<total;j++){
1638 int otuI = child_otuIndex[i][childIndex];
1639 int seqI = child_seqIndex[i][childIndex];
1641 singleTau[j] = child_singleTau[i][childIndex];
1642 aaP[otuI][nSeqsPerOTU[otuI]] = j;
1643 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
1644 nSeqsPerOTU[otuI]++;
1651 calcNewDistancesParent(0, numSeqs);
1654 catch(exception& e) {
1655 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1660 /**************************************************************************************************/
1662 void ShhherCommand::calcNewDistancesChildMPI(int startSeq, int stopSeq, vector<int>& otuIndex){
1665 vector<double> newTau(numOTUs,0);
1666 vector<double> norms(numSeqs, 0);
1669 singleTau.resize(0);
1673 for(int i=startSeq;i<stopSeq;i++){
1674 double offset = 1e8;
1675 int indexOffset = i * numOTUs;
1677 for(int j=0;j<numOTUs;j++){
1679 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1680 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1682 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1683 offset = dist[indexOffset + j];
1687 for(int j=0;j<numOTUs;j++){
1688 if(weight[j] > MIN_WEIGHT){
1689 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1690 norms[i] += newTau[j];
1697 for(int j=0;j<numOTUs;j++){
1699 newTau[j] /= norms[i];
1701 if(newTau[j] > MIN_TAU){
1702 otuIndex.push_back(j);
1703 seqIndex.push_back(i);
1704 singleTau.push_back(newTau[j]);
1710 catch(exception& e) {
1711 m->errorOut(e, "ShhherCommand", "calcNewDistancesChildMPI");
1716 /**************************************************************************************************/
1718 void ShhherCommand::calcNewDistancesChild(int startSeq, int stopSeq, vector<int>& child_otuIndex, vector<int>& child_seqIndex, vector<double>& child_singleTau){
1721 vector<double> newTau(numOTUs,0);
1722 vector<double> norms(numSeqs, 0);
1723 child_otuIndex.resize(0);
1724 child_seqIndex.resize(0);
1725 child_singleTau.resize(0);
1727 for(int i=startSeq;i<stopSeq;i++){
1728 double offset = 1e8;
1729 int indexOffset = i * numOTUs;
1732 for(int j=0;j<numOTUs;j++){
1733 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1734 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1737 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1738 offset = dist[indexOffset + j];
1742 for(int j=0;j<numOTUs;j++){
1743 if(weight[j] > MIN_WEIGHT){
1744 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1745 norms[i] += newTau[j];
1752 for(int j=0;j<numOTUs;j++){
1753 newTau[j] /= norms[i];
1755 if(newTau[j] > MIN_TAU){
1756 child_otuIndex.push_back(j);
1757 child_seqIndex.push_back(i);
1758 child_singleTau.push_back(newTau[j]);
1763 catch(exception& e) {
1764 m->errorOut(e, "ShhherCommand", "calcNewDistancesChild");
1769 /**************************************************************************************************/
1771 void ShhherCommand::calcNewDistancesParent(int startSeq, int stopSeq){
1776 vector<double> newTau(numOTUs,0);
1777 vector<double> norms(numSeqs, 0);
1778 nSeqsPerOTU.assign(numOTUs, 0);
1780 for(int i=startSeq;i<stopSeq;i++){
1781 int indexOffset = i * numOTUs;
1783 double offset = 1e8;
1785 for(int j=0;j<numOTUs;j++){
1786 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1787 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1790 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1791 offset = dist[indexOffset + j];
1795 for(int j=0;j<numOTUs;j++){
1796 if(weight[j] > MIN_WEIGHT){
1797 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1798 norms[i] += newTau[j];
1805 for(int j=0;j<numOTUs;j++){
1806 newTau[j] /= norms[i];
1809 for(int j=0;j<numOTUs;j++){
1810 if(newTau[j] > MIN_TAU){
1812 int oldTotal = total;
1816 singleTau.resize(total, 0);
1817 seqNumber.resize(total, 0);
1818 seqIndex.resize(total, 0);
1820 singleTau[oldTotal] = newTau[j];
1822 aaP[j][nSeqsPerOTU[j]] = oldTotal;
1823 aaI[j][nSeqsPerOTU[j]] = i;
1829 catch(exception& e) {
1830 m->errorOut(e, "ShhherCommand", "calcNewDistancesParent");
1835 /**************************************************************************************************/
1837 void ShhherCommand::setOTUs(){
1840 vector<double> bigTauMatrix(numOTUs * numSeqs, 0.0000);
1842 for(int i=0;i<numOTUs;i++){
1843 for(int j=0;j<nSeqsPerOTU[i];j++){
1844 int index = cumNumSeqs[i] + j;
1845 double tauValue = singleTau[seqNumber[index]];
1846 int sIndex = seqIndex[index];
1847 bigTauMatrix[sIndex * numOTUs + i] = tauValue;
1851 for(int i=0;i<numSeqs;i++){
1852 double maxTau = -1.0000;
1854 for(int j=0;j<numOTUs;j++){
1855 if(bigTauMatrix[i * numOTUs + j] > maxTau){
1856 maxTau = bigTauMatrix[i * numOTUs + j];
1861 otuData[i] = maxOTU;
1864 nSeqsPerOTU.assign(numOTUs, 0);
1866 for(int i=0;i<numSeqs;i++){
1867 int index = otuData[i];
1869 singleTau[i] = 1.0000;
1872 aaP[index][nSeqsPerOTU[index]] = i;
1873 aaI[index][nSeqsPerOTU[index]] = i;
1875 nSeqsPerOTU[index]++;
1879 catch(exception& e) {
1880 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1885 /**************************************************************************************************/
1887 void ShhherCommand::writeQualities(vector<int> otuCounts){
1890 string qualityFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.qual";
1892 ofstream qualityFile;
1893 m->openOutputFile(qualityFileName, qualityFile);
1895 qualityFile.setf(ios::fixed, ios::floatfield);
1896 qualityFile.setf(ios::showpoint);
1897 qualityFile << setprecision(6);
1899 vector<vector<int> > qualities(numOTUs);
1900 vector<double> pr(HOMOPS, 0);
1903 for(int i=0;i<numOTUs;i++){
1907 if(nSeqsPerOTU[i] > 0){
1908 qualities[i].assign(1024, -1);
1910 while(index < numFlowCells){
1911 double maxPrValue = 1e8;
1912 short maxPrIndex = -1;
1913 double count = 0.0000;
1915 pr.assign(HOMOPS, 0);
1917 for(int j=0;j<nSeqsPerOTU[i];j++){
1918 int lIndex = cumNumSeqs[i] + j;
1919 double tauValue = singleTau[seqNumber[lIndex]];
1920 int sequenceIndex = aaI[i][j];
1921 short intensity = flowDataIntI[sequenceIndex * numFlowCells + index];
1925 for(int s=0;s<HOMOPS;s++){
1926 pr[s] += tauValue * singleLookUp[s * NUMBINS + intensity];
1930 maxPrIndex = uniqueFlowgrams[centroids[i] * numFlowCells + index];
1931 maxPrValue = pr[maxPrIndex];
1933 if(count > MIN_COUNT){
1935 double norm = 0.0000;
1937 for(int s=0;s<HOMOPS;s++){
1938 norm += exp(-(pr[s] - maxPrValue));
1941 for(int s=1;s<=maxPrIndex;s++){
1943 double temp = 0.0000;
1945 U += exp(-(pr[s-1]-maxPrValue))/norm;
1953 temp = floor(-10 * temp);
1954 value = (int)floor(temp);
1955 if(value > 100){ value = 100; }
1957 qualities[i][base] = (int)value;
1967 if(otuCounts[i] > 0){
1968 qualityFile << '>' << seqNameVector[mapUniqueToSeq[i]] << endl;
1970 int j=4; //need to get past the first four bases
1971 while(qualities[i][j] != -1){
1972 qualityFile << qualities[i][j] << ' ';
1975 qualityFile << endl;
1978 qualityFile.close();
1981 catch(exception& e) {
1982 m->errorOut(e, "ShhherCommand", "writeQualities");
1987 /**************************************************************************************************/
1989 void ShhherCommand::writeSequences(vector<int> otuCounts){
1992 string bases = "TACG";
1994 string fastaFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.fasta";
1996 m->openOutputFile(fastaFileName, fastaFile);
1998 vector<string> names(numOTUs, "");
2000 for(int i=0;i<numOTUs;i++){
2001 int index = centroids[i];
2003 if(otuCounts[i] > 0){
2004 fastaFile << '>' << seqNameVector[aaI[i][0]] << endl;
2006 for(int j=8;j<numFlowCells;j++){
2008 char base = bases[j % 4];
2009 for(int k=0;k<uniqueFlowgrams[index * numFlowCells + j];k++){
2018 catch(exception& e) {
2019 m->errorOut(e, "ShhherCommand", "writeSequences");
2024 /**************************************************************************************************/
2026 void ShhherCommand::writeNames(vector<int> otuCounts){
2028 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.final.names";
2030 m->openOutputFile(nameFileName, nameFile);
2032 for(int i=0;i<numOTUs;i++){
2033 if(otuCounts[i] > 0){
2034 nameFile << seqNameVector[aaI[i][0]] << '\t' << seqNameVector[aaI[i][0]];
2036 for(int j=1;j<nSeqsPerOTU[i];j++){
2037 nameFile << ',' << seqNameVector[aaI[i][j]];
2045 catch(exception& e) {
2046 m->errorOut(e, "ShhherCommand", "writeNames");
2051 /**************************************************************************************************/
2053 void ShhherCommand::writeGroups(){
2055 string fileRoot = flowFileName.substr(0,flowFileName.find_last_of('.'));
2056 string groupFileName = fileRoot + ".pn.groups";
2058 m->openOutputFile(groupFileName, groupFile);
2060 for(int i=0;i<numSeqs;i++){
2061 groupFile << seqNameVector[i] << '\t' << fileRoot << endl;
2065 catch(exception& e) {
2066 m->errorOut(e, "ShhherCommand", "writeGroups");
2071 /**************************************************************************************************/
2073 void ShhherCommand::writeClusters(vector<int> otuCounts){
2075 string otuCountsFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".pn.counts";
2076 ofstream otuCountsFile;
2077 m->openOutputFile(otuCountsFileName, otuCountsFile);
2079 string bases = "TACG";
2081 for(int i=0;i<numOTUs;i++){
2082 //output the translated version of the centroid sequence for the otu
2083 if(otuCounts[i] > 0){
2084 int index = centroids[i];
2086 otuCountsFile << "ideal\t";
2087 for(int j=8;j<numFlowCells;j++){
2088 char base = bases[j % 4];
2089 for(int s=0;s<uniqueFlowgrams[index * numFlowCells + j];s++){
2090 otuCountsFile << base;
2093 otuCountsFile << endl;
2095 for(int j=0;j<nSeqsPerOTU[i];j++){
2096 int sequence = aaI[i][j];
2097 otuCountsFile << seqNameVector[sequence] << '\t';
2099 for(int k=8;k<lengths[sequence];k++){
2100 char base = bases[k % 4];
2101 int freq = int(0.01 * (double)flowDataIntI[sequence * numFlowCells + k] + 0.5);
2103 for(int s=0;s<freq;s++){
2104 otuCountsFile << base;
2107 otuCountsFile << endl;
2109 otuCountsFile << endl;
2112 otuCountsFile.close();
2114 catch(exception& e) {
2115 m->errorOut(e, "ShhherCommand", "writeClusters");
2120 //**********************************************************************************************************************