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
29 //**********************************************************************************************************************
30 vector<string> ShhherCommand::setParameters(){
32 CommandParameter pflow("flow", "InputTypes", "", "", "none", "fileflow", "none",false,false); parameters.push_back(pflow);
33 CommandParameter pfile("file", "InputTypes", "", "", "none", "fileflow", "none",false,false); parameters.push_back(pfile);
34 CommandParameter plookup("lookup", "InputTypes", "", "", "none", "none", "none",false,true); parameters.push_back(plookup);
35 CommandParameter pcutoff("cutoff", "Number", "", "0.01", "", "", "",false,false); parameters.push_back(pcutoff);
36 CommandParameter pprocessors("processors", "Number", "", "1", "", "", "",false,false); parameters.push_back(pprocessors);
37 CommandParameter pmaxiter("maxiter", "Number", "", "1000", "", "", "",false,false); parameters.push_back(pmaxiter);
38 CommandParameter psigma("sigma", "Number", "", "60", "", "", "",false,false); parameters.push_back(psigma);
39 CommandParameter pmindelta("mindelta", "Number", "", "0.000001", "", "", "",false,false); parameters.push_back(pmindelta);
40 CommandParameter porder("order", "String", "", "", "", "", "",false,false); parameters.push_back(porder);
41 CommandParameter pinputdir("inputdir", "String", "", "", "", "", "",false,false); parameters.push_back(pinputdir);
42 CommandParameter poutputdir("outputdir", "String", "", "", "", "", "",false,false); parameters.push_back(poutputdir);
44 vector<string> myArray;
45 for (int i = 0; i < parameters.size(); i++) { myArray.push_back(parameters[i].name); }
49 m->errorOut(e, "ShhherCommand", "setParameters");
53 //**********************************************************************************************************************
54 string ShhherCommand::getHelpString(){
56 string helpString = "";
57 helpString += "The shhh.seqs command reads a file containing flowgrams and creates a file of corrected sequences.\n";
61 m->errorOut(e, "ShhherCommand", "getHelpString");
65 //**********************************************************************************************************************
67 ShhherCommand::ShhherCommand(){
69 abort = true; calledHelp = true;
72 //initialize outputTypes
73 // vector<string> tempOutNames;
74 // outputTypes["pn.dist"] = tempOutNames;
78 m->errorOut(e, "ShhherCommand", "ShhherCommand");
83 //**********************************************************************************************************************
85 ShhherCommand::ShhherCommand(string option) {
89 MPI_Comm_rank(MPI_COMM_WORLD, &pid); //find out who we are
90 MPI_Comm_size(MPI_COMM_WORLD, &ncpus);
96 abort = false; calledHelp = false;
99 //allow user to run help
100 if(option == "help") { help(); abort = true; calledHelp = true; }
101 else if(option == "citation") { citation(); abort = true; calledHelp = true;}
104 vector<string> myArray = setParameters();
106 OptionParser parser(option);
107 map<string,string> parameters = parser.getParameters();
109 ValidParameters validParameter;
110 map<string,string>::iterator it;
112 //check to make sure all parameters are valid for command
113 for (it = parameters.begin(); it != parameters.end(); it++) {
114 if (validParameter.isValidParameter(it->first, myArray, it->second) != true) { abort = true; }
117 //initialize outputTypes
118 vector<string> tempOutNames;
119 // outputTypes["pn.dist"] = tempOutNames;
120 // outputTypes["fasta"] = tempOutNames;
122 //if the user changes the input directory command factory will send this info to us in the output parameter
123 string inputDir = validParameter.validFile(parameters, "inputdir", false);
124 if (inputDir == "not found"){ inputDir = ""; }
127 it = parameters.find("flow");
128 //user has given a template file
129 if(it != parameters.end()){
130 path = m->hasPath(it->second);
131 //if the user has not given a path then, add inputdir. else leave path alone.
132 if (path == "") { parameters["flow"] = inputDir + it->second; }
135 it = parameters.find("lookup");
136 //user has given a template file
137 if(it != parameters.end()){
138 path = m->hasPath(it->second);
139 //if the user has not given a path then, add inputdir. else leave path alone.
140 if (path == "") { parameters["lookup"] = inputDir + it->second; }
143 it = parameters.find("file");
144 //user has given a template file
145 if(it != parameters.end()){
146 path = m->hasPath(it->second);
147 //if the user has not given a path then, add inputdir. else leave path alone.
148 if (path == "") { parameters["file"] = inputDir + it->second; }
153 //check for required parameters
154 flowFileName = validParameter.validFile(parameters, "flow", true);
155 flowFilesFileName = validParameter.validFile(parameters, "file", true);
156 if (flowFileName == "not found" && flowFilesFileName == "not found") {
157 m->mothurOut("values for either flow or file must be provided for the shhh.seqs command.");
158 m->mothurOutEndLine();
161 else if (flowFileName == "not open" || flowFilesFileName == "not open") { abort = true; }
163 if(flowFileName != "not found"){
164 compositeFASTAFileName = "";
165 compositeNamesFileName = "";
170 compositeFASTAFileName = flowFilesFileName.substr(0, flowFilesFileName.length()-10) + "shhh.fasta";
171 m->openOutputFile(compositeFASTAFileName, temp);
174 compositeNamesFileName = flowFilesFileName.substr(0, flowFilesFileName.length()-10) + "shhh.names";
175 m->openOutputFile(compositeNamesFileName, temp);
179 //if the user changes the output directory command factory will send this info to us in the output parameter
180 outputDir = validParameter.validFile(parameters, "outputdir", false); if (outputDir == "not found"){
182 outputDir += m->hasPath(flowFileName); //if user entered a file with a path then preserve it
186 //check for optional parameter and set defaults
187 // ...at some point should added some additional type checking...
189 temp = validParameter.validFile(parameters, "lookup", true);
190 if (temp == "not found") {
191 lookupFileName = "LookUp_Titanium.pat";
195 ableToOpen = m->openInputFile(lookupFileName, in, "noerror");
198 //if you can't open it, try input location
199 if (ableToOpen == 1) {
200 if (inputDir != "") { //default path is set
201 string tryPath = inputDir + lookupFileName;
202 m->mothurOut("Unable to open " + lookupFileName + ". Trying input directory " + tryPath); m->mothurOutEndLine();
204 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
206 lookupFileName = tryPath;
210 //if you can't open it, try default location
211 if (ableToOpen == 1) {
212 if (m->getDefaultPath() != "") { //default path is set
213 string tryPath = m->getDefaultPath() + m->getSimpleName(lookupFileName);
214 m->mothurOut("Unable to open " + lookupFileName + ". Trying default " + tryPath); m->mothurOutEndLine();
216 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
218 lookupFileName = tryPath;
222 //if you can't open it its not in current working directory or inputDir, try mothur excutable location
223 if (ableToOpen == 1) {
224 string exepath = m->argv;
225 string tempPath = exepath;
226 for (int i = 0; i < exepath.length(); i++) { tempPath[i] = tolower(exepath[i]); }
227 exepath = exepath.substr(0, (tempPath.find_last_of('m')));
229 string tryPath = m->getFullPathName(exepath) + m->getSimpleName(lookupFileName);
230 m->mothurOut("Unable to open " + lookupFileName + ". Trying mothur's executable location " + tryPath); m->mothurOutEndLine();
232 ableToOpen = m->openInputFile(tryPath, in2, "noerror");
234 lookupFileName = tryPath;
237 if (ableToOpen == 1) { m->mothurOut("Unable to open " + lookupFileName + "."); m->mothurOutEndLine(); abort=true; }
239 else if(temp == "not open") {
241 lookupFileName = validParameter.validFile(parameters, "lookup", false);
243 //if you can't open it its not inputDir, try mothur excutable location
244 string exepath = m->argv;
245 string tempPath = exepath;
246 for (int i = 0; i < exepath.length(); i++) { tempPath[i] = tolower(exepath[i]); }
247 exepath = exepath.substr(0, (tempPath.find_last_of('m')));
249 string tryPath = m->getFullPathName(exepath) + lookupFileName;
250 m->mothurOut("Unable to open " + lookupFileName + ". Trying mothur's executable location " + tryPath); m->mothurOutEndLine();
252 int ableToOpen = m->openInputFile(tryPath, in2, "noerror");
254 lookupFileName = tryPath;
256 if (ableToOpen == 1) { m->mothurOut("Unable to open " + lookupFileName + "."); m->mothurOutEndLine(); abort=true; }
257 }else { lookupFileName = temp; }
259 temp = validParameter.validFile(parameters, "processors", false); if (temp == "not found"){ temp = m->getProcessors(); }
260 m->setProcessors(temp);
261 convert(temp, processors);
263 temp = validParameter.validFile(parameters, "cutoff", false); if (temp == "not found"){ temp = "0.01"; }
264 convert(temp, cutoff);
266 temp = validParameter.validFile(parameters, "mindelta", false); if (temp == "not found"){ temp = "0.000001"; }
267 convert(temp, minDelta);
269 temp = validParameter.validFile(parameters, "maxiter", false); if (temp == "not found"){ temp = "1000"; }
270 convert(temp, maxIters);
272 temp = validParameter.validFile(parameters, "sigma", false);if (temp == "not found") { temp = "60"; }
273 convert(temp, sigma);
275 flowOrder = validParameter.validFile(parameters, "order", false);
276 if (flowOrder == "not found"){ flowOrder = "TACG"; }
277 else if(flowOrder.length() != 4){
278 m->mothurOut("The value of the order option must be four bases long\n");
288 catch(exception& e) {
289 m->errorOut(e, "ShhherCommand", "ShhherCommand");
293 //**********************************************************************************************************************
295 int ShhherCommand::execute(){
297 if (abort == true) { if (calledHelp) { return 0; } return 2; }
304 for(int i=1;i<ncpus;i++){
305 MPI_Send(&abort, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
307 if(abort == 1){ return 0; }
311 m->mothurOut("\nGetting preliminary data...\n");
315 vector<string> flowFileVector;
316 if(flowFilesFileName != "not found"){
319 ifstream flowFilesFile;
320 m->openInputFile(flowFilesFileName, flowFilesFile);
321 while(flowFilesFile){
322 flowFilesFile >> fName;
323 flowFileVector.push_back(fName);
324 m->gobble(flowFilesFile);
328 flowFileVector.push_back(flowFileName);
330 int numFiles = flowFileVector.size();
332 for(int i=1;i<ncpus;i++){
333 MPI_Send(&numFiles, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
336 for(int i=0;i<numFiles;i++){
337 double begClock = clock();
338 unsigned long int begTime = time(NULL);
340 flowFileName = flowFileVector[i];
342 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
343 m->mothurOut("Reading flowgrams...\n");
346 m->mothurOut("Identifying unique flowgrams...\n");
349 m->mothurOut("Calculating distances between flowgrams...\n");
351 strcpy(fileName, flowFileName.c_str());
353 for(int i=1;i<ncpus;i++){
354 MPI_Send(&fileName[0], 1024, MPI_CHAR, i, tag, MPI_COMM_WORLD);
356 MPI_Send(&numSeqs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
357 MPI_Send(&numUniques, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
358 MPI_Send(&numFlowCells, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
359 MPI_Send(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, i, tag, MPI_COMM_WORLD);
360 MPI_Send(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
361 MPI_Send(&mapUniqueToSeq[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
362 MPI_Send(&mapSeqToUnique[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
363 MPI_Send(&lengths[0], numSeqs, MPI_INT, i, tag, MPI_COMM_WORLD);
364 MPI_Send(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
365 MPI_Send(&cutoff, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
368 string distFileName = flowDistMPI(0, int(sqrt(1.0/float(ncpus)) * numUniques));
371 for(int i=1;i<ncpus;i++){
372 MPI_Recv(&done, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
374 m->appendFiles((distFileName + ".temp." + toString(i)), distFileName);
375 remove((distFileName + ".temp." + toString(i)).c_str());
378 string namesFileName = createNamesFile();
380 m->mothurOut("\nClustering flowgrams...\n");
381 string listFileName = cluster(distFileName, namesFileName);
383 getOTUData(listFileName);
385 remove(distFileName.c_str());
386 remove(namesFileName.c_str());
387 remove(listFileName.c_str());
391 for(int i=1;i<ncpus;i++){
392 MPI_Send(&numOTUs, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
393 MPI_Send(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
394 MPI_Send(&uniqueFlowgrams[0], numFlowCells * numUniques, MPI_SHORT, i, tag, MPI_COMM_WORLD);
395 MPI_Send(&sigma, 1, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
402 int numOTUsOnCPU = numOTUs / ncpus;
403 int numSeqsOnCPU = numSeqs / ncpus;
404 m->mothurOut("\nDenoising flowgrams...\n");
405 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
407 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
409 double cycClock = clock();
410 unsigned long int cycTime = time(NULL);
413 int total = singleTau.size();
414 for(int i=1;i<ncpus;i++){
415 MPI_Send(&total, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
416 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
417 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
419 MPI_Send(&singleTau[0], total, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
420 MPI_Send(&seqNumber[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
421 MPI_Send(&seqIndex[0], total, MPI_INT, i, tag, MPI_COMM_WORLD);
422 MPI_Send(&nSeqsPerOTU[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
423 MPI_Send(&cumNumSeqs[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
426 calcCentroidsDriver(0, numOTUsOnCPU);
428 for(int i=1;i<ncpus;i++){
429 int otuStart = i * numOTUs / ncpus;
430 int otuStop = (i + 1) * numOTUs / ncpus;
432 vector<int> tempCentroids(numOTUs, 0);
433 vector<short> tempChange(numOTUs, 0);
435 MPI_Recv(&tempCentroids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
436 MPI_Recv(&tempChange[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD, &status);
438 for(int j=otuStart;j<otuStop;j++){
439 centroids[j] = tempCentroids[j];
440 change[j] = tempChange[j];
444 maxDelta = getNewWeights();
445 double nLL = getLikelihood();
448 for(int i=1;i<ncpus;i++){
449 MPI_Send(¢roids[0], numOTUs, MPI_INT, i, tag, MPI_COMM_WORLD);
450 MPI_Send(&weight[0], numOTUs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD);
451 MPI_Send(&change[0], numOTUs, MPI_SHORT, i, tag, MPI_COMM_WORLD);
454 calcNewDistancesParent(0, numSeqsOnCPU);
456 total = singleTau.size();
458 for(int i=1;i<ncpus;i++){
460 int seqStart = i * numSeqs / ncpus;
461 int seqStop = (i + 1) * numSeqs / ncpus;
463 MPI_Recv(&childTotal, 1, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
465 vector<int> childSeqIndex(childTotal, 0);
466 vector<double> childSingleTau(childTotal, 0);
467 vector<double> childDist(numSeqs * numOTUs, 0);
468 vector<int> otuIndex(childTotal, 0);
470 MPI_Recv(&childSeqIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
471 MPI_Recv(&childSingleTau[0], childTotal, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
472 MPI_Recv(&childDist[0], numOTUs * numSeqs, MPI_DOUBLE, i, tag, MPI_COMM_WORLD, &status);
473 MPI_Recv(&otuIndex[0], childTotal, MPI_INT, i, tag, MPI_COMM_WORLD, &status);
475 int oldTotal = total;
477 singleTau.resize(total, 0);
478 seqIndex.resize(total, 0);
479 seqNumber.resize(total, 0);
483 for(int j=oldTotal;j<total;j++){
484 int otuI = otuIndex[childIndex];
485 int seqI = childSeqIndex[childIndex];
487 singleTau[j] = childSingleTau[childIndex];
489 aaP[otuI][nSeqsPerOTU[otuI]] = j;
490 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
495 int index = seqStart * numOTUs;
496 for(int j=seqStart;j<seqStop;j++){
497 for(int k=0;k<numOTUs;k++){
498 dist[index] = childDist[index];
506 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
508 if((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
510 for(int i=1;i<ncpus;i++){
511 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD);
516 for(int i=1;i<ncpus;i++){
517 MPI_Send(&live, 1, MPI_INT, i, tag, MPI_COMM_WORLD); //send kill command
523 m->mothurOut("\nFinalizing...\n");
527 vector<int> otuCounts(numOTUs, 0);
528 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
529 calcCentroidsDriver(0, numOTUs);
531 writeQualities(otuCounts);
532 writeSequences(otuCounts);
533 writeNames(otuCounts);
534 writeClusters(otuCounts);
538 m->mothurOut("Total time to process " + toString(flowFileName) + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
544 MPI_Recv(&abort, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
545 if(abort){ return 0; }
548 MPI_Recv(&numFiles, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
550 for(int i=0;i<numFiles;i++){
551 //Now into the pyrodist part
555 MPI_Recv(&fileName, 1024, MPI_CHAR, 0, tag, MPI_COMM_WORLD, &status);
556 MPI_Recv(&numSeqs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
557 MPI_Recv(&numUniques, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
558 MPI_Recv(&numFlowCells, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
560 flowDataIntI.resize(numSeqs * numFlowCells);
561 flowDataPrI.resize(numSeqs * numFlowCells);
562 mapUniqueToSeq.resize(numSeqs);
563 mapSeqToUnique.resize(numSeqs);
564 lengths.resize(numSeqs);
565 jointLookUp.resize(NUMBINS * NUMBINS);
567 MPI_Recv(&flowDataIntI[0], numSeqs * numFlowCells, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
568 MPI_Recv(&flowDataPrI[0], numSeqs * numFlowCells, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
569 MPI_Recv(&mapUniqueToSeq[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
570 MPI_Recv(&mapSeqToUnique[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
571 MPI_Recv(&lengths[0], numSeqs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
572 MPI_Recv(&jointLookUp[0], NUMBINS * NUMBINS, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
573 MPI_Recv(&cutoff, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
575 flowFileName = string(fileName);
576 int flowDistStart = int(sqrt(float(pid)/float(ncpus)) * numUniques);
577 int flowDistEnd = int(sqrt(float(pid+1)/float(ncpus)) * numUniques);
579 string distanceStringChild = flowDistMPI(flowDistStart, flowDistEnd);
582 MPI_Send(&done, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
584 //Now into the pyronoise part
585 MPI_Recv(&numOTUs, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
587 singleLookUp.resize(HOMOPS * NUMBINS);
588 uniqueFlowgrams.resize(numUniques * numFlowCells);
589 weight.resize(numOTUs);
590 centroids.resize(numOTUs);
591 change.resize(numOTUs);
592 dist.assign(numOTUs * numSeqs, 0);
593 nSeqsPerOTU.resize(numOTUs);
594 cumNumSeqs.resize(numOTUs);
596 MPI_Recv(&singleLookUp[0], singleLookUp.size(), MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
597 MPI_Recv(&uniqueFlowgrams[0], uniqueFlowgrams.size(), MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
598 MPI_Recv(&sigma, 1, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
600 int startOTU = pid * numOTUs / ncpus;
601 int endOTU = (pid + 1) * numOTUs / ncpus;
603 int startSeq = pid * numSeqs / ncpus;
604 int endSeq = (pid + 1) * numSeqs /ncpus;
610 MPI_Recv(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
611 singleTau.assign(total, 0.0000);
612 seqNumber.assign(total, 0);
613 seqIndex.assign(total, 0);
615 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
616 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
617 MPI_Recv(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
618 MPI_Recv(&seqNumber[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
619 MPI_Recv(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
620 MPI_Recv(&nSeqsPerOTU[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
621 MPI_Recv(&cumNumSeqs[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
623 calcCentroidsDriver(startOTU, endOTU);
625 MPI_Send(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD);
626 MPI_Send(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD);
628 MPI_Recv(¢roids[0], numOTUs, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
629 MPI_Recv(&weight[0], numOTUs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD, &status);
630 MPI_Recv(&change[0], numOTUs, MPI_SHORT, 0, tag, MPI_COMM_WORLD, &status);
632 vector<int> otuIndex(total, 0);
633 calcNewDistancesChildMPI(startSeq, endSeq, otuIndex);
634 total = otuIndex.size();
636 MPI_Send(&total, 1, MPI_INT, 0, tag, MPI_COMM_WORLD);
637 MPI_Send(&seqIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
638 MPI_Send(&singleTau[0], total, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
639 MPI_Send(&dist[0], numOTUs * numSeqs, MPI_DOUBLE, 0, tag, MPI_COMM_WORLD);
640 MPI_Send(&otuIndex[0], total, MPI_INT, 0, tag, MPI_COMM_WORLD);
642 MPI_Recv(&live, 1, MPI_INT, 0, tag, MPI_COMM_WORLD, &status);
646 MPI_Barrier(MPI_COMM_WORLD);
649 if(compositeFASTAFileName != ""){
650 outputNames.push_back(compositeFASTAFileName);
651 outputNames.push_back(compositeNamesFileName);
654 m->mothurOutEndLine();
655 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
656 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
657 m->mothurOutEndLine();
662 catch(exception& e) {
663 m->errorOut(e, "ShhherCommand", "execute");
668 /**************************************************************************************************/
670 string ShhherCommand::flowDistMPI(int startSeq, int stopSeq){
672 ostringstream outStream;
673 outStream.setf(ios::fixed, ios::floatfield);
674 outStream.setf(ios::dec, ios::basefield);
675 outStream.setf(ios::showpoint);
676 outStream.precision(6);
678 int begTime = time(NULL);
679 double begClock = clock();
681 for(int i=startSeq;i<stopSeq;i++){
682 for(int j=0;j<i;j++){
683 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
685 if(flowDistance < 1e-6){
686 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
688 else if(flowDistance <= cutoff){
689 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
693 m->mothurOut(toString(i) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
697 m->mothurOut(toString(stopSeq) + '\t' + toString(time(NULL) - begTime) + '\t' + toString((clock()-begClock)/CLOCKS_PER_SEC) + '\n');
699 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.dist";
700 if(pid != 0){ fDistFileName += ".temp." + toString(pid); }
702 ofstream distFile(fDistFileName.c_str());
703 distFile << outStream.str();
706 return fDistFileName;
708 catch(exception& e) {
709 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
715 //**********************************************************************************************************************
717 int ShhherCommand::execute(){
719 if (abort == true) { return 0; }
724 vector<string> flowFileVector;
725 if(flowFilesFileName != "not found"){
728 ifstream flowFilesFile;
729 m->openInputFile(flowFilesFileName, flowFilesFile);
730 while(flowFilesFile){
731 flowFilesFile >> fName;
732 flowFileVector.push_back(fName);
733 m->gobble(flowFilesFile);
737 flowFileVector.push_back(flowFileName);
739 int numFiles = flowFileVector.size();
742 for(int i=0;i<numFiles;i++){
743 flowFileName = flowFileVector[i];
745 m->mothurOut("\n>>>>>\tProcessing " + flowFileName + " (file " + toString(i+1) + " of " + toString(numFiles) + ")\t<<<<<\n");
746 m->mothurOut("Reading flowgrams...\n");
749 m->mothurOut("Identifying unique flowgrams...\n");
753 m->mothurOut("Calculating distances between flowgrams...\n");
754 string distFileName = createDistFile(processors);
755 string namesFileName = createNamesFile();
757 m->mothurOut("\nClustering flowgrams...\n");
758 string listFileName = cluster(distFileName, namesFileName);
759 getOTUData(listFileName);
761 remove(distFileName.c_str());
762 remove(namesFileName.c_str());
763 remove(listFileName.c_str());
770 double begClock = clock();
771 unsigned long int begTime = time(NULL);
773 m->mothurOut("\nDenoising flowgrams...\n");
774 m->mothurOut("iter\tmaxDelta\tnLL\t\tcycletime\n");
776 while((maxIters == 0 && maxDelta > minDelta) || iter < MIN_ITER || (maxDelta > minDelta && iter < maxIters)){
778 double cycClock = clock();
779 unsigned long int cycTime = time(NULL);
784 maxDelta = getNewWeights();
785 double nLL = getLikelihood();
792 m->mothurOut(toString(iter) + '\t' + toString(maxDelta) + '\t' + toString(nLL) + '\t' + toString(time(NULL) - cycTime) + '\t' + toString((clock() - cycClock)/(double)CLOCKS_PER_SEC) + '\n');
796 m->mothurOut("\nFinalizing...\n");
800 vector<int> otuCounts(numOTUs, 0);
801 for(int i=0;i<numSeqs;i++) { otuCounts[otuData[i]]++; }
803 calcCentroidsDriver(0, numOTUs);
804 writeQualities(otuCounts);
805 writeSequences(otuCounts);
806 writeNames(otuCounts);
807 writeClusters(otuCounts);
810 m->mothurOut("Total time to process " + flowFileName + ":\t" + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/(double)CLOCKS_PER_SEC) + '\n');
813 if(compositeFASTAFileName != ""){
814 outputNames.push_back(compositeFASTAFileName);
815 outputNames.push_back(compositeNamesFileName);
818 m->mothurOutEndLine();
819 m->mothurOut("Output File Names: "); m->mothurOutEndLine();
820 for (int i = 0; i < outputNames.size(); i++) { m->mothurOut(outputNames[i]); m->mothurOutEndLine(); }
821 m->mothurOutEndLine();
825 catch(exception& e) {
826 m->errorOut(e, "ShhherCommand", "execute");
831 /**************************************************************************************************/
833 void ShhherCommand::getFlowData(){
836 m->openInputFile(flowFileName, flowFile);
839 seqNameVector.clear();
841 flowDataIntI.clear();
845 int currentNumFlowCells;
849 flowFile >> numFlowCells;
850 int index = 0;//pcluster
851 while(!flowFile.eof()){
852 flowFile >> seqName >> currentNumFlowCells;
853 lengths.push_back(currentNumFlowCells);
855 seqNameVector.push_back(seqName);
856 nameMap[seqName] = index++;//pcluster
858 for(int i=0;i<numFlowCells;i++){
859 flowFile >> intensity;
860 if(intensity > 9.99) { intensity = 9.99; }
861 int intI = int(100 * intensity + 0.0001);
862 flowDataIntI.push_back(intI);
868 numSeqs = seqNameVector.size();
870 for(int i=0;i<numSeqs;i++){
871 int iNumFlowCells = i * numFlowCells;
872 for(int j=lengths[i];j<numFlowCells;j++){
873 flowDataIntI[iNumFlowCells + j] = 0;
878 catch(exception& e) {
879 m->errorOut(e, "ShhherCommand", "getFlowData");
884 /**************************************************************************************************/
886 void ShhherCommand::getSingleLookUp(){
888 // these are the -log probabilities that a signal corresponds to a particular homopolymer length
889 singleLookUp.assign(HOMOPS * NUMBINS, 0);
893 m->openInputFile(lookupFileName, lookUpFile);
895 for(int i=0;i<HOMOPS;i++){
897 lookUpFile >> logFracFreq;
899 for(int j=0;j<NUMBINS;j++) {
900 lookUpFile >> singleLookUp[index];
906 catch(exception& e) {
907 m->errorOut(e, "ShhherCommand", "getSingleLookUp");
912 /**************************************************************************************************/
914 void ShhherCommand::getJointLookUp(){
917 // the most likely joint probability (-log) that two intenities have the same polymer length
918 jointLookUp.resize(NUMBINS * NUMBINS, 0);
920 for(int i=0;i<NUMBINS;i++){
921 for(int j=0;j<NUMBINS;j++){
923 double minSum = 100000000;
925 for(int k=0;k<HOMOPS;k++){
926 double sum = singleLookUp[k * NUMBINS + i] + singleLookUp[k * NUMBINS + j];
928 if(sum < minSum) { minSum = sum; }
930 jointLookUp[i * NUMBINS + j] = minSum;
934 catch(exception& e) {
935 m->errorOut(e, "ShhherCommand", "getJointLookUp");
940 /**************************************************************************************************/
942 double ShhherCommand::getProbIntensity(int intIntensity){
945 double minNegLogProb = 100000000;
948 for(int i=0;i<HOMOPS;i++){//loop signal strength
949 float negLogProb = singleLookUp[i * NUMBINS + intIntensity];
950 if(negLogProb < minNegLogProb) { minNegLogProb = negLogProb; }
953 return minNegLogProb;
955 catch(exception& e) {
956 m->errorOut(e, "ShhherCommand", "getProbIntensity");
961 /**************************************************************************************************/
963 void ShhherCommand::getUniques(){
968 uniqueFlowgrams.assign(numFlowCells * numSeqs, -1);
969 uniqueCount.assign(numSeqs, 0); // anWeights
970 uniqueLengths.assign(numSeqs, 0);
971 mapSeqToUnique.assign(numSeqs, -1);
972 mapUniqueToSeq.assign(numSeqs, -1);
974 vector<short> uniqueFlowDataIntI(numFlowCells * numSeqs, -1);
976 for(int i=0;i<numSeqs;i++){
979 vector<short> current(numFlowCells);
980 for(int j=0;j<numFlowCells;j++){
981 current[j] = short(((flowDataIntI[i * numFlowCells + j] + 50.0)/100.0));
984 for(int j=0;j<numUniques;j++){
985 int offset = j * numFlowCells;
989 if(lengths[i] < uniqueLengths[j]) { shorterLength = lengths[i]; }
990 else { shorterLength = uniqueLengths[j]; }
992 for(int k=0;k<shorterLength;k++){
993 if(current[k] != uniqueFlowgrams[offset + k]){
1000 mapSeqToUnique[i] = j;
1003 if(lengths[i] > uniqueLengths[j]) { uniqueLengths[j] = lengths[i]; }
1009 if(index == numUniques){
1010 uniqueLengths[numUniques] = lengths[i];
1011 uniqueCount[numUniques] = 1;
1012 mapSeqToUnique[i] = numUniques;//anMap
1013 mapUniqueToSeq[numUniques] = i;//anF
1015 for(int k=0;k<numFlowCells;k++){
1016 uniqueFlowgrams[numUniques * numFlowCells + k] = current[k];
1017 uniqueFlowDataIntI[numUniques * numFlowCells + k] = flowDataIntI[i * numFlowCells + k];
1023 uniqueFlowDataIntI.resize(numFlowCells * numUniques);
1024 uniqueLengths.resize(numUniques);
1026 flowDataPrI.resize(numSeqs * numFlowCells, 0);
1027 for(int i=0;i<flowDataPrI.size();i++) { flowDataPrI[i] = getProbIntensity(flowDataIntI[i]); }
1029 catch(exception& e) {
1030 m->errorOut(e, "ShhherCommand", "getUniques");
1035 /**************************************************************************************************/
1037 float ShhherCommand::calcPairwiseDist(int seqA, int seqB){
1039 int minLength = lengths[mapSeqToUnique[seqA]];
1040 if(lengths[seqB] < minLength){ minLength = lengths[mapSeqToUnique[seqB]]; }
1042 int ANumFlowCells = seqA * numFlowCells;
1043 int BNumFlowCells = seqB * numFlowCells;
1047 for(int i=0;i<minLength;i++){
1048 int flowAIntI = flowDataIntI[ANumFlowCells + i];
1049 float flowAPrI = flowDataPrI[ANumFlowCells + i];
1051 int flowBIntI = flowDataIntI[BNumFlowCells + i];
1052 float flowBPrI = flowDataPrI[BNumFlowCells + i];
1053 dist += jointLookUp[flowAIntI * NUMBINS + flowBIntI] - flowAPrI - flowBPrI;
1056 dist /= (float) minLength;
1059 catch(exception& e) {
1060 m->errorOut(e, "ShhherCommand", "calcPairwiseDist");
1065 /**************************************************************************************************/
1067 void ShhherCommand::flowDistParentFork(string distFileName, int startSeq, int stopSeq){
1070 ostringstream outStream;
1071 outStream.setf(ios::fixed, ios::floatfield);
1072 outStream.setf(ios::dec, ios::basefield);
1073 outStream.setf(ios::showpoint);
1074 outStream.precision(6);
1076 int begTime = time(NULL);
1077 double begClock = clock();
1079 for(int i=startSeq;i<stopSeq;i++){
1080 for(int j=0;j<i;j++){
1081 float flowDistance = calcPairwiseDist(mapUniqueToSeq[i], mapUniqueToSeq[j]);
1083 if(flowDistance < 1e-6){
1084 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << 0.000000 << endl;
1086 else if(flowDistance <= cutoff){
1087 outStream << mapUniqueToSeq[i] << '\t' << mapUniqueToSeq[j] << '\t' << flowDistance << endl;
1091 m->mothurOut(toString(i) + "\t" + toString(time(NULL) - begTime));
1092 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1093 m->mothurOutEndLine();
1096 m->mothurOut(toString(stopSeq-1) + "\t" + toString(time(NULL) - begTime));
1097 m->mothurOut("\t" + toString((clock()-begClock)/CLOCKS_PER_SEC));
1098 m->mothurOutEndLine();
1100 ofstream distFile(distFileName.c_str());
1101 distFile << outStream.str();
1104 catch(exception& e) {
1105 m->errorOut(e, "ShhherCommand", "flowDistParentFork");
1110 /**************************************************************************************************/
1112 string ShhherCommand::createDistFile(int processors){
1114 string fDistFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.dist";
1116 unsigned long int begTime = time(NULL);
1117 double begClock = clock();
1122 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1123 if(processors == 1) { flowDistParentFork(fDistFileName, 0, numUniques); }
1124 else{ //you have multiple processors
1126 if (numSeqs < processors){ processors = 1; }
1128 vector<int> start(processors, 0);
1129 vector<int> end(processors, 0);
1131 for (int i = 0; i < processors; i++) {
1132 start[i] = int(sqrt(float(i)/float(processors)) * numUniques);
1133 end[i] = int(sqrt(float(i+1)/float(processors)) * numUniques);
1137 vector<int> processIDs;
1139 //loop through and create all the processes you want
1140 while (process != processors) {
1144 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1146 }else if (pid == 0){
1147 flowDistParentFork(fDistFileName + toString(getpid()) + ".temp", start[process], end[process]);
1150 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1152 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1157 //parent does its part
1158 flowDistParentFork(fDistFileName, start[0], end[0]);
1160 //force parent to wait until all the processes are done
1161 for (int i=0;i<processIDs.size();i++) {
1162 int temp = processIDs[i];
1166 //append and remove temp files
1167 for (int i=0;i<processIDs.size();i++) {
1168 m->appendFiles((fDistFileName + toString(processIDs[i]) + ".temp"), fDistFileName);
1169 remove((fDistFileName + toString(processIDs[i]) + ".temp").c_str());
1175 flowDistParentFork(fDistFileName, 0, numUniques);
1178 m->mothurOutEndLine();
1180 m->mothurOut("Total time: " + toString(time(NULL) - begTime) + '\t' + toString((clock() - begClock)/CLOCKS_PER_SEC) + '\n');
1183 return fDistFileName;
1185 catch(exception& e) {
1186 m->errorOut(e, "ShhherCommand", "createDistFile");
1192 /**************************************************************************************************/
1194 string ShhherCommand::createNamesFile(){
1197 vector<string> duplicateNames(numUniques, "");
1198 for(int i=0;i<numSeqs;i++){
1199 duplicateNames[mapSeqToUnique[i]] += seqNameVector[i] + ',';
1202 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.names";
1205 m->openOutputFile(nameFileName, nameFile);
1207 for(int i=0;i<numUniques;i++){
1208 // nameFile << seqNameVector[mapUniqueToSeq[i]] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1209 nameFile << mapUniqueToSeq[i] << '\t' << duplicateNames[i].substr(0, duplicateNames[i].find_last_of(',')) << endl;
1213 return nameFileName;
1215 catch(exception& e) {
1216 m->errorOut(e, "ShhherCommand", "createNamesFile");
1221 //**********************************************************************************************************************
1223 string ShhherCommand::cluster(string distFileName, string namesFileName){
1226 ReadMatrix* read = new ReadColumnMatrix(distFileName);
1227 read->setCutoff(cutoff);
1229 NameAssignment* clusterNameMap = new NameAssignment(namesFileName);
1230 clusterNameMap->readMap();
1231 read->read(clusterNameMap);
1233 ListVector* list = read->getListVector();
1234 SparseMatrix* matrix = read->getMatrix();
1237 delete clusterNameMap;
1239 RAbundVector* rabund = new RAbundVector(list->getRAbundVector());
1241 Cluster* cluster = new CompleteLinkage(rabund, list, matrix, cutoff, "furthest");
1242 string tag = cluster->getTag();
1244 double clusterCutoff = cutoff;
1245 while (matrix->getSmallDist() <= clusterCutoff && matrix->getNNodes() > 0){
1246 cluster->update(clusterCutoff);
1249 list->setLabel(toString(cutoff));
1251 string listFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.list";
1253 m->openOutputFile(listFileName, listFile);
1254 list->print(listFile);
1257 delete matrix; delete cluster; delete rabund; delete list;
1259 return listFileName;
1261 catch(exception& e) {
1262 m->errorOut(e, "ShhherCommand", "cluster");
1267 /**************************************************************************************************/
1269 void ShhherCommand::getOTUData(string listFileName){
1273 m->openInputFile(listFileName, listFile);
1276 listFile >> label >> numOTUs;
1278 otuData.assign(numSeqs, 0);
1279 cumNumSeqs.assign(numOTUs, 0);
1280 nSeqsPerOTU.assign(numOTUs, 0);
1281 aaP.clear();aaP.resize(numOTUs);
1287 string singleOTU = "";
1289 for(int i=0;i<numOTUs;i++){
1291 listFile >> singleOTU;
1293 istringstream otuString(singleOTU);
1297 string seqName = "";
1299 for(int j=0;j<singleOTU.length();j++){
1300 char letter = otuString.get();
1306 map<string,int>::iterator nmIt = nameMap.find(seqName);
1307 int index = nmIt->second;
1309 nameMap.erase(nmIt);
1313 aaP[i].push_back(index);
1318 map<string,int>::iterator nmIt = nameMap.find(seqName);
1320 int index = nmIt->second;
1321 nameMap.erase(nmIt);
1325 aaP[i].push_back(index);
1330 sort(aaP[i].begin(), aaP[i].end());
1331 for(int j=0;j<nSeqsPerOTU[i];j++){
1332 seqNumber.push_back(aaP[i][j]);
1334 for(int j=nSeqsPerOTU[i];j<numSeqs;j++){
1335 aaP[i].push_back(0);
1341 for(int i=1;i<numOTUs;i++){
1342 cumNumSeqs[i] = cumNumSeqs[i-1] + nSeqsPerOTU[i-1];
1345 seqIndex = seqNumber;
1350 catch(exception& e) {
1351 m->errorOut(e, "ShhherCommand", "getOTUData");
1356 /**************************************************************************************************/
1358 void ShhherCommand::initPyroCluster(){
1360 dist.assign(numSeqs * numOTUs, 0);
1361 change.assign(numOTUs, 1);
1362 centroids.assign(numOTUs, -1);
1363 weight.assign(numOTUs, 0);
1364 singleTau.assign(numSeqs, 1.0);
1366 nSeqsBreaks.assign(processors+1, 0);
1367 nOTUsBreaks.assign(processors+1, 0);
1370 for(int i=0;i<processors;i++){
1371 nSeqsBreaks[i+1] = nSeqsBreaks[i] + (int)((double) numSeqs / (double) processors);
1372 nOTUsBreaks[i+1] = nOTUsBreaks[i] + (int)((double) numOTUs / (double) processors);
1374 nSeqsBreaks[processors] = numSeqs;
1375 nOTUsBreaks[processors] = numOTUs;
1377 catch(exception& e) {
1378 m->errorOut(e, "ShhherCommand", "initPyroCluster");
1383 /**************************************************************************************************/
1385 void ShhherCommand::fill(){
1388 for(int i=0;i<numOTUs;i++){
1389 cumNumSeqs[i] = index;
1390 for(int j=0;j<nSeqsPerOTU[i];j++){
1391 seqNumber[index] = aaP[i][j];
1392 seqIndex[index] = aaI[i][j];
1398 catch(exception& e) {
1399 m->errorOut(e, "ShhherCommand", "fill");
1404 /**************************************************************************************************/
1406 void ShhherCommand::calcCentroids(){
1409 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1411 if(processors == 1) {
1412 calcCentroidsDriver(0, numOTUs);
1414 else{ //you have multiple processors
1415 if (numOTUs < processors){ processors = 1; }
1418 vector<int> processIDs;
1420 //loop through and create all the processes you want
1421 while (process != processors) {
1425 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1427 }else if (pid == 0){
1428 calcCentroidsDriver(nOTUsBreaks[process], nOTUsBreaks[process+1]);
1431 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1433 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1438 //parent does its part
1439 calcCentroidsDriver(nOTUsBreaks[0], nOTUsBreaks[1]);
1441 //force parent to wait until all the processes are done
1442 for (int i=0;i<processIDs.size();i++) {
1443 int temp = processIDs[i];
1449 calcCentroidsDriver(0, numOTUs);
1452 catch(exception& e) {
1453 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1458 /**************************************************************************************************/
1460 void ShhherCommand::calcCentroidsDriver(int start, int finish){
1462 //this function gets the most likely homopolymer length at a flow position for a group of sequences
1468 for(int i=start;i<finish;i++){
1472 int minFlowGram = 100000000;
1473 double minFlowValue = 1e8;
1474 change[i] = 0; //FALSE
1476 for(int j=0;j<nSeqsPerOTU[i];j++){
1477 count += singleTau[seqNumber[cumNumSeqs[i] + j]];
1480 if(nSeqsPerOTU[i] > 0 && count > MIN_COUNT){
1481 vector<double> adF(nSeqsPerOTU[i]);
1482 vector<int> anL(nSeqsPerOTU[i]);
1484 for(int j=0;j<nSeqsPerOTU[i];j++){
1485 int index = cumNumSeqs[i] + j;
1486 int nI = seqIndex[index];
1487 int nIU = mapSeqToUnique[nI];
1490 for(k=0;k<position;k++){
1496 anL[position] = nIU;
1497 adF[position] = 0.0000;
1502 for(int j=0;j<nSeqsPerOTU[i];j++){
1503 int index = cumNumSeqs[i] + j;
1504 int nI = seqIndex[index];
1506 double tauValue = singleTau[seqNumber[index]];
1508 for(int k=0;k<position;k++){
1509 double dist = getDistToCentroid(anL[k], nI, lengths[nI]);
1510 adF[k] += dist * tauValue;
1514 for(int j=0;j<position;j++){
1515 if(adF[j] < minFlowValue){
1517 minFlowValue = adF[j];
1521 if(centroids[i] != anL[minFlowGram]){
1523 centroids[i] = anL[minFlowGram];
1526 else if(centroids[i] != -1){
1532 catch(exception& e) {
1533 m->errorOut(e, "ShhherCommand", "calcCentroidsDriver");
1538 /**************************************************************************************************/
1540 double ShhherCommand::getDistToCentroid(int cent, int flow, int length){
1543 int flowAValue = cent * numFlowCells;
1544 int flowBValue = flow * numFlowCells;
1548 for(int i=0;i<length;i++){
1549 dist += singleLookUp[uniqueFlowgrams[flowAValue] * NUMBINS + flowDataIntI[flowBValue]];
1554 return dist / (double)length;
1556 catch(exception& e) {
1557 m->errorOut(e, "ShhherCommand", "getDistToCentroid");
1562 /**************************************************************************************************/
1564 double ShhherCommand::getNewWeights(){
1567 double maxChange = 0;
1569 for(int i=0;i<numOTUs;i++){
1571 double difference = weight[i];
1574 for(int j=0;j<nSeqsPerOTU[i];j++){
1575 int index = cumNumSeqs[i] + j;
1576 double tauValue = singleTau[seqNumber[index]];
1577 weight[i] += tauValue;
1580 difference = fabs(weight[i] - difference);
1581 if(difference > maxChange){ maxChange = difference; }
1585 catch(exception& e) {
1586 m->errorOut(e, "ShhherCommand", "getNewWeights");
1591 /**************************************************************************************************/
1593 double ShhherCommand::getLikelihood(){
1597 vector<long double> P(numSeqs, 0);
1600 for(int i=0;i<numOTUs;i++){
1601 if(weight[i] > MIN_WEIGHT){
1607 for(int i=0;i<numOTUs;i++){
1608 for(int j=0;j<nSeqsPerOTU[i];j++){
1609 int index = cumNumSeqs[i] + j;
1610 int nI = seqIndex[index];
1611 double singleDist = dist[seqNumber[index]];
1613 P[nI] += weight[i] * exp(-singleDist * sigma);
1617 for(int i=0;i<numSeqs;i++){
1618 if(P[i] == 0){ P[i] = DBL_EPSILON; }
1623 nLL = nLL -(double)numSeqs * log(sigma);
1627 catch(exception& e) {
1628 m->errorOut(e, "ShhherCommand", "getNewWeights");
1633 /**************************************************************************************************/
1635 void ShhherCommand::checkCentroids(){
1637 vector<int> unique(numOTUs, 1);
1639 for(int i=0;i<numOTUs;i++){
1640 if(centroids[i] == -1 || weight[i] < MIN_WEIGHT){
1645 for(int i=0;i<numOTUs;i++){
1647 for(int j=i+1;j<numOTUs;j++){
1650 if(centroids[j] == centroids[i]){
1654 weight[i] += weight[j];
1662 catch(exception& e) {
1663 m->errorOut(e, "ShhherCommand", "checkCentroids");
1668 /**************************************************************************************************/
1670 void ShhherCommand::calcNewDistances(){
1673 #if defined (__APPLE__) || (__MACH__) || (linux) || (__linux)
1675 if(processors == 1) {
1676 calcNewDistancesParent(0, numSeqs);
1678 else{ //you have multiple processors
1679 if (numSeqs < processors){ processors = 1; }
1681 vector<vector<int> > child_otuIndex(processors);
1682 vector<vector<int> > child_seqIndex(processors);
1683 vector<vector<double> > child_singleTau(processors);
1684 vector<int> totals(processors);
1687 vector<int> processIDs;
1689 //loop through and create all the processes you want
1690 while (process != processors) {
1694 processIDs.push_back(pid); //create map from line number to pid so you can append files in correct order later
1696 }else if (pid == 0){
1697 calcNewDistancesChild(nSeqsBreaks[process], nSeqsBreaks[process+1], child_otuIndex[process], child_seqIndex[process], child_singleTau[process]);
1698 totals[process] = child_otuIndex[process].size();
1702 m->mothurOut("[ERROR]: unable to spawn the necessary processes. Error code: " + toString(pid)); m->mothurOutEndLine();
1704 for (int i=0;i<processIDs.size();i++) { int temp = processIDs[i]; kill (temp, SIGINT); }
1709 //parent does its part
1710 calcNewDistancesParent(nSeqsBreaks[0], nSeqsBreaks[1]);
1711 int total = seqIndex.size();
1713 //force parent to wait until all the processes are done
1714 for (int i=0;i<processIDs.size();i++) {
1715 int temp = processIDs[i];
1719 for(int i=1;i<processors;i++){
1720 int oldTotal = total;
1723 singleTau.resize(total, 0);
1724 seqIndex.resize(total, 0);
1725 seqNumber.resize(total, 0);
1729 for(int j=oldTotal;j<total;j++){
1730 int otuI = child_otuIndex[i][childIndex];
1731 int seqI = child_seqIndex[i][childIndex];
1733 singleTau[j] = child_singleTau[i][childIndex];
1734 aaP[otuI][nSeqsPerOTU[otuI]] = j;
1735 aaI[otuI][nSeqsPerOTU[otuI]] = seqI;
1736 nSeqsPerOTU[otuI]++;
1743 calcNewDistancesParent(0, numSeqs);
1746 catch(exception& e) {
1747 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1752 /**************************************************************************************************/
1754 void ShhherCommand::calcNewDistancesChildMPI(int startSeq, int stopSeq, vector<int>& otuIndex){
1757 vector<double> newTau(numOTUs,0);
1758 vector<double> norms(numSeqs, 0);
1765 for(int i=startSeq;i<stopSeq;i++){
1766 double offset = 1e8;
1767 int indexOffset = i * numOTUs;
1769 for(int j=0;j<numOTUs;j++){
1771 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1772 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1774 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1775 offset = dist[indexOffset + j];
1779 for(int j=0;j<numOTUs;j++){
1780 if(weight[j] > MIN_WEIGHT){
1781 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1782 norms[i] += newTau[j];
1789 for(int j=0;j<numOTUs;j++){
1791 newTau[j] /= norms[i];
1793 if(newTau[j] > MIN_TAU){
1794 otuIndex.push_back(j);
1795 seqIndex.push_back(i);
1796 singleTau.push_back(newTau[j]);
1802 catch(exception& e) {
1803 m->errorOut(e, "ShhherCommand", "calcNewDistancesChildMPI");
1808 /**************************************************************************************************/
1810 void ShhherCommand::calcNewDistancesChild(int startSeq, int stopSeq, vector<int>& child_otuIndex, vector<int>& child_seqIndex, vector<double>& child_singleTau){
1813 vector<double> newTau(numOTUs,0);
1814 vector<double> norms(numSeqs, 0);
1815 child_otuIndex.resize(0);
1816 child_seqIndex.resize(0);
1817 child_singleTau.resize(0);
1819 for(int i=startSeq;i<stopSeq;i++){
1820 double offset = 1e8;
1821 int indexOffset = i * numOTUs;
1824 for(int j=0;j<numOTUs;j++){
1825 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1826 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1829 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1830 offset = dist[indexOffset + j];
1834 for(int j=0;j<numOTUs;j++){
1835 if(weight[j] > MIN_WEIGHT){
1836 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1837 norms[i] += newTau[j];
1844 for(int j=0;j<numOTUs;j++){
1845 newTau[j] /= norms[i];
1847 if(newTau[j] > MIN_TAU){
1848 child_otuIndex.push_back(j);
1849 child_seqIndex.push_back(i);
1850 child_singleTau.push_back(newTau[j]);
1855 catch(exception& e) {
1856 m->errorOut(e, "ShhherCommand", "calcNewDistancesChild");
1861 /**************************************************************************************************/
1863 void ShhherCommand::calcNewDistancesParent(int startSeq, int stopSeq){
1868 vector<double> newTau(numOTUs,0);
1869 vector<double> norms(numSeqs, 0);
1870 nSeqsPerOTU.assign(numOTUs, 0);
1872 for(int i=startSeq;i<stopSeq;i++){
1873 int indexOffset = i * numOTUs;
1875 double offset = 1e8;
1877 for(int j=0;j<numOTUs;j++){
1878 if(weight[j] > MIN_WEIGHT && change[j] == 1){
1879 dist[indexOffset + j] = getDistToCentroid(centroids[j], i, lengths[i]);
1882 if(weight[j] > MIN_WEIGHT && dist[indexOffset + j] < offset){
1883 offset = dist[indexOffset + j];
1887 for(int j=0;j<numOTUs;j++){
1888 if(weight[j] > MIN_WEIGHT){
1889 newTau[j] = exp(sigma * (-dist[indexOffset + j] + offset)) * weight[j];
1890 norms[i] += newTau[j];
1897 for(int j=0;j<numOTUs;j++){
1898 newTau[j] /= norms[i];
1901 for(int j=0;j<numOTUs;j++){
1902 if(newTau[j] > MIN_TAU){
1904 int oldTotal = total;
1908 singleTau.resize(total, 0);
1909 seqNumber.resize(total, 0);
1910 seqIndex.resize(total, 0);
1912 singleTau[oldTotal] = newTau[j];
1914 aaP[j][nSeqsPerOTU[j]] = oldTotal;
1915 aaI[j][nSeqsPerOTU[j]] = i;
1921 catch(exception& e) {
1922 m->errorOut(e, "ShhherCommand", "calcNewDistancesParent");
1927 /**************************************************************************************************/
1929 void ShhherCommand::setOTUs(){
1932 vector<double> bigTauMatrix(numOTUs * numSeqs, 0.0000);
1934 for(int i=0;i<numOTUs;i++){
1935 for(int j=0;j<nSeqsPerOTU[i];j++){
1936 int index = cumNumSeqs[i] + j;
1937 double tauValue = singleTau[seqNumber[index]];
1938 int sIndex = seqIndex[index];
1939 bigTauMatrix[sIndex * numOTUs + i] = tauValue;
1943 for(int i=0;i<numSeqs;i++){
1944 double maxTau = -1.0000;
1946 for(int j=0;j<numOTUs;j++){
1947 if(bigTauMatrix[i * numOTUs + j] > maxTau){
1948 maxTau = bigTauMatrix[i * numOTUs + j];
1953 otuData[i] = maxOTU;
1956 nSeqsPerOTU.assign(numOTUs, 0);
1958 for(int i=0;i<numSeqs;i++){
1959 int index = otuData[i];
1961 singleTau[i] = 1.0000;
1964 aaP[index][nSeqsPerOTU[index]] = i;
1965 aaI[index][nSeqsPerOTU[index]] = i;
1967 nSeqsPerOTU[index]++;
1971 catch(exception& e) {
1972 m->errorOut(e, "ShhherCommand", "calcNewDistances");
1977 /**************************************************************************************************/
1979 void ShhherCommand::writeQualities(vector<int> otuCounts){
1982 string qualityFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.qual";
1984 ofstream qualityFile;
1985 m->openOutputFile(qualityFileName, qualityFile);
1987 qualityFile.setf(ios::fixed, ios::floatfield);
1988 qualityFile.setf(ios::showpoint);
1989 qualityFile << setprecision(6);
1991 vector<vector<int> > qualities(numOTUs);
1992 vector<double> pr(HOMOPS, 0);
1995 for(int i=0;i<numOTUs;i++){
1999 if(nSeqsPerOTU[i] > 0){
2000 qualities[i].assign(1024, -1);
2002 while(index < numFlowCells){
2003 double maxPrValue = 1e8;
2004 short maxPrIndex = -1;
2005 double count = 0.0000;
2007 pr.assign(HOMOPS, 0);
2009 for(int j=0;j<nSeqsPerOTU[i];j++){
2010 int lIndex = cumNumSeqs[i] + j;
2011 double tauValue = singleTau[seqNumber[lIndex]];
2012 int sequenceIndex = aaI[i][j];
2013 short intensity = flowDataIntI[sequenceIndex * numFlowCells + index];
2017 for(int s=0;s<HOMOPS;s++){
2018 pr[s] += tauValue * singleLookUp[s * NUMBINS + intensity];
2022 maxPrIndex = uniqueFlowgrams[centroids[i] * numFlowCells + index];
2023 maxPrValue = pr[maxPrIndex];
2025 if(count > MIN_COUNT){
2027 double norm = 0.0000;
2029 for(int s=0;s<HOMOPS;s++){
2030 norm += exp(-(pr[s] - maxPrValue));
2033 for(int s=1;s<=maxPrIndex;s++){
2035 double temp = 0.0000;
2037 U += exp(-(pr[s-1]-maxPrValue))/norm;
2045 temp = floor(-10 * temp);
2046 value = (int)floor(temp);
2047 if(value > 100){ value = 100; }
2049 qualities[i][base] = (int)value;
2059 if(otuCounts[i] > 0){
2060 qualityFile << '>' << seqNameVector[mapUniqueToSeq[i]] << endl;
2062 int j=4; //need to get past the first four bases
2063 while(qualities[i][j] != -1){
2064 qualityFile << qualities[i][j] << ' ';
2067 qualityFile << endl;
2070 qualityFile.close();
2071 outputNames.push_back(qualityFileName);
2074 catch(exception& e) {
2075 m->errorOut(e, "ShhherCommand", "writeQualities");
2080 /**************************************************************************************************/
2082 void ShhherCommand::writeSequences(vector<int> otuCounts){
2085 string fastaFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.fasta";
2087 m->openOutputFile(fastaFileName, fastaFile);
2089 vector<string> names(numOTUs, "");
2091 for(int i=0;i<numOTUs;i++){
2092 int index = centroids[i];
2094 if(otuCounts[i] > 0){
2095 fastaFile << '>' << seqNameVector[aaI[i][0]] << endl;
2099 for(int j=0;j<numFlowCells;j++){
2101 char base = flowOrder[j % 4];
2102 for(int k=0;k<uniqueFlowgrams[index * numFlowCells + j];k++){
2107 fastaFile << newSeq.substr(4) << endl;
2112 outputNames.push_back(fastaFileName);
2114 if(compositeFASTAFileName != ""){
2115 m->appendFiles(fastaFileName, compositeFASTAFileName);
2118 catch(exception& e) {
2119 m->errorOut(e, "ShhherCommand", "writeSequences");
2124 /**************************************************************************************************/
2126 void ShhherCommand::writeNames(vector<int> otuCounts){
2128 string nameFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.names";
2130 m->openOutputFile(nameFileName, nameFile);
2132 for(int i=0;i<numOTUs;i++){
2133 if(otuCounts[i] > 0){
2134 nameFile << seqNameVector[aaI[i][0]] << '\t' << seqNameVector[aaI[i][0]];
2136 for(int j=1;j<nSeqsPerOTU[i];j++){
2137 nameFile << ',' << seqNameVector[aaI[i][j]];
2144 outputNames.push_back(nameFileName);
2147 if(compositeNamesFileName != ""){
2148 m->appendFiles(nameFileName, compositeNamesFileName);
2151 catch(exception& e) {
2152 m->errorOut(e, "ShhherCommand", "writeNames");
2157 /**************************************************************************************************/
2159 void ShhherCommand::writeGroups(){
2161 string fileRoot = flowFileName.substr(0,flowFileName.find_last_of('.'));
2162 string groupFileName = fileRoot + ".shhh.groups";
2164 m->openOutputFile(groupFileName, groupFile);
2166 for(int i=0;i<numSeqs;i++){
2167 groupFile << seqNameVector[i] << '\t' << fileRoot << endl;
2170 outputNames.push_back(groupFileName);
2173 catch(exception& e) {
2174 m->errorOut(e, "ShhherCommand", "writeGroups");
2179 /**************************************************************************************************/
2181 void ShhherCommand::writeClusters(vector<int> otuCounts){
2183 string otuCountsFileName = flowFileName.substr(0,flowFileName.find_last_of('.')) + ".shhh.counts";
2184 ofstream otuCountsFile;
2185 m->openOutputFile(otuCountsFileName, otuCountsFile);
2187 string bases = flowOrder;
2189 for(int i=0;i<numOTUs;i++){
2190 //output the translated version of the centroid sequence for the otu
2191 if(otuCounts[i] > 0){
2192 int index = centroids[i];
2194 otuCountsFile << "ideal\t";
2195 for(int j=8;j<numFlowCells;j++){
2196 char base = bases[j % 4];
2197 for(int s=0;s<uniqueFlowgrams[index * numFlowCells + j];s++){
2198 otuCountsFile << base;
2201 otuCountsFile << endl;
2203 for(int j=0;j<nSeqsPerOTU[i];j++){
2204 int sequence = aaI[i][j];
2205 otuCountsFile << seqNameVector[sequence] << '\t';
2209 for(int k=0;k<lengths[sequence];k++){
2210 char base = bases[k % 4];
2211 int freq = int(0.01 * (double)flowDataIntI[sequence * numFlowCells + k] + 0.5);
2213 for(int s=0;s<freq;s++){
2215 //otuCountsFile << base;
2218 otuCountsFile << newSeq.substr(4) << endl;
2220 otuCountsFile << endl;
2223 otuCountsFile.close();
2224 outputNames.push_back(otuCountsFileName);
2227 catch(exception& e) {
2228 m->errorOut(e, "ShhherCommand", "writeClusters");
2233 //**********************************************************************************************************************